safer reduction
MPI_SUM can lead to invalid/undefined values for enum
This commit is contained in:
parent
712b22b678
commit
27a8b63a3d
|
@ -573,7 +573,7 @@ subroutine formResidual(da_local,x_local, &
|
|||
call utilities_constitutiveResponse(status,P_current,&
|
||||
P_av,C_volAvg,devNull, &
|
||||
F,params%Delta_t,params%rotation_BC)
|
||||
call MPI_Allreduce(MPI_IN_PLACE,status,1_MPI_INTEGER_KIND,MPI_INTEGER,MPI_SUM,MPI_COMM_WORLD,err_MPI)
|
||||
call MPI_Allreduce(MPI_IN_PLACE,status,1_MPI_INTEGER_KIND,MPI_INTEGER,MPI_MAX,MPI_COMM_WORLD,err_MPI)
|
||||
call parallelization_chkerr(err_MPI)
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
|
|
@ -518,7 +518,7 @@ subroutine formResidual(residual_subdomain, F, &
|
|||
call utilities_constitutiveResponse(status,P, &
|
||||
P_av,C_volAvg,C_minMaxAvg, &
|
||||
F,params%Delta_t,params%rotation_BC)
|
||||
call MPI_Allreduce(MPI_IN_PLACE,status,1_MPI_INTEGER_KIND,MPI_INTEGER,MPI_SUM,MPI_COMM_WORLD,err_MPI)
|
||||
call MPI_Allreduce(MPI_IN_PLACE,status,1_MPI_INTEGER_KIND,MPI_INTEGER,MPI_MAX,MPI_COMM_WORLD,err_MPI)
|
||||
call parallelization_chkerr(err_MPI)
|
||||
err_div = utilities_divergenceRMS(P)
|
||||
end associate
|
||||
|
|
|
@ -612,7 +612,7 @@ subroutine formResidual(residual_subdomain, FandF_tau, &
|
|||
#endif
|
||||
P_av,C_volAvg,C_minMaxAvg, &
|
||||
F - r_F_tau/num%beta,params%Delta_t,params%rotation_BC)
|
||||
call MPI_Allreduce(MPI_IN_PLACE,status,1_MPI_INTEGER_KIND,MPI_INTEGER,MPI_SUM,MPI_COMM_WORLD,err_MPI)
|
||||
call MPI_Allreduce(MPI_IN_PLACE,status,1_MPI_INTEGER_KIND,MPI_INTEGER,MPI_MAX,MPI_COMM_WORLD,err_MPI)
|
||||
#ifdef __GFORTRAN__
|
||||
err_div = utilities_divergenceRMS(r_F)
|
||||
#else
|
||||
|
|
|
@ -461,7 +461,7 @@ subroutine FEM_mechanical_formResidual(dm_local,xx_local,f_local,dummy,err_PETSc
|
|||
!--------------------------------------------------------------------------------------------------
|
||||
! evaluate constitutive response
|
||||
call utilities_constitutiveResponse(status,params%Delta_t,P_av,ForwardData)
|
||||
call MPI_Allreduce(MPI_IN_PLACE,status,1_MPI_INTEGER_KIND,MPI_INTEGER,MPI_SUM,MPI_COMM_WORLD,err_MPI)
|
||||
call MPI_Allreduce(MPI_IN_PLACE,status,1_MPI_INTEGER_KIND,MPI_INTEGER,MPI_MAX,MPI_COMM_WORLD,err_MPI)
|
||||
call parallelization_chkerr(err_MPI)
|
||||
ForwardData = .false.
|
||||
|
||||
|
|
Loading…
Reference in New Issue