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