ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF.H
Go to the documentation of this file.
1 #ifndef ERF_H_
2 #define ERF_H_
3 
4 #include <string>
5 #include <limits>
6 #include <memory>
7 
8 #ifdef _OPENMP
9 #include <omp.h>
10 #endif
11 
12 #include <AMReX_AmrCore.H>
13 #include <AMReX_BCRec.H>
14 #include <AMReX_Print.H>
15 
16 #include <AMReX_ParallelDescriptor.H>
17 #include <AMReX_ParmParse.H>
18 #include <AMReX_MultiFabUtil.H>
19 #include <AMReX_FillPatchUtil.H>
20 #include <AMReX_VisMF.H>
21 #include <AMReX_PhysBCFunct.H>
22 #include <AMReX_YAFluxRegister.H>
23 #include <AMReX_ErrorList.H>
24 #include <AMReX_EBFabFactory.H>
25 #include <AMReX_EBMultiFabUtil.H>
26 
27 #include <ERF_EB.H>
28 
29 #ifdef ERF_USE_FFT
30 #include <AMReX_FFT_Poisson.H>
31 #endif
32 
33 #ifdef AMREX_MEM_PROFILING
34 #include <AMReX_MemProfiler.H>
35 #endif
36 
37 #include "ERF_ProbCommon.H"
38 
39 #include <ERF_IndexDefines.H>
40 #include <ERF_DataStruct.H>
41 #include <ERF_TurbPertStruct.H>
42 #include <ERF_InputSoundingData.H>
43 #include <ERF_InputSpongeData.H>
44 #include <ERF_ABLMost.H>
45 #include <ERF_Derive.H>
46 #include <ERF_ReadBndryPlanes.H>
47 #include <ERF_WriteBndryPlanes.H>
48 #include <ERF_MRI.H>
49 #include <ERF_PhysBCFunct.H>
50 #include <ERF_FillPatcher.H>
51 #include <ERF_SampleData.H>
52 #include <ERF_ForestDrag.H>
53 
54 #ifdef ERF_USE_PARTICLES
55 #include "ERF_ParticleData.H"
56 #endif
57 
60 #include "ERF_LandSurface.H"
61 #ifdef ERF_USE_WINDFARM
62 #include "ERF_WindFarm.H"
63 #endif
64 
65 #ifdef ERF_USE_RRTMGP
66 #include "ERF_Radiation.H"
67 #endif
68 
69 #include <iostream>
70 
71 #ifdef AMREX_LAZY
72 #include <AMReX_Lazy.H>
73 #endif
74 
75 #ifdef ERF_USE_MULTIBLOCK
77 #endif
78 
79 /**
80  * Enum of possible initialization types
81 */
82 AMREX_ENUM(InitType,
83  None, Input_Sounding, Ideal, Real, Metgrid, Uniform
84 );
85 
86 /**
87  * Enum of possible interpolation types between coarse/fine
88 */
89 AMREX_ENUM(StateInterpType,
90  FullState, Perturbational
91 );
92 
93 /**
94  * Enum of possible plotfile types
95 */
96 AMREX_ENUM(PlotFileType,
97  None, Amrex, Netcdf
98 );
99 
100 
101 #if 0
102 /**
103  * Enum of possible coarse/fine interpolation options
104 */
105 namespace InterpType {
106  enum {
107  PCInterp = 0,
108  NodeBilinear,
109  CellConservativeLinear,
110  CellBilinear,
111  CellQuadratic,
112  CellConservativeProtected,
113  CellConservativeQuartic,
114  FaceLinear,
115  FaceConserativeLinear
116  };
117 }
118 #endif
119 
120 /**
121  * Main class in ERF code, instantiated from main.cpp
122 */
123 
124 class ERF
125  : public amrex::AmrCore
126 {
127 public:
128 
129  ////////////////
130  // public member functions
131 
132  // constructor - reads in parameters from inputs file
133  // - sizes multilevel arrays and data structures
134  ERF ();
135  ~ERF () override;
136 
137  void ERF_shared ();
138 
139  // Declare a default move constructor so we ensure the destructor is
140  // not called when we return an object of this class by value
141  ERF (ERF&&) noexcept = delete;
142 
143  // Declare a default move assignment operator
144  ERF& operator=(ERF&& other) noexcept = delete;
145 
146  // Delete the copy constructor
147  ERF (const ERF& other) = delete;
148  //
149  // Delete the copy assignment operator
150  ERF& operator=(const ERF& other) = delete;
151 
152  // Advance solution to final time
153  void Evolve ();
154 
155  // Tag cells for refinement
156  void ErrorEst (int lev, amrex::TagBoxArray& tags, amrex::Real time, int ngrow) override;
157 
158  // Initialize multilevel data
159  void InitData ();
160 
161  // Initialize multilevel data before MultiBlock
162  void InitData_pre ();
163 
164  // Initialize multilevel data after MultiBlock
165  void InitData_post ();
166 
167  void WriteMyEBSurface ();
168 
169  // Compute the divergence -- whether EB, no-terrain, flat terrain or general terrain
170  void compute_divergence (int lev, amrex::MultiFab& rhs, amrex::Array<amrex::MultiFab const*,AMREX_SPACEDIM> rho0_u_const,
171  amrex::Geometry const& geom_at_lev);
172 
173  // Project the velocities to be divergence-free -- this is only relevant if anelastic == 1
174  void project_velocities (int lev, amrex::Real dt, amrex::Vector<amrex::MultiFab >& vars,
175  amrex::MultiFab& p);
176 
177  // Project the velocities to be divergence-free with a thin body
178  void project_velocities_tb (int lev, amrex::Real dt, amrex::Vector<amrex::MultiFab >& vars,
179  amrex::MultiFab& p);
180 
181  // Calculate wall distance by solving a Poisson equation
182  void poisson_wall_dist (int lev);
183 
184 #ifdef ERF_USE_FFT
185  void solve_with_fft (int lev, amrex::MultiFab& rhs, amrex::MultiFab& p,
186  amrex::Array<amrex::MultiFab,AMREX_SPACEDIM>& fluxes);
187 #endif
188  void solve_with_EB_mlmg (int lev, amrex::Vector<amrex::MultiFab>& rhs, amrex::Vector<amrex::MultiFab>& p,
189  amrex::Vector<amrex::Array<amrex::MultiFab,AMREX_SPACEDIM>>& fluxes);
190 
191  void solve_with_gmres (int lev, amrex::Vector<amrex::MultiFab>& rhs, amrex::Vector<amrex::MultiFab>& p,
192  amrex::Vector<amrex::Array<amrex::MultiFab,AMREX_SPACEDIM>>& fluxes);
193  void solve_with_mlmg (int lev, amrex::Vector<amrex::MultiFab>& rhs, amrex::Vector<amrex::MultiFab>& p,
194  amrex::Vector<amrex::Array<amrex::MultiFab,AMREX_SPACEDIM>>& fluxes);
195 
196  // Define the projection bc's based on the domain bc types
197  amrex::Array<amrex::LinOpBCType,AMREX_SPACEDIM>
198  get_projection_bc (amrex::Orientation::Side side) const noexcept;
199  bool projection_has_dirichlet (amrex::Array<amrex::LinOpBCType,AMREX_SPACEDIM> bcs) const;
200 
201  // Init (NOT restart or regrid)
202  void init_only (int lev, amrex::Real time);
203 
204  // Restart
205  void restart ();
206 
207  // Is it time to write a plotfile or checkpoint?
208  bool writeNow (const amrex::Real cur_time, const amrex::Real dt, const int nstep,
209  const int plot_int, const amrex::Real plot_per);
210 
211  // Called after every level 0 timestep
212  void post_timestep (int nstep, amrex::Real time, amrex::Real dt_lev);
213 
214  // Diagnostics
215  void sum_integrated_quantities (amrex::Real time);
216  void write_1D_profiles (amrex::Real time);
217  void write_1D_profiles_stag (amrex::Real time);
218 
219  amrex::Real cloud_fraction (amrex::Real time);
220 
221  // Fill the physical boundary conditions for cell-centered velocity (diagnostic only)
222  void FillBdyCCVels (amrex::Vector<amrex::MultiFab>& mf_cc_vel);
223 
224  void sample_points (int lev, amrex::Real time, amrex::IntVect cell, amrex::MultiFab& mf);
225  void sample_lines (int lev, amrex::Real time, amrex::IntVect cell, amrex::MultiFab& mf);
226 
227  void derive_diag_profiles (amrex::Real time,
228  amrex::Gpu::HostVector<amrex::Real>& h_avg_u , amrex::Gpu::HostVector<amrex::Real>& h_avg_v,
229  amrex::Gpu::HostVector<amrex::Real>& h_avg_w , amrex::Gpu::HostVector<amrex::Real>& h_avg_rho,
230  amrex::Gpu::HostVector<amrex::Real>& h_avg_th , amrex::Gpu::HostVector<amrex::Real>& h_avg_ksgs,
231  amrex::Gpu::HostVector<amrex::Real>& h_avg_Kmv , amrex::Gpu::HostVector<amrex::Real>& h_avg_Khv,
232  amrex::Gpu::HostVector<amrex::Real>& h_avg_qv , amrex::Gpu::HostVector<amrex::Real>& h_avg_qc,
233  amrex::Gpu::HostVector<amrex::Real>& h_avg_qr ,
234  amrex::Gpu::HostVector<amrex::Real>& h_avg_wqv , amrex::Gpu::HostVector<amrex::Real>& h_avg_wqc,
235  amrex::Gpu::HostVector<amrex::Real>& h_avg_wqr ,
236  amrex::Gpu::HostVector<amrex::Real>& h_avg_qi , amrex::Gpu::HostVector<amrex::Real>& h_avg_qs,
237  amrex::Gpu::HostVector<amrex::Real>& h_avg_qg ,
238  amrex::Gpu::HostVector<amrex::Real>& h_avg_uu , amrex::Gpu::HostVector<amrex::Real>& h_avg_uv,
239  amrex::Gpu::HostVector<amrex::Real>& h_avg_uw,
240  amrex::Gpu::HostVector<amrex::Real>& h_avg_vv , amrex::Gpu::HostVector<amrex::Real>& h_avg_vw,
241  amrex::Gpu::HostVector<amrex::Real>& h_avg_ww,
242  amrex::Gpu::HostVector<amrex::Real>& h_avg_uth , amrex::Gpu::HostVector<amrex::Real>& h_avg_vth,
243  amrex::Gpu::HostVector<amrex::Real>& h_avg_wth, amrex::Gpu::HostVector<amrex::Real>& h_avg_thth,
244  amrex::Gpu::HostVector<amrex::Real>& h_avg_ku, amrex::Gpu::HostVector<amrex::Real>& h_avg_kv,
245  amrex::Gpu::HostVector<amrex::Real>& h_avg_kw,
246  amrex::Gpu::HostVector<amrex::Real>& h_avg_p,
247  amrex::Gpu::HostVector<amrex::Real>& h_avg_pu, amrex::Gpu::HostVector<amrex::Real>& h_avg_pv,
248  amrex::Gpu::HostVector<amrex::Real>& h_avg_pw, amrex::Gpu::HostVector<amrex::Real>& h_avg_wthv);
249  void derive_diag_profiles_stag (amrex::Real time,
250  amrex::Gpu::HostVector<amrex::Real>& h_avg_u , amrex::Gpu::HostVector<amrex::Real>& h_avg_v,
251  amrex::Gpu::HostVector<amrex::Real>& h_avg_w , amrex::Gpu::HostVector<amrex::Real>& h_avg_rho,
252  amrex::Gpu::HostVector<amrex::Real>& h_avg_th , amrex::Gpu::HostVector<amrex::Real>& h_avg_ksgs,
253  amrex::Gpu::HostVector<amrex::Real>& h_avg_Kmv , amrex::Gpu::HostVector<amrex::Real>& h_avg_Khv,
254  amrex::Gpu::HostVector<amrex::Real>& h_avg_qv , amrex::Gpu::HostVector<amrex::Real>& h_avg_qc,
255  amrex::Gpu::HostVector<amrex::Real>& h_avg_qr ,
256  amrex::Gpu::HostVector<amrex::Real>& h_avg_wqv , amrex::Gpu::HostVector<amrex::Real>& h_avg_wqc,
257  amrex::Gpu::HostVector<amrex::Real>& h_avg_wqr ,
258  amrex::Gpu::HostVector<amrex::Real>& h_avg_qi , amrex::Gpu::HostVector<amrex::Real>& h_avg_qs,
259  amrex::Gpu::HostVector<amrex::Real>& h_avg_qg ,
260  amrex::Gpu::HostVector<amrex::Real>& h_avg_uu , amrex::Gpu::HostVector<amrex::Real>& h_avg_uv,
261  amrex::Gpu::HostVector<amrex::Real>& h_avg_uw,
262  amrex::Gpu::HostVector<amrex::Real>& h_avg_vv , amrex::Gpu::HostVector<amrex::Real>& h_avg_vw,
263  amrex::Gpu::HostVector<amrex::Real>& h_avg_ww,
264  amrex::Gpu::HostVector<amrex::Real>& h_avg_uth , amrex::Gpu::HostVector<amrex::Real>& h_avg_vth,
265  amrex::Gpu::HostVector<amrex::Real>& h_avg_wth, amrex::Gpu::HostVector<amrex::Real>& h_avg_thth,
266  amrex::Gpu::HostVector<amrex::Real>& h_avg_ku, amrex::Gpu::HostVector<amrex::Real>& h_avg_kv,
267  amrex::Gpu::HostVector<amrex::Real>& h_avg_kw,
268  amrex::Gpu::HostVector<amrex::Real>& h_avg_p,
269  amrex::Gpu::HostVector<amrex::Real>& h_avg_pu, amrex::Gpu::HostVector<amrex::Real>& h_avg_pv,
270  amrex::Gpu::HostVector<amrex::Real>& h_avg_pw, amrex::Gpu::HostVector<amrex::Real>& h_avg_wthv);
271 
272  void derive_stress_profiles (amrex::Gpu::HostVector<amrex::Real>& h_avg_tau11, amrex::Gpu::HostVector<amrex::Real>& h_avg_tau12,
273  amrex::Gpu::HostVector<amrex::Real>& h_avg_tau13, amrex::Gpu::HostVector<amrex::Real>& h_avg_tau22,
274  amrex::Gpu::HostVector<amrex::Real>& h_avg_tau23, amrex::Gpu::HostVector<amrex::Real>& h_avg_tau33,
275  amrex::Gpu::HostVector<amrex::Real>& h_avg_hfx3, amrex::Gpu::HostVector<amrex::Real>& h_avg_q1fx3,
276  amrex::Gpu::HostVector<amrex::Real>& h_avg_q2fx3, amrex::Gpu::HostVector<amrex::Real>& h_avg_diss);
277  void derive_stress_profiles_stag (amrex::Gpu::HostVector<amrex::Real>& h_avg_tau11, amrex::Gpu::HostVector<amrex::Real>& h_avg_tau12,
278  amrex::Gpu::HostVector<amrex::Real>& h_avg_tau13, amrex::Gpu::HostVector<amrex::Real>& h_avg_tau22,
279  amrex::Gpu::HostVector<amrex::Real>& h_avg_tau23, amrex::Gpu::HostVector<amrex::Real>& h_avg_tau33,
280  amrex::Gpu::HostVector<amrex::Real>& h_avg_hfx3, amrex::Gpu::HostVector<amrex::Real>& h_avg_q1fx3,
281  amrex::Gpu::HostVector<amrex::Real>& h_avg_q2fx3, amrex::Gpu::HostVector<amrex::Real>& h_avg_diss);
282 
283  // Perform the volume-weighted sum
284  amrex::Real
285  volWgtSumMF (int lev, const amrex::MultiFab& mf, int comp,
286  const amrex::MultiFab& mapfac,
287  bool local, bool finemask);
288 
289  // Decide if it is time to take an action
290  static bool is_it_time_for_action (int nstep, amrex::Real time, amrex::Real dt,
291  int action_interval, amrex::Real action_per);
292 
293  // Make a new level using provided BoxArray and DistributionMapping and
294  // fill with interpolated coarse level data.
295  // overrides the pure virtual function in AmrCore
296  void MakeNewLevelFromCoarse (int lev, amrex::Real time, const amrex::BoxArray& ba,
297  const amrex::DistributionMapping& dm) override;
298 
299  // Remake an existing level using provided BoxArray and DistributionMapping and
300  // fill with existing fine and coarse data.
301  // overrides the pure virtual function in AmrCore
302  void RemakeLevel (int lev, amrex::Real time, const amrex::BoxArray& ba,
303  const amrex::DistributionMapping& dm) override;
304 
305  // Delete level data
306  // overrides the pure virtual function in AmrCore
307  void ClearLevel (int lev) override;
308 
309  // Make a new level from scratch using provided BoxArray and DistributionMapping.
310  // Only used during initialization.
311  // overrides the pure virtual function in AmrCore
312  void MakeNewLevelFromScratch (int lev, amrex::Real time, const amrex::BoxArray& ba,
313  const amrex::DistributionMapping& dm) override;
314 
315  // compute dt from CFL considerations
316  amrex::Real estTimeStep (int lev, long& dt_fast_ratio) const;
317 
318 #ifdef ERF_USE_WW3_COUPLING
319  //amrex::Print() << " About to call send_to_ww3 from ERF.H" << std::endl;
320  void send_to_ww3(int lev);
321  //amrex::Print() << " About to call read_waves from ERF.H" << std::endl;
322  void read_waves(int lev);
323  //void send_to_ww3 (int lev);
324  //void read_waves (int lev);
325 
326  //void send_to_ww3 (int lev);
327 #endif
328 
329  // Interface for advancing the data at one level by one "slow" timestep
330  void advance_dycore (int level,
331  amrex::Vector<amrex::MultiFab>& state_old,
332  amrex::Vector<amrex::MultiFab>& state_new,
333  amrex::MultiFab& xvel_old, amrex::MultiFab& yvel_old, amrex::MultiFab& zvel_old,
334  amrex::MultiFab& xvel_new, amrex::MultiFab& yvel_new, amrex::MultiFab& zvel_new,
335  amrex::MultiFab& source, amrex::MultiFab& xmom_src,
336  amrex::MultiFab& ymom_src, amrex::MultiFab& zmom_src,
337  amrex::Geometry fine_geom,
338  amrex::Real dt, amrex::Real time);
339 
340  void advance_microphysics (int lev,
341  amrex::MultiFab& cons_in,
342  const amrex::Real& dt_advance,
343  const int& iteration,
344  const amrex::Real& time);
345 
346  void advance_lsm (int lev,
347  amrex::MultiFab& /*cons_in*/,
348  const amrex::Real& dt_advance);
349 
350 #if defined(ERF_USE_RRTMGP)
351  void advance_radiation (int lev,
352  amrex::MultiFab& cons_in,
353  const amrex::Real& dt_advance);
354 #endif
355 
356  amrex::MultiFab& build_fine_mask (int lev);
357 
358  void MakeHorizontalAverages ();
359  void MakeDiagnosticAverage (amrex::Vector<amrex::Real>& h_havg, amrex::MultiFab& S, int n);
360  void derive_upwp (amrex::Vector<amrex::Real>& h_havg);
361 
362  // write plotfile to disk
363  void WritePlotFile (int which, PlotFileType plotfile_type, amrex::Vector<std::string> plot_var_names);
364 
365  void WriteMultiLevelPlotfileWithTerrain (const std::string &plotfilename,
366  int nlevels,
367  const amrex::Vector<const amrex::MultiFab*> &mf,
368  const amrex::Vector<const amrex::MultiFab*> &mf_nd,
369  const amrex::Vector<std::string> &varnames,
370  const amrex::Vector<amrex::Geometry>& my_geom,
371  amrex::Real time,
372  const amrex::Vector<int> &level_steps,
373  const amrex::Vector<amrex::IntVect>& my_ref_ratio,
374  const std::string &versionName = "HyperCLaw-V1.1",
375  const std::string &levelPrefix = "Level_",
376  const std::string &mfPrefix = "Cell",
377  const amrex::Vector<std::string>& extra_dirs = amrex::Vector<std::string>()) const;
378 
379 
380  void WriteGenericPlotfileHeaderWithTerrain (std::ostream &HeaderFile,
381  int nlevels,
382  const amrex::Vector<amrex::BoxArray> &bArray,
383  const amrex::Vector<std::string> &varnames,
384  const amrex::Vector<amrex::Geometry>& my_geom,
385  amrex::Real time,
386  const amrex::Vector<int> &level_steps,
387  const amrex::Vector<amrex::IntVect>& my_ref_ratio,
388  const std::string &versionName,
389  const std::string &levelPrefix,
390  const std::string &mfPrefix) const;
391 
392  void erf_enforce_hse (int lev,
393  amrex::MultiFab& dens, amrex::MultiFab& pres, amrex::MultiFab& pi, amrex::MultiFab& th,
394  std::unique_ptr<amrex::MultiFab>& z_cc);
395 
396 #ifdef ERF_USE_NETCDF
397  //! Write a timestep to 1D vertical column output for coupling
398  void writeToNCColumnFile (int lev,
399  const std::string& colfile_name, amrex::Real xloc, amrex::Real yloc,
400  amrex::Real time);
401 #endif //ERF_USE_NETCDF
402 
403  void init_from_input_sounding (int lev);
404 
405  void input_sponge (int lev);
406 
407  void init_from_hse (int lev);
408 
409  void init_thin_body (int lev, const amrex::BoxArray& ba, const amrex::DistributionMapping& dm);
410 
411 #ifdef ERF_USE_MULTIBLOCK
412  // constructor used when ERF is created by a multiblock driver
413  // calls AmrCore constructor -> AmrMesh constructor
414  ERF (const amrex::RealBox& rb, int max_level_in,
415  const amrex::Vector<int>& n_cell_in, int coord,
416  const amrex::Vector<amrex::IntVect>& ref_ratio,
417  const amrex::Array<int,AMREX_SPACEDIM>& is_per,
418  std::string prefix);
419 
420  // Advance a block specified number of time steps
421  void Evolve_MB (int MBstep, int max_block_step);
422 
423  // get the current time values and dt
424  amrex::Real get_t_old (int lev) { return t_old[lev]; }
425  amrex::Real get_t_new (int lev) { return t_new[lev]; }
426  amrex::Real get_dt (int lev) { return dt[lev]; }
427 
428  // Set parmparse prefix for MultiBlock
429  void SetParmParsePrefix (std::string name) { pp_prefix = name; }
430 
431  // Set 'this' from multiblock container
432  void SetMultiBlockPointer (MultiBlockContainer *mbc) { m_mbc = mbc; }
433 
434  // Public data copy for MB
435  std::vector<amrex::Box> domain_p;
436  amrex::Vector<amrex::Vector<amrex::MultiFab> > vars_new;
437  amrex::Vector<amrex::Vector<amrex::MultiFab> > vars_old;
438 
439  // Velocity time averaged field
440  amrex::Vector<std::unique_ptr<amrex::MultiFab>> vel_t_avg;
441  amrex::Vector<amrex::Real> t_avg_cnt;
442 #endif
443 
444  std::string pp_prefix {"erf"};
445 
446  void fill_from_bndryregs (const amrex::Vector<amrex::MultiFab*>& mfs,
447  amrex::Real time);
448 
449 #ifdef ERF_USE_NETCDF
450  void fill_from_realbdy (const amrex::Vector<amrex::MultiFab*>& mfs,
451  amrex::Real time,
452  bool cons_only,
453  int icomp_cons,
454  int ncomp_cons,
455  amrex::IntVect ngvect_cons,
456  amrex::IntVect ngvect_vels);
457 #endif
458 
459 #ifdef ERF_USE_NETCDF
460  void init_from_wrfinput (int lev);
461  void init_from_metgrid (int lev);
462 #endif // ERF_USE_NETCDF
463 
464 #ifdef ERF_USE_WINDFARM
465  void init_windfarm(int lev);
466  void advance_windfarm (const amrex::Geometry& a_geom,
467  const amrex::Real& dt_advance,
468  amrex::MultiFab& cons_in,
469  amrex::MultiFab& U_old, amrex::MultiFab& V_old, amrex::MultiFab& W_old,
470  amrex::MultiFab& mf_vars_windfarm,
471  const amrex::MultiFab& mf_Nturb,
472  const amrex::MultiFab& mf_SMark,
473  const amrex::Real& time);
474 #endif
475 
476  void MakeEBGeometry ();
477  void make_eb_box ();
479  void redistribute_term (int lev, int ncomp,
480  amrex::MultiFab& result,
481  amrex::MultiFab& result_tmp,
482  amrex::MultiFab const& state,
483  amrex::BCRec const* bc,
484  amrex::Real const dt);
485  void redistribute_term (amrex::MFIter const& mfi, int lev, int ncomp,
486  amrex::MultiFab& result,
487  amrex::MultiFab& result_tmp,
488  amrex::MultiFab const& state,
489  amrex::BCRec const* bc,
490  amrex::Real const dt);
491 
492  // more flexible version of AverageDown() that lets you average down across multiple levels
493  void AverageDownTo (int crse_lev, int scomp, int ncomp); // NOLINT
494 
495  // write checkpoint file to disk
496  void WriteCheckpointFile () const;
497 
498  // read checkpoint file from disk
499  void ReadCheckpointFile ();
500 
501  // read checkpoint file from disk -- called after instantiating m_most
502  void ReadCheckpointFileMOST ();
503 
504  void init_zphys (int lev, amrex::Real time);
505  void remake_zphys (int lev, amrex::Real time, std::unique_ptr<amrex::MultiFab>& temp_zphys_nd);
506  void update_terrain_arrays (int lev);
507 
508 private:
509 
510  ///////////////////////////
511  // private member functions
512  ///////////////////////////
513 
514  // read in some parameters from inputs file
515  void ReadParameters ();
516  void ParameterSanityChecks ();
517 
518  // set covered coarse cells/faces to be the average of overlying fine cells/faces
519  void AverageDown ();
520 
521  void update_diffusive_arrays (int lev, const amrex::BoxArray& ba, const amrex::DistributionMapping& dm);
522 
523  void Construct_ERFFillPatchers (int lev);
524 
525  void Define_ERFFillPatchers (int lev);
526 
527  void init1DArrays ();
528 
529  void init_bcs ();
530 
531  void init_custom (int lev);
532 
533  void init_uniform (int lev);
534 
535  void init_stuff (int lev, const amrex::BoxArray& ba, const amrex::DistributionMapping& dm,
536  amrex::Vector<amrex::MultiFab>& lev_new, amrex::Vector<amrex::MultiFab>& lev_old,
537  amrex::MultiFab& tmp_base_state,
538  std::unique_ptr<amrex::MultiFab>& tmp_zphys_nd);
539 
540  // Initialize the Turbulent perturbation
541  void turbPert_update (const int lev, const amrex::Real dt);
542  void turbPert_amplitude (const int lev);
543 
544  // Initialize the integrator object
545  void initialize_integrator (int lev, amrex::MultiFab& cons_mf, amrex::MultiFab& vel_mf);
546 
547  // Create the physbcs objects
548  void make_physbcs (int lev);
549 
550  // Initialize the microphysics object
551  void initializeMicrophysics (const int&);
552 
553 #ifdef ERF_USE_WINDFARM
554  // Initialize the windfarm object
555  void initializeWindFarm (const int&);
556 #endif
557 
558  // Compute a vector of new MultiFabs by copying from valid region and filling ghost cells
559  //
560  // NOTE: FillPatch takes in an empty MF, and returns cell-centered + velocities (not momenta)
561  //
562  // This one works only at level = 0 (base state does not change)
563  void FillPatch (int lev, amrex::Real time,
564  const amrex::Vector<amrex::MultiFab*>& mfs_vel,
565  bool cons_only=false);
566 
567  // This one works only at level > 0 (base state does change)
568  void FillPatch (int lev, amrex::Real time,
569  const amrex::Vector<amrex::MultiFab*>& mfs_vel,
570  const amrex::Vector<amrex::MultiFab*>& mfs_mom,
571  const amrex::MultiFab& old_base_state,
572  const amrex::MultiFab& new_base_state,
573  bool fillset=true, bool cons_only=false);
574 
575  // Compute new multifabs by copying data from valid region and filling ghost cells.
576  // Unlike FillPatch, FillIntermediatePatch will use the supplied multifabs instead of fine level data.
577  // This is to support filling boundary cells at an intermediate time between old/new times
578  // on the fine level when valid data at a specific time is already available (such as
579  // at each RK stage when integrating between initial and final times at a given level).
580  //
581  // NOTE: FillIntermediatePatch takes in updated momenta, and returns both updated velocity and momenta
582  //
583  void FillIntermediatePatch (int lev, amrex::Real time,
584  const amrex::Vector<amrex::MultiFab*>& mfs_vel,
585  const amrex::Vector<amrex::MultiFab*>& mfs_mom,
586  int ng_cons, int ng_vel, bool cons_only, int icomp_cons, int ncomp_cons,
587  bool allow_most_bcs = true);
588 
589  // Fill all multifabs (and all components) in a vector of multifabs corresponding to the
590  // grid variables defined in vars_old and vars_new just as FillCoarsePatch.
591  void FillCoarsePatch (int lev, amrex::Real time);
592 
593  // advance a level by dt
594  // includes a recursive call for finer levels
595  void timeStep (int lev, amrex::Real time, int iteration);
596 
597  // advance a single level for a single time step
598  void Advance (int lev, amrex::Real time, amrex::Real dt_lev, int iteration, int ncycle);
599 
600  //! Initialize HSE
601  void initHSE ();
602  void initHSE (int lev);
603 
604  //! Initialize Rayleigh damping profiles
605  void initRayleigh ();
606 
607  //! Initialize sponge profiles
608  void initSponge ();
609 
610  //! Set Rayleigh mean profiles from input sounding
611  void setRayleighRefFromSounding (bool restarting);
612 
613  //! Set sponge mean profiles from input sounding
614  void setSpongeRefFromSounding (bool restarting);
615 
616  // a wrapper for estTimeStep()
617  void ComputeDt (int step = -1);
618 
619  // get plotfile name
620  [[nodiscard]] std::string PlotFileName (int lev) const;
621 
622  // set plotfile variables names
623  static amrex::Vector<std::string> PlotFileVarNames (amrex::Vector<std::string> plot_var_names) ;
624 
625  // set which variables and derived quantities go into plotfiles
626  void setPlotVariables (const std::string& pp_plot_var_names, amrex::Vector<std::string>& plot_var_names);
627  // append variables to plot
628  void appendPlotVariables (const std::string& pp_plot_var_names, amrex::Vector<std::string>& plot_var_names);
629 
630 #ifdef ERF_USE_NETCDF
631  //! Write plotfile using NETCDF
632  void writeNCPlotFile (int lev, int which, const std::string& dir,
633  const amrex::Vector<const amrex::MultiFab*> &mf,
634  const amrex::Vector<std::string> &plot_var_names,
635  const amrex::Vector<int>& level_steps, amrex::Real time) const;
636 
637  //! Write checkpointFile using NetCdf
638  void WriteNCCheckpointFile () const;
639 
640  //! Read checkpointFile for restart
641  void ReadNCCheckpointFile ();
642 
643  //! Write MultiFab in NetCDF format
644  static void WriteNCMultiFab (const amrex::FabArray<amrex::FArrayBox> &fab,
645  const std::string& name,
646  bool set_ghost = false) ;
647 
648  //! Read MultiFab in NetCDF format
649  void ReadNCMultiFab (amrex::FabArray<amrex::FArrayBox> &mf,
650  const std::string &name,
651  int coordinatorProc = amrex::ParallelDescriptor::IOProcessorNumber(),
652  int allow_empty_mf = 0);
653 
654  //! Create 1D vertical column output for coupling
655  void createNCColumnFile (int lev,
656  const std::string& colfile_name, amrex::Real xloc, amrex::Real yloc);
657 
658  // Copy from the NC*fabs into the MultiFabs holding the boundary data
659  void init_from_wrfbdy (amrex::Vector<amrex::FArrayBox*> x_vel_lateral,
660  amrex::Vector<amrex::FArrayBox*> y_vel_lateral,
661  amrex::Vector<amrex::FArrayBox*> z_vel_lateral,
662  amrex::Vector<amrex::FArrayBox*> T_lateral);
663 
664  amrex::Real start_bdy_time;
665 
666  // *** *** FArrayBox's for holding the SURFACE data
667  // amrex::IArrayBox NC_IVGTYP_fab; // Vegetation type (IVGTYP); Discrete numbers;
668  // amrex::FArrayBox NC_z0_fab; // Surface Roughness, z0 = z0 (IVGTYP)
669  // amrex::FArrayBox NC_PSFC_fab; // Surface pressure
670 
671  // TODO: Clarify the relation between SST and TSK
672  // amrex::FArrayBox NC_SST_fab; // Sea Surface Temperature; Defined even for land area
673  // amrex::FArrayBox NC_TSK_fab; // Surface Skin Temperature; Appears to be same as SST...
674 
675  // Vectors (over time) of Vector (over variables) of FArrayBoxs for holding the data read from the wrfbdy NetCDF file
676  amrex::Vector<amrex::Vector<amrex::FArrayBox>> bdy_data_xlo;
677  amrex::Vector<amrex::Vector<amrex::FArrayBox>> bdy_data_xhi;
678  amrex::Vector<amrex::Vector<amrex::FArrayBox>> bdy_data_ylo;
679  amrex::Vector<amrex::Vector<amrex::FArrayBox>> bdy_data_yhi;
680 
681  amrex::Real bdy_time_interval;
682  amrex::Vector<std::unique_ptr<amrex::MultiFab>> lat_m, lon_m;
683 #endif // ERF_USE_NETCDF
684 
685  // Struct for working with the sounding data we take as an input
687 
688  // Struct for working with the sponge data we take as an input
690 
691  // Vector (6 planes) of DeviceVectors (ncell in plane) for Dirichlet BC data
692  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real>> xvel_bc_data;
693  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real>> yvel_bc_data;
694  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real>> zvel_bc_data;
695 
696  // Function to read and populate above vectors (if input file exists)
697  void init_Dirichlet_bc_data (const std::string input_file);
698 
699  // Read the file passed to amr.restart and use it as an initial condition for
700  // the current simulation. Supports a different number of components and
701  // ghost cells.
703 
704  // Initialize the new-time data at a level from the initial_data MultiFab
705  void InitializeLevelFromData (int lev, const amrex::MultiFab& initial_data);
706 
707  // utility to skip to next line in Header
708  static void GotoNextLine (std::istream& is);
709 
710  // Single level functions called by advance()
711  void post_update (amrex::MultiFab& state_mf, amrex::Real time, const amrex::Geometry& geom);
712  void fill_rhs (amrex::MultiFab& rhs_mf, const amrex::MultiFab& state_mf, amrex::Real time, const amrex::Geometry& geom);
713 
714  ////////////////
715  // private data members
716 
717  std::unique_ptr<ProblemBase> prob = nullptr; // problem-specific functions
718 
719  amrex::Vector<int> num_boxes_at_level; // how many boxes specified at each level by tagging criteria
720  amrex::Vector<int> num_files_at_level; // how many wrfinput files specified at each level
721  amrex::Vector<amrex::Vector<amrex::Box>> boxes_at_level; // the boxes specified at each level by tagging criteria
722 
723  amrex::Vector<int> istep; // which step?
724  amrex::Vector<int> nsubsteps; // how many substeps on each level?
725 
726  // keep track of old time, new time, and time step at each level
727  amrex::Vector<amrex::Real> t_new;
728  amrex::Vector<amrex::Real> t_old;
729  amrex::Vector<amrex::Real> dt;
730  amrex::Vector<long> dt_mri_ratio;
731 
732  // array of multifabs to store the solution at each level of refinement
733  // after advancing a level we use "swap".
734 #ifndef ERF_USE_MULTIBLOCK
735  amrex::Vector<amrex::Vector<amrex::MultiFab> > vars_new;
736  amrex::Vector<amrex::Vector<amrex::MultiFab> > vars_old;
737 
738  // Velocity time averaged field
739  amrex::Vector<std::unique_ptr<amrex::MultiFab>> vel_t_avg;
740  amrex::Vector<amrex::Real> t_avg_cnt;
741 #endif
742  amrex::Vector<std::unique_ptr<MRISplitIntegrator<amrex::Vector<amrex::MultiFab> > > > mri_integrator_mem;
743 
744  amrex::Vector<amrex::MultiFab> pp_inc;
745 
746  // Vector over levels of routines to impose physical boundary conditions
747  amrex::Vector<std::unique_ptr<ERFPhysBCFunct_cons>> physbcs_cons;
748  amrex::Vector<std::unique_ptr<ERFPhysBCFunct_u>> physbcs_u;
749  amrex::Vector<std::unique_ptr<ERFPhysBCFunct_v>> physbcs_v;
750  amrex::Vector<std::unique_ptr<ERFPhysBCFunct_w>> physbcs_w;
751  amrex::Vector<std::unique_ptr<ERFPhysBCFunct_base>> physbcs_base;
752 
753  // Store primitive variable for MOST BC
754  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Theta_prim;
755  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Qv_prim;
756  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Qr_prim;
757 
758  // Scratch space for time integrator
759  amrex::Vector<amrex::MultiFab> rU_old;
760  amrex::Vector<amrex::MultiFab> rU_new;
761  amrex::Vector<amrex::MultiFab> rV_old;
762  amrex::Vector<amrex::MultiFab> rV_new;
763  amrex::Vector<amrex::MultiFab> rW_old;
764  amrex::Vector<amrex::MultiFab> rW_new;
765 
766  // amrex::Vector<amrex::MultiFab> xmom_crse_rhs;
767  // amrex::Vector<amrex::MultiFab> ymom_crse_rhs;
768  amrex::Vector<amrex::MultiFab> zmom_crse_rhs;
769 
770  std::unique_ptr<Microphysics> micro;
771  amrex::Vector<amrex::Vector<amrex::MultiFab*>> qmoist; // (lev,ncomp) This has up to 8 components: qt, qv, qc, qi, qp, qr, qs, qg
772 
773  // Variables for wind farm parametrization models
774 
775 #ifdef ERF_USE_WINDFARM
776  std::unique_ptr<WindFarm> windfarm;
777  amrex::Vector<amrex::MultiFab> Nturb;
778  amrex::Vector<amrex::MultiFab> vars_windfarm; // Fitch: Vabs, Vabsdt, dudt, dvdt, dTKEdt
779  // EWP: dudt, dvdt, dTKEdt
780 
781  amrex::Vector<amrex::MultiFab> SMark; // A multifab that holds an integer corresponding
782  // to the number of the wind turbine to sample
783  // velocity upstream of the turbine
784 #endif
785 
787  amrex::Vector<amrex::Vector<amrex::MultiFab*>> lsm_data; // (lev,ncomp) Components: theta, q1, q2
788  amrex::Vector<amrex::Vector<amrex::MultiFab*>> lsm_flux; // (lev,ncomp) Components: theta, q1, q2
789 
790 #if defined(ERF_USE_RRTMGP)
791  amrex::Vector<std::unique_ptr<Radiation>> rad;
792  amrex::Vector<std::unique_ptr<amrex::MultiFab>> qheating_rates; // radiation heating rate source terms
793 
794  // Containers for additional SLM inputs
795  amrex::Vector<std::unique_ptr<amrex::MultiFab>> sw_lw_fluxes; // Direct SW (visible, NIR), Diffuse SW (visible, NIR), LW flux
796  amrex::Vector<std::unique_ptr<amrex::MultiFab>> solar_zenith; // Solar zenith angle
797 
798  bool plot_rad = false;
799 #endif
800 
801  // Fillpatcher classes for coarse-fine boundaries
802  int cf_width{0};
803  int cf_set_width{0};
804  amrex::Vector<ERFFillPatcher> FPr_c;
805  amrex::Vector<ERFFillPatcher> FPr_u;
806  amrex::Vector<ERFFillPatcher> FPr_v;
807  amrex::Vector<ERFFillPatcher> FPr_w;
808 
809  // Diffusive stresses and Smag
810  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Tau11_lev, Tau22_lev, Tau33_lev;
811  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Tau12_lev, Tau21_lev;
812  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Tau13_lev, Tau31_lev;
813  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Tau23_lev, Tau32_lev;
814  amrex::Vector<std::unique_ptr<amrex::MultiFab>> eddyDiffs_lev;
815  amrex::Vector<std::unique_ptr<amrex::MultiFab>> SmnSmn_lev;
816 
817  // Sea Surface Temps and Land Masks (lev, ntimes)
818  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab>>> sst_lev;
819  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::iMultiFab>>> lmask_lev;
820 
821  // Other SFS terms
822  amrex::Vector<std::unique_ptr<amrex::MultiFab>> SFS_hfx1_lev, SFS_hfx2_lev, SFS_hfx3_lev;
823  amrex::Vector<std::unique_ptr<amrex::MultiFab>> SFS_diss_lev;
824  amrex::Vector<std::unique_ptr<amrex::MultiFab>> SFS_q1fx1_lev, SFS_q1fx2_lev, SFS_q1fx3_lev;
825  amrex::Vector<std::unique_ptr<amrex::MultiFab>> SFS_q2fx3_lev;
826 
827  // Grid stretching
828  amrex::Vector<amrex::Vector<amrex::Real>> zlevels_stag; // nominal height levels
829 
830  // Data structures for terrain-fitted coordinates
831  amrex::Vector<std::unique_ptr<amrex::MultiFab>> z_phys_nd;
832  amrex::Vector<std::unique_ptr<amrex::MultiFab>> z_phys_cc;
833 
834  amrex::Vector<std::unique_ptr<amrex::MultiFab>> detJ_cc;
835  amrex::Vector<std::unique_ptr<amrex::MultiFab>> ax;
836  amrex::Vector<std::unique_ptr<amrex::MultiFab>> ay;
837  amrex::Vector<std::unique_ptr<amrex::MultiFab>> az;
838 
839  amrex::Vector<std::unique_ptr<amrex::MultiFab>> z_phys_nd_src;
840  amrex::Vector<std::unique_ptr<amrex::MultiFab>> detJ_cc_src;
841  amrex::Vector<std::unique_ptr<amrex::MultiFab>> ax_src;
842  amrex::Vector<std::unique_ptr<amrex::MultiFab>> ay_src;
843  amrex::Vector<std::unique_ptr<amrex::MultiFab>> az_src;
844 
845  amrex::Vector<std::unique_ptr<amrex::MultiFab>> z_phys_nd_new;
846  amrex::Vector<std::unique_ptr<amrex::MultiFab>> detJ_cc_new;
847  amrex::Vector<std::unique_ptr<amrex::MultiFab>> ax_new;
848  amrex::Vector<std::unique_ptr<amrex::MultiFab>> ay_new;
849  amrex::Vector<std::unique_ptr<amrex::MultiFab>> az_new;
850 
851  amrex::Vector<std::unique_ptr<amrex::MultiFab>> z_t_rk;
852 
853  // Data structures for terrain immersed forcing
854  amrex::Vector<std::unique_ptr<amrex::MultiFab>> terrain_blanking;
855 
856  // Wall distance function
857  amrex::Vector<std::unique_ptr<amrex::MultiFab>> walldist;
858 
859  // Map scale factors
860  amrex::Vector<std::unique_ptr<amrex::MultiFab>> mapfac_m;
861  amrex::Vector<std::unique_ptr<amrex::MultiFab>> mapfac_u;
862  amrex::Vector<std::unique_ptr<amrex::MultiFab>> mapfac_v;
863 
864  amrex::Vector<amrex::Vector<amrex::Real>> stretched_dz_h;
865  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real>> stretched_dz_d;
866 
867  amrex::Vector<amrex::MultiFab> base_state;
868  amrex::Vector<amrex::MultiFab> base_state_new;
869 
870  // Wave coupling data
871  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Hwave;
872  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Lwave;
873  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Hwave_onegrid;
874  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Lwave_onegrid;
875  bool finished_wave = false;
876 
877  // array of flux registers
878  amrex::Vector<amrex::YAFluxRegister*> advflux_reg;
879 
880  // A BCRec is essentially a 2*DIM integer array storing the boundary
881  // condition type at each lo/hi walls in each direction. We have one BCRec
882  // for each component of the cell-centered variables and each velocity component.
883  amrex::Vector <amrex::BCRec> domain_bcs_type;
884  amrex::Gpu::DeviceVector<amrex::BCRec> domain_bcs_type_d;
885 
886  // We store these so that we can print them out in the job_info file
887  amrex::Array<std::string,2*AMREX_SPACEDIM> domain_bc_type;
888 
889  // These hold the Dirichlet values at walls which need them ...
890  amrex::Array<amrex::Array<amrex::Real, AMREX_SPACEDIM*2>, AMREX_SPACEDIM+NBCVAR_max> m_bc_extdir_vals;
891 
892  // These hold the Neumann values at walls which need them ...
893  amrex::Array<amrex::Array<amrex::Real, AMREX_SPACEDIM*2>, AMREX_SPACEDIM+NBCVAR_max> m_bc_neumann_vals;
894 
895  // These are the "physical" boundary condition types (e.g. "inflow")
896  amrex::GpuArray<ERF_BC, AMREX_SPACEDIM*2> phys_bc_type;
897 
898  // These are the masks for the thin immersed body
899  amrex::Vector<std::unique_ptr<amrex::iMultiFab>> xflux_imask;
900  amrex::Vector<std::unique_ptr<amrex::iMultiFab>> yflux_imask;
901  amrex::Vector<std::unique_ptr<amrex::iMultiFab>> zflux_imask;
902  //amrex::Vector<std::unique_ptr<amrex::iMultiFab>> overset_imask;
903 
904  // These are the body forces that result from the thin immersed body
905  amrex::Vector<std::unique_ptr<amrex::MultiFab>> thin_xforce;
906  amrex::Vector<std::unique_ptr<amrex::MultiFab>> thin_yforce;
907  amrex::Vector<std::unique_ptr<amrex::MultiFab>> thin_zforce;
908 
911 
914 
915  ////////////////
916  // runtime parameters
917 
918  // maximum number of steps and stop time
919  int max_step = std::numeric_limits<int>::max();
920  amrex::Real start_time = 0.0;
921  amrex::Real stop_time = std::numeric_limits<amrex::Real>::max();
922 
923  // if >= 0 we restart from a checkpoint
924  std::string restart_chkfile = "";
925 
926  // Time step controls
927  static amrex::Real cfl;
928  static amrex::Real sub_cfl;
929  static amrex::Real init_shrink;
930  static amrex::Real change_max;
931  static amrex::Real dt_max_initial;
932  static amrex::Real dt_max;
933 
934  // Fixed dt for level 0 timesteps (only used if positive)
935  amrex::Vector<amrex::Real> fixed_dt;
936  amrex::Vector<amrex::Real> fixed_fast_dt;
937  static int fixed_mri_dt_ratio;
938 
939  // how often each level regrids the higher levels of refinement
940  // (after a level advances that many time steps)
941  int regrid_int = -1;
942 
943  // plotfile prefix and frequency
944  std::string plot_file_1 {"plt_1_"};
945  std::string plot_file_2 {"plt_2_"};
947  int m_plot_int_1 = -1;
948  int m_plot_int_2 = -1;
949  amrex::Real m_plot_per_1 = -1.0;
950  amrex::Real m_plot_per_2 = -1.0;
951  bool m_plot_face_vels = false;
952 
953  bool plot_lsm = false;
954 
955  // other sampling output control
956  int profile_int = -1;
957  bool destag_profiles = true;
958 
959  // Checkpoint type, prefix and frequency
960  std::string check_file {"chk"};
961  std::string check_type {"native"};
962  std::string restart_type {"native"};
963  int m_check_int = -1;
964  amrex::Real m_check_per = -1.0;
965 
966  amrex::Vector<std::string> plot_var_names_1;
967  amrex::Vector<std::string> plot_var_names_2;
968  const amrex::Vector<std::string> cons_names {"density", "rhotheta", "rhoKE", "rhoadv_0",
969  "rhoQ1", "rhoQ2", "rhoQ3",
970  "rhoQ4", "rhoQ5", "rhoQ6"};
971 
972  // Note that the order of variable names here must match the order in ERF_Derive.cpp
973  const amrex::Vector<std::string> derived_names {"soundspeed", "temp", "theta", "KE", "scalar",
974  "vorticity_x","vorticity_y","vorticity_z",
975  "magvel", "divU",
976  "pres_hse", "dens_hse", "theta_hse", "pressure", "pert_pres", "pert_dens",
977  "eq_pot_temp", "num_turb", "SMark0", "SMark1",
978  "dpdx", "dpdy", "pres_hse_x", "pres_hse_y",
979  "z_phys", "detJ" , "mapfac", "lat_m", "lon_m",
980  // Time averaged velocity
981  "u_t_avg", "v_t_avg", "w_t_avg", "umag_t_avg",
982  // eddy viscosity
983  "nut",
984  // eddy diffusivity of momentum
985  "Kmv","Kmh",
986  // eddy diffusivity of heat
987  "Khv","Khh",
988  // turbulence lengthscale
989  "Lturb",
990  // wall distance
991  "walldist",
992  // dissipation
993  "diss",
994  // moisture vars
995  "qt", "qv", "qc", "qi", "qp", "qrain", "qsnow", "qgraup", "qsat",
996  "rain_accum", "snow_accum", "graup_accum",
997  // Terrain IB mask
998  "terrain_IB_mask",
999  // EB variables
1000  "volfrac"
1001 #ifdef ERF_COMPUTE_ERROR
1002  // error vars
1003  ,"xvel_err", "yvel_err", "zvel_err", "pp_err"
1004 #endif
1005 #ifdef ERF_USE_RRTMGP
1006  ,"qsrc_sw", "qsrc_lw"
1007 #endif
1008  };
1009 
1010  // algorithm choices
1012 
1013  // Turbulent perturbation structure
1015 
1016 #ifdef ERF_USE_PARTICLES
1017  // Particle container with all particle species
1018  ParticleData particleData;
1019 
1020  // variables and functions for tracers particles
1021  bool m_use_tracer_particles; /*!< tracer particles that advect with flow */
1022  bool m_use_hydro_particles; /*!< tracer particles that fall with gravity */
1023 
1024  /*! Read tracer and hydro particles parameters */
1025  void readTracersParams();
1026 
1027  /*! Initialize tracer and hydro particles */
1028  void initializeTracers ( amrex::ParGDBBase*,
1029  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& );
1030 
1031  /*! Restart tracer and hydro particles */
1032  void restartTracers ( amrex::ParGDBBase*, const std::string& );
1033 
1034  /*! Evolve tracers and hydro particles */
1035  void evolveTracers( int,
1036  amrex::Real,
1037  amrex::Vector<amrex::Vector<amrex::MultiFab>>&,
1038  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& );
1039 
1040 #endif
1041 
1042 #ifdef ERF_USE_MULTIBLOCK
1043  MultiBlockContainer *m_mbc = nullptr;
1044 #endif
1045 
1046  static int verbose;
1047  static int mg_verbose;
1048  static bool use_fft;
1049 
1050  // Diagnostic output interval
1051  static int sum_interval;
1052  static int pert_interval;
1053  static amrex::Real sum_per;
1054 
1055  // Write in native AMReX or NetCDF format for each plotfile
1056  static PlotFileType plotfile_type_1;
1057  static PlotFileType plotfile_type_2;
1058 
1059  static InitType init_type;
1060  static StateInterpType interpolation_type;
1061 
1062  // sponge_type: "input_sponge"
1063  static std::string sponge_type;
1064 
1065  // use_real_bcs: only true if 1) ( (init_type == InitType::Real) or (init_type == InitType::Metgrid) )
1066  // AND 2) we want to use the bc's from the WRF bdy file
1067  static bool use_real_bcs;
1068 
1069  // NetCDF initialization (wrfinput) file
1070  static amrex::Vector<amrex::Vector<std::string>> nc_init_file;
1071 
1072  // NetCDF initialization (wrfbdy/met_em) file
1073  static std::string nc_bdy_file;
1074  int real_width{0};
1076 
1077  // Flag to trigger initialization from input_sounding like WRF's ideal.exe
1078  // used with init_type == InitType::Input_Sounding
1079  static bool init_sounding_ideal;
1080 
1081  // Options for vertical interpolation of met_em*.nc data.
1084  bool metgrid_debug_dry{false};
1085  bool metgrid_debug_psfc{false};
1086  bool metgrid_debug_msf{false};
1090  bool metgrid_use_sfc{true};
1091  bool metgrid_retain_sfc{false};
1092  amrex::Real metgrid_proximity{500.0};
1095 
1096  // 1D CDF output (for ingestion in AMR-Wind)
1097  static int output_1d_column;
1098  static int column_interval;
1099  static amrex::Real column_per;
1100  static amrex::Real column_loc_x;
1101  static amrex::Real column_loc_y;
1102  static std::string column_file_name;
1103 
1104  // 2D BndryRegister output (for ingestion in AMR-Wind)
1107  static amrex::Real bndry_output_planes_per;
1109 
1110  // 2D BndryRegister input
1112 
1113  static int ng_dens_hse;
1114  static int ng_pres_hse;
1115 
1116  // Custom source terms
1117  amrex::Vector< amrex::Vector<amrex::Real> > h_rhotheta_src;
1118  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real> > d_rhotheta_src;
1119 
1120  amrex::Vector< amrex::Vector<amrex::Real> > h_rhoqt_src;
1121  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real> > d_rhoqt_src;
1122 
1123  amrex::Vector< amrex::Vector<amrex::Real> > h_w_subsid;
1124  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real> > d_w_subsid;
1125 
1126  amrex::Vector< amrex::Vector<amrex::Real> > h_u_geos;
1127  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real> > d_u_geos;
1128 
1129  amrex::Vector< amrex::Vector<amrex::Real> > h_v_geos;
1130  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real> > d_v_geos;
1131 
1132  // Function to read and populate above host vectors (if input file exists)
1133  void init_geo_wind_profile (const std::string input_file,
1134  amrex::Vector<amrex::Real>& u_geos,
1135  amrex::Gpu::DeviceVector<amrex::Real>& u_geos_d,
1136  amrex::Vector<amrex::Real>& v_geos,
1137  amrex::Gpu::DeviceVector<amrex::Real>& v_geos_d,
1138  const amrex::Geometry& lgeom,
1139  const amrex::Vector<amrex::Real>& zlev_stag);
1140 
1141  // This is a vector over levels of vectors across quantities of Vectors
1142  amrex::Vector<amrex::Vector< amrex::Vector<amrex::Real> > > h_rayleigh_ptrs;
1143 
1144  // This is a vector over levels of vectors across quantities of Vectors
1145  amrex::Vector<amrex::Vector< amrex::Vector<amrex::Real> > > h_sponge_ptrs;
1146 
1147  // This is a vector over levels of vectors across quantities of DeviceVectors
1148  amrex::Vector<amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real> > > d_rayleigh_ptrs;
1149 
1150  // This is a vector over levels of vectors across quantities of DeviceVectors
1151  amrex::Vector<amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real> > > d_sponge_ptrs;
1152 
1153  amrex::Vector<amrex::Real> h_havg_density;
1154  amrex::Vector<amrex::Real> h_havg_temperature;
1155  amrex::Vector<amrex::Real> h_havg_pressure;
1156  amrex::Vector<amrex::Real> h_havg_qv;
1157  amrex::Vector<amrex::Real> h_havg_qc;
1158 
1159  amrex::Gpu::DeviceVector<amrex::Real> d_havg_density;
1160  amrex::Gpu::DeviceVector<amrex::Real> d_havg_temperature;
1161  amrex::Gpu::DeviceVector<amrex::Real> d_havg_pressure;
1162  amrex::Gpu::DeviceVector<amrex::Real> d_havg_qv;
1163  amrex::Gpu::DeviceVector<amrex::Real> d_havg_qc;
1164 
1165  void refinement_criteria_setup ();
1166 
1167  std::unique_ptr<WriteBndryPlanes> m_w2d = nullptr;
1168  std::unique_ptr<ReadBndryPlanes> m_r2d = nullptr;
1169  std::unique_ptr<ABLMost> m_most = nullptr;
1170  amrex::Vector<std::unique_ptr<ForestDrag>> m_forest_drag;
1171 
1172  //
1173  // Holds info for dynamically generated tagging criteria
1174  //
1175  static amrex::Vector<amrex::AMRErrorTag> ref_tags;
1176 
1177  //
1178  // Build a mask that zeroes out values on a coarse level underlying
1179  // grids on the next finest level
1180  //
1181  amrex::MultiFab fine_mask;
1182 
1183  amrex::Vector<amrex::Real> dz_min;
1184 
1185  static AMREX_FORCE_INLINE
1186  int
1187  ComputeGhostCells (const AdvChoice& advChoice, bool use_num_diff)
1188  {
1189  if (use_num_diff)
1190  {
1191  return 3;
1192  } else {
1193  if (
1199  || (advChoice.dryscal_vert_adv_type == AdvType::Centered_6th) )
1200  { return 3; }
1201  else if (
1203  || (advChoice.dycore_vert_adv_type == AdvType::Upwind_5th)
1207  || (advChoice.dryscal_vert_adv_type == AdvType::Upwind_5th) )
1208  { return 3; }
1209  else if (
1211  || (advChoice.dryscal_vert_adv_type == AdvType::Weno_5)
1212  || (advChoice.moistscal_horiz_adv_type == AdvType::Weno_5)
1213  || (advChoice.moistscal_vert_adv_type == AdvType::Weno_5)
1214  || (advChoice.moistscal_horiz_adv_type == AdvType::Weno_5Z)
1215  || (advChoice.moistscal_vert_adv_type == AdvType::Weno_5Z)
1216  || (advChoice.dryscal_horiz_adv_type == AdvType::Weno_5Z)
1217  || (advChoice.dryscal_vert_adv_type == AdvType::Weno_5Z) )
1218  { return 3; }
1219  else if (
1221  || (advChoice.dryscal_vert_adv_type == AdvType::Weno_7)
1222  || (advChoice.moistscal_horiz_adv_type == AdvType::Weno_7)
1223  || (advChoice.moistscal_vert_adv_type == AdvType::Weno_7)
1224  || (advChoice.moistscal_horiz_adv_type == AdvType::Weno_7Z)
1225  || (advChoice.moistscal_vert_adv_type == AdvType::Weno_7Z)
1226  || (advChoice.dryscal_horiz_adv_type == AdvType::Weno_7Z)
1227  || (advChoice.dryscal_vert_adv_type == AdvType::Weno_7Z) )
1228  { return 4; }
1229  else
1230  { return 2; }
1231  }
1232  }
1233 
1234  AMREX_FORCE_INLINE
1235  amrex::YAFluxRegister* getAdvFluxReg (int lev)
1236  {
1237  return advflux_reg[lev];
1238  }
1239 
1240  AMREX_FORCE_INLINE
1241  std::ostream&
1242  DataLog (int i)
1243  {
1244  return *datalog[i];
1245  }
1246 
1247  AMREX_FORCE_INLINE
1248  int
1249  NumDataLogs () noexcept
1250  {
1251  return datalog.size();
1252  }
1253 
1254  AMREX_FORCE_INLINE
1255  std::ostream&
1257  {
1258  return *sampleptlog[i];
1259  }
1260 
1261  AMREX_FORCE_INLINE
1262  int
1264  {
1265  return sampleptlog.size();
1266  }
1267 
1268  AMREX_FORCE_INLINE
1269  std::ostream&
1271  {
1272  return *samplelinelog[i];
1273  }
1274 
1275  AMREX_FORCE_INLINE
1276  int
1277  NumSampleLineLogs () noexcept
1278  {
1279  return samplelinelog.size();
1280  }
1281 
1282  amrex::IntVect&
1283  SamplePoint (int i)
1284  {
1285  return samplepoint[i];
1286  }
1287 
1288  AMREX_FORCE_INLINE
1289  int
1290  NumSamplePoints () noexcept
1291  {
1292  return samplepoint.size();
1293  }
1294 
1295  amrex::IntVect&
1296  SampleLine (int i)
1297  {
1298  return sampleline[i];
1299  }
1300 
1301  AMREX_FORCE_INLINE
1302  int
1303  NumSampleLines () noexcept
1304  {
1305  return sampleline.size();
1306  }
1307 
1308  static amrex::Real startCPUTime;
1309  static amrex::Real previousCPUTimeUsed;
1310 
1311  static amrex::Real
1313  {
1314  int numCores = amrex::ParallelDescriptor::NProcs();
1315 #ifdef _OPENMP
1316  numCores = numCores * omp_get_max_threads();
1317 #endif
1318 
1319  amrex::Real T =
1320  numCores * (amrex::ParallelDescriptor::second() - startCPUTime) +
1322 
1323  return T;
1324  }
1325 
1326  void setRecordDataInfo (int i, const std::string& filename) // NOLINT
1327  {
1328  if (amrex::ParallelDescriptor::IOProcessor())
1329  {
1330  datalog[i] = std::make_unique<std::fstream>();
1331  datalog[i]->open(filename.c_str(),std::ios::out|std::ios::app);
1332  if (!datalog[i]->good()) {
1333  amrex::FileOpenFailed(filename);
1334  }
1335  }
1336  amrex::ParallelDescriptor::Barrier("ERF::setRecordDataInfo");
1337  }
1338 
1339  void setRecordSamplePointInfo (int i, int lev, amrex::IntVect& cell, const std::string& filename) // NOLINT
1340  {
1341  amrex::MultiFab dummy(grids[lev],dmap[lev],1,0);
1342  for (amrex::MFIter mfi(dummy); mfi.isValid(); ++mfi)
1343  {
1344  const amrex::Box& bx = mfi.validbox();
1345  if (bx.contains(cell)) {
1346  sampleptlog[i] = std::make_unique<std::fstream>();
1347  sampleptlog[i]->open(filename.c_str(),std::ios::out|std::ios::app);
1348  if (!sampleptlog[i]->good()) {
1349  amrex::FileOpenFailed(filename);
1350  }
1351  }
1352  }
1353  amrex::ParallelDescriptor::Barrier("ERF::setRecordSamplePointInfo");
1354  }
1355 
1356  void setRecordSampleLineInfo (int i, int lev, amrex::IntVect& cell, const std::string& filename) // NOLINT
1357  {
1358  amrex::MultiFab dummy(grids[lev],dmap[lev],1,0);
1359  for (amrex::MFIter mfi(dummy); mfi.isValid(); ++mfi)
1360  {
1361  const amrex::Box& bx = mfi.validbox();
1362  if (bx.contains(cell)) {
1363  samplelinelog[i] = std::make_unique<std::fstream>();
1364  samplelinelog[i]->open(filename.c_str(),std::ios::out|std::ios::app);
1365  if (!samplelinelog[i]->good()) {
1366  amrex::FileOpenFailed(filename);
1367  }
1368  }
1369  }
1370  amrex::ParallelDescriptor::Barrier("ERF::setRecordSampleLineInfo");
1371  }
1372 
1373  // Data sampler for line and plane output
1375  amrex::Real sampler_per = -1.0;
1376  std::unique_ptr<SampleData> data_sampler = nullptr;
1377 
1378  amrex::Vector<std::unique_ptr<std::fstream> > datalog;
1379  amrex::Vector<std::string> datalogname;
1380 
1381  amrex::Vector<std::unique_ptr<std::fstream> > sampleptlog;
1382  amrex::Vector<std::string> sampleptlogname;
1383  amrex::Vector<amrex::IntVect> samplepoint;
1384 
1385  amrex::Vector<std::unique_ptr<std::fstream> > samplelinelog;
1386  amrex::Vector<std::string> samplelinelogname;
1387  amrex::Vector<amrex::IntVect> sampleline;
1388 
1389  //! The filename of the ith datalog file.
1390  [[nodiscard]] std::string DataLogName (int i) const noexcept { return datalogname[i]; }
1391 
1392  //! The filename of the ith sampleptlog file.
1393  [[nodiscard]] std::string SamplePointLogName (int i) const noexcept { return sampleptlogname[i]; }
1394 
1395  //! The filename of the ith samplelinelog file.
1396  [[nodiscard]] std::string SampleLineLogName (int i) const noexcept { return samplelinelogname[i]; }
1397 
1398  // amrex::Vector<std::unique_ptr<amrex::EBFArrayBoxFactory> > m_factory;
1399  amrex::Vector<std::unique_ptr<amrex::FabFactory<amrex::FArrayBox> > > m_factory;
1400 
1401  [[nodiscard]] amrex::FabFactory<amrex::FArrayBox> const&
1402  Factory (int lev) const noexcept { return *m_factory[lev]; }
1403  [[nodiscard]] amrex::EBFArrayBoxFactory const&
1404  EBFactory (int lev) const noexcept {
1405  return static_cast<amrex::EBFArrayBoxFactory const&>(*m_factory[lev]);
1406  }
1407 
1408  [[nodiscard]] static int nghost_eb_basic ()
1409  { return 5; }
1410 
1411  // We need 5 for doing StateRedistribution; otherwise 4 would be enough
1412  [[nodiscard]] static int nghost_eb_volume ()
1413  { return 5; }
1414 
1415  [[nodiscard]] static int nghost_eb_full ()
1416  { return 4; }
1417 
1418 #ifdef ERF_USE_FFT
1419  amrex::Vector<std::unique_ptr<amrex::FFT::Poisson<amrex::MultiFab>>> m_3D_poisson;
1420  amrex::Vector<std::unique_ptr<amrex::FFT::PoissonHybrid<amrex::MultiFab>>> m_2D_poisson;
1421 #endif
1422 
1423 public:
1424  void writeJobInfo (const std::string& dir) const;
1425  static void writeBuildInfo (std::ostream& os);
1426 
1427  static void print_banner(MPI_Comm /*comm*/, std::ostream& /*out*/);
1428  static void print_usage(MPI_Comm /*comm*/, std::ostream& /*out*/);
1429  static void print_error(MPI_Comm /*comm*/, const std::string& msg);
1430  static void print_summary(std::ostream&);
1431  static void print_tpls(std::ostream& /*out*/);
1432 };
1433 
1434 #endif
AMREX_ENUM(InitType, None, Input_Sounding, Ideal, Real, Metgrid, Uniform)
Contains the Eulerian microphysics class.
#define NBCVAR_max
Definition: ERF_IndexDefines.H:29
@ Centered_6th
Contains the Lagrangian microphysics class.
Definition: ERF.H:126
void MakeHorizontalAverages()
Definition: ERF.cpp:1747
amrex::Vector< amrex::MultiFab > rU_new
Definition: ERF.H:760
int last_check_file_step
Definition: ERF.H:912
amrex::Vector< std::unique_ptr< amrex::MultiFab > > walldist
Definition: ERF.H:857
void initRayleigh()
Initialize Rayleigh damping profiles.
Definition: ERF_InitRayleigh.cpp:14
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Tau23_lev
Definition: ERF.H:813
bool metgrid_basic_linear
Definition: ERF.H:1088
amrex::Vector< amrex::Vector< amrex::Box > > boxes_at_level
Definition: ERF.H:721
std::unique_ptr< ABLMost > m_most
Definition: ERF.H:1169
amrex::Vector< std::string > samplelinelogname
Definition: ERF.H:1386
void MakeEBGeometry()
int max_step
Definition: ERF.H:919
bool metgrid_debug_msf
Definition: ERF.H:1086
void setRayleighRefFromSounding(bool restarting)
Set Rayleigh mean profiles from input sounding.
Definition: ERF_InitRayleigh.cpp:55
int last_plot_file_step_2
Definition: ERF.H:910
amrex::Vector< std::unique_ptr< MRISplitIntegrator< amrex::Vector< amrex::MultiFab > > > > mri_integrator_mem
Definition: ERF.H:742
void Evolve()
Definition: ERF.cpp:390
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > d_rhotheta_src
Definition: ERF.H:1118
amrex::Vector< amrex::MultiFab > pp_inc
Definition: ERF.H:744
void make_eb_box()
amrex::Vector< ERFFillPatcher > FPr_u
Definition: ERF.H:805
amrex::Real cloud_fraction(amrex::Real time)
Definition: ERF_WriteScalarProfiles.cpp:173
amrex::Vector< amrex::IntVect > sampleline
Definition: ERF.H:1387
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SFS_q1fx3_lev
Definition: ERF.H:824
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Hwave_onegrid
Definition: ERF.H:873
amrex::Vector< amrex::Vector< amrex::MultiFab > > vars_new
Definition: ERF.H:735
ERF(ERF &&) noexcept=delete
AMREX_FORCE_INLINE int NumSampleLineLogs() noexcept
Definition: ERF.H:1277
static void print_tpls(std::ostream &)
Definition: ERF_ConsoleIO.cpp:137
amrex::Vector< amrex::Real > dz_min
Definition: ERF.H:1183
std::string DataLogName(int i) const noexcept
The filename of the ith datalog file.
Definition: ERF.H:1390
static PlotFileType plotfile_type_1
Definition: ERF.H:1056
amrex::Vector< std::unique_ptr< amrex::MultiFab > > thin_yforce
Definition: ERF.H:906
void ErrorEst(int lev, amrex::TagBoxArray &tags, amrex::Real time, int ngrow) override
Definition: ERF_Tagging.cpp:15
std::string plot_file_2
Definition: ERF.H:945
void FillBdyCCVels(amrex::Vector< amrex::MultiFab > &mf_cc_vel)
Definition: ERF_FillBdyCCVels.cpp:11
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SFS_hfx3_lev
Definition: ERF.H:822
amrex::Vector< amrex::Vector< amrex::Real > > h_w_subsid
Definition: ERF.H:1123
static void print_banner(MPI_Comm, std::ostream &)
Definition: ERF_ConsoleIO.cpp:60
amrex::Vector< std::unique_ptr< std::fstream > > datalog
Definition: ERF.H:1378
static amrex::Real sum_per
Definition: ERF.H:1053
amrex::Vector< ERFFillPatcher > FPr_v
Definition: ERF.H:806
int cf_set_width
Definition: ERF.H:803
amrex::Vector< std::unique_ptr< amrex::MultiFab > > mapfac_u
Definition: ERF.H:861
static amrex::Real previousCPUTimeUsed
Definition: ERF.H:1309
void setPlotVariables(const std::string &pp_plot_var_names, amrex::Vector< std::string > &plot_var_names)
Definition: ERF_Plotfile.cpp:17
amrex::Gpu::DeviceVector< amrex::Real > d_havg_temperature
Definition: ERF.H:1160
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Tau31_lev
Definition: ERF.H:812
amrex::Vector< amrex::Vector< amrex::Vector< amrex::Real > > > h_sponge_ptrs
Definition: ERF.H:1145
void fill_from_bndryregs(const amrex::Vector< amrex::MultiFab * > &mfs, amrex::Real time)
Definition: ERF_BoundaryConditionsBndryReg.cpp:13
amrex::Vector< std::unique_ptr< amrex::MultiFab > > mapfac_m
Definition: ERF.H:860
void setRecordDataInfo(int i, const std::string &filename)
Definition: ERF.H:1326
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SFS_hfx1_lev
Definition: ERF.H:822
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_extdir_vals
Definition: ERF.H:890
void project_velocities_tb(int lev, amrex::Real dt, amrex::Vector< amrex::MultiFab > &vars, amrex::MultiFab &p)
Definition: ERF_PoissonSolve_tb.cpp:20
void init_from_input_sounding(int lev)
Definition: ERF_InitFromInputSounding.cpp:50
amrex::Gpu::DeviceVector< amrex::Real > d_havg_qv
Definition: ERF.H:1162
static amrex::Real column_loc_y
Definition: ERF.H:1101
static int mg_verbose
Definition: ERF.H:1047
void ReadParameters()
Definition: ERF.cpp:1398
static amrex::Vector< std::string > PlotFileVarNames(amrex::Vector< std::string > plot_var_names)
Definition: ERF_Plotfile.cpp:174
void InitializeFromFile()
amrex::Vector< std::unique_ptr< ERFPhysBCFunct_cons > > physbcs_cons
Definition: ERF.H:747
static bool init_sounding_ideal
Definition: ERF.H:1079
amrex::Vector< std::unique_ptr< amrex::MultiFab > > ay_new
Definition: ERF.H:848
ERF()
Definition: ERF.cpp:89
void init_Dirichlet_bc_data(const std::string input_file)
Definition: ERF_InitBCs.cpp:600
~ERF() override
std::string restart_type
Definition: ERF.H:962
amrex::Vector< std::unique_ptr< amrex::MultiFab > > z_phys_nd_src
Definition: ERF.H:839
amrex::Vector< std::unique_ptr< amrex::MultiFab > > detJ_cc
Definition: ERF.H:834
amrex::Real m_plot_per_1
Definition: ERF.H:949
void WriteMyEBSurface()
Definition: ERF_EBWriteSurface.cpp:5
amrex::Vector< std::unique_ptr< amrex::MultiFab > > z_phys_cc
Definition: ERF.H:832
amrex::Vector< std::unique_ptr< amrex::MultiFab > > eddyDiffs_lev
Definition: ERF.H:814
static SolverChoice solverChoice
Definition: ERF.H:1011
std::string check_type
Definition: ERF.H:961
amrex::Vector< ERFFillPatcher > FPr_c
Definition: ERF.H:804
static bool use_fft
Definition: ERF.H:1048
bool m_plot_face_vels
Definition: ERF.H:951
amrex::Vector< amrex::MultiFab > base_state_new
Definition: ERF.H:868
amrex::Vector< std::unique_ptr< amrex::MultiFab > > az
Definition: ERF.H:837
int regrid_int
Definition: ERF.H:941
amrex::Vector< amrex::Real > fixed_dt
Definition: ERF.H:935
void write_1D_profiles_stag(amrex::Real time)
Definition: ERF_Write1DProfiles_stag.cpp:23
void WriteGenericPlotfileHeaderWithTerrain(std::ostream &HeaderFile, int nlevels, const amrex::Vector< amrex::BoxArray > &bArray, const amrex::Vector< std::string > &varnames, const amrex::Vector< amrex::Geometry > &my_geom, amrex::Real time, const amrex::Vector< int > &level_steps, const amrex::Vector< amrex::IntVect > &my_ref_ratio, const std::string &versionName, const std::string &levelPrefix, const std::string &mfPrefix) const
Definition: ERF_Plotfile.cpp:1691
amrex::Vector< amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > > d_sponge_ptrs
Definition: ERF.H:1151
amrex::Vector< amrex::Vector< amrex::Real > > h_rhoqt_src
Definition: ERF.H:1120
amrex::Vector< long > dt_mri_ratio
Definition: ERF.H:730
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vel_t_avg
Definition: ERF.H:739
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SFS_q2fx3_lev
Definition: ERF.H:825
static amrex::Real dt_max
Definition: ERF.H:932
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::iMultiFab > > > lmask_lev
Definition: ERF.H:819
amrex::Real stop_time
Definition: ERF.H:921
AMREX_FORCE_INLINE int NumSamplePointLogs() noexcept
Definition: ERF.H:1263
amrex::Vector< amrex::Real > h_havg_pressure
Definition: ERF.H:1155
void update_diffusive_arrays(int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Definition: ERF_MakeNewArrays.cpp:380
static int verbose
Definition: ERF.H:1046
amrex::Vector< amrex::Real > h_havg_qc
Definition: ERF.H:1157
amrex::Vector< std::unique_ptr< ERFPhysBCFunct_w > > physbcs_w
Definition: ERF.H:750
static int nghost_eb_basic()
Definition: ERF.H:1408
void Advance(int lev, amrex::Real time, amrex::Real dt_lev, int iteration, int ncycle)
Definition: ERF_Advance.cpp:20
static std::string column_file_name
Definition: ERF.H:1102
amrex::Vector< std::unique_ptr< std::fstream > > samplelinelog
Definition: ERF.H:1385
amrex::Vector< std::unique_ptr< amrex::MultiFab > > terrain_blanking
Definition: ERF.H:854
std::unique_ptr< Microphysics > micro
Definition: ERF.H:770
amrex::Vector< amrex::MultiFab > base_state
Definition: ERF.H:867
AMREX_FORCE_INLINE amrex::YAFluxRegister * getAdvFluxReg(int lev)
Definition: ERF.H:1235
amrex::Vector< amrex::Real > h_havg_density
Definition: ERF.H:1153
amrex::Vector< std::unique_ptr< amrex::MultiFab > > z_phys_nd_new
Definition: ERF.H:845
amrex::Vector< amrex::Real > fixed_fast_dt
Definition: ERF.H:936
bool metgrid_retain_sfc
Definition: ERF.H:1091
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Qv_prim
Definition: ERF.H:755
static int sum_interval
Definition: ERF.H:1051
static int pert_interval
Definition: ERF.H:1052
void restart()
Definition: ERF.cpp:1266
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Tau13_lev
Definition: ERF.H:812
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SFS_q1fx2_lev
Definition: ERF.H:824
void ReadCheckpointFileMOST()
Definition: ERF_Checkpoint.cpp:627
void FillIntermediatePatch(int lev, amrex::Real time, const amrex::Vector< amrex::MultiFab * > &mfs_vel, const amrex::Vector< amrex::MultiFab * > &mfs_mom, int ng_cons, int ng_vel, bool cons_only, int icomp_cons, int ncomp_cons, bool allow_most_bcs=true)
Definition: ERF_FillIntermediatePatch.cpp:28
amrex::IntVect & SampleLine(int i)
Definition: ERF.H:1296
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Tau21_lev
Definition: ERF.H:811
amrex::Vector< amrex::MultiFab > rV_new
Definition: ERF.H:762
std::string PlotFileName(int lev) const
void write_1D_profiles(amrex::Real time)
Definition: ERF_Write1DProfiles.cpp:15
void initialize_integrator(int lev, amrex::MultiFab &cons_mf, amrex::MultiFab &vel_mf)
Definition: ERF_MakeNewArrays.cpp:589
amrex::Vector< std::unique_ptr< amrex::MultiFab > > thin_zforce
Definition: ERF.H:907
static InitType init_type
Definition: ERF.H:1059
amrex::Vector< amrex::BCRec > domain_bcs_type
Definition: ERF.H:883
amrex::Array< amrex::LinOpBCType, AMREX_SPACEDIM > get_projection_bc(amrex::Orientation::Side side) const noexcept
Definition: ERF_SolveWithMLMG.cpp:17
int m_plot_int_1
Definition: ERF.H:947
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Qr_prim
Definition: ERF.H:756
void post_timestep(int nstep, amrex::Real time, amrex::Real dt_lev)
Definition: ERF.cpp:477
std::string pp_prefix
Definition: ERF.H:444
std::string SampleLineLogName(int i) const noexcept
The filename of the ith samplelinelog file.
Definition: ERF.H:1396
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > sst_lev
Definition: ERF.H:818
amrex::Vector< std::unique_ptr< amrex::MultiFab > > thin_xforce
Definition: ERF.H:905
AMREX_FORCE_INLINE int NumSamplePoints() noexcept
Definition: ERF.H:1290
bool metgrid_use_sfc
Definition: ERF.H:1090
void erf_enforce_hse(int lev, amrex::MultiFab &dens, amrex::MultiFab &pres, amrex::MultiFab &pi, amrex::MultiFab &th, std::unique_ptr< amrex::MultiFab > &z_cc)
Definition: ERF_Init1D.cpp:159
static amrex::Real bndry_output_planes_per
Definition: ERF.H:1107
amrex::Real m_check_per
Definition: ERF.H:964
void init_custom(int lev)
Definition: ERF_InitCustom.cpp:26
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Tau33_lev
Definition: ERF.H:810
std::unique_ptr< ProblemBase > prob
Definition: ERF.H:717
amrex::Vector< int > num_files_at_level
Definition: ERF.H:720
void WriteMultiLevelPlotfileWithTerrain(const std::string &plotfilename, int nlevels, const amrex::Vector< const amrex::MultiFab * > &mf, const amrex::Vector< const amrex::MultiFab * > &mf_nd, const amrex::Vector< std::string > &varnames, const amrex::Vector< amrex::Geometry > &my_geom, amrex::Real time, const amrex::Vector< int > &level_steps, const amrex::Vector< amrex::IntVect > &my_ref_ratio, const std::string &versionName="HyperCLaw-V1.1", const std::string &levelPrefix="Level_", const std::string &mfPrefix="Cell", const amrex::Vector< std::string > &extra_dirs=amrex::Vector< std::string >()) const
Definition: ERF_Plotfile.cpp:1604
void init_bcs()
Definition: ERF_InitBCs.cpp:20
int profile_int
Definition: ERF.H:956
bool metgrid_debug_quiescent
Definition: ERF.H:1082
amrex::Vector< std::unique_ptr< ERFPhysBCFunct_u > > physbcs_u
Definition: ERF.H:748
amrex::Vector< amrex::Real > t_new
Definition: ERF.H:727
bool destag_profiles
Definition: ERF.H:957
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > stretched_dz_d
Definition: ERF.H:865
amrex::Vector< amrex::Real > t_avg_cnt
Definition: ERF.H:740
amrex::Vector< std::string > plot_var_names_1
Definition: ERF.H:966
void remake_zphys(int lev, amrex::Real time, std::unique_ptr< amrex::MultiFab > &temp_zphys_nd)
Definition: ERF_MakeNewArrays.cpp:543
int m_check_int
Definition: ERF.H:963
amrex::Real estTimeStep(int lev, long &dt_fast_ratio) const
Definition: ERF_ComputeTimestep.cpp:54
static std::string sponge_type
Definition: ERF.H:1063
static amrex::Real startCPUTime
Definition: ERF.H:1308
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SFS_diss_lev
Definition: ERF.H:823
void redistribute_term(amrex::MFIter const &mfi, int lev, int ncomp, amrex::MultiFab &result, amrex::MultiFab &result_tmp, amrex::MultiFab const &state, amrex::BCRec const *bc, amrex::Real const dt)
amrex::Vector< amrex::MultiFab > rU_old
Definition: ERF.H:759
amrex::Vector< amrex::Real > t_old
Definition: ERF.H:728
AMREX_FORCE_INLINE int NumSampleLines() noexcept
Definition: ERF.H:1303
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Theta_prim
Definition: ERF.H:754
void redistribute_term(int lev, int ncomp, amrex::MultiFab &result, amrex::MultiFab &result_tmp, amrex::MultiFab const &state, amrex::BCRec const *bc, amrex::Real const dt)
void init1DArrays()
amrex::Vector< std::unique_ptr< amrex::MultiFab > > z_t_rk
Definition: ERF.H:851
void compute_divergence(int lev, amrex::MultiFab &rhs, amrex::Array< amrex::MultiFab const *, AMREX_SPACEDIM > rho0_u_const, amrex::Geometry const &geom_at_lev)
Definition: ERF_ComputeDivergence.cpp:10
void fill_rhs(amrex::MultiFab &rhs_mf, const amrex::MultiFab &state_mf, amrex::Real time, const amrex::Geometry &geom)
void appendPlotVariables(const std::string &pp_plot_var_names, amrex::Vector< std::string > &plot_var_names)
Definition: ERF_Plotfile.cpp:128
amrex::Vector< std::unique_ptr< ERFPhysBCFunct_v > > physbcs_v
Definition: ERF.H:749
amrex::Vector< std::string > plot_var_names_2
Definition: ERF.H:967
static amrex::Real column_per
Definition: ERF.H:1099
static int output_bndry_planes
Definition: ERF.H:1105
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > d_v_geos
Definition: ERF.H:1130
int last_plot_file_step_1
Definition: ERF.H:909
void update_terrain_arrays(int lev)
Definition: ERF_MakeNewArrays.cpp:578
static std::string nc_bdy_file
Definition: ERF.H:1073
bool metgrid_interp_theta
Definition: ERF.H:1087
static amrex::Vector< amrex::Vector< std::string > > nc_init_file
Definition: ERF.H:1070
void init_only(int lev, amrex::Real time)
Definition: ERF.cpp:1286
static int input_bndry_planes
Definition: ERF.H:1111
static StateInterpType interpolation_type
Definition: ERF.H:1060
amrex::Gpu::DeviceVector< amrex::Real > d_havg_qc
Definition: ERF.H:1163
void AverageDown()
Definition: ERF_AverageDown.cpp:16
amrex::Vector< amrex::Vector< amrex::Real > > h_v_geos
Definition: ERF.H:1129
bool projection_has_dirichlet(amrex::Array< amrex::LinOpBCType, AMREX_SPACEDIM > bcs) const
Definition: ERF_PoissonSolve_tb.cpp:8
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Lwave_onegrid
Definition: ERF.H:874
InputSoundingData input_sounding_data
Definition: ERF.H:686
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > d_rhoqt_src
Definition: ERF.H:1121
std::unique_ptr< SampleData > data_sampler
Definition: ERF.H:1376
amrex::MultiFab fine_mask
Definition: ERF.H:1181
amrex::Gpu::DeviceVector< amrex::Real > d_havg_density
Definition: ERF.H:1159
void init_from_hse(int lev)
Definition: ERF_InitFromHSE.cpp:32
static bool use_real_bcs
Definition: ERF.H:1067
amrex::Vector< amrex::YAFluxRegister * > advflux_reg
Definition: ERF.H:878
int metgrid_force_sfc_k
Definition: ERF.H:1094
amrex::Vector< amrex::Vector< amrex::Real > > h_rhotheta_src
Definition: ERF.H:1117
static int ng_pres_hse
Definition: ERF.H:1114
static amrex::Real bndry_output_planes_start_time
Definition: ERF.H:1108
static amrex::Real cfl
Definition: ERF.H:927
void solve_with_mlmg(int lev, amrex::Vector< amrex::MultiFab > &rhs, amrex::Vector< amrex::MultiFab > &p, amrex::Vector< amrex::Array< amrex::MultiFab, AMREX_SPACEDIM >> &fluxes)
Definition: ERF_SolveWithMLMG.cpp:40
amrex::Gpu::DeviceVector< amrex::BCRec > domain_bcs_type_d
Definition: ERF.H:884
amrex::Vector< std::unique_ptr< ForestDrag > > m_forest_drag
Definition: ERF.H:1170
void derive_diag_profiles_stag(amrex::Real time, amrex::Gpu::HostVector< amrex::Real > &h_avg_u, amrex::Gpu::HostVector< amrex::Real > &h_avg_v, amrex::Gpu::HostVector< amrex::Real > &h_avg_w, amrex::Gpu::HostVector< amrex::Real > &h_avg_rho, amrex::Gpu::HostVector< amrex::Real > &h_avg_th, amrex::Gpu::HostVector< amrex::Real > &h_avg_ksgs, amrex::Gpu::HostVector< amrex::Real > &h_avg_Kmv, amrex::Gpu::HostVector< amrex::Real > &h_avg_Khv, amrex::Gpu::HostVector< amrex::Real > &h_avg_qv, amrex::Gpu::HostVector< amrex::Real > &h_avg_qc, amrex::Gpu::HostVector< amrex::Real > &h_avg_qr, amrex::Gpu::HostVector< amrex::Real > &h_avg_wqv, amrex::Gpu::HostVector< amrex::Real > &h_avg_wqc, amrex::Gpu::HostVector< amrex::Real > &h_avg_wqr, amrex::Gpu::HostVector< amrex::Real > &h_avg_qi, amrex::Gpu::HostVector< amrex::Real > &h_avg_qs, amrex::Gpu::HostVector< amrex::Real > &h_avg_qg, amrex::Gpu::HostVector< amrex::Real > &h_avg_uu, amrex::Gpu::HostVector< amrex::Real > &h_avg_uv, amrex::Gpu::HostVector< amrex::Real > &h_avg_uw, amrex::Gpu::HostVector< amrex::Real > &h_avg_vv, amrex::Gpu::HostVector< amrex::Real > &h_avg_vw, amrex::Gpu::HostVector< amrex::Real > &h_avg_ww, amrex::Gpu::HostVector< amrex::Real > &h_avg_uth, amrex::Gpu::HostVector< amrex::Real > &h_avg_vth, amrex::Gpu::HostVector< amrex::Real > &h_avg_wth, amrex::Gpu::HostVector< amrex::Real > &h_avg_thth, amrex::Gpu::HostVector< amrex::Real > &h_avg_ku, amrex::Gpu::HostVector< amrex::Real > &h_avg_kv, amrex::Gpu::HostVector< amrex::Real > &h_avg_kw, amrex::Gpu::HostVector< amrex::Real > &h_avg_p, amrex::Gpu::HostVector< amrex::Real > &h_avg_pu, amrex::Gpu::HostVector< amrex::Real > &h_avg_pv, amrex::Gpu::HostVector< amrex::Real > &h_avg_pw, amrex::Gpu::HostVector< amrex::Real > &h_avg_wthv)
Definition: ERF_Write1DProfiles_stag.cpp:298
amrex::Real sampler_per
Definition: ERF.H:1375
InputSpongeData input_sponge_data
Definition: ERF.H:689
std::string restart_chkfile
Definition: ERF.H:924
bool metgrid_use_below_sfc
Definition: ERF.H:1089
void init_zphys(int lev, amrex::Real time)
Definition: ERF_MakeNewArrays.cpp:471
amrex::Vector< std::string > sampleptlogname
Definition: ERF.H:1382
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > xvel_bc_data
Definition: ERF.H:692
void InitData_pre()
Definition: ERF.cpp:631
amrex::IntVect & SamplePoint(int i)
Definition: ERF.H:1283
amrex::Vector< amrex::Vector< amrex::Real > > h_u_geos
Definition: ERF.H:1126
void InitializeLevelFromData(int lev, const amrex::MultiFab &initial_data)
void sum_integrated_quantities(amrex::Real time)
Definition: ERF_WriteScalarProfiles.cpp:14
amrex::Vector< std::string > datalogname
Definition: ERF.H:1379
void initHSE()
Initialize HSE.
Definition: ERF_Init1D.cpp:140
static amrex::Real column_loc_x
Definition: ERF.H:1100
amrex::Vector< std::unique_ptr< amrex::MultiFab > > ax
Definition: ERF.H:835
amrex::Vector< std::unique_ptr< amrex::MultiFab > > z_phys_nd
Definition: ERF.H:831
void MakeDiagnosticAverage(amrex::Vector< amrex::Real > &h_havg, amrex::MultiFab &S, int n)
Definition: ERF.cpp:1853
amrex::Real metgrid_proximity
Definition: ERF.H:1092
amrex::Vector< amrex::Real > h_havg_temperature
Definition: ERF.H:1154
amrex::Vector< std::unique_ptr< std::fstream > > sampleptlog
Definition: ERF.H:1381
static PlotFileType plotfile_type_2
Definition: ERF.H:1057
void poisson_wall_dist(int lev)
Definition: ERF_PoissonWallDist.cpp:20
void solve_with_gmres(int lev, amrex::Vector< amrex::MultiFab > &rhs, amrex::Vector< amrex::MultiFab > &p, amrex::Vector< amrex::Array< amrex::MultiFab, AMREX_SPACEDIM >> &fluxes)
Definition: ERF_SolveWithGMRES.cpp:12
void WritePlotFile(int which, PlotFileType plotfile_type, amrex::Vector< std::string > plot_var_names)
Definition: ERF_Plotfile.cpp:186
amrex::Vector< std::unique_ptr< amrex::MultiFab > > detJ_cc_src
Definition: ERF.H:840
amrex::Gpu::DeviceVector< amrex::Real > d_havg_pressure
Definition: ERF.H:1161
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SmnSmn_lev
Definition: ERF.H:815
const amrex::Vector< std::string > derived_names
Definition: ERF.H:973
amrex::Real start_time
Definition: ERF.H:920
amrex::Vector< std::unique_ptr< amrex::MultiFab > > ay_src
Definition: ERF.H:842
void sample_points(int lev, amrex::Real time, amrex::IntVect cell, amrex::MultiFab &mf)
Definition: ERF_WriteScalarProfiles.cpp:249
AMREX_FORCE_INLINE std::ostream & DataLog(int i)
Definition: ERF.H:1242
void writeJobInfo(const std::string &dir) const
Definition: ERF_WriteJobInfo.cpp:9
void solve_with_EB_mlmg(int lev, amrex::Vector< amrex::MultiFab > &rhs, amrex::Vector< amrex::MultiFab > &p, amrex::Vector< amrex::Array< amrex::MultiFab, AMREX_SPACEDIM >> &fluxes)
Definition: ERF_SolveWithEBMLMG.cpp:19
void ComputeDt(int step=-1)
Definition: ERF_ComputeTimestep.cpp:11
amrex::Vector< int > nsubsteps
Definition: ERF.H:724
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > yflux_imask
Definition: ERF.H:900
amrex::Vector< amrex::MultiFab > rW_new
Definition: ERF.H:764
amrex::Vector< std::unique_ptr< amrex::MultiFab > > az_new
Definition: ERF.H:849
std::unique_ptr< WriteBndryPlanes > m_w2d
Definition: ERF.H:1167
std::string plot_file_1
Definition: ERF.H:944
AMREX_FORCE_INLINE std::ostream & SampleLineLog(int i)
Definition: ERF.H:1270
amrex::Vector< amrex::Vector< amrex::MultiFab * > > lsm_flux
Definition: ERF.H:788
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Tau32_lev
Definition: ERF.H:813
std::string SamplePointLogName(int i) const noexcept
The filename of the ith sampleptlog file.
Definition: ERF.H:1393
void turbPert_update(const int lev, const amrex::Real dt)
Definition: ERF_InitTurbPert.cpp:12
void InitData_post()
Definition: ERF.cpp:677
void refinement_criteria_setup()
Definition: ERF_Tagging.cpp:127
static int nghost_eb_volume()
Definition: ERF.H:1412
bool metgrid_debug_dry
Definition: ERF.H:1084
static int bndry_output_planes_interval
Definition: ERF.H:1106
void init_geo_wind_profile(const std::string input_file, amrex::Vector< amrex::Real > &u_geos, amrex::Gpu::DeviceVector< amrex::Real > &u_geos_d, amrex::Vector< amrex::Real > &v_geos, amrex::Gpu::DeviceVector< amrex::Real > &v_geos_d, const amrex::Geometry &lgeom, const amrex::Vector< amrex::Real > &zlev_stag)
Definition: ERF_InitGeowind.cpp:10
void ERF_shared()
Definition: ERF.cpp:104
static int nghost_eb_full()
Definition: ERF.H:1415
void RemakeLevel(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Definition: ERF_MakeNewLevel.cpp:314
static void GotoNextLine(std::istream &is)
Definition: ERF_Checkpoint.cpp:15
void MakeNewLevelFromScratch(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Definition: ERF_MakeNewLevel.cpp:23
void init_stuff(int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, amrex::Vector< amrex::MultiFab > &lev_new, amrex::Vector< amrex::MultiFab > &lev_old, amrex::MultiFab &tmp_base_state, std::unique_ptr< amrex::MultiFab > &tmp_zphys_nd)
Definition: ERF_MakeNewArrays.cpp:23
amrex::Vector< amrex::MultiFab > zmom_crse_rhs
Definition: ERF.H:768
bool metgrid_debug_isothermal
Definition: ERF.H:1083
void initSponge()
Initialize sponge profiles.
Definition: ERF_InitSponge.cpp:35
void derive_diag_profiles(amrex::Real time, amrex::Gpu::HostVector< amrex::Real > &h_avg_u, amrex::Gpu::HostVector< amrex::Real > &h_avg_v, amrex::Gpu::HostVector< amrex::Real > &h_avg_w, amrex::Gpu::HostVector< amrex::Real > &h_avg_rho, amrex::Gpu::HostVector< amrex::Real > &h_avg_th, amrex::Gpu::HostVector< amrex::Real > &h_avg_ksgs, amrex::Gpu::HostVector< amrex::Real > &h_avg_Kmv, amrex::Gpu::HostVector< amrex::Real > &h_avg_Khv, amrex::Gpu::HostVector< amrex::Real > &h_avg_qv, amrex::Gpu::HostVector< amrex::Real > &h_avg_qc, amrex::Gpu::HostVector< amrex::Real > &h_avg_qr, amrex::Gpu::HostVector< amrex::Real > &h_avg_wqv, amrex::Gpu::HostVector< amrex::Real > &h_avg_wqc, amrex::Gpu::HostVector< amrex::Real > &h_avg_wqr, amrex::Gpu::HostVector< amrex::Real > &h_avg_qi, amrex::Gpu::HostVector< amrex::Real > &h_avg_qs, amrex::Gpu::HostVector< amrex::Real > &h_avg_qg, amrex::Gpu::HostVector< amrex::Real > &h_avg_uu, amrex::Gpu::HostVector< amrex::Real > &h_avg_uv, amrex::Gpu::HostVector< amrex::Real > &h_avg_uw, amrex::Gpu::HostVector< amrex::Real > &h_avg_vv, amrex::Gpu::HostVector< amrex::Real > &h_avg_vw, amrex::Gpu::HostVector< amrex::Real > &h_avg_ww, amrex::Gpu::HostVector< amrex::Real > &h_avg_uth, amrex::Gpu::HostVector< amrex::Real > &h_avg_vth, amrex::Gpu::HostVector< amrex::Real > &h_avg_wth, amrex::Gpu::HostVector< amrex::Real > &h_avg_thth, amrex::Gpu::HostVector< amrex::Real > &h_avg_ku, amrex::Gpu::HostVector< amrex::Real > &h_avg_kv, amrex::Gpu::HostVector< amrex::Real > &h_avg_kw, amrex::Gpu::HostVector< amrex::Real > &h_avg_p, amrex::Gpu::HostVector< amrex::Real > &h_avg_pu, amrex::Gpu::HostVector< amrex::Real > &h_avg_pv, amrex::Gpu::HostVector< amrex::Real > &h_avg_pw, amrex::Gpu::HostVector< amrex::Real > &h_avg_wthv)
Definition: ERF_Write1DProfiles.cpp:192
int real_width
Definition: ERF.H:1074
void MakeNewLevelFromCoarse(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Definition: ERF_MakeNewLevel.cpp:209
amrex::Vector< std::unique_ptr< amrex::FabFactory< amrex::FArrayBox > > > m_factory
Definition: ERF.H:1399
amrex::Vector< std::unique_ptr< amrex::MultiFab > > ax_src
Definition: ERF.H:841
void input_sponge(int lev)
Definition: ERF_InitSponge.cpp:17
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > d_u_geos
Definition: ERF.H:1127
void Define_ERFFillPatchers(int lev)
Definition: ERF.cpp:1925
AMREX_FORCE_INLINE int NumDataLogs() noexcept
Definition: ERF.H:1249
TurbulentPerturbation turbPert
Definition: ERF.H:1014
amrex::Vector< amrex::MultiFab > rW_old
Definition: ERF.H:763
void advance_lsm(int lev, amrex::MultiFab &, const amrex::Real &dt_advance)
Definition: ERF_AdvanceLSM.cpp:5
void FillCoarsePatch(int lev, amrex::Real time)
Definition: ERF_FillCoarsePatch.cpp:21
void ClearLevel(int lev) override
Definition: ERF_MakeNewLevel.cpp:473
static amrex::Vector< amrex::AMRErrorTag > ref_tags
Definition: ERF.H:1175
void make_physbcs(int lev)
Definition: ERF_MakeNewArrays.cpp:611
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > d_w_subsid
Definition: ERF.H:1124
amrex::Vector< ERFFillPatcher > FPr_w
Definition: ERF.H:807
int real_set_width
Definition: ERF.H:1075
void FillPatch(int lev, amrex::Real time, const amrex::Vector< amrex::MultiFab * > &mfs_vel, const amrex::Vector< amrex::MultiFab * > &mfs_mom, const amrex::MultiFab &old_base_state, const amrex::MultiFab &new_base_state, bool fillset=true, bool cons_only=false)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SFS_hfx2_lev
Definition: ERF.H:822
amrex::Vector< amrex::Vector< amrex::Real > > zlevels_stag
Definition: ERF.H:828
amrex::EBFArrayBoxFactory const & EBFactory(int lev) const noexcept
Definition: ERF.H:1404
amrex::Vector< amrex::Vector< amrex::MultiFab * > > lsm_data
Definition: ERF.H:787
amrex::Vector< amrex::Vector< amrex::Real > > stretched_dz_h
Definition: ERF.H:864
void WriteCheckpointFile() const
Definition: ERF_Checkpoint.cpp:25
int sampler_interval
Definition: ERF.H:1374
static int output_1d_column
Definition: ERF.H:1097
bool metgrid_debug_psfc
Definition: ERF.H:1085
amrex::Vector< std::unique_ptr< amrex::MultiFab > > az_src
Definition: ERF.H:843
static int fixed_mri_dt_ratio
Definition: ERF.H:937
int m_plot_int_2
Definition: ERF.H:948
amrex::Vector< amrex::Real > dt
Definition: ERF.H:729
static amrex::Real init_shrink
Definition: ERF.H:929
static amrex::Real dt_max_initial
Definition: ERF.H:931
static bool is_it_time_for_action(int nstep, amrex::Real time, amrex::Real dt, int action_interval, amrex::Real action_per)
Definition: ERF_WriteScalarProfiles.cpp:466
amrex::Vector< amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > > d_rayleigh_ptrs
Definition: ERF.H:1148
void InitData()
Definition: ERF.cpp:622
void advance_microphysics(int lev, amrex::MultiFab &cons_in, const amrex::Real &dt_advance, const int &iteration, const amrex::Real &time)
Definition: ERF_AdvanceMicrophysics.cpp:5
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Lwave
Definition: ERF.H:872
void project_velocities(int lev, amrex::Real dt, amrex::Vector< amrex::MultiFab > &vars, amrex::MultiFab &p)
Definition: ERF_PoissonSolve.cpp:10
void ParameterSanityChecks()
Definition: ERF.cpp:1681
void derive_stress_profiles(amrex::Gpu::HostVector< amrex::Real > &h_avg_tau11, amrex::Gpu::HostVector< amrex::Real > &h_avg_tau12, amrex::Gpu::HostVector< amrex::Real > &h_avg_tau13, amrex::Gpu::HostVector< amrex::Real > &h_avg_tau22, amrex::Gpu::HostVector< amrex::Real > &h_avg_tau23, amrex::Gpu::HostVector< amrex::Real > &h_avg_tau33, amrex::Gpu::HostVector< amrex::Real > &h_avg_hfx3, amrex::Gpu::HostVector< amrex::Real > &h_avg_q1fx3, amrex::Gpu::HostVector< amrex::Real > &h_avg_q2fx3, amrex::Gpu::HostVector< amrex::Real > &h_avg_diss)
Definition: ERF_Write1DProfiles.cpp:482
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Tau12_lev
Definition: ERF.H:811
int cf_width
Definition: ERF.H:802
amrex::Vector< amrex::Vector< amrex::Vector< amrex::Real > > > h_rayleigh_ptrs
Definition: ERF.H:1142
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > zflux_imask
Definition: ERF.H:901
bool m_expand_plotvars_to_unif_rr
Definition: ERF.H:946
int plot_file_on_restart
Definition: ERF.H:913
void Construct_ERFFillPatchers(int lev)
Definition: ERF.cpp:1899
void post_update(amrex::MultiFab &state_mf, amrex::Real time, const amrex::Geometry &geom)
amrex::Vector< int > num_boxes_at_level
Definition: ERF.H:719
amrex::FabFactory< amrex::FArrayBox > const & Factory(int lev) const noexcept
Definition: ERF.H:1402
static amrex::Real sub_cfl
Definition: ERF.H:928
static void print_error(MPI_Comm, const std::string &msg)
Definition: ERF_ConsoleIO.cpp:43
static int ng_dens_hse
Definition: ERF.H:1113
std::unique_ptr< ReadBndryPlanes > m_r2d
Definition: ERF.H:1168
amrex::Vector< std::unique_ptr< amrex::MultiFab > > ay
Definition: ERF.H:836
amrex::Vector< std::unique_ptr< amrex::MultiFab > > ax_new
Definition: ERF.H:847
amrex::Vector< amrex::Vector< amrex::MultiFab * > > qmoist
Definition: ERF.H:771
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SFS_q1fx1_lev
Definition: ERF.H:824
static amrex::Real getCPUTime()
Definition: ERF.H:1312
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > zvel_bc_data
Definition: ERF.H:694
void derive_stress_profiles_stag(amrex::Gpu::HostVector< amrex::Real > &h_avg_tau11, amrex::Gpu::HostVector< amrex::Real > &h_avg_tau12, amrex::Gpu::HostVector< amrex::Real > &h_avg_tau13, amrex::Gpu::HostVector< amrex::Real > &h_avg_tau22, amrex::Gpu::HostVector< amrex::Real > &h_avg_tau23, amrex::Gpu::HostVector< amrex::Real > &h_avg_tau33, amrex::Gpu::HostVector< amrex::Real > &h_avg_hfx3, amrex::Gpu::HostVector< amrex::Real > &h_avg_q1fx3, amrex::Gpu::HostVector< amrex::Real > &h_avg_q2fx3, amrex::Gpu::HostVector< amrex::Real > &h_avg_diss)
Definition: ERF_Write1DProfiles_stag.cpp:612
LandSurface lsm
Definition: ERF.H:786
static amrex::Real change_max
Definition: ERF.H:930
void AverageDownTo(int crse_lev, int scomp, int ncomp)
Definition: ERF_AverageDown.cpp:36
void setRecordSampleLineInfo(int i, int lev, amrex::IntVect &cell, const std::string &filename)
Definition: ERF.H:1356
void setSpongeRefFromSounding(bool restarting)
Set sponge mean profiles from input sounding.
Definition: ERF_InitSponge.cpp:65
amrex::Vector< std::unique_ptr< amrex::MultiFab > > detJ_cc_new
Definition: ERF.H:846
amrex::Vector< std::unique_ptr< amrex::MultiFab > > mapfac_v
Definition: ERF.H:862
amrex::Array< std::string, 2 *AMREX_SPACEDIM > domain_bc_type
Definition: ERF.H:887
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > yvel_bc_data
Definition: ERF.H:693
void init_uniform(int lev)
Definition: ERF_InitUniform.cpp:17
static AMREX_FORCE_INLINE int ComputeGhostCells(const AdvChoice &advChoice, bool use_num_diff)
Definition: ERF.H:1187
static void writeBuildInfo(std::ostream &os)
Definition: ERF_WriteJobInfo.cpp:137
amrex::GpuArray< ERF_BC, AMREX_SPACEDIM *2 > phys_bc_type
Definition: ERF.H:896
std::string check_file
Definition: ERF.H:960
amrex::Vector< amrex::IntVect > samplepoint
Definition: ERF.H:1383
amrex::Vector< amrex::Real > h_havg_qv
Definition: ERF.H:1156
amrex::Real volWgtSumMF(int lev, const amrex::MultiFab &mf, int comp, const amrex::MultiFab &mapfac, bool local, bool finemask)
Definition: ERF_WriteScalarProfiles.cpp:378
static void print_usage(MPI_Comm, std::ostream &)
Definition: ERF_ConsoleIO.cpp:26
amrex::Vector< amrex::MultiFab > rV_old
Definition: ERF.H:761
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Tau11_lev
Definition: ERF.H:810
void make_eb_regular()
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Hwave
Definition: ERF.H:871
amrex::Vector< int > istep
Definition: ERF.H:723
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > xflux_imask
Definition: ERF.H:899
void derive_upwp(amrex::Vector< amrex::Real > &h_havg)
int metgrid_order
Definition: ERF.H:1093
bool finished_wave
Definition: ERF.H:875
void ReadCheckpointFile()
Definition: ERF_Checkpoint.cpp:273
bool writeNow(const amrex::Real cur_time, const amrex::Real dt, const int nstep, const int plot_int, const amrex::Real plot_per)
Definition: ERF.cpp:1976
amrex::Vector< amrex::Vector< amrex::MultiFab > > vars_old
Definition: ERF.H:736
amrex::MultiFab & build_fine_mask(int lev)
Definition: ERF_WriteScalarProfiles.cpp:429
amrex::Vector< std::unique_ptr< ERFPhysBCFunct_base > > physbcs_base
Definition: ERF.H:751
void init_thin_body(int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Definition: ERF_MakeNewLevel.cpp:512
AMREX_FORCE_INLINE std::ostream & SamplePointLog(int i)
Definition: ERF.H:1256
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_neumann_vals
Definition: ERF.H:893
void advance_dycore(int level, amrex::Vector< amrex::MultiFab > &state_old, amrex::Vector< amrex::MultiFab > &state_new, amrex::MultiFab &xvel_old, amrex::MultiFab &yvel_old, amrex::MultiFab &zvel_old, amrex::MultiFab &xvel_new, amrex::MultiFab &yvel_new, amrex::MultiFab &zvel_new, amrex::MultiFab &source, amrex::MultiFab &xmom_src, amrex::MultiFab &ymom_src, amrex::MultiFab &zmom_src, amrex::Geometry fine_geom, amrex::Real dt, amrex::Real time)
Definition: ERF_AdvanceDycore.cpp:37
void setRecordSamplePointInfo(int i, int lev, amrex::IntVect &cell, const std::string &filename)
Definition: ERF.H:1339
static int column_interval
Definition: ERF.H:1098
static void print_summary(std::ostream &)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Tau22_lev
Definition: ERF.H:810
void turbPert_amplitude(const int lev)
Definition: ERF_InitTurbPert.cpp:41
void sample_lines(int lev, amrex::Real time, amrex::IntVect cell, amrex::MultiFab &mf)
Definition: ERF_WriteScalarProfiles.cpp:287
void initializeMicrophysics(const int &)
Definition: ERF.cpp:1231
bool plot_lsm
Definition: ERF.H:953
const amrex::Vector< std::string > cons_names
Definition: ERF.H:968
void FillPatch(int lev, amrex::Real time, const amrex::Vector< amrex::MultiFab * > &mfs_vel, bool cons_only=false)
void timeStep(int lev, amrex::Real time, int iteration)
Definition: ERF_TimeStep.cpp:16
amrex::Real m_plot_per_2
Definition: ERF.H:950
Definition: ERF_LandSurface.H:14
Definition: ERF_MultiBlockContainer.H:8
@ pres
Definition: ERF_Kessler.H:33
@ T
Definition: ERF_IndexDefines.H:99
Definition: ERF_ConsoleIO.cpp:12
Definition: ERF_AdvStruct.H:19
AdvType moistscal_horiz_adv_type
Definition: ERF_AdvStruct.H:286
AdvType dycore_vert_adv_type
Definition: ERF_AdvStruct.H:283
AdvType moistscal_vert_adv_type
Definition: ERF_AdvStruct.H:287
AdvType dryscal_horiz_adv_type
Definition: ERF_AdvStruct.H:284
AdvType dycore_horiz_adv_type
Definition: ERF_AdvStruct.H:282
AdvType dryscal_vert_adv_type
Definition: ERF_AdvStruct.H:285
Definition: ERF_InputSoundingData.H:22
Definition: ERF_InputSpongeData.H:19
Definition: ERF_DataStruct.H:82
Definition: ERF_TurbPertStruct.H:18