PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
timeintegration.cpp
Go to the documentation of this file.
1
8// here goes headers from our PetIBM
10
11namespace petibm
12{
13namespace timeintegration
14{
15PetscErrorCode TimeIntegrationBase::printInfo() const
16{
17 PetscFunctionBeginUser;
18
19 PetscErrorCode ierr;
20
21 std::string info;
22
23 info += (std::string(80, '=') + "\n");
24 info += ("Time Integration [" + name + "]\n");
25 info += (std::string(80, '=') + "\n");
26
27 info += ("\tScheme: " + scheme + "\n\n");
28 info +=
29 ("\tCoefficient of implicit term: " + std::to_string(implicitCoeff) +
30 "\n\n");
31 info += "\tCoefficients of Explicit terms: [";
32 for (auto it : explicitCoeffs) info += (std::to_string(it) + ", ");
33 info += "]\n\n";
34
35 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s", info.c_str()); CHKERRQ(ierr);
36
37 PetscFunctionReturn(0);
38} // printInfo
39
40PetscErrorCode createTimeIntegration(const std::string &name,
41 const YAML::Node &node,
42 type::TimeIntegration &integration)
43{
44 PetscFunctionBeginUser;
45
46 std::string scheme;
47
48 if (!node["parameters"].IsDefined())
49 {
50 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
51 "Could not find the key \"parameters\" in the YAML node "
52 "passed to the function \"createTimeIntegration\"!");
53 }
54
55 if (!node["parameters"][name].IsDefined())
56 {
57 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
58 "Could not find the key \"%s\" under sub-node \"parameters\" "
59 "in the YAML node passed to the function "
60 "\"createTimeIntegration\"!",
61 name.c_str());
62 }
63
64 scheme = node["parameters"][name].as<std::string>();
65
66 // set up coefficients.
67 if (scheme == "EULER_EXPLICIT")
68 integration = std::make_shared<Euler_Explicit>(name);
69 else if (scheme == "EULER_IMPLICIT")
70 integration = std::make_shared<Euler_Implicit>(name);
71 else if (scheme == "ADAMS_BASHFORTH_2")
72 integration = std::make_shared<Adams_Bashforth_2>(name);
73 else if (scheme == "CRANK_NICOLSON")
74 integration = std::make_shared<Crank_Nicolson>(name);
75 else
76 {
77 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
78 "The time integration scheme \"%s\" does not exist.\n",
79 scheme.c_str());
80 }
81
82 PetscFunctionReturn(0);
83} // createTimeIntegration
84
85} // end of namespace timeintegration
86} // end of namespace petibm
const type::RealVec1D explicitCoeffs
Coefficients of explicit terms.
PetscErrorCode printInfo() const
Print information to standard output.
const std::string scheme
Name of the scheme.
const PetscReal implicitCoeff
Coefficient of implicit term.
const std::string name
Name of current instance.
PetscErrorCode createTimeIntegration(const std::string &name, const YAML::Node &node, type::TimeIntegration &integration)
factory function for type::TimeIntegration.
std::shared_ptr< timeintegration::TimeIntegrationBase > TimeIntegration
Definition of type::TimeIntegration.
A toolbox for building flow solvers.
Definition: bodypack.h:52
Definition of TimeIntegration related classes.