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 condition type. Analytic reads from prob.cpp. Netcdf is from file */
30enum class IC_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
84/** \brief masks */
85enum class MaskType {
87};
88
90 public:
91 /** \brief read in and initialize parameters */
93 {
94 amrex::ParmParse pp(pp_prefix);
95
96 pp.queryAdd("flat_bathymetry", flat_bathymetry);
97
98 // Which horizontal advection scheme for tracers
99 static std::string tracer_hadv_string = "upstream3";
100 pp.queryAdd("tracer_horizontal_advection_scheme",tracer_hadv_string);
101 if (tracer_hadv_string == "centered4")
103 else if (tracer_hadv_string == "upstream3")
105 else
106 amrex::Error("Advection scheme unknown.");
107
108 // Which horizontal advection scheme
109 static std::string uv_hadv_string = "upstream3";
110 pp.queryAdd("uv_horizontal_advection_scheme",uv_hadv_string);
111 if (uv_hadv_string == "upstream3")
113 else if (uv_hadv_string == "centered2")
115 else
116 amrex::Error("UV advection scheme unknown.");
117
118 pp.queryAdd("rdrag", rdrag);
119 pp.queryAdd("rdrag2", rdrag2);
120 pp.queryAdd("Zos", Zos);
121 pp.queryAdd("Zob", Zob);
122 pp.queryAdd("Cdb_max", Cdb_max);
123 pp.queryAdd("Cdb_min", Cdb_min);
124
125 // Include salinity?
126 pp.queryAdd("use_salt", use_salt);
127
128 // Include Coriolis forcing?
129 pp.queryAdd("use_coriolis", use_coriolis);
130
131 // Include prestep / lagged predictor / corrections
132 pp.queryAdd("use_prestep", use_prestep);
133
134 //This affect forcing and some of the mixing terms for velocity
135 pp.queryAdd("use_uv3dmix", use_uv3dmix);
136
137 //This accounts for the main 2d loop but may not include coupling and copying properly
138 pp.queryAdd("use_barotropic", use_barotropic);
139
140 // Whether to do rivers. By default, rivers are temp and salt sources. Rivers
141 // have to be momentum sources.
142 pp.queryAdd("do_rivers", do_rivers);
143
144 pp.queryAdd("do_rivers_temp", do_rivers_temp);
145 pp.queryAdd("do_rivers_salt", do_rivers_salt);
146 pp.queryAdd("do_rivers_scalar", do_rivers_scalar);
147
148 do_rivers_cons.resize(NCONS);
149 // If we aren't doing rivers, set all of these flags to false
153
154 pp.queryAdd("init_l1ad_h", init_l1ad_h);
155 pp.queryAdd("init_l1ad_T", init_l1ad_T);
156
157 pp.queryAdd("init_ana_h", init_ana_h);
158 pp.queryAdd("init_ana_T", init_ana_T);
159
160 pp.queryAdd("init_l0int_h", init_l0int_h);
161 pp.queryAdd("init_l0int_T", init_l0int_T);
162
163 static std::string eos_type_string = "linear";
164 pp.queryAdd("eos_type",eos_type_string);
165 if (eos_type_string == "linear" || eos_type_string == "Linear" ||
166 eos_type_string == "lin" || eos_type_string == "Lin") {
168 pp.queryAdd("Tcoef",Tcoef);
169 pp.queryAdd("Scoef",Scoef);
170 } else if (eos_type_string == "nonlinear" || eos_type_string == "Nonlinear" ||
171 eos_type_string == "non-linear" || eos_type_string == "Non-linear" ||
172 eos_type_string == "nonlin" || eos_type_string == "Nonlin") {
174 } else {
175 amrex::Abort("Dont know this eos_type");
176 }
177 pp.queryAdd("R0",R0);
178 pp.queryAdd("S0",S0);
179 pp.queryAdd("T0",T0);
180 pp.queryAdd("rho0", rho0);
181
182 pp.queryAdd("bulk_fluxes",bulk_fluxes);
183 if (bulk_fluxes) {
184 do_salt_flux = true;
185 do_temp_flux = true;
186 }
187 pp.queryAdd("air_pressure",Pair);
188 pp.queryAdd("air_temperature",Tair);
189 pp.queryAdd("air_humidity",Hair);
190 pp.queryAdd("surface_radiation_flux",srflux);
191 pp.queryAdd("srflx_from_netcdf",srflx_from_netcdf);
192 pp.queryAdd("cloud",cloud);
193 pp.queryAdd("rain",rain);
194 pp.queryAdd("blk_ZQ",blk_ZQ);
195 pp.queryAdd("blk_ZT",blk_ZT);
196 pp.queryAdd("blk_ZW",blk_ZW);
197 pp.queryAdd("eminusp",eminusp);
198 pp.queryAdd("eminusp_correct_ssh",eminusp_correct_ssh);
199
200 // Control flags for loading atmospheric variables from NetCDF
201 pp.queryAdd("Tair_from_netcdf",Tair_from_netcdf);
202 pp.queryAdd("qair_from_netcdf",qair_from_netcdf);
203 pp.queryAdd("qair_is_percent",qair_is_percent);
204 pp.queryAdd("Pair_from_netcdf",Pair_from_netcdf);
205 pp.queryAdd("rain_from_netcdf",rain_from_netcdf);
206 pp.queryAdd("cloud_from_netcdf",cloud_from_netcdf);
207 pp.queryAdd("EminusP_from_netcdf",EminusP_from_netcdf);
208
209 if (eminusp_correct_ssh and !eminusp) {
210 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)");
211 }
212 if (eminusp and !bulk_fluxes) {
213 amrex::Abort("Evaporation minus precipitation (E-P) requires bulk flux parametrizations (remora.bulk_fluxes=true)");
214 }
215
216 //Read in linear eos parameters
217 //Grid stretching
218 pp.queryAdd("theta_s",theta_s);
219 pp.queryAdd("theta_b",theta_b);
220 pp.queryAdd("tcline",tcline);
221
222 //coriolis factor
223 pp.queryAdd("coriolis_f0",coriolis_f0);
224 pp.queryAdd("coriolis_beta",coriolis_beta);
225
226 pp.queryAdd("Akv_bak",Akv_bak);
227 pp.queryAdd("Akt_bak",Akt_bak);
228
229
230 static std::string grid_scale_type_string = "constant";
231 pp.queryAdd("grid_scale_type",grid_scale_type_string);
232
233 if (amrex::toLower(grid_scale_type_string) == "constant") {
235 } else if (amrex::toLower(grid_scale_type_string) == "custom") {
236 amrex::Warning("Initialization of grid scale from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
238 } else if (amrex::toLower(grid_scale_type_string) == "analytic") {
240 } else {
241 amrex::Error("Don't know this grid_scale_type");
242 }
243
244 static std::string ic_type_string = "analytic";
245 bool found_ic_bc = pp.queryAdd("ic_bc_type", ic_type_string);
246 pp.queryAdd("ic_type", ic_type_string);
247
248 if (found_ic_bc) {
249 amrex::Warning("remora.ic_bc_type is now called remora.ic_type, and will eventually be deprecated");
250 }
251
252 if ( amrex::toLower(ic_type_string) == "custom") {
253 amrex::Warning("Problem initialization from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
255 } else if ( amrex::toLower(ic_type_string) == "analytic") {
257 } else if ( amrex::toLower(ic_type_string) == "netcdf") {
259 } else if ( amrex::toLower(ic_type_string) == "real") {
260 amrex::Warning("Problem initialization from NetCDF (remora.ic_type) is now called 'netcdf'. 'real' will be deprecated");
262 } else {
263 amrex::Error("Don't know this ic_type");
264 }
265
266 // Which type of refinement
267 static std::string coupling_type_string = "TwoWay";
268 pp.queryAdd("coupling_type",coupling_type_string);
269 if (amrex::toLower(coupling_type_string) == "twoway" ||
270 amrex::toLower(coupling_type_string) == "two_way") {
272 } else if (amrex::toLower(coupling_type_string) == "oneway" ||
273 amrex::toLower(coupling_type_string) == "one_way") {
275 } else {
276 amrex::Abort("Dont know this coupling_type");
277 }
278
279 // Which type of coriolis forcing
280 if (use_coriolis) {
281 static std::string coriolis_type_string = "beta_plane";
282 pp.queryAdd("coriolis_type",coriolis_type_string);
283 if ( amrex::toLower(coriolis_type_string) == "custom") {
284 amrex::Warning("Coriolis initialization from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
286 } else if ( amrex::toLower(coriolis_type_string) == "analytic") {
288 } else if ((amrex::toLower(coriolis_type_string) == "beta_plane") ||
289 (amrex::toLower(coriolis_type_string) == "betaplane")) {
291 } else if ( (amrex::toLower(coriolis_type_string) == "netcdf")) {
293 } else if ( (amrex::toLower(coriolis_type_string) == "real")) {
294 amrex::Warning("Coriolis initialization from NetCDF is now called 'netcdf'. 'real' will be deprecated");
296 } else {
297 amrex::Abort("Don't know this coriolis_type");
298 }
299 }
300
301 static std::string smflux_type_string = "analytic";
302 int smflux_specified = pp.queryAdd("smflux_type",smflux_type_string);
303 if ( amrex::toLower(smflux_type_string) == "custom") {
304 amrex::Warning("Surface momentum flux initialization from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
306 } else if ( amrex::toLower(smflux_type_string) == "analytic") {
308 } else if ( amrex::toLower(smflux_type_string) == "netcdf") {
310 } else {
311 amrex::Abort("Don't know this smflux_type");
312 }
313
314 static std::string wind_type_string = "analytic";
315 int wind_specified = pp.queryAdd("wind_type",wind_type_string);
316 if ( amrex::toLower(wind_type_string) == "custom") {
317 amrex::Warning("Surface wind initialization from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
319 } else if ( amrex::toLower(wind_type_string) == "analytic") {
321 } else if ( amrex::toLower(wind_type_string) == "netcdf") {
323 } else {
324 amrex::Abort("Don't know this smflux_type");
325 }
326
327 if (wind_specified && smflux_specified) {
328 amrex::Abort("Cannot specify both wind and surface momentum flux");
329 }
330
331 static std::string mask_type_string = "none";
332 // If initial condition type is netcdf, default to netcdf masks
333 if (ic_type == IC_Type::netcdf) {
334 mask_type_string = "netcdf";
335 }
336 pp.queryAdd("mask_type", mask_type_string);
337 if (amrex::toLower(mask_type_string) == "none") {
339 } else if (amrex::toLower(mask_type_string) == "analytic") {
341 } else if (amrex::toLower(mask_type_string) == "netcdf") {
343 } else {
344 amrex::Abort("Don't know this mask_type");
345 }
346
347 static std::string bottom_stress_type_string = "linear";
348 pp.queryAdd("bottom_stress_type", bottom_stress_type_string);
349 if (amrex::toLower(bottom_stress_type_string) == "linear") {
351 } else if (amrex::toLower(bottom_stress_type_string) == "quadratic") {
353 } else if (amrex::toLower(bottom_stress_type_string) == "logarithmic") {
355 } else {
356 amrex::Abort("Don't know this bottom_stress_type");
357 }
358
359 amrex::Real tnu2_salt = amrex::Real(0.0);
360 amrex::Real tnu2_temp = amrex::Real(0.0);
361 amrex::Real tnu2_scalar = amrex::Real(0.0);
362 static std::string horiz_mixing_type_string = "analytic";
363 pp.queryAdd("horizontal_mixing_type", horiz_mixing_type_string);
364 if (amrex::toLower(horiz_mixing_type_string) == "analytical" ||
365 amrex::toLower(horiz_mixing_type_string) == "analytic") {
367 } else if (amrex::toLower(horiz_mixing_type_string) == "constant") {
369 } else {
370 amrex::Abort("Don't know this horizontal mixing type");
371 }
372 pp.queryAdd("visc2",visc2);
373 pp.queryAdd("tnu2_salt",tnu2_salt);
374 pp.queryAdd("tnu2_temp",tnu2_temp);
375 pp.queryAdd("tnu2_scalar",tnu2_scalar);
376 tnu2.resize(NCONS);
377 if (use_salt) {
378 tnu2[0] = tnu2_temp;
379 tnu2[1] = tnu2_salt;
380 tnu2[2] = tnu2_scalar;
381 } else {
382 tnu2[0] = tnu2_temp;
383 tnu2[1] = tnu2_scalar;
384 }
385
386 pp.queryAdd("Akk_bak", Akk_bak);
387 pp.queryAdd("Akp_bak", Akp_bak);
388 static std::string vert_mixing_type_string = "analytic";
389 static std::string gls_stability_type_string = "Canuto_A";
390 pp.queryAdd("vertical_mixing_type", vert_mixing_type_string);
391 pp.queryAdd("gls_stability_type", gls_stability_type_string);
392 if (amrex::toLower(vert_mixing_type_string) == "analytical" ||
393 amrex::toLower(vert_mixing_type_string) == "analytic") {
395 } else if (amrex::toLower(vert_mixing_type_string) == "gls") {
397 if (amrex::toLower(gls_stability_type_string) == "canuto_a") {
399 }
400 else if (amrex::toLower(gls_stability_type_string) == "canuto_b") {
402 }
403 else if (amrex::toLower(gls_stability_type_string) == "galperin") {
405 }
406 else {
407 amrex::Abort("Don't know this GLS stability type");
408 }
409 } else {
410 amrex::Abort("Don't know this vertical mixing type");
411 }
412 // Read in GLS params
414 pp.queryAdd("gls_P", gls_p);
415 pp.queryAdd("gls_M", gls_m);
416 pp.queryAdd("gls_N", gls_n);
417 pp.queryAdd("gls_Kmin", gls_Kmin);
418 pp.queryAdd("gls_Pmin", gls_Pmin);
419
420 pp.queryAdd("gls_cmu0", gls_cmu0);
421 pp.queryAdd("gls_c1", gls_c1);
422 pp.queryAdd("gls_c2", gls_c2);
423 pp.queryAdd("gls_c3m", gls_c3m);
424 pp.queryAdd("gls_c3p", gls_c3p);
425 pp.queryAdd("gls_sigk", gls_sigk);
426 pp.queryAdd("gls_sigp", gls_sigp);
428 gls_Gh0 = amrex::Real(0.0329); // 0.0329 GOTM, 0.0673 Burchard
429 gls_Ghcri = amrex::Real(0.03);
430 gls_L1 = amrex::Real(0.107);
431 gls_L2 = amrex::Real(0.0032);
432 gls_L3 = amrex::Real(0.0864);
433 gls_L4 = amrex::Real(0.12);
434 gls_L5 = amrex::Real(11.9);
435 gls_L6 = amrex::Real(0.4);
436 gls_L7 = amrex::Real(0.0);
437 gls_L8 = amrex::Real(0.48);
439 gls_Gh0 = amrex::Real(0.0444); // 0.044 GOTM, 0.0673 Burchard
440 gls_Ghcri = amrex::Real(0.0414);
441 gls_L1 = amrex::Real(0.127);
442 gls_L2 = amrex::Real(0.00336);
443 gls_L3 = amrex::Real(0.0906);
444 gls_L4 = amrex::Real(0.101);
445 gls_L5 = amrex::Real(11.2);
446 gls_L6 = amrex::Real(0.4);
447 gls_L7 = amrex::Real(0.0);
448 gls_L8 = amrex::Real(0.318);
449 } else {
450 gls_Gh0 = amrex::Real(0.028);
451 gls_Ghcri = amrex::Real(0.02);
452 }
453 }
454
455 // Read and compute inverse nudging coeffs from inputs given in days,
456 // and store in a vector corresponding to BdyVars enum
457 amrex::Real tnudg = amrex::Real(0.0);
458 amrex::Real znudg = amrex::Real(0.0);
459 amrex::Real m2nudg = amrex::Real(0.0);
460 amrex::Real m3nudg = amrex::Real(0.0);
461 pp.queryAdd("tnudg",tnudg);
462 pp.queryAdd("znudg",znudg);
463 pp.queryAdd("m2nudg",m2nudg);
464 pp.queryAdd("m3nudg",m3nudg);
465 pp.queryAdd("obcfac",obcfac);
466
468 nudg_coeff[BdyVars::u ] = (m3nudg > 0.0) ? amrex::Real(1.0) / (m3nudg * amrex::Real(86400.0)) : 0.0;//BdyVars::u
469 nudg_coeff[BdyVars::v ] = (m3nudg > 0.0) ? amrex::Real(1.0) / (m3nudg * amrex::Real(86400.0)) : 0.0;//BdyVars::v
470 nudg_coeff[BdyVars::t ] = ( tnudg > 0.0) ? amrex::Real(1.0) / ( tnudg * amrex::Real(86400.0)) : 0.0;//BdyVars::t
471 nudg_coeff[BdyVars::s ] = ( tnudg > 0.0) ? amrex::Real(1.0) / ( tnudg * amrex::Real(86400.0)) : 0.0;//BdyVars::s
472 nudg_coeff[BdyVars::ubar ] = (m2nudg > 0.0) ? amrex::Real(1.0) / (m2nudg * amrex::Real(86400.0)) : 0.0;//BdyVars::ubar
473 nudg_coeff[BdyVars::vbar ] = (m2nudg > 0.0) ? amrex::Real(1.0) / (m2nudg * amrex::Real(86400.0)) : 0.0;//BdyVars::vbar
474 nudg_coeff[BdyVars::zeta ] = ( znudg > 0.0) ? amrex::Real(1.0) / ( znudg * amrex::Real(86400.0)) : 0.0;//BdyVars::zeta
475
476 pp.queryAdd("do_m3_clim_nudg", do_m3_clim_nudg);
477 pp.queryAdd("do_m2_clim_nudg", do_m2_clim_nudg);
478 pp.queryAdd("do_temp_clim_nudg", do_temp_clim_nudg);
479 pp.queryAdd("do_salt_clim_nudg", do_salt_clim_nudg);
480
482 do_any_clim_nudg = true;
483 }
484#ifndef REMORA_USE_NETCDF
485 if (do_any_clim_nudg) {
486 amrex::Abort("Climatology nudging requires building with NetCDF");
487 }
488#endif
489 }
490
491 void display()
492 {
493 amrex::Print() << "SOLVER CHOICE: " << std::endl;
494 amrex::Print() << "use_salt : " << use_salt << std::endl;
495 amrex::Print() << "use_coriolis : " << use_coriolis << std::endl;
496 amrex::Print() << "use_prestep : " << use_prestep << std::endl;
497 amrex::Print() << "use_uv3dmix : " << use_uv3dmix << std::endl;
498 amrex::Print() << "use_barotropic : " << use_barotropic << std::endl;
499 amrex::Print() << "flat_bathymetry : " << flat_bathymetry << std::endl;
500 amrex::Print() << "spatial_order : " << spatial_order << std::endl;
501
502 if (ic_type == IC_Type::analytic) {
503 amrex::Print() << "Using analytic initial onditions" << std::endl;
504 }
505 else if (ic_type == IC_Type::netcdf) {
506 amrex::Print() << "Using NetCDF initial conditions" << std::endl;
507 }
508
510 amrex::Print() << "Horizontal advection scheme for tracers: " << "Centered 4" << std::endl;
511 }
513 amrex::Print() << "Horizontal advection scheme for tracers: " << "Upstream 3" << std::endl;
514 }
515 else {
516 amrex::Error("Invalid horizontal advection scheme for tracers.");
517 }
518
520 amrex::Print() << "Horizontal advection scheme for momenta: " << "Centered 2" << std::endl;
521 }
523 amrex::Print() << "Horizontal advection scheme for momenta: " << "Upstream 3" << std::endl;
524 }
525 else {
526 amrex::Error("Invalid horizontal advection scheme for momenta.");
527 }
528
530 amrex::Print() << "Using two-way coupling " << std::endl;
531 } else if (coupling_type == CouplingType::one_way) {
532 amrex::Print() << "Using one-way coupling " << std::endl;
533 }
534
535 if (use_coriolis) {
537 amrex::Print() << "Using analytic coriolis forcing " << std::endl;
538 } else if (coriolis_type == Cor_Type::beta_plane) {
539 amrex::Print() << "Using beta plane coriolis forcing " << std::endl;
540 } else if (coriolis_type == Cor_Type::netcdf) {
541 amrex::Print() << "Using coriolis forcing loaded from file " << std::endl;
542 }
543 }
544 }
545
546 // Default prefix
547 std::string pp_prefix {"remora"};
548
549 bool flat_bathymetry = false;
550
551 bool use_salt = true;
552
553 // Specify what additional physics/forcing modules we use
554 bool use_coriolis = false;
555
556 // Specify whether terms are used for debugging purposes
557 bool use_prestep = true;
558 bool use_uv3dmix = true;
559 bool use_baroclinic = true;
560 bool use_barotropic = true;
561
562 bool bulk_fluxes = false;
563 bool do_temp_flux = false;
564 bool do_salt_flux = false;
565
566 // Control flags for loading atmospheric variables from NetCDF
567 bool Tair_from_netcdf = false;
568 bool qair_from_netcdf = false;
569 bool Pair_from_netcdf = false;
570 bool srflx_from_netcdf = false;
571 bool rain_from_netcdf = false;
572 bool cloud_from_netcdf = false;
574 bool qair_is_percent = false;
575
576 bool do_rivers = false;
577 bool do_rivers_temp = true;
578 bool do_rivers_salt = true;
579 bool do_rivers_scalar = false;
580 amrex::Vector<int> do_rivers_cons;
581
582 bool init_l1ad_h = false;
583 bool init_l1ad_T = false;
584
585 bool init_ana_h = false;
586 bool init_ana_T = false;
587
588 bool init_l0int_h = true;
589 bool init_l0int_T = true;
590
593
594 // Coupling options are "OneWay" or "TwoWay"
596
597 // IC and BC Type: "analytic" or "netcdf"
599
600 // Coriolis forcing type
602
603 // Surface momentum flux type
605
606 // Surface wind speed type
608
609 // EOS type
611
612 // Bottom stress type
614
615 // Land/sea mask type
617
618 // Mixing type and parameters
622
623 // Type for grid scale (pm and pn)
625
626 // Stretching and depth parameters which may need to be read from inputs
627 amrex::Real theta_s = amrex::Real(3.0);
628 amrex::Real theta_b = amrex::Real(0.0);
629 amrex::Real tcline = amrex::Real(150.0);
630
631 // Linear drag coefficient [m/s]
632 amrex::Real rdrag = amrex::Real(3e-4);
633 // Quadratic drag coefficient [dimensionless]
634 amrex::Real rdrag2 = amrex::Real(3e-3);
635
636 // Momentum stress scales [m]
637 amrex::Real Zob = amrex::Real(2e-2);
638 amrex::Real Zos = amrex::Real(2e-2);
639
640 amrex::Real Cdb_max = amrex::Real(0.5);
641 amrex::Real Cdb_min = amrex::Real(1e-6);
642
643 // Linear equation of state parameters
644 amrex::Real R0 = amrex::Real(1028); // background density value (Kg/m3) used in Linear Equation of State
645 amrex::Real S0 = amrex::Real(35.0); // background salinity (nondimensional) constant
646 amrex::Real T0 = amrex::Real(5.0); // background potential temperature (Celsius) constant
647 amrex::Real Tcoef = amrex::Real(1.7e-4); // linear equation of state parameter (1/Celsius)
648 amrex::Real Scoef = amrex::Real(0.0); // linear equation of state parameter (nondimensional)
649 amrex::Real rho0 = amrex::Real(1025.0); // Mean density (Kg/m3) used when Boussinesq approx is inferred
650
651 // Coriolis forcing
652 amrex::Real coriolis_f0 = amrex::Real(0.0); // f-plane constant (1/s)
653 amrex::Real coriolis_beta = amrex::Real(0.0); // beta-plane constant (1/s/m)
654
655 // Air pressure
656 amrex::Real Pair = amrex::Real(1013.48);
657 // Air temperature
658 amrex::Real Tair = amrex::Real(23.567);
659 // Relative humidity (air)
660 amrex::Real Hair = amrex::Real(0.776);
661 // Cloud cover fraction (0=clear sky, 1=overcast)
662 amrex::Real cloud = amrex::Real(0.0);
663 // Precipitation rate (kg/m2/s)
664 amrex::Real rain = amrex::Real(0.0);
665 // Height (m) of atmospheric measurements for Bulk fluxes parametrization
666 amrex::Real blk_ZQ = amrex::Real(10.0); // air humidity
667 amrex::Real blk_ZT = amrex::Real(10.0); // air temperature
668 amrex::Real blk_ZW = amrex::Real(10.0); // winds
669
670 bool eminusp = false;
672
673 // Surface radiation flux
674 amrex::Real srflux = amrex::Real(0.0);
675
676 // Spatial discretization
678
679 // Horizontal mixing parameters
680 amrex::Real visc2 = amrex::Real(0.0);
681 amrex::Vector<amrex::Real> tnu2;
682
683 // GLS params
684 amrex::Real gls_p = amrex::Real(3.0);
685 amrex::Real gls_m = amrex::Real(1.5);
686 amrex::Real gls_n = amrex::Real(-1.0);
687 amrex::Real gls_Kmin = amrex::Real(7.6e-6);
688 amrex::Real gls_Pmin = amrex::Real(1.0e-12);
689
690 amrex::Real gls_cmu0 = amrex::Real(0.5477);
691 amrex::Real gls_c1 = amrex::Real(1.44);
692 amrex::Real gls_c2 = amrex::Real(1.92);
693 amrex::Real gls_c3m = amrex::Real(-0.4);
694 amrex::Real gls_c3p = amrex::Real(1.0);
695 amrex::Real gls_sigk = amrex::Real(1.0);
696 amrex::Real gls_sigp = amrex::Real(1.3);
697
698 // Turbulence closure
699 amrex::Real Akk_bak = amrex::Real(5.0e-6);
700 amrex::Real Akp_bak = amrex::Real(5.0e-6);
701 amrex::Real Akv_bak = amrex::Real(5.0e-6);
702 amrex::Real Akt_bak = amrex::Real(1.0e-6); // Note: this is a vector with one component per tracer in ROMS
703
704 // Params for stability functions.
705 amrex::Real gls_Gh0;
706 amrex::Real gls_Ghcri;
707 amrex::Real gls_Ghmin = amrex::Real(-0.28);
708 amrex::Real gls_E2 = amrex::Real(1.33);
709 // Params only for Canuto stability
710 amrex::Real gls_L1;
711 amrex::Real gls_L2;
712 amrex::Real gls_L3;
713 amrex::Real gls_L4;
714 amrex::Real gls_L5;
715 amrex::Real gls_L6;
716 amrex::Real gls_L7;
717 amrex::Real gls_L8;
718
719 // Params for some GLS and also Mellor-Yamada
720 amrex::Real my_A1 = amrex::Real(0.92);
721 amrex::Real my_A2 = amrex::Real(0.74);
722 amrex::Real my_B1 = amrex::Real(16.6);
723 amrex::Real my_B2 = amrex::Real(10.1);
724 amrex::Real my_C1 = amrex::Real(0.08);
725 amrex::Real my_C2 = amrex::Real(0.7);
726 amrex::Real my_C3 = amrex::Real(0.2);
727 amrex::Real my_E1 = amrex::Real(1.8);
728 amrex::Real my_E2 = amrex::Real(1.33);
729 amrex::Real my_Gh0 = amrex::Real(0.0233);
730 amrex::Real my_Sq = amrex::Real(0.2);
731 amrex::Real my_dtfac = amrex::Real(0.05);
732 amrex::Real my_lmax = amrex::Real(0.53);
733 amrex::Real my_qmin = amrex::Real(1.0E-8);
734
735 // Nudging time scales in 1/s
736 amrex::Vector<amrex::Real> nudg_coeff;
737
738 // Factor between passive (outflow) and active (inflow) open boundary
739 // conditions.
740 amrex::Real obcfac = amrex::Real(0.0);
741
742 // Whether to do climatoogy nudging
743 bool do_m2_clim_nudg = false;
744 bool do_m3_clim_nudg = false;
745 bool do_temp_clim_nudg = false;
746 bool do_salt_clim_nudg = false;
747 bool do_any_clim_nudg = false;
748
750};
751
752#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.
MaskType
masks
Cor_Type
Coriolis factor.
WindType
surface wind
GLS_StabilityType
stability function for GLS
CouplingType
Type of coupling between levels in AMR.
VertMixingType
vertical mixing type
IC_Type
Type of initial condition type. Analytic reads from prob.cpp. Netcdf is from file.
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
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