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