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