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