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 */
53
54/** \brief How to scale scaled_to_grid coefficients on AMR levels */
57};
58
59/** \brief stability function for GLS */
63
64/** \brief equation of state */
65enum class EOSType {
67};
68
69/** \brief bottom stress formulation */
73
74/** \brief initialization for pm and pn */
75enum class GridScaleType {
77};
78
79/** \brief surface momentum flux */
80enum class SMFluxType {
82};
83
84/** \brief surface wind */
85enum class WindType {
87};
88
89/** \brief masks */
90enum class MaskType {
92};
93
94/** \brief harmonic mixing; which surfaces to calculate along */
97};
98
100 public:
101 /** \brief read in and initialize parameters */
102 void init_params(int ncons)
103 {
104 amrex::ParmParse pp(pp_prefix);
105
106 pp.queryAdd("flat_bathymetry", flat_bathymetry);
107
108 // Which horizontal advection scheme for tracers
109 static std::string tracer_hadv_string = "upstream3";
110 pp.queryAdd("tracer_horizontal_advection_scheme",tracer_hadv_string);
111 if (tracer_hadv_string == "centered4")
113 else if (tracer_hadv_string == "upstream3")
115 else
116 amrex::Error("Advection scheme unknown.");
117
118 // Which horizontal advection scheme
119 static std::string uv_hadv_string = "upstream3";
120 pp.queryAdd("uv_horizontal_advection_scheme",uv_hadv_string);
121 if (uv_hadv_string == "upstream3")
123 else if (uv_hadv_string == "centered2")
125 else
126 amrex::Error("UV advection scheme unknown.");
127
128 pp.queryAdd("rdrag", rdrag);
129 pp.queryAdd("rdrag2", rdrag2);
130 pp.queryAdd("Zos", Zos);
131 pp.queryAdd("Zob", Zob);
132 pp.queryAdd("Cdb_max", Cdb_max);
133 pp.queryAdd("Cdb_min", Cdb_min);
134
135 // Include salinity?
136 pp.queryAdd("use_salt", use_salt);
137
138 // Include Coriolis forcing?
139 pp.queryAdd("use_coriolis", use_coriolis);
140
141 // Include prestep / lagged predictor / corrections
142 pp.queryAdd("use_prestep", use_prestep);
143
144 //This affect forcing and some of the mixing terms for velocity
145 pp.queryAdd("use_uv3dmix", use_uv3dmix);
146
147 //This accounts for the main 2d loop but may not include coupling and copying properly
148 pp.queryAdd("use_barotropic", use_barotropic);
149
150 // Whether to do rivers. By default, rivers are temp and salt sources. Rivers
151 // have to be momentum sources.
152 pp.queryAdd("do_rivers", do_rivers);
153
154 pp.queryAdd("do_rivers_temp", do_rivers_temp);
155 pp.queryAdd("do_rivers_salt", do_rivers_salt);
156 pp.queryAdd("do_rivers_scalar", do_rivers_scalar);
157
158 do_rivers_cons.assign(ncons, 0);
159 // If we aren't doing rivers, set all of these flags to false
162 for (int icomp = Tracer_comp; icomp < ncons; ++icomp) {
163 do_rivers_cons[icomp] = do_rivers ? do_rivers_scalar : false;
164 }
165
166 pp.queryAdd("init_l1ad_T", init_l1ad_T);
167
168 pp.queryAdd("init_ana_T", init_ana_T);
169
170 pp.queryAdd("init_l0int_T", init_l0int_T);
171
172 static std::string eos_type_string = "linear";
173 pp.queryAdd("eos_type",eos_type_string);
174 if (eos_type_string == "linear" || eos_type_string == "Linear" ||
175 eos_type_string == "lin" || eos_type_string == "Lin") {
177 pp.queryAdd("Tcoef",Tcoef);
178 pp.queryAdd("Scoef",Scoef);
179 } else if (eos_type_string == "nonlinear" || eos_type_string == "Nonlinear" ||
180 eos_type_string == "non-linear" || eos_type_string == "Non-linear" ||
181 eos_type_string == "nonlin" || eos_type_string == "Nonlin") {
183 } else {
184 amrex::Abort("Dont know this eos_type");
185 }
186 pp.queryAdd("R0",R0);
187 pp.queryAdd("S0",S0);
188 pp.queryAdd("T0",T0);
189 pp.queryAdd("rho0", rho0);
190
191 pp.queryAdd("bulk_fluxes",bulk_fluxes);
192 if (bulk_fluxes) {
193 do_salt_flux = true;
194 do_temp_flux = true;
195 }
196 // Outputs forcing variables if true
197 pp.queryAdd("output_forcing", output_forcing);
198 pp.queryAdd("air_pressure",Pair);
199 pp.queryAdd("air_temperature",Tair);
200 pp.queryAdd("air_humidity",Hair);
201 pp.queryAdd("surface_radiation_flux",srflux);
202 pp.queryAdd("srflx_from_netcdf",srflx_from_netcdf);
203 pp.queryAdd("longwave_down", longwave_down);
204 pp.queryAdd("longwave_down_from_netcdf",longwave_down_from_netcdf);
205 pp.queryAdd("longwave_netcdf_is_net", longwave_netcdf_is_net);
206 pp.queryAdd("longwave_netcdf_varname", longwave_netcdf_varname);
207 pp.queryAdd("cloud",cloud);
208 pp.queryAdd("rain",rain);
209 pp.queryAdd("blk_ZQ",blk_ZQ);
210 pp.queryAdd("blk_ZT",blk_ZT);
211 pp.queryAdd("blk_ZW",blk_ZW);
212 pp.queryAdd("eminusp",eminusp);
213 pp.queryAdd("eminusp_correct_ssh",eminusp_correct_ssh);
214 // Control flags for loading atmospheric variables from NetCDF
215 pp.queryAdd("Tair_from_netcdf",Tair_from_netcdf);
216 pp.queryAdd("qair_from_netcdf",qair_from_netcdf);
217 pp.queryAdd("qair_is_percent",qair_is_percent);
218 pp.queryAdd("Pair_from_netcdf",Pair_from_netcdf);
219 pp.queryAdd("rain_from_netcdf",rain_from_netcdf);
220 pp.queryAdd("cloud_from_netcdf",cloud_from_netcdf);
221 pp.queryAdd("EminusP_from_netcdf",EminusP_from_netcdf);
222
224 // If user provides net longwave from NetCDF, force use of file-based longwave.
225 longwave_down = true;
227 }
228
230 amrex::Warning("remora.longwave_netcdf_varname is set but remora.longwave_down_from_netcdf is false; the value will be ignored.");
231 }
232
233 if (eminusp_correct_ssh and !eminusp) {
234 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)");
235 }
236 if (eminusp and !bulk_fluxes) {
237 amrex::Abort("Evaporation minus precipitation (E-P) requires bulk flux parametrizations (remora.bulk_fluxes=true)");
238 }
239
240 //Read in linear eos parameters
241 //Grid stretching
242 pp.queryAdd("theta_s",theta_s);
243 pp.queryAdd("theta_b",theta_b);
244 pp.queryAdd("tcline",tcline);
245
246 //coriolis factor
247 pp.queryAdd("coriolis_f0",coriolis_f0);
248 pp.queryAdd("coriolis_beta",coriolis_beta);
249
250 pp.queryAdd("Akv_bak",Akv_bak);
251 pp.queryAdd("Akt_bak",Akt_bak);
252
253
254 static std::string grid_scale_type_string = "constant";
255 pp.queryAdd("grid_scale_type",grid_scale_type_string);
256
257 if (amrex::toLower(grid_scale_type_string) == "constant") {
259 } else if (amrex::toLower(grid_scale_type_string) == "custom") {
260 amrex::Warning("Initialization of grid scale from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
262 } else if (amrex::toLower(grid_scale_type_string) == "analytic") {
264 } else {
265 amrex::Error("Don't know this grid_scale_type");
266 }
267
268 static std::string ic_type_string = "analytic";
269 bool found_ic_bc = pp.queryAdd("ic_bc_type", ic_type_string);
270 pp.queryAdd("ic_type", ic_type_string);
271
272 if (found_ic_bc) {
273 amrex::Warning("remora.ic_bc_type is now called remora.ic_type, and will eventually be deprecated");
274 }
275
276 if ( amrex::toLower(ic_type_string) == "custom") {
277 amrex::Warning("Problem initialization from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
279 } else if ( amrex::toLower(ic_type_string) == "analytic") {
281 } else if ( amrex::toLower(ic_type_string) == "netcdf") {
283 } else if ( amrex::toLower(ic_type_string) == "real") {
284 amrex::Warning("Problem initialization from NetCDF (remora.ic_type) is now called 'netcdf'. 'real' will be deprecated");
286 } else {
287 amrex::Error("Don't know this ic_type");
288 }
289
290 // Which type of refinement
291 static std::string coupling_type_string = "TwoWay";
292 pp.queryAdd("coupling_type",coupling_type_string);
293 if (amrex::toLower(coupling_type_string) == "twoway" ||
294 amrex::toLower(coupling_type_string) == "two_way") {
296 } else if (amrex::toLower(coupling_type_string) == "oneway" ||
297 amrex::toLower(coupling_type_string) == "one_way") {
299 } else {
300 amrex::Abort("Dont know this coupling_type");
301 }
302
303 // Which type of coriolis forcing
304 if (use_coriolis) {
305 static std::string coriolis_type_string = "beta_plane";
306 pp.queryAdd("coriolis_type",coriolis_type_string);
307 if ( amrex::toLower(coriolis_type_string) == "custom") {
308 amrex::Warning("Coriolis initialization from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
310 } else if ( amrex::toLower(coriolis_type_string) == "analytic") {
312 } else if ((amrex::toLower(coriolis_type_string) == "beta_plane") ||
313 (amrex::toLower(coriolis_type_string) == "betaplane")) {
315 } else if ( (amrex::toLower(coriolis_type_string) == "netcdf")) {
317 } else if ( (amrex::toLower(coriolis_type_string) == "real")) {
318 amrex::Warning("Coriolis initialization from NetCDF is now called 'netcdf'. 'real' will be deprecated");
320 } else {
321 amrex::Abort("Don't know this coriolis_type");
322 }
323 }
324
325 static std::string smflux_type_string = "analytic";
326 int smflux_specified = pp.queryAdd("smflux_type",smflux_type_string);
327 if ( amrex::toLower(smflux_type_string) == "custom") {
328 amrex::Warning("Surface momentum flux initialization from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
330 } else if ( amrex::toLower(smflux_type_string) == "analytic") {
332 } else if ( amrex::toLower(smflux_type_string) == "netcdf") {
334 } else {
335 amrex::Abort("Don't know this smflux_type");
336 }
337
338 static std::string wind_type_string = "analytic";
339 int wind_specified = pp.queryAdd("wind_type",wind_type_string);
340 if ( amrex::toLower(wind_type_string) == "custom") {
341 amrex::Warning("Surface wind initialization from prob.cpp is now called 'analytic'. 'custom' will be deprecated");
343 } else if ( amrex::toLower(wind_type_string) == "analytic") {
345 } else if ( amrex::toLower(wind_type_string) == "netcdf") {
347 } else {
348 amrex::Abort("Don't know this smflux_type");
349 }
350
351 if (wind_specified && smflux_specified) {
352 amrex::Abort("Cannot specify both wind and surface momentum flux");
353 }
354
355 static std::string mask_type_string = "none";
356 // If initial condition type is netcdf, default to netcdf masks
357 if (ic_type == IC_Type::netcdf) {
358 mask_type_string = "netcdf";
359 }
360 pp.queryAdd("mask_type", mask_type_string);
361 if (amrex::toLower(mask_type_string) == "none") {
363 } else if (amrex::toLower(mask_type_string) == "analytic") {
365 } else if (amrex::toLower(mask_type_string) == "netcdf") {
367 } else {
368 amrex::Abort("Don't know this mask_type");
369 }
370
371 static std::string bottom_stress_type_string = "linear";
372 pp.queryAdd("bottom_stress_type", bottom_stress_type_string);
373 if (amrex::toLower(bottom_stress_type_string) == "linear") {
375 } else if (amrex::toLower(bottom_stress_type_string) == "quadratic") {
377 } else if (amrex::toLower(bottom_stress_type_string) == "logarithmic") {
379 } else {
380 amrex::Abort("Don't know this bottom_stress_type");
381 }
382
383 amrex::Real tnu2_salt = amrex::Real(0.0);
384 amrex::Real tnu2_temp = amrex::Real(0.0);
385 amrex::Real tnu2_scalar = amrex::Real(0.0);
386 static std::string horiz_mixing_type_string = "analytic";
387 pp.queryAdd("horizontal_mixing_type", horiz_mixing_type_string);
388 if (amrex::toLower(horiz_mixing_type_string) == "analytical" ||
389 amrex::toLower(horiz_mixing_type_string) == "analytic") {
391 } else if (amrex::toLower(horiz_mixing_type_string) == "constant") {
393 } else if (amrex::toLower(horiz_mixing_type_string) == "scaled_to_grid") {
395 } else {
396 amrex::Abort("Don't know this horizontal mixing type");
397 }
398 pp.queryAdd("visc2",visc2);
399 pp.queryAdd("tnu2_salt",tnu2_salt);
400 pp.queryAdd("tnu2_temp",tnu2_temp);
401 pp.queryAdd("tnu2_scalar",tnu2_scalar);
402
403 // For scaled_to_grid runs with AMR refinement: optionally scale the coefficients
404 // by the horizontal refinement ratio (linear in grid size).
405 static std::string scaled_to_grid_amr_scaling_string = "none";
406 pp.queryAdd("scaled_to_grid_amr_scaling", scaled_to_grid_amr_scaling_string);
407 if (amrex::toLower(scaled_to_grid_amr_scaling_string) == "none") {
409 } else if (amrex::toLower(scaled_to_grid_amr_scaling_string) == "linear") {
411 } else {
412 amrex::Abort("Don't know this scaled_to_grid_amr_scaling option");
413 }
414
415 tnu2.assign(ncons, tnu2_scalar);
416 if (ncons > Temp_comp) {
417 tnu2[Temp_comp] = tnu2_temp;
418 }
419 if (ncons > Salt_comp) {
420 tnu2[Salt_comp] = tnu2_salt;
421 }
422
423 static std::string harmonic_mixing_type_string = "s";
424 pp.queryAdd("harmonic_mixing_type", harmonic_mixing_type_string);
425 if (amrex::toLower(harmonic_mixing_type_string) == "s") {
427 } else if (amrex::toLower(harmonic_mixing_type_string) == "geopotential" ||
428 amrex::toLower(harmonic_mixing_type_string) == "geo") {
430 } else {
431 amrex::Abort("Don't know this harmonic_mixing_type");
432 }
433
434 pp.queryAdd("Akk_bak", Akk_bak);
435 pp.queryAdd("Akp_bak", Akp_bak);
436 static std::string vert_mixing_type_string = "analytic";
437 static std::string gls_stability_type_string = "Canuto_A";
438 pp.queryAdd("vertical_mixing_type", vert_mixing_type_string);
439 pp.queryAdd("gls_stability_type", gls_stability_type_string);
440 if (amrex::toLower(vert_mixing_type_string) == "analytical" ||
441 amrex::toLower(vert_mixing_type_string) == "analytic") {
443 } else if (amrex::toLower(vert_mixing_type_string) == "gls") {
445 if (amrex::toLower(gls_stability_type_string) == "canuto_a") {
447 }
448 else if (amrex::toLower(gls_stability_type_string) == "canuto_b") {
450 }
451 else if (amrex::toLower(gls_stability_type_string) == "galperin") {
453 }
454 else {
455 amrex::Abort("Don't know this GLS stability type");
456 }
457 } else {
458 amrex::Abort("Don't know this vertical mixing type");
459 }
460 // Read in GLS params
462 pp.queryAdd("gls_P", gls_p);
463 pp.queryAdd("gls_M", gls_m);
464 pp.queryAdd("gls_N", gls_n);
465 pp.queryAdd("gls_Kmin", gls_Kmin);
466 pp.queryAdd("gls_Pmin", gls_Pmin);
467
468 pp.queryAdd("gls_cmu0", gls_cmu0);
469 pp.queryAdd("gls_c1", gls_c1);
470 pp.queryAdd("gls_c2", gls_c2);
471 pp.queryAdd("gls_c3m", gls_c3m);
472 pp.queryAdd("gls_c3p", gls_c3p);
473 pp.queryAdd("gls_sigk", gls_sigk);
474 pp.queryAdd("gls_sigp", gls_sigp);
476 gls_Gh0 = amrex::Real(0.0329); // 0.0329 GOTM, 0.0673 Burchard
477 gls_Ghcri = amrex::Real(0.03);
478 gls_L1 = amrex::Real(0.107);
479 gls_L2 = amrex::Real(0.0032);
480 gls_L3 = amrex::Real(0.0864);
481 gls_L4 = amrex::Real(0.12);
482 gls_L5 = amrex::Real(11.9);
483 gls_L6 = amrex::Real(0.4);
484 gls_L7 = amrex::Real(0.0);
485 gls_L8 = amrex::Real(0.48);
487 gls_Gh0 = amrex::Real(0.0444); // 0.044 GOTM, 0.0673 Burchard
488 gls_Ghcri = amrex::Real(0.0414);
489 gls_L1 = amrex::Real(0.127);
490 gls_L2 = amrex::Real(0.00336);
491 gls_L3 = amrex::Real(0.0906);
492 gls_L4 = amrex::Real(0.101);
493 gls_L5 = amrex::Real(11.2);
494 gls_L6 = amrex::Real(0.4);
495 gls_L7 = amrex::Real(0.0);
496 gls_L8 = amrex::Real(0.318);
497 } else {
498 gls_Gh0 = amrex::Real(0.028);
499 gls_Ghcri = amrex::Real(0.02);
500 }
501 }
502
503 // Read and compute inverse nudging coeffs from inputs given in days,
504 // and store in a vector corresponding to BdyVars enum
505 amrex::Real tnudg = amrex::Real(0.0);
506 amrex::Real znudg = amrex::Real(0.0);
507 amrex::Real m2nudg = amrex::Real(0.0);
508 amrex::Real m3nudg = amrex::Real(0.0);
509 pp.queryAdd("tnudg",tnudg);
510 pp.queryAdd("znudg",znudg);
511 pp.queryAdd("m2nudg",m2nudg);
512 pp.queryAdd("m3nudg",m3nudg);
513 pp.queryAdd("obcfac",obcfac);
514
516 nudg_coeff[BdyVars::u ] = (m3nudg > 0.0) ? amrex::Real(1.0) / (m3nudg * amrex::Real(86400.0)) : 0.0;//BdyVars::u
517 nudg_coeff[BdyVars::v ] = (m3nudg > 0.0) ? amrex::Real(1.0) / (m3nudg * amrex::Real(86400.0)) : 0.0;//BdyVars::v
518 nudg_coeff[BdyVars::t ] = ( tnudg > 0.0) ? amrex::Real(1.0) / ( tnudg * amrex::Real(86400.0)) : 0.0;//BdyVars::t
519 nudg_coeff[BdyVars::s ] = ( tnudg > 0.0) ? amrex::Real(1.0) / ( tnudg * amrex::Real(86400.0)) : 0.0;//BdyVars::s
520 nudg_coeff[BdyVars::ubar ] = (m2nudg > 0.0) ? amrex::Real(1.0) / (m2nudg * amrex::Real(86400.0)) : 0.0;//BdyVars::ubar
521 nudg_coeff[BdyVars::vbar ] = (m2nudg > 0.0) ? amrex::Real(1.0) / (m2nudg * amrex::Real(86400.0)) : 0.0;//BdyVars::vbar
522 nudg_coeff[BdyVars::zeta ] = ( znudg > 0.0) ? amrex::Real(1.0) / ( znudg * amrex::Real(86400.0)) : 0.0;//BdyVars::zeta
523
524 pp.queryAdd("do_m3_clim_nudg", do_m3_clim_nudg);
525 pp.queryAdd("do_m2_clim_nudg", do_m2_clim_nudg);
526 pp.queryAdd("do_temp_clim_nudg", do_temp_clim_nudg);
527 pp.queryAdd("do_salt_clim_nudg", do_salt_clim_nudg);
528
530 do_any_clim_nudg = true;
531 }
532#ifndef REMORA_USE_NETCDF
533 if (do_any_clim_nudg) {
534 amrex::Abort("Climatology nudging requires building with NetCDF");
535 }
536#endif
537 }
538
539 void display()
540 {
541 amrex::Print() << "SOLVER CHOICE: " << std::endl;
542 amrex::Print() << "use_salt : " << use_salt << std::endl;
543 amrex::Print() << "use_coriolis : " << use_coriolis << std::endl;
544 amrex::Print() << "use_prestep : " << use_prestep << std::endl;
545 amrex::Print() << "use_uv3dmix : " << use_uv3dmix << std::endl;
546 amrex::Print() << "use_barotropic : " << use_barotropic << std::endl;
547 amrex::Print() << "flat_bathymetry : " << flat_bathymetry << std::endl;
548 amrex::Print() << "spatial_order : " << spatial_order << std::endl;
549
550 if (ic_type == IC_Type::analytic) {
551 amrex::Print() << "Using analytic initial onditions" << std::endl;
552 }
553 else if (ic_type == IC_Type::netcdf) {
554 amrex::Print() << "Using NetCDF initial conditions" << std::endl;
555 }
556
558 amrex::Print() << "Horizontal advection scheme for tracers: " << "Centered 4" << std::endl;
559 }
561 amrex::Print() << "Horizontal advection scheme for tracers: " << "Upstream 3" << std::endl;
562 }
563 else {
564 amrex::Error("Invalid horizontal advection scheme for tracers.");
565 }
566
568 amrex::Print() << "Horizontal advection scheme for momenta: " << "Centered 2" << std::endl;
569 }
571 amrex::Print() << "Horizontal advection scheme for momenta: " << "Upstream 3" << std::endl;
572 }
573 else {
574 amrex::Error("Invalid horizontal advection scheme for momenta.");
575 }
576
578 amrex::Print() << "Using two-way coupling " << std::endl;
579 } else if (coupling_type == CouplingType::one_way) {
580 amrex::Print() << "Using one-way coupling " << std::endl;
581 }
582
583 if (use_coriolis) {
585 amrex::Print() << "Using analytic coriolis forcing " << std::endl;
586 } else if (coriolis_type == Cor_Type::beta_plane) {
587 amrex::Print() << "Using beta plane coriolis forcing " << std::endl;
588 } else if (coriolis_type == Cor_Type::netcdf) {
589 amrex::Print() << "Using coriolis forcing loaded from file " << std::endl;
590 }
591 }
592 }
593
594 // Default prefix
595 std::string pp_prefix {"remora"};
596
597 bool flat_bathymetry = false;
598
599 bool use_salt = true;
600
601 // Specify what additional physics/forcing modules we use
602 bool use_coriolis = false;
603
604 // Specify whether terms are used for debugging purposes
605 bool use_prestep = true;
606 bool use_uv3dmix = true;
607 bool use_baroclinic = true;
608 bool use_barotropic = true;
609
610 bool bulk_fluxes = false;
611
612 bool output_forcing = false;
613 bool do_temp_flux = false;
614 bool do_salt_flux = false;
615 bool longwave_down = false;
616
617 // Control flags for loading atmospheric variables from NetCDF
618 bool Tair_from_netcdf = false;
619 bool qair_from_netcdf = false;
620 bool Pair_from_netcdf = false;
621 bool srflx_from_netcdf = false;
624 bool rain_from_netcdf = false;
625 bool cloud_from_netcdf = false;
627 bool qair_is_percent = false;
628 std::string longwave_netcdf_varname = "lwrad";
629
630
631 bool do_rivers = false;
632 bool do_rivers_temp = true;
633 bool do_rivers_salt = true;
634 bool do_rivers_scalar = false;
635 amrex::Vector<int> do_rivers_cons;
636
637 bool init_l1ad_T = false;
638
639 bool init_ana_T = false;
640
641 bool init_l0int_T = true;
642
645
646 // Coupling options are "OneWay" or "TwoWay"
648
649 // IC and BC Type: "analytic" or "netcdf"
651
652 // Coriolis forcing type
654
655 // Surface momentum flux type
657
658 // Surface wind speed type
660
661 // EOS type
663
664 // Bottom stress type
666
667 // Land/sea mask type
669
670 // Mixing type and parameters
676
677 // Type for grid scale (pm and pn)
679
680 // Stretching and depth parameters which may need to be read from inputs
681 amrex::Real theta_s = amrex::Real(3.0);
682 amrex::Real theta_b = amrex::Real(0.0);
683 amrex::Real tcline = amrex::Real(150.0);
684
685 // Linear drag coefficient [m/s]
686 amrex::Real rdrag = amrex::Real(3e-4);
687 // Quadratic drag coefficient [dimensionless]
688 amrex::Real rdrag2 = amrex::Real(3e-3);
689
690 // Momentum stress scales [m]
691 amrex::Real Zob = amrex::Real(2e-2);
692 amrex::Real Zos = amrex::Real(2e-2);
693
694 amrex::Real Cdb_max = amrex::Real(0.5);
695 amrex::Real Cdb_min = amrex::Real(1e-6);
696
697 // Linear equation of state parameters
698 amrex::Real R0 = amrex::Real(1028); // background density value (Kg/m3) used in Linear Equation of State
699 amrex::Real S0 = amrex::Real(35.0); // background salinity (nondimensional) constant
700 amrex::Real T0 = amrex::Real(5.0); // background potential temperature (Celsius) constant
701 amrex::Real Tcoef = amrex::Real(1.7e-4); // linear equation of state parameter (1/Celsius)
702 amrex::Real Scoef = amrex::Real(0.0); // linear equation of state parameter (nondimensional)
703 amrex::Real rho0 = amrex::Real(1025.0); // Mean density (Kg/m3) used when Boussinesq approx is inferred
704
705 // Coriolis forcing
706 amrex::Real coriolis_f0 = amrex::Real(0.0); // f-plane constant (1/s)
707 amrex::Real coriolis_beta = amrex::Real(0.0); // beta-plane constant (1/s/m)
708
709 // Air pressure
710 amrex::Real Pair = amrex::Real(1013.48);
711 // Air temperature
712 amrex::Real Tair = amrex::Real(23.567);
713 // Relative humidity (air)
714 amrex::Real Hair = amrex::Real(0.776);
715 // Cloud cover fraction (0=clear sky, 1=overcast)
716 amrex::Real cloud = amrex::Real(0.0);
717 // Precipitation rate (kg/m2/s)
718 amrex::Real rain = amrex::Real(0.0);
719 // Height (m) of atmospheric measurements for Bulk fluxes parametrization
720 amrex::Real blk_ZQ = amrex::Real(10.0); // air humidity
721 amrex::Real blk_ZT = amrex::Real(10.0); // air temperature
722 amrex::Real blk_ZW = amrex::Real(10.0); // winds
723
724 bool eminusp = false;
726
727 // Surface radiation flux
728 amrex::Real srflux = amrex::Real(0.0);
729
730 // Spatial discretization
732
733 // Horizontal mixing parameters
734 amrex::Real visc2 = amrex::Real(0.0);
735 amrex::Vector<amrex::Real> tnu2;
736
737 // GLS params
738 amrex::Real gls_p = amrex::Real(3.0);
739 amrex::Real gls_m = amrex::Real(1.5);
740 amrex::Real gls_n = amrex::Real(-1.0);
741 amrex::Real gls_Kmin = amrex::Real(7.6e-6);
742 amrex::Real gls_Pmin = amrex::Real(1.0e-12);
743
744 amrex::Real gls_cmu0 = amrex::Real(0.5477);
745 amrex::Real gls_c1 = amrex::Real(1.44);
746 amrex::Real gls_c2 = amrex::Real(1.92);
747 amrex::Real gls_c3m = amrex::Real(-0.4);
748 amrex::Real gls_c3p = amrex::Real(1.0);
749 amrex::Real gls_sigk = amrex::Real(1.0);
750 amrex::Real gls_sigp = amrex::Real(1.3);
751
752 // Turbulence closure
753 amrex::Real Akk_bak = amrex::Real(5.0e-6);
754 amrex::Real Akp_bak = amrex::Real(5.0e-6);
755 amrex::Real Akv_bak = amrex::Real(5.0e-6);
756 amrex::Real Akt_bak = amrex::Real(1.0e-6); // Note: this is a vector with one component per tracer in ROMS
757
758 // Params for stability functions.
759 amrex::Real gls_Gh0;
760 amrex::Real gls_Ghcri;
761 amrex::Real gls_Ghmin = amrex::Real(-0.28);
762 amrex::Real gls_E2 = amrex::Real(1.33);
763 // Params only for Canuto stability
764 amrex::Real gls_L1;
765 amrex::Real gls_L2;
766 amrex::Real gls_L3;
767 amrex::Real gls_L4;
768 amrex::Real gls_L5;
769 amrex::Real gls_L6;
770 amrex::Real gls_L7;
771 amrex::Real gls_L8;
772
773 // Params for some GLS and also Mellor-Yamada
774 amrex::Real my_A1 = amrex::Real(0.92);
775 amrex::Real my_A2 = amrex::Real(0.74);
776 amrex::Real my_B1 = amrex::Real(16.6);
777 amrex::Real my_B2 = amrex::Real(10.1);
778 amrex::Real my_C1 = amrex::Real(0.08);
779 amrex::Real my_C2 = amrex::Real(0.7);
780 amrex::Real my_C3 = amrex::Real(0.2);
781 amrex::Real my_E1 = amrex::Real(1.8);
782 amrex::Real my_E2 = amrex::Real(1.33);
783 amrex::Real my_Gh0 = amrex::Real(0.0233);
784 amrex::Real my_Sq = amrex::Real(0.2);
785 amrex::Real my_dtfac = amrex::Real(0.05);
786 amrex::Real my_lmax = amrex::Real(0.53);
787 amrex::Real my_qmin = amrex::Real(1.0E-8);
788
789 // Nudging time scales in 1/s
790 amrex::Vector<amrex::Real> nudg_coeff;
791
792 // Factor between passive (outflow) and active (inflow) open boundary
793 // conditions.
794 amrex::Real obcfac = amrex::Real(0.0);
795
796 // Whether to do climatoogy nudging
797 bool do_m2_clim_nudg = false;
798 bool do_m3_clim_nudg = false;
799 bool do_temp_clim_nudg = false;
800 bool do_salt_clim_nudg = false;
801 bool do_any_clim_nudg = false;
802
804};
805
806#endif
BottomStressType
bottom stress formulation
GridScaleType
initialization for pm and pn
HarmonicMixingType
harmonic mixing; which surfaces to calculate along
SMFluxType
surface momentum flux
AdvectionScheme
Horizontal advection schemes.
PlotfileType
plotfile format
HorizMixingType
horizontal viscosity/diffusion type
Coord
Coordinates.
MaskType
masks
ScaledToGridAMRScaling
How to scale scaled_to_grid coefficients on AMR levels.
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
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
ScaledToGridAMRScaling scaled_to_grid_amr_scaling
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
void init_params(int ncons)
read in and initialize parameters
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
HarmonicMixingType harmonic_mixing_type
amrex::Real Akp_bak
amrex::Real my_Sq
CouplingType coupling_type
amrex::Real gls_E2