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