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