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