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