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