Skip to content

Commit 3d70877

Browse files
authored
Merge pull request #15533 from cxp484/master
FDS Source: Tabulate CVODE warning messages in verbose mode.
2 parents 97e7341 + 8f21a99 commit 3d70877

File tree

3 files changed

+88
-7
lines changed

3 files changed

+88
-7
lines changed

Source/chem.f90

Lines changed: 11 additions & 7 deletions
Original file line numberDiff line numberDiff line change
@@ -807,6 +807,7 @@ SUBROUTINE CVODE_SERIAL(CC,ZZ_0, TMP_IN, TMP_UNMIX, PR_IN, ZETA0, TAU_MIX, CELL_
807807
TMP_OUT, CHEM_TIME, WRITE_SUBSTEPS, CVODE_CALL_OPTION)
808808
USE PHYSICAL_FUNCTIONS, ONLY : MOLAR_CONC_TO_MASS_FRAC, CALC_EQUIV_RATIO, GET_ENTHALPY, GET_MOLECULAR_WEIGHT
809809
USE COMP_FUNCTIONS, ONLY: GET_FILE_NUMBER
810+
USE CHEMCONS, ONLY: CVODE_WARNING_CELLS,CVODE_ERR_CODE_MIN,CVODE_ERR_CODE_MAX
810811
USE GLOBAL_CONSTANTS
811812
USE FCVODE_MOD ! FORTRAN INTERFACE TO CVODE
812813
USE FSUNDIALS_CONTEXT_MOD ! FORTRAN INTERFACE TO SUNCONTEXT
@@ -1015,15 +1016,18 @@ SUBROUTINE CVODE_SERIAL(CC,ZZ_0, TMP_IN, TMP_UNMIX, PR_IN, ZETA0, TAU_MIX, CELL_
10151016
ENDIF
10161017

10171018
IF (IERR_C .NE. CV_SUCCESS) THEN
1018-
IF (IERR_C == CV_TOO_MUCH_WORK) THEN
1019-
WRITE(LU_ERR,'(A, 2E18.8, I8, A)')" WARN: CVODE took all internal substeps. CUR_CFD_TIME, DT, MAXTRY=", CUR_CFD_TIME, &
1020-
(TEND-TCUR), MAXTRY, ". If the warning persists, reduce the timestep."
1021-
ELSE
1022-
WRITE(LU_ERR,'(A, I4, A, 2E18.8, A)')" WARN: CVODE didn't finish ODE solution with message code:", IERR_C, &
1023-
" and CUR_CFD_TIME, DT=", CUR_CFD_TIME, (TEND-TCUR), ". If the warning persists, reduce the timestep."
1019+
IF(IERR_C>=CVODE_ERR_CODE_MIN .AND. IERR_C<=CVODE_ERR_CODE_MAX) THEN
1020+
CVODE_WARNING_CELLS(IERR_C) = CVODE_WARNING_CELLS(IERR_C) + 1
10241021
ENDIF
1025-
10261022
IF (DEBUG) THEN
1023+
IF (IERR_C == CV_TOO_MUCH_WORK) THEN
1024+
WRITE(LU_ERR,'(A, 2E18.8, I8, A)')" WARN: CVODE took all internal substeps. CUR_CFD_TIME, DT, MAXTRY=", CUR_CFD_TIME, &
1025+
(TEND-TCUR), MAXTRY, ". If the warning persists, reduce the timestep."
1026+
ELSE
1027+
WRITE(LU_ERR,'(A, I4, A, 2E18.8, A)')" WARN: CVODE didn't finish ODE solution with message code:", IERR_C, &
1028+
" and CUR_CFD_TIME, DT=", CUR_CFD_TIME, (TEND-TCUR), ". If the warning persists, reduce the timestep."
1029+
ENDIF
1030+
10271031
CALL MOLAR_CONC_TO_MASS_FRAC(CC(1:N_TRACKED_SPECIES), ZZ(1:N_TRACKED_SPECIES))
10281032
CALL CALC_EQUIV_RATIO(ZZ(1:N_TRACKED_SPECIES), EQUIV)
10291033
DO NS = 1, N_TRACKED_SPECIES

Source/cons.f90

Lines changed: 12 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -912,6 +912,10 @@ MODULE CHEMCONS
912912
INTEGER :: MAX_CVODE_SUBSTEPS=100000
913913
INTEGER :: CVODE_MAX_TRY=4
914914
INTEGER :: CVODE_ORDER=0
915+
INTEGER :: CVODE_ERR_CODE_MIN=-100
916+
INTEGER :: CVODE_ERR_CODE_MAX=100
917+
INTEGER :: CVODE_WARNING_CELLS(-100:100)! Index of the array is error code and value is Cell count
918+
CHARACTER(LEN=100) :: CVODE_WARN_MESSAGES(-100:100)
915919

916920
! FOR WRITING CVODE SUBSTEPS
917921
LOGICAL :: WRITE_CVODE_SUBSTEPS = .FALSE.
@@ -932,4 +936,12 @@ MODULE CHEMCONS
932936
INTEGER :: N_IGNITION_ZONES = 0
933937
TYPE(IGNITION_ZONE_TYPE), DIMENSION(MAX_IGNITION_ZONES) :: IGNITION_ZONES !< Coordinates of ignition zones
934938

939+
CONTAINS
940+
SUBROUTINE INIT_CVODE_WARN_MESSAGES()
941+
CVODE_WARN_MESSAGES = 'CVODE didn''t finish ODE solution with this code.'
942+
CVODE_WARN_MESSAGES(-1) = 'CVODE took all internal substeps.'
943+
CVODE_WARN_MESSAGES(-3) = 'Minimum step size was reached.'
944+
CVODE_WARN_MESSAGES(-4) = 'Convergence test failure.'
945+
END SUBROUTINE INIT_CVODE_WARN_MESSAGES
946+
935947
END MODULE CHEMCONS

Source/fire.f90

Lines changed: 65 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -41,6 +41,7 @@ SUBROUTINE COMBUSTION_LOAD_BALANCED(T,DT)
4141

4242
USE SOOT_ROUTINES, ONLY: SOOT_SURFACE_OXIDATION
4343
USE COMP_FUNCTIONS, ONLY: CURRENT_TIME
44+
USE CHEMCONS, ONLY: CVODE_WARNING_CELLS,INIT_CVODE_WARN_MESSAGES
4445
REAL(EB), INTENT(IN) :: T,DT
4546
INTEGER :: NM,ICC,JCC
4647
REAL(EB) :: TNOW
@@ -55,12 +56,14 @@ SUBROUTINE COMBUSTION_LOAD_BALANCED(T,DT)
5556
COMBUSTION_INIT = .TRUE.
5657
NVAR_TO_SEND = N_TRACKED_SPECIES +7 ! NS, TEMP, RHO, PRES, MU, DELTA, VOL, IGN_ZN
5758
NVAR_TO_RECEIVE = N_TRACKED_SPECIES +4 ! NS, Q_OUT, MIX_TIME_OUT, CHI_R_OUT, CHEM_SUBIT_TMP_OUT
59+
CALL INIT_CVODE_WARN_MESSAGES()
5860
ENDIF
5961

6062
DO NM=LOWER_MESH_INDEX,UPPER_MESH_INDEX
6163
CALL POINT_TO_MESH(NM)
6264
Q = 0._EB
6365
CHI_R = 0._EB
66+
CVODE_WARNING_CELLS = 0
6467
IF (CC_IBM) THEN
6568
DO ICC=1,MESHES(NM)%N_CUTCELL_MESH
6669
DO JCC=1,CUT_CELL(ICC)%NCELL
@@ -109,6 +112,7 @@ SUBROUTINE COMBUSTION_GENERAL_LOAD_BALANCED(T,DT)
109112
REAL(EB) :: MIX_TIME_OUT, Q_OUT, CHI_R_OUT, MYTEMP, MYRHO, MYMU, DELTA, VOL, TNOW2
110113
INTEGER :: IGN_ZN
111114

115+
112116
Q_EXISTS = .FALSE.
113117

114118
!------
@@ -329,6 +333,13 @@ SUBROUTINE COMBUSTION_GENERAL_LOAD_BALANCED(T,DT)
329333
ENDDO
330334
ENDIF
331335

336+
! Output accumulated CVODE warnings (instead of cell by cell) in verbose mode
337+
#ifdef WITH_SUNDIALS
338+
IF (COMBUSTION_ODE_SOLVER == CVODE_SOLVER .AND. VERBOSE) THEN
339+
CALL DUMP_CVODE_WARNING_SUMMARY(NCHEM_ACTIVE_CELLS_AND_CC)
340+
ENDIF
341+
#endif
342+
332343
END SUBROUTINE COMBUSTION_GENERAL_LOAD_BALANCED
333344

334345

@@ -1269,6 +1280,60 @@ SUBROUTINE CALC_AFT_REAC_AND_PROD(ZZ,ZZ_REAC,ZZ_PROD)
12691280

12701281
END SUBROUTINE CALC_AFT_REAC_AND_PROD
12711282

1283+
!> \brief Dump CVODE warning summary
1284+
SUBROUTINE DUMP_CVODE_WARNING_SUMMARY(NCHEM_ACTIVE_CELLS_AND_CC)
1285+
USE CHEMCONS, ONLY: CVODE_WARNING_CELLS,CVODE_ERR_CODE_MIN,CVODE_ERR_CODE_MAX,CVODE_WARN_MESSAGES
1286+
1287+
INTEGER, INTENT(IN) :: NCHEM_ACTIVE_CELLS_AND_CC
1288+
INTEGER :: NCHEM_ACTIVE_CELLS_AND_CC_GLOBAL,IERR, CVODE_ERR_CODE, TOT_WARN_CELLS
1289+
REAL(EB) :: PCNT
1290+
CHARACTER(LEN=16) :: CELL_STR, PCNT_STR
1291+
CHARACTER(LEN=40) :: CELL_INFO
1292+
1293+
IF (N_MPI_PROCESSES > 1) THEN
1294+
CALL MPI_ALLREDUCE(MPI_IN_PLACE, CVODE_WARNING_CELLS, SIZE(CVODE_WARNING_CELLS), &
1295+
MPI_INTEGER, MPI_SUM, MPI_COMM_WORLD, IERR)
1296+
1297+
CALL MPI_ALLREDUCE(NCHEM_ACTIVE_CELLS_AND_CC, NCHEM_ACTIVE_CELLS_AND_CC_GLOBAL, &
1298+
1, MPI_INTEGER, MPI_SUM, MPI_COMM_WORLD, IERR)
1299+
ELSE
1300+
NCHEM_ACTIVE_CELLS_AND_CC_GLOBAL = NCHEM_ACTIVE_CELLS_AND_CC
1301+
END IF
1302+
TOT_WARN_CELLS = SUM(CVODE_WARNING_CELLS)
1303+
1304+
IF (TOT_WARN_CELLS > 0 .AND. NCHEM_ACTIVE_CELLS_AND_CC_GLOBAL > 0 .AND. MY_RANK == 0) THEN
1305+
1306+
WRITE(LU_ERR, '(A)') '--------------- CVODE Warning Summary ----------------'
1307+
WRITE(LU_ERR, '(A)') ' Code Cells (% Chem) Message'
1308+
1309+
DO CVODE_ERR_CODE = CVODE_ERR_CODE_MIN, CVODE_ERR_CODE_MAX
1310+
IF (CVODE_WARNING_CELLS(CVODE_ERR_CODE) > 0) THEN
1311+
1312+
PCNT = REAL(CVODE_WARNING_CELLS(CVODE_ERR_CODE)) / &
1313+
REAL(NCHEM_ACTIVE_CELLS_AND_CC_GLOBAL) * 100.0_EB
1314+
1315+
! Convert to strings and left-align
1316+
WRITE(CELL_STR, '(I0)') CVODE_WARNING_CELLS(CVODE_ERR_CODE)
1317+
WRITE(PCNT_STR, '(F6.2)') PCNT
1318+
CELL_STR = ADJUSTL(CELL_STR)
1319+
PCNT_STR = ADJUSTL(PCNT_STR)
1320+
1321+
! Combine count + percentage into one left-aligned string
1322+
WRITE(CELL_INFO, '(A, " (", A, " %)")') TRIM(CELL_STR), TRIM(PCNT_STR)
1323+
CELL_INFO = ADJUSTL(CELL_INFO)
1324+
1325+
! Print one line of summary
1326+
WRITE(LU_ERR, '(1X, I4, 5X, A, 5X, A)') CVODE_ERR_CODE, &
1327+
TRIM(CELL_INFO), TRIM(CVODE_WARN_MESSAGES(CVODE_ERR_CODE))
1328+
END IF
1329+
END DO
1330+
1331+
WRITE(LU_ERR, '(A)') '------------------------------------------------------'
1332+
1333+
END IF
1334+
1335+
END SUBROUTINE DUMP_CVODE_WARNING_SUMMARY
1336+
12721337
#endif
12731338

12741339
SUBROUTINE CHECK_AUTO_IGNITION(EXTINCT,TMP_IN,AIT,IIC,JJC,KKC,REAC_INDEX)

0 commit comments

Comments
 (0)