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