cuIBM
A GPU-based Immersed Boundary Method code
parseSimulationFile.cu
Go to the documentation of this file.
1 
8 #include <fstream>
9 
10 #include "io.h"
11 #include "utilities/parameterDB.h"
12 #include "yaml-cpp/yaml.h"
13 
14 
19 namespace io
20 {
21 
30 {
31  if (s == "EULER_EXPLICIT")
32  return EULER_EXPLICIT;
33  else if (s == "EULER_IMPLICIT")
34  return EULER_IMPLICIT;
35  else if (s == "ADAMS_BASHFORTH_2")
36  return ADAMS_BASHFORTH_2;
37  else if (s == "RUNGE_KUTTA_3")
38  return RUNGE_KUTTA_3;
39  else if (s == "CRANK_NICOLSON")
40  return CRANK_NICOLSON;
41  else
42  {
43  printf("Error: Unknown timeScheme '%s'!\n", s.c_str());
44  exit(-1);
45  }
46 } // timeSchemeFromString
47 
48 
57 {
58  if (s == "NONE")
59  return NONE;
60  else if (s == "DIAGONAL")
61  return DIAGONAL;
62  else if (s == "SMOOTHED_AGGREGATION")
63  return SMOOTHED_AGGREGATION;
64  else if (s == "AINV")
65  return AINV;
66  else
67  {
68  printf("Error: Unknown preconditioner '%s'!\n", s.c_str());
69  exit(-1);
70  };
71 } // preconditionerTypeFromString
72 
73 
82 {
83  if (s == "NAVIER_STOKES")
84  return NAVIER_STOKES;
85  else if (s == "SAIKI_BIRINGEN")
86  return SAIKI_BIRINGEN;
87  else if (s == "DIRECT_FORCING")
88  return DIRECT_FORCING;
89  else if (s == "TAIRA_COLONIUS")
90  return TAIRA_COLONIUS;
91  else if (s == "FADLUN_ET_AL")
92  return FADLUN_ET_AL;
93  else if (s == "DIFFUSION")
94  return DIFFUSION;
95  else if (s == "DF_MODIFIED")
96  return DF_MODIFIED;
97  else if (s == "FEA_MODIFIED")
98  return FEA_MODIFIED;
99  else if (s == "DF_IMPROVED")
100  return DF_IMPROVED;
101  else
102  {
103  printf("Error: Unknown ibmScheme '%s'!\n", s.c_str());
104  exit(-1);
105  }
106 } // ibmSchemeFromString
107 
108 
117 {
118  if (s == "CONSTANT")
119  return CONSTANT;
120  else if (s == "LINEAR")
121  return LINEAR;
122  else if (s == "QUADRATIC")
123  return QUADRATIC;
124  else
125  {
126  printf("Error: Unknown interpolationType '%s'!\n", s.c_str());
127  exit(-1);
128  };
129 } // interpolationTypeFromString
130 
131 
138 void parseSimulation(const YAML::Node &node, parameterDB &DB)
139 {
140  // write to DB
141  std::string dbKey = "simulation";
142  DB[dbKey]["dt"].set<real>(node["dt"].as<real>(0.02));
143  DB[dbKey]["scaleCV"].set<real>(node["scaleCV"].as<real>(2.0));
144  DB[dbKey]["startStep"].set<int>(node["startStep"].as<int>(0));
145  DB[dbKey]["nt"].set<int>(node["nt"].as<int>(100));
146  DB[dbKey]["nsave"].set<int>(node["nsave"].as<int>(100));
147  DB[dbKey]["ibmScheme"].set<ibmScheme>(
148  ibmSchemeFromString(node["ibmScheme"].as<std::string>("NAVIER_STOKES")));
149  DB[dbKey]["convTimeScheme"].set<timeScheme>(
150  timeSchemeFromString(node["timeScheme"][0].as<std::string>("EULER_EXPLICIT")));
151  DB[dbKey]["diffTimeScheme"].set<timeScheme>(
152  timeSchemeFromString(node["timeScheme"][1].as<std::string>("EULER_IMPLICIT")));
153  DB[dbKey]["interpolationType"].set<interpolationType>(
154  interpolationTypeFromString(node["interpolationType"].as<std::string>("LINEAR")));
155 
156  const YAML::Node &solvers = node["linearSolvers"];
157  for (unsigned int i=0; i<solvers.size(); i++)
158  {
159  std::string system = solvers[i]["system"].as<std::string>();
160  std::string dbKey = system + "Solve";
161  DB[dbKey]["solver"].set<std::string>(solvers[i]["solver"].as<std::string>("CG"));
162  if (DB[dbKey]["solver"].get<std::string>() == "GMRES")
163  DB[dbKey]["restart"].set<int>(solvers[i]["restart"].as<int>(50));
164  DB[dbKey]["preconditioner"].set<preconditionerType>(
165  preconditionerTypeFromString(solvers[i]["preconditioner"].as<std::string>("DIAGONAL")));
166  DB[dbKey]["rTol"].set<real>(solvers[i]["relTolerance"].as<real>(1.0E-05));
167  DB[dbKey]["aTol"].set<real>(solvers[i]["absTolerance"].as<real>(1.0E-50));
168  DB[dbKey]["maxIterations"].set<int>(solvers[i]["maxIterations"].as<int>(10000));
169  DB[dbKey]["monitor"].set<int>(solvers[i]["monitor"].as<bool>(false));
170  }
171 } // parseSimulation
172 
173 
180 void parseSimulationFile(std::string &simFile, parameterDB &DB)
181 {
182  printf("Parsing YAML file with simulation parameters ...\n");
183  YAML::Node nodes = YAML::LoadFile(simFile);
184  for(unsigned int i=0; i<nodes.size(); i++)
185  parseSimulation(nodes[i], DB);
186 } // parseSimulationFile
187 
188 } // End of namespace io
Fadlun et al. modified.
Definition: types.h:83
constant
Definition: types.h:94
double real
Is a float or a double depending on the machine precision.
Definition: types.h:116
ibmScheme ibmSchemeFromString(std::string s)
Converts a string to an IBM scheme.
fully discrete direct forcing method
Definition: types.h:78
ibmScheme
Specifies the immersed boundary method used to solve the flow.
Definition: types.h:74
implicit Euler method (first order)
Definition: types.h:63
Declaration of the functions of the namespace io.
timeScheme
Specifies the numerical scheme used for time-integration.
Definition: types.h:60
preconditionerType preconditionerTypeFromString(std::string s)
Converts a string to a preconditioner type.
void parseSimulationFile(std::string &simFile, parameterDB &DB)
Parses simParams.yaml and stores the simulation parameters.
preconditionerType
Specifies the type of preconditioner.
Definition: types.h:104
std::map< std::string, componentParameter > parameterDB
Map from a string to a componentParameter.
Definition: parameterDB.h:64
smoothed aggregation preconditioner
Definition: types.h:108
Contains functions related to I/O tasks.
no immersed bodies - Perot (1993)
Definition: types.h:76
Direct Forcing Improved.
Definition: types.h:84
Fadlun et al (2000)
Definition: types.h:79
Saiki & Biringen (1996)
Definition: types.h:77
Declaration of the class property.
quadratic
Definition: types.h:96
interpolationType
Specifies the type of interpolation.
Definition: types.h:92
Taira & Colonius (2007)
Definition: types.h:80
no preconditioner
Definition: types.h:106
approximate inverse preconditioner
Definition: types.h:109
second-order Adams-Bashforth scheme
Definition: types.h:64
third-order low storage Runge-Kutta method
Definition: types.h:65
Crank-Nicolson scheme (second order)
Definition: types.h:66
linear
Definition: types.h:95
timeScheme timeSchemeFromString(std::string s)
Converts a string to a time-integration scheme type.
diagonal preconditioner
Definition: types.h:107
void parseSimulation(const YAML::Node &node, parameterDB &DB)
Fills the database with the simulation parameters.
Diffusion.
Definition: types.h:81
explicit Euler method (first order)
Definition: types.h:62
Direct Forcing modified.
Definition: types.h:82
interpolationType interpolationTypeFromString(std::string s)
Converts a string to a interpolation type.