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