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