-
Notifications
You must be signed in to change notification settings - Fork 44
/
Copy pathPeleLMeX.H
2084 lines (1791 loc) · 66.7 KB
/
PeleLMeX.H
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
#ifndef PeleLMeX_H
#define PeleLMeX_H
// PeleLMeX includes
#include "PeleLMeX_UserKeys.H"
#include "PeleLMeX_Index.H"
#include "PeleLMeX_Derive.H"
#include "pelelmex_prob_parm.H"
#include "PeleLMeX_DiffusionOp.H"
#include "DiagBase.H"
#include "PeleLMeX_FlowControllerData.H"
#include "PeleLMeX_BPatch.H"
#ifdef PELE_USE_EFIELD
#include "PrecondOp.H"
#endif
// PelePhysics lib
#include <mechanism.H>
#include <PelePhysics.H>
#include <PMFData.H>
#include <turbinflow.H>
#include <ReactorBase.H>
#include <Utilities.H>
// AMReX-Hydro lib
#include <hydro_MacProjector.H>
#include <hydro_NodalProjector.H>
#include <hydro_godunov_ppm.H>
// AMReX lib
#include <AMReX_AmrCore.H>
#include <AMReX_MultiFabUtil.H>
#include <AMReX_BC_TYPES.H>
#include <AMReX_ErrorList.H>
#include <AMReX_VisMF.H>
#ifdef PELE_USE_RADIATION
#include <PeleLMRad.H>
#endif
// Forward declarations
#ifdef PELE_USE_SPRAY
class SprayParticleContainer;
#endif
#ifdef PELE_USE_SOOT
class SootModel;
#endif
class BPatch;
const std::string PrettyLine = " " + std::string(78, '=') + "\n";
class PeleLM : public amrex::AmrCore
{
public:
enum TimeStamp { AmrOldTime, AmrHalfTime, AmrNewTime };
// constructor
PeleLM();
// destructor
~PeleLM() override;
// Setup function
void Setup();
// Init function
void Init();
// Simulation function
void Evolve();
// Unit Test function
void Evaluate();
// Time advance function
void Advance(int is_init);
//-----------------------------------------------------------------------------
// Virtual AmrCore functions
void regrid(int lbase, amrex::Real time, bool initial = false) override;
void MakeNewLevelFromScratch(
int lev,
amrex::Real time,
const amrex::BoxArray& ba,
const amrex::DistributionMapping& dm) override;
void ErrorEst(
int lev, amrex::TagBoxArray& tags, amrex::Real time, int ng) override;
void MakeNewLevelFromCoarse(
int lev,
amrex::Real time,
const amrex::BoxArray& ba,
const amrex::DistributionMapping& dm) override;
void RemakeLevel(
int lev,
amrex::Real time,
const amrex::BoxArray& ba,
const amrex::DistributionMapping& dm) override;
void ClearLevel(int lev) override;
//-----------------------------------------------------------------------------
//-----------------------------------------------------------------------------
// SETUP
// ReadParameters
void readParameters();
// ReadIOParameters
void readIOParameters();
// ReadProblemParameters
void readProbParm();
// ReadGridFile
void
readGridFile(std::string grid_file, amrex::Vector<amrex::BoxArray>& input_ba);
#ifdef PELE_USE_SPRAY
// ReadSprayParameters
void SprayReadParameters();
// SpraySetup
void SpraySetup();
// SprayCreateData
void SprayCreateData();
// SprayInit
void SprayInit();
#endif
// VariablesSetup
void variablesSetup();
// DerivedSetup
void derivedSetup();
// EvaluateSetup
void evaluateSetup();
// TaggingSetup
void taggingSetup();
// ResizeArray
void resizeArray();
//-----------------------------------------------------------------------------
//-----------------------------------------------------------------------------
/**
* \brief Data container for the state data and long lived derived data
*/
struct LevelData
{
LevelData() = default;
LevelData(
amrex::BoxArray const& ba,
amrex::DistributionMapping const& dm,
amrex::FabFactory<amrex::FArrayBox> const& factory,
int a_incompressible,
int a_has_divu,
int a_nAux,
int a_nGrowState,
int a_use_soret,
int a_do_les);
// cell-centered state multifabs
amrex::MultiFab
state; // Stqte data: velocity, density, rhoYs, rhoH, Temp, RhoRT
amrex::MultiFab auxiliaries; // Auxiliary variables (passive scalars and
// others) (dim:m_nAux)
amrex::MultiFab gp; // pressure gradient (dim:AMREX_SPACEDIM)
amrex::MultiFab divu; // Velocity divergence constraint
// node-centered state multifabs
amrex::MultiFab press; // nodal pressure (dim:1)
// cell-centered transport multifabs
amrex::MultiFab visc_cc; // Viscosity (dim:1)
amrex::Array<amrex::MultiFab, AMREX_SPACEDIM>
visc_turb_fc; // Turbulent Viscosity (dim:)
amrex::Array<amrex::MultiFab, AMREX_SPACEDIM>
lambda_turb_fc; // Turbulent thermal conductivity (dim:)
amrex::MultiFab diff_cc; // Diffusivity (dim:NUM_SPECIES+2)
#ifdef PELE_USE_EFIELD
amrex::MultiFab diffE_cc; // Electron diffusivity (dim:1)
amrex::MultiFab mobE_cc; // Electron mobility (dim:1)
amrex::MultiFab mob_cc; // Species mobility (dim:NUM_IONS)
#endif
};
/**
* \brief Data container for reaction data
*/
struct LevelDataReact
{
LevelDataReact() = default;
LevelDataReact(
const amrex::BoxArray& ba,
const amrex::DistributionMapping& dm,
const amrex::FabFactory<amrex::FArrayBox>& factory);
amrex::MultiFab I_R; // Species reaction rates
amrex::MultiFab functC; // Implicit integrator function call count
#ifdef PELE_USE_EFIELD
amrex::MultiFab I_RnE; // Electron number density reaction term
#endif
};
//-----------------------------------------------------------------------------
//-----------------------------------------------------------------------------
/**
* \brief Advance function data container for diffusion term
*/
struct AdvanceDiffData
{
AdvanceDiffData() = default;
AdvanceDiffData(
int a_finestLevel,
const amrex::Vector<amrex::BoxArray>& ba,
const amrex::Vector<amrex::DistributionMapping>& dm,
const amrex::Vector<std::unique_ptr<amrex::FabFactory<amrex::FArrayBox>>>&
factory,
int a_nGrowAdv,
int a_use_wbar,
int a_use_soret,
int is_init = 0);
amrex::Vector<amrex::MultiFab> Dn; // Diffusion term t^n
amrex::Vector<amrex::MultiFab> Dnp1; // Diffusion term t^(n+1,k)
amrex::Vector<amrex::MultiFab> Dhat; // Diffusion term t^(n+1,k+1)
amrex::Vector<amrex::MultiFab> Dwbar; // Diffusion term of Wbar
amrex::Vector<amrex::Array<amrex::MultiFab, AMREX_SPACEDIM>>
wbar_fluxes; // Wbar flux correction
amrex::Vector<amrex::MultiFab> DT; // Diffusion term of T (soret flux)
amrex::Vector<amrex::Array<amrex::MultiFab, AMREX_SPACEDIM>>
soret_fluxes; // Soret flux correction
};
/**
* \brief Advance function data container for advection term
*/
struct AdvanceAdvData
{
AdvanceAdvData() = default;
AdvanceAdvData(
int a_finestLevel,
const amrex::Vector<amrex::BoxArray>& ba,
const amrex::Vector<amrex::DistributionMapping>& dm,
const amrex::Vector<std::unique_ptr<amrex::FabFactory<amrex::FArrayBox>>>&
factory,
int a_incompressible,
int a_nGrowAdv,
int a_nGrowMAC);
amrex::Vector<amrex::Array<amrex::MultiFab, AMREX_SPACEDIM>>
umac; // MAC face velocity
amrex::Vector<amrex::MultiFab> AofS; // Scalar advection term at t^(n+1/2)
amrex::Vector<amrex::MultiFab> chi; // Thermodynamic constraint
amrex::Vector<amrex::MultiFab>
Forcing; // Scalar forcing for both advection and diffusion
amrex::Vector<amrex::MultiFab> mac_divu; // divu used in MAC projection
#ifdef PELE_USE_EFIELD
amrex::Vector<amrex::Array<amrex::MultiFab, AMREX_SPACEDIM>>
uDrift; // ions drift face velocity
#endif
};
//-----------------------------------------------------------------------------
//-----------------------------------------------------------------------------
// INIT
/**
* \brief Top level initialization function
*/
void initData();
/**
* \brief Fill the initial level data container
* using the user defined init function
*/
void initLevelData(int lev);
/**
* \brief Fill the initial level data container
* using the data from a pltfile
* \param a_pltFile path to PeleLMeX plot file
*/
void initLevelDataFromPlt(int a_lev, const std::string& a_dataPltFile);
/**
* \brief Project the initial solution velocity field,
* initializing the reaction term and providing a first
* estimate of the time step size
*/
void projectInitSolution();
/**
* \brief Iterate the advance function to provide an
* initial pressure field.
*/
void initialIterations();
// initFromGridFile
void InitFromGridFile(amrex::Real time);
// regridFromFile
void regridFromGridFile(int lbase, amrex::Real time, bool initial);
/**
* \brief Check the consistency of the run
* parameters.
*/
void checkRunParams();
void checkSetupParams();
//-----------------------------------------------------------------------------
// REGRID / LOAD BALANCE
/**
* \brief Compute load balancing cost estimate on a given level
* filling the PeleLM class LayoutData with default method
* \param a_lev target level
*/
void computeCosts(int a_lev);
/**
* \brief Compute load balancing cost estimate on a given level
* \param a_lev target level
* \param a_costs LayoutData holding the cost for each box
* \param a_costMethod cost method
*/
void computeCosts(
int a_lev, amrex::LayoutData<amrex::Real>& a_costs, int a_costMethod);
/**
* \brief Create/update the DMap used for chemistry on all levels
*/
void loadBalanceChem();
/**
* \brief Create/update the DMap used for chemistry on a given level
* \param a_lev level of interest
*/
void loadBalanceChemLev(int a_lev);
//-----------------------------------------------------------------------------
//-----------------------------------------------------------------------------
// NODAL PROJECTION
/**
* \brief Initial velocity projection. The pressure and
* pressure gradient are reset to zero after the projection.
* Sigma is 1/rho instead of dt/rho.
*/
void initialProjection();
/**
* \brief Initial velocity projection considering zero
* velocity divergence (no cell-centered RHS) to establish
* the isostatic pressure when using gravity
* Sigma is 1/rho instead of dt/rho.
*/
void initialPressProjection();
/**
* \brief Advance function velocity projection.
* \param is_init flag initial iteration calls
* \param a_rhoTime time used to fill sigma
* \param a_dt time time step size
*/
void velocityProjection(
int is_init, const PeleLM::TimeStamp& a_rhoTime, const amrex::Real& a_dt);
/**
* \brief Actual nodal projection call function.
* \param a_vel vector of velocity MultiFabs
* \param a_sigma vector of variable coefficient
* \param rhs_cc vector of cell-centered projection RHS (can be empty)
* \param rhs_cc vector of node-centered projection RHS (can be empty)
* \param increment_gp flag incremental projection (where vel is U^{np1*} -
* U^{n}) \param scaling_factor used for constant coefficient projection
*/
void doNodalProject(
const amrex::Vector<amrex::MultiFab*>& a_vel,
const amrex::Vector<amrex::MultiFab*>& a_sigma,
const amrex::Vector<amrex::MultiFab*>& rhs_cc,
const amrex::Vector<const amrex::MultiFab*>& rhs_nd,
int incremental,
amrex::Real scaling_factor);
/**
* \brief For 2D-RZ, scale multifab components by radius
* including 1 ghost cell on Dirichlet BC
* \param a_lev level index
* \param a_mf MultiFab to act upon
*/
void scaleProj_RZ(int a_lev, amrex::MultiFab& a_mf);
/**
* \brief For 2D-RZ, unscale multifab components by radius
* \param a_lev level index
* \param a_mf MultiFab to act upon
*/
void unscaleProj_RZ(int a_lev, amrex::MultiFab& a_mf);
//-----------------------------------------------------------------------------
//-----------------------------------------------------------------------------
// MAC PROJECTION
/**
* \brief Re-create the MAC projector (reset the std::unique_ptr)
* if the circonstances require it
*/
void resetMacProjector();
/**
* \brief Use Godunov to obtain face-centered velocities
* from cell-centered ones
* \param advData container for the face-centered velocity umac
*/
void predictVelocity(std::unique_ptr<AdvanceAdvData>& advData);
/**
* \brief Create MAC projection RHS by fillpatching divU
* between old and new time
* \param advData container for mac_divu
*/
void createMACRHS(std::unique_ptr<AdvanceAdvData>& advData);
/**
* \brief Create MAC projection RHS by fillpatching divU
* between old and new time
* \param a_sdcIter current SDC iteration index
* \param a_time time stamp to access thermodynamic pressure
* \param advData container for chi and incremented mac_divu
*/
void addChiIncrement(
int a_sdcIter,
const PeleLM::TimeStamp& a_time,
std::unique_ptr<AdvanceAdvData>& advData);
/**
* \brief Perform MAC projection to obtained divergence constrained
* face-centered (MAC) velocities
* \param a_time time stamp to access the density used in the projection
* \param advData container for the face-centered velocities umac
* \param a_divu velocity divergence constraint (can be empty)
*/
void macProject(
const PeleLM::TimeStamp& a_time,
std::unique_ptr<AdvanceAdvData>& advData,
const amrex::Vector<amrex::MultiFab*>& a_divu);
/**
* \brief Fill the projected face velocities ghost faces, ensuring
* that the ghost cell velocities are constrained by the ghost cell
* divergence constraint
* \param a_lev target level
* \param a_nGrow number of ghost faces
* \param crse_geom pointer to the geometry of the next coarser level
* \param fine_geom pointer to the geometry of the current level
* \param u_mac_crse face-centered umac of the next coarser level
* \param u_mac_fine face-centered umac of the current level (those updated by
* the function) \param divu pointer to the current level divU \param
* crse_ratio refinement ratio with the next coarser level
*/
void create_constrained_umac_grown(
int a_nGrow,
const amrex::Geometry* crse_geom,
const amrex::Geometry* fine_geom,
amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM> u_mac_crse,
amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM> u_mac_fine,
const amrex::IntVect& crse_ratio);
/**
* \brief Setup the linear operator boundary condition for the
* MAC projection in a given direction
* \param a_side an AMReX orientation object
*/
amrex::Array<amrex::LinOpBCType, AMREX_SPACEDIM>
getMACProjectionBC(amrex::Orientation::Side a_side);
//-----------------------------------------------------------------------------
//-----------------------------------------------------------------------------
// DIFFUSION
/**
* \brief Compute the cell-centered viscosity on all levels
* filling the specified container of Old or New LevelData
* \param a_time either Old or New stamp
*/
void calcViscosity(const PeleLM::TimeStamp& a_time);
/**
* \brief Compute the face-centered turbulent viscosity on all levels
* filling the specified container of Old or New LevelData
* \param a_time either Old or New stamp
*/
void calcTurbViscosity(const PeleLM::TimeStamp& a_time);
/**
* \brief Compute the cell-centereds diffusivity on all levels
* filling the specified container of Old or New LevelData
* \param a_time either Old or New stamp
*/
void calcDiffusivity(const PeleLM::TimeStamp& a_time);
// get edge-centered diffusivity on a per level / per comp basis
/**
* \brief Compute face-averaged diffusivity (or else) from
* cell-centered data, accounting for Dirichlet BCs, allowing
* to call user-defined zero_visc() function and possibly adding
* a turbulent contribution
* \param lev target level
* \param beta_comp first component to use in a_diff_cc MultiFab
* \param ncomp number of components of the outgoing data
* \param doZeroVisc zero_visc() call trigger
* \param a_bcrec boundary conditions for the required components
* \param a_diff_cc cell-centered data
* \param addTurbContrib turbulent contribution increment trigger
*/
amrex::Array<amrex::MultiFab, AMREX_SPACEDIM> getDiffusivity(
int lev,
int beta_comp,
int ncomp,
int doZeroVisc,
amrex::Vector<amrex::BCRec> bcrec,
const amrex::MultiFab& beta_cc,
int addTurbContrib = 0);
/**
* \brief Compute the explicit face-centered diffusion fluxes on all levels
* for species and enthalpy (NUM_SPECIES+2 components) using Old or
* New state/coefficient data from the LevelData container
* \param a_time either Old or New stamp
* \param a_fluxes outgoing fluxes container, the first NUM_SPECIES components
* filled \param a_EBfluxes optional EB fluxes container \param a_wbarfluxes
* optional Wbar fluxes container \param a_soretfluxes optional Soret fluxes
* container
*/
void computeDifferentialDiffusionFluxes(
const PeleLM::TimeStamp& a_time,
const amrex::Vector<amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>>&
a_fluxes,
const amrex::Vector<amrex::MultiFab*>& a_EBfluxes,
const amrex::Vector<amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>>&
a_wbarfluxes,
const amrex::Vector<amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>>&
a_soretfluxes);
void correctIsothermalBoundary(
const TimeStamp& a_time,
const amrex::Vector<amrex::MultiFab*>& a_spec_boundary,
const amrex::Vector<amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>>&
a_wbarfluxes,
const amrex::Vector<amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>>&
a_soretfluxes);
/**
* \brief Compute the explicit cell-centered diffusion term on all levels
* by computing the divergence of the diffusion fluxes at Old or New time,
* filling the proper container in diffData
* \param a_time either Old or New stamp
* \param diffData container for the returning diffusion terms
* \param is_init initial computation flag
*/
void computeDifferentialDiffusionTerms(
const PeleLM::TimeStamp& a_time,
std::unique_ptr<AdvanceDiffData>& diffData,
int is_init = 0);
/**
* \brief Add the Wbar contribution to the face-centered species diffusion
* fluxes on all levels optionally returning that contribution in a separate
* container \param a_fluxes diffusion fluxes to be updated \param
* a_wbarfluxes optional container to return the Wbar contribution \param
* a_spec species rhoYs state data on all levels \param a_rho density state
* data on all levels \param a_beta cell-centered species diffusivity on all
* levels
*/
void addWbarTerm(
const amrex::Vector<amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>>&
a_spfluxes,
const amrex::Vector<amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>>&
a_spwbarfluxes,
amrex::Vector<amrex::MultiFab const*> const& a_spec,
amrex::Vector<amrex::MultiFab const*> const& a_rho,
amrex::Vector<amrex::MultiFab const*> const& a_beta,
amrex::Vector<amrex::MultiFab const*> const& a_boundary);
/**
* \brief Add the Soret contribution to the face-centered species diffusion
* fluxes on all levels optionally returning that contribution in a separate
* container \param a_fluxes diffusion fluxes to be updated \param
* a_soretfluxes optional container to return the Soret contribution \param
* a_spec species rhoYs state data on all levels \param a_rho density state
* data on all levels \param a_temp temperature state data on all levels
* \param a_beta cell-centered species diffusivity on all levels
*/
void addSoretTerm(
const amrex::Vector<amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>>&
a_spfluxes,
const amrex::Vector<amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>>&
a_spsoretfluxes,
amrex::Vector<amrex::MultiFab const*> const& a_temp,
amrex::Vector<amrex::MultiFab const*> const& a_beta);
/**
* \brief Modify the species diffusion fluxes by computing a correction
* velocity, ensuring they sum up to zero. \param a_fluxes diffusion fluxes to
* be updated \param a_spec species rhoYs state data on all levels
*/
template <typename EOSType>
void adjustSpeciesFluxes(
const amrex::Vector<amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>>&
a_spfluxes,
amrex::Vector<amrex::MultiFab const*> const& a_spec);
/**
* \brief Compute the enthalpy flux due to species differential diffusion,
* filling the NUM_SPECIES+1 component of a_fluxes
* \param a_fluxes outgoing fluxes container, filled on the NUM_SPECIES+1
* component \param a_temp temperature state data on all levels
*/
void computeSpeciesEnthalpyFlux(
const amrex::Vector<amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>>&
a_fluxes,
amrex::Vector<amrex::MultiFab const*> const& a_temp);
/**
* \brief Implicit diffusion solves
* \param advData container for the external forcing of the linear solve
* \param diffData container for the returning implicit diffusion term
*/
void differentialDiffusionUpdate(
std::unique_ptr<AdvanceAdvData>& advData,
std::unique_ptr<AdvanceDiffData>& diffData);
/**
* \brief Setup the implicit deltaT linear solve
* \param a_rhs outgoing RHS for the linear solve
* \param a_Tsave outgoing temporary container for temperature while solving
* for deltaT \param a_rhoCp outgoing linear solve coefficient \param advData
* container for external forcing used to build RHS \param diffData container
* for diffusion forcing used to build RHS
*/
void deltaTIter_prepare(
const amrex::Vector<amrex::MultiFab*>& a_rhs,
const amrex::Vector<amrex::MultiFab*>& a_Tsave,
const amrex::Vector<amrex::MultiFab*>& a_rhoCp,
std::unique_ptr<AdvanceAdvData>& advData,
std::unique_ptr<AdvanceDiffData>& diffData);
/**
* \brief Restore state, update enthalpy and compute fluxes after deltaT
* linear solve within the deltaT iteration scheme \param a_dtiter deltaT
* iteration index \param a_fluxes diffusion fluxes, with updated Fourier and
* differential species fluxes \param a_ebfluxes optional temperature fluxes
* at the EB \param a_Tsave temporary container to restore state temperature
* \param diffData container for diffusion forcing used to build RHS
* \param a_deltaT_norm norm of deltaT residual
*/
void deltaTIter_update(
int a_dtiter,
const amrex::Vector<amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>>&
a_fluxes,
const amrex::Vector<amrex::MultiFab*>& a_ebfluxes,
const amrex::Vector<amrex::MultiFab const*>& a_Tsave,
std::unique_ptr<AdvanceDiffData>& diffData,
amrex::Real& a_deltaT_norm);
/**
* \brief Compute the viscous force term in the momentum equation using AMReX
* TensorOp
* \param a_time state time for velocity, density and viscosity, either Old or
* New \param a_divtau outgoing divTau term \param use_density flag to scale
* divTau by 1/rho \param scale prefactor of divTau, used to weight by num.
* scheme property (e.g. 0.5 for Crank Nicholson)
*/
void computeDivTau(
const PeleLM::TimeStamp& a_time,
const amrex::Vector<amrex::MultiFab*>& a_divtau,
int use_density,
amrex::Real scale = 1.0);
/**
* \brief Implicit diffusion of the of velocity field using AMReX
* TensorOp. Always act on New time data and include Crank Nicholson factor
*/
void diffuseVelocity();
/**
* \brief Return LinOpBC for the scalar diffusion operator
* \param a_side AMReX orientation
* \param a_bc BCRec for the scalar of interest
*/
amrex::Array<amrex::LinOpBCType, AMREX_SPACEDIM> getDiffusionLinOpBC(
amrex::Orientation::Side a_side, const amrex::BCRec& a_bc);
/**
* \brief Return AMReX_SPACEDIM LinOpBC vector for the velocity tensor
* operator \param a_side AMReX orientation \param a_bc BCRec for all the
* velocity components
*/
amrex::Vector<amrex::Array<amrex::LinOpBCType, AMREX_SPACEDIM>>
getDiffusionTensorOpBC(
amrex::Orientation::Side a_side, const amrex::Vector<amrex::BCRec> a_bc);
//-----------------------------------------------------------------------------
//-----------------------------------------------------------------------------
// REACTION
/**
* \brief Top-level reaction solve function, performing the chemistry
* integration on all the AMR levels. \param advData container for Forcing,
* the advection/diffusion component of the advance function in the SDC
* context
*/
void advanceChemistry(std::unique_ptr<AdvanceAdvData>& advData);
/**
* \brief Performing the chemistry integration on a given level, for a fixed
* step length and relying on AmrCore BoxArray/DMap, not masking fine-covered
* regions \param lev level of interest \param a_dt integration length \param
* a_extForcing advection/diffusion forcing
*/
void advanceChemistry(
int lev, const amrex::Real& a_dt, amrex::MultiFab& a_extForcing);
/**
* \brief Performing the chemistry integration on a given level, for a fixed
* step length and using the BoxArray/DMap masking fine-covered region and
* specifically load balanced for reaction
* \param lev level of interest
* \param a_dt integration length
* \param a_extForcing advection/diffusion forcing
*/
void advanceChemistryBAChem(
int lev, const amrex::Real& a_dt, amrex::MultiFab& a_extForcing);
/**
* \brief Top-level instantaneous reaction rate function, acting on all levels
* \param a_I_R outgoing multi-level container inst. RR container
* \param a_time state time used to compute RR, either Old or New
*/
void computeInstantaneousReactionRate(
const amrex::Vector<amrex::MultiFab*>& a_I_R,
const PeleLM::TimeStamp& a_time);
/**
* \brief Compute instantaneous reaction rate on a given level
* \param lev level of interest
* \param a_time state time used to compute RR, either Old or New
* \param a_I_R outgoing MultiFab for inst. RR container
*/
void computeInstantaneousReactionRate(
int lev, const PeleLM::TimeStamp& a_time, amrex::MultiFab* a_I_R);
/**
* \brief Compute the heat release rate on a given level, using the
* LevelDataReact I_R data
* \param lev level of interest
* \param a_HR outgoing MultiFab
*/
void getHeatRelease(int lev, amrex::MultiFab* a_HR);
//-----------------------------------------------------------------------------
//-----------------------------------------------------------------------------
// Time step size
/**
* \brief Compute the time step across all levels and processes.
* \param is_init flag to trigger init-only dt factors
* \param a_time state time to work on, either Old or New
*/
amrex::Real computeDt(int is_init, const PeleLM::TimeStamp& a_time);
/**
* \brief Compute the CFL step size
* \param a_time state time to work on, either Old or New
*/
amrex::Real estConvectiveDt(const PeleLM::TimeStamp& a_time);
/**
* \brief Compute the divU step size
* \param a_time state time to work on, either Old or New
*/
amrex::Real estDivUDt(const PeleLM::TimeStamp& a_time);
void checkDt(const PeleLM::TimeStamp& a_time, const amrex::Real& a_dt);
//-----------------------------------------------------------------------------
//-----------------------------------------------------------------------------
// FORCES
void getVelForces(
const TimeStamp& a_time,
const amrex::Vector<amrex::MultiFab*>& a_divTau,
const amrex::Vector<amrex::MultiFab*>& a_velForce,
int nGrowForce,
int add_gradP);
void getVelForces(
const TimeStamp& a_time,
int lev,
amrex::MultiFab* a_divTau,
amrex::MultiFab* a_velForce,
int add_gradP);
void getVelForces(
int lev,
const amrex::Box& bx,
const amrex::Real& a_time,
amrex::Array4<amrex::Real> const& force,
amrex::Array4<const amrex::Real> const& vel,
amrex::Array4<const amrex::Real> const& rho,
amrex::Array4<const amrex::Real> const& rhoY,
amrex::Array4<const amrex::Real> const& rhoh,
amrex::Array4<const amrex::Real> const& temp,
amrex::Array4<const amrex::Real> const& extMom,
amrex::Array4<const amrex::Real> const& extRho);
void addSpark(const PeleLM::TimeStamp& a_time);
//-----------------------------------------------------------------------------
#ifdef PELE_USE_SPRAY
//-----------------------------------------------------------------------------
// Sprays
void setupVirtualParticles(const int level);
void removeVirtualParticles(const int level);
void setupGhostParticles(const int ngrow, const int level);
void removeGhostParticles(const int level);
amrex::Real SprayEstDt();
void SprayMKD(const amrex::Real time, const amrex::Real dt);
void
SprayMKDLevel(const int level, const amrex::Real time, const amrex::Real dt);
void SprayInjectRedist();
void SprayPostRegrid();
void SpraySetState(const amrex::Real& a_flow_dt);
void SprayAddSource(const int level);
static int num_spray_src;
static bool do_spray_particles;
static std::unique_ptr<SprayParticleContainer> SprayPC;
static std::unique_ptr<SprayParticleContainer> VirtPC;
static std::unique_ptr<SprayParticleContainer> GhostPC;
// State used for spray interpolation
amrex::Vector<std::unique_ptr<amrex::MultiFab>> m_spraysource;
amrex::Vector<std::unique_ptr<amrex::MultiFab>> m_spraystate;
//-----------------------------------------------------------------------------
#endif
#ifdef PELE_USE_SOOT
//-----------------------------------------------------------------------------
// Soot
// SetSootIndx
void setSootIndx();
// CleanupSootModel
void cleanupSootModel();
// AddSootDerivePlotVars
// void addSootDerivePlotVars(PeleLMDeriveList& derive_lst);
// Compute the soot source terms
void computeSootSource(
const PeleLM::TimeStamp& a_timestamp, const amrex::Real a_dt);
// Clip soot moments
void clipSootMoments();
//-----------------------------------------------------------------------------
#endif
#ifdef PELE_USE_RADIATION
// Initialize the radiation module
void RadInit();
// Compute the radiation source term
void computeRadSource(const PeleLM::TimeStamp& a_timestamp);
#endif
//-----------------------------------------------------------------------------
// External Sources
bool m_user_defined_ext_sources;
bool m_ext_sources_SDC;
void getExternalSources(
int is_initIter,
const PeleLM::TimeStamp& a_timestamp_old,
const PeleLM::TimeStamp& a_timestamp_new);
//-----------------------------------------------------------------------------
// EOS
/**
* \brief Compute thermodynamic pressure (RhoRT) from state on all levels
* \param a_time state time to work on, either Old or New
*/
void setThermoPress(const PeleLM::TimeStamp& a_time);
/**
* \brief Compute thermodynamic pressure (RhoRT) from state on a single level
* \param lev level of interest
* \param a_time state time to work on, either Old or New
*/
void setThermoPress(int lev, const PeleLM::TimeStamp& a_time);
/**
* \brief Enforce rho = \sum rhoYs on a given level
* \param lev level of interest
* \param a_time state time to work on, either Old or New
*/
void setRhoToSumRhoY(int lev, const PeleLM::TimeStamp& a_time);
/**
* \brief Compute temperature from rhoY and rhoH on all levels
* \param a_time state time to work on, either Old or New
*/
void setTemperature(const PeleLM::TimeStamp& a_time);
/**
* \brief Compute temperature from rhoY and rhoH on a given level
* \param lev level of interest
* \param a_time state time to work on, either Old or New
*/
void setTemperature(int lev, const PeleLM::TimeStamp& a_time);
/**
* \brief Compute the divergence constraint from reaction/diffusion
* terms
* \param is_init flag for init-only checks
* \param computeDiff flag to trigger computing the diffusion terms or
* using the one already precompute and available in diffData
* \param a_time state time to work with
* \param diffData container for the diffusion terms
*/
void calcDivU(
int is_init,
int computeDiff,
int do_avgDown,
const PeleLM::TimeStamp& a_time,
std::unique_ptr<AdvanceDiffData>& diffData);
/**
* \brief Compute the pressure drift term using RhoRT and
* state data at a_time, acting on all levels
* \param a_time state time to work on, either Old or New
* \param a_dPdt outgoing container
*/
void calc_dPdt(
const PeleLM::TimeStamp& a_time,
const amrex::Vector<amrex::MultiFab*>& a_dPdt);
/**
* \brief Compute the pressure drift term using RhoRT and
* state data at a_time on a given level
* \param lev level of interest
* \param a_time state time to work on, either Old or New
* \param a_dPdt outgoing container
*/
void
calc_dPdt(int lev, const PeleLM::TimeStamp& a_time, amrex::MultiFab* a_dPdt);
/**
* \brief For the closed chamber algorithm, compute averaged
* mac_divU and theta, update divU and the spatially averaged
* thermodynamic pressure. Return the mean divU.
* \param advData container including mac_divu
*/
amrex::Real adjustPandDivU(std::unique_ptr<AdvanceAdvData>& advData);
//-----------------------------------------------------------------------------
//-----------------------------------------------------------------------------
// I/O
void WritePlotFile();
bool writePlotNow() const;
bool checkMessage(const std::string& a_action) const;
void WriteCheckPointFile();
void ReadCheckPointFile();
bool writeCheckNow() const;
void WriteJobInfo(const std::string& path) const;
void WriteHeader(const std::string& name, bool is_checkpoint) const;
void WriteDebugPlotFile(
const amrex::Vector<const amrex::MultiFab*>& a_MF,
const std::string& pltname);
//-----------------------------------------------------------------------------
//-----------------------------------------------------------------------------
// Diagnostics
void createDiagnostics();
void updateDiagnostics();
void doDiagnostics();
//-----------------------------------------------------------------------------