REMORA
Energy Research and Forecasting: An Atmospheric Modeling Code
SolverChoice Struct Reference

#include <DataStruct.H>

Collaboration diagram for SolverChoice:

Public Member Functions

void init_params ()
 
void display ()
 

Public Attributes

std::string pp_prefix {"remora"}
 
bool flat_bathymetry = true
 
bool use_salt = true
 
bool use_coriolis = false
 
bool use_prestep = true
 
bool use_uv3dmix = true
 
bool use_baroclinic = true
 
bool use_barotropic = true
 
AdvectionScheme Hadv_scheme
 
CouplingType coupling_type
 
IC_BC_Type ic_bc_type
 
Cor_Type coriolis_type
 
amrex::Real theta_s = 3.0
 
amrex::Real theta_b = 0.0
 
amrex::Real tcline = 150.0
 
amrex::Real rdrag = 3e-4
 
amrex::Real R0 = 1028
 
amrex::Real S0 = 35.0
 
amrex::Real T0 = 5.0
 
amrex::Real Tcoef = 1.7e-4
 
amrex::Real Scoef = 0.0
 
amrex::Real rho0 = 1025.0
 
amrex::Real coriolis_f0 = 0.0
 
amrex::Real coriolis_beta = 0.0
 
amrex::Real g = 9.81
 
int spatial_order = 2
 

Member Function Documentation

◆ display()

void SolverChoice::display ( )
inline
142  {
143  amrex::Print() << "SOLVER CHOICE: " << std::endl;
144  amrex::Print() << "use_salt : " << use_salt << std::endl;
145  amrex::Print() << "use_coriolis : " << use_coriolis << std::endl;
146  amrex::Print() << "use_prestep : " << use_prestep << std::endl;
147  amrex::Print() << "use_uv3dmix : " << use_uv3dmix << std::endl;
148  amrex::Print() << "use_barotropic : " << use_barotropic << std::endl;
149  amrex::Print() << "flat_bathymetry : " << flat_bathymetry << std::endl;
150  amrex::Print() << "spatial_order : " << spatial_order << std::endl;
151 
153  amrex::Print() << "Using custom initial and boundary conditions (No mesoscale forcing!)" << std::endl;
154  }
155  else if (ic_bc_type == IC_BC_Type::Real) {
156  amrex::Print() << "Using REAL initial and boundary conditions (Mesoscale forcing!)" << std::endl;
157  }
158 
160  amrex::Print() << "Horizontal advection scheme: " << "Centered 4" << std::endl;
161  }
163  amrex::Print() << "Horizontal advection scheme: " << "Upstream 3" << std::endl;
164  }
165  else {
166  amrex::Error("Invalid horizontal advection scheme.");
167  }
168 
170  amrex::Print() << "Using two-way coupling " << std::endl;
171  } else if (coupling_type == CouplingType::OneWay) {
172  amrex::Print() << "Using one-way coupling " << std::endl;
173  }
174 
175  if (use_coriolis) {
177  amrex::Print() << "Using custom coriolis forcing " << std::endl;
178  } else if (coriolis_type == Cor_Type::Beta_Plane) {
179  amrex::Print() << "Using beta plane coriolis forcing " << std::endl;
180  } else if (coriolis_type == Cor_Type::Real) {
181  amrex::Print() << "Using coriolis forcing loaded from file " << std::endl;
182  }
183  }
184  }
bool use_uv3dmix
Definition: DataStruct.H:198
bool use_barotropic
Definition: DataStruct.H:200
IC_BC_Type ic_bc_type
Definition: DataStruct.H:208
AdvectionScheme Hadv_scheme
Definition: DataStruct.H:202
bool flat_bathymetry
Definition: DataStruct.H:189
Cor_Type coriolis_type
Definition: DataStruct.H:211
int spatial_order
Definition: DataStruct.H:237
bool use_coriolis
Definition: DataStruct.H:194
bool use_prestep
Definition: DataStruct.H:197
CouplingType coupling_type
Definition: DataStruct.H:205
bool use_salt
Definition: DataStruct.H:191

◆ init_params()

void SolverChoice::init_params ( )
inline
44  {
45  amrex::ParmParse pp(pp_prefix);
46 
47  pp.query("flat_bathymetry", flat_bathymetry);
48 
49  // Which horizontal advection scheme
50  static std::string hadv_string = "upstream3";
51  pp.query("horizontal_advection_scheme",hadv_string);
52  if (hadv_string == "centered4")
54  else if (hadv_string == "upstream3")
56  else
57  amrex::Error("Advection scheme unknown.");
58 
59  pp.query("rdrag", rdrag);
60 
61  // Order of spatial discretization
62  pp.query("spatial_order", spatial_order);
63 
64  // Include salinity?
65  pp.query("use_salt", use_salt);
66 
67  // Include Coriolis forcing?
68  pp.query("use_coriolis", use_coriolis);
69 
70  // Include prestep / lagged predictor / corrections
71  pp.query("use_prestep", use_prestep);
72 
73  //This affect forcing and some of the mixing terms for velocity
74  pp.query("use_uv3dmix", use_uv3dmix);
75 
76  //This accounts for the main 2d loop but may not include coupling and copying properly
77  pp.query("use_barotropic", use_barotropic);
78 
79  //Read in linear eos parameters
80  pp.query("R0",R0);
81  pp.query("S0",S0);
82  pp.query("T0",T0);
83  pp.query("Tcoef",Tcoef);
84  pp.query("Scoef",Scoef);
85  pp.query("rho0", rho0);
86 
87  //Grid stretching
88  pp.query("theta_s",theta_s);
89  pp.query("theta_b",theta_b);
90  pp.query("tcline",tcline);
91 
92  //coriolis factor
93  pp.query("coriolis_f0",coriolis_f0);
94  pp.query("coriolis_beta",coriolis_beta);
95 
96  pp.query("ggrav",g);
97 
98  static std::string ic_bc_type_string = "custom";
99  pp.query("ic_bc_type", ic_bc_type_string);
100 
101  if ( !ic_bc_type_string.compare("Custom") ||
102  !ic_bc_type_string.compare("custom") ) {
104  } else if ( !ic_bc_type_string.compare("Real") ||
105  !ic_bc_type_string.compare("real") ) {
107  } else {
108  amrex::Error("Don't know this ic_bc_type");
109  }
110 
111  // Which type of refinement
112  static std::string coupling_type_string = "TwoWay";
113  pp.query("coupling_type",coupling_type_string);
114  if (coupling_type_string == "TwoWay") {
116  } else if (coupling_type_string == "OneWay") {
118  } else {
119  amrex::Abort("Dont know this coupling_type");
120  }
121 
122  // Which type of coriolis forcing
123  if (use_coriolis) {
124  static std::string coriolis_type_string = "beta_plane";
125  pp.query("coriolis_type",coriolis_type_string);
126  if ( (coriolis_type_string == "Custom") || (coriolis_type_string == "custom") ) {
128  } else if ((coriolis_type_string == "Beta_Plane") ||
129  (coriolis_type_string == "beta_Plane") ||
130  (coriolis_type_string == "beta_plane") ||
131  (coriolis_type_string == "Beta_plane")) {
133  } else if ( (coriolis_type_string == "Real") || (coriolis_type_string == "real") ) {
135  } else {
136  amrex::Abort("Don't know this coriolis_type");
137  }
138  }
139  }
amrex::Real coriolis_beta
Definition: DataStruct.H:231
amrex::Real coriolis_f0
Definition: DataStruct.H:230
amrex::Real Tcoef
Definition: DataStruct.H:225
std::string pp_prefix
Definition: DataStruct.H:187
amrex::Real theta_b
Definition: DataStruct.H:215
amrex::Real theta_s
Definition: DataStruct.H:214
amrex::Real R0
Definition: DataStruct.H:222
amrex::Real tcline
Definition: DataStruct.H:216
amrex::Real rho0
Definition: DataStruct.H:227
amrex::Real rdrag
Definition: DataStruct.H:219
amrex::Real T0
Definition: DataStruct.H:224
amrex::Real Scoef
Definition: DataStruct.H:226
amrex::Real g
Definition: DataStruct.H:234
amrex::Real S0
Definition: DataStruct.H:223

Member Data Documentation

◆ coriolis_beta

amrex::Real SolverChoice::coriolis_beta = 0.0

Referenced by init_params().

◆ coriolis_f0

amrex::Real SolverChoice::coriolis_f0 = 0.0

Referenced by init_params().

◆ coriolis_type

Cor_Type SolverChoice::coriolis_type

Referenced by display(), and init_params().

◆ coupling_type

CouplingType SolverChoice::coupling_type

Referenced by display(), and init_params().

◆ flat_bathymetry

bool SolverChoice::flat_bathymetry = true

Referenced by display(), and init_params().

◆ g

amrex::Real SolverChoice::g = 9.81

Referenced by init_params().

◆ Hadv_scheme

AdvectionScheme SolverChoice::Hadv_scheme

Referenced by display(), and init_params().

◆ ic_bc_type

IC_BC_Type SolverChoice::ic_bc_type

Referenced by display(), and init_params().

◆ pp_prefix

std::string SolverChoice::pp_prefix {"remora"}

Referenced by init_params().

◆ R0

amrex::Real SolverChoice::R0 = 1028

Referenced by init_params().

◆ rdrag

amrex::Real SolverChoice::rdrag = 3e-4

Referenced by init_params().

◆ rho0

amrex::Real SolverChoice::rho0 = 1025.0

Referenced by init_params().

◆ S0

amrex::Real SolverChoice::S0 = 35.0

Referenced by init_params().

◆ Scoef

amrex::Real SolverChoice::Scoef = 0.0

Referenced by init_params().

◆ spatial_order

int SolverChoice::spatial_order = 2

Referenced by display(), and init_params().

◆ T0

amrex::Real SolverChoice::T0 = 5.0

Referenced by init_params().

◆ tcline

amrex::Real SolverChoice::tcline = 150.0

Referenced by init_params().

◆ Tcoef

amrex::Real SolverChoice::Tcoef = 1.7e-4

Referenced by init_params().

◆ theta_b

amrex::Real SolverChoice::theta_b = 0.0

Referenced by init_params().

◆ theta_s

amrex::Real SolverChoice::theta_s = 3.0

Referenced by init_params().

◆ use_baroclinic

bool SolverChoice::use_baroclinic = true

◆ use_barotropic

bool SolverChoice::use_barotropic = true

Referenced by display(), and init_params().

◆ use_coriolis

bool SolverChoice::use_coriolis = false

Referenced by display(), and init_params().

◆ use_prestep

bool SolverChoice::use_prestep = true

Referenced by display(), and init_params().

◆ use_salt

bool SolverChoice::use_salt = true

Referenced by display(), and init_params().

◆ use_uv3dmix

bool SolverChoice::use_uv3dmix = true

Referenced by display(), and init_params().


The documentation for this struct was generated from the following file: