REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA.H
Go to the documentation of this file.
1#ifndef REMORA_H_
2#define REMORA_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_InterpFaceRegister.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
25#ifdef AMREX_MEM_PROFILING
26#include <AMReX_MemProfiler.H>
27#endif
28
29#include <REMORA_Math.H>
30#include <REMORA_PhysBCFunct.H>
31#include <REMORA_FillPatcher.H>
32
33#include <REMORA_IndexDefines.H>
35#include <REMORA_DataStruct.H>
36#include <REMORA_Derive.H>
37#include "REMORA_prob_common.H"
38
39#ifdef REMORA_USE_PARTICLES
40#include "REMORA_ParticleData.H"
41#endif
42
43#ifdef REMORA_USE_NETCDF
44#include "REMORA_NCInterface.H"
45#include "REMORA_NCTimeSeries.H"
48#endif
49
50#ifdef REMORA_USE_MOAB
51#include "REMORA_MOAB.H"
52#endif
54#include <iostream>
55
56#ifdef AMREX_LAZY
57#include <AMReX_Lazy.H>
58#endif
59
60#ifndef AMREX_USE_MPI
61using amrex::MPI_COMM_WORLD;
62using amrex::MPI_Comm;
63#endif
64
76
77// Forward declarations to resolve circular dependency
78class ProblemBase;
79#ifdef REMORA_USE_NETCDF
80class NCTimeSeries;
83#endif
84
85/*!
86 * \brief Class that stores all relevant simulation state data with methods for time stepping
87 */
88class REMORA
89 : public amrex::AmrCore
90{
91public:
92
93 REMORA ();
94 virtual ~REMORA ();
95
96 /** \brief Advance solution to final time */
97 void Evolve ();
98
99 /** \brief Tag cells for refinement */
100 virtual void ErrorEst (int lev, amrex::TagBoxArray& tags, amrex::Real time, int ngrow) override;
101
102 /** \brief Initialize multilevel data */
103 void InitData ();
104
105 /** \brief Init (NOT restart or regrid) */
106 void init_only (int lev, amrex::Real time);
107
108 /** \briefRestart */
109 void restart ();
110
111 /** \brief Called after every level 0 timestep */
112 void post_timestep (int nstep, amrex::Real time, amrex::Real dt_lev);
113
114 // Diagnostics
115
116 /** \brief Integrate conserved quantities for diagnostics */
117 void sum_integrated_quantities (amrex::Real time);
118
119 /** \brief Perform the volume-weighted sum */
120 amrex::Real
121 volWgtSumMF (int lev,
122 const amrex::MultiFab& mf, int comp, bool local, bool finemask);
123
124 /** \brief Decide if it is time to take an action */
125 bool is_it_time_for_action (int nstep, amrex::Real time, amrex::Real dt,
126 int action_interval, amrex::Real action_per);
127
128 /** \brief Make a new level using provided BoxArray and DistributionMapping and
129 * fill with interpolated coarse level data.
130 * Overrides the pure virtual function in AmrCore
131 */
132 virtual void MakeNewLevelFromCoarse (int lev, amrex::Real time, const amrex::BoxArray& ba,
133 const amrex::DistributionMapping& dm) override;
134
135 /** \brief Remake an existing level using provided BoxArray and DistributionMapping and
136 * fill with existing fine and coarse data.
137 * Overrides the pure virtual function in AmrCore
138 */
139 virtual void RemakeLevel (int lev, amrex::Real time, const amrex::BoxArray& ba,
140 const amrex::DistributionMapping& dm) override;
141
142 /** \brief Delete level data
143 * Overrides the pure virtual function in AmrCore
144 */
145 virtual void ClearLevel (int lev) override;
146
147 /** \brief Make a new level from scratch using provided BoxArray and DistributionMapping.
148 * Only used during initialization.
149 * Overrides the pure virtual function in AmrCore
150 */
151 virtual void MakeNewLevelFromScratch (int lev, amrex::Real time, const amrex::BoxArray& ba,
152 const amrex::DistributionMapping& dm) override;
153
154 /** \brief Set pm and pn arrays on level lev. Only works if using analytic initialization */
155 void set_grid_scale (int lev);
156
157 /** \brief Set zeta components to be equal to time-averaged Zt_avg1 */
158 void set_zeta_to_Ztavg (int lev);
159
160 /** \brief Calculate u-, v-, and psi-point masks based on rho-point masks after analytic initialization */
161 void calculate_nodal_masks (int lev);
162
163 /** \brief Copy maskr to all z levels */
164 void fill_3d_masks (int lev);
165
166 /** \brief Set psi-point mask to be consistent with rho-point mask */
167 void update_mskp (int lev);
168
169 /** \brief compute dt from CFL considerations */
170 amrex::Real estTimeStep (int lev) const;
171
172 /** \brief Interface for advancing the data at one level by one "slow" timestep */
173 void remora_advance(int level,
174 amrex::MultiFab& cons_old, amrex::MultiFab& cons_new,
175 amrex::MultiFab& xvel_old, amrex::MultiFab& yvel_old, amrex::MultiFab& zvel_old,
176 amrex::MultiFab& xvel_new, amrex::MultiFab& yvel_new, amrex::MultiFab& zvel_new,
177 amrex::MultiFab& source,
178 const amrex::Geometry fine_geom,
179 const amrex::Real dt, const amrex::Real time);
180
181 /** \brief Make mask to zero out covered cells (for mesh refinement) */
182 amrex::MultiFab& build_fine_mask(int lev);
183
184 /** \brief main driver for writing AMReX plotfiles */
185 void WritePlotFile (int istep);
186
187 /** \brief write out particular data to an AMReX plotfile */
188 void WriteMultiLevelPlotfileWithBathymetry (const std::string &plotfilename,
189 int nlevels,
190 const amrex::Vector<const amrex::MultiFab*> &mf,
191 const amrex::Vector<const amrex::MultiFab*> &mf_nd,
192 const amrex::Vector<const amrex::MultiFab*> &mf_u,
193 const amrex::Vector<const amrex::MultiFab*> &mf_v,
194 const amrex::Vector<const amrex::MultiFab*> &mf_w,
195 const amrex::Vector<const amrex::MultiFab*> &mf_2d_rho,
196 const amrex::Vector<const amrex::MultiFab*> &mf_2d_u,
197 const amrex::Vector<const amrex::MultiFab*> &mf_2d_v,
198 const amrex::Vector<std::string> &varnames_3d,
199 const amrex::Vector<std::string> &varnames_2d_rho,
200 const amrex::Vector<std::string> &varnames_2d_u,
201 const amrex::Vector<std::string> &varnames_2d_v,
202 const amrex::Vector<amrex::Geometry>& my_geom,
203 amrex::Real time,
204 const amrex::Vector<int> &level_steps,
205 const amrex::Vector<amrex::IntVect>& rr,
206 const std::string &versionName = "HyperCLaw-V1.1",
207 const std::string &levelPrefix = "Level_",
208 const std::string &mfPrefix = "Cell",
209 const amrex::Vector<std::string>& extra_dirs = amrex::Vector<std::string>()) const;
210
211
212
213 /** \brief write out header data for an AMReX plotfile */
214 void WriteGenericPlotfileHeaderWithBathymetry (std::ostream &HeaderFile,
215 int nlevels,
216 const amrex::Vector<amrex::BoxArray> &bArray,
217 const amrex::Vector<std::string> &varnames_3d,
218 const amrex::Vector<std::string> &varnames_2d_rho,
219 const amrex::Vector<std::string> &varnames_2d_u,
220 const amrex::Vector<std::string> &varnames_2d_v,
221 const amrex::Vector<amrex::Geometry>& my_geom,
222 amrex::Real time,
223 const amrex::Vector<int> &level_steps,
224 const amrex::Vector<amrex::IntVect>& rr,
225 const std::string &versionName,
226 const std::string &levelPrefix,
227 const std::string &mfPrefix) const;
228
229 /** \brief default prefix for input file parameters */
230 std::string pp_prefix {"remora"};
231
232 /** \brief multilevel data container for last step's scalar data: temperature, salinity, passive tracer */
233 amrex::Vector<amrex::MultiFab*> cons_old;
234 /** \brief multilevel data container for last step's x velocities (u in ROMS) */
235 amrex::Vector<amrex::MultiFab*> xvel_old;
236 /** \brief multilevel data container for last step's y velocities (v in ROMS) */
237 amrex::Vector<amrex::MultiFab*> yvel_old;
238 /** \brief multilevel data container for last step's z velocities (largely unused; W stored separately) */
239 amrex::Vector<amrex::MultiFab*> zvel_old;
240
241 /** \brief multilevel data container for current step's scalar data: temperature, salinity, passive tracer */
242 amrex::Vector<amrex::MultiFab*> cons_new;
243 /** \brief multilevel data container for current step's x velocities (u in ROMS) */
244 amrex::Vector<amrex::MultiFab*> xvel_new;
245 /** \brief multilevel data container for current step's y velocities (v in ROMS) */
246 amrex::Vector<amrex::MultiFab*> yvel_new;
247 /** \brief multilevel data container for current step's z velocities (largely unused; W stored separately) */
248 amrex::Vector<amrex::MultiFab*> zvel_new;
249
250 // Program state data is represented by vectors of pointers to AMReX Multifabs.
251 // There is one pointer per level
252
253 /** \brief Bathymetry data (2D, positive valued, h in ROMS) **/
254 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_h;
255
256 /** \brief Bathymetry data on the whole domain at each potential level **/
257 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_h_full_domain;
258
259 /** \brief Width of cells in the vertical (z-) direction (3D, Hz in ROMS) */
260 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Hz;
261 /** \brief u-volume flux (3D) */
262 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Huon;
263 /** \brief v-volume flux (3D) */
264 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Hvom;
265 /** \brief u velocity RHS (3D, includes horizontal and vertical advection) */
266 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_ru;
267 /** \brief v velocity RHS (3D, includes horizontal and vertical advection) */
268 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rv;
269 /** \brief u velocity RHS (2D, includes horizontal and vertical advection) */
270 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_ru2d;
271 /** \brief v velocity RHS (2D, includes horizontal and vertical advection) */
272 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rv2d;
273 /** \brief u velocity RHS, integrated, including advection and bottom/surface stresses (2D) */
274 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rufrc;
275 /** \brief v velocity RHS, integrated, including advection and bottom/surface stresses (2D) */
276 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rvfrc;
277 /** \brief Vertical viscosity coefficient (3D) */
278 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Akv;
279 /** \brief Vertical diffusion coefficient (3D) */
280 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Akt;
281 /** \brief Harmonic viscosity defined on the psi points (corners of horizontal grid cells) */
282 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_visc2_p;
283 /** \brief Harmonic viscosity defined on the rho points (centers) */
284 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_visc2_r;
285 /** \brief Harmonic diffusivity for temperature / salinity */
286 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_diff2;
287
288 /** \brief z coordinates at rho points (cell centers) */
289 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_z_r;
290
291 /** \brief z coordinates at w points (faces between z-cells) */
292 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_z_w;
293 /** \brief Scaled vertical coordinate (range [0,1]) that transforms to z, defined at rho points (cell centers) */
294 amrex::Gpu::DeviceVector<amrex::Real> s_r;
295 /** \brief Scaled vertical coordinate (range [0,1]) that transforms to z, defined at w-points (cell faces) */
296 amrex::Gpu::DeviceVector<amrex::Real> s_w;
297
298 //amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_s_r;
299 //amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_s_w;
300 //amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Cs_r;
301 //amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Cs_w;
302
303 /** \brief Stretching coefficients at rho points */
304 amrex::Gpu::DeviceVector<amrex::Real> Cs_r;
305 /** \brief Stretching coefficients at w points */
306 amrex::Gpu::DeviceVector<amrex::Real> Cs_w;
307
308 /** \brief z coordinates at psi points (cell nodes) **/
309 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_z_phys_nd;
310
311 /** \brief Average of the free surface, zeta (2D) */
312 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Zt_avg1;
313
314 /** \brief Surface stress in the u direction */
315 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_sustr;
316 /** \brief Surface stress in the v direction */
317 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_svstr;
318
319 /** \brief Wind in the u direction, defined at rho-points */
320 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_uwind;
321 /** \brief Wind in the v direction, defined at rho-points */
322 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_vwind;
323 /** \brief Air temperature [°C], defined at rho-points */
324 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Tair;
325 /** \brief Specific humidity [kg/kg], defined at rho-points */
326 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_qair;
327 /** \brief Air pressure [mb], defined at rho-points */
328 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Pair;
329
330 /** \brief Shortwave radiation flux [W/m²], defined at rho-points */
331 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_srflx;
332 /** \brief longwave radiation */
333 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_lrflx;
334 /** \brief Downward longwave radiation */
335 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_longwave_down;
336 /** \brief latent heat flux */
337 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_lhflx;
338 /** \brief sensible heat flux */
339 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_shflx;
340
341 /** \brief Surface tracer flux; working arrays */
342 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_stflx;
343 /** \brief Surface tracer flux; input arrays */
344 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_stflux;
345 /** \brief Bottom tracer flux; working arrays */
346 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_btflx;
347 /** \brief Bottom tracer flux; input arrays */
348 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_btflux;
349
350 /** \brief precipitation rate [kg/m^2/s] */
351 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rain;
352 /** \brief evaporation rate [kg/m^2/s] */
353 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_evap;
354 /** \brief cloud cover fraction [0-1], defined at rho-points */
355 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_cloud;
356 /** \brief evaporation minus precipitation [kg/m^2/s], defined at rho-points */
357 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_EminusP;
358
359 /** \brief Linear drag coefficient [m/s], defined at rho points */
360 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rdrag;
361 /** \brief Quadratic drag coefficient [unitless], defined at rho points */
362 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rdrag2;
363 /** \brief Bottom roughness length [m], defined at rho points */
364 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_ZoBot;
365
366 /** \brief Bottom stress in the u direction */
367 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_bustr;
368 /** \brief Bottom stress in the v direction */
369 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_bvstr;
370
371 /** \brief time average of barotropic x velocity flux (2D) */
372 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_DU_avg1;
373 /** \brief correct time average of barotropic x velocity flux for coupling (2D) */
374 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_DU_avg2;
375 /** \brief time average of barotropic y velocity flux */
376 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_DV_avg1;
377 /** \brief correct time average of barotropic y velocity flux for coupling (2D) */
378 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_DV_avg2;
379 /** \brief barotropic x velocity for the RHS (2D) */
380 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rubar;
381 /** \brief barotropic y velocity for the RHS (2D) */
382 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rvbar;
383 /** \brief free surface height for the RHS (2D) */
384 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rzeta;
385 /** \brief barotropic x velocity (2D) */
386 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_ubar;
387 /** \brief barotropic y velocity (2D) */
388 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_vbar;
389 /** \brief free surface height (2D) */
390 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_zeta;
391
392 /** \brief land/sea mask at cell centers (2D) */
393 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_mskr;
394 /** \brief land/sea mask at x-faces (2D) */
395 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_msku;
396 /** \brief land/sea mask at y-faces (2D) */
397 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_mskv;
398 /** \brief land/sea mask at cell corners (2D) */
399 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_mskp;
400 /** \brief land/sea mask at cell centers, copied to all z levels (3D) */
401 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_mskr3d;
402
403 /** \brief horizontal scaling factor: 1 / dx (2D) */
404 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_pm;
405 /** \brief horizontal scaling factor: 1 / dy (2D) */
406 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_pn;
407
408 /** \brief coriolis factor (2D) */
409 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_fcor;
410
411 /** \brief x_grid on rho points (2D) */
412 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_xr;
413 /** \brief y_grid on rho points (2D) */
414 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_yr;
415
416 /** \brief x_grid on u-points (2D) */
417 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_xu;
418 /** \brief y_grid on u-points (2D) */
419 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_yu;
420
421 /** \brief x_grid on v-points (2D) */
422 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_xv;
423 /** \brief y_grid on v-points (2D) */
424 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_yv;
425
426 /** \brief x_grid on psi-points (2D) */
427 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_xp;
428 /** \brief y_grid on psi-points (2D) */
429 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_yp;
430
431 /** \brief additional scratch space for calculations on temp, salt, etc */
432 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_sstore;
433
434 /** \brief density perturbation */
435 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rhoS;
436 /** \brief vertically-averaged density */
437 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rhoA;
438 /** \brief Brunt-Vaisala frequency (3D) */
439 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_bvf;
440 /** \brief Thermal expansion coefficient (3D) */
441 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_alpha;
442 /** \brief Saline contraction coefficient (3D) */
443 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_beta;
444
445 /** \brief Weights for calculating avg1 in 2D advance */
446 amrex::Vector<amrex::Real> vec_weight1;
447 /** \brief Weights for calculating avg2 in 2D advance */
448 amrex::Vector<amrex::Real> vec_weight2;
449
450 // GLS
451 /** \brief Turbulent kinetic energy */
452 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_tke;
453 /** \brief Turbulent generic length scale */
454 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_gls;
455 /** \brief Vertical mixing turbulent length scale */
456 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Lscale;
457 /** \brief Turbulent kinetic energy vertical diffusion coefficient */
458 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Akk;
459 /** \brief Turbulent length scale vertical diffusion coefficient */
460 amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Akp;
461
462 /** \brief Climatology nudging coefficients */
463 amrex::Vector<amrex::Vector<std::unique_ptr<amrex::MultiFab>>> vec_nudg_coeff;
464
465 /** \brief advance a single level for a single time step */
466 void Advance (int lev, amrex::Real time, amrex::Real dt_lev, int iteration, int ncycle);
467
468 /** \brief Set everything up for a step on a level */
469 void setup_step(int lev, amrex::Real time, amrex::Real dt_lev);
470
471 /** \brief 3D advance on a single level */
472 void advance_3d_ml (int lev, amrex::Real dt_lev);
473
474 /** \brief 2D advance, one predictor/corrector step */
475 void advance_2d_onestep (int lev, amrex::Real dt_lev, amrex::Real dtfast_lev, int my_iif, int nfast_counter);
476
477 /** \brief Perform a 2D predictor (predictor_2d_step=True) or corrector (predictor_2d_step=False) step */
478 void advance_2d (int lev,
479 amrex::MultiFab const* mf_rhoS,
480 amrex::MultiFab const* mf_rhoA,
481 amrex::MultiFab * mf_ru2d,
482 amrex::MultiFab * mf_rv2d,
483 amrex::MultiFab * mf_rufrc,
484 amrex::MultiFab * mf_rvfrc,
485 amrex::MultiFab * mf_Zt_avg1,
486 std::unique_ptr<amrex::MultiFab>& mf_DU_avg1,
487 std::unique_ptr<amrex::MultiFab>& mf_DU_avg2,
488 std::unique_ptr<amrex::MultiFab>& mf_DV_avg1,
489 std::unique_ptr<amrex::MultiFab>& mf_DV_avg2,
490 std::unique_ptr<amrex::MultiFab>& mf_rubar,
491 std::unique_ptr<amrex::MultiFab>& mf_rvbar,
492 std::unique_ptr<amrex::MultiFab>& mf_rzeta,
493 std::unique_ptr<amrex::MultiFab>& mf_ubar,
494 std::unique_ptr<amrex::MultiFab>& mf_vbar,
495 amrex::MultiFab * mf_zeta,
496 amrex::MultiFab const* mf_h,
497 amrex::MultiFab const* mf_pm,
498 amrex::MultiFab const* mf_pn,
499 amrex::MultiFab const* mf_fcor,
500 amrex::MultiFab const* mf_visc2_p,
501 amrex::MultiFab const* mf_visc2_r,
502 amrex::MultiFab const* mf_mskr,
503 amrex::MultiFab const* mf_msku,
504 amrex::MultiFab const* mf_mskv,
505 amrex::MultiFab const* mf_mskp,
506 amrex::Real dtfast_lev,
507 bool predictor_2d_step,
508 bool first_2d_step, int my_iif,
509 int & next_indx1);
510
511 /** \brief Advance the 3D variables */
512 void advance_3d (int lev,
513 amrex::MultiFab& mf_cons,
514 amrex::MultiFab& mf_u , amrex::MultiFab& mf_v,
515 amrex::MultiFab* mf_sstore,
516 amrex::MultiFab* mf_ru , amrex::MultiFab* mf_rv,
517 std::unique_ptr<amrex::MultiFab>& mf_DU_avg1,
518 std::unique_ptr<amrex::MultiFab>& mf_DU_avg2,
519 std::unique_ptr<amrex::MultiFab>& mf_DV_avg1,
520 std::unique_ptr<amrex::MultiFab>& mf_DV_avg2,
521 std::unique_ptr<amrex::MultiFab>& mf_ubar,
522 std::unique_ptr<amrex::MultiFab>& mf_vbar,
523 std::unique_ptr<amrex::MultiFab>& mf_Akv,
524 std::unique_ptr<amrex::MultiFab>& mf_Akt,
525 std::unique_ptr<amrex::MultiFab>& mf_Hz,
526 std::unique_ptr<amrex::MultiFab>& mf_Huon,
527 std::unique_ptr<amrex::MultiFab>& mf_Hvom,
528 std::unique_ptr<amrex::MultiFab>& mf_z_w,
529 amrex::MultiFab const* mf_h,
530 amrex::MultiFab const* mf_pm,
531 amrex::MultiFab const* mf_pn,
532 amrex::MultiFab const* mf_mskr,
533 amrex::MultiFab const* mf_msku,
534 amrex::MultiFab const* mf_mskv,
535 const int N,
536 const amrex::Real dt_lev);
537
538 /** \brief Calculate bulk temperature, salinity, wind fluxes */
539 void bulk_fluxes (int lev,
540 amrex::MultiFab* mf_cons,
541 amrex::MultiFab* mf_uwind,
542 amrex::MultiFab* mf_vwind,
543 amrex::MultiFab* mf_Tair,
544 amrex::MultiFab* mf_qair,
545 amrex::MultiFab* mf_Pair,
546 amrex::MultiFab* mf_srflx,
547 amrex::MultiFab* mf_longwave_down,
548 amrex::MultiFab* mf_evap,
549 amrex::MultiFab* mf_sustr,
550 amrex::MultiFab* mf_svstr,
551 amrex::MultiFab* mf_stflux,
552 amrex::MultiFab* mf_lrflx,
553 amrex::MultiFab* mf_lhflx,
554 amrex::MultiFab* mf_shflx,
555 const int N);
556
557 /** \brief Wrapper function for prestep */
558 void prestep (int lev,
559 amrex::MultiFab& mf_uold, amrex::MultiFab& mf_vold,
560 amrex::MultiFab& mf_u, amrex::MultiFab& mf_v,
561 amrex::MultiFab* mf_ru,
562 amrex::MultiFab* mf_rv,
563 amrex::MultiFab& S_old,
564 amrex::MultiFab& S_new,
565 amrex::MultiFab& mf_W, amrex::MultiFab& mf_DC,
566 /* MF mf_FC? */
567 const amrex::MultiFab* mf_z_r,
568 const amrex::MultiFab* mf_z_w,
569 const amrex::MultiFab* mf_h,
570 const amrex::MultiFab* mf_pm,
571 const amrex::MultiFab* mf_pn,
572 const amrex::MultiFab* mf_sustr,
573 const amrex::MultiFab* mf_svstr,
574 const amrex::MultiFab* mf_bustr,
575 const amrex::MultiFab* mf_bvstr,
576 const amrex::MultiFab* mf_msku,
577 const amrex::MultiFab* mf_mskv,
578 const int iic, const int nfirst,
579 const int nnew, int nstp, int nrhs,
580 int N, const amrex::Real dt_lev);
581
582 /** \brief Prestep advection calculations for the tracers */
583 void prestep_t_advection (int lev, const amrex::Box& tbx,
584 const amrex::Box& gbx,
585 const amrex::Array4<amrex::Real >& tempold,
586 const amrex::Array4<amrex::Real >& tempcache,
587 const amrex::Array4<amrex::Real >& Hz,
588 const amrex::Array4<amrex::Real >& Huon,
589 const amrex::Array4<amrex::Real >& Hvom,
590 const amrex::Array4<amrex::Real >& W,
591 const amrex::Array4<amrex::Real >& DC,
592 const amrex::Array4<amrex::Real >& FC,
593 const amrex::Array4<amrex::Real >& sstore,
594 const amrex::Array4<amrex::Real const>& z_w,
595 const amrex::Array4<amrex::Real const>& h,
596 const amrex::Array4<amrex::Real const>& pm,
597 const amrex::Array4<amrex::Real const>& pn,
598 const amrex::Array4<amrex::Real const>& msku,
599 const amrex::Array4<amrex::Real const>& mskv,
600 const amrex::Array4<int const>& river_pos,
601 const amrex::Array4<amrex::Real const>& river_source,
602 int iic, int ntfirst, int nrhs, int N,
603 const amrex::Real dt_lev);
604
605 /** \brief RHS terms for tracer */
606 void rhs_t_3d (int lev,
607 const amrex::Box& bx,
608 const amrex::Array4<amrex::Real >& t,
609 const amrex::Array4<amrex::Real const>& tempstore,
610 const amrex::Array4<amrex::Real const>& Huon,
611 const amrex::Array4<amrex::Real const>& Hvom,
612 const amrex::Array4<amrex::Real const>& Hz,
613 const amrex::Array4<amrex::Real const>& pn,
614 const amrex::Array4<amrex::Real const>& pm,
615 const amrex::Array4<amrex::Real const>& W,
616 const amrex::Array4<amrex::Real >& FC,
617 const amrex::Array4<amrex::Real const>& mskr,
618 const amrex::Array4<amrex::Real const>& msku,
619 const amrex::Array4<amrex::Real const>& mskv,
620 const amrex::Array4<int const>& river_pos,
621 const amrex::Array4<amrex::Real const>& river_source,
622 int nrhs, int nnew, int N, const amrex::Real dt_lev);
623
624 /** \brief RHS terms for 3D momentum */
625 void rhs_uv_3d (int lev,
626 const amrex::Box& xbx,
627 const amrex::Box& ybx,
628 const amrex::Array4<amrex::Real const>& uold,
629 const amrex::Array4<amrex::Real const>& vold,
630 const amrex::Array4<amrex::Real >& ru,
631 const amrex::Array4<amrex::Real >& rv,
632 const amrex::Array4<amrex::Real >& rufrc,
633 const amrex::Array4<amrex::Real >& rvfrc,
634 const amrex::Array4<amrex::Real const>& sustr,
635 const amrex::Array4<amrex::Real const>& svstr,
636 const amrex::Array4<amrex::Real const>& bustr,
637 const amrex::Array4<amrex::Real const>& bvstr,
638 const amrex::Array4<amrex::Real const>& Huon,
639 const amrex::Array4<amrex::Real const>& Hvom,
640 const amrex::Array4<amrex::Real const>& pm,
641 const amrex::Array4<amrex::Real const>& pn,
642 const amrex::Array4<amrex::Real const>& W,
643 const amrex::Array4<amrex::Real >& FC,
644 int nrhs, int N);
645
646 /** \brief RHS terms for 2D momentum */
647 void rhs_uv_2d (int lev,
648 const amrex::Box& xbx,
649 const amrex::Box& ybx,
650 const amrex::Array4<amrex::Real const>& uold,
651 const amrex::Array4<amrex::Real const>& vold,
652 const amrex::Array4<amrex::Real >& ru,
653 const amrex::Array4<amrex::Real >& rv,
654 const amrex::Array4<amrex::Real const>& Duon,
655 const amrex::Array4<amrex::Real const>& Dvom,
656 const int nrhs);
657
658 /** \brief Wrapper around equation of state calculation */
659 void rho_eos (const amrex::Box& bx,
660 const amrex::Array4<amrex::Real const>& state,
661 const amrex::Array4<amrex::Real >& rho,
662 const amrex::Array4<amrex::Real >& rhoA,
663 const amrex::Array4<amrex::Real >& rhoS,
664 const amrex::Array4<amrex::Real >& bvf,
665 const amrex::Array4<amrex::Real >& alpha,
666 const amrex::Array4<amrex::Real >& beta,
667 const amrex::Array4<amrex::Real const>& Hz,
668 const amrex::Array4<amrex::Real const>& z_w,
669 const amrex::Array4<amrex::Real const>& z_r,
670 const amrex::Array4<amrex::Real const>& h,
671 const amrex::Array4<amrex::Real const>& mskr,
672 const int N);
673
674 /** \brief Calculate density and related quantities from linear equation of state */
675 void lin_eos (const amrex::Box& bx,
676 const amrex::Array4<amrex::Real const>& state,
677 const amrex::Array4<amrex::Real >& rho,
678 const amrex::Array4<amrex::Real >& rhoA,
679 const amrex::Array4<amrex::Real >& rhoS,
680 const amrex::Array4<amrex::Real >& bvf,
681 const amrex::Array4<amrex::Real const>& Hz,
682 const amrex::Array4<amrex::Real const>& z_w,
683 const amrex::Array4<amrex::Real const>& z_r,
684 const amrex::Array4<amrex::Real const>& h,
685 const amrex::Array4<amrex::Real const>& mskr,
686 const int N);
687
688 /** \brief Calculate density and related quantities from nonlinear equation of state */
689 void nonlin_eos (const amrex::Box& bx,
690 const amrex::Array4<amrex::Real const>& state,
691 const amrex::Array4<amrex::Real >& rho,
692 const amrex::Array4<amrex::Real >& rhoA,
693 const amrex::Array4<amrex::Real >& rhoS,
694 const amrex::Array4<amrex::Real >& bvf,
695 const amrex::Array4<amrex::Real >& alpha,
696 const amrex::Array4<amrex::Real >& beta,
697 const amrex::Array4<amrex::Real const>& Hz,
698 const amrex::Array4<amrex::Real const>& z_w,
699 const amrex::Array4<amrex::Real const>& z_r,
700 const amrex::Array4<amrex::Real const>& h,
701 const amrex::Array4<amrex::Real const>& mskr,
702 const int N);
703
704 /** \brief Calculate pressure gradient */
705 void prsgrd (const amrex::Box& bx,
706 const amrex::Box& gbx,
707 const amrex::Box& utbx,
708 const amrex::Box& vtbx,
709 const amrex::Array4<amrex::Real >& ru,
710 const amrex::Array4<amrex::Real >& rv,
711 const amrex::Array4<amrex::Real const>& pn,
712 const amrex::Array4<amrex::Real const>& pm,
713 const amrex::Array4<amrex::Real const>& rho,
714 const amrex::Array4<amrex::Real >& FC,
715 const amrex::Array4<amrex::Real const>& Hz,
716 const amrex::Array4<amrex::Real const>& z_r,
717 const amrex::Array4<amrex::Real const>& z_w,
718 const amrex::Array4<amrex::Real const>& msku,
719 const amrex::Array4<amrex::Real const>& mskv,
720 const int nrhs, const int N);
721
722 /** \brief Update velocities or tracers with diffusion/viscosity
723 * as the last part of the prestep */
724 void prestep_diffusion (const amrex::Box& bx,
725 const amrex::Box& gbx,
726 const int ioff, const int joff,
727 const amrex::Array4<amrex::Real >& vel,
728 const amrex::Array4<amrex::Real const>& vel_old,
729 const amrex::Array4<amrex::Real >& rvel,
730 const amrex::Array4<amrex::Real const>& Hz,
731 const amrex::Array4<amrex::Real const>& Akv,
732 const amrex::Array4<amrex::Real >& FC,
733 const amrex::Array4<amrex::Real const>& sstr,
734 const amrex::Array4<amrex::Real const>& bstr,
735 const amrex::Array4<amrex::Real const>& z_r,
736 const amrex::Array4<amrex::Real const>& pm,
737 const amrex::Array4<amrex::Real const>& pn,
738 const int iic, const int ntfirst, const int nnew, int nstp, int nrhs, int N,
739 const amrex::Real lambda, const amrex::Real dt_lev);
740
741 /** \brief Calculate effects of vertical viscosity or diffusivity */
742 void vert_visc_3d (const amrex::Box& bx,
743 const int ioff, const int joff,
744 const amrex::Array4<amrex::Real >& phi,
745 const amrex::Array4<amrex::Real const>& Hz,
746 const amrex::Array4<amrex::Real >& Hzk,
747 const amrex::Array4<amrex::Real >& AK,
748 const amrex::Array4<amrex::Real const>& Akv,
749 const amrex::Array4<amrex::Real >& BC,
750 const amrex::Array4<amrex::Real >& DC,
751 const amrex::Array4<amrex::Real >& FC,
752 const amrex::Array4<amrex::Real >& CF,
753 const int nnew, const int N,
754 const amrex::Real dt_lev);
755
756 /** \brief Correct mass flux */
757 void update_massflux_3d (int lev, const amrex::Box& bx,
758 const int ioff, const int joff,
759 const amrex::Array4<amrex::Real >& phi,
760 const amrex::Array4<amrex::Real >& phibar,
761 const amrex::Array4<amrex::Real >& Hphi,
762 const amrex::Array4<amrex::Real const>& Hz,
763 const amrex::Array4<amrex::Real const>& pm_or_pn,
764 const amrex::Array4<amrex::Real const>& Dphi1,
765 const amrex::Array4<amrex::Real const>& Dphi2,
766 const amrex::Array4<amrex::Real >& DC,
767 const amrex::Array4<amrex::Real >& FC,
768 const amrex::Array4<amrex::Real const>& msk,
769 const int nnew);
770
771 /** \brief Adjust 3D momentum variables based on vertical mean momentum */
772 void vert_mean_3d (const amrex::Box& bx,
773 const int ioff, const int joff,
774 const amrex::Array4<amrex::Real >& phi,
775 const amrex::Array4<amrex::Real const>& Hz,
776 const amrex::Array4<amrex::Real const>& Dphi_avg1,
777 const amrex::Array4<amrex::Real >& DC,
778 const amrex::Array4<amrex::Real >& CF,
779 const amrex::Array4<amrex::Real const>& pm_or_pn,
780 const amrex::Array4<amrex::Real const>& msk,
781 const int nnew, const int N);
782
783 /** \brief Harmonic viscosity */
784 void uv3dmix (const amrex::Box& xbx,
785 const amrex::Box& ybx,
786 const amrex::Array4<amrex::Real >& u,
787 const amrex::Array4<amrex::Real >& v,
788 const amrex::Array4<amrex::Real const>& uold,
789 const amrex::Array4<amrex::Real const>& vold,
790 const amrex::Array4<amrex::Real >& rufrc,
791 const amrex::Array4<amrex::Real >& rvfrc,
792 const amrex::Array4<amrex::Real const>& visc2_p,
793 const amrex::Array4<amrex::Real const>& visc2_r,
794 const amrex::Array4<amrex::Real const>& Hz,
795 const amrex::Array4<amrex::Real const>& pm,
796 const amrex::Array4<amrex::Real const>& pn,
797 const amrex::Array4<amrex::Real const>& mskp,
798 int nrhs, int nnew,
799 const amrex::Real dt_lev);
800
801 /** \brief Wrapper for harmonic diffusivity for tracers */
802 void t3dmix2 (const amrex::Box& bx,
803 const amrex::Array4<amrex::Real >& state,
804 const amrex::Array4<amrex::Real >& state_rhs,
805 const amrex::Array4<amrex::Real const>& diff2,
806 const amrex::Array4<amrex::Real const>& Hz,
807 const amrex::Array4<amrex::Real const>& z_r,
808 const amrex::Array4<amrex::Real const>& pm,
809 const amrex::Array4<amrex::Real const>& pn,
810 const amrex::Array4<amrex::Real const>& msku,
811 const amrex::Array4<amrex::Real const>& mskv,
812 const amrex::Real dt_lev,
813 const int ncomp, const int N);
814
815 /** \brief Harmonic diffusivity for tracers along S-coordinate level surfaces */
816 void t3dmix2_s(const amrex::Box& bx,
817 const amrex::Array4<amrex::Real >& state,
818 const amrex::Array4<amrex::Real >& state_rhs,
819 const amrex::Array4<amrex::Real const>& diff2,
820 const amrex::Array4<amrex::Real const>& Hz,
821 const amrex::Array4<amrex::Real const>& pm,
822 const amrex::Array4<amrex::Real const>& pn,
823 const amrex::Array4<amrex::Real const>& msku,
824 const amrex::Array4<amrex::Real const>& mskv,
825 const amrex::Real dt_lev, const int ncomp);
826
827 /** \brief Harmonic diffusivity for tracers along geopotential surfaces */
828 void t3dmix2_geo (const amrex::Box& bx,
829 const amrex::Array4<amrex::Real >& state,
830 const amrex::Array4<amrex::Real >& state_rhs,
831 const amrex::Array4<amrex::Real const>& diff2,
832 const amrex::Array4<amrex::Real const>& Hz,
833 const amrex::Array4<amrex::Real const>& z_r,
834 const amrex::Array4<amrex::Real const>& pm,
835 const amrex::Array4<amrex::Real const>& pn,
836 const amrex::Array4<amrex::Real const>& msku,
837 const amrex::Array4<amrex::Real const>& mskv,
838 const amrex::Real dt_lev,
839 const int ncomp, const int N);
840
841 /** \brief Calculate Coriolis terms */
842 void coriolis (const amrex::Box& xbx,
843 const amrex::Box& ybx,
844 const amrex::Array4<amrex::Real const>& uold,
845 const amrex::Array4<amrex::Real const>& vold,
846 const amrex::Array4<amrex::Real >& ru,
847 const amrex::Array4<amrex::Real >& rv,
848 const amrex::Array4<amrex::Real const>& Hz,
849 const amrex::Array4<amrex::Real const>& fomn,
850 int nrhs, int nr);
851
852 /** \brief Apply climatology nudging */
853 void apply_clim_nudg (const amrex::Box& bx,
854 int ioff, int joff,
855 const amrex::Array4<amrex::Real >& var,
856 const amrex::Array4<amrex::Real const>& var_old,
857 const amrex::Array4<amrex::Real const>& var_clim,
858 const amrex::Array4<amrex::Real const>& clim_coeff,
859 const amrex::Array4<amrex::Real const>& Hz,
860 const amrex::Array4<amrex::Real const>& pm,
861 const amrex::Array4<amrex::Real const>& pn,
862 const amrex::Real dt_lev = amrex::Real(0.0));
863
864 /** \brief Set 2D momentum arrays from 3D momentum */
865 void set_2darrays (int lev);
866
867 /** \brief Set Zt_avg1 to zeta */
868 void set_zeta_average (int lev);
869
870 /** \brief Initialize zeta from file or analytic */
871 void set_zeta (int lev);
872
873 /** \brief Initialize bathymetry from file or analytic */
874 void set_bathymetry (int lev);
875
876 /** \brief Copy over bathymetry data that has been averaged down from high resolution
877 * input netcdf file */
878 void set_bathymetry_averaged_down (int lev);
879
880 /** \brief Initialize Coriolis factor from file or analytic */
881 void set_coriolis (int lev);
882
883 /** \brief Initialize land-sea masks from file or analytic */
884 void set_masks (int lev);
885
886 /** \brief Calculate vertical stretched coordinates */
887 void stretch_transform (int lev);
888
889 /** \brief Initialize vertical mixing coefficients from file or analytic */
890 void init_set_vmix (int lev);
891 /** \brief Set vertical mixing coefficients from analytic */
892 void set_analytic_vmix (int lev);
893 /** \brief Initialize GLS variables */
894 void init_gls_vmix (int lev, SolverChoice solver_choice);
895 /** \brief Prestep for GLS calculation */
896 void gls_prestep (int lev, amrex::MultiFab* mf_gls, amrex::MultiFab* mf_tke,
897 amrex::MultiFab& mf_W,
898 amrex::MultiFab* mf_msku, amrex::MultiFab* mf_mskv,
899 const int nstp, const int nnew, const int iic, const int ntfirst,
900 const int N, const amrex::Real dt_lev);
901 /** \brief Corrector step for GLS calculation */
902 void gls_corrector (int lev, amrex::MultiFab* mf_gls, amrex::MultiFab* mf_tke,
903 amrex::MultiFab& mf_W, amrex::MultiFab* mf_Akv, amrex::MultiFab* mf_Akt,
904 amrex::MultiFab* mf_Akk, amrex::MultiFab* mf_Akp,
905 amrex::MultiFab* mf_mskr,
906 amrex::MultiFab* mf_msku, amrex::MultiFab* mf_mskv,
907 const int nstp, const int nnew,
908 const int N, const amrex::Real dt_lev);
909
910 /** \brief Scale RHS momentum variables by 1/cell area, needed before FillPatch to different levels*/
911 void scale_rhs_vars ();
912 /** \brief Scale RHS momentum variables by cell area, needed after FillPatch to different levels*/
913 void scale_rhs_vars_inv ();
914
915 /** \brief Initialize or calculate surface momentum flux from file or analytic */
916 void set_smflux (int lev);
917
918 /** \brief Initialize or calculate wind speed from file or analytic */
919 void set_wind (int lev);
920
921 /** \brief Initialize horizontal mixing coefficients */
922 void set_hmixcoef (int lev);
923
924 /** \brief Initialize flat bathymetry to value from problo */
925 void init_flat_bathymetry (int lev);
926
927 /** \brief Initialize or calculate bottom drag */
928 void set_drag (int lev);
929
930 /** \brief Set weights for averaging 3D variables to 2D */
931 void set_weights (int lev);
932
933 /** \brief Fill the physical boundary conditions for cell-centered velocity (diagnostic only) */
934 void FillBdyCCVels (int lev, amrex::MultiFab& mf_cc_vel);
935
936 /** \brief Fill a new MultiFab by copying in phi from valid region and filling ghost cells */
937 void FillPatch (int lev, amrex::Real time,
938 amrex::MultiFab& mf_to_be_filled,
939 amrex::Vector<amrex::MultiFab*> const& mfs,
940 const int bccomp,
941 const int bdy_var_type = BdyVars::null,
942 const int icomp=0,
943 const bool fill_all=true,
944 const bool fill_set=true,
945 const int n_not_fill=0,
946 const int icomp_calc=0,
947 const amrex::Real dt = amrex::Real(0.0),
948 const amrex::MultiFab& mf_calc = amrex::MultiFab());
949
950 /** \brief Fill a new MultiFab by copying in phi from valid region and filling ghost cells without applying boundary conditions */
951 void FillPatchNoBC (int lev, amrex::Real time,
952 amrex::MultiFab& mf_to_be_filled,
953 amrex::Vector<amrex::MultiFab*> const& mfs,
954 const int bdy_var_type = BdyVars::null,
955 const int icomp=0,
956 const bool fill_all=true,
957 const bool fill_set=true);
958 /** \brief Fill boundary data from netcdf file */
959 void fill_from_bdyfiles (int lev,
960 amrex::MultiFab& mf_to_fill,
961 const amrex::MultiFab& mf_mask,
962 const amrex::Real time,
963 const int bccomp,
964 const int bdy_var_type,
965 const int icomp_to_fill,
966 const int icomp_calc = 0,
967 const amrex::MultiFab& mf_calc = amrex::MultiFab(),
968 const amrex::Real = amrex::Real(0.0));
969
970 /** \brief fill an entire multifab by interpolating from the coarser level
971 * using the piecewise constant interpolater */
972 void FillCoarsePatchPC (int lev, amrex::Real time, amrex::MultiFab* mf_fine,
973 amrex::MultiFab* mf_crse,
974 const int bccomp,
975 const int bdy_var_type = BdyVars::null,
976 const int icomp = 0,
977 const bool fill_all = true,
978 const int n_not_fill=0,
979 const int icomp_calc=0,
980 const amrex::Real dt = amrex::Real(0.0),
981 const amrex::MultiFab& mf_calc = amrex::MultiFab());
982
983 /** \brief fill an entire multifab by interpolating from the coarser level */
984 void FillCoarsePatch (int lev, amrex::Real time, amrex::MultiFab* mf_fine,
985 amrex::MultiFab* mf_crse,
986 const int bccomp,
987 const int bdy_var_type = BdyVars::null,
988 const int icomp = 0,
989 const bool fill_all = true,
990 const int n_not_fill=0,
991 const int icomp_calc=0,
992 const amrex::Real dt = amrex::Real(0.0),
993 const amrex::MultiFab& mf_calc = amrex::MultiFab());
994
995 /** \brief fill an entire multifab by interpolating from the coarser level,
996 * explicitly specifying interpolator to use */
997 void FillCoarsePatchMap (int lev, amrex::Real time, amrex::MultiFab* mf_fine,
998 amrex::MultiFab* mf_crse,
999 const int bccomp,
1000 const int bdy_var_type = BdyVars::null,
1001 const int icomp = 0,
1002 const bool fill_all = true,
1003 const int n_not_fill=0,
1004 const int icomp_calc=0,
1005 const amrex::Real dt = amrex::Real(0.0),
1006 const amrex::MultiFab& mf_calc = amrex::MultiFab(),
1007 amrex::Interpolater* mapper = nullptr);
1008
1009 /** \brief initialize and calculate stretch coefficients */
1010 void init_stretch_coeffs ();
1011 /** \brief calculate vertical stretch coefficients */
1012 void calc_stretch_coeffs ();
1013
1014 /** \brief Calculate Coriolis parameters from beta plane parametrization */
1015 void init_beta_plane_coriolis (int lev);
1016
1017 /** \brief Full domain bathymetry data initialization from analytic */
1019
1020 /** \brief Problem initialization from NetCDF file */
1021 void init_data_from_netcdf (int lev);
1022 /** \brief Boundary data initialization from NetCDF file */
1023 void init_bdry_from_netcdf (int lev);
1024 /** \brief Mask data initialization from NetCDF file */
1025 void init_masks_from_netcdf (int lev);
1026 /** \brief Bathymetry data initialization from NetCDF file */
1027 void init_bathymetry_from_netcdf (int lev);
1028 /** \brief Full domain bathymetry data initialization from NetCDF file */
1030 /** \brief Sea-surface height data initialization from NetCDF file */
1031 void init_zeta_from_netcdf (int lev);
1032 /** \brief Coriolis parameter data initialization from NetCDF file */
1033 void init_coriolis_from_netcdf (int lev);
1034 /** \brief Climatology nudging coefficient initialization from NetCDF file */
1035 void init_clim_nudg_coeff_from_netcdf (int lev);
1036 /** \brief Wrapper to initialize climatology nudging coefficient */
1037 void init_clim_nudg_coeff (int lev);
1038 void init_riv_pos_from_netcdf (int lev);
1039
1040 /** \brief Allocate multifabs for storing full-domain bathymetry data */
1042
1043 /** \brief Convert data in a multifab from inverse days to inverse seconds */
1044 void convert_inv_days_to_inv_s (amrex::MultiFab*);
1045
1046 /** \brief Mask data arrays before writing output */
1047 void mask_arrays_for_write (int lev, amrex::Real fill_value, amrex::Real fill_where);
1048
1049 /** \brief Average down from level lev+1 to lev in mf, including grow cells */
1050 void average_down_with_grow_cells(int lev, amrex::Vector<std::unique_ptr<amrex::MultiFab>>& mf);
1051
1052 int xvel_bc () const noexcept { return BCVars::xvel_bc(ncons); }
1053 int yvel_bc () const noexcept { return BCVars::yvel_bc(ncons); }
1054 int zvel_bc () const noexcept { return BCVars::zvel_bc(ncons); }
1055 int ubar_bc () const noexcept { return BCVars::ubar_bc(ncons); }
1056 int vbar_bc () const noexcept { return BCVars::vbar_bc(ncons); }
1057 int zeta_bc () const noexcept { return BCVars::zeta_bc(ncons); }
1058 int tke_bc () const noexcept { return BCVars::tke_bc(ncons); }
1060 int foextrap_bc () const noexcept { return BCVars::foextrap_bc(ncons); }
1061 int u2d_simple_bc () const noexcept { return BCVars::u2d_simple_bc(ncons); }
1062 int v2d_simple_bc () const noexcept { return BCVars::v2d_simple_bc(ncons); }
1063 int num_bc_vars () const noexcept { return BCVars::NumTypes(ncons); }
1064
1065
1066#ifdef REMORA_USE_NETCDF
1068#endif
1069
1070private:
1071
1072 ///////////////////////////
1073 // private member functions
1074 ///////////////////////////
1075
1076 /** \brief read in some parameters from inputs file */
1077 void ReadParameters();
1078
1079 /** \brief Build runtime scalar names after nscalar is known */
1080 void init_scalar_metadata();
1081
1082 /** \brief set covered coarse cells to be the average of overlying fine cells */
1083 void AverageDown ();
1084
1085 /** \brief Read in boundary parameters from input file and set up data structures */
1086 void init_bcs();
1087
1088 /** \brief Initialize initial problem data from analytic functions */
1089 void init_analytic(int lev);
1090
1091 /** \brief Allocate MultiFabs for state and evolution variables */
1092 void init_stuff (int lev, const amrex::BoxArray& ba, const amrex::DistributionMapping& dm);
1093
1094 /** \brief Allocate MultiFabs for masks */
1095 void init_masks (int lev, const amrex::BoxArray& ba, const amrex::DistributionMapping& dm);
1096
1097 /** \brief Resize variable containers to accommodate data on levels 0 to max_lev */
1098 void resize_stuff (int lev);
1099
1100 /** \brief more flexible version of AverageDown() that lets you average down across multiple levels */
1101 void AverageDownTo (int crse_lev);
1102
1103 /** \brief Construct FillPatchers */
1104 void Construct_REMORAFillPatchers (int lev);
1105
1106 /** \brief Define FillPatchers */
1107 void Define_REMORAFillPatchers (int lev);
1108
1109 /** \brief utility to copy in data from old and/or new state into another multifab */
1110 TimeInterpolatedData GetDataAtTime (int lev, amrex::Real time);
1111
1112 /** \brief advance a level by dt, includes a recursive call for finer levels */
1113 void timeStep (int lev, amrex::Real time, int iteration);
1114
1115 /** \brief advance all levels by dt, loops over finer levels */
1116 void timeStepML (amrex::Real time, int iteration);
1117
1118 /** \brief a wrapper for estTimeStep() */
1119 void ComputeDt ();
1120
1121 /** \brief get plotfile name */
1122 std::string PlotFileName (int lev) const;
1123
1124 /* set which variables and derived quantities go into plotfiles */
1125 void set2DPlotVariables (const std::string& pp_plot_var_names_2d);
1126 void set3DPlotVariables (const std::string& pp_plot_var_names_3d);
1127
1128 /* append variables to plot */
1129 void append2DPlotVariables (const std::string& pp_plot_var_names_2d);
1130 void append3DPlotVariables (const std::string& pp_plot_var_names_3d);
1131
1132#ifdef REMORA_USE_NETCDF
1133 /** \brief Write plotfile using NetCDF (wrapper) */
1134 void WriteNCPlotFile (int istep, amrex::MultiFab const* plotMF);
1135
1136 /** \brief Write a particular NetCDF plotfile */
1137 void WriteNCPlotFile_which (int lev, int which_subdomain,
1138 amrex::MultiFab const* plotMF,
1139 bool write_header, ncutils::NCFile& ncf,
1140 bool is_history);
1141
1142 /** \brief Write MultiFab in NetCDF format */
1143 void WriteNCMultiFab (const amrex::FabArray<amrex::FArrayBox> &fab,
1144 const std::string& name,
1145 bool set_ghost = false) const;
1146
1147 /** \brief Read MultiFab in NetCDF format */
1148 void ReadNCMultiFab (amrex::FabArray<amrex::FArrayBox> &fab,
1149 const std::string &name,
1150 int coordinatorProc = amrex::ParallelDescriptor::IOProcessorNumber(),
1151 int allow_empty_mf = 0);
1152
1153 /** \brief Start time in the time series of boundary data */
1154 amrex::Real start_bdy_time;
1155 /** \brief Interval between boundary data times */
1157
1158 /** \brief Whether to output NetCDF files as a single history file with several time steps */
1160
1161 /** \brief Data container for u-component surface momentum flux read from file */
1162 std::unique_ptr<NCTimeSeries> sustr_data_from_file;
1163 /** \brief Data container for v-component surface momentum flux read from file */
1164 std::unique_ptr<NCTimeSeries> svstr_data_from_file;
1165 /** \brief Data container for u-direction wind read from file */
1166 std::unique_ptr<NCTimeSeries> Uwind_data_from_file;
1167 /** \brief Data container for v-direction wind read from file */
1168 std::unique_ptr<NCTimeSeries> Vwind_data_from_file;
1169 /** \brief Data container for air temperature read from file */
1170 std::unique_ptr<NCTimeSeries> Tair_data_from_file;
1171 /** \brief Data container for specific humidity read from file */
1172 std::unique_ptr<NCTimeSeries> qair_data_from_file;
1173 /** \brief Data container for air pressure read from file */
1174 std::unique_ptr<NCTimeSeries> Pair_data_from_file;
1175 /** \brief Data container for shortwave radiation flux read from file */
1176 std::unique_ptr<NCTimeSeries> srflx_data_from_file;
1177 /** \brief Data container for downward longwave radiation flux read from file */
1178 std::unique_ptr<NCTimeSeries> longwave_down_data_from_file;
1179 /** \brief Data container for precipitation rate read from file */
1180 std::unique_ptr<NCTimeSeries> rain_data_from_file;
1181 /** \brief Data container for cloud cover fraction read from file */
1182 std::unique_ptr<NCTimeSeries> cloud_data_from_file;
1183 /** \brief Data container for evaporation minus precipitation read from file */
1184 std::unique_ptr<NCTimeSeries> EminusP_data_from_file;
1185
1186 /** \brief Data container for ubar climatology data read from file */
1187 std::unique_ptr<NCTimeSeries> ubar_clim_data_from_file;
1188 /** \brief Data container for vbar climatology data read from file */
1189 std::unique_ptr<NCTimeSeries> vbar_clim_data_from_file;
1190 /** \brief Data container for u-velocity climatology data read from file */
1191 std::unique_ptr<NCTimeSeries> u_clim_data_from_file;
1192 /** \brief Data container for v-velocity climatology data read from file */
1193 std::unique_ptr<NCTimeSeries> v_clim_data_from_file;
1194 /** \brief Data container for temperature climatology data read from file */
1195 std::unique_ptr<NCTimeSeries> temp_clim_data_from_file;
1196 /** \brief Data container for salinity climatology data read from file */
1197 std::unique_ptr<NCTimeSeries> salt_clim_data_from_file;
1198
1199 /** \brief Vector of data containers for scalar data in rivers */
1200 amrex::Vector<std::unique_ptr<NCTimeSeriesRiver>> river_source_cons;
1201 /** \brief Data container for momentum transport in rivers */
1202 std::unique_ptr<NCTimeSeriesRiver> river_source_transport;
1203 /** \brief Data container for vertically integrated momentum transport in rivers */
1204 std::unique_ptr<NCTimeSeriesRiver> river_source_transportbar;
1205
1206 /** \brief Vector over BdyVars of boundary series data containers */
1207 amrex::Vector<amrex::Vector<std::unique_ptr<NCTimeSeriesBoundary>>> boundary_series;
1208
1209#endif // REMORA_USE_NETCDF
1210 /** \brief iMultiFab for river positions; contents are indices of rivers */
1211 amrex::Vector<std::unique_ptr<amrex::iMultiFab>> vec_river_position;
1212 /** \brief Vector over rivers of river direction: 0: u-face; 1: v-face; 2: w-face */
1213 amrex::Gpu::DeviceVector<int> river_direction;
1214
1215
1216#ifdef REMORA_USE_MOAB
1217 /** \brief Initialize connection to MOAB */
1218 void InitMOABMesh();
1219#endif
1220 /** \brief Nudging width at coarse-fine interface */
1221 int cf_width{0};
1222 /** \brief Width for fixing values at coarse-fine interface */
1224 // Fillpatcher classes for coarse-fine boundaries
1225 /** \brief Vector over levels of FillPatchers for scalars */
1226 amrex::Vector<REMORAFillPatcher> FPr_c;
1227 /** \brief Vector over levels of FillPatchers for u (3D) */
1228 amrex::Vector<REMORAFillPatcher> FPr_u;
1229 /** \brief Vector over levels of FillPatchers for v (3D) */
1230 amrex::Vector<REMORAFillPatcher> FPr_v;
1231 /** \brief Vector over levels of FillPatchers for w */
1232 amrex::Vector<REMORAFillPatcher> FPr_w;
1233
1234 // Fillpatcher classes for 2d velocity variables
1235 /** \brief Vector over levels of FillPatchers for ubar (2D) */
1236 amrex::Vector<REMORAFillPatcher> FPr_ubar;
1237 /** \brief Vector over levels of FillPatchers for vbar (2D) */
1238 amrex::Vector<REMORAFillPatcher> FPr_vbar;
1239
1240 /** \brief write checkpoint file to disk */
1241 void WriteCheckpointFile ();
1242
1243 /** \brief read checkpoint file from disk */
1244 void ReadCheckpointFile ();
1245
1246 /** \brief Read the file passed to remora.restart and use it as an initial condition for the current simulation */
1248
1249 /** \brief Initialize the new-time data at a level from the initial_data MultiFab */
1250 void InitializeLevelFromData (int lev, const amrex::MultiFab& initial_data);
1251
1252 /** \brief utility to skip to next line in Header */
1253 static void GotoNextLine (std::istream& is);
1254
1255 ////////////////
1256 // private data members
1257
1258 /** \brief Pointer to container of analytical functions for problem definition */
1259 std::unique_ptr<ProblemBase> prob = nullptr;
1260
1261 /** \brief how many boxes specified at each level by tagging criteria */
1262 amrex::Vector<int> num_boxes_at_level;
1263 /** \brief how many netcdf input files specified at each level */
1264 amrex::Vector<int> num_files_at_level;
1265 /** \brief the boxes specified at each level by tagging criteria */
1266 amrex::Vector<amrex::Vector<amrex::Box>> boxes_at_level;
1267
1268 /** \brief which step? */
1269 amrex::Vector<int> istep;
1270 /** \brief How many substeps on each level? */
1271 amrex::Vector<int> nsubsteps;
1272 /** \brief new time at each level */
1273 amrex::Vector<amrex::Real> t_new;
1274 /** \brief old time at each level */
1275 amrex::Vector<amrex::Real> t_old;
1276 /** \brief time step at each level */
1277 amrex::Vector<amrex::Real> dt;
1278
1279 /** \brief whether to set boundary conditions by variable rather than just by side */
1281
1282 /** \brief Vector (over level) of functors to apply physical boundary conditions */
1283 amrex::Vector<std::unique_ptr<REMORAPhysBCFunct>> physbcs;
1284
1285 /** \brief array of flux registers for refluxing in multilevel */
1286 amrex::Vector<std::unique_ptr<amrex::YAFluxRegister>> advflux_reg;
1287
1288 ////////////////
1289 // boundary data members
1290
1291 // A BCRec is essentially a 2*DIM integer array storing the boundary
1292 // condition type at each lo/hi walls in each direction. We have one BCRec
1293 // for each component of the cell-centered variables and each velocity component.
1294 /** \brief vector (over BCVars) of BCRecs */
1295 amrex::Vector <amrex::BCRec> domain_bcs_type;
1296 /** \brief GPU vector (over BCVars) of BCRecs */
1297 amrex::Gpu::DeviceVector<amrex::BCRec> domain_bcs_type_d;
1298
1299 /** \brief Array of strings describing domain boundary conditions */
1300 amrex::Array<std::string,2*AMREX_SPACEDIM> domain_bc_type;
1301
1302 /** \brief Array holding the Dirichlet values at walls which need them */
1303 amrex::Vector<amrex::GpuArray<amrex::Real, AMREX_SPACEDIM*2>> m_bc_extdir_vals;
1304
1305 /** \brief Array holding the "physical" boundary condition types (e.g. "inflow") */
1306 amrex::Vector<amrex::GpuArray<REMORA_BC, AMREX_SPACEDIM*2>> phys_bc_type;
1307
1308 /** \brief These are flags that indicate whether we need to read in boundary data from file */
1309 amrex::GpuArray<amrex::GpuArray<bool, AMREX_SPACEDIM*2>,BdyVars::NumTypes+1> phys_bc_need_data;
1310
1311 /** \brief Container to connect boundary data being read in boundary condition containers.
1312 *
1313 * Needed for phys_bc_need_data to work correctly
1314 */
1315 amrex::Vector<int> bdy_index;
1316
1317 /** \brief Step when we last output a plotfile */
1319 /** \brief Simulation time when we last output a plotfile */
1321
1322 /** \brief Step when we last output a checkpoint file */
1324 /** \brief Simulation time when we last output a checkpoint file */
1326 /** \brief Whether to output a plotfile on restart from checkpoint */
1328
1329 ////////////////
1330 // runtime parameters
1331
1332 /** \brief maximum number of steps */
1333 int max_step = std::numeric_limits<int>::max();
1334 /** \brief Time to stop */
1335 amrex::Real stop_time = std::numeric_limits<amrex::Real>::max();
1336
1337 /** \brief Time of the start of the simulation, in seconds */
1338 amrex::Real start_time = 0.0;
1339
1340 /** \brief If set, restart from this checkpoint file */
1341 std::string restart_chkfile = "";
1342
1343 /** \brief Number of passive scalars carried in the state */
1344 int nscalar = 1;
1345 /** \brief Number of conserved scalars in the state (temperature + salt + passive scalars) */
1347
1348 // Time step controls
1349 /** \brief CFL condition */
1350 static amrex::Real cfl;
1351 /** \brief Fraction maximum change in subsequent time steps */
1352 static amrex::Real change_max;
1353 /** \brief User specified fixed baroclinic time step */
1354 static amrex::Real fixed_dt;
1355 /** \brief User specified fixed barotropic time step */
1356 static amrex::Real fixed_fast_dt;
1357 /** \brief User specified, number of barotropic steps per baroclinic step */
1359 /** \brief Number of fast steps to take */
1361
1362 /** \brief Whether to substep fine levels in time */
1363 int do_substep = 0;
1364
1365 /** \brief how often each level regrids the higher levels of refinement (after a level advances that many time steps) */
1366 int regrid_int = 2;
1367
1368 // I/O controls
1369 /** \brief Plotfile prefix */
1370 std::string plot_file_name {"plt_"};
1371 /** \brief Plotfile output interval in iterations */
1372 int plot_int = -1;
1373 /** \brief Plotfile output interval in seconds */
1374 amrex::Real plot_int_time = amrex::Real(-1.0);
1375 /** \brief Checkpoint file prefix */
1376 std::string check_file {"chk"};
1377 /** \brief Checkpoint output interval in iterations */
1378 int check_int = -1;
1379 /** \brief Checkpoint output interval in seconds */
1380 amrex::Real check_int_time = amrex::Real(-1.0);
1381
1382 /** \brief Whether to chunk netcdf history file */
1384 /** \brief Number of time steps per netcdf history file.
1385 *
1386 -1 is the default and means code will automatically compute
1387 number of steps per file based on grid size to keep files < 2 GB */
1389 /** \brief Counter for which time index we are writing to in the netcdf history file */
1391
1392 /** \brief Names of 3D variables to output to AMReX plotfile */
1393 amrex::Vector<std::string> plot_var_names_3d;
1394 /** \brief Names of 2D variables to output to AMReX plotfile */
1395 amrex::Vector<std::string> plot_var_names_2d;
1396 /** \brief Names of scalars for plotfile output */
1397 amrex::Vector<std::string> cons_names {"temp", "salt", "tracer"};
1398 // Note that the order of variable names here must match the order in PlotFile.cpp
1399 /** \brief Names of derived fields for plotfiles */
1400 const amrex::Vector<std::string> derived_names {"vorticity"};
1401
1402 /** \brief Container for algorithmic choices */
1404
1405 /** \brief whether plotfile variables should be expanded to a uniform refinement ratio */
1407
1408 /** \brief fill value for masked arrays in amrex plotfiles */
1409 amrex::Real plotfile_fill_value = amrex::Real(0.0);
1410 /** \brief fill value for masked arrays in netcdf output */
1411 amrex::Real netcdf_fill_value = amrex::Real(1.0e37);
1412
1413
1414#ifdef REMORA_USE_PARTICLES
1415 /** \brief Particle container with all particle species */
1416 ParticleData particleData;
1417
1418 // variables and functions for tracers particles
1419 /** \brief tracer particles that advect with flow */
1420 bool m_use_tracer_particles;
1421 /** \brief tracer particles that fall with gravity */
1422 bool m_use_hydro_particles;
1423
1424 /** \brief Read tracer and hydro particles parameters */
1425 void readTracersParams();
1426
1427 /** \brief Initialize tracer and hydro particles */
1428 void initializeTracers ( amrex::ParGDBBase*,
1429 const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& );
1430
1431 /** \brief Evolve tracers and hydro particles */
1432 void evolveTracers( int lev,
1433 amrex::Real dt,
1434 amrex::Vector<amrex::MultiFab const*>& flow_vel,
1435 const amrex::Vector<std::unique_ptr<amrex::MultiFab>>& );
1436
1437#endif
1438 /** \brief Verbosity level of output */
1439 static int verbose;
1440
1441 /** \brief Diagnostic sum output interval in number of steps */
1442 static int sum_interval;
1443 /** \brief Diagnostic sum output interval in time */
1444 static amrex::Real sum_per;
1445
1446 /** \brief Minimum number of digits in plotfile name or chunked history file */
1448
1449 /** \brief Whether to write the staggered velocities (not averaged to cell centers) */
1451
1452 /** \brief Native or NetCDF plotfile output */
1454
1455 /** \brief NetCDF initialization file */
1456 static amrex::Vector<amrex::Vector<std::string>> nc_init_file;
1457 /** \brief NetCDF grid file */
1458 static amrex::Vector<amrex::Vector<std::string>> nc_grid_file;
1459 /** \brief NetCDF boundary data */
1460 static amrex::Vector<std::string> nc_bdry_file;
1461
1462 /** \brief Which level the high resolution bathymetry is at */
1464 /** \brief Grid file for high resolution bathymetry */
1466 /** \brief Box for the full domain at nc_hires_grid_level */
1468 /** \brief Cumulative refinement ratio between level 0 and level i */
1469 amrex::Vector<amrex::IntVect> cum_ref_ratios;
1470
1471 /** \brief NetCDF forcing file(s) */
1472 amrex::Vector<std::string> nc_frc_file;
1473 /** \brief NetCDF river file(s) */
1474 amrex::Vector<std::string> nc_riv_file;
1475
1476 /** \brief NetCDF climatology history file(s) */
1477 amrex::Vector<std::string> nc_clim_his_file;
1478 /** \brief NetCDF climatology coefficient file */
1480
1481 /** \brief Default name of time field for boundary data */
1482 std::string bdry_time_varname = "ocean_time";
1483 /** \brief Name of time fields for boundary data */
1484 amrex::Vector<std::string> bdry_time_name_byvar;
1485
1486 /** \brief Name of time field for ubar climatology data */
1487 std::string clim_ubar_time_varname = "ocean_time";
1488 /** \brief Name of time field for vbar climatology data */
1489 std::string clim_vbar_time_varname = "ocean_time";
1490 /** \brief Name of time field for u climatology data */
1491 std::string clim_u_time_varname = "ocean_time";
1492 /** \brief Name of time field for v climatology data */
1493 std::string clim_v_time_varname = "ocean_time";
1494 /** \brief Name of time field for salinity climatology data */
1495 std::string clim_salt_time_varname = "ocean_time";
1496 /** \brief Name of time field for temperature climatology data */
1497 std::string clim_temp_time_varname = "ocean_time";
1498 /** \brief Name of time field for river time */
1499 std::string riv_time_varname = "river_time";
1500 /** \brief Name of time field for forcing data */
1501 std::string frc_time_varname = "";
1502
1503 /** \brief Set refinement criteria */
1505
1506 /** \brief Holds info for dynamically generated tagging criteria */
1507 static amrex::Vector<amrex::AMRErrorTag> ref_tags;
1508
1509 /** \brief Mask that zeroes out values on a coarse level underlying grids on the next finest level */
1510 amrex::MultiFab fine_mask;
1511
1512 /** \brief Helper function to determine number of ghost cells */
1513 AMREX_FORCE_INLINE
1514 int
1515 ComputeGhostCells(const int& spatial_order) {
1516 int nGhostCells = 0;
1517 switch (spatial_order) {
1518 case 2:
1519 nGhostCells = 2; // We need this many to compute the eddy viscosity in the ghost cells
1520 break;
1521 case 3:
1522 nGhostCells = 2;
1523 break;
1524 case 4:
1525 nGhostCells = 2;
1526 break;
1527 case 5:
1528 nGhostCells = 3;
1529 break;
1530 case 6:
1531 nGhostCells = 3;
1532 break;
1533 default:
1534 amrex::Error("Must specify spatial order to be 2,3,4,5 or 6");
1535 }
1536
1537 return nGhostCells;
1538 }
1539
1540 /** \brief Helper function for IO stream */
1541 AMREX_FORCE_INLINE
1542 std::ostream&
1543 DataLog (int i)
1544 {
1545 return *datalog[i];
1546 }
1547
1548 AMREX_FORCE_INLINE
1549 int
1550 NumDataLogs () noexcept
1551 {
1552 return datalog.size();
1553 }
1554
1555 /** \brief Variable for CPU timing */
1556 static amrex::Real startCPUTime;
1557 /** \brief Accumulator variable for CPU time used thusfar */
1558 static amrex::Real previousCPUTimeUsed;
1559
1560 /** \brief Get CPU time used */
1561 amrex::Real
1563 {
1564 int numCores = amrex::ParallelDescriptor::NProcs();
1565#ifdef _OPENMP
1566 numCores = numCores * omp_get_max_threads();
1567#endif
1568
1569 amrex::Real T =
1570 numCores * (amrex::Real(amrex::ParallelDescriptor::second()) - startCPUTime) +
1572
1573 return T;
1574 }
1575
1576 void setRecordDataInfo (int i, const std::string& filename)
1577 {
1578 if (amrex::ParallelDescriptor::IOProcessor())
1579 {
1580 datalog[i] = std::make_unique<std::fstream>();
1581 datalog[i]->open(filename.c_str(),std::ios::out|std::ios::app);
1582 if (!datalog[i]->good()) {
1583 amrex::FileOpenFailed(filename);
1584 }
1585 }
1586 amrex::ParallelDescriptor::Barrier("REMORA::setRecordDataInfo");
1587 }
1588
1589 amrex::Vector<std::unique_ptr<std::fstream> > datalog;
1590 amrex::Vector<std::string> datalogname;
1591
1592 /** \brief The filename of the ith datalog file. */
1593 const std::string DataLogName (int i) const noexcept { return datalogname[i]; }
1594
1595public:
1596 /** \brief Write job info to stdout */
1597 void writeJobInfo (const std::string& dir) const;
1598 /** \brief Write build info to os */
1599 static void writeBuildInfo (std::ostream& os);
1600
1601 static void print_banner(MPI_Comm /*comm*/, std::ostream& /*out*/);
1602 static void print_usage(MPI_Comm /*comm*/, std::ostream& /*out*/);
1603 static void print_error(MPI_Comm /*comm*/, const std::string& msg);
1604 static void print_summary(std::ostream&);
1605 static void print_tpls(std::ostream& /*out*/);
1606
1607 /** \brief Accessor method for t_old to expose to outside classes */
1608 amrex::Real get_t_old(int lev) const;
1609
1610 /** \brief Evaluate stability function psi for wind speed */
1611 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
1612 static amrex::Real bulk_psiu(amrex::Real ZoL) {
1613 // Compute stability function, psi
1614 using namespace amrex;
1615
1616 Real psiu;
1617 if (ZoL < 0.0_rt) {
1618 // Unstable conditions.
1619 Real x = std::pow(1.0_rt-15.0_rt*ZoL,0.25_rt);
1620 Real psik = 2.0_rt*std::log(0.5_rt*(1.0_rt+x))+
1621 std::log(0.5_rt*(1.0_rt+x*x))-
1622 2.0_rt*std::atan(x)+0.5_rt*PI;
1623 // For very unstable conditions, use free-convection (Fairall).
1624 Real sqrt3 = std::sqrt(3.0_rt);
1625 Real y = std::pow(1.0_rt-10.15_rt*ZoL,1.0_rt/3.0_rt);
1626 Real psic = 1.5_rt*std::log(1.0_rt/3.0_rt*(1.0_rt+y+y*y))-
1627 sqrt3*std::atan((1.0_rt+2.0_rt*y)/sqrt3)+PI/sqrt3;
1628 // Match Kansas and free-convection forms with weighting Fw.
1629 Real Zol2 = ZoL*ZoL;
1630 Real Fw = Zol2/(1.0_rt+Zol2);
1631
1632 psiu = (1.0_rt-Fw)*psik+Fw*psic;
1633 } else {
1634 // Stable conditions
1635 Real cff=std::min(50.0_rt,0.35_rt*ZoL);
1636 psiu = -((1.0_rt+ZoL)+0.6667_rt*(ZoL-14.28_rt)/
1637 std::exp(cff)+8.525_rt);
1638 }
1639 return psiu;
1640 }
1641
1642 /** \brief Evaluate stability function psi for moisture and heat */
1643 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
1644 static amrex::Real bulk_psit(amrex::Real ZoL) {
1645 // Compute stability function, psi
1646 using namespace amrex;
1647
1648 Real psit;
1649 if (ZoL < 0.0_rt) {
1650 // Unstable conditions.
1651 Real x = std::pow(1.0_rt-15.0_rt*ZoL,0.5_rt);
1652 Real psik = 2.0_rt*std::log(0.5_rt*(1.0_rt+x));
1653 // For very unstable conditions, use free-convection (Fairall).
1654 Real sqrt3=std::sqrt(3.0_rt);
1655 Real y=std::pow(1.0_rt-34.15_rt*ZoL,1.0_rt/3.0_rt);
1656 Real psic=1.5_rt*std::log(1.0_rt/3.0_rt*(1.0_rt+y+y*y))-
1657 sqrt3*std::atan((1.0_rt+2.0_rt*y)/sqrt3)+PI/sqrt3;
1658
1659 // Match Kansas and free-convection forms with weighting Fw.
1660 Real ZoL2=ZoL*ZoL;
1661 Real Fw=ZoL2/(1.0_rt+ZoL2);
1662 psit = (1.0_rt-Fw)*psik+Fw*psic;
1663 } else {
1664 // Stable conditions.
1665 Real cff=std::min(50.0_rt,0.35_rt*ZoL);
1666 psit=-(std::pow(1.0_rt+2.0_rt*ZoL,1.5_rt)+
1667 0.6667_rt*(ZoL-14.28_rt)/std::exp(cff)+8.525_rt);
1668 }
1669 return psit;
1670 }
1671};
1672
1673#endif
constexpr amrex::Real PI
PlotfileType
plotfile format
#define Tracer_comp
A class to hold and interpolate time series data read from a NetCDF file.
A class to hold and interpolate time series data read from a NetCDF file.
Class that stores all relevant simulation state data with methods for time stepping.
Definition REMORA.H:90
static PlotfileType plotfile_type
Native or NetCDF plotfile output.
Definition REMORA.H:1453
std::string nc_grid_file_hires
Grid file for high resolution bathymetry.
Definition REMORA.H:1465
void FillCoarsePatchPC(int lev, amrex::Real time, amrex::MultiFab *mf_fine, amrex::MultiFab *mf_crse, const int bccomp, const int bdy_var_type=BdyVars::null, const int icomp=0, const bool fill_all=true, const int n_not_fill=0, const int icomp_calc=0, const amrex::Real dt=amrex::Real(0.0), const amrex::MultiFab &mf_calc=amrex::MultiFab())
fill an entire multifab by interpolating from the coarser level using the piecewise constant interpol...
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_EminusP
evaporation minus precipitation [kg/m^2/s], defined at rho-points
Definition REMORA.H:357
amrex::Vector< std::string > nc_riv_file
NetCDF river file(s)
Definition REMORA.H:1474
void allocate_bathymetry_full_domain()
Allocate multifabs for storing full-domain bathymetry data.
std::string riv_time_varname
Name of time field for river time.
Definition REMORA.H:1499
int foextrap_periodic_bc() const noexcept
Definition REMORA.H:1059
static void GotoNextLine(std::istream &is)
utility to skip to next line in Header
int nfast
Number of fast steps to take.
Definition REMORA.H:1360
amrex::Vector< std::string > nc_clim_his_file
NetCDF climatology history file(s)
Definition REMORA.H:1477
int ncons
Number of conserved scalars in the state (temperature + salt + passive scalars)
Definition REMORA.H:1346
void WriteNCMultiFab(const amrex::FabArray< amrex::FArrayBox > &fab, const std::string &name, bool set_ghost=false) const
Write MultiFab in NetCDF format.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv2d
v velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:272
void prsgrd(const amrex::Box &bx, const amrex::Box &gbx, const amrex::Box &utbx, const amrex::Box &vtbx, const amrex::Array4< amrex::Real > &ru, const amrex::Array4< amrex::Real > &rv, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &rho, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &z_w, const amrex::Array4< amrex::Real const > &msku, const amrex::Array4< amrex::Real const > &mskv, const int nrhs, const int N)
Calculate pressure gradient.
static amrex::Real fixed_dt
User specified fixed baroclinic time step.
Definition REMORA.H:1354
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_evap
evaporation rate [kg/m^2/s]
Definition REMORA.H:353
amrex::Real last_plot_file_time
Simulation time when we last output a plotfile.
Definition REMORA.H:1320
void update_massflux_3d(int lev, const amrex::Box &bx, const int ioff, const int joff, const amrex::Array4< amrex::Real > &phi, const amrex::Array4< amrex::Real > &phibar, const amrex::Array4< amrex::Real > &Hphi, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pm_or_pn, const amrex::Array4< amrex::Real const > &Dphi1, const amrex::Array4< amrex::Real const > &Dphi2, const amrex::Array4< amrex::Real > &DC, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real const > &msk, const int nnew)
Correct mass flux.
void scale_rhs_vars()
Scale RHS momentum variables by 1/cell area, needed before FillPatch to different levels.
int zvel_bc() const noexcept
Definition REMORA.H:1054
void init_bathymetry_from_netcdf(int lev)
Bathymetry data initialization from NetCDF file.
void init_bcs()
Read in boundary parameters from input file and set up data structures.
int xvel_bc() const noexcept
Definition REMORA.H:1052
amrex::Vector< amrex::BCRec > domain_bcs_type
vector (over BCVars) of BCRecs
Definition REMORA.H:1295
bool set_bcs_by_var
whether to set boundary conditions by variable rather than just by side
Definition REMORA.H:1280
void calculate_nodal_masks(int lev)
Calculate u-, v-, and psi-point masks based on rho-point masks after analytic initialization.
std::unique_ptr< NCTimeSeries > qair_data_from_file
Data container for specific humidity read from file.
Definition REMORA.H:1172
static amrex::Real previousCPUTimeUsed
Accumulator variable for CPU time used thusfar.
Definition REMORA.H:1558
amrex::Vector< std::string > cons_names
Names of scalars for plotfile output.
Definition REMORA.H:1397
amrex::Vector< std::unique_ptr< amrex::YAFluxRegister > > advflux_reg
array of flux registers for refluxing in multilevel
Definition REMORA.H:1286
std::unique_ptr< NCTimeSeries > sustr_data_from_file
Data container for u-component surface momentum flux read from file.
Definition REMORA.H:1162
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_fcor
coriolis factor (2D)
Definition REMORA.H:409
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_btflux
Bottom tracer flux; input arrays.
Definition REMORA.H:348
void init_gls_vmix(int lev, SolverChoice solver_choice)
Initialize GLS variables.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rubar
barotropic x velocity for the RHS (2D)
Definition REMORA.H:380
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_h
Bathymetry data (2D, positive valued, h in ROMS)
Definition REMORA.H:254
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pm
horizontal scaling factor: 1 / dx (2D)
Definition REMORA.H:404
amrex::MultiFab fine_mask
Mask that zeroes out values on a coarse level underlying grids on the next finest level.
Definition REMORA.H:1510
void set2DPlotVariables(const std::string &pp_plot_var_names_2d)
amrex::Vector< REMORAFillPatcher > FPr_v
Vector over levels of FillPatchers for v (3D)
Definition REMORA.H:1230
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ZoBot
Bottom roughness length [m], defined at rho points.
Definition REMORA.H:364
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DU_avg2
correct time average of barotropic x velocity flux for coupling (2D)
Definition REMORA.H:374
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_lrflx
longwave radiation
Definition REMORA.H:333
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yv
y_grid on v-points (2D)
Definition REMORA.H:424
static void print_error(MPI_Comm, const std::string &msg)
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive tracer
Definition REMORA.H:242
static bool write_history_file
Whether to output NetCDF files as a single history file with several time steps.
Definition REMORA.H:1159
virtual void MakeNewLevelFromCoarse(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Make a new level using provided BoxArray and DistributionMapping and fill with interpolated coarse le...
void FillCoarsePatch(int lev, amrex::Real time, amrex::MultiFab *mf_fine, amrex::MultiFab *mf_crse, const int bccomp, const int bdy_var_type=BdyVars::null, const int icomp=0, const bool fill_all=true, const int n_not_fill=0, const int icomp_calc=0, const amrex::Real dt=amrex::Real(0.0), const amrex::MultiFab &mf_calc=amrex::MultiFab())
fill an entire multifab by interpolating from the coarser level
int tke_bc() const noexcept
Definition REMORA.H:1058
void stretch_transform(int lev)
Calculate vertical stretched coordinates.
std::unique_ptr< NCTimeSeries > rain_data_from_file
Data container for precipitation rate read from file.
Definition REMORA.H:1180
amrex::Gpu::DeviceVector< amrex::Real > s_w
Scaled vertical coordinate (range [0,1]) that transforms to z, defined at w-points (cell faces)
Definition REMORA.H:296
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vwind
Wind in the v direction, defined at rho-points.
Definition REMORA.H:322
std::unique_ptr< ProblemBase > prob
Pointer to container of analytical functions for problem definition.
Definition REMORA.H:1259
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskr
land/sea mask at cell centers (2D)
Definition REMORA.H:393
void Construct_REMORAFillPatchers(int lev)
Construct FillPatchers.
Definition REMORA.cpp:374
void advance_3d(int lev, amrex::MultiFab &mf_cons, amrex::MultiFab &mf_u, amrex::MultiFab &mf_v, amrex::MultiFab *mf_sstore, amrex::MultiFab *mf_ru, amrex::MultiFab *mf_rv, std::unique_ptr< amrex::MultiFab > &mf_DU_avg1, std::unique_ptr< amrex::MultiFab > &mf_DU_avg2, std::unique_ptr< amrex::MultiFab > &mf_DV_avg1, std::unique_ptr< amrex::MultiFab > &mf_DV_avg2, std::unique_ptr< amrex::MultiFab > &mf_ubar, std::unique_ptr< amrex::MultiFab > &mf_vbar, std::unique_ptr< amrex::MultiFab > &mf_Akv, std::unique_ptr< amrex::MultiFab > &mf_Akt, std::unique_ptr< amrex::MultiFab > &mf_Hz, std::unique_ptr< amrex::MultiFab > &mf_Huon, std::unique_ptr< amrex::MultiFab > &mf_Hvom, std::unique_ptr< amrex::MultiFab > &mf_z_w, amrex::MultiFab const *mf_h, amrex::MultiFab const *mf_pm, amrex::MultiFab const *mf_pn, amrex::MultiFab const *mf_mskr, amrex::MultiFab const *mf_msku, amrex::MultiFab const *mf_mskv, const int N, const amrex::Real dt_lev)
Advance the 3D variables.
static int sum_interval
Diagnostic sum output interval in number of steps.
Definition REMORA.H:1442
int history_count
Counter for which time index we are writing to in the netcdf history file.
Definition REMORA.H:1390
amrex::Real stop_time
Time to stop.
Definition REMORA.H:1335
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rain
precipitation rate [kg/m^2/s]
Definition REMORA.H:351
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_tke
Turbulent kinetic energy.
Definition REMORA.H:452
void rho_eos(const amrex::Box &bx, const amrex::Array4< amrex::Real const > &state, const amrex::Array4< amrex::Real > &rho, const amrex::Array4< amrex::Real > &rhoA, const amrex::Array4< amrex::Real > &rhoS, const amrex::Array4< amrex::Real > &bvf, const amrex::Array4< amrex::Real > &alpha, const amrex::Array4< amrex::Real > &beta, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &z_w, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &h, const amrex::Array4< amrex::Real const > &mskr, const int N)
Wrapper around equation of state calculation.
static void print_summary(std::ostream &)
int do_substep
Whether to substep fine levels in time.
Definition REMORA.H:1363
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_stflx
Surface tracer flux; working arrays.
Definition REMORA.H:342
void Evolve()
Advance solution to final time.
Definition REMORA.cpp:163
void t3dmix2(const amrex::Box &bx, const amrex::Array4< amrex::Real > &state, const amrex::Array4< amrex::Real > &state_rhs, const amrex::Array4< amrex::Real const > &diff2, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &msku, const amrex::Array4< amrex::Real const > &mskv, const amrex::Real dt_lev, const int ncomp, const int N)
Wrapper for harmonic diffusivity for tracers.
std::string bdry_time_varname
Default name of time field for boundary data.
Definition REMORA.H:1482
amrex::Real plotfile_fill_value
fill value for masked arrays in amrex plotfiles
Definition REMORA.H:1409
void ReadCheckpointFile()
read checkpoint file from disk
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_gls
Turbulent generic length scale.
Definition REMORA.H:454
bool chunk_history_file
Whether to chunk netcdf history file.
Definition REMORA.H:1383
void writeJobInfo(const std::string &dir) const
Write job info to stdout.
int num_bc_vars() const noexcept
Definition REMORA.H:1063
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sustr
Surface stress in the u direction.
Definition REMORA.H:315
amrex::Real get_t_old(int lev) const
Accessor method for t_old to expose to outside classes.
Definition REMORA.cpp:1653
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yp
y_grid on psi-points (2D)
Definition REMORA.H:429
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xr
x_grid on rho points (2D)
Definition REMORA.H:412
int yvel_bc() const noexcept
Definition REMORA.H:1053
std::unique_ptr< NCTimeSeries > longwave_down_data_from_file
Data container for downward longwave radiation flux read from file.
Definition REMORA.H:1178
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru2d
u velocity RHS (2D, includes horizontal and vertical advection)
Definition REMORA.H:270
virtual void ClearLevel(int lev) override
Delete level data Overrides the pure virtual function in AmrCore.
amrex::Vector< std::string > datalogname
Definition REMORA.H:1590
amrex::Vector< amrex::MultiFab * > zvel_new
multilevel data container for current step's z velocities (largely unused; W stored separately)
Definition REMORA.H:248
AMREX_FORCE_INLINE int ComputeGhostCells(const int &spatial_order)
Helper function to determine number of ghost cells.
Definition REMORA.H:1515
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sstore
additional scratch space for calculations on temp, salt, etc
Definition REMORA.H:432
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xv
x_grid on v-points (2D)
Definition REMORA.H:422
void init_only(int lev, amrex::Real time)
Init (NOT restart or regrid)
Definition REMORA.cpp:1052
void init_set_vmix(int lev)
Initialize vertical mixing coefficients from file or analytic.
Definition REMORA.cpp:607
std::unique_ptr< NCTimeSeries > v_clim_data_from_file
Data container for v-velocity climatology data read from file.
Definition REMORA.H:1193
std::string clim_u_time_varname
Name of time field for u climatology data.
Definition REMORA.H:1491
void rhs_uv_3d(int lev, const amrex::Box &xbx, const amrex::Box &ybx, const amrex::Array4< amrex::Real const > &uold, const amrex::Array4< amrex::Real const > &vold, const amrex::Array4< amrex::Real > &ru, const amrex::Array4< amrex::Real > &rv, const amrex::Array4< amrex::Real > &rufrc, const amrex::Array4< amrex::Real > &rvfrc, const amrex::Array4< amrex::Real const > &sustr, const amrex::Array4< amrex::Real const > &svstr, const amrex::Array4< amrex::Real const > &bustr, const amrex::Array4< amrex::Real const > &bvstr, const amrex::Array4< amrex::Real const > &Huon, const amrex::Array4< amrex::Real const > &Hvom, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &W, const amrex::Array4< amrex::Real > &FC, int nrhs, int N)
RHS terms for 3D momentum.
void set_grid_scale(int lev)
Set pm and pn arrays on level lev. Only works if using analytic initialization.
void set_coriolis(int lev)
Initialize Coriolis factor from file or analytic.
Definition REMORA.cpp:578
int foextrap_bc() const noexcept
Definition REMORA.H:1060
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Lscale
Vertical mixing turbulent length scale.
Definition REMORA.H:456
amrex::Vector< REMORAFillPatcher > FPr_u
Vector over levels of FillPatchers for u (3D)
Definition REMORA.H:1228
std::string clim_temp_time_varname
Name of time field for temperature climatology data.
Definition REMORA.H:1497
amrex::Vector< amrex::Vector< amrex::Box > > boxes_at_level
the boxes specified at each level by tagging criteria
Definition REMORA.H:1266
void prestep_diffusion(const amrex::Box &bx, const amrex::Box &gbx, const int ioff, const int joff, const amrex::Array4< amrex::Real > &vel, const amrex::Array4< amrex::Real const > &vel_old, const amrex::Array4< amrex::Real > &rvel, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &Akv, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real const > &sstr, const amrex::Array4< amrex::Real const > &bstr, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const int iic, const int ntfirst, const int nnew, int nstp, int nrhs, int N, const amrex::Real lambda, const amrex::Real dt_lev)
Update velocities or tracers with diffusion/viscosity as the last part of the prestep.
void vert_mean_3d(const amrex::Box &bx, const int ioff, const int joff, const amrex::Array4< amrex::Real > &phi, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &Dphi_avg1, const amrex::Array4< amrex::Real > &DC, const amrex::Array4< amrex::Real > &CF, const amrex::Array4< amrex::Real const > &pm_or_pn, const amrex::Array4< amrex::Real const > &msk, const int nnew, const int N)
Adjust 3D momentum variables based on vertical mean momentum.
static amrex::Vector< amrex::AMRErrorTag > ref_tags
Holds info for dynamically generated tagging criteria.
Definition REMORA.H:1507
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Hz
Width of cells in the vertical (z-) direction (3D, Hz in ROMS)
Definition REMORA.H:260
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akt
Vertical diffusion coefficient (3D)
Definition REMORA.H:280
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_msku
land/sea mask at x-faces (2D)
Definition REMORA.H:395
std::unique_ptr< NCTimeSeriesRiver > river_source_transportbar
Data container for vertically integrated momentum transport in rivers.
Definition REMORA.H:1204
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rvfrc
v velocity RHS, integrated, including advection and bottom/surface stresses (2D)
Definition REMORA.H:276
int v2d_simple_bc() const noexcept
Definition REMORA.H:1062
std::string clim_ubar_time_varname
Name of time field for ubar climatology data.
Definition REMORA.H:1487
std::unique_ptr< NCTimeSeries > u_clim_data_from_file
Data container for u-velocity climatology data read from file.
Definition REMORA.H:1191
std::string check_file
Checkpoint file prefix.
Definition REMORA.H:1376
void FillPatchNoBC(int lev, amrex::Real time, amrex::MultiFab &mf_to_be_filled, amrex::Vector< amrex::MultiFab * > const &mfs, const int bdy_var_type=BdyVars::null, const int icomp=0, const bool fill_all=true, const bool fill_set=true)
Fill a new MultiFab by copying in phi from valid region and filling ghost cells without applying boun...
static amrex::Real startCPUTime
Variable for CPU timing.
Definition REMORA.H:1556
amrex::Vector< amrex::MultiFab * > xvel_old
multilevel data container for last step's x velocities (u in ROMS)
Definition REMORA.H:235
amrex::Real start_time
Time of the start of the simulation, in seconds.
Definition REMORA.H:1338
void init_data_from_netcdf(int lev)
Problem initialization from NetCDF file.
void init_masks_from_netcdf(int lev)
Mask data initialization from NetCDF file.
amrex::Vector< amrex::MultiFab * > yvel_new
multilevel data container for current step's y velocities (v in ROMS)
Definition REMORA.H:246
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_uwind
Wind in the u direction, defined at rho-points.
Definition REMORA.H:320
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rufrc
u velocity RHS, integrated, including advection and bottom/surface stresses (2D)
Definition REMORA.H:274
static amrex::Real fixed_fast_dt
User specified fixed barotropic time step.
Definition REMORA.H:1356
int regrid_int
how often each level regrids the higher levels of refinement (after a level advances that many time s...
Definition REMORA.H:1366
amrex::Real check_int_time
Checkpoint output interval in seconds.
Definition REMORA.H:1380
amrex::Real bdy_time_interval
Interval between boundary data times.
Definition REMORA.H:1156
void gls_prestep(int lev, amrex::MultiFab *mf_gls, amrex::MultiFab *mf_tke, amrex::MultiFab &mf_W, amrex::MultiFab *mf_msku, amrex::MultiFab *mf_mskv, const int nstp, const int nnew, const int iic, const int ntfirst, const int N, const amrex::Real dt_lev)
Prestep for GLS calculation.
void init_scalar_metadata()
Build runtime scalar names after nscalar is known.
Definition REMORA.cpp:150
amrex::Real start_bdy_time
Start time in the time series of boundary data.
Definition REMORA.H:1154
int zeta_bc() const noexcept
Definition REMORA.H:1057
amrex::Vector< amrex::GpuArray< REMORA_BC, AMREX_SPACEDIM *2 > > phys_bc_type
Array holding the "physical" boundary condition types (e.g. "inflow")
Definition REMORA.H:1306
amrex::Vector< amrex::Real > vec_weight2
Weights for calculating avg2 in 2D advance.
Definition REMORA.H:448
amrex::Vector< int > bdy_index
Container to connect boundary data being read in boundary condition containers.
Definition REMORA.H:1315
amrex::Gpu::DeviceVector< amrex::Real > s_r
Scaled vertical coordinate (range [0,1]) that transforms to z, defined at rho points (cell centers)
Definition REMORA.H:294
void Define_REMORAFillPatchers(int lev)
Define FillPatchers.
Definition REMORA.cpp:423
amrex::Vector< amrex::IntVect > cum_ref_ratios
Cumulative refinement ratio between level 0 and level i.
Definition REMORA.H:1469
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_shflx
sensible heat flux
Definition REMORA.H:339
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_visc2_p
Harmonic viscosity defined on the psi points (corners of horizontal grid cells)
Definition REMORA.H:282
void set_drag(int lev)
Initialize or calculate bottom drag.
amrex::MultiFab & build_fine_mask(int lev)
Make mask to zero out covered cells (for mesh refinement)
amrex::Real plot_int_time
Plotfile output interval in seconds.
Definition REMORA.H:1374
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rvbar
barotropic y velocity for the RHS (2D)
Definition REMORA.H:382
amrex::Vector< int > num_files_at_level
how many netcdf input files specified at each level
Definition REMORA.H:1264
void advance_2d_onestep(int lev, amrex::Real dt_lev, amrex::Real dtfast_lev, int my_iif, int nfast_counter)
2D advance, one predictor/corrector step
amrex::Vector< REMORAFillPatcher > FPr_vbar
Vector over levels of FillPatchers for vbar (2D)
Definition REMORA.H:1238
void AverageDownTo(int crse_lev)
more flexible version of AverageDown() that lets you average down across multiple levels
Definition REMORA.cpp:1608
int steps_per_history_file
Number of time steps per netcdf history file.
Definition REMORA.H:1388
void post_timestep(int nstep, amrex::Real time, amrex::Real dt_lev)
Called after every level 0 timestep.
Definition REMORA.cpp:237
int max_step
maximum number of steps
Definition REMORA.H:1333
amrex::Vector< amrex::MultiFab * > zvel_old
multilevel data container for last step's z velocities (largely unused; W stored separately)
Definition REMORA.H:239
std::unique_ptr< NCTimeSeries > svstr_data_from_file
Data container for v-component surface momentum flux read from file.
Definition REMORA.H:1164
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_z_r
z coordinates at rho points (cell centers)
Definition REMORA.H:289
amrex::Vector< std::string > nc_frc_file
NetCDF forcing file(s)
Definition REMORA.H:1472
amrex::Vector< int > num_boxes_at_level
how many boxes specified at each level by tagging criteria
Definition REMORA.H:1262
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
Definition REMORA.H:244
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_lhflx
latent heat flux
Definition REMORA.H:337
void refinement_criteria_setup()
Set refinement criteria.
int u2d_simple_bc() const noexcept
Definition REMORA.H:1061
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskp
land/sea mask at cell corners (2D)
Definition REMORA.H:399
void remora_advance(int level, amrex::MultiFab &cons_old, amrex::MultiFab &cons_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, const amrex::Geometry fine_geom, const amrex::Real dt, const amrex::Real time)
Interface for advancing the data at one level by one "slow" timestep.
int last_check_file_step
Step when we last output a checkpoint file.
Definition REMORA.H:1323
void prestep(int lev, amrex::MultiFab &mf_uold, amrex::MultiFab &mf_vold, amrex::MultiFab &mf_u, amrex::MultiFab &mf_v, amrex::MultiFab *mf_ru, amrex::MultiFab *mf_rv, amrex::MultiFab &S_old, amrex::MultiFab &S_new, amrex::MultiFab &mf_W, amrex::MultiFab &mf_DC, const amrex::MultiFab *mf_z_r, const amrex::MultiFab *mf_z_w, const amrex::MultiFab *mf_h, const amrex::MultiFab *mf_pm, const amrex::MultiFab *mf_pn, const amrex::MultiFab *mf_sustr, const amrex::MultiFab *mf_svstr, const amrex::MultiFab *mf_bustr, const amrex::MultiFab *mf_bvstr, const amrex::MultiFab *mf_msku, const amrex::MultiFab *mf_mskv, const int iic, const int nfirst, const int nnew, int nstp, int nrhs, int N, const amrex::Real dt_lev)
Wrapper function for prestep.
void init_beta_plane_coriolis(int lev)
Calculate Coriolis parameters from beta plane parametrization.
std::string clim_vbar_time_varname
Name of time field for vbar climatology data.
Definition REMORA.H:1489
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_bvf
Brunt-Vaisala frequency (3D)
Definition REMORA.H:439
virtual void RemakeLevel(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Remake an existing level using provided BoxArray and DistributionMapping and fill with existing fine ...
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskv
land/sea mask at y-faces (2D)
Definition REMORA.H:397
void init_masks(int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Allocate MultiFabs for masks.
amrex::Vector< int > nsubsteps
How many substeps on each level?
Definition REMORA.H:1271
amrex::Vector< std::unique_ptr< REMORAPhysBCFunct > > physbcs
Vector (over level) of functors to apply physical boundary conditions.
Definition REMORA.H:1283
void ComputeDt()
a wrapper for estTimeStep()
void fill_3d_masks(int lev)
Copy maskr to all z levels.
static void writeBuildInfo(std::ostream &os)
Write build info to os.
std::unique_ptr< NCTimeSeries > EminusP_data_from_file
Data container for evaporation minus precipitation read from file.
Definition REMORA.H:1184
virtual void ErrorEst(int lev, amrex::TagBoxArray &tags, amrex::Real time, int ngrow) override
Tag cells for refinement.
int plot_int
Plotfile output interval in iterations.
Definition REMORA.H:1372
void FillPatch(int lev, amrex::Real time, amrex::MultiFab &mf_to_be_filled, amrex::Vector< amrex::MultiFab * > const &mfs, const int bccomp, const int bdy_var_type=BdyVars::null, const int icomp=0, const bool fill_all=true, const bool fill_set=true, const int n_not_fill=0, const int icomp_calc=0, const amrex::Real dt=amrex::Real(0.0), const amrex::MultiFab &mf_calc=amrex::MultiFab())
Fill a new MultiFab by copying in phi from valid region and filling ghost cells.
std::unique_ptr< NCTimeSeries > cloud_data_from_file
Data container for cloud cover fraction read from file.
Definition REMORA.H:1182
amrex::Real getCPUTime() const
Get CPU time used.
Definition REMORA.H:1562
void InitData()
Initialize multilevel data.
Definition REMORA.cpp:269
void set3DPlotVariables(const std::string &pp_plot_var_names_3d)
amrex::Vector< int > istep
which step?
Definition REMORA.H:1269
void apply_clim_nudg(const amrex::Box &bx, int ioff, int joff, const amrex::Array4< amrex::Real > &var, const amrex::Array4< amrex::Real const > &var_old, const amrex::Array4< amrex::Real const > &var_clim, const amrex::Array4< amrex::Real const > &clim_coeff, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Real dt_lev=amrex::Real(0.0))
Apply climatology nudging.
void uv3dmix(const amrex::Box &xbx, const amrex::Box &ybx, const amrex::Array4< amrex::Real > &u, const amrex::Array4< amrex::Real > &v, const amrex::Array4< amrex::Real const > &uold, const amrex::Array4< amrex::Real const > &vold, const amrex::Array4< amrex::Real > &rufrc, const amrex::Array4< amrex::Real > &rvfrc, const amrex::Array4< amrex::Real const > &visc2_p, const amrex::Array4< amrex::Real const > &visc2_r, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &mskp, int nrhs, int nnew, const amrex::Real dt_lev)
Harmonic viscosity.
void WriteCheckpointFile()
write checkpoint file to disk
std::string nc_clim_coeff_file
NetCDF climatology coefficient file.
Definition REMORA.H:1479
void lin_eos(const amrex::Box &bx, const amrex::Array4< amrex::Real const > &state, const amrex::Array4< amrex::Real > &rho, const amrex::Array4< amrex::Real > &rhoA, const amrex::Array4< amrex::Real > &rhoS, const amrex::Array4< amrex::Real > &bvf, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &z_w, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &h, const amrex::Array4< amrex::Real const > &mskr, const int N)
Calculate density and related quantities from linear equation of state.
void setRecordDataInfo(int i, const std::string &filename)
Definition REMORA.H:1576
void advance_3d_ml(int lev, amrex::Real dt_lev)
3D advance on a single level
const std::string DataLogName(int i) const noexcept
The filename of the ith datalog file.
Definition REMORA.H:1593
void WriteMultiLevelPlotfileWithBathymetry(const std::string &plotfilename, int nlevels, const amrex::Vector< const amrex::MultiFab * > &mf, const amrex::Vector< const amrex::MultiFab * > &mf_nd, const amrex::Vector< const amrex::MultiFab * > &mf_u, const amrex::Vector< const amrex::MultiFab * > &mf_v, const amrex::Vector< const amrex::MultiFab * > &mf_w, const amrex::Vector< const amrex::MultiFab * > &mf_2d_rho, const amrex::Vector< const amrex::MultiFab * > &mf_2d_u, const amrex::Vector< const amrex::MultiFab * > &mf_2d_v, const amrex::Vector< std::string > &varnames_3d, const amrex::Vector< std::string > &varnames_2d_rho, const amrex::Vector< std::string > &varnames_2d_u, const amrex::Vector< std::string > &varnames_2d_v, const amrex::Vector< amrex::Geometry > &my_geom, amrex::Real time, const amrex::Vector< int > &level_steps, const amrex::Vector< amrex::IntVect > &rr, 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
write out particular data to an AMReX plotfile
void set_analytic_vmix(int lev)
Set vertical mixing coefficients from analytic.
Definition REMORA.cpp:624
void mask_arrays_for_write(int lev, amrex::Real fill_value, amrex::Real fill_where)
Mask data arrays before writing output.
void init_flat_bathymetry(int lev)
Initialize flat bathymetry to value from problo.
Definition REMORA.cpp:951
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rhoS
density perturbation
Definition REMORA.H:435
std::unique_ptr< NCTimeSeries > temp_clim_data_from_file
Data container for temperature climatology data read from file.
Definition REMORA.H:1195
amrex::Vector< std::string > bdry_time_name_byvar
Name of time fields for boundary data.
Definition REMORA.H:1484
void rhs_uv_2d(int lev, const amrex::Box &xbx, const amrex::Box &ybx, const amrex::Array4< amrex::Real const > &uold, const amrex::Array4< amrex::Real const > &vold, const amrex::Array4< amrex::Real > &ru, const amrex::Array4< amrex::Real > &rv, const amrex::Array4< amrex::Real const > &Duon, const amrex::Array4< amrex::Real const > &Dvom, const int nrhs)
RHS terms for 2D momentum.
static int file_min_digits
Minimum number of digits in plotfile name or chunked history file.
Definition REMORA.H:1447
void init_riv_pos_from_netcdf(int lev)
static amrex::Vector< std::string > nc_bdry_file
NetCDF boundary data.
Definition REMORA.H:51
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_visc2_r
Harmonic viscosity defined on the rho points (centers)
Definition REMORA.H:284
void update_mskp(int lev)
Set psi-point mask to be consistent with rho-point mask.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_svstr
Surface stress in the v direction.
Definition REMORA.H:317
std::unique_ptr< NCTimeSeries > srflx_data_from_file
Data container for shortwave radiation flux read from file.
Definition REMORA.H:1176
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Huon
u-volume flux (3D)
Definition REMORA.H:262
REMORA()
Definition REMORA.cpp:64
void FillBdyCCVels(int lev, amrex::MultiFab &mf_cc_vel)
Fill the physical boundary conditions for cell-centered velocity (diagnostic only)
void set_zeta(int lev)
Initialize zeta from file or analytic.
Definition REMORA.cpp:483
static amrex::Real change_max
Fraction maximum change in subsequent time steps.
Definition REMORA.H:1352
void calc_stretch_coeffs()
calculate vertical stretch coefficients
void init_zeta_from_netcdf(int lev)
Sea-surface height data initialization from NetCDF file.
void set_zeta_average(int lev)
Set Zt_avg1 to zeta.
void init_coriolis_from_netcdf(int lev)
Coriolis parameter data initialization from NetCDF file.
static void print_usage(MPI_Comm, std::ostream &)
void fill_from_bdyfiles(int lev, amrex::MultiFab &mf_to_fill, const amrex::MultiFab &mf_mask, const amrex::Real time, const int bccomp, const int bdy_var_type, const int icomp_to_fill, const int icomp_calc=0, const amrex::MultiFab &mf_calc=amrex::MultiFab(), const amrex::Real=amrex::Real(0.0))
Fill boundary data from netcdf file.
std::string pp_prefix
default prefix for input file parameters
Definition REMORA.H:230
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_h_full_domain
Bathymetry data on the whole domain at each potential level.
Definition REMORA.H:257
void set_bathymetry(int lev)
Initialize bathymetry from file or analytic.
Definition REMORA.cpp:511
amrex::Vector< amrex::MultiFab * > yvel_old
multilevel data container for last step's y velocities (v in ROMS)
Definition REMORA.H:237
std::unique_ptr< NCTimeSeries > ubar_clim_data_from_file
Data container for ubar climatology data read from file.
Definition REMORA.H:1187
void init_stuff(int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Allocate MultiFabs for state and evolution variables.
TimeInterpolatedData GetDataAtTime(int lev, amrex::Real time)
utility to copy in data from old and/or new state into another multifab
void Advance(int lev, amrex::Real time, amrex::Real dt_lev, int iteration, int ncycle)
advance a single level for a single time step
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rhoA
vertically-averaged density
Definition REMORA.H:437
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DV_avg1
time average of barotropic y velocity flux
Definition REMORA.H:376
amrex::Vector< REMORAFillPatcher > FPr_c
Vector over levels of FillPatchers for scalars.
Definition REMORA.H:1226
void WriteNCPlotFile(int istep, amrex::MultiFab const *plotMF)
Write plotfile using NetCDF (wrapper)
std::unique_ptr< NCTimeSeries > Tair_data_from_file
Data container for air temperature read from file.
Definition REMORA.H:1170
void rhs_t_3d(int lev, const amrex::Box &bx, const amrex::Array4< amrex::Real > &t, const amrex::Array4< amrex::Real const > &tempstore, const amrex::Array4< amrex::Real const > &Huon, const amrex::Array4< amrex::Real const > &Hvom, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &W, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real const > &mskr, const amrex::Array4< amrex::Real const > &msku, const amrex::Array4< amrex::Real const > &mskv, const amrex::Array4< int const > &river_pos, const amrex::Array4< amrex::Real const > &river_source, int nrhs, int nnew, int N, const amrex::Real dt_lev)
RHS terms for tracer.
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real bulk_psiu(amrex::Real ZoL)
Evaluate stability function psi for wind speed.
Definition REMORA.H:1612
static void print_banner(MPI_Comm, std::ostream &)
int nscalar
Number of passive scalars carried in the state.
Definition REMORA.H:1344
std::string clim_v_time_varname
Name of time field for v climatology data.
Definition REMORA.H:1493
void scale_rhs_vars_inv()
Scale RHS momentum variables by cell area, needed after FillPatch to different levels.
amrex::Vector< REMORAFillPatcher > FPr_w
Vector over levels of FillPatchers for w.
Definition REMORA.H:1232
void average_down_with_grow_cells(int lev, amrex::Vector< std::unique_ptr< amrex::MultiFab > > &mf)
Average down from level lev+1 to lev in mf, including grow cells.
Definition REMORA.cpp:1636
std::string clim_salt_time_varname
Name of time field for salinity climatology data.
Definition REMORA.H:1495
amrex::Gpu::DeviceVector< amrex::Real > Cs_r
Stretching coefficients at rho points.
Definition REMORA.H:304
std::unique_ptr< NCTimeSeries > Uwind_data_from_file
Data container for u-direction wind read from file.
Definition REMORA.H:1166
std::unique_ptr< NCTimeSeries > Pair_data_from_file
Data container for air pressure read from file.
Definition REMORA.H:1174
amrex::Vector< amrex::Real > t_new
new time at each level
Definition REMORA.H:1273
void init_stretch_coeffs()
initialize and calculate stretch coefficients
void init_bdry_from_netcdf(int lev)
Boundary data initialization from NetCDF file.
static SolverChoice solverChoice
Container for algorithmic choices.
Definition REMORA.H:1403
void resize_stuff(int lev)
Resize variable containers to accommodate data on levels 0 to max_lev.
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > vec_river_position
iMultiFab for river positions; contents are indices of rivers
Definition REMORA.H:1211
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akk
Turbulent kinetic energy vertical diffusion coefficient.
Definition REMORA.H:458
void set_masks(int lev)
Initialize land-sea masks from file or analytic.
Definition REMORA.cpp:640
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rdrag2
Quadratic drag coefficient [unitless], defined at rho points.
Definition REMORA.H:362
amrex::Vector< amrex::Vector< std::unique_ptr< NCTimeSeriesBoundary > > > boundary_series
Vector over BdyVars of boundary series data containers.
Definition REMORA.H:1207
int cf_set_width
Width for fixing values at coarse-fine interface.
Definition REMORA.H:1223
void ReadParameters()
read in some parameters from inputs file
Definition REMORA.cpp:1328
static void print_tpls(std::ostream &)
void WriteGenericPlotfileHeaderWithBathymetry(std::ostream &HeaderFile, int nlevels, const amrex::Vector< amrex::BoxArray > &bArray, const amrex::Vector< std::string > &varnames_3d, const amrex::Vector< std::string > &varnames_2d_rho, const amrex::Vector< std::string > &varnames_2d_u, const amrex::Vector< std::string > &varnames_2d_v, const amrex::Vector< amrex::Geometry > &my_geom, amrex::Real time, const amrex::Vector< int > &level_steps, const amrex::Vector< amrex::IntVect > &rr, const std::string &versionName, const std::string &levelPrefix, const std::string &mfPrefix) const
write out header data for an AMReX plotfile
void vert_visc_3d(const amrex::Box &bx, const int ioff, const int joff, const amrex::Array4< amrex::Real > &phi, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real > &Hzk, const amrex::Array4< amrex::Real > &AK, const amrex::Array4< amrex::Real const > &Akv, const amrex::Array4< amrex::Real > &BC, const amrex::Array4< amrex::Real > &DC, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real > &CF, const int nnew, const int N, const amrex::Real dt_lev)
Calculate effects of vertical viscosity or diffusivity.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru
u velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:266
void sum_integrated_quantities(amrex::Real time)
Integrate conserved quantities for diagnostics.
static int total_nc_plot_file_step
Definition REMORA.H:1067
amrex::Vector< amrex::Real > vec_weight1
Weights for calculating avg1 in 2D advance.
Definition REMORA.H:446
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xp
x_grid on psi-points (2D)
Definition REMORA.H:427
static int plot_staggered_vels
Whether to write the staggered velocities (not averaged to cell centers)
Definition REMORA.H:1450
static amrex::Vector< amrex::Vector< std::string > > nc_grid_file
NetCDF grid file.
Definition REMORA.H:53
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_longwave_down
Downward longwave radiation.
Definition REMORA.H:335
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_zeta
free surface height (2D)
Definition REMORA.H:390
int ubar_bc() const noexcept
Definition REMORA.H:1055
amrex::Gpu::DeviceVector< int > river_direction
Vector over rivers of river direction: 0: u-face; 1: v-face; 2: w-face.
Definition REMORA.H:1213
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vbar
barotropic y velocity (2D)
Definition REMORA.H:388
void ReadNCMultiFab(amrex::FabArray< amrex::FArrayBox > &fab, const std::string &name, int coordinatorProc=amrex::ParallelDescriptor::IOProcessorNumber(), int allow_empty_mf=0)
Read MultiFab in NetCDF format.
amrex::Box nc_hires_grid_box
Box for the full domain at nc_hires_grid_level.
Definition REMORA.H:1467
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DU_avg1
time average of barotropic x velocity flux (2D)
Definition REMORA.H:372
amrex::Gpu::DeviceVector< amrex::Real > Cs_w
Stretching coefficients at w points.
Definition REMORA.H:306
void set_zeta_to_Ztavg(int lev)
Set zeta components to be equal to time-averaged Zt_avg1.
bool expand_plotvars_to_unif_rr
whether plotfile variables should be expanded to a uniform refinement ratio
Definition REMORA.H:1406
void advance_2d(int lev, amrex::MultiFab const *mf_rhoS, amrex::MultiFab const *mf_rhoA, amrex::MultiFab *mf_ru2d, amrex::MultiFab *mf_rv2d, amrex::MultiFab *mf_rufrc, amrex::MultiFab *mf_rvfrc, amrex::MultiFab *mf_Zt_avg1, std::unique_ptr< amrex::MultiFab > &mf_DU_avg1, std::unique_ptr< amrex::MultiFab > &mf_DU_avg2, std::unique_ptr< amrex::MultiFab > &mf_DV_avg1, std::unique_ptr< amrex::MultiFab > &mf_DV_avg2, std::unique_ptr< amrex::MultiFab > &mf_rubar, std::unique_ptr< amrex::MultiFab > &mf_rvbar, std::unique_ptr< amrex::MultiFab > &mf_rzeta, std::unique_ptr< amrex::MultiFab > &mf_ubar, std::unique_ptr< amrex::MultiFab > &mf_vbar, amrex::MultiFab *mf_zeta, amrex::MultiFab const *mf_h, amrex::MultiFab const *mf_pm, amrex::MultiFab const *mf_pn, amrex::MultiFab const *mf_fcor, amrex::MultiFab const *mf_visc2_p, amrex::MultiFab const *mf_visc2_r, amrex::MultiFab const *mf_mskr, amrex::MultiFab const *mf_msku, amrex::MultiFab const *mf_mskv, amrex::MultiFab const *mf_mskp, amrex::Real dtfast_lev, bool predictor_2d_step, bool first_2d_step, int my_iif, int &next_indx1)
Perform a 2D predictor (predictor_2d_step=True) or corrector (predictor_2d_step=False) step.
int plot_file_on_restart
Whether to output a plotfile on restart from checkpoint.
Definition REMORA.H:1327
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_alpha
Thermal expansion coefficient (3D)
Definition REMORA.H:441
void set_2darrays(int lev)
Set 2D momentum arrays from 3D momentum.
void FillCoarsePatchMap(int lev, amrex::Real time, amrex::MultiFab *mf_fine, amrex::MultiFab *mf_crse, const int bccomp, const int bdy_var_type=BdyVars::null, const int icomp=0, const bool fill_all=true, const int n_not_fill=0, const int icomp_calc=0, const amrex::Real dt=amrex::Real(0.0), const amrex::MultiFab &mf_calc=amrex::MultiFab(), amrex::Interpolater *mapper=nullptr)
fill an entire multifab by interpolating from the coarser level, explicitly specifying interpolator t...
void init_analytic(int lev)
Initialize initial problem data from analytic functions.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ubar
barotropic x velocity (2D)
Definition REMORA.H:386
void InitializeLevelFromData(int lev, const amrex::MultiFab &initial_data)
Initialize the new-time data at a level from the initial_data MultiFab.
amrex::Vector< amrex::MultiFab * > cons_old
multilevel data container for last step's scalar data: temperature, salinity, passive tracer
Definition REMORA.H:233
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_bustr
Bottom stress in the u direction.
Definition REMORA.H:367
AMREX_FORCE_INLINE int NumDataLogs() noexcept
Definition REMORA.H:1550
amrex::Real volWgtSumMF(int lev, const amrex::MultiFab &mf, int comp, bool local, bool finemask)
Perform the volume-weighted sum.
std::string frc_time_varname
Name of time field for forcing data.
Definition REMORA.H:1501
void set_wind(int lev)
Initialize or calculate wind speed from file or analytic.
Definition REMORA.cpp:980
void init_clim_nudg_coeff_from_netcdf(int lev)
Climatology nudging coefficient initialization from NetCDF file.
amrex::Vector< REMORAFillPatcher > FPr_ubar
Vector over levels of FillPatchers for ubar (2D)
Definition REMORA.H:1236
bool is_it_time_for_action(int nstep, amrex::Real time, amrex::Real dt, int action_interval, amrex::Real action_per)
Decide if it is time to take an action.
int vbar_bc() const noexcept
Definition REMORA.H:1056
std::unique_ptr< NCTimeSeries > Vwind_data_from_file
Data container for v-direction wind read from file.
Definition REMORA.H:1168
std::string PlotFileName(int lev) const
get plotfile name
void init_bathymetry_full_domain_from_netcdf()
Full domain bathymetry data initialization from NetCDF file.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DV_avg2
correct time average of barotropic y velocity flux for coupling (2D)
Definition REMORA.H:378
void set_hmixcoef(int lev)
Initialize horizontal mixing coefficients.
Definition REMORA.cpp:666
amrex::Array< std::string, 2 *AMREX_SPACEDIM > domain_bc_type
Array of strings describing domain boundary conditions.
Definition REMORA.H:1300
amrex::Vector< std::unique_ptr< NCTimeSeriesRiver > > river_source_cons
Vector of data containers for scalar data in rivers.
Definition REMORA.H:1200
amrex::Real estTimeStep(int lev) const
compute dt from CFL considerations
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yu
y_grid on u-points (2D)
Definition REMORA.H:419
void WriteNCPlotFile_which(int lev, int which_subdomain, amrex::MultiFab const *plotMF, bool write_header, ncutils::NCFile &ncf, bool is_history)
Write a particular NetCDF plotfile.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_bvstr
Bottom stress in the v direction.
Definition REMORA.H:369
void timeStep(int lev, amrex::Real time, int iteration)
advance a level by dt, includes a recursive call for finer levels
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rzeta
free surface height for the RHS (2D)
Definition REMORA.H:384
std::unique_ptr< NCTimeSeriesRiver > river_source_transport
Data container for momentum transport in rivers.
Definition REMORA.H:1202
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_z_phys_nd
z coordinates at psi points (cell nodes)
Definition REMORA.H:309
void AverageDown()
set covered coarse cells to be the average of overlying fine cells
Definition REMORA.cpp:1595
static int fixed_ndtfast_ratio
User specified, number of barotropic steps per baroclinic step.
Definition REMORA.H:1358
amrex::Real netcdf_fill_value
fill value for masked arrays in netcdf output
Definition REMORA.H:1411
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xu
x_grid on u-points (2D)
Definition REMORA.H:417
void timeStepML(amrex::Real time, int iteration)
advance all levels by dt, loops over finer levels
amrex::GpuArray< amrex::GpuArray< bool, AMREX_SPACEDIM *2 >, BdyVars::NumTypes+1 > phys_bc_need_data
These are flags that indicate whether we need to read in boundary data from file.
Definition REMORA.H:1309
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > vec_nudg_coeff
Climatology nudging coefficients.
Definition REMORA.H:463
void set_weights(int lev)
Set weights for averaging 3D variables to 2D.
void prestep_t_advection(int lev, const amrex::Box &tbx, const amrex::Box &gbx, const amrex::Array4< amrex::Real > &tempold, const amrex::Array4< amrex::Real > &tempcache, const amrex::Array4< amrex::Real > &Hz, const amrex::Array4< amrex::Real > &Huon, const amrex::Array4< amrex::Real > &Hvom, const amrex::Array4< amrex::Real > &W, const amrex::Array4< amrex::Real > &DC, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real > &sstore, const amrex::Array4< amrex::Real const > &z_w, const amrex::Array4< amrex::Real const > &h, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &msku, const amrex::Array4< amrex::Real const > &mskv, const amrex::Array4< int const > &river_pos, const amrex::Array4< amrex::Real const > &river_source, int iic, int ntfirst, int nrhs, int N, const amrex::Real dt_lev)
Prestep advection calculations for the tracers.
AMREX_FORCE_INLINE std::ostream & DataLog(int i)
Helper function for IO stream.
Definition REMORA.H:1543
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pn
horizontal scaling factor: 1 / dy (2D)
Definition REMORA.H:406
amrex::Vector< std::unique_ptr< std::fstream > > datalog
Definition REMORA.H:1589
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akv
Vertical viscosity coefficient (3D)
Definition REMORA.H:278
static amrex::Real cfl
CFL condition.
Definition REMORA.H:1350
void append3DPlotVariables(const std::string &pp_plot_var_names_3d)
virtual void MakeNewLevelFromScratch(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Make a new level from scratch using provided BoxArray and DistributionMapping. Only used during initi...
amrex::Vector< std::string > plot_var_names_3d
Names of 3D variables to output to AMReX plotfile.
Definition REMORA.H:1393
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_stflux
Surface tracer flux; input arrays.
Definition REMORA.H:344
static int verbose
Verbosity level of output.
Definition REMORA.H:1439
void setup_step(int lev, amrex::Real time, amrex::Real dt_lev)
Set everything up for a step on a level.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rdrag
Linear drag coefficient [m/s], defined at rho points.
Definition REMORA.H:360
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_cloud
cloud cover fraction [0-1], defined at rho-points
Definition REMORA.H:355
void nonlin_eos(const amrex::Box &bx, const amrex::Array4< amrex::Real const > &state, const amrex::Array4< amrex::Real > &rho, const amrex::Array4< amrex::Real > &rhoA, const amrex::Array4< amrex::Real > &rhoS, const amrex::Array4< amrex::Real > &bvf, const amrex::Array4< amrex::Real > &alpha, const amrex::Array4< amrex::Real > &beta, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &z_w, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &h, const amrex::Array4< amrex::Real const > &mskr, const int N)
Calculate density and related quantities from nonlinear equation of state.
void t3dmix2_geo(const amrex::Box &bx, const amrex::Array4< amrex::Real > &state, const amrex::Array4< amrex::Real > &state_rhs, const amrex::Array4< amrex::Real const > &diff2, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &msku, const amrex::Array4< amrex::Real const > &mskv, const amrex::Real dt_lev, const int ncomp, const int N)
Harmonic diffusivity for tracers along geopotential surfaces.
std::string plot_file_name
Plotfile prefix.
Definition REMORA.H:1370
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Zt_avg1
Average of the free surface, zeta (2D)
Definition REMORA.H:312
std::unique_ptr< NCTimeSeries > salt_clim_data_from_file
Data container for salinity climatology data read from file.
Definition REMORA.H:1197
int check_int
Checkpoint output interval in iterations.
Definition REMORA.H:1378
void set_smflux(int lev)
Initialize or calculate surface momentum flux from file or analytic.
Definition REMORA.cpp:961
void WritePlotFile(int istep)
main driver for writing AMReX plotfiles
std::string restart_chkfile
If set, restart from this checkpoint file.
Definition REMORA.H:1341
void init_clim_nudg_coeff(int lev)
Wrapper to initialize climatology nudging coefficient.
void init_bathymetry_full_domain_from_analytic()
Full domain bathymetry data initialization from analytic.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv
v velocity RHS (3D, includes horizontal and vertical advection)
Definition REMORA.H:268
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_btflx
Bottom tracer flux; working arrays.
Definition REMORA.H:346
int cf_width
Nudging width at coarse-fine interface.
Definition REMORA.H:1221
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real bulk_psit(amrex::Real ZoL)
Evaluate stability function psi for moisture and heat.
Definition REMORA.H:1644
static amrex::Vector< amrex::Vector< std::string > > nc_init_file
NetCDF initialization file.
Definition REMORA.H:52
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_beta
Saline contraction coefficient (3D)
Definition REMORA.H:443
int last_plot_file_step
Step when we last output a plotfile.
Definition REMORA.H:1318
amrex::Vector< std::string > plot_var_names_2d
Names of 2D variables to output to AMReX plotfile.
Definition REMORA.H:1395
void bulk_fluxes(int lev, amrex::MultiFab *mf_cons, amrex::MultiFab *mf_uwind, amrex::MultiFab *mf_vwind, amrex::MultiFab *mf_Tair, amrex::MultiFab *mf_qair, amrex::MultiFab *mf_Pair, amrex::MultiFab *mf_srflx, amrex::MultiFab *mf_longwave_down, amrex::MultiFab *mf_evap, amrex::MultiFab *mf_sustr, amrex::MultiFab *mf_svstr, amrex::MultiFab *mf_stflux, amrex::MultiFab *mf_lrflx, amrex::MultiFab *mf_lhflx, amrex::MultiFab *mf_shflx, const int N)
Calculate bulk temperature, salinity, wind fluxes.
void t3dmix2_s(const amrex::Box &bx, const amrex::Array4< amrex::Real > &state, const amrex::Array4< amrex::Real > &state_rhs, const amrex::Array4< amrex::Real const > &diff2, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &msku, const amrex::Array4< amrex::Real const > &mskv, const amrex::Real dt_lev, const int ncomp)
Harmonic diffusivity for tracers along S-coordinate level surfaces.
const amrex::Vector< std::string > derived_names
Names of derived fields for plotfiles.
Definition REMORA.H:1400
amrex::Vector< amrex::Real > t_old
old time at each level
Definition REMORA.H:1275
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskr3d
land/sea mask at cell centers, copied to all z levels (3D)
Definition REMORA.H:401
amrex::Vector< amrex::GpuArray< amrex::Real, AMREX_SPACEDIM *2 > > m_bc_extdir_vals
Array holding the Dirichlet values at walls which need them.
Definition REMORA.H:1303
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_srflx
Shortwave radiation flux [W/m²], defined at rho-points.
Definition REMORA.H:331
amrex::Real last_check_file_time
Simulation time when we last output a checkpoint file.
Definition REMORA.H:1325
void append2DPlotVariables(const std::string &pp_plot_var_names_2d)
void set_bathymetry_averaged_down(int lev)
Copy over bathymetry data that has been averaged down from high resolution input netcdf file.
Definition REMORA.cpp:562
amrex::Vector< amrex::Real > dt
time step at each level
Definition REMORA.H:1277
static amrex::Real sum_per
Diagnostic sum output interval in time.
Definition REMORA.H:1444
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Pair
Air pressure [mb], defined at rho-points.
Definition REMORA.H:328
void InitializeFromFile()
Read the file passed to remora.restart and use it as an initial condition for the current simulation.
void coriolis(const amrex::Box &xbx, const amrex::Box &ybx, const amrex::Array4< amrex::Real const > &uold, const amrex::Array4< amrex::Real const > &vold, const amrex::Array4< amrex::Real > &ru, const amrex::Array4< amrex::Real > &rv, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &fomn, int nrhs, int nr)
Calculate Coriolis terms.
amrex::Gpu::DeviceVector< amrex::BCRec > domain_bcs_type_d
GPU vector (over BCVars) of BCRecs.
Definition REMORA.H:1297
virtual ~REMORA()
Definition REMORA.cpp:145
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_qair
Specific humidity [kg/kg], defined at rho-points.
Definition REMORA.H:326
int hires_grid_level
Which level the high resolution bathymetry is at.
Definition REMORA.H:1463
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_z_w
z coordinates at w points (faces between z-cells)
Definition REMORA.H:292
void restart()
Definition REMORA.cpp:470
std::unique_ptr< NCTimeSeries > vbar_clim_data_from_file
Data container for vbar climatology data read from file.
Definition REMORA.H:1189
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Tair
Air temperature [°C], defined at rho-points.
Definition REMORA.H:324
void convert_inv_days_to_inv_s(amrex::MultiFab *)
Convert data in a multifab from inverse days to inverse seconds.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yr
y_grid on rho points (2D)
Definition REMORA.H:414
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Hvom
v-volume flux (3D)
Definition REMORA.H:264
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akp
Turbulent length scale vertical diffusion coefficient.
Definition REMORA.H:460
void gls_corrector(int lev, amrex::MultiFab *mf_gls, amrex::MultiFab *mf_tke, amrex::MultiFab &mf_W, amrex::MultiFab *mf_Akv, amrex::MultiFab *mf_Akt, amrex::MultiFab *mf_Akk, amrex::MultiFab *mf_Akp, amrex::MultiFab *mf_mskr, amrex::MultiFab *mf_msku, amrex::MultiFab *mf_mskv, const int nstp, const int nnew, const int N, const amrex::Real dt_lev)
Corrector step for GLS calculation.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_diff2
Harmonic diffusivity for temperature / salinity.
Definition REMORA.H:286
int NumTypes(int ncons) noexcept
int xvel_bc(int ncons) noexcept
int zvel_bc(int ncons) noexcept
int yvel_bc(int ncons) noexcept
int vbar_bc(int ncons) noexcept
int foextrap_bc(int ncons) noexcept
int v2d_simple_bc(int ncons) noexcept
int u2d_simple_bc(int ncons) noexcept
int zeta_bc(int ncons) noexcept
int ubar_bc(int ncons) noexcept
int tke_bc(int ncons) noexcept
int foextrap_periodic_bc(int ncons) noexcept
@ CellConservativeQuartic
Definition REMORA.H:73
@ CellBilinear
Definition REMORA.H:70
@ PCInterp
Definition REMORA.H:67
@ NodeBilinear
Definition REMORA.H:68
@ CellConservativeLinear
Definition REMORA.H:69
@ CellQuadratic
Definition REMORA.H:71
@ CellConservativeProtected
Definition REMORA.H:72