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 // Outputs forcing variables if true
188 pp.queryAdd("output_forcing", output_forcing);
189 pp.queryAdd("air_pressure",Pair);
190 pp.queryAdd("air_temperature",Tair);
191 pp.queryAdd("air_humidity",Hair);
192 pp.queryAdd("surface_radiation_flux",srflux);
193 pp.queryAdd("srflx_from_netcdf",srflx_from_netcdf);
194 pp.queryAdd("longwave_down", longwave_down);
195 pp.queryAdd("longwave_down_from_netcdf",longwave_down_from_netcdf);
196 pp.queryAdd("longwave_netcdf_is_net", longwave_netcdf_is_net);
197 pp.queryAdd("longwave_netcdf_varname", longwave_netcdf_varname);
198 pp.queryAdd("cloud",cloud);
199 pp.queryAdd("rain",rain);
200 pp.queryAdd("blk_ZQ",blk_ZQ);
201 pp.queryAdd("blk_ZT",blk_ZT);
202 pp.queryAdd("blk_ZW",blk_ZW);
203 pp.queryAdd("eminusp",eminusp);
204 pp.queryAdd("eminusp_correct_ssh",eminusp_correct_ssh);
205 // Control flags for loading atmospheric variables from NetCDF
206 pp.queryAdd("Tair_from_netcdf",Tair_from_netcdf);
207 pp.queryAdd("qair_from_netcdf",qair_from_netcdf);
208 pp.queryAdd("qair_is_percent",qair_is_percent);
209 pp.queryAdd("Pair_from_netcdf",Pair_from_netcdf);
210 pp.queryAdd("rain_from_netcdf",rain_from_netcdf);
211 pp.queryAdd("cloud_from_netcdf",cloud_from_netcdf);
212 pp.queryAdd("EminusP_from_netcdf",EminusP_from_netcdf);
213
215 // If user provides net longwave from NetCDF, force use of file-based longwave.
216 longwave_down = true;
218 }
219
221 amrex::Warning("remora.longwave_netcdf_varname is set but remora.longwave_down_from_netcdf is false; the value will be ignored.");
222 }
223
224 if (eminusp_correct_ssh and !eminusp) {
225 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)");
226 }
227 if (eminusp and !bulk_fluxes) {
228 amrex::Abort("Evaporation minus precipitation (E-P) requires bulk flux parametrizations (remora.bulk_fluxes=true)");
229 }
230
231 //Read in linear eos parameters
232 //Grid stretching
233 pp.queryAdd("theta_s",theta_s);
234 pp.queryAdd("theta_b",theta_b);
235 pp.queryAdd("tcline",tcline);
236
237 //coriolis factor
238 pp.queryAdd("coriolis_f0",coriolis_f0);
239 pp.queryAdd("coriolis_beta",coriolis_beta);
240
241 pp.queryAdd("Akv_bak",Akv_bak);
242 pp.queryAdd("Akt_bak",Akt_bak);
243
244
245 static std::string grid_scale_type_string = "constant";
246 pp.queryAdd("grid_scale_type",grid_scale_type_string);
247
248 if (amrex::toLower(grid_scale_type_string) == "constant") {
250 } else if (amrex::toLower(grid_scale_type_string) == "custom") {
251 amrex::Warning("Initialization of grid scale from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
253 } else if (amrex::toLower(grid_scale_type_string) == "analytic") {
255 } else {
256 amrex::Error("Don't know this grid_scale_type");
257 }
258
259 static std::string ic_type_string = "analytic";
260 bool found_ic_bc = pp.queryAdd("ic_bc_type", ic_type_string);
261 pp.queryAdd("ic_type", ic_type_string);
262
263 if (found_ic_bc) {
264 amrex::Warning("remora.ic_bc_type is now called remora.ic_type, and will eventually be deprecated");
265 }
266
267 if ( amrex::toLower(ic_type_string) == "custom") {
268 amrex::Warning("Problem initialization from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
270 } else if ( amrex::toLower(ic_type_string) == "analytic") {
272 } else if ( amrex::toLower(ic_type_string) == "netcdf") {
274 } else if ( amrex::toLower(ic_type_string) == "real") {
275 amrex::Warning("Problem initialization from NetCDF (remora.ic_type) is now called 'netcdf'. 'real' will be deprecated");
277 } else {
278 amrex::Error("Don't know this ic_type");
279 }
280
281 // Which type of refinement
282 static std::string coupling_type_string = "TwoWay";
283 pp.queryAdd("coupling_type",coupling_type_string);
284 if (amrex::toLower(coupling_type_string) == "twoway" ||
285 amrex::toLower(coupling_type_string) == "two_way") {
287 } else if (amrex::toLower(coupling_type_string) == "oneway" ||
288 amrex::toLower(coupling_type_string) == "one_way") {
290 } else {
291 amrex::Abort("Dont know this coupling_type");
292 }
293
294 // Which type of coriolis forcing
295 if (use_coriolis) {
296 static std::string coriolis_type_string = "beta_plane";
297 pp.queryAdd("coriolis_type",coriolis_type_string);
298 if ( amrex::toLower(coriolis_type_string) == "custom") {
299 amrex::Warning("Coriolis initialization from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
301 } else if ( amrex::toLower(coriolis_type_string) == "analytic") {
303 } else if ((amrex::toLower(coriolis_type_string) == "beta_plane") ||
304 (amrex::toLower(coriolis_type_string) == "betaplane")) {
306 } else if ( (amrex::toLower(coriolis_type_string) == "netcdf")) {
308 } else if ( (amrex::toLower(coriolis_type_string) == "real")) {
309 amrex::Warning("Coriolis initialization from NetCDF is now called 'netcdf'. 'real' will be deprecated");
311 } else {
312 amrex::Abort("Don't know this coriolis_type");
313 }
314 }
315
316 static std::string smflux_type_string = "analytic";
317 int smflux_specified = pp.queryAdd("smflux_type",smflux_type_string);
318 if ( amrex::toLower(smflux_type_string) == "custom") {
319 amrex::Warning("Surface momentum flux initialization from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
321 } else if ( amrex::toLower(smflux_type_string) == "analytic") {
323 } else if ( amrex::toLower(smflux_type_string) == "netcdf") {
325 } else {
326 amrex::Abort("Don't know this smflux_type");
327 }
328
329 static std::string wind_type_string = "analytic";
330 int wind_specified = pp.queryAdd("wind_type",wind_type_string);
331 if ( amrex::toLower(wind_type_string) == "custom") {
332 amrex::Warning("Surface wind initialization from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
334 } else if ( amrex::toLower(wind_type_string) == "analytic") {
336 } else if ( amrex::toLower(wind_type_string) == "netcdf") {
338 } else {
339 amrex::Abort("Don't know this smflux_type");
340 }
341
342 if (wind_specified && smflux_specified) {
343 amrex::Abort("Cannot specify both wind and surface momentum flux");
344 }
345
346 static std::string mask_type_string = "none";
347 // If initial condition type is netcdf, default to netcdf masks
348 if (ic_type == IC_Type::netcdf) {
349 mask_type_string = "netcdf";
350 }
351 pp.queryAdd("mask_type", mask_type_string);
352 if (amrex::toLower(mask_type_string) == "none") {
354 } else if (amrex::toLower(mask_type_string) == "analytic") {
356 } else if (amrex::toLower(mask_type_string) == "netcdf") {
358 } else {
359 amrex::Abort("Don't know this mask_type");
360 }
361
362 static std::string bottom_stress_type_string = "linear";
363 pp.queryAdd("bottom_stress_type", bottom_stress_type_string);
364 if (amrex::toLower(bottom_stress_type_string) == "linear") {
366 } else if (amrex::toLower(bottom_stress_type_string) == "quadratic") {
368 } else if (amrex::toLower(bottom_stress_type_string) == "logarithmic") {
370 } else {
371 amrex::Abort("Don't know this bottom_stress_type");
372 }
373
374 amrex::Real tnu2_salt = amrex::Real(0.0);
375 amrex::Real tnu2_temp = amrex::Real(0.0);
376 amrex::Real tnu2_scalar = amrex::Real(0.0);
377 static std::string horiz_mixing_type_string = "analytic";
378 pp.queryAdd("horizontal_mixing_type", horiz_mixing_type_string);
379 if (amrex::toLower(horiz_mixing_type_string) == "analytical" ||
380 amrex::toLower(horiz_mixing_type_string) == "analytic") {
382 } else if (amrex::toLower(horiz_mixing_type_string) == "constant") {
384 } else {
385 amrex::Abort("Don't know this horizontal mixing type");
386 }
387 pp.queryAdd("visc2",visc2);
388 pp.queryAdd("tnu2_salt",tnu2_salt);
389 pp.queryAdd("tnu2_temp",tnu2_temp);
390 pp.queryAdd("tnu2_scalar",tnu2_scalar);
391 tnu2.resize(NCONS);
392 if (use_salt) {
393 tnu2[0] = tnu2_temp;
394 tnu2[1] = tnu2_salt;
395 tnu2[2] = tnu2_scalar;
396 } else {
397 tnu2[0] = tnu2_temp;
398 tnu2[1] = tnu2_scalar;
399 }
400
401 pp.queryAdd("Akk_bak", Akk_bak);
402 pp.queryAdd("Akp_bak", Akp_bak);
403 static std::string vert_mixing_type_string = "analytic";
404 static std::string gls_stability_type_string = "Canuto_A";
405 pp.queryAdd("vertical_mixing_type", vert_mixing_type_string);
406 pp.queryAdd("gls_stability_type", gls_stability_type_string);
407 if (amrex::toLower(vert_mixing_type_string) == "analytical" ||
408 amrex::toLower(vert_mixing_type_string) == "analytic") {
410 } else if (amrex::toLower(vert_mixing_type_string) == "gls") {
412 if (amrex::toLower(gls_stability_type_string) == "canuto_a") {
414 }
415 else if (amrex::toLower(gls_stability_type_string) == "canuto_b") {
417 }
418 else if (amrex::toLower(gls_stability_type_string) == "galperin") {
420 }
421 else {
422 amrex::Abort("Don't know this GLS stability type");
423 }
424 } else {
425 amrex::Abort("Don't know this vertical mixing type");
426 }
427 // Read in GLS params
429 pp.queryAdd("gls_P", gls_p);
430 pp.queryAdd("gls_M", gls_m);
431 pp.queryAdd("gls_N", gls_n);
432 pp.queryAdd("gls_Kmin", gls_Kmin);
433 pp.queryAdd("gls_Pmin", gls_Pmin);
434
435 pp.queryAdd("gls_cmu0", gls_cmu0);
436 pp.queryAdd("gls_c1", gls_c1);
437 pp.queryAdd("gls_c2", gls_c2);
438 pp.queryAdd("gls_c3m", gls_c3m);
439 pp.queryAdd("gls_c3p", gls_c3p);
440 pp.queryAdd("gls_sigk", gls_sigk);
441 pp.queryAdd("gls_sigp", gls_sigp);
443 gls_Gh0 = amrex::Real(0.0329); // 0.0329 GOTM, 0.0673 Burchard
444 gls_Ghcri = amrex::Real(0.03);
445 gls_L1 = amrex::Real(0.107);
446 gls_L2 = amrex::Real(0.0032);
447 gls_L3 = amrex::Real(0.0864);
448 gls_L4 = amrex::Real(0.12);
449 gls_L5 = amrex::Real(11.9);
450 gls_L6 = amrex::Real(0.4);
451 gls_L7 = amrex::Real(0.0);
452 gls_L8 = amrex::Real(0.48);
454 gls_Gh0 = amrex::Real(0.0444); // 0.044 GOTM, 0.0673 Burchard
455 gls_Ghcri = amrex::Real(0.0414);
456 gls_L1 = amrex::Real(0.127);
457 gls_L2 = amrex::Real(0.00336);
458 gls_L3 = amrex::Real(0.0906);
459 gls_L4 = amrex::Real(0.101);
460 gls_L5 = amrex::Real(11.2);
461 gls_L6 = amrex::Real(0.4);
462 gls_L7 = amrex::Real(0.0);
463 gls_L8 = amrex::Real(0.318);
464 } else {
465 gls_Gh0 = amrex::Real(0.028);
466 gls_Ghcri = amrex::Real(0.02);
467 }
468 }
469
470 // Read and compute inverse nudging coeffs from inputs given in days,
471 // and store in a vector corresponding to BdyVars enum
472 amrex::Real tnudg = amrex::Real(0.0);
473 amrex::Real znudg = amrex::Real(0.0);
474 amrex::Real m2nudg = amrex::Real(0.0);
475 amrex::Real m3nudg = amrex::Real(0.0);
476 pp.queryAdd("tnudg",tnudg);
477 pp.queryAdd("znudg",znudg);
478 pp.queryAdd("m2nudg",m2nudg);
479 pp.queryAdd("m3nudg",m3nudg);
480 pp.queryAdd("obcfac",obcfac);
481
483 nudg_coeff[BdyVars::u ] = (m3nudg > 0.0) ? amrex::Real(1.0) / (m3nudg * amrex::Real(86400.0)) : 0.0;//BdyVars::u
484 nudg_coeff[BdyVars::v ] = (m3nudg > 0.0) ? amrex::Real(1.0) / (m3nudg * amrex::Real(86400.0)) : 0.0;//BdyVars::v
485 nudg_coeff[BdyVars::t ] = ( tnudg > 0.0) ? amrex::Real(1.0) / ( tnudg * amrex::Real(86400.0)) : 0.0;//BdyVars::t
486 nudg_coeff[BdyVars::s ] = ( tnudg > 0.0) ? amrex::Real(1.0) / ( tnudg * amrex::Real(86400.0)) : 0.0;//BdyVars::s
487 nudg_coeff[BdyVars::ubar ] = (m2nudg > 0.0) ? amrex::Real(1.0) / (m2nudg * amrex::Real(86400.0)) : 0.0;//BdyVars::ubar
488 nudg_coeff[BdyVars::vbar ] = (m2nudg > 0.0) ? amrex::Real(1.0) / (m2nudg * amrex::Real(86400.0)) : 0.0;//BdyVars::vbar
489 nudg_coeff[BdyVars::zeta ] = ( znudg > 0.0) ? amrex::Real(1.0) / ( znudg * amrex::Real(86400.0)) : 0.0;//BdyVars::zeta
490
491 pp.queryAdd("do_m3_clim_nudg", do_m3_clim_nudg);
492 pp.queryAdd("do_m2_clim_nudg", do_m2_clim_nudg);
493 pp.queryAdd("do_temp_clim_nudg", do_temp_clim_nudg);
494 pp.queryAdd("do_salt_clim_nudg", do_salt_clim_nudg);
495
497 do_any_clim_nudg = true;
498 }
499#ifndef REMORA_USE_NETCDF
500 if (do_any_clim_nudg) {
501 amrex::Abort("Climatology nudging requires building with NetCDF");
502 }
503#endif
504 }
505
506 void display()
507 {
508 amrex::Print() << "SOLVER CHOICE: " << std::endl;
509 amrex::Print() << "use_salt : " << use_salt << std::endl;
510 amrex::Print() << "use_coriolis : " << use_coriolis << std::endl;
511 amrex::Print() << "use_prestep : " << use_prestep << std::endl;
512 amrex::Print() << "use_uv3dmix : " << use_uv3dmix << std::endl;
513 amrex::Print() << "use_barotropic : " << use_barotropic << std::endl;
514 amrex::Print() << "flat_bathymetry : " << flat_bathymetry << std::endl;
515 amrex::Print() << "spatial_order : " << spatial_order << std::endl;
516
517 if (ic_type == IC_Type::analytic) {
518 amrex::Print() << "Using analytic initial onditions" << std::endl;
519 }
520 else if (ic_type == IC_Type::netcdf) {
521 amrex::Print() << "Using NetCDF initial conditions" << std::endl;
522 }
523
525 amrex::Print() << "Horizontal advection scheme for tracers: " << "Centered 4" << std::endl;
526 }
528 amrex::Print() << "Horizontal advection scheme for tracers: " << "Upstream 3" << std::endl;
529 }
530 else {
531 amrex::Error("Invalid horizontal advection scheme for tracers.");
532 }
533
535 amrex::Print() << "Horizontal advection scheme for momenta: " << "Centered 2" << std::endl;
536 }
538 amrex::Print() << "Horizontal advection scheme for momenta: " << "Upstream 3" << std::endl;
539 }
540 else {
541 amrex::Error("Invalid horizontal advection scheme for momenta.");
542 }
543
545 amrex::Print() << "Using two-way coupling " << std::endl;
546 } else if (coupling_type == CouplingType::one_way) {
547 amrex::Print() << "Using one-way coupling " << std::endl;
548 }
549
550 if (use_coriolis) {
552 amrex::Print() << "Using analytic coriolis forcing " << std::endl;
553 } else if (coriolis_type == Cor_Type::beta_plane) {
554 amrex::Print() << "Using beta plane coriolis forcing " << std::endl;
555 } else if (coriolis_type == Cor_Type::netcdf) {
556 amrex::Print() << "Using coriolis forcing loaded from file " << std::endl;
557 }
558 }
559 }
560
561 // Default prefix
562 std::string pp_prefix {"remora"};
563
564 bool flat_bathymetry = false;
565
566 bool use_salt = true;
567
568 // Specify what additional physics/forcing modules we use
569 bool use_coriolis = false;
570
571 // Specify whether terms are used for debugging purposes
572 bool use_prestep = true;
573 bool use_uv3dmix = true;
574 bool use_baroclinic = true;
575 bool use_barotropic = true;
576
577 bool bulk_fluxes = false;
578
579 bool output_forcing = false;
580 bool do_temp_flux = false;
581 bool do_salt_flux = false;
582 bool longwave_down = false;
583
584 // Control flags for loading atmospheric variables from NetCDF
585 bool Tair_from_netcdf = false;
586 bool qair_from_netcdf = false;
587 bool Pair_from_netcdf = false;
588 bool srflx_from_netcdf = false;
591 bool rain_from_netcdf = false;
592 bool cloud_from_netcdf = false;
594 bool qair_is_percent = false;
595 std::string longwave_netcdf_varname = "lwrad";
596
597
598 bool do_rivers = false;
599 bool do_rivers_temp = true;
600 bool do_rivers_salt = true;
601 bool do_rivers_scalar = false;
602 amrex::Vector<int> do_rivers_cons;
603
604 bool init_l1ad_h = false;
605 bool init_l1ad_T = false;
606
607 bool init_ana_h = false;
608 bool init_ana_T = false;
609
610 bool init_l0int_h = true;
611 bool init_l0int_T = true;
612
615
616 // Coupling options are "OneWay" or "TwoWay"
618
619 // IC and BC Type: "analytic" or "netcdf"
621
622 // Coriolis forcing type
624
625 // Surface momentum flux type
627
628 // Surface wind speed type
630
631 // EOS type
633
634 // Bottom stress type
636
637 // Land/sea mask type
639
640 // Mixing type and parameters
644
645 // Type for grid scale (pm and pn)
647
648 // Stretching and depth parameters which may need to be read from inputs
649 amrex::Real theta_s = amrex::Real(3.0);
650 amrex::Real theta_b = amrex::Real(0.0);
651 amrex::Real tcline = amrex::Real(150.0);
652
653 // Linear drag coefficient [m/s]
654 amrex::Real rdrag = amrex::Real(3e-4);
655 // Quadratic drag coefficient [dimensionless]
656 amrex::Real rdrag2 = amrex::Real(3e-3);
657
658 // Momentum stress scales [m]
659 amrex::Real Zob = amrex::Real(2e-2);
660 amrex::Real Zos = amrex::Real(2e-2);
661
662 amrex::Real Cdb_max = amrex::Real(0.5);
663 amrex::Real Cdb_min = amrex::Real(1e-6);
664
665 // Linear equation of state parameters
666 amrex::Real R0 = amrex::Real(1028); // background density value (Kg/m3) used in Linear Equation of State
667 amrex::Real S0 = amrex::Real(35.0); // background salinity (nondimensional) constant
668 amrex::Real T0 = amrex::Real(5.0); // background potential temperature (Celsius) constant
669 amrex::Real Tcoef = amrex::Real(1.7e-4); // linear equation of state parameter (1/Celsius)
670 amrex::Real Scoef = amrex::Real(0.0); // linear equation of state parameter (nondimensional)
671 amrex::Real rho0 = amrex::Real(1025.0); // Mean density (Kg/m3) used when Boussinesq approx is inferred
672
673 // Coriolis forcing
674 amrex::Real coriolis_f0 = amrex::Real(0.0); // f-plane constant (1/s)
675 amrex::Real coriolis_beta = amrex::Real(0.0); // beta-plane constant (1/s/m)
676
677 // Air pressure
678 amrex::Real Pair = amrex::Real(1013.48);
679 // Air temperature
680 amrex::Real Tair = amrex::Real(23.567);
681 // Relative humidity (air)
682 amrex::Real Hair = amrex::Real(0.776);
683 // Cloud cover fraction (0=clear sky, 1=overcast)
684 amrex::Real cloud = amrex::Real(0.0);
685 // Precipitation rate (kg/m2/s)
686 amrex::Real rain = amrex::Real(0.0);
687 // Height (m) of atmospheric measurements for Bulk fluxes parametrization
688 amrex::Real blk_ZQ = amrex::Real(10.0); // air humidity
689 amrex::Real blk_ZT = amrex::Real(10.0); // air temperature
690 amrex::Real blk_ZW = amrex::Real(10.0); // winds
691
692 bool eminusp = false;
694
695 // Surface radiation flux
696 amrex::Real srflux = amrex::Real(0.0);
697
698 // Spatial discretization
700
701 // Horizontal mixing parameters
702 amrex::Real visc2 = amrex::Real(0.0);
703 amrex::Vector<amrex::Real> tnu2;
704
705 // GLS params
706 amrex::Real gls_p = amrex::Real(3.0);
707 amrex::Real gls_m = amrex::Real(1.5);
708 amrex::Real gls_n = amrex::Real(-1.0);
709 amrex::Real gls_Kmin = amrex::Real(7.6e-6);
710 amrex::Real gls_Pmin = amrex::Real(1.0e-12);
711
712 amrex::Real gls_cmu0 = amrex::Real(0.5477);
713 amrex::Real gls_c1 = amrex::Real(1.44);
714 amrex::Real gls_c2 = amrex::Real(1.92);
715 amrex::Real gls_c3m = amrex::Real(-0.4);
716 amrex::Real gls_c3p = amrex::Real(1.0);
717 amrex::Real gls_sigk = amrex::Real(1.0);
718 amrex::Real gls_sigp = amrex::Real(1.3);
719
720 // Turbulence closure
721 amrex::Real Akk_bak = amrex::Real(5.0e-6);
722 amrex::Real Akp_bak = amrex::Real(5.0e-6);
723 amrex::Real Akv_bak = amrex::Real(5.0e-6);
724 amrex::Real Akt_bak = amrex::Real(1.0e-6); // Note: this is a vector with one component per tracer in ROMS
725
726 // Params for stability functions.
727 amrex::Real gls_Gh0;
728 amrex::Real gls_Ghcri;
729 amrex::Real gls_Ghmin = amrex::Real(-0.28);
730 amrex::Real gls_E2 = amrex::Real(1.33);
731 // Params only for Canuto stability
732 amrex::Real gls_L1;
733 amrex::Real gls_L2;
734 amrex::Real gls_L3;
735 amrex::Real gls_L4;
736 amrex::Real gls_L5;
737 amrex::Real gls_L6;
738 amrex::Real gls_L7;
739 amrex::Real gls_L8;
740
741 // Params for some GLS and also Mellor-Yamada
742 amrex::Real my_A1 = amrex::Real(0.92);
743 amrex::Real my_A2 = amrex::Real(0.74);
744 amrex::Real my_B1 = amrex::Real(16.6);
745 amrex::Real my_B2 = amrex::Real(10.1);
746 amrex::Real my_C1 = amrex::Real(0.08);
747 amrex::Real my_C2 = amrex::Real(0.7);
748 amrex::Real my_C3 = amrex::Real(0.2);
749 amrex::Real my_E1 = amrex::Real(1.8);
750 amrex::Real my_E2 = amrex::Real(1.33);
751 amrex::Real my_Gh0 = amrex::Real(0.0233);
752 amrex::Real my_Sq = amrex::Real(0.2);
753 amrex::Real my_dtfac = amrex::Real(0.05);
754 amrex::Real my_lmax = amrex::Real(0.53);
755 amrex::Real my_qmin = amrex::Real(1.0E-8);
756
757 // Nudging time scales in 1/s
758 amrex::Vector<amrex::Real> nudg_coeff;
759
760 // Factor between passive (outflow) and active (inflow) open boundary
761 // conditions.
762 amrex::Real obcfac = amrex::Real(0.0);
763
764 // Whether to do climatoogy nudging
765 bool do_m2_clim_nudg = false;
766 bool do_m3_clim_nudg = false;
767 bool do_temp_clim_nudg = false;
768 bool do_salt_clim_nudg = false;
769 bool do_any_clim_nudg = false;
770
772};
773
774#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 Tracer_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
std::string longwave_netcdf_varname
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