propagate non-binary status information

This commit is contained in:
Martin Diehl 2024-01-10 00:36:22 +01:00
parent 89ebad19d5
commit 0384f93d46
No known key found for this signature in database
GPG Key ID: 1FD50837275A0A9B
3 changed files with 11 additions and 13 deletions

View File

@ -57,7 +57,7 @@ module grid_thermal_spectral
integer :: totalIter = 0 !< total iteration in current increment integer :: totalIter = 0 !< total iteration in current increment
real(pREAL), dimension(3,3) :: K_ref real(pREAL), dimension(3,3) :: K_ref
real(pREAL) :: mu_ref, Delta_t_ real(pREAL) :: mu_ref, Delta_t_
logical :: broken integer(kind(STATUS_OK)) :: status
public :: & public :: &
grid_thermal_spectral_init, & grid_thermal_spectral_init, &
@ -208,7 +208,7 @@ function grid_thermal_spectral_solution(Delta_t) result(solution)
call SNESGetConvergedReason(SNES_thermal,reason,err_PETSc) call SNESGetConvergedReason(SNES_thermal,reason,err_PETSc)
CHKERRQ(err_PETSc) CHKERRQ(err_PETSc)
solution%converged = reason > 0 .and. .not. broken solution%converged = reason > 0 .and. status == STATUS_OK
solution%iterationsNeeded = merge(totalIter,num%itmax,solution%converged) solution%iterationsNeeded = merge(totalIter,num%itmax,solution%converged)
call SNESGetDM(SNES_thermal,DM_thermal,err_PETSc) call SNESGetDM(SNES_thermal,DM_thermal,err_PETSc)
@ -323,11 +323,9 @@ subroutine formResidual(residual_subdomain,x_scal,r,dummy,err_PETSc)
integer :: i, j, k, ce integer :: i, j, k, ce
real(pREAL), dimension(3,cells(1),cells(2),cells3) :: vectorField real(pREAL), dimension(3,cells(1),cells(2),cells3) :: vectorField
integer(kind(STATUS_OK)) :: status
call homogenization_thermal_response(status,Delta_t_,1,product(cells(1:2))*cells3) call homogenization_thermal_response(status,Delta_t_,1,product(cells(1:2))*cells3)
broken = status /= STATUS_OK
associate(T => x_scal) associate(T => x_scal)
vectorField = utilities_ScalarGradient(T) vectorField = utilities_ScalarGradient(T)

View File

@ -121,21 +121,19 @@ end subroutine FEM_utilities_init
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief calculates constitutive response !> @brief calculates constitutive response
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
subroutine utilities_constitutiveResponse(broken, Delta_t,P_av,forwardData) subroutine utilities_constitutiveResponse(status, Delta_t,P_av,forwardData)
logical, intent(out) :: broken integer(kind(STATUS_OK)), intent(out) :: status
real(pREAL), intent(in) :: Delta_t !< loading time real(pREAL), intent(in) :: Delta_t !< loading time
logical, intent(in) :: forwardData !< age results logical, intent(in) :: forwardData !< age results
real(pREAL),intent(out), dimension(3,3) :: P_av !< average PK stress real(pREAL),intent(out), dimension(3,3) :: P_av !< average PK stress
integer(MPI_INTEGER_KIND) :: err_MPI integer(MPI_INTEGER_KIND) :: err_MPI
integer(kind(STATUS_OK)) :: status
print'(/,1x,a)', '... evaluating constitutive response ......................................' print'(/,1x,a)', '... evaluating constitutive response ......................................'
call homogenization_mechanical_response(status,Delta_t,1,mesh_maxNips*mesh_NcpElems) ! calculate P field call homogenization_mechanical_response(status,Delta_t,1,mesh_maxNips*mesh_NcpElems) ! calculate P field
broken = status /= STATUS_OK
cutBack = .false. cutBack = .false.
P_av = sum(homogenization_P,dim=3) * wgt P_av = sum(homogenization_P,dim=3) * wgt

View File

@ -25,6 +25,7 @@ module mesh_mechanical_FEM
use FEM_quadrature use FEM_quadrature
use homogenization use homogenization
use math use math
use constants
#if (PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>14) && !defined(PETSC_HAVE_MPI_F90MODULE_VISIBILITY) #if (PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>14) && !defined(PETSC_HAVE_MPI_F90MODULE_VISIBILITY)
implicit none(type,external) implicit none(type,external)
@ -68,7 +69,8 @@ module mesh_mechanical_FEM
character(len=pSTRLEN) :: incInfo character(len=pSTRLEN) :: incInfo
real(pREAL), dimension(3,3) :: & real(pREAL), dimension(3,3) :: &
P_av = 0.0_pREAL P_av = 0.0_pREAL
logical :: ForwardData, broken logical :: ForwardData
integer(kind(STATUS_OK)) :: status
real(pREAL), parameter :: eps = 1.0e-18_pREAL real(pREAL), parameter :: eps = 1.0e-18_pREAL
external :: & ! ToDo: write interfaces external :: & ! ToDo: write interfaces
@ -311,7 +313,7 @@ subroutine FEM_mechanical_init(mechBC,num_mesh)
call DMPlexVecSetClosure(mechanical_mesh,section,solution_local,cell,px_scal,5,err_PETSc) call DMPlexVecSetClosure(mechanical_mesh,section,solution_local,cell,px_scal,5,err_PETSc)
CHKERRQ(err_PETSc) CHKERRQ(err_PETSc)
end do end do
call utilities_constitutiveResponse(broken,0.0_pREAL,devNull,.true.) call utilities_constitutiveResponse(status,0.0_pREAL,devNull,.true.)
end subroutine FEM_mechanical_init end subroutine FEM_mechanical_init
@ -458,8 +460,8 @@ subroutine FEM_mechanical_formResidual(dm_local,xx_local,f_local,dummy,err_PETSc
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! evaluate constitutive response ! evaluate constitutive response
call utilities_constitutiveResponse(broken,params%Delta_t,P_av,ForwardData) call utilities_constitutiveResponse(status,params%Delta_t,P_av,ForwardData)
call MPI_Allreduce(MPI_IN_PLACE,broken,1_MPI_INTEGER_KIND,MPI_LOGICAL,MPI_LOR,MPI_COMM_WORLD,err_MPI) call MPI_Allreduce(MPI_IN_PLACE,status,1_MPI_INTEGER_KIND,MPI_INTEGER,MPI_SUM,MPI_COMM_WORLD,err_MPI)
call parallelization_chkerr(err_MPI) call parallelization_chkerr(err_MPI)
ForwardData = .false. ForwardData = .false.
@ -747,7 +749,7 @@ subroutine FEM_mechanical_converged(snes_local,PETScIter,xnorm,snorm,fnorm,reaso
divTol = max(maxval(abs(P_av(1:dimPlex,1:dimPlex)))*num%eps_struct_rtol,num%eps_struct_atol) divTol = max(maxval(abs(P_av(1:dimPlex,1:dimPlex)))*num%eps_struct_rtol,num%eps_struct_atol)
call SNESConvergedDefault(snes_local,PETScIter,xnorm,snorm,fnorm/divTol,reason,dummy,err_PETSc) call SNESConvergedDefault(snes_local,PETScIter,xnorm,snorm,fnorm/divTol,reason,dummy,err_PETSc)
CHKERRQ(err_PETSc) CHKERRQ(err_PETSc)
if (broken) reason = SNES_DIVERGED_FUNCTION_DOMAIN if (status /= STATUS_OK) reason = SNES_DIVERGED_FUNCTION_DOMAIN
print'(/,1x,a,a,i0,a,f0.3)', trim(incInfo), & print'(/,1x,a,a,i0,a,f0.3)', trim(incInfo), &
' @ Iteration ',PETScIter,' mechanical residual norm = ',fnorm/divTol ' @ Iteration ',PETScIter,' mechanical residual norm = ',fnorm/divTol
print'(/,1x,a,/,2(3(2x,f12.4,1x)/),3(2x,f12.4,1x))', & print'(/,1x,a,/,2(3(2x,f12.4,1x)/),3(2x,f12.4,1x))', &