Skip to content

Commit 8800c84

Browse files
committed
FDS Source: add time ramp RAMP_ZETA_0 for INITIAL_UNMIXED_FRACTION on COMB
1 parent f72fe51 commit 8800c84

File tree

3 files changed

+22
-13
lines changed

3 files changed

+22
-13
lines changed

Source/cons.f90

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -439,6 +439,7 @@ MODULE GLOBAL_CONSTANTS
439439
INTEGER :: N_PASSIVE_SCALARS=0 !< Number of passive scalars
440440
INTEGER :: N_TOTAL_SCALARS=0 !< Number of total scalars, tracked and passive
441441
INTEGER :: N_FIXED_CHEMISTRY_SUBSTEPS=-1 !< Number of chemistry substeps in combustion routine
442+
INTEGER :: ZETA_0_RAMP_INDEX=0 !< Ramp index for initial unmixed fraction
442443

443444
LOGICAL :: OUTPUT_CHEM_IT=.FALSE.
444445
LOGICAL :: REAC_SOURCE_CHECK=.FALSE.

Source/fire.f90

Lines changed: 15 additions & 12 deletions
Original file line numberDiff line numberDiff line change
@@ -96,12 +96,13 @@ SUBROUTINE COMBUSTION_GENERAL_LOAD_BALANCED(T,DT)
9696
USE PHYSICAL_FUNCTIONS, ONLY: GET_SPECIFIC_GAS_CONSTANT,GET_MASS_FRACTION_ALL,GET_SPECIFIC_HEAT,GET_MOLECULAR_WEIGHT, &
9797
GET_SENSIBLE_ENTHALPY_Z,IS_REALIZABLE
9898
USE CHEMCONS, ONLY: DO_CHEM_LOAD_BALANCE
99+
USE MATH_FUNCTIONS, ONLY: EVALUATE_RAMP
99100

100101
INTEGER :: NM,I,J,K,NC, CHEM_SUBIT_TMP,CHEM_SUBIT_TMP_OUT,ICC,JCC,INDX,NRECEIVE_CELLS
101102
INTEGER :: NCHEM_ACTIVE_CELLS_AND_CC, NCHEM_ACTIVE_CELLS_AND_CC_GLOBAL, IERR, NCHEM_ACTIVE_CELLS, NCHEM_ACTIVE_CC
102103
REAL(EB), INTENT(IN) :: T,DT
103104
REAL(EB) :: ZZ_GET(1:N_TRACKED_SPECIES),DZZ(1:N_TRACKED_SPECIES),&
104-
REAC_SOURCE_TERM_TMP(N_TRACKED_SPECIES),Q_REAC_TMP(N_REACTIONS),VCELL,PRES
105+
REAC_SOURCE_TERM_TMP(N_TRACKED_SPECIES),Q_REAC_TMP(N_REACTIONS),VCELL,PRES,ZETA_0,TSI
105106
LOGICAL :: Q_EXISTS
106107
TYPE (CC_CUTCELL_TYPE), POINTER :: CC
107108
LOGICAL :: IS_CHEM_ACTIVE
@@ -175,7 +176,10 @@ SUBROUTINE COMBUSTION_GENERAL_LOAD_BALANCED(T,DT)
175176
! For 1 MPI process (serial) distribution is not needed.
176177
!------
177178

178-
DO_CHEM_LOAD_BALANCE_IF : IF(N_MPI_PROCESSES > 1 .AND. DO_CHEM_LOAD_BALANCE) THEN
179+
TSI = T-T_BEGIN
180+
ZETA_0 = EVALUATE_RAMP(TSI,ZETA_0_RAMP_INDEX)*INITIAL_UNMIXED_FRACTION
181+
182+
CHEM_LOAD_BALANCE_IF : IF(N_MPI_PROCESSES > 1 .AND. DO_CHEM_LOAD_BALANCE) THEN
179183
IF (.NOT. ALLOCATED(NCHEM_ACTIVE_CELLS_AND_CC_RANK_WISE)) ALLOCATE(NCHEM_ACTIVE_CELLS_AND_CC_RANK_WISE(N_MPI_PROCESSES))
180184
NCHEM_ACTIVE_CELLS_AND_CC_RANK_WISE = 0
181185
CALL MPI_ALLGATHER(NCHEM_ACTIVE_CELLS_AND_CC,1,MPI_INTEGER,NCHEM_ACTIVE_CELLS_AND_CC_RANK_WISE,1,MPI_INTEGER,MPI_COMM_WORLD,IERR)
@@ -205,7 +209,7 @@ SUBROUTINE COMBUSTION_GENERAL_LOAD_BALANCED(T,DT)
205209
! Call combustion integration routine for Cartesian cell (I,J,K)
206210
CALL COMBUSTION_MODEL( T,DT,ZZ_GET,Q_OUT,MIX_TIME_OUT,CHI_R_OUT,&
207211
CHEM_SUBIT_TMP_OUT,REAC_SOURCE_TERM_TMP,Q_REAC_TMP,&
208-
MYTEMP,MYRHO,PRES,MYMU,DELTA,VOL)
212+
MYTEMP,MYRHO,PRES,MYMU,DELTA,VOL,ZETA_0)
209213
!***************************************************************************************
210214

211215
RESULTS_TO_SEND_ARRAY(1:N_TRACKED_SPECIES,INDX) = ZZ_GET
@@ -216,7 +220,6 @@ SUBROUTINE COMBUSTION_GENERAL_LOAD_BALANCED(T,DT)
216220

217221
ENDDO
218222

219-
220223
IF (COMBUSTION_ODE_SOLVER == CVODE_SOLVER) THEN
221224
T_CHEM_ODE = T_CHEM_ODE+CURRENT_TIME()-TNOW2
222225
ENDIF
@@ -230,7 +233,7 @@ SUBROUTINE COMBUSTION_GENERAL_LOAD_BALANCED(T,DT)
230233
IF (STOP_STATUS/=NO_STOP) RETURN
231234
ENDIF !NCHEM_ACTIVE_CELLS_AND_CC_GLOBAL >0
232235

233-
ELSE DO_CHEM_LOAD_BALANCE_IF
236+
ELSE CHEM_LOAD_BALANCE_IF
234237

235238
IF (NCHEM_ACTIVE_CELLS_AND_CC >0) THEN
236239
! Serial chemistry:
@@ -254,7 +257,7 @@ SUBROUTINE COMBUSTION_GENERAL_LOAD_BALANCED(T,DT)
254257
CALL COMBUSTION_MODEL( T,DT,ZZ_GET,Q(I,J,K),MIX_TIME(I,J,K),CHI_R(I,J,K),&
255258
CHEM_SUBIT_TMP,REAC_SOURCE_TERM_TMP,Q_REAC_TMP,&
256259
TMP(I,J,K),RHO(I,J,K),PRES, MU(I,J,K),&
257-
LES_FILTER_WIDTH(I,J,K),DX(I)*DY(J)*DZ(K),IIC=I,JJC=J,KKC=K )
260+
LES_FILTER_WIDTH(I,J,K),DX(I)*DY(J)*DZ(K),ZETA_0,IIC=I,JJC=J,KKC=K )
258261
!***************************************************************************************
259262
!IF (STOP_STATUS/=NO_STOP) RETURN
260263
IF (OUTPUT_CHEM_IT) CHEM_SUBIT(I,J,K) = CHEM_SUBIT_TMP
@@ -278,7 +281,7 @@ SUBROUTINE COMBUSTION_GENERAL_LOAD_BALANCED(T,DT)
278281
CC%CHI_R(JCC),&
279282
CHEM_SUBIT_TMP,REAC_SOURCE_TERM_TMP,Q_REAC_TMP,&
280283
CC%TMP(JCC),CC%RHO(JCC),PRES,MU(I,J,K),&
281-
LES_FILTER_WIDTH(I,J,K),CC%VOLUME(JCC),IIC=I,JJC=J,KKC=K)
284+
LES_FILTER_WIDTH(I,J,K),CC%VOLUME(JCC),ZETA_0,IIC=I,JJC=J,KKC=K)
282285
!***************************************************************************************
283286
CALL SET_SPECIES_SOURCE_TERM_CUTCELL(DT, ICC, JCC, ZZ_GET, DZZ, REAC_SOURCE_TERM_TMP, Q_REAC_TMP)
284287
END DO
@@ -291,7 +294,7 @@ SUBROUTINE COMBUSTION_GENERAL_LOAD_BALANCED(T,DT)
291294
T_CHEM_ODE = T_CHEM_ODE+CURRENT_TIME()-TNOW2
292295
ENDIF
293296
ENDIF !NCHEM_ACTIVE_CELLS_AND_CC >0
294-
ENDIF DO_CHEM_LOAD_BALANCE_IF
297+
ENDIF CHEM_LOAD_BALANCE_IF
295298

296299
! This volume refactoring is needed for RADIATION_FVM (CHI_R, Q) and plotting slices:
297300
IF(CC_IBM) THEN
@@ -712,13 +715,13 @@ END SUBROUTINE GATHER_CELLS_FROM_MPI_PROCESSES
712715

713716

714717
SUBROUTINE COMBUSTION_MODEL(T,DT,ZZ_GET,Q_OUT,MIX_TIME_OUT,CHI_R_OUT,CHEM_SUBIT_OUT,REAC_SOURCE_TERM_OUT,Q_REAC_OUT,&
715-
TMP_IN,RHO_IN,PRES_IN,MU_IN,DELTA,CELL_VOLUME,IIC,JJC,KKC)
718+
TMP_IN,RHO_IN,PRES_IN,MU_IN,DELTA,CELL_VOLUME,ZETA_0_IN,IIC,JJC,KKC)
716719
USE MATH_FUNCTIONS, ONLY: EVALUATE_RAMP
717720
USE PHYSICAL_FUNCTIONS, ONLY: GET_REALIZABLE_MF
718721
USE COMP_FUNCTIONS, ONLY: SHUTDOWN
719722
USE CHEMCONS, ONLY: ODE_MIN_ATOL
720723
INTEGER, INTENT(IN), OPTIONAL :: IIC,JJC,KKC
721-
REAL(EB), INTENT(IN) :: T,DT,RHO_IN,PRES_IN,MU_IN,DELTA,CELL_VOLUME
724+
REAL(EB), INTENT(IN) :: T,DT,RHO_IN,PRES_IN,MU_IN,DELTA,CELL_VOLUME,ZETA_0_IN
722725
REAL(EB), INTENT(OUT) :: Q_OUT,MIX_TIME_OUT,CHI_R_OUT,REAC_SOURCE_TERM_OUT(N_TRACKED_SPECIES),Q_REAC_OUT(N_REACTIONS)
723726
INTEGER, INTENT(OUT) :: CHEM_SUBIT_OUT
724727
REAL(EB), INTENT(INOUT) :: ZZ_GET(1:N_TRACKED_SPECIES)
@@ -730,7 +733,7 @@ SUBROUTINE COMBUSTION_MODEL(T,DT,ZZ_GET,Q_OUT,MIX_TIME_OUT,CHI_R_OUT,CHEM_SUBIT_
730733
Q_REAC_SUM(1:N_REACTIONS),Q_SUM_CHI_R,CHI_R_SUM,TIME_RAMP_FACTOR,&
731734
TOTAL_MIXED_MASS_1,TOTAL_MIXED_MASS_2,TOTAL_MIXED_MASS_4,TOTAL_MIXED_MASS,&
732735
ZETA_1,ZETA_2,ZETA_4,D_F,TMP_IN,C_U,DT_SUB_OLD,ERR_EST(N_TRACKED_SPECIES),ERR_TOL(N_TRACKED_SPECIES),ERR_TINY,&
733-
ZZ_TEMP(1:N_TRACKED_SPECIES), ATOL(1:N_TRACKED_SPECIES)
736+
ZZ_TEMP(1:N_TRACKED_SPECIES),ATOL(1:N_TRACKED_SPECIES)
734737
INTEGER :: NR,NS,ITER,TVI,RICH_ITER,TIME_ITER,RICH_ITER_MAX
735738
INTEGER, PARAMETER :: TV_ITER_MIN=5
736739
LOGICAL :: TV_FLUCT(1:N_TRACKED_SPECIES),EXTINCT,NO_REACTIONS,NO_REAC_2,NO_REAC_4
@@ -761,7 +764,7 @@ SUBROUTINE COMBUSTION_MODEL(T,DT,ZZ_GET,Q_OUT,MIX_TIME_OUT,CHI_R_OUT,CHEM_SUBIT_
761764
END SELECT
762765
ENDIF
763766

764-
ZETA_0 = INITIAL_UNMIXED_FRACTION
767+
ZETA_0 = ZETA_0_IN
765768
CELL_MASS = RHO_IN*CELL_VOLUME
766769

767770
DT_SUB_MIN = DT/REAL(MAX_CHEMISTRY_SUBSTEPS,EB)

Source/read.f90

Lines changed: 6 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -4492,11 +4492,12 @@ END SUBROUTINE PROC_SPEC_2
44924492
SUBROUTINE READ_COMB
44934493
USE CHEMCONS, ONLY: ODE_MIN_ATOL, EQUIV_RATIO_CHECK, MIN_EQUIV_RATIO, MAX_EQUIV_RATIO, DO_CHEM_LOAD_BALANCE
44944494
REAL(EB) :: ODE_REL_ERROR
4495+
CHARACTER(LABEL_LENGTH) :: RAMP_ZETA_0='null'
44954496
NAMELIST /COMB/ CHECK_REALIZABILITY,COMPUTE_ADIABATIC_FLAME_TEMPERATURE, DO_CHEM_LOAD_BALANCE, EQUIV_RATIO_CHECK, &
44964497
EXTINCTION_MODEL,FINITE_RATE_MIN_TEMP, FIXED_MIX_TIME,FREE_BURN_TEMPERATURE,&
44974498
INITIAL_UNMIXED_FRACTION, MAX_CHEMISTRY_SUBSTEPS, MAX_EQUIV_RATIO, MIN_EQUIV_RATIO, &
44984499
N_FIXED_CHEMISTRY_SUBSTEPS, ODE_MIN_ATOL,ODE_REL_ERROR,ODE_SOLVER,SUPPRESSION,TAU_CHEM, &
4499-
TAU_FLAME,ZZ_MIN_GLOBAL
4500+
TAU_FLAME,RAMP_ZETA_0,ZZ_MIN_GLOBAL
45004501

45014502
ODE_SOLVER = 'null'
45024503
ODE_REL_ERROR = -1._EB
@@ -4563,6 +4564,10 @@ SUBROUTINE READ_COMB
45634564
IF (ODE_MIN_ATOL < 0.D0) ODE_MIN_ATOL = DBLE(0.1_EB*ZZ_MIN_GLOBAL)
45644565
IF (ODE_REL_ERROR > 0._EB) GLOBAL_ODE_REL_ERROR = ODE_REL_ERROR
45654566

4567+
! Set index for ZETA_0 ramp, if needed
4568+
4569+
IF (RAMP_ZETA_0 /='null') CALL GET_RAMP_INDEX(RAMP_ZETA_0,'TIME',ZETA_0_RAMP_INDEX)
4570+
45664571
END SUBROUTINE READ_COMB
45674572

45684573

0 commit comments

Comments
 (0)