ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF.H
Go to the documentation of this file.
1 #ifndef ERF_H_
2 #define ERF_H_
3 
4 #include <string>
5 #include <limits>
6 #include <memory>
7 
8 #ifdef _OPENMP
9 #include <omp.h>
10 #endif
11 
12 #include <AMReX_AmrCore.H>
13 #include <AMReX_BCRec.H>
14 #include <AMReX_Print.H>
15 
16 #include <AMReX_ParallelDescriptor.H>
17 #include <AMReX_ParmParse.H>
18 #include <AMReX_MultiFabUtil.H>
19 #include <AMReX_FillPatchUtil.H>
20 #include <AMReX_VisMF.H>
21 #include <AMReX_PhysBCFunct.H>
22 #include <AMReX_YAFluxRegister.H>
23 #include <AMReX_ErrorList.H>
24 #include <AMReX_EBFabFactory.H>
25 #include <AMReX_EBMultiFabUtil.H>
26 
27 #include <ERF_EB.H>
28 
29 #ifdef ERF_USE_FFT
30 #include <AMReX_FFT_Poisson.H>
31 #endif
32 
33 #ifdef AMREX_MEM_PROFILING
34 #include <AMReX_MemProfiler.H>
35 #endif
36 
37 #include "ERF_ProbCommon.H"
38 
39 #include <ERF_IndexDefines.H>
40 #include <ERF_DataStruct.H>
41 #include <ERF_TurbPertStruct.H>
42 #include <ERF_InputSoundingData.H>
43 #include <ERF_InputSpongeData.H>
44 #include <ERF_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 
62 #ifdef ERF_USE_WINDFARM
63 #include "ERF_WindFarm.H"
64 #endif
65 
66 #include <ERF_RadiationModels.H>
67 
68 #ifdef ERF_USE_SHOC
69 #include "ERF_ShocInterface.H"
70 #endif
71 
72 #ifdef ERF_USE_P3
73 #include "ERF_P3Interface.H"
74 #endif
75 
76 #include <iostream>
77 
78 #ifdef AMREX_LAZY
79 #include <AMReX_Lazy.H>
80 #endif
81 
82 #ifndef AMREX_USE_MPI
83 using amrex::MPI_COMM_WORLD;
84 using amrex::MPI_Comm;
85 #endif
86 
87 #ifdef ERF_USE_MULTIBLOCK
88 class MultiBlockContainer;
89 #endif
90 
91 /**
92  * Enum of possible interpolation types between coarse/fine
93 */
94 AMREX_ENUM(StateInterpType,
95  FullState, Perturbational
96 );
97 
98 /**
99  * Enum of possible plotfile types
100 */
101 AMREX_ENUM(PlotFileType,
102  None, Amrex, Netcdf
103 );
104 
105 /**
106  * Main class in ERF code, instantiated from main.cpp
107 */
108 
109 class ERF
110  : public amrex::AmrCore
111 {
112 public:
113 
114  ////////////////
115  // public member functions
116 
117  // constructor - reads in parameters from inputs file
118  // - sizes multilevel arrays and data structures
119  ERF ();
120  ~ERF () override;
121 
122  void ERF_shared ();
123 
124  // Declare a default move constructor so we ensure the destructor is
125  // not called when we return an object of this class by value
126  ERF (ERF&&) noexcept = delete;
127 
128  // Declare a default move assignment operator
129  ERF& operator=(ERF&& other) noexcept = delete;
130 
131  // Delete the copy constructor
132  ERF (const ERF& other) = delete;
133  //
134  // Delete the copy assignment operator
135  ERF& operator=(const ERF& other) = delete;
136 
137  // Advance solution to final time
138  void Evolve ();
139 
140  // Tag cells for refinement
141  void ErrorEst (int lev, amrex::TagBoxArray& tags, amrex::Real time, int ngrow) override;
142 
143  // Refinement criteria for hurricanes
144  void HurricaneTracker (int lev,
145  amrex::Real time,
146  const amrex::MultiFab& cc_vel,
147  const amrex::Real velmag_threshold,
148  amrex::TagBoxArray* tags = nullptr);
149 
150  bool FindInitialEye (int lev,
151  const amrex::MultiFab& cc_vel,
152  const amrex::Real velmag_threshold,
153  amrex::Real& eye_x, amrex::Real& eye_y);
154 
155  // Hurricane track line output as VTKPolyline
156  amrex::Vector<std::array<amrex::Real, 2>> hurricane_track_xy;
157  amrex::Vector<std::array<amrex::Real, 2>> hurricane_eye_track_xy,
163  amrex::Vector<amrex::Vector<amrex::MultiFab>> forecast_state_1,
166  amrex::Vector<amrex::MultiFab> surface_state_1,
169 
170  std::string MakeVTKFilename(int nstep);
171  std::string MakeVTKFilename_TrackerCircle(int nstep);
172  std::string MakeVTKFilename_EyeTracker_xy(int nstep);
173  std::string MakeFilename_EyeTracker_latlon(int nstep);
174  std::string MakeFilename_EyeTracker_maxvel(int nstep);
175  std::string MakeFilename_EyeTracker_minpressure(int nstep);
176  void WriteVTKPolyline(const std::string& filename,
177  amrex::Vector<std::array<amrex::Real, 2>>& points_xy);
178 
179  void WriteLinePlot(const std::string& filename,
180  amrex::Vector<std::array<amrex::Real, 2>>& points_xy);
181 
182  // Initialize multilevel data
183  void InitData ();
184 
185  // Initialize multilevel data before MultiBlock
186  void InitData_pre ();
187 
188  // Initialize multilevel data after MultiBlock
189  void InitData_post ();
190 
191  // Interpolate 2D arrays lat_m, lon_m, sinPhi, cosPhi
192  void Interp2DArrays (int lev, const amrex::BoxArray& my_ba2d, const amrex::DistributionMapping& my_dm);
193 
195 
196  // Compute the divergence -- whether EB, no-terrain, flat terrain or general terrain
197  void compute_divergence (int lev, amrex::MultiFab& rhs, amrex::Array<amrex::MultiFab const*,AMREX_SPACEDIM> rho0_u_const,
198  amrex::Geometry const& geom_at_lev);
199 
200  // Project the velocities to be divergence-free -- this is only relevant if anelastic == 1
201  void project_initial_velocity (int lev, amrex::Real time, amrex::Real dt);
202  void project_momenta (int lev, amrex::Real l_time, amrex::Real l_dt, amrex::Vector<amrex::MultiFab >& vars);
203 
204  // Project the velocities to be divergence-free with a thin body
205  void project_velocity_tb (int lev, amrex::Real dt, amrex::Vector<amrex::MultiFab >& vars);
206 
207  // Calculate wall distance by solving a Poisson equation
208  void poisson_wall_dist (int lev);
209 
210  void make_subdomains(const amrex::BoxList& ba, amrex::Vector<amrex::BoxArray>& bins);
211 
212 #ifdef ERF_USE_FFT
213  void solve_with_fft (int lev, const amrex::Box& subdomain,
214  amrex::MultiFab& rhs, amrex::MultiFab& p,
215  amrex::Array<amrex::MultiFab,AMREX_SPACEDIM>& fluxes);
216 #endif
217  void solve_with_gmres (int lev, const amrex::Box& subdomain,
218  amrex::MultiFab& rhs, amrex::MultiFab& p,
219  amrex::Array<amrex::MultiFab,AMREX_SPACEDIM>& fluxes,
220  amrex::MultiFab& ax_sub, amrex::MultiFab& ay_sub,
221  amrex::MultiFab& az_sub,
222  amrex::MultiFab& /*dJ_sub*/, amrex::MultiFab& znd_sub);
223 
224  void ImposeBCsOnPhi (int lev, amrex::MultiFab& phi, const amrex::Box& subdomain);
225 
226  // Init (NOT restart or regrid)
227  void init_only (int lev, amrex::Real time);
228 
229  // Restart
230  void restart ();
231 
232  // Check for NaNs and other bad stuff
233  void check_state_for_nans(amrex::MultiFab const& S);
234  void check_vels_for_nans(amrex::MultiFab const& xvel, amrex::MultiFab const& yvel, amrex::MultiFab const& zvel);
235  void check_for_negative_theta (amrex::MultiFab& S);
236  void check_for_low_temp (amrex::MultiFab& S);
237 
238  // Is it time to write a plotfile or checkpoint?
239  bool writeNow (const amrex::Real cur_time, const int nstep, const int plot_int,
240  const amrex::Real plot_per, const amrex::Real dt_0, amrex::Real& last_file_time);
241 
242  // Called after every level 0 timestep
243  void post_timestep (int nstep, amrex::Real time, amrex::Real dt_lev);
244 
245  // Diagnostics
249 
250  void write_1D_profiles (amrex::Real time);
252 
254 
255  // Fill the physical boundary conditions for cell-centered velocity (diagnostic only)
256  void FillBdyCCVels (amrex::Vector<amrex::MultiFab>& mf_cc_vel, int levc=0);
257 
258  void sample_points (int lev, amrex::Real time, amrex::IntVect cell, amrex::MultiFab& mf);
259  void sample_lines (int lev, amrex::Real time, amrex::IntVect cell, amrex::MultiFab& mf);
260 
262  amrex::Gpu::HostVector<amrex::Real>& h_avg_u , amrex::Gpu::HostVector<amrex::Real>& h_avg_v,
263  amrex::Gpu::HostVector<amrex::Real>& h_avg_w , amrex::Gpu::HostVector<amrex::Real>& h_avg_rho,
264  amrex::Gpu::HostVector<amrex::Real>& h_avg_th , amrex::Gpu::HostVector<amrex::Real>& h_avg_ksgs,
265  amrex::Gpu::HostVector<amrex::Real>& h_avg_Kmv , amrex::Gpu::HostVector<amrex::Real>& h_avg_Khv,
266  amrex::Gpu::HostVector<amrex::Real>& h_avg_qv , amrex::Gpu::HostVector<amrex::Real>& h_avg_qc,
267  amrex::Gpu::HostVector<amrex::Real>& h_avg_qr ,
268  amrex::Gpu::HostVector<amrex::Real>& h_avg_wqv , amrex::Gpu::HostVector<amrex::Real>& h_avg_wqc,
269  amrex::Gpu::HostVector<amrex::Real>& h_avg_wqr ,
270  amrex::Gpu::HostVector<amrex::Real>& h_avg_qi , amrex::Gpu::HostVector<amrex::Real>& h_avg_qs,
271  amrex::Gpu::HostVector<amrex::Real>& h_avg_qg ,
272  amrex::Gpu::HostVector<amrex::Real>& h_avg_uu , amrex::Gpu::HostVector<amrex::Real>& h_avg_uv,
273  amrex::Gpu::HostVector<amrex::Real>& h_avg_uw,
274  amrex::Gpu::HostVector<amrex::Real>& h_avg_vv , amrex::Gpu::HostVector<amrex::Real>& h_avg_vw,
275  amrex::Gpu::HostVector<amrex::Real>& h_avg_ww,
276  amrex::Gpu::HostVector<amrex::Real>& h_avg_uth , amrex::Gpu::HostVector<amrex::Real>& h_avg_vth,
277  amrex::Gpu::HostVector<amrex::Real>& h_avg_wth, amrex::Gpu::HostVector<amrex::Real>& h_avg_thth,
278  amrex::Gpu::HostVector<amrex::Real>& h_avg_ku, amrex::Gpu::HostVector<amrex::Real>& h_avg_kv,
279  amrex::Gpu::HostVector<amrex::Real>& h_avg_kw,
280  amrex::Gpu::HostVector<amrex::Real>& h_avg_p,
281  amrex::Gpu::HostVector<amrex::Real>& h_avg_pu, amrex::Gpu::HostVector<amrex::Real>& h_avg_pv,
282  amrex::Gpu::HostVector<amrex::Real>& h_avg_pw, amrex::Gpu::HostVector<amrex::Real>& h_avg_wthv);
284  amrex::Gpu::HostVector<amrex::Real>& h_avg_u , amrex::Gpu::HostVector<amrex::Real>& h_avg_v,
285  amrex::Gpu::HostVector<amrex::Real>& h_avg_w , amrex::Gpu::HostVector<amrex::Real>& h_avg_rho,
286  amrex::Gpu::HostVector<amrex::Real>& h_avg_th , amrex::Gpu::HostVector<amrex::Real>& h_avg_ksgs,
287  amrex::Gpu::HostVector<amrex::Real>& h_avg_Kmv , amrex::Gpu::HostVector<amrex::Real>& h_avg_Khv,
288  amrex::Gpu::HostVector<amrex::Real>& h_avg_qv , amrex::Gpu::HostVector<amrex::Real>& h_avg_qc,
289  amrex::Gpu::HostVector<amrex::Real>& h_avg_qr ,
290  amrex::Gpu::HostVector<amrex::Real>& h_avg_wqv , amrex::Gpu::HostVector<amrex::Real>& h_avg_wqc,
291  amrex::Gpu::HostVector<amrex::Real>& h_avg_wqr ,
292  amrex::Gpu::HostVector<amrex::Real>& h_avg_qi , amrex::Gpu::HostVector<amrex::Real>& h_avg_qs,
293  amrex::Gpu::HostVector<amrex::Real>& h_avg_qg ,
294  amrex::Gpu::HostVector<amrex::Real>& h_avg_uu , amrex::Gpu::HostVector<amrex::Real>& h_avg_uv,
295  amrex::Gpu::HostVector<amrex::Real>& h_avg_uw,
296  amrex::Gpu::HostVector<amrex::Real>& h_avg_vv , amrex::Gpu::HostVector<amrex::Real>& h_avg_vw,
297  amrex::Gpu::HostVector<amrex::Real>& h_avg_ww,
298  amrex::Gpu::HostVector<amrex::Real>& h_avg_uth , amrex::Gpu::HostVector<amrex::Real>& h_avg_vth,
299  amrex::Gpu::HostVector<amrex::Real>& h_avg_wth, amrex::Gpu::HostVector<amrex::Real>& h_avg_thth,
300  amrex::Gpu::HostVector<amrex::Real>& h_avg_ku, amrex::Gpu::HostVector<amrex::Real>& h_avg_kv,
301  amrex::Gpu::HostVector<amrex::Real>& h_avg_kw,
302  amrex::Gpu::HostVector<amrex::Real>& h_avg_p,
303  amrex::Gpu::HostVector<amrex::Real>& h_avg_pu, amrex::Gpu::HostVector<amrex::Real>& h_avg_pv,
304  amrex::Gpu::HostVector<amrex::Real>& h_avg_pw, amrex::Gpu::HostVector<amrex::Real>& h_avg_wthv);
305 
306  void derive_stress_profiles (amrex::Gpu::HostVector<amrex::Real>& h_avg_tau11, amrex::Gpu::HostVector<amrex::Real>& h_avg_tau12,
307  amrex::Gpu::HostVector<amrex::Real>& h_avg_tau13, amrex::Gpu::HostVector<amrex::Real>& h_avg_tau22,
308  amrex::Gpu::HostVector<amrex::Real>& h_avg_tau23, amrex::Gpu::HostVector<amrex::Real>& h_avg_tau33,
309  amrex::Gpu::HostVector<amrex::Real>& h_avg_hfx3, amrex::Gpu::HostVector<amrex::Real>& h_avg_q1fx3,
310  amrex::Gpu::HostVector<amrex::Real>& h_avg_q2fx3, amrex::Gpu::HostVector<amrex::Real>& h_avg_diss);
311  void derive_stress_profiles_stag (amrex::Gpu::HostVector<amrex::Real>& h_avg_tau11, amrex::Gpu::HostVector<amrex::Real>& h_avg_tau12,
312  amrex::Gpu::HostVector<amrex::Real>& h_avg_tau13, amrex::Gpu::HostVector<amrex::Real>& h_avg_tau22,
313  amrex::Gpu::HostVector<amrex::Real>& h_avg_tau23, amrex::Gpu::HostVector<amrex::Real>& h_avg_tau33,
314  amrex::Gpu::HostVector<amrex::Real>& h_avg_hfx3, amrex::Gpu::HostVector<amrex::Real>& h_avg_q1fx3,
315  amrex::Gpu::HostVector<amrex::Real>& h_avg_q2fx3, amrex::Gpu::HostVector<amrex::Real>& h_avg_diss);
316 
317  // Perform the volume-weighted sum
319  volWgtSumMF (int lev, const amrex::MultiFab& mf, int comp,
320  const amrex::MultiFab& dJ, const amrex::MultiFab& mfx, const amrex::MultiFab& mfy,
321  bool finemask, bool local = true);
322 
323  void
324  volWgtColumnSum (int lev, const amrex::MultiFab& mf, int comp,
325  amrex::MultiFab& mf_2d, const amrex::MultiFab& dJ);
326 
327  // Decide if it is time to take an action
328  static bool is_it_time_for_action (int nstep, amrex::Real time, amrex::Real dt,
329  int action_interval, amrex::Real action_per);
330 
331  // Make a new level using provided BoxArray and DistributionMapping and
332  // fill with interpolated coarse level data.
333  // overrides the pure virtual function in AmrCore
334  void MakeNewLevelFromCoarse (int lev, amrex::Real time, const amrex::BoxArray& ba,
335  const amrex::DistributionMapping& dm) override;
336 
337  // Remake an existing level using provided BoxArray and DistributionMapping and
338  // fill with existing fine and coarse data.
339  // overrides the pure virtual function in AmrCore
340  void RemakeLevel (int lev, amrex::Real time, const amrex::BoxArray& ba,
341  const amrex::DistributionMapping& dm) override;
342 
343  // Delete level data
344  // overrides the pure virtual function in AmrCore
345  void ClearLevel (int lev) override;
346 
347  // Make a new level from scratch using provided BoxArray and DistributionMapping.
348  // Only used during initialization.
349  // overrides the pure virtual function in AmrCore
350  void MakeNewLevelFromScratch (int lev, amrex::Real time, const amrex::BoxArray& ba,
351  const amrex::DistributionMapping& dm) override;
352 
353  // compute dt from CFL considerations
354  amrex::Real estTimeStep (int lev, long& dt_fast_ratio) const;
355 
356 #ifdef ERF_USE_WW3_COUPLING
357  //amrex::Print() << " About to call send_to_ww3 from ERF.H" << std::endl;
358  void send_to_ww3(int lev);
359  //amrex::Print() << " About to call read_waves from ERF.H" << std::endl;
360  void read_waves(int lev);
361  //void send_to_ww3 (int lev);
362  //void read_waves (int lev);
363 
364  //void send_to_ww3 (int lev);
365 #endif
366 
367  // Interface for advancing the data at one level by one "slow" timestep
368  void advance_dycore (int level,
369  amrex::Vector<amrex::MultiFab>& state_old,
370  amrex::Vector<amrex::MultiFab>& state_new,
371  amrex::MultiFab& xvel_old, amrex::MultiFab& yvel_old, amrex::MultiFab& zvel_old,
372  amrex::MultiFab& xvel_new, amrex::MultiFab& yvel_new, amrex::MultiFab& zvel_new,
373  amrex::MultiFab& source, amrex::MultiFab& xmom_src,
374  amrex::MultiFab& ymom_src, amrex::MultiFab& zmom_src,
375  amrex::MultiFab& buoyancy, amrex::Geometry fine_geom,
376  amrex::Real dt, amrex::Real time);
377 
378  void advance_microphysics (int lev,
379  amrex::MultiFab& cons_in,
380  const amrex::Real& dt_advance,
381  const int& iteration,
382  const amrex::Real& time);
383 
384  void advance_lsm (int lev,
385  amrex::MultiFab& cons_in,
386  amrex::MultiFab& xvel_in,
387  amrex::MultiFab& yvel_in,
388  const amrex::Real& dt_advance);
389 
390  void advance_radiation (int lev,
391  amrex::MultiFab& cons_in,
392  const amrex::Real& dt_advance);
393 
394 #ifdef ERF_USE_SHOC
395  void compute_shoc_tendencies (int lev,
396  amrex::MultiFab* cons,
397  amrex::MultiFab* xvel,
398  amrex::MultiFab* yvel,
399  amrex::MultiFab* zvel,
400  amrex::Real* w_subsid,
401  amrex::MultiFab* tau13,
402  amrex::MultiFab* tau23,
403  amrex::MultiFab* hfx3,
404  amrex::MultiFab* qfx3,
405  amrex::MultiFab* eddyDiffs,
406  amrex::MultiFab* z_phys_nd,
407  const amrex::Real& dt_advance);
408 #endif
409 
410 #ifdef ERF_USE_P3
411  void compute_p3_tendencies (int lev,
412  amrex::MultiFab& cons_in,
413  amrex::MultiFab& source,
414  const amrex::Real& dt_advance);
415 #endif
416 
417  void build_fine_mask (int lev, amrex::MultiFab& fine_mask);
418 
419  void MakeHorizontalAverages ();
420  void MakeDiagnosticAverage (amrex::Vector<amrex::Real>& h_havg, amrex::MultiFab& S, int n);
421  void derive_upwp (amrex::Vector<amrex::Real>& h_havg);
422 
423  // Write plotfile to disk
424  void Write3DPlotFile (int which, PlotFileType plotfile_type, amrex::Vector<std::string> plot_var_names);
425  void Write2DPlotFile (int which, PlotFileType plotfile_type, amrex::Vector<std::string> plot_var_names);
426 
427  // Write subvolume of data into a different "plotfile"
428  void WriteSubvolume (int isub, amrex::Vector<std::string> subvol_var_names);
429 
430  void WriteMultiLevelPlotfileWithTerrain (const std::string &plotfilename,
431  int nlevels,
432  const amrex::Vector<const amrex::MultiFab*> &mf,
433  const amrex::Vector<const amrex::MultiFab*> &mf_nd,
434  const amrex::Vector<std::string> &varnames,
435  const amrex::Vector<amrex::Geometry>& my_geom,
436  amrex::Real time,
437  const amrex::Vector<int> &level_steps,
438  const amrex::Vector<amrex::IntVect>& my_ref_ratio,
439  const std::string &versionName = "HyperCLaw-V1.1",
440  const std::string &levelPrefix = "Level_",
441  const std::string &mfPrefix = "Cell",
442  const amrex::Vector<std::string>& extra_dirs = amrex::Vector<std::string>()) const;
443 
444 
445  void WriteGenericPlotfileHeaderWithTerrain (std::ostream &HeaderFile,
446  int nlevels,
447  const amrex::Vector<amrex::BoxArray> &bArray,
448  const amrex::Vector<std::string> &varnames,
449  const amrex::Vector<amrex::Geometry>& my_geom,
450  amrex::Real time,
451  const amrex::Vector<int> &level_steps,
452  const amrex::Vector<amrex::IntVect>& my_ref_ratio,
453  const std::string &versionName,
454  const std::string &levelPrefix,
455  const std::string &mfPrefix) const;
456 
457  void erf_enforce_hse (int lev,
458  amrex::MultiFab& dens, amrex::MultiFab& pres, amrex::MultiFab& pi,
459  amrex::MultiFab& th, amrex::MultiFab& qv,
460  std::unique_ptr<amrex::MultiFab>& z_cc);
461 
462 #ifdef ERF_USE_NETCDF
463  //! Write a timestep to 1D vertical column output for coupling
464  void writeToNCColumnFile (int lev,
465  const std::string& colfile_name, amrex::Real xloc, amrex::Real yloc,
466  amrex::Real time);
467 #endif //ERF_USE_NETCDF
468 
469  void init_from_input_sounding (int lev);
470 
471  void init_immersed_forcing (int lev);
472 
473  void input_sponge (int lev);
474 
475  void init_from_hse (int lev);
476 
477  void init_thin_body (int lev, const amrex::BoxArray& ba, const amrex::DistributionMapping& dm);
478 
479  void FillForecastStateMultiFabs(const int lev,
480  const std::string& filename,
481  const std::unique_ptr<amrex::MultiFab>& z_phys_nd,
482  amrex::Vector<amrex::Vector<amrex::MultiFab>>& forecast_state);
483 
484  void FillSurfaceStateMultiFabs(const int lev,
485  const std::string& filename,
486  amrex::Vector<amrex::MultiFab>& surface_state);
487 
488 
489  void WeatherDataInterpolation(const int nlevs,
490  const amrex::Real time,
491  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& z_phys_nd,
492  bool regrid_forces_file_read);
493 
494  void SurfaceDataInterpolation(const int nlevs,
495  const amrex::Real time,
496  amrex::Vector<std::unique_ptr<amrex::MultiFab>>& z_phys_nd,
497  bool regrid_forces_file_read);
498 
499  void create_random_perturbations (const int lev,
500  amrex::MultiFab& cons_pert,
501  amrex::MultiFab& xvel_pert,
502  amrex::MultiFab& yvel_pert,
503  amrex::MultiFab& zvel_pert);
504 
505  void apply_gaussian_smoothing_to_perturbations (const int lev,
506  amrex::MultiFab& cons_pert,
507  amrex::MultiFab& xvel_pert,
508  amrex::MultiFab& yvel_pert,
509  amrex::MultiFab& zvel_pert);
510 
511 #ifdef ERF_USE_MULTIBLOCK
512  // constructor used when ERF is created by a multiblock driver
513  // calls AmrCore constructor -> AmrMesh constructor
514  ERF (const amrex::RealBox& rb, int max_level_in,
515  const amrex::Vector<int>& n_cell_in, int coord,
516  const amrex::Vector<amrex::IntVect>& ref_ratio,
517  const amrex::Array<int,AMREX_SPACEDIM>& is_per,
518  std::string prefix);
519 
520  // Advance a block specified number of time steps
521  void Evolve_MB (int MBstep, int max_block_step);
522 
523  // get the current time values and dt
524  amrex::Real get_t_old (int lev) { return t_old[lev]; }
525  amrex::Real get_t_new (int lev) { return t_new[lev]; }
526  amrex::Real get_dt (int lev) { return dt[lev]; }
527 
528  // Set parmparse prefix for MultiBlock
529  void SetParmParsePrefix (std::string name) { pp_prefix = name; }
530 
531  // Set 'this' from multiblock container
532  void SetMultiBlockPointer (MultiBlockContainer *mbc) { m_mbc = mbc; }
533 
534  // Public data copy for MB
535  std::vector<amrex::Box> domain_p;
536 
537  // Primary solution data data containers
538  amrex::Vector<amrex::Vector<amrex::MultiFab> > vars_new;
539  amrex::Vector<amrex::Vector<amrex::MultiFab> > vars_old;
540 
541  // Velocity time averaged field
542  amrex::Vector<std::unique_ptr<amrex::MultiFab>> vel_t_avg;
543  amrex::Vector<amrex::Real> t_avg_cnt;
544 #endif
545 
546  std::string pp_prefix {"erf"};
547 
548  void fill_from_bndryregs (const amrex::Vector<amrex::MultiFab*>& mfs,
549  amrex::Real time);
550 
551 #ifdef ERF_USE_NETCDF
552  // Version in which we always set the velocity values but don't set theta
553  void fill_from_realbdy (const amrex::Vector<amrex::MultiFab*>& mfs,
554  amrex::Real time,
555  bool cons_only,
556  int icomp_cons,
557  int ncomp_cons,
558  amrex::IntVect ngvect_cons,
559  amrex::IntVect ngvect_vels);
560 
561  // Version in which we only set the boundary value if inflowing, includes theta
562  void fill_from_realbdy_upwind (const amrex::Vector<amrex::MultiFab*>& mfs,
563  amrex::Real time,
564  bool cons_only,
565  int icomp_cons,
566  int ncomp_cons,
567  amrex::IntVect ngvect_cons,
568  amrex::IntVect ngvect_vels);
569 #endif
570 
571 #ifdef ERF_USE_NETCDF
572  void init_from_wrfinput (int lev,
573  amrex::MultiFab& mf_C1H,
574  amrex::MultiFab& mf_C2H,
575  amrex::MultiFab& mf_MUB,
576  amrex::MultiFab& mf_PSFC);
577  void init_from_metgrid (int lev);
578  void init_from_ncfile (int lev);
579 #endif // ERF_USE_NETCDF
580 
581 #ifdef ERF_USE_WINDFARM
582  void init_windfarm(int lev);
583  void advance_windfarm (const amrex::Geometry& a_geom,
584  const amrex::Real& dt_advance,
585  amrex::MultiFab& cons_in,
586  amrex::MultiFab& U_old,
587  amrex::MultiFab& V_old,
588  amrex::MultiFab& W_old,
589  amrex::MultiFab& mf_vars_windfarm,
590  const amrex::MultiFab& mf_Nturb,
591  const amrex::MultiFab& mf_SMark,
592  const amrex::Real& time);
593 #endif
594 
595  void MakeEBGeometry ();
596  void make_eb_box ();
598 
599  // more flexible version of AverageDown() that lets you average down across multiple levels
600  void AverageDownTo (int crse_lev, int scomp, int ncomp); // NOLINT
601 
602  // write checkpoint file to disk
603  void WriteCheckpointFile () const;
604 
605  // read checkpoint file from disk
606  void ReadCheckpointFile ();
607 
608  // read checkpoint file from disk -- called after instantiating m_SurfaceLayer
610 
611  void init_zphys (int lev, amrex::Real time);
612  void remake_zphys (int lev, amrex::Real time, std::unique_ptr<amrex::MultiFab>& temp_zphys_nd);
613  void update_terrain_arrays (int lev);
614 
615 private:
616 
617  ///////////////////////////
618  // private member functions
619  ///////////////////////////
620 
621  // read in some parameters from inputs file
622  void ReadParameters ();
623  void ParameterSanityChecks ();
624 
625  // set covered coarse cells/faces to be the average of overlying fine cells/faces
626  void AverageDown ();
627 
628  void update_diffusive_arrays (int lev, const amrex::BoxArray& ba, const amrex::DistributionMapping& dm);
629 
630  void Construct_ERFFillPatchers (int lev);
631 
632  void Define_ERFFillPatchers (int lev);
633 
634  void init1DArrays ();
635 
636  void init_bcs ();
637  void init_phys_bcs (bool& rho_read, bool& read_prim_theta);
638 
639  void init_custom (int lev);
640 
641  void init_uniform (int lev);
642 
643  void init_stuff (int lev, const amrex::BoxArray& ba, const amrex::DistributionMapping& dm,
644  amrex::Vector<amrex::MultiFab>& lev_new, amrex::Vector<amrex::MultiFab>& lev_old,
645  amrex::MultiFab& tmp_base_state,
646  std::unique_ptr<amrex::MultiFab>& tmp_zphys_nd);
647 
648  // Initialize the Turbulent perturbation
649  void turbPert_update (const int lev, const amrex::Real dt);
650  void turbPert_amplitude (const int lev);
651 
652  // Initialize the integrator object
653  void initialize_integrator (int lev, amrex::MultiFab& cons_mf, amrex::MultiFab& vel_mf);
654 
655  // Create the physbcs objects
656  void make_physbcs (int lev);
657 
658  // Initialize the microphysics object
659  void initializeMicrophysics (const int&);
660 
661 #ifdef ERF_USE_WINDFARM
662  // Initialize the windfarm object
663  void initializeWindFarm (const int&);
664 #endif
665 
666  // Compute a vector of new MultiFabs by copying from valid region and filling ghost cells
667  //
668  // NOTE: FillPatch takes in an empty MF, and returns cell-centered + velocities (not momenta)
669  //
670  // This one works only at level = 0 (base state does not change)
671  void FillPatchCrseLevel (int lev, amrex::Real time,
672  const amrex::Vector<amrex::MultiFab*>& mfs_vel,
673  bool cons_only=false);
674 
675  // This one works only at level > 0 (base state does change)
676  void FillPatchFineLevel (int lev, amrex::Real time,
677  const amrex::Vector<amrex::MultiFab*>& mfs_vel,
678  const amrex::Vector<amrex::MultiFab*>& mfs_mom,
679  const amrex::MultiFab& old_base_state,
680  const amrex::MultiFab& new_base_state,
681  bool fillset=true, bool cons_only=false);
682 
683  // Compute new multifabs by copying data from valid region and filling ghost cells.
684  // Unlike FillPatch, FillIntermediatePatch will use the supplied multifabs instead of fine level data.
685  // This is to support filling boundary cells at an intermediate time between old/new times
686  // on the fine level when valid data at a specific time is already available (such as
687  // at each RK stage when integrating between initial and final times at a given level).
688  //
689  // NOTE: FillIntermediatePatch takes in updated momenta, and returns both updated velocity and momenta
690  //
691  void FillIntermediatePatch (int lev, amrex::Real time,
692  const amrex::Vector<amrex::MultiFab*>& mfs_vel,
693  const amrex::Vector<amrex::MultiFab*>& mfs_mom,
694  int ng_cons, int ng_vel, bool cons_only, int icomp_cons, int ncomp_cons);
695 
696  // Fill all multifabs (and all components) in a vector of multifabs corresponding to the
697  // grid variables defined in vars_old and vars_new just as FillCoarsePatch.
698  void FillCoarsePatch (int lev, amrex::Real time);
699 
700  // advance a level by dt
701  // includes a recursive call for finer levels
702  void timeStep (int lev, amrex::Real time, int iteration);
703 
704  // advance a single level for a single time step
705  void Advance (int lev, amrex::Real time, amrex::Real dt_lev, int iteration, int ncycle);
706 
707  //! Initialize HSE
708  void initHSE ();
709  void initHSE (int lev);
710 
711  //! Initialize Rayleigh damping profiles
712  void initRayleigh ();
713 
714  //! Initialize sponge profiles
715  void initSponge ();
716 
717  //! Set Rayleigh mean profiles from input sounding
718  void setRayleighRefFromSounding (bool restarting);
719 
720  //! Set sponge mean profiles from input sounding
721  void setSpongeRefFromSounding (bool restarting);
722 
723  // a wrapper for estTimeStep()
724  void ComputeDt (int step = -1);
725 
726  // get plotfile name
727  [[nodiscard]] std::string PlotFileName (int lev) const;
728 
729  // set plotfile variable names
730  static amrex::Vector<std::string> PlotFileVarNames (amrex::Vector<std::string> plot_var_names);
731 
732  // set which variables and derived quantities go into plotfiles
733  void setPlotVariables (const std::string& pp_plot_var_names, amrex::Vector<std::string>& plot_var_names);
734  void setPlotVariables2D (const std::string& pp_plot_var_names, amrex::Vector<std::string>& plot_var_names);
735  // append variables to plot
736  void appendPlotVariables (const std::string& pp_plot_var_names, amrex::Vector<std::string>& plot_var_names);
737 
738  void setSubVolVariables (const std::string& pp_subvol_var_names, amrex::Vector<std::string>& subvol_var_names);
739 
740 #ifdef ERF_USE_NETCDF
741  //! Create 1D vertical column output for coupling
742  void createNCColumnFile (int lev,
743  const std::string& colfile_name, amrex::Real xloc, amrex::Real yloc);
744 
745  // Copy from the NC*fabs into the MultiFabs holding the boundary data
746  void init_from_wrfbdy (amrex::Vector<amrex::FArrayBox*> x_vel_lateral,
747  amrex::Vector<amrex::FArrayBox*> y_vel_lateral,
748  amrex::Vector<amrex::FArrayBox*> z_vel_lateral,
749  amrex::Vector<amrex::FArrayBox*> T_lateral);
750 
751  static amrex::Real start_bdy_time;
752  static amrex::Real start_low_time;
753 
754  static amrex::Real bdy_time_interval;
755  static amrex::Real low_time_interval;
756 
757  // *** *** FArrayBox's for holding the SURFACE data
758  // amrex::IArrayBox NC_IVGTYP_fab; // Vegetation type (IVGTYP); Discrete numbers;
759  // amrex::FArrayBox NC_z0_fab; // Surface Roughness, z0 = z0 (IVGTYP)
760  // amrex::FArrayBox NC_PSFC_fab; // Surface pressure
761 
762  // TODO: Clarify the relation between SST and TSK
763  // amrex::FArrayBox NC_SST_fab; // Sea Surface Temperature; Defined even for land area
764  // amrex::FArrayBox NC_TSK_fab; // Surface Skin Temperature; Appears to be same as SST...
765 
766  // Vectors (over time) of Vector (over variables) of FArrayBoxs for holding the data read from the wrfbdy NetCDF file
767  amrex::Vector<amrex::Vector<amrex::FArrayBox>> bdy_data_xlo;
768  amrex::Vector<amrex::Vector<amrex::FArrayBox>> bdy_data_xhi;
769  amrex::Vector<amrex::Vector<amrex::FArrayBox>> bdy_data_ylo;
770  amrex::Vector<amrex::Vector<amrex::FArrayBox>> bdy_data_yhi;
771 
772  amrex::Vector<amrex::Vector<amrex::FArrayBox>> low_data_zlo;
773 
774  // Maximum value of terrain read in from wrfinput at level 0
775  amrex::Real z_top = 0.;
776 
777 #endif // ERF_USE_NETCDF
778 
779  amrex::Vector<std::unique_ptr<amrex::MultiFab>> lat_m, lon_m; // Latitude and Longitude on the grid
780 
781  amrex::Vector<std::unique_ptr<amrex::MultiFab>> sinPhi_m, cosPhi_m; // Coriolis factors on the grid
782 
783  // Struct for working with the sounding data we take as an input
785 
786  // Struct for working with the sponge data we take as an input
788 
789  // Vector (6 planes) of DeviceVectors (ncell in plane) for Dirichlet BC data
790  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real>> xvel_bc_data;
791  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real>> yvel_bc_data;
792  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real>> zvel_bc_data;
793  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real>> th_bc_data;
794 
795  // Function to read and populate above vectors (if input file exists)
796  void init_Dirichlet_bc_data (const std::string input_file);
797 
798  // Read the file passed to amr.restart and use it as an initial condition for
799  // the current simulation. Supports a different number of components and
800  // ghost cells.
802 
803  // Initialize the new-time data at a level from the initial_data MultiFab
804  void InitializeLevelFromData (int lev, const amrex::MultiFab& initial_data);
805 
806  // utility to skip to next line in Header
807  static void GotoNextLine (std::istream& is);
808 
809  // Single level functions called by advance()
810  void post_update (amrex::MultiFab& state_mf, amrex::Real time, const amrex::Geometry& geom);
811  void fill_rhs (amrex::MultiFab& rhs_mf, const amrex::MultiFab& state_mf, amrex::Real time, const amrex::Geometry& geom);
812 
813  ////////////////
814  // private data members
815 
816  std::unique_ptr<ProblemBase> prob = nullptr; // problem-specific functions
817 
818  amrex::Vector<int> num_boxes_at_level; // how many boxes specified at each level by tagging criteria
819  amrex::Vector<int> num_files_at_level; // how many wrfinput files specified at each level
820  amrex::Vector<amrex::Vector<amrex::Box>> boxes_at_level; // the boxes specified at each level by tagging criteria
821 
822  amrex::Vector<int> istep; // which step?
823  amrex::Vector<int> nsubsteps; // how many substeps on each level?
824 
825  // keep track of old time, new time, and time step at each level
826  amrex::Vector<amrex::Real> t_new;
827  amrex::Vector<amrex::Real> t_old;
828  amrex::Vector<amrex::Real> dt;
829  amrex::Vector<long> dt_mri_ratio;
830 
831 #ifndef ERF_USE_MULTIBLOCK
832  // Array of multifabs to store the solution at each level of refinement;
833  // after advancing a level we use "swap".
834  amrex::Vector<amrex::Vector<amrex::MultiFab> > vars_new;
835  amrex::Vector<amrex::Vector<amrex::MultiFab> > vars_old;
836 
837  // Pressure gradient at each level
838  amrex::Vector<amrex::Vector<amrex::MultiFab> > gradp;
839 
840  // Velocity time averaged field
841  amrex::Vector<std::unique_ptr<amrex::MultiFab>> vel_t_avg;
842  amrex::Vector<amrex::Real> t_avg_cnt;
843 #endif
844  amrex::Vector<std::unique_ptr<MRISplitIntegrator<amrex::Vector<amrex::MultiFab> > > > mri_integrator_mem;
845 
846  // Used for anelastic or when we impose an initial projection
847  amrex::Vector<amrex::MultiFab> pp_inc;
848 
849  // Used only for fast substepping
850  amrex::Vector<amrex::MultiFab> lagged_delta_rt;
851  amrex::Vector<amrex::MultiFab> avg_xmom;
852  amrex::Vector<amrex::MultiFab> avg_ymom;
853  amrex::Vector<amrex::MultiFab> avg_zmom;
854 
855  // Vector over levels of routines to impose physical boundary conditions
856  amrex::Vector<std::unique_ptr<ERFPhysBCFunct_cons>> physbcs_cons;
857  amrex::Vector<std::unique_ptr<ERFPhysBCFunct_u>> physbcs_u;
858  amrex::Vector<std::unique_ptr<ERFPhysBCFunct_v>> physbcs_v;
859  amrex::Vector<std::unique_ptr<ERFPhysBCFunct_w>> physbcs_w;
860  amrex::Vector<std::unique_ptr<ERFPhysBCFunct_base>> physbcs_base;
861 
862  // Store primitive variable for MOST BC
863  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Theta_prim;
864  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Qv_prim;
865  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Qr_prim;
866 
867  // Scratch space for time integrator
868  amrex::Vector<amrex::MultiFab> rU_old;
869  amrex::Vector<amrex::MultiFab> rU_new;
870  amrex::Vector<amrex::MultiFab> rV_old;
871  amrex::Vector<amrex::MultiFab> rV_new;
872  amrex::Vector<amrex::MultiFab> rW_old;
873  amrex::Vector<amrex::MultiFab> rW_new;
874 
875  // amrex::Vector<amrex::MultiFab> xmom_crse_rhs;
876  // amrex::Vector<amrex::MultiFab> ymom_crse_rhs;
877  amrex::Vector<amrex::MultiFab> zmom_crse_rhs;
878 
879  std::unique_ptr<Microphysics> micro;
880  amrex::Vector<amrex::Vector<amrex::MultiFab*>> qmoist; // (lev,ncomp) rain_accum, snow_accum, graup_accum
881 
882  // Variables for wind farm parametrization models
883 
884 #ifdef ERF_USE_WINDFARM
885  std::unique_ptr<WindFarm> windfarm;
886  amrex::Vector<amrex::MultiFab> Nturb;
887  amrex::Vector<amrex::MultiFab> vars_windfarm; // Fitch: Vabs, Vabsdt, dudt, dvdt, dTKEdt
888  // EWP: dudt, dvdt, dTKEdt
889 
890  amrex::Vector<amrex::MultiFab> SMark; // A multifab that holds an integer corresponding
891  // to the number of the wind turbine to sample
892  // velocity upstream of the turbine
893 #endif
894 
896  amrex::Vector<std::string> lsm_data_name; // (ncomp)
897  amrex::Vector<amrex::Vector<amrex::MultiFab*>> lsm_data; // (lev,ncomp)
898  amrex::Vector<std::string> lsm_flux_name; // (ncomp)
899  amrex::Vector<amrex::Vector<amrex::MultiFab*>> lsm_flux; // (lev,ncomp)
900 
901  amrex::Vector<std::unique_ptr<IRadiation>> rad; // Radiation model at each level
902  amrex::Vector<std::unique_ptr<amrex::MultiFab>> qheating_rates; // radiation heating rate source terms (SW, LW)
903  amrex::Vector<std::unique_ptr<amrex::MultiFab>> rad_fluxes; // radiation fluxes (SW up/dn, LW up/dn)
904 #ifdef ERF_USE_SHOC
905  amrex::Vector<std::unique_ptr<SHOCInterface>> shoc_interface; // SHOC model at each level
906 #endif
907 
908 #ifdef ERF_USE_P3
909  amrex::Vector<std::unique_ptr<P3Interface>> p3_interface; // P3 model at each level
910 #endif
911 
912  // Containers for additional SLM inputs
913  amrex::Vector<std::unique_ptr<amrex::MultiFab>> sw_lw_fluxes; // Direct SW (visible, NIR), Diffuse SW (visible, NIR), LW flux
914  amrex::Vector<std::unique_ptr<amrex::MultiFab>> solar_zenith; // Solar zenith angle
915 
916  bool plot_rad = false;
917  int rad_datalog_int = -1;
918 
919  // Fillpatcher classes for coarse-fine boundaries
920  int cf_width{0};
921  int cf_set_width{0};
922  amrex::Vector<ERFFillPatcher> FPr_c;
923  amrex::Vector<ERFFillPatcher> FPr_u;
924  amrex::Vector<ERFFillPatcher> FPr_v;
925  amrex::Vector<ERFFillPatcher> FPr_w;
926 
927  // Diffusive stresses and Smag
928  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab>>> Tau;
929  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab>>> Tau_corr; // correction to tau for implicit
930  amrex::Vector<std::unique_ptr<amrex::MultiFab>> eddyDiffs_lev;
931  amrex::Vector<std::unique_ptr<amrex::MultiFab>> SmnSmn_lev;
932 
933  // Sea Surface Temps, Skin Temperature and Land Masks (lev, ntimes)
934  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab>>> sst_lev;
935  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab>>> tsk_lev;
936  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::iMultiFab>>> lmask_lev;
937 
938  // Land grid types and urban fraction
939  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::iMultiFab>>> land_type_lev;
940  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::iMultiFab>>> soil_type_lev;
941  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab>>> urb_frac_lev;
942 
943  // Other SFS terms
944  amrex::Vector<std::unique_ptr<amrex::MultiFab>> SFS_hfx1_lev, SFS_hfx2_lev, SFS_hfx3_lev;
945  amrex::Vector<std::unique_ptr<amrex::MultiFab>> SFS_diss_lev;
946  amrex::Vector<std::unique_ptr<amrex::MultiFab>> SFS_q1fx1_lev, SFS_q1fx2_lev, SFS_q1fx3_lev;
947  amrex::Vector<std::unique_ptr<amrex::MultiFab>> SFS_q2fx3_lev;
948 
949  // Grid stretching
950  amrex::Vector<amrex::Vector<amrex::Real>> zlevels_stag; // nominal height levels
951 
952  // Data structures for terrain-fitted coordinates
953  amrex::Vector<std::unique_ptr<amrex::MultiFab>> z_phys_nd;
954  amrex::Vector<std::unique_ptr<amrex::MultiFab>> z_phys_cc;
955 
956  amrex::Vector<std::unique_ptr<amrex::MultiFab>> detJ_cc;
957  amrex::Vector<std::unique_ptr<amrex::MultiFab>> ax;
958  amrex::Vector<std::unique_ptr<amrex::MultiFab>> ay;
959  amrex::Vector<std::unique_ptr<amrex::MultiFab>> az;
960 
961  amrex::Vector<std::unique_ptr<amrex::MultiFab>> z_phys_nd_src;
962  amrex::Vector<std::unique_ptr<amrex::MultiFab>> z_phys_cc_src;
963  amrex::Vector<std::unique_ptr<amrex::MultiFab>> detJ_cc_src;
964  amrex::Vector<std::unique_ptr<amrex::MultiFab>> ax_src;
965  amrex::Vector<std::unique_ptr<amrex::MultiFab>> ay_src;
966  amrex::Vector<std::unique_ptr<amrex::MultiFab>> az_src;
967 
968  amrex::Vector<std::unique_ptr<amrex::MultiFab>> z_phys_nd_new;
969  amrex::Vector<std::unique_ptr<amrex::MultiFab>> detJ_cc_new;
970 
971  amrex::Vector<std::unique_ptr<amrex::MultiFab>> z_t_rk;
972 
973  // Data structures for terrain immersed forcing
974  amrex::Vector<std::unique_ptr<amrex::MultiFab>> terrain_blanking;
975 
976  // Wall distance function
977  amrex::Vector<std::unique_ptr<amrex::MultiFab>> walldist;
978 
979  // Map scale factors -- vector across levels of vector across types
980  amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab>>> mapfac;
981 
982  // Mask to be used to identify cells covered by finer level (if level > 0)
983  amrex::Vector<std::unique_ptr<amrex::MultiFab>> fine_mask;
984 
985  amrex::Vector<amrex::Vector<amrex::Real>> stretched_dz_h;
986  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real>> stretched_dz_d;
987 
988  amrex::Vector<amrex::MultiFab> base_state;
989  amrex::Vector<amrex::MultiFab> base_state_new;
990 
991  // Wave coupling data
992  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Hwave;
993  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Lwave;
994  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Hwave_onegrid;
995  amrex::Vector<std::unique_ptr<amrex::MultiFab>> Lwave_onegrid;
996  bool finished_wave = false;
997 
998  // array of flux registers
999  amrex::Vector<amrex::YAFluxRegister*> advflux_reg;
1000 
1001  // A BCRec is essentially a 2*DIM integer array storing the boundary
1002  // condition type at each lo/hi walls in each direction. We have one BCRec
1003  // for each component of the cell-centered variables and each velocity component.
1004  amrex::Vector <amrex::BCRec> domain_bcs_type;
1005  amrex::Gpu::DeviceVector<amrex::BCRec> domain_bcs_type_d;
1006 
1007  // We store these so that we can print them out in the job_info file
1008  amrex::Array<std::string,2*AMREX_SPACEDIM> domain_bc_type;
1009 
1010  // These hold the Dirichlet values at walls which need them ...
1011  amrex::Array<amrex::Array<amrex::Real, AMREX_SPACEDIM*2>, AMREX_SPACEDIM+NBCVAR_max> m_bc_extdir_vals;
1012 
1013  // These hold the Neumann values at walls which need them ...
1014  amrex::Array<amrex::Array<amrex::Real, AMREX_SPACEDIM*2>, AMREX_SPACEDIM+NBCVAR_max> m_bc_neumann_vals;
1015 
1016  // These are the "physical" boundary condition types (e.g. "inflow")
1017  amrex::GpuArray<ERF_BC, AMREX_SPACEDIM*2> phys_bc_type;
1018 
1019  // These are the masks for the thin immersed body
1020  amrex::Vector<std::unique_ptr<amrex::iMultiFab>> xflux_imask;
1021  amrex::Vector<std::unique_ptr<amrex::iMultiFab>> yflux_imask;
1022  amrex::Vector<std::unique_ptr<amrex::iMultiFab>> zflux_imask;
1023  //amrex::Vector<std::unique_ptr<amrex::iMultiFab>> overset_imask;
1024 
1025  // These are the body forces that result from the thin immersed body
1026  amrex::Vector<std::unique_ptr<amrex::MultiFab>> thin_xforce;
1027  amrex::Vector<std::unique_ptr<amrex::MultiFab>> thin_yforce;
1028  amrex::Vector<std::unique_ptr<amrex::MultiFab>> thin_zforce;
1029 
1035 
1041 
1042  amrex::Vector<int> last_subvol_step;
1043  amrex::Vector<amrex::Real> last_subvol_time;
1044 
1046 
1047  // Output control
1048  const int datwidth = 14;
1049  const int datprecision = 6;
1050  const int timeprecision = 13; // e.g., 1-yr LES: 31,536,000 s with dt ~ 0.01 ==> min prec = 10
1051  // epoch time for 2030-01-01 00:00:00 is 1,893,456,000 s with dt ~ 0.01 ==> min prec = 12
1052 
1053  ////////////////
1054  // runtime parameters
1055 
1056  // Maximum number of steps
1057  int max_step = -1;
1058 
1059  // Start and stop times
1062 
1063  bool use_datetime = false;
1064  const std::string datetime_format = "%Y-%m-%d %H:%M:%S"; // ISO 8601 standard
1065 
1066  // if >= 0 we restart from a checkpoint
1067  std::string restart_chkfile = "";
1068 
1069  // Time step controls
1076 
1077  // Fixed dt at each level (only used if positive)
1078  amrex::Vector<amrex::Real> fixed_dt;
1079  amrex::Vector<amrex::Real> fixed_fast_dt;
1081 
1082  // how often each level regrids the higher levels of refinement
1083  // (after a level advances that many time steps)
1084  int regrid_int = -1;
1085 
1086  // Should we regrid level 0 after restart if Ngrids < Nprocs or
1087  // max_grid_size has changed?
1089 
1090  // plotfile prefix and frequency
1091  std::string plot3d_file_1 {"plt_1_"};
1092  std::string plot3d_file_2 {"plt_2_"};
1093  std::string plot2d_file_1 {"plt2d_1_"};
1094  std::string plot2d_file_2 {"plt2d_2_"};
1095  std::string subvol_file {"subvol"};
1097  int m_plot3d_int_1 = -1;
1098  int m_plot3d_int_2 = -1;
1099  int m_plot2d_int_1 = -1;
1100  int m_plot2d_int_2 = -1;
1101 
1102  amrex::Vector<int> m_subvol_int;
1103  amrex::Vector<amrex::Real> m_subvol_per;
1104 
1109  bool m_plot_face_vels = false;
1110 
1111  bool plot_lsm = false;
1112 
1113  // other sampling output control
1114  int profile_int = -1;
1115  bool destag_profiles = true;
1116 
1117  // Checkpoint type, prefix and frequency
1118  std::string check_file {"chk"};
1119  int m_check_int = -1;
1121 
1122  amrex::Vector<std::string> subvol3d_var_names;
1123 
1124  amrex::Vector<std::string> plot3d_var_names_1;
1125  amrex::Vector<std::string> plot3d_var_names_2;
1126  amrex::Vector<std::string> plot2d_var_names_1;
1127  amrex::Vector<std::string> plot2d_var_names_2;
1128  const amrex::Vector<std::string> cons_names {"density", "rhotheta", "rhoKE", "rhoadv_0",
1129  "rhoQ1", "rhoQ2", "rhoQ3",
1130  "rhoQ4", "rhoQ5", "rhoQ6"};
1131 
1132  // **************************************************************************************
1133  // NOTE: The order of variable names here **MUST MATCH THE ORDER** in IO/ERF_Plotfile.cpp
1134  // **************************************************************************************
1135  const amrex::Vector<std::string> derived_names {"soundspeed", "temp", "theta", "KE", "scalar",
1136  "vorticity_x","vorticity_y","vorticity_z",
1137  "magvel", "divU",
1138  "pres_hse", "dens_hse", "theta_hse", "pressure", "pert_pres", "pert_dens",
1139  "eq_pot_temp",
1140 #ifdef ERF_USE_WINDFARM
1141  "num_turb", "SMark0", "SMark1",
1142 #endif
1143  "dpdx", "dpdy", "dpdz", "pres_hse_x", "pres_hse_y",
1144  "z_phys", "detJ", "mapfac", "lat_m", "lon_m",
1145  // Time averaged velocity
1146  "u_t_avg", "v_t_avg", "w_t_avg", "umag_t_avg",
1147  // eddy viscosity
1148  "nut",
1149  // eddy diffusivity of momentum
1150  "Kmv","Kmh",
1151  // eddy diffusivity of heat
1152  "Khv","Khh",
1153  // turbulence lengthscale
1154  "Lturb",
1155  // wall distance
1156  "walldist",
1157  // dissipation
1158  "diss",
1159  // moisture vars
1160  "moist_density", "qv", "qc", "qi", "qrain", "qsnow", "qgraup", "qt", "qn", "qp", "qsat",
1161  "rain_accum", "snow_accum", "graup_accum",
1162  "rel_humidity", "condensation_rate",
1163  "reflectivity", "max_reflectivity",
1164  // Terrain IB mask
1165  "terrain_IB_mask",
1166  // EB variables
1167  "volfrac"
1168 #ifdef ERF_COMPUTE_ERROR
1169  // error vars
1170  ,"xvel_err", "yvel_err", "zvel_err", "pp_err"
1171 #endif
1172  ,"qsrc_sw", "qsrc_lw"
1173  };
1174 
1175  // **************************************************************************************
1176  // NOTE: The order of variable names here **MUST MATCH THE ORDER** in IO/ERF_Plotfile.cpp
1177  // **************************************************************************************
1178  const amrex::Vector<std::string> derived_names_2d {
1179  "z_surf", "landmask", "mapfac", "lat_m", "lon_m",
1180  "u_star", "w_star", "t_star", "q_star", "Olen", "pblh",
1181  "t_surf", "q_surf", "z0", "OLR", "sens_flux", "laten_flux",
1182  "surf_pres", "integrated_qv"
1183  };
1184 
1185  // **************************************************************************************
1186  // NOTE: The order of variable names here **MUST MATCH THE ORDER** in IO/ERF_WriteSubVolume.cpp
1187  // **************************************************************************************
1188  const amrex::Vector<std::string> derived_subvol_names {"soundspeed", "temp", "theta", "KE", "scalar"};
1189 
1190  // algorithm choices
1192 
1193  // Turbulent perturbation structure
1195 
1196 #ifdef ERF_USE_PARTICLES
1197  // Particle container with all particle species
1198  ParticleData particleData;
1199 
1200  // variables and functions for tracers particles
1201  bool m_use_tracer_particles; /*!< tracer particles that advect with flow with optional sedimentation */
1202 
1203  /*! Read tracer and hydro particles parameters */
1204  void readTracersParams();
1205 
1206  /*! Initialize tracer and hydro particles */
1207  void initializeTracers ( amrex::ParGDBBase*,
1208  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>&,
1209  const amrex::Real time);
1210 
1211  /*! Restart tracer and hydro particles */
1212  void restartTracers ( amrex::ParGDBBase*, const std::string& );
1213 
1214  /*! Evolve tracers and hydro particles */
1215  void evolveTracers( int,
1216  amrex::Real,
1217  amrex::Vector<amrex::Vector<amrex::MultiFab>>&,
1218  const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& );
1219 
1220 #endif
1221 
1222 #ifdef ERF_USE_MULTIBLOCK
1223  MultiBlockContainer *m_mbc = nullptr;
1224 #endif
1225 
1226  static int verbose;
1227  static int mg_verbose;
1228  static bool use_fft;
1229 
1230  static int check_for_nans;
1231 
1232  // Diagnostic output interval
1233  static int sum_interval;
1234  static int pert_interval;
1236 
1237  // Write in native AMReX or NetCDF format for each plotfile
1238  static PlotFileType plotfile3d_type_1;
1239  static PlotFileType plotfile3d_type_2;
1240  static PlotFileType plotfile2d_type_1;
1241  static PlotFileType plotfile2d_type_2;
1242 
1245 
1246  static StateInterpType interpolation_type;
1247 
1248  // sponge_type: "input_sponge"
1249  static std::string sponge_type;
1250 
1251  // NetCDF initialization (wrfinput) file
1252  static amrex::Vector<amrex::Vector<std::string>> nc_init_file;
1253  static amrex::Vector<amrex::Vector<int>> have_read_nc_init_file;
1254 
1255  // NetCDF initialization (wrfbdy/met_em) file
1256  static std::string nc_bdy_file;
1257  int real_width{0};
1259  bool real_extrap_w{true};
1260 
1261  // NetCDF initialization of zlo boundary values
1262  static std::string nc_low_file;
1263 
1264  // Options for vertical interpolation of met_em*.nc data.
1267  bool metgrid_debug_dry{false};
1268  bool metgrid_debug_psfc{false};
1269  bool metgrid_debug_msf{false};
1273  bool metgrid_use_sfc{true};
1274  bool metgrid_retain_sfc{false};
1278 
1279  amrex::Vector<amrex::BoxArray> ba1d;
1280  amrex::Vector<amrex::BoxArray> ba2d;
1281  std::unique_ptr<amrex::MultiFab> mf_C1H;
1282  std::unique_ptr<amrex::MultiFab> mf_C2H;
1283  std::unique_ptr<amrex::MultiFab> mf_MUB;
1284 
1285  amrex::Vector<std::unique_ptr<amrex::MultiFab>> mf_PSFC;
1286 
1287  // 1D CDF output (for ingestion in AMR-Wind)
1288  static int output_1d_column;
1289  static int column_interval;
1293  static std::string column_file_name;
1294 
1295  // 2D BndryRegister output (for ingestion in AMR-Wind)
1300 
1301  // 2D BndryRegister input
1303 
1304  static int ng_dens_hse;
1305  static int ng_pres_hse;
1306 
1307  // Custom source terms
1308  amrex::Vector<std::unique_ptr<amrex::MultiFab>> rhotheta_src;
1309  amrex::Vector<std::unique_ptr<amrex::MultiFab>> rhoqt_src;
1310 
1311  amrex::Vector< amrex::Vector<amrex::Real> > h_w_subsid;
1312  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real> > d_w_subsid;
1313 
1314  amrex::Vector< amrex::Vector<amrex::Real> > h_u_geos;
1315  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real> > d_u_geos;
1316 
1317  amrex::Vector< amrex::Vector<amrex::Real> > h_v_geos;
1318  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real> > d_v_geos;
1319 
1320  // Function to read and populate above host vectors (if input file exists)
1321  void init_geo_wind_profile (const std::string input_file,
1322  amrex::Vector<amrex::Real>& u_geos,
1323  amrex::Gpu::DeviceVector<amrex::Real>& u_geos_d,
1324  amrex::Vector<amrex::Real>& v_geos,
1325  amrex::Gpu::DeviceVector<amrex::Real>& v_geos_d,
1326  const amrex::Geometry& lgeom,
1327  const amrex::Vector<amrex::Real>& zlev_stag);
1328 
1329  // This is a vector over levels of vectors across quantities of Vectors
1330  amrex::Vector<amrex::Vector<amrex::Vector<amrex::Real> > > h_rayleigh_ptrs;
1331  amrex::Vector<amrex::Vector<amrex::Vector<amrex::Real> > > h_sponge_ptrs;
1332 
1333  // These are vectors over levels of Vectors
1334  amrex::Vector<amrex::Vector<amrex::Real> > h_sinesq_ptrs;
1335  amrex::Vector<amrex::Vector<amrex::Real> > h_sinesq_stag_ptrs;
1336 
1337  // This is a vector over levels of vectors across quantities of DeviceVectors
1338  amrex::Vector<amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real> > > d_rayleigh_ptrs;
1339  amrex::Vector<amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real> > > d_sponge_ptrs;
1340 
1341  // These are vectors over levels of DeviceVectors
1342  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real> > d_sinesq_ptrs;
1343  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real> > d_sinesq_stag_ptrs;
1344 
1345  amrex::Vector<amrex::Real> h_havg_density;
1346  amrex::Vector<amrex::Real> h_havg_temperature;
1347  amrex::Vector<amrex::Real> h_havg_pressure;
1348  amrex::Vector<amrex::Real> h_havg_qv;
1349  amrex::Vector<amrex::Real> h_havg_qc;
1350 
1351  amrex::Gpu::DeviceVector<amrex::Real> d_havg_density;
1352  amrex::Gpu::DeviceVector<amrex::Real> d_havg_temperature;
1353  amrex::Gpu::DeviceVector<amrex::Real> d_havg_pressure;
1354  amrex::Gpu::DeviceVector<amrex::Real> d_havg_qv;
1355  amrex::Gpu::DeviceVector<amrex::Real> d_havg_qc;
1356 
1357  void refinement_criteria_setup ();
1358 
1359  std::unique_ptr<WriteBndryPlanes> m_w2d = nullptr;
1360  std::unique_ptr<ReadBndryPlanes> m_r2d = nullptr;
1361  std::unique_ptr<SurfaceLayer> m_SurfaceLayer = nullptr;
1362  amrex::Vector<std::unique_ptr<ForestDrag>> m_forest_drag;
1363 
1364  //
1365  // Holds info for dynamically generated tagging criteria
1366  //
1367  static amrex::Vector<amrex::AMRErrorTag> ref_tags;
1368 
1369  amrex::Vector<amrex::Vector<amrex::BoxArray>> subdomains;
1370 
1371  amrex::Vector<amrex::Real> dz_min;
1372 
1373  static AMREX_FORCE_INLINE
1374  int
1376  {
1377  int ngrow = 0;
1378 
1379  if (sc.use_num_diff)
1380  {
1381  ngrow = 3;
1382  } else {
1383  if (
1390  { ngrow = 3; }
1391  else if (
1398  { ngrow = 3; }
1399  else if (
1408  { ngrow = 3; }
1409  else if (
1418  { ngrow = 4; }
1419  else
1420  {
1421  if (sc.terrain_type == TerrainType::EB){
1422  ngrow = 4;
1423  } else {
1424  ngrow = 2;
1425  }
1426  }
1427  }
1428 
1429  return ngrow;
1430  }
1431 
1432  AMREX_FORCE_INLINE
1433  amrex::YAFluxRegister* getAdvFluxReg (int lev)
1434  {
1435  return advflux_reg[lev];
1436  }
1437 
1438  AMREX_FORCE_INLINE
1439  std::ostream&
1440  DataLog (int i)
1441  {
1442  return *datalog[i];
1443  }
1444 
1445  AMREX_FORCE_INLINE
1446  std::ostream&
1447  DerDataLog (int i)
1448  {
1449  return *der_datalog[i];
1450  }
1451 
1452  AMREX_FORCE_INLINE
1453  int
1454  NumDataLogs () noexcept
1455  {
1456  return datalog.size();
1457  }
1458 
1459  AMREX_FORCE_INLINE
1460  int
1461  NumDerDataLogs () noexcept
1462  {
1463  return der_datalog.size();
1464  }
1465 
1466 
1467  AMREX_FORCE_INLINE
1468  std::ostream&
1470  {
1471  return *sampleptlog[i];
1472  }
1473 
1474  AMREX_FORCE_INLINE
1475  int
1477  {
1478  return sampleptlog.size();
1479  }
1480 
1481  AMREX_FORCE_INLINE
1482  std::ostream&
1484  {
1485  return *samplelinelog[i];
1486  }
1487 
1488  AMREX_FORCE_INLINE
1489  int
1490  NumSampleLineLogs () noexcept
1491  {
1492  return samplelinelog.size();
1493  }
1494 
1495  amrex::IntVect&
1496  SamplePoint (int i)
1497  {
1498  return samplepoint[i];
1499  }
1500 
1501  AMREX_FORCE_INLINE
1502  int
1503  NumSamplePoints () noexcept
1504  {
1505  return samplepoint.size();
1506  }
1507 
1508  amrex::IntVect&
1509  SampleLine (int i)
1510  {
1511  return sampleline[i];
1512  }
1513 
1514  AMREX_FORCE_INLINE
1515  int
1516  NumSampleLines () noexcept
1517  {
1518  return sampleline.size();
1519  }
1520 
1523 
1524  static amrex::Real
1526  {
1527  int numCores = amrex::ParallelDescriptor::NProcs();
1528 #ifdef _OPENMP
1529  numCores = numCores * omp_get_max_threads();
1530 #endif
1531 
1532  amrex::Real T =
1533  numCores * (amrex::ParallelDescriptor::second() - startCPUTime) +
1535 
1536  return T;
1537  }
1538 
1539  void setRecordDataInfo (int i, const std::string& filename) // NOLINT
1540  {
1541  if (amrex::ParallelDescriptor::IOProcessor())
1542  {
1543  datalog[i] = std::make_unique<std::fstream>();
1544  datalog[i]->open(filename.c_str(),std::ios::out|std::ios::app);
1545  if (!datalog[i]->good()) {
1546  amrex::FileOpenFailed(filename);
1547  }
1548  }
1549  amrex::ParallelDescriptor::Barrier("ERF::setRecordDataInfo");
1550  }
1551 
1552  void setRecordDerDataInfo (int i, const std::string& filename) // NOLINT
1553  {
1554  if (amrex::ParallelDescriptor::IOProcessor())
1555  {
1556  der_datalog[i] = std::make_unique<std::fstream>();
1557  der_datalog[i]->open(filename.c_str(),std::ios::out|std::ios::app);
1558  if (!der_datalog[i]->good()) {
1559  amrex::FileOpenFailed(filename);
1560  }
1561  }
1562  amrex::ParallelDescriptor::Barrier("ERF::setRecordDerDataInfo");
1563  }
1564 
1565  void setRecordEnergyDataInfo (int i, const std::string& filename) // NOLINT
1566  {
1567  if (amrex::ParallelDescriptor::IOProcessor())
1568  {
1569  tot_e_datalog[i] = std::make_unique<std::fstream>();
1570  tot_e_datalog[i]->open(filename.c_str(),std::ios::out|std::ios::app);
1571  if (!tot_e_datalog[i]->good()) {
1572  amrex::FileOpenFailed(filename);
1573  }
1574  }
1575  amrex::ParallelDescriptor::Barrier("ERF::setRecordEnergyDataInfo");
1576  }
1577 
1578  void setRecordSamplePointInfo (int i, int lev, amrex::IntVect& cell, const std::string& filename) // NOLINT
1579  {
1580  amrex::MultiFab dummy(grids[lev],dmap[lev],1,0);
1581  for (amrex::MFIter mfi(dummy); mfi.isValid(); ++mfi)
1582  {
1583  const amrex::Box& bx = mfi.validbox();
1584  if (bx.contains(cell)) {
1585  sampleptlog[i] = std::make_unique<std::fstream>();
1586  sampleptlog[i]->open(filename.c_str(),std::ios::out|std::ios::app);
1587  if (!sampleptlog[i]->good()) {
1588  amrex::FileOpenFailed(filename);
1589  }
1590  }
1591  }
1592  amrex::ParallelDescriptor::Barrier("ERF::setRecordSamplePointInfo");
1593  }
1594 
1595  void setRecordSampleLineInfo (int i, int lev, amrex::IntVect& cell, const std::string& filename) // NOLINT
1596  {
1597  amrex::MultiFab dummy(grids[lev],dmap[lev],1,0);
1598  for (amrex::MFIter mfi(dummy); mfi.isValid(); ++mfi)
1599  {
1600  const amrex::Box& bx = mfi.validbox();
1601  if (bx.contains(cell)) {
1602  samplelinelog[i] = std::make_unique<std::fstream>();
1603  samplelinelog[i]->open(filename.c_str(),std::ios::out|std::ios::app);
1604  if (!samplelinelog[i]->good()) {
1605  amrex::FileOpenFailed(filename);
1606  }
1607  }
1608  }
1609  amrex::ParallelDescriptor::Barrier("ERF::setRecordSampleLineInfo");
1610  }
1611 
1612  // Data sampler for line and plane output
1617  std::unique_ptr<LineSampler> line_sampler = nullptr;
1618  std::unique_ptr<PlaneSampler> plane_sampler = nullptr;
1619 
1620  amrex::Vector<std::unique_ptr<std::fstream> > datalog;
1621  amrex::Vector<std::unique_ptr<std::fstream> > der_datalog;
1622  amrex::Vector<std::unique_ptr<std::fstream> > tot_e_datalog;
1623  amrex::Vector<std::string> datalogname;
1624  amrex::Vector<std::string> der_datalogname;
1625  amrex::Vector<std::string> tot_e_datalogname;
1626 
1627  amrex::Vector<std::unique_ptr<std::fstream> > sampleptlog;
1628  amrex::Vector<std::string> sampleptlogname;
1629  amrex::Vector<amrex::IntVect> samplepoint;
1630 
1631  amrex::Vector<std::unique_ptr<std::fstream> > samplelinelog;
1632  amrex::Vector<std::string> samplelinelogname;
1633  amrex::Vector<amrex::IntVect> sampleline;
1634 
1635  //! The filename of the ith datalog file.
1636  [[nodiscard]] std::string DataLogName (int i) const noexcept { return datalogname[i]; }
1637  [[nodiscard]] std::string DerDataLogName (int i) const noexcept { return der_datalogname[i]; }
1638 
1639  //! The filename of the ith sampleptlog file.
1640  [[nodiscard]] std::string SamplePointLogName (int i) const noexcept { return sampleptlogname[i]; }
1641 
1642  //! The filename of the ith samplelinelog file.
1643  [[nodiscard]] std::string SampleLineLogName (int i) const noexcept { return samplelinelogname[i]; }
1644 
1645  // array of EB objects
1646  amrex::Vector<std::unique_ptr<eb_>> eb;
1647 
1648  [[nodiscard]] eb_ const& get_eb (int lev) const noexcept {
1649  AMREX_ASSERT(lev >= 0 && lev < eb.size() && eb[lev] != nullptr);
1650  return *eb[lev];
1651  }
1652 
1653  [[nodiscard]] amrex::EBFArrayBoxFactory const&
1654  EBFactory (int lev) const noexcept {
1655  return *(eb[lev]->get_const_factory());
1656  }
1657 
1658  [[nodiscard]] static int nghost_eb_basic ()
1659  { return 5; }
1660 
1661  // We need 5 for doing StateRedistribution; otherwise 4 would be enough
1662  [[nodiscard]] static int nghost_eb_volume ()
1663  { return 5; }
1664 
1665  [[nodiscard]] static int nghost_eb_full ()
1666  { return 4; }
1667 
1668 #ifdef ERF_USE_FFT
1669  amrex::Vector<std::unique_ptr<amrex::FFT::Poisson<amrex::MultiFab>>> m_3D_poisson;
1670  amrex::Vector<std::unique_ptr<amrex::FFT::PoissonHybrid<amrex::MultiFab>>> m_2D_poisson;
1671 #endif
1672 
1673 public:
1674  void writeJobInfo (const std::string& dir) const;
1675  static void writeBuildInfo (std::ostream& os);
1676 
1677  static void print_banner(MPI_Comm /*comm*/, std::ostream& /*out*/);
1678  static void print_usage(MPI_Comm /*comm*/, std::ostream& /*out*/);
1679  static void print_error(MPI_Comm /*comm*/, const std::string& msg);
1680  static void print_summary(std::ostream&);
1681  static void print_tpls(std::ostream& /*out*/);
1682 };
1683 
1684 #endif
AMREX_ENUM(StateInterpType, FullState, Perturbational)
@ tau23
Definition: ERF_DataStruct.H:32
@ tau13
Definition: ERF_DataStruct.H:32
Contains the Eulerian microphysics class.
struct @19 out
#define NBCVAR_max
Definition: ERF_IndexDefines.H:29
@ Centered_6th
Contains the Lagrangian microphysics class.
amrex::Real Real
Definition: ERF_ShocInterface.H:19
Definition: ERF.H:111
void SurfaceDataInterpolation(const int nlevs, const amrex::Real time, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &z_phys_nd, bool regrid_forces_file_read)
Definition: ERF_SurfaceDataInterpolation.cpp:142
void MakeHorizontalAverages()
Definition: ERF.cpp:2792
amrex::Vector< amrex::MultiFab > rU_new
Definition: ERF.H:869
static int last_check_file_step
Definition: ERF.H:1034
amrex::Vector< std::unique_ptr< amrex::MultiFab > > walldist
Definition: ERF.H:977
amrex::Vector< amrex::Vector< amrex::Vector< amrex::Real > > > h_sponge_ptrs
Definition: ERF.H:1331
void initRayleigh()
Initialize Rayleigh damping profiles.
Definition: ERF_InitRayleigh.cpp:14
static amrex::Real start_time
Definition: ERF.H:1060
bool metgrid_basic_linear
Definition: ERF.H:1271
amrex::Vector< amrex::Vector< amrex::Box > > boxes_at_level
Definition: ERF.H:820
amrex::Vector< std::string > samplelinelogname
Definition: ERF.H:1632
void MakeEBGeometry()
int max_step
Definition: ERF.H:1057
bool metgrid_debug_msf
Definition: ERF.H:1269
AMREX_FORCE_INLINE std::ostream & DerDataLog(int i)
Definition: ERF.H:1447
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > mapfac
Definition: ERF.H:980
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:94
amrex::Vector< std::unique_ptr< MRISplitIntegrator< amrex::Vector< amrex::MultiFab > > > > mri_integrator_mem
Definition: ERF.H:844
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > d_sinesq_stag_ptrs
Definition: ERF.H:1343
void check_vels_for_nans(amrex::MultiFab const &xvel, amrex::MultiFab const &yvel, amrex::MultiFab const &zvel)
Definition: ERF.cpp:3065
void Evolve()
Definition: ERF.cpp:576
amrex::Vector< amrex::MultiFab > avg_xmom
Definition: ERF.H:851
amrex::Vector< amrex::MultiFab > pp_inc
Definition: ERF.H:847
static amrex::Real last_plot2d_file_time_2
Definition: ERF.H:1039
void make_eb_box()
amrex::Vector< ERFFillPatcher > FPr_u
Definition: ERF.H:923
amrex::Real cloud_fraction(amrex::Real time)
Definition: ERF_WriteScalarProfiles.cpp:451
amrex::Vector< amrex::IntVect > sampleline
Definition: ERF.H:1633
amrex::Vector< std::string > subvol3d_var_names
Definition: ERF.H:1122
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SFS_q1fx3_lev
Definition: ERF.H:946
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Hwave_onegrid
Definition: ERF.H:994
amrex::Vector< amrex::Vector< amrex::MultiFab > > vars_new
Definition: ERF.H:834
ERF(ERF &&) noexcept=delete
void project_initial_velocity(int lev, amrex::Real time, amrex::Real dt)
Definition: ERF_PoissonSolve.cpp:31
AMREX_FORCE_INLINE int NumSampleLineLogs() noexcept
Definition: ERF.H:1490
static void print_tpls(std::ostream &)
Definition: ERF_ConsoleIO.cpp:137
amrex::Vector< amrex::Real > dz_min
Definition: ERF.H:1371
amrex::Real m_plot2d_per_1
Definition: ERF.H:1107
amrex::Vector< amrex::MultiFab > lagged_delta_rt
Definition: ERF.H:850
amrex::Real plane_sampling_per
Definition: ERF.H:1616
std::string DataLogName(int i) const noexcept
The filename of the ith datalog file.
Definition: ERF.H:1636
std::string plot2d_file_2
Definition: ERF.H:1094
amrex::Vector< std::unique_ptr< amrex::MultiFab > > thin_yforce
Definition: ERF.H:1027
void ErrorEst(int lev, amrex::TagBoxArray &tags, amrex::Real time, int ngrow) override
Definition: ERF_Tagging.cpp:24
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SFS_hfx3_lev
Definition: ERF.H:944
amrex::Vector< std::unique_ptr< amrex::MultiFab > > sw_lw_fluxes
Definition: ERF.H:913
amrex::Vector< amrex::Vector< amrex::Real > > h_w_subsid
Definition: ERF.H:1311
static void print_banner(MPI_Comm, std::ostream &)
Definition: ERF_ConsoleIO.cpp:60
void volWgtColumnSum(int lev, const amrex::MultiFab &mf, int comp, amrex::MultiFab &mf_2d, const amrex::MultiFab &dJ)
Definition: ERF_VolWgtSum.cpp:82
static amrex::Real last_plot2d_file_time_1
Definition: ERF.H:1038
amrex::Vector< std::unique_ptr< std::fstream > > datalog
Definition: ERF.H:1620
static amrex::Real sum_per
Definition: ERF.H:1235
std::string MakeFilename_EyeTracker_maxvel(int nstep)
Definition: ERF_TrackerOutput.cpp:66
amrex::Vector< ERFFillPatcher > FPr_v
Definition: ERF.H:924
const amrex::Vector< std::string > derived_names_2d
Definition: ERF.H:1178
int cf_set_width
Definition: ERF.H:921
static amrex::Real previousCPUTimeUsed
Definition: ERF.H:1522
void setPlotVariables(const std::string &pp_plot_var_names, amrex::Vector< std::string > &plot_var_names)
Definition: ERF_Plotfile.cpp:25
amrex::Gpu::DeviceVector< amrex::Real > d_havg_temperature
Definition: ERF.H:1352
static int last_plot2d_file_step_2
Definition: ERF.H:1033
void fill_from_bndryregs(const amrex::Vector< amrex::MultiFab * > &mfs, amrex::Real time)
Definition: ERF_BoundaryConditionsBndryReg.cpp:13
const int timeprecision
Definition: ERF.H:1050
void setRecordDataInfo(int i, const std::string &filename)
Definition: ERF.H:1539
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SFS_hfx1_lev
Definition: ERF.H:944
amrex::Vector< std::array< amrex::Real, 2 > > hurricane_eye_track_xy
Definition: ERF.H:157
amrex::Vector< amrex::BoxArray > ba2d
Definition: ERF.H:1280
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_extdir_vals
Definition: ERF.H:1011
void init_from_input_sounding(int lev)
Definition: ERF_InitFromInputSounding.cpp:52
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:164
amrex::Vector< amrex::Vector< amrex::MultiFab > > gradp
Definition: ERF.H:838
std::string plot3d_file_1
Definition: ERF.H:1091
void FillPatchCrseLevel(int lev, amrex::Real time, const amrex::Vector< amrex::MultiFab * > &mfs_vel, bool cons_only=false)
Definition: ERF_FillPatch.cpp:288
eb_ const & get_eb(int lev) const noexcept
Definition: ERF.H:1648
amrex::Gpu::DeviceVector< amrex::Real > d_havg_qv
Definition: ERF.H:1354
static amrex::Real column_loc_y
Definition: ERF.H:1292
static bool plot_file_on_restart
Definition: ERF.H:1045
static int mg_verbose
Definition: ERF.H:1227
void ReadParameters()
Definition: ERF.cpp:2226
static amrex::Vector< std::string > PlotFileVarNames(amrex::Vector< std::string > plot_var_names)
Definition: ERF_Plotfile.cpp:306
amrex::Vector< amrex::Vector< amrex::MultiFab > > forecast_state_interp
Definition: ERF.H:165
bool FindInitialEye(int lev, const amrex::MultiFab &cc_vel, const amrex::Real velmag_threshold, amrex::Real &eye_x, amrex::Real &eye_y)
Definition: ERF_Tagging.cpp:548
void InitializeFromFile()
amrex::Vector< std::unique_ptr< ERFPhysBCFunct_cons > > physbcs_cons
Definition: ERF.H:856
void project_velocity_tb(int lev, amrex::Real dt, amrex::Vector< amrex::MultiFab > &vars)
Definition: ERF_PoissonSolve_tb.cpp:22
amrex::Vector< std::unique_ptr< amrex::MultiFab > > mf_PSFC
Definition: ERF.H:1285
ERF()
Definition: ERF.cpp:134
void init_Dirichlet_bc_data(const std::string input_file)
Definition: ERF_InitBCs.cpp:675
~ERF() override
amrex::Vector< std::unique_ptr< amrex::MultiFab > > z_phys_nd_src
Definition: ERF.H:961
amrex::Vector< std::unique_ptr< amrex::MultiFab > > detJ_cc
Definition: ERF.H:956
amrex::Vector< std::string > lsm_flux_name
Definition: ERF.H:898
void WriteMyEBSurface()
amrex::Vector< std::unique_ptr< amrex::MultiFab > > z_phys_cc
Definition: ERF.H:954
amrex::Vector< std::unique_ptr< eb_ > > eb
Definition: ERF.H:1646
amrex::Vector< std::unique_ptr< amrex::MultiFab > > eddyDiffs_lev
Definition: ERF.H:930
static SolverChoice solverChoice
Definition: ERF.H:1191
amrex::Vector< ERFFillPatcher > FPr_c
Definition: ERF.H:922
bool plot_rad
Definition: ERF.H:916
amrex::Vector< std::array< amrex::Real, 2 > > hurricane_tracker_circle
Definition: ERF.H:161
static bool use_fft
Definition: ERF.H:1228
bool m_plot_face_vels
Definition: ERF.H:1109
amrex::Vector< std::array< amrex::Real, 2 > > hurricane_maxvel_vs_time
Definition: ERF.H:159
amrex::Vector< amrex::MultiFab > base_state_new
Definition: ERF.H:989
std::string plot3d_file_2
Definition: ERF.H:1092
amrex::Vector< std::unique_ptr< amrex::MultiFab > > az
Definition: ERF.H:959
int regrid_int
Definition: ERF.H:1084
amrex::Vector< amrex::Real > fixed_dt
Definition: ERF.H:1078
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:1876
void solve_with_gmres(int lev, const amrex::Box &subdomain, amrex::MultiFab &rhs, amrex::MultiFab &p, amrex::Array< amrex::MultiFab, AMREX_SPACEDIM > &fluxes, amrex::MultiFab &ax_sub, amrex::MultiFab &ay_sub, amrex::MultiFab &az_sub, amrex::MultiFab &, amrex::MultiFab &znd_sub)
Definition: ERF_SolveWithGMRES.cpp:12
amrex::Vector< amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > > d_sponge_ptrs
Definition: ERF.H:1339
amrex::Vector< long > dt_mri_ratio
Definition: ERF.H:829
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > Tau
Definition: ERF.H:928
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vel_t_avg
Definition: ERF.H:841
void sum_energy_quantities(amrex::Real time)
Definition: ERF_WriteScalarProfiles.cpp:312
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SFS_q2fx3_lev
Definition: ERF.H:947
static amrex::Real dt_max
Definition: ERF.H:1075
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::iMultiFab > > > lmask_lev
Definition: ERF.H:936
AMREX_FORCE_INLINE int NumSamplePointLogs() noexcept
Definition: ERF.H:1476
amrex::Vector< amrex::Real > h_havg_pressure
Definition: ERF.H:1347
void update_diffusive_arrays(int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Definition: ERF_MakeNewArrays.cpp:537
static int verbose
Definition: ERF.H:1226
amrex::Vector< amrex::Real > h_havg_qc
Definition: ERF.H:1349
amrex::Vector< std::unique_ptr< ERFPhysBCFunct_w > > physbcs_w
Definition: ERF.H:859
static int nghost_eb_basic()
Definition: ERF.H:1658
void Advance(int lev, amrex::Real time, amrex::Real dt_lev, int iteration, int ncycle)
Definition: ERF_Advance.cpp:20
void FillForecastStateMultiFabs(const int lev, const std::string &filename, const std::unique_ptr< amrex::MultiFab > &z_phys_nd, amrex::Vector< amrex::Vector< amrex::MultiFab >> &forecast_state)
Definition: ERF_WeatherDataInterpolation.cpp:64
static std::string column_file_name
Definition: ERF.H:1293
amrex::Vector< std::unique_ptr< std::fstream > > samplelinelog
Definition: ERF.H:1631
static amrex::Real last_plot3d_file_time_2
Definition: ERF.H:1037
amrex::Vector< std::unique_ptr< amrex::MultiFab > > terrain_blanking
Definition: ERF.H:974
std::unique_ptr< Microphysics > micro
Definition: ERF.H:879
int m_plot2d_int_2
Definition: ERF.H:1100
amrex::Vector< amrex::MultiFab > base_state
Definition: ERF.H:988
AMREX_FORCE_INLINE amrex::YAFluxRegister * getAdvFluxReg(int lev)
Definition: ERF.H:1433
int m_plot3d_int_1
Definition: ERF.H:1097
amrex::Vector< amrex::Real > h_havg_density
Definition: ERF.H:1345
amrex::Vector< std::unique_ptr< amrex::MultiFab > > z_phys_nd_new
Definition: ERF.H:968
amrex::Vector< amrex::Real > fixed_fast_dt
Definition: ERF.H:1079
void project_momenta(int lev, amrex::Real l_time, amrex::Real l_dt, amrex::Vector< amrex::MultiFab > &vars)
Definition: ERF_PoissonSolve.cpp:100
const amrex::Vector< std::string > derived_subvol_names
Definition: ERF.H:1188
amrex::Vector< std::array< amrex::Real, 2 > > hurricane_minpressure_vs_time
Definition: ERF.H:160
bool metgrid_retain_sfc
Definition: ERF.H:1274
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Qv_prim
Definition: ERF.H:864
static int sum_interval
Definition: ERF.H:1233
static int pert_interval
Definition: ERF.H:1234
amrex::Real line_sampling_per
Definition: ERF.H:1615
void restart()
Definition: ERF.cpp:2035
static int last_plot3d_file_step_2
Definition: ERF.H:1031
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SFS_q1fx2_lev
Definition: ERF.H:946
amrex::IntVect & SampleLine(int i)
Definition: ERF.H:1509
void build_fine_mask(int lev, amrex::MultiFab &fine_mask)
Definition: ERF_VolWgtSum.cpp:125
int file_name_digits
Definition: ERF.H:1243
amrex::Vector< amrex::MultiFab > rV_new
Definition: ERF.H:871
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:835
amrex::Vector< std::unique_ptr< amrex::MultiFab > > thin_zforce
Definition: ERF.H:1028
amrex::Vector< std::string > plot3d_var_names_2
Definition: ERF.H:1125
amrex::Vector< amrex::BCRec > domain_bcs_type
Definition: ERF.H:1004
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Qr_prim
Definition: ERF.H:865
void post_timestep(int nstep, amrex::Real time, amrex::Real dt_lev)
Definition: ERF.cpp:734
std::unique_ptr< amrex::MultiFab > mf_MUB
Definition: ERF.H:1283
std::string pp_prefix
Definition: ERF.H:546
std::string SampleLineLogName(int i) const noexcept
The filename of the ith samplelinelog file.
Definition: ERF.H:1643
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > sst_lev
Definition: ERF.H:934
amrex::Vector< std::string > plot2d_var_names_1
Definition: ERF.H:1126
amrex::Vector< std::unique_ptr< amrex::MultiFab > > thin_xforce
Definition: ERF.H:1026
AMREX_FORCE_INLINE int NumSamplePoints() noexcept
Definition: ERF.H:1503
std::unique_ptr< amrex::MultiFab > mf_C2H
Definition: ERF.H:1282
bool metgrid_use_sfc
Definition: ERF.H:1273
amrex::Real m_plot2d_per_2
Definition: ERF.H:1108
static amrex::Real bndry_output_planes_per
Definition: ERF.H:1298
void setPlotVariables2D(const std::string &pp_plot_var_names, amrex::Vector< std::string > &plot_var_names)
Definition: ERF_Plotfile.cpp:187
amrex::Real m_check_per
Definition: ERF.H:1120
void init_custom(int lev)
Definition: ERF_InitCustomPertState.cpp:26
amrex::Vector< std::unique_ptr< IRadiation > > rad
Definition: ERF.H:901
std::unique_ptr< ProblemBase > prob
Definition: ERF.H:816
amrex::Vector< int > num_files_at_level
Definition: ERF.H:819
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:1789
void init_bcs()
Definition: ERF_InitBCs.cpp:287
int profile_int
Definition: ERF.H:1114
bool metgrid_debug_quiescent
Definition: ERF.H:1265
amrex::Vector< std::unique_ptr< ERFPhysBCFunct_u > > physbcs_u
Definition: ERF.H:857
amrex::Vector< amrex::Real > t_new
Definition: ERF.H:826
bool destag_profiles
Definition: ERF.H:1115
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > stretched_dz_d
Definition: ERF.H:986
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > th_bc_data
Definition: ERF.H:793
amrex::Vector< amrex::Real > t_avg_cnt
Definition: ERF.H:842
void FillPatchFineLevel(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)
Definition: ERF_FillPatch.cpp:20
std::string DerDataLogName(int i) const noexcept
Definition: ERF.H:1637
void remake_zphys(int lev, amrex::Real time, std::unique_ptr< amrex::MultiFab > &temp_zphys_nd)
Definition: ERF_MakeNewArrays.cpp:758
int m_check_int
Definition: ERF.H:1119
amrex::Vector< std::unique_ptr< amrex::MultiFab > > solar_zenith
Definition: ERF.H:914
amrex::Real estTimeStep(int lev, long &dt_fast_ratio) const
Definition: ERF_ComputeTimestep.cpp:54
static std::string sponge_type
Definition: ERF.H:1249
static amrex::Real startCPUTime
Definition: ERF.H:1521
amrex::Vector< amrex::MultiFab > surface_state_1
Definition: ERF.H:166
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SFS_diss_lev
Definition: ERF.H:945
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > d_sinesq_ptrs
Definition: ERF.H:1342
amrex::Vector< amrex::MultiFab > rU_old
Definition: ERF.H:868
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:827
AMREX_FORCE_INLINE int NumSampleLines() noexcept
Definition: ERF.H:1516
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Theta_prim
Definition: ERF.H:863
void init1DArrays()
amrex::Vector< std::unique_ptr< amrex::MultiFab > > z_t_rk
Definition: ERF.H:971
void HurricaneTracker(int lev, amrex::Real time, const amrex::MultiFab &cc_vel, const amrex::Real velmag_threshold, amrex::TagBoxArray *tags=nullptr)
Definition: ERF_Tagging.cpp:648
static int check_for_nans
Definition: ERF.H:1230
void make_subdomains(const amrex::BoxList &ba, amrex::Vector< amrex::BoxArray > &bins)
Definition: ERF_MakeSubdomains.cpp:6
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_TrackerOutput.cpp:11
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:229
amrex::Vector< std::unique_ptr< ERFPhysBCFunct_v > > physbcs_v
Definition: ERF.H:858
static amrex::Real column_per
Definition: ERF.H:1290
amrex::Vector< std::string > tot_e_datalogname
Definition: ERF.H:1625
static amrex::Real stop_time
Definition: ERF.H:1061
static int output_bndry_planes
Definition: ERF.H:1296
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > d_v_geos
Definition: ERF.H:1318
void update_terrain_arrays(int lev)
Definition: ERF_MakeNewArrays.cpp:818
static std::string nc_bdy_file
Definition: ERF.H:1256
bool metgrid_interp_theta
Definition: ERF.H:1270
static amrex::Vector< amrex::Vector< std::string > > nc_init_file
Definition: ERF.H:1252
void init_only(int lev, amrex::Real time)
Definition: ERF.cpp:2097
static int input_bndry_planes
Definition: ERF.H:1302
static StateInterpType interpolation_type
Definition: ERF.H:1246
amrex::Gpu::DeviceVector< amrex::Real > d_havg_qc
Definition: ERF.H:1355
void AverageDown()
Definition: ERF_AverageDown.cpp:16
amrex::Vector< amrex::Vector< amrex::Real > > h_v_geos
Definition: ERF.H:1317
bool regrid_level_0_on_restart
Definition: ERF.H:1088
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Lwave_onegrid
Definition: ERF.H:995
InputSoundingData input_sounding_data
Definition: ERF.H:784
amrex::Vector< amrex::Vector< amrex::Real > > h_sinesq_ptrs
Definition: ERF.H:1334
void Write2DPlotFile(int which, PlotFileType plotfile_type, amrex::Vector< std::string > plot_var_names)
Definition: ERF_Plotfile.cpp:1965
amrex::Gpu::DeviceVector< amrex::Real > d_havg_density
Definition: ERF.H:1351
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rhotheta_src
Definition: ERF.H:1308
void init_from_hse(int lev)
Definition: ERF_InitFromHSE.cpp:32
const std::string datetime_format
Definition: ERF.H:1064
amrex::Vector< amrex::YAFluxRegister * > advflux_reg
Definition: ERF.H:999
int metgrid_force_sfc_k
Definition: ERF.H:1277
void setSubVolVariables(const std::string &pp_subvol_var_names, amrex::Vector< std::string > &subvol_var_names)
Definition: ERF_WriteSubvolume.cpp:9
static int ng_pres_hse
Definition: ERF.H:1305
static amrex::Real bndry_output_planes_start_time
Definition: ERF.H:1299
bool real_extrap_w
Definition: ERF.H:1259
static amrex::Real cfl
Definition: ERF.H:1070
amrex::Vector< std::unique_ptr< amrex::MultiFab > > fine_mask
Definition: ERF.H:983
amrex::Gpu::DeviceVector< amrex::BCRec > domain_bcs_type_d
Definition: ERF.H:1005
amrex::Vector< std::unique_ptr< ForestDrag > > m_forest_drag
Definition: ERF.H:1362
const int datwidth
Definition: ERF.H:1048
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::Vector< amrex::BoxArray > ba1d
Definition: ERF.H:1279
InputSpongeData input_sponge_data
Definition: ERF.H:787
amrex::Vector< amrex::Vector< amrex::BoxArray > > subdomains
Definition: ERF.H:1369
void check_state_for_nans(amrex::MultiFab const &S)
Definition: ERF.cpp:3042
std::string restart_chkfile
Definition: ERF.H:1067
bool metgrid_use_below_sfc
Definition: ERF.H:1272
void init_zphys(int lev, amrex::Real time)
Definition: ERF_MakeNewArrays.cpp:662
amrex::Vector< std::string > sampleptlogname
Definition: ERF.H:1628
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > xvel_bc_data
Definition: ERF.H:790
void InitData_pre()
Definition: ERF.cpp:956
amrex::IntVect & SamplePoint(int i)
Definition: ERF.H:1496
bool use_datetime
Definition: ERF.H:1063
amrex::Vector< std::array< amrex::Real, 2 > > hurricane_track_xy
Definition: ERF.H:156
amrex::Vector< amrex::Vector< amrex::Real > > h_u_geos
Definition: ERF.H:1314
void InitializeLevelFromData(int lev, const amrex::MultiFab &initial_data)
std::string subvol_file
Definition: ERF.H:1095
amrex::Vector< amrex::Real > m_subvol_per
Definition: ERF.H:1103
int rad_datalog_int
Definition: ERF.H:917
void sum_derived_quantities(amrex::Real time)
Definition: ERF_WriteScalarProfiles.cpp:177
void sum_integrated_quantities(amrex::Real time)
Definition: ERF_WriteScalarProfiles.cpp:15
amrex::Vector< std::string > datalogname
Definition: ERF.H:1623
void initHSE()
Initialize HSE.
Definition: ERF_Init1D.cpp:146
static amrex::Real column_loc_x
Definition: ERF.H:1291
amrex::Vector< amrex::MultiFab > surface_state_2
Definition: ERF.H:167
amrex::Vector< std::unique_ptr< amrex::MultiFab > > ax
Definition: ERF.H:957
amrex::Vector< std::unique_ptr< amrex::MultiFab > > z_phys_nd
Definition: ERF.H:953
void MakeDiagnosticAverage(amrex::Vector< amrex::Real > &h_havg, amrex::MultiFab &S, int n)
Definition: ERF.cpp:2898
amrex::Real metgrid_proximity
Definition: ERF.H:1275
amrex::Real volWgtSumMF(int lev, const amrex::MultiFab &mf, int comp, const amrex::MultiFab &dJ, const amrex::MultiFab &mfx, const amrex::MultiFab &mfy, bool finemask, bool local=true)
Definition: ERF_VolWgtSum.cpp:20
void setRecordDerDataInfo(int i, const std::string &filename)
Definition: ERF.H:1552
amrex::Vector< amrex::Real > h_havg_temperature
Definition: ERF.H:1346
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:1627
void poisson_wall_dist(int lev)
Definition: ERF_PoissonWallDist.cpp:22
amrex::Vector< std::array< amrex::Real, 2 > > hurricane_eye_track_latlon
Definition: ERF.H:158
amrex::Vector< std::unique_ptr< amrex::MultiFab > > detJ_cc_src
Definition: ERF.H:963
amrex::Gpu::DeviceVector< amrex::Real > d_havg_pressure
Definition: ERF.H:1353
std::string plot2d_file_1
Definition: ERF.H:1093
amrex::Vector< std::string > der_datalogname
Definition: ERF.H:1624
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SmnSmn_lev
Definition: ERF.H:931
const amrex::Vector< std::string > derived_names
Definition: ERF.H:1135
std::string MakeVTKFilename_TrackerCircle(int nstep)
Definition: ERF_TrackerOutput.cpp:24
amrex::Vector< std::unique_ptr< amrex::MultiFab > > ay_src
Definition: ERF.H:965
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:1440
void writeJobInfo(const std::string &dir) const
Definition: ERF_WriteJobInfo.cpp:10
std::string MakeVTKFilename_EyeTracker_xy(int nstep)
Definition: ERF_TrackerOutput.cpp:38
void ComputeDt(int step=-1)
Definition: ERF_ComputeTimestep.cpp:11
amrex::Vector< int > nsubsteps
Definition: ERF.H:823
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > yflux_imask
Definition: ERF.H:1021
amrex::Vector< amrex::MultiFab > rW_new
Definition: ERF.H:873
amrex::Vector< amrex::MultiFab > weather_forecast_data_2
Definition: ERF.H:162
amrex::Vector< std::unique_ptr< amrex::MultiFab > > lon_m
Definition: ERF.H:779
std::unique_ptr< WriteBndryPlanes > m_w2d
Definition: ERF.H:1359
AMREX_FORCE_INLINE std::ostream & SampleLineLog(int i)
Definition: ERF.H:1483
amrex::Vector< amrex::Vector< amrex::MultiFab * > > lsm_flux
Definition: ERF.H:899
amrex::Vector< std::string > plot3d_var_names_1
Definition: ERF.H:1124
std::string SamplePointLogName(int i) const noexcept
The filename of the ith sampleptlog file.
Definition: ERF.H:1640
void turbPert_update(const int lev, const amrex::Real dt)
Definition: ERF_InitTurbPert.cpp:12
void InitData_post()
Definition: ERF.cpp:980
void refinement_criteria_setup()
Definition: ERF_Tagging.cpp:320
void apply_gaussian_smoothing_to_perturbations(const int lev, amrex::MultiFab &cons_pert, amrex::MultiFab &xvel_pert, amrex::MultiFab &yvel_pert, amrex::MultiFab &zvel_pert)
Definition: ERF_InitCustomPertState.cpp:135
static int nghost_eb_volume()
Definition: ERF.H:1662
bool metgrid_debug_dry
Definition: ERF.H:1267
static AMREX_FORCE_INLINE int ComputeGhostCells(const SolverChoice &sc)
Definition: ERF.H:1375
static int bndry_output_planes_interval
Definition: ERF.H:1297
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:149
amrex::Vector< std::string > plot2d_var_names_2
Definition: ERF.H:1127
void WeatherDataInterpolation(const int nlevs, const amrex::Real time, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &z_phys_nd, bool regrid_forces_file_read)
Definition: ERF_WeatherDataInterpolation.cpp:347
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > Tau_corr
Definition: ERF.H:929
void WriteSubvolume(int isub, amrex::Vector< std::string > subvol_var_names)
Definition: ERF_WriteSubvolume.cpp:145
static int nghost_eb_full()
Definition: ERF.H:1665
void RemakeLevel(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Definition: ERF_MakeNewLevel.cpp:505
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:877
bool metgrid_debug_isothermal
Definition: ERF.H:1266
amrex::Vector< std::string > lsm_data_name
Definition: ERF.H:896
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rhoqt_src
Definition: ERF.H:1309
void initSponge()
Initialize sponge profiles.
Definition: ERF_InitSponge.cpp:35
std::unique_ptr< PlaneSampler > plane_sampler
Definition: ERF.H:1618
void check_for_low_temp(amrex::MultiFab &S)
Definition: ERF.cpp:3092
amrex::Vector< std::unique_ptr< std::fstream > > der_datalog
Definition: ERF.H:1621
static amrex::Vector< amrex::Vector< int > > have_read_nc_init_file
Definition: ERF.H:1253
amrex::Real m_plot3d_per_2
Definition: ERF.H:1106
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:1622
amrex::Vector< std::unique_ptr< amrex::MultiFab > > sinPhi_m
Definition: ERF.H:781
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:1257
void MakeNewLevelFromCoarse(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Definition: ERF_MakeNewLevel.cpp:272
amrex::Vector< int > last_subvol_step
Definition: ERF.H:1042
amrex::Vector< std::unique_ptr< amrex::MultiFab > > ax_src
Definition: ERF.H:964
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > urb_frac_lev
Definition: ERF.H:941
amrex::Vector< std::unique_ptr< amrex::MultiFab > > z_phys_cc_src
Definition: ERF.H:962
void input_sponge(int lev)
Definition: ERF_InitSponge.cpp:17
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > d_u_geos
Definition: ERF.H:1315
void Define_ERFFillPatchers(int lev)
Definition: ERF.cpp:2970
AMREX_FORCE_INLINE int NumDataLogs() noexcept
Definition: ERF.H:1454
static PlotFileType plotfile3d_type_2
Definition: ERF.H:1239
void setRecordEnergyDataInfo(int i, const std::string &filename)
Definition: ERF.H:1565
TurbulentPerturbation turbPert
Definition: ERF.H:1194
amrex::Vector< amrex::MultiFab > rW_old
Definition: ERF.H:872
amrex::Vector< amrex::Vector< amrex::MultiFab > > forecast_state_2
Definition: ERF.H:164
void check_for_negative_theta(amrex::MultiFab &S)
Definition: ERF.cpp:3127
void ImposeBCsOnPhi(int lev, amrex::MultiFab &phi, const amrex::Box &subdomain)
Definition: ERF_ImposeBCsOnPhi.cpp:12
void Interp2DArrays(int lev, const amrex::BoxArray &my_ba2d, const amrex::DistributionMapping &my_dm)
Definition: ERF.cpp:1851
void FillCoarsePatch(int lev, amrex::Real time)
Definition: ERF_FillCoarsePatch.cpp:21
void ClearLevel(int lev) override
Definition: ERF_MakeNewLevel.cpp:768
std::string MakeFilename_EyeTracker_minpressure(int nstep)
Definition: ERF_TrackerOutput.cpp:80
std::unique_ptr< SurfaceLayer > m_SurfaceLayer
Definition: ERF.H:1361
static PlotFileType plotfile2d_type_2
Definition: ERF.H:1241
int plane_sampling_interval
Definition: ERF.H:1614
static amrex::Vector< amrex::AMRErrorTag > ref_tags
Definition: ERF.H:1367
void make_physbcs(int lev)
Definition: ERF_MakeNewArrays.cpp:857
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > d_w_subsid
Definition: ERF.H:1312
amrex::Vector< ERFFillPatcher > FPr_w
Definition: ERF.H:925
int real_set_width
Definition: ERF.H:1258
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SFS_hfx2_lev
Definition: ERF.H:944
void WriteLinePlot(const std::string &filename, amrex::Vector< std::array< amrex::Real, 2 >> &points_xy)
Definition: ERF_Write1DProfiles.cpp:574
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::iMultiFab > > > soil_type_lev
Definition: ERF.H:940
bool writeNow(const amrex::Real cur_time, const int nstep, const int plot_int, const amrex::Real plot_per, const amrex::Real dt_0, amrex::Real &last_file_time)
Definition: ERF.cpp:3021
amrex::Vector< amrex::Vector< amrex::Real > > zlevels_stag
Definition: ERF.H:950
AMREX_FORCE_INLINE int NumDerDataLogs() noexcept
Definition: ERF.H:1461
amrex::EBFArrayBoxFactory const & EBFactory(int lev) const noexcept
Definition: ERF.H:1654
amrex::Vector< amrex::Vector< amrex::MultiFab * > > lsm_data
Definition: ERF.H:897
void init_immersed_forcing(int lev)
Definition: ERF_InitImmersedForcing.cpp:15
amrex::Vector< amrex::Vector< amrex::Real > > stretched_dz_h
Definition: ERF.H:985
bool use_real_time_in_pltname
Definition: ERF.H:1244
int m_plot2d_int_1
Definition: ERF.H:1099
void WriteCheckpointFile() const
Definition: ERF_Checkpoint.cpp:26
static int output_1d_column
Definition: ERF.H:1288
bool metgrid_debug_psfc
Definition: ERF.H:1268
amrex::Vector< std::unique_ptr< amrex::MultiFab > > az_src
Definition: ERF.H:966
static int fixed_mri_dt_ratio
Definition: ERF.H:1080
void create_random_perturbations(const int lev, amrex::MultiFab &cons_pert, amrex::MultiFab &xvel_pert, amrex::MultiFab &yvel_pert, amrex::MultiFab &zvel_pert)
Definition: ERF_InitCustomPertState.cpp:113
amrex::Vector< amrex::Real > dt
Definition: ERF.H:828
static amrex::Real init_shrink
Definition: ERF.H:1072
void WriteVTKPolyline(const std::string &filename, amrex::Vector< std::array< amrex::Real, 2 >> &points_xy)
Definition: ERF_TrackerOutput.cpp:94
void Write3DPlotFile(int which, PlotFileType plotfile_type, amrex::Vector< std::string > plot_var_names)
Definition: ERF_Plotfile.cpp:318
void advance_radiation(int lev, amrex::MultiFab &cons_in, const amrex::Real &dt_advance)
Definition: ERF_AdvanceRadiation.cpp:5
static amrex::Real dt_max_initial
Definition: ERF.H:1074
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:653
static std::string nc_low_file
Definition: ERF.H:1262
amrex::Vector< amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > > d_rayleigh_ptrs
Definition: ERF.H:1338
void InitData()
Definition: ERF.cpp:947
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:993
amrex::Vector< std::unique_ptr< amrex::MultiFab > > cosPhi_m
Definition: ERF.H:781
void ParameterSanityChecks()
Definition: ERF.cpp:2725
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::iMultiFab > > > land_type_lev
Definition: ERF.H:939
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
amrex::Vector< amrex::Vector< amrex::MultiFab > > forecast_state_1
Definition: ERF.H:163
int cf_width
Definition: ERF.H:920
static int last_plot2d_file_step_1
Definition: ERF.H:1032
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > zflux_imask
Definition: ERF.H:1022
bool m_expand_plotvars_to_unif_rr
Definition: ERF.H:1096
amrex::Real m_plot3d_per_1
Definition: ERF.H:1105
void Construct_ERFFillPatchers(int lev)
Definition: ERF.cpp:2944
void post_update(amrex::MultiFab &state_mf, amrex::Real time, const amrex::Geometry &geom)
amrex::Vector< int > num_boxes_at_level
Definition: ERF.H:818
void init_phys_bcs(bool &rho_read, bool &read_prim_theta)
Definition: ERF_InitBCs.cpp:20
void FillSurfaceStateMultiFabs(const int lev, const std::string &filename, amrex::Vector< amrex::MultiFab > &surface_state)
Definition: ERF_SurfaceDataInterpolation.cpp:19
static amrex::Real sub_cfl
Definition: ERF.H:1071
static void print_error(MPI_Comm, const std::string &msg)
Definition: ERF_ConsoleIO.cpp:43
static int ng_dens_hse
Definition: ERF.H:1304
std::unique_ptr< ReadBndryPlanes > m_r2d
Definition: ERF.H:1360
amrex::Vector< std::unique_ptr< amrex::MultiFab > > ay
Definition: ERF.H:958
amrex::Vector< amrex::Vector< amrex::MultiFab * > > qmoist
Definition: ERF.H:880
amrex::Vector< std::unique_ptr< amrex::MultiFab > > SFS_q1fx1_lev
Definition: ERF.H:946
static amrex::Real getCPUTime()
Definition: ERF.H:1525
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > zvel_bc_data
Definition: ERF.H:792
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:895
static amrex::Real change_max
Definition: ERF.H:1073
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:1595
void setSpongeRefFromSounding(bool restarting)
Set sponge mean profiles from input sounding.
Definition: ERF_InitSponge.cpp:65
amrex::Vector< amrex::MultiFab > avg_ymom
Definition: ERF.H:852
int line_sampling_interval
Definition: ERF.H:1613
amrex::Vector< amrex::Vector< amrex::Real > > h_sinesq_stag_ptrs
Definition: ERF.H:1335
amrex::Vector< std::unique_ptr< amrex::MultiFab > > detJ_cc_new
Definition: ERF.H:969
amrex::Array< std::string, 2 *AMREX_SPACEDIM > domain_bc_type
Definition: ERF.H:1008
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real > > yvel_bc_data
Definition: ERF.H:791
void init_uniform(int lev)
Definition: ERF_InitUniform.cpp:17
static void writeBuildInfo(std::ostream &os)
Definition: ERF_WriteJobInfo.cpp:144
amrex::GpuArray< ERF_BC, AMREX_SPACEDIM *2 > phys_bc_type
Definition: ERF.H:1017
std::string check_file
Definition: ERF.H:1118
amrex::Vector< amrex::IntVect > samplepoint
Definition: ERF.H:1629
amrex::Vector< amrex::MultiFab > surface_state_interp
Definition: ERF.H:168
amrex::Vector< amrex::Real > h_havg_qv
Definition: ERF.H:1348
amrex::Vector< amrex::Real > last_subvol_time
Definition: ERF.H:1043
amrex::Vector< amrex::MultiFab > weather_forecast_data_1
Definition: ERF.H:162
static void print_usage(MPI_Comm, std::ostream &)
Definition: ERF_ConsoleIO.cpp:26
amrex::Vector< amrex::MultiFab > rV_old
Definition: ERF.H:870
amrex::Vector< std::unique_ptr< amrex::MultiFab > > lat_m
Definition: ERF.H:779
amrex::Vector< amrex::Vector< amrex::Vector< amrex::Real > > > h_rayleigh_ptrs
Definition: ERF.H:1330
void make_eb_regular()
static amrex::Real last_check_file_time
Definition: ERF.H:1040
std::unique_ptr< amrex::MultiFab > mf_C1H
Definition: ERF.H:1281
static int last_plot3d_file_step_1
Definition: ERF.H:1030
static amrex::Real last_plot3d_file_time_1
Definition: ERF.H:1036
std::unique_ptr< LineSampler > line_sampler
Definition: ERF.H:1617
amrex::Vector< std::unique_ptr< amrex::MultiFab > > Hwave
Definition: ERF.H:992
amrex::Vector< int > istep
Definition: ERF.H:822
const int datprecision
Definition: ERF.H:1049
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > xflux_imask
Definition: ERF.H:1020
void derive_upwp(amrex::Vector< amrex::Real > &h_havg)
static PlotFileType plotfile2d_type_1
Definition: ERF.H:1240
int metgrid_order
Definition: ERF.H:1276
bool finished_wave
Definition: ERF.H:996
void ReadCheckpointFile()
Definition: ERF_Checkpoint.cpp:458
amrex::Vector< amrex::Vector< amrex::MultiFab > > vars_old
Definition: ERF.H:835
amrex::Vector< std::unique_ptr< ERFPhysBCFunct_base > > physbcs_base
Definition: ERF.H:860
void init_thin_body(int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Definition: ERF_MakeNewLevel.cpp:837
static PlotFileType plotfile3d_type_1
Definition: ERF.H:1238
amrex::Vector< int > m_subvol_int
Definition: ERF.H:1102
AMREX_FORCE_INLINE std::ostream & SamplePointLog(int i)
Definition: ERF.H:1469
std::string MakeFilename_EyeTracker_latlon(int nstep)
Definition: ERF_TrackerOutput.cpp:52
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NBCVAR_max > m_bc_neumann_vals
Definition: ERF.H:1014
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > tsk_lev
Definition: ERF.H:935
void setRecordSamplePointInfo(int i, int lev, amrex::IntVect &cell, const std::string &filename)
Definition: ERF.H:1578
static int column_interval
Definition: ERF.H:1289
static void print_summary(std::ostream &)
int m_plot3d_int_2
Definition: ERF.H:1098
amrex::Vector< std::unique_ptr< amrex::MultiFab > > qheating_rates
Definition: ERF.H:902
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
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rad_fluxes
Definition: ERF.H:903
void initializeMicrophysics(const int &)
Definition: ERF.cpp:2001
bool plot_lsm
Definition: ERF.H:1111
amrex::Vector< amrex::MultiFab > avg_zmom
Definition: ERF.H:853
const amrex::Vector< std::string > cons_names
Definition: ERF.H:1128
void ReadCheckpointFileSurfaceLayer()
Definition: ERF_Checkpoint.cpp:1058
void timeStep(int lev, amrex::Real time, int iteration)
Definition: ERF_TimeStep.cpp:17
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
@ xvel
Definition: ERF_IndexDefines.H:141
@ cons
Definition: ERF_IndexDefines.H:140
@ zvel
Definition: ERF_IndexDefines.H:143
@ yvel
Definition: ERF_IndexDefines.H:142
Definition: ERF_ConsoleIO.cpp:12
real(c_double), parameter, private pi
Definition: ERF_module_mp_morr_two_moment.F90:100
integer, private isub
Definition: ERF_module_mp_morr_two_moment.F90:164
AdvType moistscal_horiz_adv_type
Definition: ERF_AdvStruct.H:423
AdvType dycore_vert_adv_type
Definition: ERF_AdvStruct.H:420
AdvType moistscal_vert_adv_type
Definition: ERF_AdvStruct.H:424
AdvType dryscal_horiz_adv_type
Definition: ERF_AdvStruct.H:421
AdvType dycore_horiz_adv_type
Definition: ERF_AdvStruct.H:419
AdvType dryscal_vert_adv_type
Definition: ERF_AdvStruct.H:422
Definition: ERF_InputSoundingData.H:22
Definition: ERF_InputSpongeData.H:19
Definition: ERF_DataStruct.H:129
bool use_num_diff
Definition: ERF_DataStruct.H:1164
AdvChoice advChoice
Definition: ERF_DataStruct.H:1056
static TerrainType terrain_type
Definition: ERF_DataStruct.H:1036
Definition: ERF_TurbPertStruct.H:22