REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_DataStruct.H
Go to the documentation of this file.
1#ifndef _REMORA_DATA_STRUCT_H_
2#define _REMORA_DATA_STRUCT_H_
3
4#include <string>
5#include <iostream>
6
7#include <AMReX_ParmParse.H>
8#include <AMReX_Print.H>
9#include <AMReX_Gpu.H>
10
11#include <REMORA_Constants.H>
12#include "REMORA_IndexDefines.H"
13
14/** \brief Type of coupling between levels in AMR */
15enum struct CouplingType {
17};
18
19/** \brief Coordinates */
20enum class Coord {
21 x, y, z
22};
23
24/** \brief Horizontal advection schemes */
28
29/** \brief Type of initial/boundary condition type. Analytic reads from prob.cpp. Netcdf is from file */
30enum class IC_BC_Type {
32};
33
34/** \brief Coriolis factor */
35enum class Cor_Type {
37};
38
39/** \brief plotfile format */
40enum class PlotfileType {
42};
43
44/** \brief vertical mixing type */
45enum class VertMixingType {
47};
48
49/** \brief horizontal viscosity/diffusion type */
50enum class HorizMixingType {
52};
53
54/** \brief stability function for GLS */
58
59/** \brief equation of state */
60enum class EOSType {
62};
63
64/** \brief bottom stress formulation */
68
69/** \brief initialization for pm and pn */
70enum class GridScaleType {
72};
73
74/** \brief surface momentum flux */
75enum class SMFluxType {
77};
78
79/** \brief surface wind */
80enum class WindType {
82};
83
85 public:
86 /** \brief read in and initialize parameters */
88 {
89 amrex::ParmParse pp(pp_prefix);
90
91 pp.query("flat_bathymetry", flat_bathymetry);
92
93 // Which horizontal advection scheme for tracers
94 static std::string tracer_hadv_string = "upstream3";
95 pp.query("tracer_horizontal_advection_scheme",tracer_hadv_string);
96 if (tracer_hadv_string == "centered4")
98 else if (tracer_hadv_string == "upstream3")
100 else
101 amrex::Error("Advection scheme unknown.");
102
103 // Which horizontal advection scheme
104 static std::string uv_hadv_string = "upstream3";
105 pp.query("uv_horizontal_advection_scheme",uv_hadv_string);
106 if (uv_hadv_string == "upstream3")
108 else if (uv_hadv_string == "centered2")
110 else
111 amrex::Error("UV advection scheme unknown.");
112
113 pp.query("rdrag", rdrag);
114 pp.query("rdrag2", rdrag2);
115 pp.query("Zos", Zos);
116 pp.query("Zob", Zob);
117 pp.query("Cdb_max", Cdb_max);
118 pp.query("Cdb_min", Cdb_min);
119
120 // Include salinity?
121 pp.query("use_salt", use_salt);
122
123 // Include Coriolis forcing?
124 pp.query("use_coriolis", use_coriolis);
125
126 // Include prestep / lagged predictor / corrections
127 pp.query("use_prestep", use_prestep);
128
129 //This affect forcing and some of the mixing terms for velocity
130 pp.query("use_uv3dmix", use_uv3dmix);
131
132 //This accounts for the main 2d loop but may not include coupling and copying properly
133 pp.query("use_barotropic", use_barotropic);
134
135 // Whether to do rivers. By default, rivers are temp and salt sources. Rivers
136 // have to be momentum sources.
137 pp.query("do_rivers", do_rivers);
138
139 pp.query("do_rivers_temp", do_rivers_temp);
140 pp.query("do_rivers_salt", do_rivers_salt);
141 pp.query("do_rivers_scalar", do_rivers_scalar);
142
143 do_rivers_cons.resize(NCONS);
144 // If we aren't doing rivers, set all of these flags to false
148
149 pp.query("init_l1ad_h", init_l1ad_h);
150 pp.query("init_l1ad_T", init_l1ad_T);
151
152 pp.query("init_ana_h", init_ana_h);
153 pp.query("init_ana_T", init_ana_T);
154
155 pp.query("init_l0int_h", init_l0int_h);
156 pp.query("init_l0int_T", init_l0int_T);
157
158 static std::string eos_type_string = "linear";
159 pp.query("eos_type",eos_type_string);
160 if (eos_type_string == "linear" || eos_type_string == "Linear" ||
161 eos_type_string == "lin" || eos_type_string == "Lin") {
163 pp.query("Tcoef",Tcoef);
164 pp.query("Scoef",Scoef);
165 } else if (eos_type_string == "nonlinear" || eos_type_string == "Nonlinear" ||
166 eos_type_string == "non-linear" || eos_type_string == "Non-linear" ||
167 eos_type_string == "nonlin" || eos_type_string == "Nonlin") {
169 } else {
170 amrex::Abort("Dont know this eos_type");
171 }
172 pp.query("R0",R0);
173 pp.query("S0",S0);
174 pp.query("T0",T0);
175 pp.query("rho0", rho0);
176
177 pp.query("bulk_fluxes",bulk_fluxes);
178 if (bulk_fluxes) {
179 do_salt_flux = true;
180 do_temp_flux = true;
181 }
182 pp.query("air_pressure",Pair);
183 pp.query("air_temperature",Tair);
184 pp.query("air_humidity",Hair);
185 pp.query("surface_radiation_flux",srflux);
186 pp.query("cloud",cloud);
187 pp.query("rain",rain);
188 pp.query("blk_ZQ",blk_ZQ);
189 pp.query("blk_ZT",blk_ZT);
190 pp.query("blk_ZW",blk_ZW);
191 pp.query("eminusp",eminusp);
192 pp.query("eminusp_correct_ssh",eminusp_correct_ssh);
193
194 if (eminusp_correct_ssh and !eminusp) {
195 amrex::Abort("If evaporation minus precipitation (E-P) sea surface height correct in is on, E-P must be on as well (remora.eminusp=true)");
196 }
197 if (eminusp and !bulk_fluxes) {
198 amrex::Abort("Evaporation minus precipitation (E-P) requires bulk flux parametrizations (remora.bulk_fluxes=true)");
199 }
200
201 //Read in linear eos parameters
202 //Grid stretching
203 pp.query("theta_s",theta_s);
204 pp.query("theta_b",theta_b);
205 pp.query("tcline",tcline);
206
207 //coriolis factor
208 pp.query("coriolis_f0",coriolis_f0);
209 pp.query("coriolis_beta",coriolis_beta);
210
211 pp.query("Akv_bak",Akv_bak);
212 pp.query("Akt_bak",Akt_bak);
213
214
215 static std::string grid_scale_type_string = "constant";
216 pp.query("grid_scale_type",grid_scale_type_string);
217
218 if (amrex::toLower(grid_scale_type_string) == "constant") {
220 } else if (amrex::toLower(grid_scale_type_string) == "custom") {
221 amrex::Warning("Initialization of grid scale from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
223 } else if (amrex::toLower(grid_scale_type_string) == "analytic") {
225 } else {
226 amrex::Error("Don't know this grid_scale_type");
227 }
228
229 static std::string ic_bc_type_string = "analytic";
230 pp.query("ic_bc_type", ic_bc_type_string);
231
232 if ( amrex::toLower(ic_bc_type_string) == "custom") {
233 amrex::Warning("Problem initialization from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
235 } else if ( amrex::toLower(ic_bc_type_string) == "analytic") {
237 } else if ( amrex::toLower(ic_bc_type_string) == "netcdf") {
239 } else if ( amrex::toLower(ic_bc_type_string) == "real") {
240 amrex::Warning("Problem initialization from NetCDF (remora.ic_bc_type) is now called 'netcdf'. 'real' will be deprecated");
242 } else {
243 amrex::Error("Don't know this ic_bc_type");
244 }
245
246 // Which type of refinement
247 static std::string coupling_type_string = "TwoWay";
248 pp.query("coupling_type",coupling_type_string);
249 if (amrex::toLower(coupling_type_string) == "twoway" ||
250 amrex::toLower(coupling_type_string) == "two_way") {
252 } else if (amrex::toLower(coupling_type_string) == "oneway" ||
253 amrex::toLower(coupling_type_string) == "one_way") {
255 } else {
256 amrex::Abort("Dont know this coupling_type");
257 }
258
259 // Which type of coriolis forcing
260 if (use_coriolis) {
261 static std::string coriolis_type_string = "beta_plane";
262 pp.query("coriolis_type",coriolis_type_string);
263 if ( amrex::toLower(coriolis_type_string) == "custom") {
264 amrex::Warning("Coriolis initialization from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
266 } else if ( amrex::toLower(coriolis_type_string) == "analytic") {
268 } else if ((amrex::toLower(coriolis_type_string) == "beta_plane") ||
269 (amrex::toLower(coriolis_type_string) == "betaplane")) {
271 } else if ( (amrex::toLower(coriolis_type_string) == "netcdf")) {
273 } else if ( (amrex::toLower(coriolis_type_string) == "real")) {
274 amrex::Warning("Coriolis initialization from NetCDF is now called 'netcdf'. 'real' will be deprecated");
276 } else {
277 amrex::Abort("Don't know this coriolis_type");
278 }
279 }
280
281 static std::string smflux_type_string = "analytic";
282 int smflux_specified = pp.query("smflux_type",smflux_type_string);
283 if ( amrex::toLower(smflux_type_string) == "custom") {
284 amrex::Warning("Surface momentum flux initialization from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
286 } else if ( amrex::toLower(smflux_type_string) == "analytic") {
288 } else if ( amrex::toLower(smflux_type_string) == "netcdf") {
290 } else {
291 amrex::Abort("Don't know this smflux_type");
292 }
293
294 static std::string wind_type_string = "analytic";
295 int wind_specified = pp.query("wind_type",wind_type_string);
296 if ( amrex::toLower(wind_type_string) == "custom") {
297 amrex::Warning("Surface wind initialization from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
299 } else if ( amrex::toLower(wind_type_string) == "analytic") {
301 } else if ( amrex::toLower(wind_type_string) == "netcdf") {
303 } else {
304 amrex::Abort("Don't know this smflux_type");
305 }
306
307 if (wind_specified && smflux_specified) {
308 amrex::Abort("Cannot specify both wind and surface momentum flux");
309 }
310
311 static std::string bottom_stress_type_string = "linear";
312 pp.query("bottom_stress_type", bottom_stress_type_string);
313 if (amrex::toLower(bottom_stress_type_string) == "linear") {
315 } else if (amrex::toLower(bottom_stress_type_string) == "quadratic") {
317 } else if (amrex::toLower(bottom_stress_type_string) == "logarithmic") {
319 } else {
320 amrex::Abort("Don't know this bottom_stress_type");
321 }
322
323 amrex::Real tnu2_salt = amrex::Real(0.0);
324 amrex::Real tnu2_temp = amrex::Real(0.0);
325 amrex::Real tnu2_scalar = amrex::Real(0.0);
326 static std::string horiz_mixing_type_string = "analytic";
327 pp.query("horizontal_mixing_type", horiz_mixing_type_string);
328 if (amrex::toLower(horiz_mixing_type_string) == "analytical" ||
329 amrex::toLower(horiz_mixing_type_string) == "analytic") {
331 } else if (amrex::toLower(horiz_mixing_type_string) == "constant") {
333 } else {
334 amrex::Abort("Don't know this horizontal mixing type");
335 }
336 pp.query("visc2",visc2);
337 pp.query("tnu2_salt",tnu2_salt);
338 pp.query("tnu2_temp",tnu2_temp);
339 pp.query("tnu2_scalar",tnu2_scalar);
340 tnu2.resize(NCONS);
341 if (use_salt) {
342 tnu2[0] = tnu2_temp;
343 tnu2[1] = tnu2_salt;
344 tnu2[2] = tnu2_scalar;
345 } else {
346 tnu2[0] = tnu2_temp;
347 tnu2[1] = tnu2_scalar;
348 }
349
350 pp.query("Akk_bak", Akk_bak);
351 pp.query("Akp_bak", Akp_bak);
352 static std::string vert_mixing_type_string = "analytic";
353 static std::string gls_stability_type_string = "Canuto_A";
354 pp.query("vertical_mixing_type", vert_mixing_type_string);
355 pp.query("gls_stability_type", gls_stability_type_string);
356 if (amrex::toLower(vert_mixing_type_string) == "analytical" ||
357 amrex::toLower(vert_mixing_type_string) == "analytic") {
359 } else if (amrex::toLower(vert_mixing_type_string) == "gls") {
361 if (amrex::toLower(gls_stability_type_string) == "canuto_a") {
363 }
364 else if (amrex::toLower(gls_stability_type_string) == "canuto_b") {
366 }
367 else if (amrex::toLower(gls_stability_type_string) == "galperin") {
369 }
370 else {
371 amrex::Abort("Don't know this GLS stability type");
372 }
373 } else {
374 amrex::Abort("Don't know this vertical mixing type");
375 }
376 // Read in GLS params
378 pp.query("gls_P", gls_p);
379 pp.query("gls_M", gls_m);
380 pp.query("gls_N", gls_n);
381 pp.query("gls_Kmin", gls_Kmin);
382 pp.query("gls_Pmin", gls_Pmin);
383
384 pp.query("gls_cmu0", gls_cmu0);
385 pp.query("gls_c1", gls_c1);
386 pp.query("gls_c2", gls_c2);
387 pp.query("gls_c3m", gls_c3m);
388 pp.query("gls_c3p", gls_c3p);
389 pp.query("gls_sigk", gls_sigk);
390 pp.query("gls_sigp", gls_sigp);
392 gls_Gh0 = amrex::Real(0.0329); // 0.0329 GOTM, 0.0673 Burchard
393 gls_Ghcri = amrex::Real(0.03);
394 gls_L1 = amrex::Real(0.107);
395 gls_L2 = amrex::Real(0.0032);
396 gls_L3 = amrex::Real(0.0864);
397 gls_L4 = amrex::Real(0.12);
398 gls_L5 = amrex::Real(11.9);
399 gls_L6 = amrex::Real(0.4);
400 gls_L7 = amrex::Real(0.0);
401 gls_L8 = amrex::Real(0.48);
403 gls_Gh0 = amrex::Real(0.0444); // 0.044 GOTM, 0.0673 Burchard
404 gls_Ghcri = amrex::Real(0.0414);
405 gls_L1 = amrex::Real(0.127);
406 gls_L2 = amrex::Real(0.00336);
407 gls_L3 = amrex::Real(0.0906);
408 gls_L4 = amrex::Real(0.101);
409 gls_L5 = amrex::Real(11.2);
410 gls_L6 = amrex::Real(0.4);
411 gls_L7 = amrex::Real(0.0);
412 gls_L8 = amrex::Real(0.318);
413 } else {
414 gls_Gh0 = amrex::Real(0.028);
415 gls_Ghcri = amrex::Real(0.02);
416 }
417 }
418
419 // Read and compute inverse nudging coeffs from inputs given in days,
420 // and store in a vector corresponding to BdyVars enum
421 amrex::Real tnudg = amrex::Real(0.0);
422 amrex::Real znudg = amrex::Real(0.0);
423 amrex::Real m2nudg = amrex::Real(0.0);
424 amrex::Real m3nudg = amrex::Real(0.0);
425 pp.query("tnudg",tnudg);
426 pp.query("znudg",znudg);
427 pp.query("m2nudg",m2nudg);
428 pp.query("m3nudg",m3nudg);
429 pp.query("obcfac",obcfac);
430
432 nudg_coeff[BdyVars::u ] = (m3nudg > 0.0) ? amrex::Real(1.0) / (m3nudg * amrex::Real(86400.0)) : 0.0;//BdyVars::u
433 nudg_coeff[BdyVars::v ] = (m3nudg > 0.0) ? amrex::Real(1.0) / (m3nudg * amrex::Real(86400.0)) : 0.0;//BdyVars::v
434 nudg_coeff[BdyVars::t ] = ( tnudg > 0.0) ? amrex::Real(1.0) / ( tnudg * amrex::Real(86400.0)) : 0.0;//BdyVars::t
435 nudg_coeff[BdyVars::s ] = ( tnudg > 0.0) ? amrex::Real(1.0) / ( tnudg * amrex::Real(86400.0)) : 0.0;//BdyVars::s
436 nudg_coeff[BdyVars::ubar ] = (m2nudg > 0.0) ? amrex::Real(1.0) / (m2nudg * amrex::Real(86400.0)) : 0.0;//BdyVars::ubar
437 nudg_coeff[BdyVars::vbar ] = (m2nudg > 0.0) ? amrex::Real(1.0) / (m2nudg * amrex::Real(86400.0)) : 0.0;//BdyVars::vbar
438 nudg_coeff[BdyVars::zeta ] = ( znudg > 0.0) ? amrex::Real(1.0) / ( znudg * amrex::Real(86400.0)) : 0.0;//BdyVars::zeta
439
440 pp.query("do_m3_clim_nudg", do_m3_clim_nudg);
441 pp.query("do_m2_clim_nudg", do_m2_clim_nudg);
442 pp.query("do_temp_clim_nudg", do_temp_clim_nudg);
443 pp.query("do_salt_clim_nudg", do_salt_clim_nudg);
444
446 do_any_clim_nudg = true;
447 }
448#ifndef REMORA_USE_NETCDF
449 if (do_any_clim_nudg) {
450 amrex::Abort("Climatology nudging requires building with NetCDF");
451 }
452#endif
453 }
454
455 void display()
456 {
457 amrex::Print() << "SOLVER CHOICE: " << std::endl;
458 amrex::Print() << "use_salt : " << use_salt << std::endl;
459 amrex::Print() << "use_coriolis : " << use_coriolis << std::endl;
460 amrex::Print() << "use_prestep : " << use_prestep << std::endl;
461 amrex::Print() << "use_uv3dmix : " << use_uv3dmix << std::endl;
462 amrex::Print() << "use_barotropic : " << use_barotropic << std::endl;
463 amrex::Print() << "flat_bathymetry : " << flat_bathymetry << std::endl;
464 amrex::Print() << "spatial_order : " << spatial_order << std::endl;
465
467 amrex::Print() << "Using analytic initial and boundary conditions" << std::endl;
468 }
469 else if (ic_bc_type == IC_BC_Type::netcdf) {
470 amrex::Print() << "Using NetCDF initial and boundary conditions" << std::endl;
471 }
472
474 amrex::Print() << "Horizontal advection scheme for tracers: " << "Centered 4" << std::endl;
475 }
477 amrex::Print() << "Horizontal advection scheme for tracers: " << "Upstream 3" << std::endl;
478 }
479 else {
480 amrex::Error("Invalid horizontal advection scheme for tracers.");
481 }
482
484 amrex::Print() << "Horizontal advection scheme for momenta: " << "Centered 2" << std::endl;
485 }
487 amrex::Print() << "Horizontal advection scheme for momenta: " << "Upstream 3" << std::endl;
488 }
489 else {
490 amrex::Error("Invalid horizontal advection scheme for momenta.");
491 }
492
494 amrex::Print() << "Using two-way coupling " << std::endl;
495 } else if (coupling_type == CouplingType::one_way) {
496 amrex::Print() << "Using one-way coupling " << std::endl;
497 }
498
499 if (use_coriolis) {
501 amrex::Print() << "Using analytic coriolis forcing " << std::endl;
502 } else if (coriolis_type == Cor_Type::beta_plane) {
503 amrex::Print() << "Using beta plane coriolis forcing " << std::endl;
504 } else if (coriolis_type == Cor_Type::netcdf) {
505 amrex::Print() << "Using coriolis forcing loaded from file " << std::endl;
506 }
507 }
508 }
509
510 // Default prefix
511 std::string pp_prefix {"remora"};
512
513 bool flat_bathymetry = false;
514
515 bool use_salt = true;
516
517 // Specify what additional physics/forcing modules we use
518 bool use_coriolis = false;
519
520 // Specify whether terms are used for debugging purposes
521 bool use_prestep = true;
522 bool use_uv3dmix = true;
523 bool use_baroclinic = true;
524 bool use_barotropic = true;
525
526 bool bulk_fluxes = false;
527 bool do_temp_flux = false;
528 bool do_salt_flux = false;
529
530 bool do_rivers = false;
531 bool do_rivers_temp = true;
532 bool do_rivers_salt = true;
533 bool do_rivers_scalar = false;
534 amrex::Vector<int> do_rivers_cons;
535
536 bool init_l1ad_h = false;
537 bool init_l1ad_T = false;
538
539 bool init_ana_h = false;
540 bool init_ana_T = false;
541
542 bool init_l0int_h = true;
543 bool init_l0int_T = true;
544
547
548 // Coupling options are "OneWay" or "TwoWay"
550
551 // IC and BC Type: "analytic" or "netcdf"
553
554 // Coriolis forcing type
556
557 // Surface momentum flux type
559
560 // Surface wind speed type
562
563 // EOS type
565
566 // Bottom stress type
568
569 // Mixing type and parameters
573
574 // Type for grid scale (pm and pn)
576
577 // Stretching and depth parameters which may need to be read from inputs
578 amrex::Real theta_s = amrex::Real(3.0);
579 amrex::Real theta_b = amrex::Real(0.0);
580 amrex::Real tcline = amrex::Real(150.0);
581
582 // Linear drag coefficient [m/s]
583 amrex::Real rdrag = amrex::Real(3e-4);
584 // Quadratic drag coefficient [dimensionless]
585 amrex::Real rdrag2 = amrex::Real(3e-3);
586
587 // Momentum stress scales [m]
588 amrex::Real Zob = amrex::Real(2e-2);
589 amrex::Real Zos = amrex::Real(2e-2);
590
591 amrex::Real Cdb_max = amrex::Real(0.5);
592 amrex::Real Cdb_min = amrex::Real(1e-6);
593
594 // Linear equation of state parameters
595 amrex::Real R0 = amrex::Real(1028); // background density value (Kg/m3) used in Linear Equation of State
596 amrex::Real S0 = amrex::Real(35.0); // background salinity (nondimensional) constant
597 amrex::Real T0 = amrex::Real(5.0); // background potential temperature (Celsius) constant
598 amrex::Real Tcoef = amrex::Real(1.7e-4); // linear equation of state parameter (1/Celsius)
599 amrex::Real Scoef = amrex::Real(0.0); // linear equation of state parameter (nondimensional)
600 amrex::Real rho0 = amrex::Real(1025.0); // Mean density (Kg/m3) used when Boussinesq approx is inferred
601
602 // Coriolis forcing
603 amrex::Real coriolis_f0 = amrex::Real(0.0); // f-plane constant (1/s)
604 amrex::Real coriolis_beta = amrex::Real(0.0); // beta-plane constant (1/s/m)
605
606 // Air pressure
607 amrex::Real Pair = amrex::Real(1013.48);
608 // Air temperature
609 amrex::Real Tair = amrex::Real(23.567);
610 // Relative humidity (air)
611 amrex::Real Hair = amrex::Real(0.776);
612 // Cloud cover fraction (0=clear sky, 1=overcast)
613 amrex::Real cloud = amrex::Real(0.0);
614 // Precipitation rate (kg/m2/s)
615 amrex::Real rain = amrex::Real(0.0);
616 // Height (m) of atmospheric measurements for Bulk fluxes parametrization
617 amrex::Real blk_ZQ = amrex::Real(10.0); // air humidity
618 amrex::Real blk_ZT = amrex::Real(10.0); // air temperature
619 amrex::Real blk_ZW = amrex::Real(10.0); // winds
620
621 bool eminusp = false;
623
624 // Surface radiation flux
625 amrex::Real srflux = amrex::Real(0.0);
626
627 // Spatial discretization
629
630 // Horizontal mixing parameters
631 amrex::Real visc2 = amrex::Real(0.0);
632 amrex::Vector<amrex::Real> tnu2;
633
634 // GLS params
635 amrex::Real gls_p = amrex::Real(3.0);
636 amrex::Real gls_m = amrex::Real(1.5);
637 amrex::Real gls_n = amrex::Real(-1.0);
638 amrex::Real gls_Kmin = amrex::Real(7.6e-6);
639 amrex::Real gls_Pmin = amrex::Real(1.0e-12);
640
641 amrex::Real gls_cmu0 = amrex::Real(0.5477);
642 amrex::Real gls_c1 = amrex::Real(1.44);
643 amrex::Real gls_c2 = amrex::Real(1.92);
644 amrex::Real gls_c3m = amrex::Real(-0.4);
645 amrex::Real gls_c3p = amrex::Real(1.0);
646 amrex::Real gls_sigk = amrex::Real(1.0);
647 amrex::Real gls_sigp = amrex::Real(1.3);
648
649 // Turbulence closure
650 amrex::Real Akk_bak = amrex::Real(5.0e-6);
651 amrex::Real Akp_bak = amrex::Real(5.0e-6);
652 amrex::Real Akv_bak = amrex::Real(5.0e-6);
653 amrex::Real Akt_bak = amrex::Real(1.0e-6); // Note: this is a vector with one component per tracer in ROMS
654
655 // Params for stability functions.
656 amrex::Real gls_Gh0;
657 amrex::Real gls_Ghcri;
658 amrex::Real gls_Ghmin = amrex::Real(-0.28);
659 amrex::Real gls_E2 = amrex::Real(1.33);
660 // Params only for Canuto stability
661 amrex::Real gls_L1;
662 amrex::Real gls_L2;
663 amrex::Real gls_L3;
664 amrex::Real gls_L4;
665 amrex::Real gls_L5;
666 amrex::Real gls_L6;
667 amrex::Real gls_L7;
668 amrex::Real gls_L8;
669
670 // Params for some GLS and also Mellor-Yamada
671 amrex::Real my_A1 = amrex::Real(0.92);
672 amrex::Real my_A2 = amrex::Real(0.74);
673 amrex::Real my_B1 = amrex::Real(16.6);
674 amrex::Real my_B2 = amrex::Real(10.1);
675 amrex::Real my_C1 = amrex::Real(0.08);
676 amrex::Real my_C2 = amrex::Real(0.7);
677 amrex::Real my_C3 = amrex::Real(0.2);
678 amrex::Real my_E1 = amrex::Real(1.8);
679 amrex::Real my_E2 = amrex::Real(1.33);
680 amrex::Real my_Gh0 = amrex::Real(0.0233);
681 amrex::Real my_Sq = amrex::Real(0.2);
682 amrex::Real my_dtfac = amrex::Real(0.05);
683 amrex::Real my_lmax = amrex::Real(0.53);
684 amrex::Real my_qmin = amrex::Real(1.0E-8);
685
686 // Nudging time scales in 1/s
687 amrex::Vector<amrex::Real> nudg_coeff;
688
689 // Factor between passive (outflow) and active (inflow) open boundary
690 // conditions.
691 amrex::Real obcfac = amrex::Real(0.0);
692
693 // Whether to do climatoogy nudging
694 bool do_m2_clim_nudg = false;
695 bool do_m3_clim_nudg = false;
696 bool do_temp_clim_nudg = false;
697 bool do_salt_clim_nudg = false;
698 bool do_any_clim_nudg = false;
699};
700
701#endif
BottomStressType
bottom stress formulation
GridScaleType
initialization for pm and pn
SMFluxType
surface momentum flux
AdvectionScheme
Horizontal advection schemes.
PlotfileType
plotfile format
HorizMixingType
horizontal viscosity/diffusion type
Coord
Coordinates.
Cor_Type
Coriolis factor.
IC_BC_Type
Type of initial/boundary condition type. Analytic reads from prob.cpp. Netcdf is from file.
WindType
surface wind
GLS_StabilityType
stability function for GLS
CouplingType
Type of coupling between levels in AMR.
VertMixingType
vertical mixing type
EOSType
equation of state
#define Temp_comp
#define Scalar_comp
#define Salt_comp
#define NCONS
GLS_StabilityType gls_stability_type
HorizMixingType horiz_mixing_type
amrex::Real Cdb_min
amrex::Real Akv_bak
amrex::Real blk_ZT
amrex::Real cloud
amrex::Real coriolis_beta
amrex::Real gls_sigp
amrex::Vector< amrex::Real > nudg_coeff
amrex::Real my_E2
amrex::Real rdrag2
amrex::Real gls_m
amrex::Real my_lmax
amrex::Real my_Gh0
amrex::Real my_B2
amrex::Real my_A1
IC_BC_Type ic_bc_type
amrex::Vector< amrex::Real > tnu2
amrex::Real gls_sigk
amrex::Real gls_L3
amrex::Real coriolis_f0
AdvectionScheme uv_Hadv_scheme
amrex::Vector< int > do_rivers_cons
amrex::Real Tcoef
amrex::Real gls_cmu0
void init_params()
read in and initialize parameters
std::string pp_prefix
amrex::Real gls_L6
amrex::Real gls_L2
AdvectionScheme tracer_Hadv_scheme
amrex::Real Akk_bak
amrex::Real gls_L1
amrex::Real blk_ZW
amrex::Real my_E1
amrex::Real my_C3
amrex::Real gls_n
amrex::Real theta_b
amrex::Real theta_s
amrex::Real Akt_bak
amrex::Real visc2
amrex::Real tcline
amrex::Real gls_c3m
BottomStressType bottom_stress_type
amrex::Real my_A2
amrex::Real gls_Gh0
amrex::Real gls_Ghmin
amrex::Real my_C1
amrex::Real gls_c1
amrex::Real gls_L5
amrex::Real gls_L8
amrex::Real my_dtfac
amrex::Real srflux
SMFluxType smflux_type
amrex::Real my_C2
amrex::Real gls_Kmin
amrex::Real rdrag
VertMixingType vert_mixing_type
amrex::Real gls_c3p
amrex::Real gls_p
amrex::Real gls_L4
amrex::Real my_B1
amrex::Real blk_ZQ
amrex::Real gls_L7
amrex::Real Cdb_max
GridScaleType grid_scale_type
amrex::Real my_qmin
amrex::Real gls_Ghcri
amrex::Real obcfac
amrex::Real gls_Pmin
amrex::Real gls_c2
amrex::Real Scoef
amrex::Real Akp_bak
amrex::Real my_Sq
CouplingType coupling_type
amrex::Real gls_E2