-
Notifications
You must be signed in to change notification settings - Fork 0
/
PeleLM.H
922 lines (778 loc) · 33.7 KB
/
PeleLM.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
#ifndef _PeleLM_H_
#define _PeleLM_H_
#include <AMReX_AuxBoundaryData.H>
#include <NavierStokesBase.H>
#include <AMReX_ForkJoin.H>
#include <mechanism.H>
#include <pelelm_prob_parm.H>
#include <PeleLM_parm.H>
#include "PelePhysics.H"
#include "PMFData.H"
#include "ReactorBase.H"
#include <list>
#include <map>
#include <utility>
#ifdef AMREX_PARTICLES
#include <AMReX_Particles.H>
#endif
//
// Note: define TEMPERATURE if you want a variable T in the
// State_Type part of the state
// whose component index is Temp and
// whose evoution equation is
// \pd (rho T)/\pd t + diver (\rho U T) = \diver k/c_p grad T
// define RADIATION only if TEMPERATURE is also defined and the evolution equation
// for T is
// \pd (rho T)/\pd t + diver (\rho U T) = \diver k/c_p grad T - 1/c_p diver q_rad
//
// Note that component Temp is T, not rho T. This was done so
// that we could use the existing diffusion operator and
// multigrid code.
//
extern "C"
{
typedef void (*LMEF)(BL_FORT_IFAB_ARG_ANYD(tag),
const int* tagval, const int* clearval,
const BL_FORT_FAB_ARG_ANYD(data),
const int* lo, const int* hi, const int* nvar,
const int* domain_lo, const int* domain_hi,
const amrex::Real* dx, const amrex::Real* xlo,
const amrex::Real* prob_lo, const amrex::Real* time,
const int* level, const amrex::Real* value);
typedef void (*LMEF_BOX) (BL_FORT_IFAB_ARG_ANYD(tag),
const int* tagval, const int* clearval,
const amrex::Real* blo, const amrex::Real* bhi,
const int* lo, const int* hi,
const int* domain_lo, const int* domain_hi,
const amrex::Real* dx, const amrex::Real* xlo,
const amrex::Real* prob_lo, const amrex::Real* time,
const int* level);
}
class LM_Error_Value
:
public amrex::ErrorRec::ErrorFunc
{
public:
LM_Error_Value()
:
lmef(0), lmef_box(0), value(), min_time(), max_time(), max_level() {}
LM_Error_Value (amrex::Real min_time, amrex::Real max_time, int max_level);
LM_Error_Value (LMEF lmef, amrex::Real value, amrex::Real min_time,
amrex::Real max_time, int max_level);
LM_Error_Value (LMEF_BOX lmef_box, const amrex::RealBox& box, amrex::Real min_time,
amrex::Real max_time, int max_level);
virtual ~LM_Error_Value () {}
virtual amrex::ErrorRec::ErrorFunc* clone () const override {
if (BoxTag()) {
return new LM_Error_Value(lmef_box,box,min_time,max_time,max_level);
}
return new LM_Error_Value(lmef,value,min_time,max_time,max_level);
}
void tagCells(int* tag, const int* tlo, const int* thi,
const int* tagval, const int* clearval,
const amrex::Real* data, const int* dlo, const int* dhi,
const int* lo, const int* hi, const int* nvar,
const int* domain_lo, const int* domain_hi,
const amrex::Real* dx, const amrex::Real* xlo,
const amrex::Real* prob_lo, const amrex::Real* time,
const int* level) const;
void tagCells1(int* tag, const int* tlo, const int* thi,
const int* tagval, const int* clearval,
const int* lo, const int* hi,
const int* domain_lo, const int* domain_hi,
const amrex::Real* dx, const amrex::Real* xlo,
const amrex::Real* prob_lo, const amrex::Real* time,
const int* level) const;
int MaxLevel() const {return max_level;}
amrex::Real MinTime() const {return min_time;}
amrex::Real MaxTime() const {return max_time;}
amrex::Real Value() const {return value;}
bool BoxTag() const {return lmef==0 && lmef_box!=0;}
protected:
LMEF lmef;
LMEF_BOX lmef_box;
amrex::Real value, min_time, max_time;
amrex::RealBox box;
int max_level;
};
// Forward declarations
#ifdef AMREX_PARTICLES
class SprayParticleContainer;
#endif
#ifdef SOOT_MODEL
class SootModel;
#endif
class PeleLM
:
public NavierStokesBase
{
public:
PeleLM ();
PeleLM (amrex::Amr& papa,
int lev,
const amrex::Geometry& level_geom,
const amrex::BoxArray& bl,
const amrex::DistributionMapping& dm,
amrex::Real time);
virtual ~PeleLM ();
////////////////////////////////////////////////////////////////////////////
// AmrLevel virtual functions //
////////////////////////////////////////////////////////////////////////////
//
// Advance grids at this level in time.
//
virtual amrex::Real advance (amrex::Real time,
amrex::Real dt,
int iteration,
int ncycle) override;
virtual void checkPoint (const std::string& dir,
std::ostream& os,
amrex::VisMF::How how = amrex::VisMF::OneFilePerCPU,
bool dump_old = true) override;
//
// Returns a MultiFab containing the derived data for this level.
// The user is responsible for deleting this pointer when done
// with it. If ngrow>0 the MultiFab is built on the appropriately
// grown amrex::BoxArray.
//
virtual std::unique_ptr<amrex::MultiFab> derive (const std::string& name,
amrex::Real time,
int ngrow) override;
//
// This version of derive() fills the dcomp'th component of mf with the derived quantity.
//
virtual void derive (const std::string& name,
amrex::Real time,
amrex::MultiFab& mf,
int dcomp) override;
//
// Init data on this level after regridding if old level
// did not exist previously.
//
virtual void init () override;
//
// Init data on this level from another NavierStokes (during regrid).
//
virtual void init (amrex::AmrLevel& old) override;
virtual void initData () override;
#ifdef AMREX_PARTICLES
//
// Read particle-related inputs
//
static void readParticleParams();
//
// Define particle related variables
//
static void defineParticles();
//
// Return names of spray source terms
//
static std::string spraySrcName(const int i);
//
// Define gas phase state MF for computing spray source terms
//
void defineSprayStateMF();
//
// Initialize particle locations and velocities (and strengths if relevant)
//
void initParticles();
//
// Create relevant particle data
//
void createParticleData();
//
// Timestamp particles
//
void particleTimestamp(int ngrow);
//
// Default verbosity of Particle class
//
static int particle_verbose;
//
// How to initialize at restart
//
void particlePostRestart(bool is_checkpoint = true);
//
// Redistribute
//
virtual void particle_redistribute(int lbase = 0, bool init_part = false) override;
//
// Setup ghost and virtual particles and do moveKickDrift
//
void particleMKD(const amrex::Real time,
const amrex::Real dt,
const int ghost_width,
const int spray_n_grow,
const int tmp_src_width,
const int where_width,
amrex::MultiFab& tmp_spray_source);
//
// Do the moveKick for active and ghost particles
//
void particleMK(const amrex::Real time,
const amrex::Real dt,
const int spray_n_grow,
const int tmp_src_width,
amrex::MultiFab& tmp_spray_source);
//
// Setup virtual particles if necessary
//
void setupVirtualParticles();
//
// Remove virtual particles if necessary
//
void removeVirtualParticles();
//
// Setup ghost particles (for finer levels) if necessary
//
void setupGhostParticles(int ngrow);
//
// Remove ghost particles (for this level) if necessary
//
void removeGhostParticles();
//
// Estimate timestep for particles
//
void particleEstTimeStep(amrex::Real& est_dt);
//
// Default cfl of particles in Particle class
//
static amrex::Real particle_cfl;
//
// Should we write particle ascii files?
//
static int write_spray_ascii_files;
//
// Should we write the gas phase spray source terms?
//
static int plot_spray_src;
/// A state array with ghost zones
amrex::MultiFab Sborder;
static int particle_mass_tran;
static int particle_mom_tran;
static int num_spray_src;
static amrex::Vector<std::string> sprayFuelNames;
#endif
#ifdef USE_ASCENT
virtual void goAscent (int nstep);
#endif
//
// Contains operations to be done only after a full coarse timestep.
//
virtual void postCoarseTimeStep (amrex::Real cumtime) override;
virtual void post_init (amrex::Real stop_time) override;
virtual void post_regrid (int lbase, int new_finest) override;
virtual void post_restart () override;
virtual void post_timestep (int iteration) override;
virtual void restart (amrex::Amr& papa,
std::istream& is,
bool bReadSpecial = false) override;
//
// Set the variables that are put in the plotfile....
//
virtual void setPlotVariables () override;
//
// Set time levels of state data.
//
virtual void setTimeLevel (amrex::Real time,
amrex::Real dt_old,
amrex::Real dt_new) override;
//
// Write plot file stuff to specified directory.
//
virtual void writePlotFile (const std::string& dir,
std::ostream& os,
amrex::VisMF::How how) override;
////////////////////////////////////////////////////////////////////////////
// PeleLM public static functions //
////////////////////////////////////////////////////////////////////////////
//
// Define data descriptors.
//
static void variableSetUp ();
static void rhoydotSetUp ();
#ifdef AMREX_PARTICLES
static void spraydotSetUp ();
#endif
#ifdef SOOT_MODEL
static void setSootIndx ();
static void sootsrcSetUp ();
void computeSootSrc (amrex::Real time, amrex::Real dt);
void clipSootMoments();
void estSootTimeStep (amrex::Real& estdt);
#endif
//
// Cleanup data descriptors at end of run.
//
static void variableCleanUp ();
static void getSpeciesNames(amrex::Vector<std::string>& spn);
static int getSpeciesIdx(const std::string& spName);
void advance_chemistry (amrex::MultiFab& mf_old,
amrex::MultiFab& mf_new,
amrex::Real dt,
const amrex::MultiFab& Force);
void compute_scalar_advection_fluxes_and_divergence (const amrex::MultiFab& Force,
const amrex::MultiFab& divu,
amrex::Real dt);
// Mixture fraction machinery
static amrex::Array<amrex::Real, 4> Beta_mix;
static amrex::Array<amrex::Real, NUM_SPECIES> spec_Bilger_fact;
static amrex::Real Zfu;
static amrex::Real Zox;
static bool mixture_fraction_ready;
////////////////////////////////////////////////////////////////////////////
// Overriding Virtual Functions in NavierStokesBase //
////////////////////////////////////////////////////////////////////////////
//
// Setup for a level timestep.
//
virtual void advance_setup (amrex::Real time,
amrex::Real dt,
int iteration,
int ncycle) override;
virtual void avgDown () override; // Average down for all the state types.
//
// Note: these two functions must be supplied.
// If divu is not included in the state, then
// they can be no-op functions
//
virtual void calc_divu (amrex::Real time,
amrex::Real dt,
amrex::MultiFab& divu) override;
virtual void calcViscosity (const amrex::Real time,
const amrex::Real dt,
const int iteration,
const int ncycle) override;
virtual void calcDiffusivity (const amrex::Real time) override;
virtual amrex::Real estTimeStep () override;
virtual void getViscosity (amrex::MultiFab* viscosity[BL_SPACEDIM],
const amrex::Real time) override;
virtual void getViscTerms (amrex::MultiFab& visc_terms,
int src_comp,
int num_comp,
amrex::Real time) override;
virtual void getForce (amrex::FArrayBox& force,
const amrex::Box& bx,
int strt_comp,
int num_comp,
const amrex::Real time,
const amrex::FArrayBox& Vel,
const amrex::FArrayBox& Scal,
int scalScomp,
const amrex::MFIter& mfi) override;
virtual void mac_sync () override;
//
// Crse/fine fixup functions.
//
virtual void reflux () override;
//
// Reset time levels for the initial iterations.
//
virtual void resetState (amrex::Real time,
amrex::Real dt_old,
amrex::Real dt_new) override;
virtual void sum_integrated_quantities () override;
virtual void scalar_advection_update (amrex::Real dt,
int first_scalar,
int last_scalar) override;
virtual void velocity_diffusion_update (amrex::Real dt) override;
virtual void errorEst (amrex::TagBoxArray& tags,
int clearval,
int tagval,
amrex::Real time,
int n_error_buf,
int ngrow) override;
////////////////////////////////////////////////////////////////////////////
// PeleLM protected static functions //
////////////////////////////////////////////////////////////////////////////
static void Initialize ();
static void Initialize_specific ();
static void init_mixture_fraction ();
static void Finalize ();
#ifdef AMREX_PARTICLES
static int do_spray_particles;
static SprayParticleContainer* theSprayPC ();
static SprayParticleContainer* theVirtPC ();
static SprayParticleContainer* theGhostPC ();
#endif
#ifdef SOOT_MODEL
static int do_soot_solve;
static int plot_soot_src;
static SootModel* soot_model;
static int sootsrc_Type;
static int first_soot;
static int NUM_SOOT_VARS;
static int num_soot_src;
#endif
//
// MultiFab for holding source terms like from soot and sprays
//
amrex::MultiFab external_sources;
void add_external_sources(amrex::Real time, amrex::Real dt);
// enum YdotAction { HT_EstimateYdotNew, HT_ImproveYdotOld, HT_LeaveYdotAlone };
// enum Solver_Status {HT_InProgress, HT_Stalled, HT_Solved};
void adjust_spec_diffusion_fluxes (amrex::MultiFab* const * flux,
const amrex::MultiFab& S,
const amrex::BCRec& bc);
void calcDiffusivity_Wbar (const amrex::Real time);
void calc_dpdt (amrex::Real time,
amrex::Real dt,
amrex::MultiFab& dpdt);
void checkTimeStep (amrex::Real a_time,
amrex::Real a_dt);
void compute_differential_diffusion_fluxes (const amrex::MultiFab& S,
const amrex::MultiFab* Scrse,
amrex::MultiFab* const * flux,
const amrex::MultiFab* const * beta,
amrex::Real dt,
amrex::Real time,
bool include_Wbar_fluxes);
void compute_differential_diffusion_terms (amrex::MultiFab& D,
amrex::MultiFab& DD,
amrex::Real time,
amrex::Real dt,
bool include_Wbar_terms);
void compute_enthalpy_fluxes (amrex::MultiFab* const* flux,
const amrex::MultiFab* const* beta,
amrex::Real time);
enum HowToFillGrow {HT_ZERO_GROW_CELLS, HT_EXTRAP_GROW_CELLS, HT_NUM_GROW_OPTIONS};
void compute_instantaneous_reaction_rates (amrex::MultiFab& R,
const amrex::MultiFab& S,
amrex::Real time,
int nGrow = 0,
HowToFillGrow how = HT_ZERO_GROW_CELLS);
amrex::Real adjust_p_and_divu_for_closed_chamber(amrex::MultiFab& mac_divu);
void compute_rhohmix (amrex::Real time,
amrex::MultiFab& rhohmix,
int dComp);
void compute_rhoRT (const amrex::MultiFab& S,
amrex::MultiFab& P,
int pComp);
void compute_Wbar_fluxes(const amrex::MultiFab &a_scalars,
amrex::Real time,
int inc_flag,
amrex::Real inc_coeff);
void define_data ();
void differential_diffusion_update (amrex::MultiFab& Force,
amrex::MultiFab& DWbar,
int FComp,
amrex::MultiFab& D,
int DComp,
amrex::MultiFab& DD);
void differential_spec_diffuse_sync (amrex::Real dt,
bool Wbar_corrector,
bool last_mac_sync_iter);
void diffuse_velocity_setup (amrex::Real dt,
amrex::MultiFab*& delta_rhs,
FluxBoxes& fb_betan,
FluxBoxes& fb_betanp1);
void flux_divergenceRD (const amrex::MultiFab &a_state,
int stateComp,
amrex::MultiFab &a_divergence,
int divComp,
const amrex::MultiFab* const* extensive_fluxes,
int fluxComp,
int nComp,
amrex::BCRec const* d_bc,
amrex::Real scale,
amrex::Real a_dt,
int areSpeciesFluxes = 0);
void flux_divergence (amrex::MultiFab& a_divergence,
int divComp,
const amrex::MultiFab* const* a_ext_fluxes,
int fluxComp,
int nComp,
amrex::Real scale);
void getDiffusivity (amrex::MultiFab* diffusivity[AMREX_SPACEDIM],
const amrex::Real time,
const int state_comp,
const int dst_comp,
const int num_comp);
void getDiffusivity_Wbar (amrex::MultiFab* diffusivity[AMREX_SPACEDIM],
const amrex::Real time);
amrex::DistributionMapping getFuncCountDM (const amrex::BoxArray& bxba, int ngrow);
PeleLM& getLevel (int lev)
{
return *(PeleLM*) &parent->getLevel(lev);
}
void initDataOtherTypes ();
void post_init_press (amrex::Real& dt_init,
amrex::Vector<int>& nc_save,
amrex::Vector<amrex::Real>& dt_save);
void set_htt_hmixTYP ();
void set_reasonable_grow_cells_for_R(amrex::Real time);
void setThermoPress(amrex::Real time);
void set_typical_values(bool restart);
void reset_typical_values(const amrex::MultiFab& S);
void update_typical_values_chem();
void activeControl(const int step,
const int restart,
const amrex::Real time,
const amrex::Real dt);
void initActiveControl();
static void parseComposition(amrex::Vector<std::string> compositionIn,
std::string compositionType,
amrex::Real *massFrac);
void state_stats (amrex::MultiFab& S);
amrex::Real getMFsum(amrex::MultiFab &a_MF, int comp);
amrex::Real getCellsCount();
void zeroBoundaryVisc (amrex::MultiFab* beta[AMREX_SPACEDIM],
const amrex::Real time,
const int state_comp,
const int dst_comp,
const int ncomp) const;
void diffuse_scalar_fj (const amrex::Vector<amrex::MultiFab*>& S_old,
const amrex::Vector<amrex::MultiFab*>& Rho_old,
amrex::Vector<amrex::MultiFab*>& S_new,
const amrex::Vector<amrex::MultiFab*>& Rho_new,
int S_comp,
int num_comp,
int Rho_comp,
amrex::Real prev_time,
amrex::Real curr_time,
amrex::Real be_cn_theta,
const amrex::MultiFab& rho_half,
int rho_flag,
amrex::MultiFab* const* fluxn,
amrex::MultiFab* const* fluxnp1,
int fluxComp,
amrex::MultiFab* delta_rhs,
int rhsComp,
const amrex::MultiFab* alpha,
int alphaComp,
const amrex::MultiFab* const* betan,
const amrex::MultiFab* const* betanp1,
int betaComp,
const amrex::Vector<amrex::Real>& visc_coef,
int visc_coef_comp,
const amrex::MultiFab& volume,
const amrex::MultiFab* const* area,
const amrex::IntVect& cratio,
const amrex::BCRec& bc,
const amrex::Geometry& geom,
bool add_hoop_stress,
bool add_old_time_divFlux = true,
const amrex::Vector<int>& is_diffusive = amrex::Vector<int>());
void diffusionFJDriver(amrex::ForkJoin& fj,
amrex::Real prev_time,
amrex::Real curr_time,
amrex::Real be_cn_theta,
int rho_flag,
const amrex::Vector<amrex::Real>& visc_coef,
int visc_coef_comp,
const amrex::IntVect& cratio,
const amrex::BCRec& bc,
const amrex::Geometry& geom,
bool add_hoop_stress,
bool add_old_time_divFlux,
const amrex::Vector<int>& is_diffusive,
bool has_coarse_data,
bool has_delta_rhs,
bool has_alpha_in,
bool has_betan,
bool has_betanp1);
//
// Functions for interpolating from cell centers to cell edges
//
enum FPLoc { HT_Edge = 0, HT_Center };
//
static FPLoc fpi_phys_loc (int p_bc);
//
static void center_to_edge_fancy (const amrex::FArrayBox& cfab,
amrex::FArrayBox& efab,
const amrex::Box& ccBox,
const amrex::Box& eBox,
int sComp,
int dComp,
int nComp,
const amrex::Box& domain,
const FPLoc& bc_lo,
const FPLoc& bc_hi);
static void init_once ();
static void init_transport (int iEG);
void reactionRateRhoY_pphys(amrex::FArrayBox& RhoYdot,
const amrex::FArrayBox& RhoY,
const amrex::FArrayBox& RhoH,
const amrex::FArrayBox& T,
const amrex::BaseFab<int>& mask,
const amrex::Box& box,
int sCompRhoY,
int sCompRhoH,
int sCompT,
int sCompRhoYdot) const;
void getHmixGivenTY_pphys (amrex::FArrayBox& hmix,
const amrex::FArrayBox& T,
const amrex::FArrayBox& Y,
const amrex::Box& box,
int sCompT,
int sCompY,
int sCompH) const;
static int getTGivenHY_pphys (amrex::FArrayBox& T,
const amrex::FArrayBox& H,
const amrex::FArrayBox& Y,
const amrex::Box& box,
int sCompH,
int sCompY,
int sCompT,
const amrex::Real& errMAX = -1);
void HfromT_pphys(const int* lo, const int* hi,
amrex::Real* H, ARLIM_P(Hlo), ARLIM_P(Hhi),
const amrex::Real* T, ARLIM_P(Tlo), ARLIM_P(Thi));
static void RhoH_to_Temp (amrex::MultiFab& S,
int nGrow = 0,
int dominmax = false);
////////////////////////////////////////////////////////////////////////////
// Private Data //
////////////////////////////////////////////////////////////////////////////
static amrex::Vector<amrex::AMRErrorTag> errtags;
amrex::Vector<std::unique_ptr<FluxBoxes> > raii_fbs;
amrex::MultiFab** EdgeState;
amrex::MultiFab** EdgeFlux;
amrex::MultiFab** SpecDiffusionFluxn;
amrex::MultiFab** SpecDiffusionFluxnp1;
amrex::MultiFab** SpecDiffusionFluxWbar;
static bool plot_reactions;
static bool plot_consumption;
static bool plot_heat_release;
std::map<std::string,std::unique_ptr<amrex::MultiFab> > auxDiag;
static std::map<std::string,amrex::Vector<std::string> > auxDiag_names;
bool updateFluxReg;
bool is_predictor;
// these refer to the old and new-time ambient pressure for level 0
static amrex::Real p_amb_old;
static amrex::Real p_amb_new;
static amrex::Real dp0dt;
static amrex::Real thetabar;
static amrex::Real lev0cellCount;
static int closed_chamber;
static amrex::Real dpdt_factor;
// Active control data
static bool ctrl_active;
static bool ctrl_use_temp;
static amrex::Real ctrl_tauControl;
static amrex::Real ctrl_cfix;
static amrex::Real ctrl_coftOld;
static amrex::Real ctrl_sest;
static amrex::Real ctrl_corr;
static amrex::Real ctrl_V_in;
static amrex::Real ctrl_V_in_old;
static amrex::Real ctrl_changeMax;
static amrex::Real ctrl_tBase;
static amrex::Real ctrl_dV;
static amrex::Real ctrl_scale;
static amrex::Real ctrl_zBase;
static amrex::Real ctrl_h;
static amrex::Real ctrl_velMax;
static amrex::Real ctrl_temperature;
static int ctrl_verbose;
static int ctrl_NavgPts;
static int ctrl_nfilled;
static int ctrl_flameDir;
static int ctrl_pseudoGravity;
static int ctrl_method;
static std::string ctrl_AChistory;
static amrex::Vector<amrex::Real> ctrl_time_pts;
static amrex::Vector<amrex::Real> ctrl_velo_pts;
static amrex::Vector<amrex::Real> ctrl_cntl_pts;
static ProbParm* prob_parm;
static ProbParm* prob_parm_d;
static ACParm* ac_parm;
static ACParm* ac_parm_d;
static pele::physics::transport::TransportParams<
pele::physics::PhysicsType::transport_type>
trans_parms;
static pele::physics::PMF::PmfData pmf_data;
// Reactor. TODO: let's make those non-static to have per level reactors
static std::string chem_integrator;
static std::unique_ptr<pele::physics::reactions::ReactorBase> m_reactor;
//
// Protected static data.
//
static int num_divu_iters;
static int init_once_done;
static int nSpecGroup;
static int RhoH;
static int do_OT_radiation;
static int do_heat_sink;
static int unity_Le;
static int use_wbar;
static int FuncCount_Type;
static int RhoYdot_Type;
#ifdef AMREX_PARTICLES
static int spraydot_Type;
#endif
static int divu_ceiling;
static amrex::Real min_rho_divu_ceiling;
static amrex::Real divu_dt_factor;
static amrex::Real htt_tempmin;
static amrex::Real htt_tempmax;
static amrex::Real htt_hmixTYP;
static amrex::Real typical_Y_val_min;
static int RhoRT;
static int have_rhort;
static int first_spec;
static int last_spec;
static bool solve_passives;
static int first_passive;
static amrex::Vector<std::string> spec_names;
static int floor_species;
static amrex::Real prandtl;
static amrex::Real schmidt;
static amrex::Real constant_thick_val;
static int do_diffuse_sync;
static int do_reflux_visc;
static int do_set_rho_to_species_sum;
static int clipSpeciesOnRegrid;
static int zeroBndryVisc;
static int do_check_divudt;
static int hack_nochem;
static int hack_nospecdiff;
static int hack_noavgdivu;
static amrex::Real trac_diff_coef;
static std::string turbFile;
static std::string fuelName;
static std::string productName;
static amrex::Vector<std::string> consumptionName;
static int calcDiffusivity_count;
static amrex::Vector<amrex::Real> typical_values;
static bool avg_down_chem;
static int reset_typical_vals_int;
static std::map<std::string,amrex::Real> typical_values_FileVals;
static int sdc_iterMAX;
static int num_mac_sync_iter;
static int syncEntireHierarchy;
static int deltaT_verbose;
static int deltaT_crashOnConvFail;
static int iter_debug;
static int mHtoTiterMAX;
static amrex::Vector<amrex::Real> mTmpData;
#ifdef AMREX_USE_EB
static std::string diffusion_redistribution_type;
#endif
static int ncells_chem;
static bool use_typ_vals_chem;
static bool def_harm_avg_cen2edge;
static int nGrowAdvForcing;
static int nGrowDivU;
};
void pelelm_cc_ext_fill (amrex::Box const& bx, amrex::FArrayBox& data,
const int dcomp, const int numcomp,
amrex::Geometry const& geom, const amrex::Real time,
const amrex::Vector<amrex::BCRec>& bcr, const int bcomp,
const int scomp);
void pelelm_press_fill (amrex::Box const& bx, amrex::FArrayBox& data,
const int dcomp, const int numcomp,
amrex::Geometry const& geom, const amrex::Real time,
const amrex::Vector<amrex::BCRec>& bcr, const int bcomp,
const int scomp);
void pelelm_face_fill (amrex::Box const& bx, amrex::FArrayBox& data,
const int dcomp, const int numcomp,
amrex::Geometry const& geom, const amrex::Real time,
const amrex::Vector<amrex::BCRec>& bcr, const int bcomp,
const int scomp);
void pelelm_dummy_fill (amrex::Box const& bx, amrex::FArrayBox& data,
const int dcomp, const int numcomp,
amrex::Geometry const& geom, const amrex::Real time,
const amrex::Vector<amrex::BCRec>& bcr, const int bcomp,
const int scomp);
#endif /*_PeleLM_H_*/