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