-
Notifications
You must be signed in to change notification settings - Fork 0
/
wall.f90
1760 lines (1457 loc) · 63.2 KB
/
wall.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
MODULE WALL_ROUTINES
! Compute the wall boundary conditions
USE PRECISION_PARAMETERS
USE GLOBAL_CONSTANTS
USE MESH_POINTERS
IMPLICIT NONE
PRIVATE
CHARACTER(255), PARAMETER :: wallid='$Id: wall.f90 9977 2012-02-03 19:52:00Z mcgratta $'
CHARACTER(255), PARAMETER :: wallrev='$Revision: 9977 $'
CHARACTER(255), PARAMETER :: walldate='$Date: 2012-02-03 11:52:00 -0800 (Fri, 03 Feb 2012) $'
PUBLIC WALL_BC,GET_REV_wall
CONTAINS
SUBROUTINE WALL_BC(T,NM)
! This is the main control routine for this module
USE COMP_FUNCTIONS, ONLY: SECOND
REAL(EB) :: TNOW
REAL(EB), INTENT(IN) :: T
INTEGER, INTENT(IN) :: NM
IF (EVACUATION_ONLY(NM)) RETURN
TNOW=SECOND()
CALL POINT_TO_MESH(NM)
CALL DIFFUSIVITY_BC
CALL THERMAL_BC(T)
IF (ANY(SPECIES_MIXTURE%DEPOSITING)) CALL CALC_DEPOSITION(NM)
CALL SPECIES_BC(T)
CALL DENSITY_BC
IF (HVAC_SOLVE) CALL HVAC_BC
TUSED(6,NM)=TUSED(6,NM)+SECOND()-TNOW
END SUBROUTINE WALL_BC
SUBROUTINE DIFFUSIVITY_BC
! Calculate the term RHODW=RHO*D at the wall
INTEGER :: IW,N,ITMP
REAL(EB), POINTER, DIMENSION(:,:,:) :: RHOP
TYPE(WALL_TYPE), POINTER :: WC=>NULL()
IF (PREDICTOR) RHOP => RHOS
IF (CORRECTOR) RHOP => RHO
WALL_LOOP: DO IW=1,N_EXTERNAL_WALL_CELLS+N_INTERNAL_WALL_CELLS
WC=>WALL(IW)
IF (WC%BOUNDARY_TYPE==NULL_BOUNDARY .OR. WC%BOUNDARY_TYPE==INTERPOLATED_BOUNDARY) CYCLE WALL_LOOP
DO N=1,N_TRACKED_SPECIES
IF (LES) THEN
WC%RHODW(N) = MU(WC%IIG,WC%JJG,WC%KKG)*RSC*WC%RHO_F/RHOP(WC%IIG,WC%JJG,WC%KKG)
ELSE
ITMP = MIN(4999,NINT(WC%TMP_F))
WC%RHODW(N) = WC%RHO_F*D_Z(ITMP,N)
ENDIF
ENDDO
ENDDO WALL_LOOP
END SUBROUTINE DIFFUSIVITY_BC
SUBROUTINE THERMAL_BC(T)
! Thermal boundary conditions for adiabatic, fixed temperature, fixed flux and interpolated boundaries.
! One dimensional heat transfer and pyrolysis is done in PYROLYSIS, which is called at the end of this routine.
USE MATH_FUNCTIONS, ONLY: EVALUATE_RAMP
USE PHYSICAL_FUNCTIONS, ONLY : GET_SPECIFIC_GAS_CONSTANT
REAL(EB) :: DT_BC,T,TSI,TMP_G,RHO_G,DTMP,TMP_OTHER,RAMP_FACTOR,QNET,FDERIV,UN,ARO,UWO,QEXTRA,RSUM_W, &
ZZ_G_ALL(MAX_SPECIES),RHO_ZZ_F(MAX_SPECIES),ZZ_GET(0:N_TRACKED_SPECIES)
INTEGER :: IOR,II,JJ,KK,SURF_INDEX,IIG,JJG,KKG,IW,NOM,IIO,JJO,KKO
REAL(EB), POINTER, DIMENSION(:,:,:) :: UU=>NULL(),VV=>NULL(),WW=>NULL(),RHOP=>NULL(),OM_RHOP=>NULL()
REAL(EB), POINTER, DIMENSION(:,:,:,:) :: ZZP=>NULL(),OM_ZZP=>NULL()
LOGICAL :: INFLOW
TYPE (SURFACE_TYPE), POINTER :: SF=>NULL()
TYPE (VENTS_TYPE), POINTER :: VT=>NULL()
TYPE (OMESH_TYPE), POINTER :: OM=>NULL()
TYPE (MESH_TYPE), POINTER :: MM=>NULL()
REAL(EB), POINTER, DIMENSION(:,:) :: PBAR_P=>NULL()
! unstrurctured geom (experimental)
INTEGER :: TRI_INDEX,I,J,K
REAL(EB) :: NXNY_REAL,NX_REAL,IC_REAL,CELL_VOLUME,CELL_AREA
REAL(EB), POINTER, DIMENSION(:,:,:) :: TMP_SUM=>NULL(),N_CELLS=>NULL()
TYPE(FACET_TYPE), POINTER :: FACE=>NULL()
TYPE(CUTCELL_LINKED_LIST_TYPE), POINTER :: CL=>NULL()
TYPE(WALL_TYPE), POINTER :: WC=>NULL()
IF (VEG_LEVEL_SET_UNCOUPLED) RETURN
IF (PREDICTOR) THEN
UU => US
VV => VS
WW => WS
RHOP => RHOS
ZZP => ZZS
PBAR_P => PBAR_S
ELSE
UU => U
VV => V
WW => W
RHOP => RHO
ZZP => ZZ
PBAR_P => PBAR
ENDIF
! Loop through all boundary cells and apply heat transfer method, except for thermally-thick cells
HEAT_FLUX_LOOP: DO IW=1,N_EXTERNAL_WALL_CELLS+N_INTERNAL_WALL_CELLS+N_VIRTUAL_WALL_CELLS
WC=>WALL(IW)
IF (WC%BOUNDARY_TYPE==NULL_BOUNDARY) CYCLE HEAT_FLUX_LOOP
II = WC%II
JJ = WC%JJ
KK = WC%KK
IIG = WC%IIG
JJG = WC%JJG
KKG = WC%KKG
IOR = WC%IOR
SURF_INDEX = WC%SURF_INDEX
! Consider special boundary conditions
IF (WC%BOUNDARY_TYPE/=SOLID_BOUNDARY) THEN
IF (WC%BOUNDARY_TYPE==INTERPOLATED_BOUNDARY) SURF_INDEX = INTERPOLATED_SURF_INDEX
IF (WC%BOUNDARY_TYPE==OPEN_BOUNDARY) SURF_INDEX = OPEN_SURF_INDEX
IF (WC%BOUNDARY_TYPE==MIRROR_BOUNDARY) SURF_INDEX = MIRROR_SURF_INDEX
ENDIF
! Choose the SURFace type
SF => SURFACE(SURF_INDEX)
! Compute surface temperature, TMP_F, and convective heat flux, QCONF, for various boundary conditions
METHOD_OF_HEAT_TRANSFER: SELECT CASE(SF%THERMAL_BC_INDEX)
CASE (NO_CONVECTION) METHOD_OF_HEAT_TRANSFER
WC%TMP_F = TMP(IIG,JJG,KKG)
CASE (INFLOW_OUTFLOW) METHOD_OF_HEAT_TRANSFER
! Base inflow/outflow decision on velocity component with same predictor/corrector attribute
INFLOW = .FALSE.
SELECT CASE(IOR)
CASE( 1)
UN = UU(II,JJ,KK)
CASE(-1)
UN = -UU(II-1,JJ,KK)
CASE( 2)
UN = VV(II,JJ,KK)
CASE(-2)
UN = -VV(II,JJ-1,KK)
CASE( 3)
UN = WW(II,JJ,KK)
CASE(-3)
UN = -WW(II,JJ,KK-1)
END SELECT
IF (UN>=0._EB) INFLOW = .TRUE.
IF (INFLOW) THEN
WC%TMP_F = TMP_0(KK)
IF (WC%VENT_INDEX>0) THEN
VT => VENTS(WC%VENT_INDEX)
IF (VT%TMP_EXTERIOR>0._EB) WC%TMP_F = VT%TMP_EXTERIOR
ENDIF
IF (N_TRACKED_SPECIES>0) WC%ZZ_F(1:N_TRACKED_SPECIES) = SPECIES_MIXTURE(1:N_TRACKED_SPECIES)%ZZ0
ELSE
WC%TMP_F = TMP(IIG,JJG,KKG)
IF (N_TRACKED_SPECIES>0) WC%ZZ_F(1:N_TRACKED_SPECIES)=ZZP(IIG,JJG,KKG,1:N_TRACKED_SPECIES)
ENDIF
TMP(II,JJ,KK) = 2._EB*WC%TMP_F - TMP(IIG,JJG,KKG)
IF (N_TRACKED_SPECIES>0) ZZP(II,JJ,KK,1:N_TRACKED_SPECIES) = &
2._EB*WC%ZZ_F(1:N_TRACKED_SPECIES)-ZZP(IIG,JJG,KKG,1:N_TRACKED_SPECIES)
WC%QCONF = 2._EB*WC%KW*(TMP(IIG,JJG,KKG)-WC%TMP_F)*WC%RDN
CASE (SPECIFIED_TEMPERATURE) METHOD_OF_HEAT_TRANSFER
TMP_G = TMP(IIG,JJG,KKG)
IF (ABS(WC%TW-T_BEGIN) <= SPACING(WC%TW) .AND. SF%RAMP_INDEX(TIME_TEMP)>=1) THEN
TSI = T
ELSE
TSI = T - WC%TW
ENDIF
IF (WC%UW<=0._EB) THEN
IF (SF%TMP_FRONT>0._EB) THEN
WC%TMP_F = TMP_0(KK) + EVALUATE_RAMP(TSI,SF%TAU(TIME_TEMP),SF%RAMP_INDEX(TIME_TEMP))*(SF%TMP_FRONT-TMP_0(KK))
ELSE
WC%TMP_F = TMP_0(KK)
ENDIF
ELSE
WC%TMP_F = TMP_G ! If gas is being drawn from the domain, set the boundary temperature to the gas temperature
ENDIF
DTMP = TMP_G - WC%TMP_F
WC%HEAT_TRANS_COEF = HEAT_TRANSFER_COEFFICIENT(IW,IIG,JJG,KKG,IOR,DTMP,SF%H_FIXED,SF%GEOMETRY,SF%CONV_LENGTH)
WC%QCONF = WC%HEAT_TRANS_COEF*DTMP
CASE (NET_FLUX_BC) METHOD_OF_HEAT_TRANSFER
IF (ABS(WC%TW-T_BEGIN)<= SPACING(WC%TW ) .AND. SF%RAMP_INDEX(TIME_HEAT)>=1) THEN
TSI = T
ELSE
TSI = T - WC%TW
ENDIF
TMP_G = TMP(IIG,JJG,KKG)
TMP_OTHER = WC%TMP_F
RAMP_FACTOR = EVALUATE_RAMP(TSI,SF%TAU(TIME_HEAT),SF%RAMP_INDEX(TIME_HEAT))
QNET = -RAMP_FACTOR*SF%NET_HEAT_FLUX*WC%AREA_ADJUST
ADLOOP: DO
DTMP = TMP_G - TMP_OTHER
IF (ABS(QNET) > 0._EB .AND. ABS(DTMP) <ZERO_P) DTMP=1._EB
WC%HEAT_TRANS_COEF = HEAT_TRANSFER_COEFFICIENT(IW,IIG,JJG,KKG,IOR,DTMP,SF%H_FIXED,SF%GEOMETRY,SF%CONV_LENGTH)
IF (RADIATION) THEN
QEXTRA = WC%HEAT_TRANS_COEF*DTMP + WC%QRADIN - WC%E_WALL * SIGMA * TMP_OTHER ** 4 - QNET
FDERIV = -WC%HEAT_TRANS_COEF - 4._EB * WC%E_WALL * SIGMA * TMP_OTHER ** 3
ELSE
QEXTRA = WC%HEAT_TRANS_COEF*DTMP - QNET
FDERIV = -WC%HEAT_TRANS_COEF
ENDIF
IF (ABS(FDERIV) > ZERO_P) TMP_OTHER = TMP_OTHER - QEXTRA / FDERIV
IF (ABS(TMP_OTHER - WC%TMP_F) / WC%TMP_F < 0.0001) THEN
WC%TMP_F = TMP_OTHER
EXIT ADLOOP
ELSE
WC%TMP_F = TMP_OTHER
CYCLE ADLOOP
ENDIF
ENDDO ADLOOP
DTMP = TMP_G - WC%TMP_F
WC%HEAT_TRANS_COEF = HEAT_TRANSFER_COEFFICIENT(IW,IIG,JJG,KKG,IOR,DTMP,SF%H_FIXED,SF%GEOMETRY,SF%CONV_LENGTH)
WC%QCONF = WC%HEAT_TRANS_COEF*DTMP
CASE (CONVECTIVE_FLUX_BC) METHOD_OF_HEAT_TRANSFER
IF (ABS(WC%TW-T_BEGIN) <= SPACING(WC%TW) .AND. SF%RAMP_INDEX(TIME_HEAT)>=1) THEN
TSI = T
ELSE
TSI = T - WC%TW
ENDIF
RAMP_FACTOR = EVALUATE_RAMP(TSI,SF%TAU(TIME_HEAT),SF%RAMP_INDEX(TIME_HEAT))
IF (SF%TMP_FRONT>0._EB) THEN
WC%TMP_F = TMPA + RAMP_FACTOR*(SF%TMP_FRONT-TMPA)
ELSE
WC%TMP_F = TMP_0(KK)
ENDIF
WC%QCONF = -RAMP_FACTOR*SF%CONVECTIVE_HEAT_FLUX*WC%AREA_ADJUST
CASE (INTERPOLATED_BC) METHOD_OF_HEAT_TRANSFER
NOM = WC%NOM
OM => OMESH(NOM)
IF (PREDICTOR) THEN
OM_RHOP => OM%RHOS
IF (N_TRACKED_SPECIES>0) OM_ZZP => OM%ZZS
ELSE
OM_RHOP => OM%RHO
IF (N_TRACKED_SPECIES>0) OM_ZZP => OM%ZZ
ENDIF
MM => MESHES(NOM)
RHO_G = RHOP(IIG,JJG,KKG)
WC%RHO_F = RHO_G ! Initialize face value of RHO with RHO_G
IF (N_TRACKED_SPECIES>0) THEN
ZZ_G_ALL(1:N_TRACKED_SPECIES) = ZZP(IIG,JJG,KKG,1:N_TRACKED_SPECIES)
RHO_ZZ_F(1:N_TRACKED_SPECIES) = RHO_G*ZZ_G_ALL(1:N_TRACKED_SPECIES) ! Initialize face value of RHO_ZZ with RHO_G*ZZ_G
ENDIF
DO KKO=WC%NOM_IB(3),WC%NOM_IB(6)
DO JJO=WC%NOM_IB(2),WC%NOM_IB(5)
DO IIO=WC%NOM_IB(1),WC%NOM_IB(4)
SELECT CASE(IOR)
CASE( 1)
ARO = MIN(1._EB , RDY(JJ)*RDZ(KK)*MM%DY(JJO)*MM%DZ(KKO)) * 2._EB*DX(II)/(MM%DX(IIO)+DX(II))
UWO = -OM%U(IIO,JJO,KKO)
CASE(-1)
ARO = MIN(1._EB , RDY(JJ)*RDZ(KK)*MM%DY(JJO)*MM%DZ(KKO)) * 2._EB*DX(II)/(MM%DX(IIO)+DX(II))
UWO = OM%U(IIO-1,JJO,KKO)
CASE( 2)
ARO = MIN(1._EB , RDX(II)*RDZ(KK)*MM%DX(IIO)*MM%DZ(KKO)) * 2._EB*DY(JJ)/(MM%DY(JJO)+DY(JJ))
UWO = -OM%V(IIO,JJO,KKO)
CASE(-2)
ARO = MIN(1._EB , RDX(II)*RDZ(KK)*MM%DX(IIO)*MM%DZ(KKO)) * 2._EB*DY(JJ)/(MM%DY(JJO)+DY(JJ))
UWO = OM%V(IIO,JJO-1,KKO)
CASE( 3)
ARO = MIN(1._EB , RDX(II)*RDY(JJ)*MM%DX(IIO)*MM%DY(JJO)) * 2._EB*DZ(KK)/(MM%DZ(KKO)+DZ(KK))
UWO = -OM%W(IIO,JJO,KKO)
CASE(-3)
ARO = MIN(1._EB , RDX(II)*RDY(JJ)*MM%DX(IIO)*MM%DY(JJO)) * 2._EB*DZ(KK)/(MM%DZ(KKO)+DZ(KK))
UWO = OM%W(IIO,JJO,KKO-1)
END SELECT
WC%RHO_F = WC%RHO_F + 0.5_EB*ARO*(OM_RHOP(IIO,JJO,KKO)-RHO_G)
IF (N_TRACKED_SPECIES>0) RHO_ZZ_F(1:N_TRACKED_SPECIES) = RHO_ZZ_F(1:N_TRACKED_SPECIES) + &
0.5_EB*ARO*(OM_RHOP(IIO,JJO,KKO)*OM_ZZP(IIO,JJO,KKO,1:N_TRACKED_SPECIES)-RHO_G*ZZ_G_ALL(1:N_TRACKED_SPECIES))
ENDDO
ENDDO
ENDDO
RHOP(II,JJ,KK) = MIN( RHOMAX , MAX( RHOMIN , 2._EB*WC%RHO_F-RHO_G ))
IF (N_TRACKED_SPECIES==0) THEN
TMP(II,JJ,KK) = PBAR_P(KK,WC%PRESSURE_ZONE_WALL)/(SPECIES_MIXTURE(0)%RCON*RHOP(II,JJ,KK))
ELSE
WC%ZZ_F(1:N_TRACKED_SPECIES) = RHO_ZZ_F(1:N_TRACKED_SPECIES)/WC%RHO_F
ZZP(II,JJ,KK,1:N_TRACKED_SPECIES) = (2._EB*RHO_ZZ_F(1:N_TRACKED_SPECIES) - RHO_G*ZZ_G_ALL(1:N_TRACKED_SPECIES))/ &
RHOP(II,JJ,KK)
ZZ_GET(1:N_TRACKED_SPECIES) = MAX(0._EB,ZZP(II,JJ,KK,1:N_TRACKED_SPECIES))
CALL GET_SPECIFIC_GAS_CONSTANT(ZZ_GET,RSUM_W)
TMP(II,JJ,KK) = PBAR_P(KK,WC%PRESSURE_ZONE_WALL)/(RSUM_W*RHOP(II,JJ,KK))
ENDIF
WC%QCONF = 0._EB
WC%TMP_F = 0.5*(TMP(IIG,JJG,KKG)+TMP(II,JJ,KK))
END SELECT METHOD_OF_HEAT_TRANSFER
ENDDO HEAT_FLUX_LOOP
! For thermally-thick boundary conditions, call the routine PYROLYSIS
IF (CORRECTOR) THEN
WALL_COUNTER = WALL_COUNTER + 1
IF (WALL_COUNTER==WALL_INCREMENT) THEN
DT_BC = T - BC_CLOCK
BC_CLOCK = T
CALL PYROLYSIS(T,DT_BC)
WALL_COUNTER = 0
ENDIF
ENDIF
! Apply boundary conditions from unstructured geometry (under construction)
GEOM_IF: IF (N_FACE>0) THEN
NXNY_REAL = REAL(IBAR*JBAR,EB)+1.E-10_EB
NX_REAL = REAL(IBAR,EB)+1.E-10_EB
N_CELLS => WORK1
N_CELLS = 0._EB
TMP_SUM => WORK2
TMP_SUM = 0._EB
GEOM_HEAT_FLUX_LOOP: DO TRI_INDEX=1,N_FACE
FACE=>FACET(TRI_INDEX)
CL=>FACE%CUTCELL_LIST
SF=>SURFACE(FACE%SURF_INDEX)
CUTCELL_LOOP: DO
IF ( .NOT. ASSOCIATED(CL) ) EXIT ! if the next index does not exist, exit the loop
IC_REAL = REAL(CL%INDEX,EB)
K = CEILING(IC_REAL/NXNY_REAL)
J = CEILING((IC_REAL-(K-1)*NXNY_REAL)/NX_REAL)
I = NINT(IC_REAL-(K-1)*NXNY_REAL-(J-1)*NX_REAL)
CELL_VOLUME = DX(I)*DY(J)*DZ(K)
CELL_AREA = CELL_VOLUME**TWTH
GEOM_METHOD_OF_HEAT_TRANSFER: SELECT CASE(SF%THERMAL_BC_INDEX)
CASE(SPECIFIED_TEMPERATURE_FROM_FILE)
TMP_SUM(I,J,K) = TMP_SUM(I,J,K) + FACE%TMP_F
N_CELLS(I,J,K) = N_CELLS(I,J,K) + 1._EB
DTMP = TMP_SUM(I,J,K)/N_CELLS(I,J,K) - TMP(I,J,K)
Q(I,J,K) = MIN(FACE%AREA,CELL_AREA)*SF%H_FIXED*DTMP/CELL_VOLUME
CASE(SPECIFIED_TEMPERATURE)
TMP_SUM(I,J,K) = TMP_SUM(I,J,K) + SF%TMP_FRONT
N_CELLS(I,J,K) = N_CELLS(I,J,K) + 1._EB
DTMP = TMP_SUM(I,J,K)/N_CELLS(I,J,K) - TMP(I,J,K)
Q(I,J,K) = MIN(FACE%AREA,CELL_AREA)*SF%H_FIXED*DTMP/CELL_VOLUME
CASE(CONVECTIVE_FLUX_BC)
Q(I,J,K) = MIN(FACE%AREA,CELL_AREA)*SF%CONVECTIVE_HEAT_FLUX/CELL_VOLUME
END SELECT GEOM_METHOD_OF_HEAT_TRANSFER
CL=>CL%NEXT ! point to the next index in the linked list
ENDDO CUTCELL_LOOP
ENDDO GEOM_HEAT_FLUX_LOOP
ENDIF GEOM_IF
END SUBROUTINE THERMAL_BC
SUBROUTINE SPECIES_BC(T)
! Compute the species mass fractions at the boundary, ZZ_F
USE PHYSICAL_FUNCTIONS, ONLY: GET_SPECIFIC_HEAT,GET_AVERAGE_SPECIFIC_HEAT,GET_SPECIFIC_GAS_CONSTANT
USE MATH_FUNCTIONS, ONLY : EVALUATE_RAMP
REAL(EB) :: T,ZZ_G,UN,DD,MFT,TSI,RADIUS,AREA_SCALING,&
RVC,RHO_NEW,M_DOT_PPP,CP,CPBAR,MW_RATIO,H_G,DELTA_H_G,ZZ_GET(0:N_TRACKED_SPECIES),CPBAR2,RSUM_F
INTEGER :: I,SURF_INDEX,IIG,JJG,KKG,IW,II,JJ,KK,N,NS,ITER
TYPE (SURFACE_TYPE), POINTER :: SF=>NULL()
TYPE (PARTICLE_CLASS_TYPE), POINTER :: PC=>NULL()
TYPE (PARTICLE_TYPE), POINTER :: LP=>NULL()
REAL(EB), POINTER, DIMENSION(:,:,:) :: RHOP=>NULL()
REAL(EB), POINTER, DIMENSION(:,:,:,:) :: ZZP=>NULL()
TYPE (WALL_TYPE), POINTER :: WC=>NULL()
IF (VEG_LEVEL_SET_UNCOUPLED) RETURN
IF (PREDICTOR) THEN
RHOP => RHOS
ZZP => ZZS
ELSE
RHOP => RHO
ZZP => ZZ
ENDIF
! Loop through the wall cells, apply mass boundary conditions
WALL_CELL_LOOP: DO IW=1,N_EXTERNAL_WALL_CELLS+N_INTERNAL_WALL_CELLS+N_VIRTUAL_WALL_CELLS
WC => WALL(IW)
IF (WC%BOUNDARY_TYPE==OPEN_BOUNDARY) CYCLE WALL_CELL_LOOP
IF (WC%BOUNDARY_TYPE==NULL_BOUNDARY) CYCLE WALL_CELL_LOOP
IF (WC%BOUNDARY_TYPE==INTERPOLATED_BOUNDARY) CYCLE WALL_CELL_LOOP
SURF_INDEX = WC%SURF_INDEX
! Special cases that over-ride the boundary condition index, SURF_INDEX
IF (WC%BOUNDARY_TYPE==MIRROR_BOUNDARY) SURF_INDEX = MIRROR_SURF_INDEX
! Set the SURFace type
SF => SURFACE(SURF_INDEX)
! Special cases
IF (N_TRACKED_SPECIES==0 .AND. .NOT. SF%SPECIES_BC_INDEX==SPECIFIED_MASS_FLUX) CYCLE WALL_CELL_LOOP
IF (N_TRACKED_SPECIES==0 .AND. SF%SPECIES_BC_INDEX==SPECIFIED_MASS_FLUX .AND. ABS(SF%MASS_FLUX(0))<=ZERO_P) CYCLE WALL_CELL_LOOP
! Get the wall indices
II = WC%II
JJ = WC%JJ
KK = WC%KK
IIG = WC%IIG
JJG = WC%JJG
KKG = WC%KKG
! Check if suppression by water is to be applied
IF (CORRECTOR .AND. SF%E_COEFFICIENT>0._EB) THEN
IF (SUM(WC%LP_MPUA(:))>0._EB .AND. T>WC%TW) THEN
WC%EW = WC%EW + SF%E_COEFFICIENT*SUM(WC%LP_MPUA(:))*DT
ENDIF
ENDIF
! Apply the different species boundary conditions to non-thermally thick solids
METHOD_OF_MASS_TRANSFER: SELECT CASE(SF%SPECIES_BC_INDEX)
CASE (INFLOW_OUTFLOW_MASS_FLUX) METHOD_OF_MASS_TRANSFER
! OPEN boundary species BC is done in THERMAL_BC under INFLOW_OUTFLOW
CASE (NO_MASS_FLUX) METHOD_OF_MASS_TRANSFER
IF (.NOT.SOLID(CELL_INDEX(IIG,JJG,KKG))) WC%ZZ_F(1:N_TRACKED_SPECIES) = ZZP(IIG,JJG,KKG,1:N_TRACKED_SPECIES)
CASE (SPECIFIED_MASS_FRACTION) METHOD_OF_MASS_TRANSFER
IF (ABS(WC%TW-T_BEGIN)< SPACING(WC%TW) .AND. ANY(SF%RAMP_INDEX>=1)) THEN
IF (PREDICTOR) TSI = T + DT
IF (CORRECTOR) TSI = T
ELSE
IF (PREDICTOR) TSI = T + DT - WC%TW
IF (CORRECTOR) TSI = T - WC%TW
ENDIF
IF (WC%UWS<0._EB) THEN
DO N=1,N_TRACKED_SPECIES
WC%ZZ_F(N) = SPECIES_MIXTURE(N)%ZZ0 + EVALUATE_RAMP(TSI,SF%TAU(N),SF%RAMP_INDEX(N))* &
(SF%MASS_FRACTION(N)-SPECIES_MIXTURE(N)%ZZ0)
ENDDO
ELSE
WC%ZZ_F(1:N_TRACKED_SPECIES) = ZZP(IIG,JJG,KKG,1:N_TRACKED_SPECIES)
ENDIF
CASE (SPECIFIED_MASS_FLUX) METHOD_OF_MASS_TRANSFER
! If the current time is before the "activation" time, TW, apply simple BCs and get out
IF (T < WC%TW) THEN
WC%MASSFLUX(0) = 0._EB
IF (N_TRACKED_SPECIES > 0) THEN
IF (.NOT.SOLID(CELL_INDEX(IIG,JJG,KKG))) THEN
WC%ZZ_F(1:N_TRACKED_SPECIES) = ZZP(IIG,JJG,KKG,1:N_TRACKED_SPECIES)
ENDIF
IF (PREDICTOR) WC%UWS = 0._EB
WC%MASSFLUX(1:N_TRACKED_SPECIES) = 0._EB
WC%MASSFLUX_ACTUAL(1:N_TRACKED_SPECIES) = 0._EB
ENDIF
CYCLE WALL_CELL_LOOP
ENDIF
! Zero out the running counter of Mass Flux Total (MFT)
MFT = 0._EB
! If the user has specified the burning rate, evaluate the ramp and other related parameters
SUM_MASSFLUX_LOOP: DO N=0,N_TRACKED_SPECIES
IF (SF%MASS_FLUX(N) > 0._EB) THEN ! Use user-specified ramp-up of mass flux
IF (ABS(WC%TW-T_BEGIN)< SPACING(WC%TW) .AND. SF%RAMP_INDEX(N)>=1) THEN
TSI = T + DT
ELSE
TSI = T + DT - WC%TW
ENDIF
WC%MASSFLUX(N) = EVALUATE_RAMP(TSI,SF%TAU(N),SF%RAMP_INDEX(N))*SF%MASS_FLUX(N)
WC%MASSFLUX_ACTUAL(N) = WC%MASSFLUX(N)
WC%MASSFLUX(N) = SF%ADJUST_BURN_RATE(N)*WC%MASSFLUX(N)
ENDIF
WC%MASSFLUX(N) = WC%MASSFLUX(N)*WC%AREA_ADJUST
MFT = MFT + WC%MASSFLUX(N)
ENDDO SUM_MASSFLUX_LOOP
IF (WC%EW>0._EB) WC%MASSFLUX(:) = WC%MASSFLUX(:)*EXP(-WC%EW)
! Add total consumed mass to various summing arrays
CONSUME_MASS: IF (CORRECTOR .AND. SF%THERMALLY_THICK) THEN
DO N=1,N_TRACKED_SPECIES
OBSTRUCTION(WC%OBST_INDEX)%MASS = OBSTRUCTION(WC%OBST_INDEX)%MASS - WC%MASSFLUX_ACTUAL(N)*DT*WC%AW
ENDDO
ENDIF CONSUME_MASS
! Compute the cell face value of the species mass fraction to get the right mass flux
DO ITER=1,5
UN = MFT/WC%RHO_F
SPECIES_LOOP: DO N=1,N_TRACKED_SPECIES
DD = 2.*WC%RHODW(N)*WC%RDN
ZZ_G = ZZP(IIG,JJG,KKG,N)
WC%ZZ_F(N) = ( WC%MASSFLUX(N) + DD*ZZ_G ) / (DD + UN*WC%RHO_F)
ENDDO SPECIES_LOOP
ZZ_GET(1:N_TRACKED_SPECIES) = MAX(0._EB,WC%ZZ_F(1:N_TRACKED_SPECIES))
CALL GET_SPECIFIC_GAS_CONSTANT(ZZ_GET,RSUM_F)
WC%RHO_F = PBAR_S(KK,WC%PRESSURE_ZONE_WALL)/(RSUM_F*WC%TMP_F)
ENDDO
IF (PREDICTOR) WC%UWS = -UN
IF (CORRECTOR) WC%UW = -UN
END SELECT METHOD_OF_MASS_TRANSFER
! Only set species mass fraction in the ghost cell if it is not solid
IF (IW<=N_EXTERNAL_WALL_CELLS .AND. N_TRACKED_SPECIES > 0 .AND. &
.NOT.SOLID(CELL_INDEX(II,JJ,KK)) .AND. .NOT.SOLID(CELL_INDEX(IIG,JJG,KKG))) &
ZZP(II,JJ,KK,1:N_TRACKED_SPECIES) = 2._EB*WC%ZZ_F(1:N_TRACKED_SPECIES) - ZZP(IIG,JJG,KKG,1:N_TRACKED_SPECIES)
ENDDO WALL_CELL_LOOP
! Add evaporating gases from virtual particles to the mesh using a volumetric source term
IF (VIRTUAL_PARTICLES .AND. CORRECTOR) THEN
PARTICLE_LOOP: DO I=1,NLP
LP => PARTICLE(I)
PC => PARTICLE_CLASS(LP%CLASS)
SURF_INDEX = PC%SURF_INDEX
IF (SURF_INDEX<1) CYCLE PARTICLE_LOOP
SF => SURFACE(SURF_INDEX)
IW = LP%WALL_INDEX
WC=>WALL(IW)
II = WC%II
JJ = WC%JJ
KK = WC%KK
IF (SF%SHRINK) THEN
RADIUS = SUM(WC%LAYER_THICKNESS)
ELSEIF (SF%THERMALLY_THICK) THEN
RADIUS = SF%THICKNESS
ELSE
RADIUS = SF%RADIUS
ENDIF
LP%R = RADIUS
IF (ABS(RADIUS)<ZERO_P) CYCLE PARTICLE_LOOP
AREA_SCALING = 1._EB
SELECT CASE(SF%GEOMETRY)
CASE(SURF_CARTESIAN)
WC%AW = 2._EB*SF%LENGTH*SF%WIDTH
CASE(SURF_CYLINDRICAL)
WC%AW = TWOPI*RADIUS*SF%LENGTH
IF (SF%THERMAL_BC_INDEX == THERMALLY_THICK) AREA_SCALING = (SF%THICKNESS/RADIUS)
CASE(SURF_SPHERICAL)
WC%AW = 4._EB*PI*RADIUS**2
IF (SF%THERMAL_BC_INDEX == THERMALLY_THICK) AREA_SCALING = (SF%THICKNESS/RADIUS)**2
END SELECT
! In PYROLYSIS, all the mass fluxes are normalized by a virtual area based on the INITIAL radius.
! Here, correct the mass flux using the CURRENT radius. Also, multiply by LP%PWT to account for split particles
AREA_SCALING = AREA_SCALING*LP%PWT
WC%MASSFLUX(0:N_TRACKED_SPECIES) = WC%MASSFLUX(0:N_TRACKED_SPECIES) *AREA_SCALING
WC%MASSFLUX_ACTUAL(0:N_TRACKED_SPECIES) = WC%MASSFLUX_ACTUAL(0:N_TRACKED_SPECIES)*AREA_SCALING
RVC = RDX(II)*RRN(II)*RDY(JJ)*RDZ(KK)
IF (N_TRACKED_SPECIES > 0) ZZ_GET(1:N_TRACKED_SPECIES) = ZZP(II,JJ,KK,1:N_TRACKED_SPECIES)
CALL GET_SPECIFIC_HEAT(ZZ_GET,CP,TMP(II,JJ,KK))
H_G = CP*TMP(II,JJ,KK)
DO NS=0,N_TRACKED_SPECIES
IF (ABS(WC%MASSFLUX(NS))<=ZERO_P) CYCLE
IF (N_TRACKED_SPECIES>0) THEN
MW_RATIO = SPECIES_MIXTURE(NS)%RCON/RSUM(II,JJ,KK)
ELSE
MW_RATIO = 1._EB
ENDIF
M_DOT_PPP = WC%MASSFLUX(NS)*WC%AW*RVC
ZZ_GET=0._EB
IF (NS>0) ZZ_GET(NS)=1._EB
CALL GET_AVERAGE_SPECIFIC_HEAT(ZZ_GET,CPBAR,TMP(II,JJ,KK))
CALL GET_AVERAGE_SPECIFIC_HEAT(ZZ_GET,CPBAR2,WC%TMP_F)
DELTA_H_G = CPBAR2*WC%TMP_F-CPBAR*TMP(II,JJ,KK)
D_LAGRANGIAN(II,JJ,KK) = D_LAGRANGIAN(II,JJ,KK) + LP%PWT*M_DOT_PPP*(MW_RATIO + DELTA_H_G/H_G)/RHO(II,JJ,KK)
RHO_NEW = RHO(II,JJ,KK) + M_DOT_PPP*DT
IF (NS>0) THEN
ZZP(II,JJ,KK,NS) = (RHO(II,JJ,KK)*ZZP(II,JJ,KK,NS) + M_DOT_PPP*DT)/RHO_NEW
ELSE
IF (N_TRACKED_SPECIES>0) ZZP(II,JJ,KK,1:N_TRACKED_SPECIES) = RHO(II,JJ,KK)*ZZP(II,JJ,KK,1:N_TRACKED_SPECIES)/RHO_NEW
ENDIF
RHO(II,JJ,KK) = RHO_NEW
ENDDO
D_LAGRANGIAN(II,JJ,KK) = D_LAGRANGIAN(II,JJ,KK) - WC%QCONF*WC%AW*RVC/(RHO(II,JJ,KK)*H_G) * LP%PWT
ENDDO PARTICLE_LOOP
ENDIF
END SUBROUTINE SPECIES_BC
SUBROUTINE DENSITY_BC
! Compute density at wall from wall temperatures and mass fractions
USE PHYSICAL_FUNCTIONS, ONLY : GET_SPECIFIC_GAS_CONSTANT
REAL(EB) :: ZZ_GET(0:N_TRACKED_SPECIES),RSUM_F
INTEGER :: IW,II,JJ,KK,IOR
REAL(EB), POINTER, DIMENSION(:,:) :: PBAR_P=>NULL()
REAL(EB), POINTER, DIMENSION(:,:,:) :: RHOP=>NULL()
REAL(EB), POINTER, DIMENSION(:,:,:,:) :: ZZP=>NULL()
TYPE (WALL_TYPE), POINTER :: WC=>NULL()
IF (VEG_LEVEL_SET_UNCOUPLED) RETURN
IF (PREDICTOR) THEN
PBAR_P => PBAR_S
RHOP => RHOS
ZZP => ZZS
ELSE
PBAR_P => PBAR
RHOP => RHO
ZZP => ZZ
ENDIF
WALL_CELL_LOOP: DO IW=1,N_EXTERNAL_WALL_CELLS+N_INTERNAL_WALL_CELLS
WC => WALL(IW)
IF (WC%BOUNDARY_TYPE==NULL_BOUNDARY) CYCLE WALL_CELL_LOOP
II = WC%II
JJ = WC%JJ
KK = WC%KK
IOR = WC%IOR
! Determine ghost cell value of RSUM=R0*Sum(Y_i/M_i)
IF (N_TRACKED_SPECIES>0) ZZ_GET(1:N_TRACKED_SPECIES) = MAX(0._EB,WC%ZZ_F(1:N_TRACKED_SPECIES))
CALL GET_SPECIFIC_GAS_CONSTANT(ZZ_GET,RSUM_F)
! Compute density at boundary cell face
IF (WC%BOUNDARY_TYPE/=INTERPOLATED_BOUNDARY) WC%RHO_F = PBAR_P(KK,WC%PRESSURE_ZONE_WALL)/(RSUM_F*WC%TMP_F)
! Set ghost cell values for open and interpolated boundaries
IF (WC%BOUNDARY_TYPE==OPEN_BOUNDARY .OR. WC%BOUNDARY_TYPE==INTERPOLATED_BOUNDARY) THEN
IF (N_TRACKED_SPECIES>0) ZZ_GET(1:N_TRACKED_SPECIES) = MAX(0._EB,ZZP(II,JJ,KK,1:N_TRACKED_SPECIES))
CALL GET_SPECIFIC_GAS_CONSTANT(ZZ_GET,RSUM(II,JJ,KK))
RHOP(II,JJ,KK) = PBAR_P(KK,WC%PRESSURE_ZONE_WALL)/(RSUM(II,JJ,KK)*TMP(II,JJ,KK))
ENDIF
ENDDO WALL_CELL_LOOP
END SUBROUTINE DENSITY_BC
SUBROUTINE HVAC_BC
! Compute density at wall from wall temperatures and mass fractions
USE HVAC_ROUTINES, ONLY : NODE_CP,NODE_TMP,NODE_RHO,DUCT_U,LEAK_PATH,NODE_ZZ
USE PHYSICAL_FUNCTIONS, ONLY : GET_SPECIFIC_GAS_CONSTANT
REAL(EB) :: ZZ_GET(0:N_TRACKED_SPECIES),ZZ_G(1:N_TRACKED_SPECIES),UN,MFT,DD,RSUM_F
REAL(EB) :: RHO_0,TMP_0,ZZ_0(1:N_TRACKED_SPECIES),UW_0,ZZ_ERR
INTEGER :: IIG,JJG,KKG,IW,II,JJ,KK,N,SURF_INDEX,COUNTER,DN,DN2,DU,IZ1,IZ2
LOGICAL :: ITER = .FALSE.
REAL(EB), POINTER, DIMENSION(:,:) :: PBAR_P=>NULL()
REAL(EB), POINTER, DIMENSION(:,:,:) :: RHOP=>NULL()
REAL(EB), POINTER, DIMENSION(:,:,:,:) :: ZZP=>NULL()
TYPE (SURFACE_TYPE), POINTER :: SF=>NULL()
TYPE (WALL_TYPE),POINTER :: WC=>NULL()
IF (PREDICTOR) THEN
RHOP => RHOS
IF (N_TRACKED_SPECIES > 0) ZZP => ZZS
PBAR_P => PBAR_S
ELSE
RHOP => RHO
IF (N_TRACKED_SPECIES > 0) ZZP => ZZ
PBAR_P => PBAR
ENDIF
! Loop over all internal and external wall cells
WALL_CELL_LOOP: DO IW=1,N_EXTERNAL_WALL_CELLS+N_INTERNAL_WALL_CELLS
WC=>WALL(IW)
SURF_INDEX = WC%SURF_INDEX
SF => SURFACE(SURF_INDEX)
IF (SF%SPECIES_BC_INDEX/=HVAC_BOUNDARY .AND. SF%THERMAL_BC_INDEX/=HVAC_BOUNDARY .AND. &
.NOT. ANY(SF%LEAK_PATH>0._EB)) CYCLE WALL_CELL_LOOP
II = WC%II
JJ = WC%JJ
KK = WC%KK
IIG = WC%IIG
JJG = WC%JJG
KKG = WC%KKG
COUNTER = 0
! Compute R*Sum(Y_i/W_i) at the wall
IF (N_TRACKED_SPECIES>0) THEN
ZZ_G = ZZP(IIG,JJG,KKG,:)
ZZ_GET(1:N_TRACKED_SPECIES) = MAX(0._EB,WC%ZZ_F(1:N_TRACKED_SPECIES))
CALL GET_SPECIFIC_GAS_CONSTANT(ZZ_GET,RSUM_F)
ELSE
RSUM_F = RSUM0
ENDIF
IF (WC%VENT_INDEX > 0 .AND. .NOT. ANY(SF%LEAK_PATH>0._EB)) THEN
IF (VENTS(WC%VENT_INDEX)%NODE_INDEX > 0) THEN
DN=VENTS(WC%VENT_INDEX)%NODE_INDEX
DU=DUCTNODE(DN)%DUCT_INDEX(1)
MFT = -DUCTNODE(DN)%DIR(1)*DUCT_U(DU,1)*NODE_RHO(DN,1)/VENTS(WC%VENT_INDEX)%FDS_AREA
IF (DUCTNODE(DN)%DIR(1)*DUCT_U(DU,1) > 0._EB) THEN
WC%TMP_F = NODE_TMP(DN,1)
ELSE
WC%TMP_F = TMP(IIG,JJG,KKG)
ENDIF
ENDIF
WC%HEAT_TRANS_COEF = 0._EB
WC%QCONF = 0._EB
ELSE !leakage
IF (WC%PRESSURE_ZONE_WALL==SF%LEAK_PATH(1)) THEN
IZ1 = SF%LEAK_PATH(1)
IZ2 = SF%LEAK_PATH(2)
ELSE
IZ1 = SF%LEAK_PATH(2)
IZ2 = SF%LEAK_PATH(1)
ENDIF
DU = LEAK_PATH(MIN(IZ1,IZ2),MAX(IZ1,IZ2))
IF (DUCTNODE(DUCT(DU)%NODE_INDEX(1))%ZONE_INDEX==WC%PRESSURE_ZONE_WALL) THEN
DN2=DUCT(DU)%NODE_INDEX(1)
ELSE
DN2=DUCT(DU)%NODE_INDEX(2)
ENDIF
IF(DUCT_U(DU,1) > 0._EB) THEN
DN=DUCT(DU)%NODE_INDEX(1)
ELSE
DN=DUCT(DU)%NODE_INDEX(2)
ENDIF
MFT = -DUCTNODE(DN2)%DIR(1)*DUCT_U(DU,1)*NODE_RHO(DN,1)/FDS_LEAK_AREA(IZ1,IZ2,1)
ENDIF
WC%RHO_F = PBAR_P(KK,WC%PRESSURE_ZONE_WALL)/(RSUM_F*WC%TMP_F)
UN = -MFT/WC%RHO_F
! Iterate to get the appropriate normal velocity and density
SPECIES_IF_1: IF (N_TRACKED_SPECIES==0) THEN
WC%MASSFLUX(0) = MFT
ELSE SPECIES_IF_1
IF (UN > 0._EB) THEN
WC%MASSFLUX(:) = 0._EB
WC%MASSFLUX(1:N_TRACKED_SPECIES) = -NODE_ZZ(DN,1:N_TRACKED_SPECIES,1)*MFT
WC%MASSFLUX(0) = -MFT - SUM(WC%MASSFLUX(1:N_TRACKED_SPECIES))
ENDIF
ENDIF SPECIES_IF_1
ITER = .TRUE.
TMP_0 = WC%TMP_F
DO WHILE (ITER)
ITER = .FALSE.
RHO_0 = WC%RHO_F
ZZ_0 = WC%ZZ_F
UW_0 = -UN
UN = -MFT/WC%RHO_F
IF (PREDICTOR) WC%UWS = -UN
IF (CORRECTOR) WC%UW = -UN
SPECIES_IF: IF (N_TRACKED_SPECIES==0) THEN
RSUM_F = RSUM0
ELSE SPECIES_IF
ZZ_ERR = 0._EB
IF (UN <= 0._EB) THEN
WC%ZZ_F = ZZ_G
ELSE
DO N=1,N_TRACKED_SPECIES
DD = 2._EB*WC%RHODW(N)*WC%RDN
WC%ZZ_F(N) = ( WC%MASSFLUX(N) + DD*ZZ_G(N) ) / ( DD + UN*WC%RHO_F )
IF (ZZ_0(N) > 1.E-10_EB) ZZ_ERR = MAX(ZZ_ERR,ABS(WC%ZZ_F(N)-ZZ_0(N))/ZZ_0(N))
ENDDO
ENDIF
IF (COUNTER > 5) WC%ZZ_F(:) = 0.4_EB*WC%ZZ_F(:)+0.6_EB*ZZ_0(:)
ZZ_GET(1:N_TRACKED_SPECIES) = MAX(0._EB,WC%ZZ_F(1:N_TRACKED_SPECIES))
CALL GET_SPECIFIC_GAS_CONSTANT(ZZ_GET,RSUM_F)
IF (ZZ_ERR > 1.E-6_EB) ITER = .TRUE.
ENDIF SPECIES_IF
WC%RHO_F = PBAR_P(KK,WC%PRESSURE_ZONE_WALL)/(RSUM_F*WC%TMP_F)
! Decide to continue iterating
IF (ABS(RHO_0 - WC%RHO_F)/RHO_0 > 1.E-6_EB ) ITER = .TRUE.
IF (ABS(UW_0)>0._EB) THEN
IF (ABS(UW_0 + UN)/ABS(UW_0) > 1.E-6_EB) ITER = .TRUE.
ENDIF
COUNTER = COUNTER + 1
IF (COUNTER > 20) ITER = .FALSE.
ENDDO
ENDDO WALL_CELL_LOOP
END SUBROUTINE HVAC_BC
SUBROUTINE PYROLYSIS(T,DT_BC)
! Loop through all the boundary cells that require a 1-D heat transfer calc
USE PHYSICAL_FUNCTIONS, ONLY: GET_MOLECULAR_WEIGHT
USE GEOMETRY_FUNCTIONS
USE MATH_FUNCTIONS, ONLY: EVALUATE_RAMP
REAL(EB) :: DTMP,QNETF,QDXKF,QDXKB,RR,TMP_G,T,RFACF,RFACB,RFACF2,RFACB2,PPCLAUS,PPSURF,DT_BC, &
DXKF,DXKB,REACTION_RATE,QRADINB,RFLUX_UP,RFLUX_DOWN,E_WALLB, &
HVRG,Z_MF_G, RSUM_G, MFLUX, MFLUX_S, VOLSUM,ZPRSUM, &
DXF, DXB,HTCB,Q_WATER_F,Q_WATER_B,TMP_F_OLD, DX_GRID, RHO_S0,DT2_BC,TOLERANCE,C_S_ADJUST_UNITS,&
MW_G,H_MASS,X_G,X_W,D_AIR,MU_AIR,U2,V2,W2,RE_L,SC_AIR,SH_FAC_WALL,SHERWOOD,VELCON,RHO_G,TMP_BACK,TMP_WGT
!REAL(EB) :: ZZ_S,RSUM_S
INTEGER :: SURF_INDEX,IIG,JJG,KKG,IIB,JJB,KKB,IWB,NWP,I,J,NR,NN,NNN,NL,II,JJ,KK,IW,IOR,N,I_OBST,NS,ITMP
REAL(EB) :: SMALLEST_CELL_SIZE(MAX_LAYERS),THICKNESS,ZZ_GET(0:N_TRACKED_SPECIES)
REAL(EB),ALLOCATABLE,DIMENSION(:) :: TMP_W_NEW
REAL(EB),ALLOCATABLE,DIMENSION(:,:) :: INT_WGT
REAL(EB), POINTER, DIMENSION(:,:,:) :: UU=>NULL(),VV=>NULL(),WW=>NULL(),RHOG=>NULL()
REAL(EB), POINTER, DIMENSION(:,:) :: PBARP
INTEGER :: N_LAYER_CELLS_NEW(MAX_LAYERS), NWP_NEW,I_GRAD,STEPCOUNT,SMIX_PTR
LOGICAL :: POINT_SHRINK, RECOMPUTE,ITERATE,E_FOUND
TYPE (WALL_TYPE), POINTER :: WC=>NULL()
TYPE (SURFACE_TYPE), POINTER :: SF=>NULL()
TYPE (MATERIAL_TYPE), POINTER :: ML=>NULL()
! Liquid evap
IF (PREDICTOR) THEN
UU => U
VV => V
WW => W
RHOG => RHO
PBARP => PBAR
ELSE
UU => US
VV => VS
WW => WS
RHOG => RHOS
PBARP => PBAR_S
ENDIF
SC_AIR = 0.6_EB ! NU_AIR/D_AIR (Incropera & DeWitt, Chap 7, External Flow)
SH_FAC_WALL = 0.037_EB*SC_AIR**ONTH
! Special adjustment of specific heat for steady state applications
C_S_ADJUST_UNITS = 1000._EB/TIME_SHRINK_FACTOR
! Loop through the thermally-thick wall cells
WALL_CELL_LOOP: DO IW=1,N_EXTERNAL_WALL_CELLS+N_INTERNAL_WALL_CELLS+N_VIRTUAL_WALL_CELLS
WC => WALL(IW)
IF (WC%BOUNDARY_TYPE/=SOLID_BOUNDARY .AND. WC%BOUNDARY_TYPE/=VIRTUAL_BOUNDARY) CYCLE WALL_CELL_LOOP
SURF_INDEX = WC%SURF_INDEX
SF => SURFACE(SURF_INDEX)
IF (SF%THERMAL_BC_INDEX /= THERMALLY_THICK) CYCLE WALL_CELL_LOOP
II = WC%II
JJ = WC%JJ
KK = WC%KK
IOR = WC%IOR
IIG = WC%IIG
JJG = WC%JJG
KKG = WC%KKG
IF (WC%BURNAWAY) THEN
WC%MASSFLUX = 0._EB
CYCLE WALL_CELL_LOOP
ENDIF
SELECT CASE(SF%GEOMETRY)
CASE(SURF_CARTESIAN)
I_GRAD = 0
CASE(SURF_CYLINDRICAL)
I_GRAD = 1
CASE(SURF_SPHERICAL)
I_GRAD = 2
END SELECT
! Compute convective heat flux at the surface
TMP_G = TMP(IIG,JJG,KKG)
IF (ASSUMED_GAS_TEMPERATURE > 0._EB) TMP_G = ASSUMED_GAS_TEMPERATURE ! This is just for diagnostic calcs
TMP_F_OLD = WC%TMP_F
DTMP = TMP_G - TMP_F_OLD
WC%HEAT_TRANS_COEF = HEAT_TRANSFER_COEFFICIENT(IW,IIG,JJG,KKG,IOR,DTMP,SF%H_FIXED,SF%GEOMETRY,SF%CONV_LENGTH)
WC%QCONF = WC%HEAT_TRANS_COEF*DTMP
! Compute back side emissivity
E_WALLB = SF%EMISSIVITY_BACK
IF (E_WALLB < 0._EB .AND. SF%BACKING /= INSULATED) THEN
E_WALLB = 0._EB
VOLSUM = 0._EB
IF (SF%SHRINK) THEN
NWP = SUM(WC%N_LAYER_CELLS)
ELSE
NWP = SF%N_CELLS
ENDIF
DO N=1,SF%N_MATL
IF (WC%RHO_S(NWP,N)<=ZERO_P) CYCLE
ML => MATERIAL(SF%MATL_INDEX(N))
VOLSUM = VOLSUM + WC%RHO_S(NWP,N)/ML%RHO_S
E_WALLB = E_WALLB + WC%RHO_S(NWP,N)*ML%EMISSIVITY/ML%RHO_S
ENDDO
IF (VOLSUM > 0._EB) E_WALLB = E_WALLB/VOLSUM
ENDIF
! Get heat losses from convection and radiation out of back of surface
SELECT CASE(SF%BACKING)
CASE(VOID) ! Non-insulated backing to an ambient void
IF (SF%TMP_BACK>0._EB) THEN
TMP_BACK = SF%TMP_BACK
ELSE
TMP_BACK = TMP_0(KK)
ENDIF
DTMP = TMP_BACK - WC%TMP_B
HTCB = HEAT_TRANSFER_COEFFICIENT(IW,-1,-1,-1,IOR,DTMP,SF%H_FIXED,SF%GEOMETRY,SF%CONV_LENGTH)
QRADINB = E_WALLB*SIGMA*TMP_BACK**4
Q_WATER_B = 0._EB
CASE(INSULATED)
HTCB = 0._EB
QRADINB = 0._EB
E_WALLB = 0._EB
Q_WATER_B = 0._EB
TMP_BACK = TMPA
CASE(EXPOSED)
IWB = WC%BACK_INDEX
Q_WATER_B = 0._EB
IF (WALL(IWB)%BOUNDARY_TYPE==SOLID_BOUNDARY) THEN
IIB = WALL(IWB)%IIG
JJB = WALL(IWB)%JJG
KKB = WALL(IWB)%KKG
TMP_BACK = TMP(IIB,JJB,KKB)
DTMP = TMP_BACK - WC%TMP_B
HTCB = HEAT_TRANSFER_COEFFICIENT(IWB,IIB,JJB,KKB,IOR,DTMP,SF%H_FIXED,SF%GEOMETRY,SF%CONV_LENGTH)
WALL(IWB)%HEAT_TRANS_COEF = HTCB
QRADINB = WALL(IWB)%QRADIN
IF (NLP>0) Q_WATER_B = -SUM(WALL(IWB)%LP_CPUA(:))
ELSE
TMP_BACK = TMP_0(KK)
DTMP = TMP_BACK - WC%TMP_B
HTCB = HEAT_TRANSFER_COEFFICIENT(IW,-1,-1,-1,IOR,DTMP,SF%H_FIXED,SF%GEOMETRY,SF%CONV_LENGTH)