PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
parser.h
Go to the documentation of this file.
1
8#pragma once
9
10#include <string>
11
12#include <petscsys.h>
13#include <symengine/lambda_double.h>
14#include <yaml-cpp/yaml.h>
15
16#include <petibm/type.h>
17
18namespace petibm
19{
25namespace parser
26{
54PetscErrorCode getSettings(YAML::Node &node);
55
66PetscErrorCode readYAMLFile(const std::string &filePath, YAML::Node &node);
67
83PetscErrorCode parseMesh(const YAML::Node &meshNode, PetscInt &dim,
85 type::IntVec1D &nTotal, type::RealVec2D &dL);
86
101PetscErrorCode parseOneAxis(const YAML::Node &axis, PetscInt &dir,
102 PetscReal &bg, PetscReal &ed, PetscInt &nTotal,
103 type::RealVec1D &dL);
104
118PetscErrorCode parseSubDomains(const YAML::Node &subs, const PetscReal bg,
119 PetscInt &nTotal, PetscReal &ed,
120 type::RealVec1D &dL);
121
133PetscErrorCode parseOneSubDomain(const YAML::Node &sub, const PetscReal bg,
134 PetscInt &n, PetscReal &ed,
135 type::RealVec1D &dL);
136
149PetscErrorCode parseBCs(const YAML::Node &node, type::IntVec2D &bcTypes,
150 type::RealVec2D &bcValues);
151
163PetscErrorCode parseICs(const YAML::Node &node, std::vector<SymEngine::LambdaRealDoubleVisitor> &ics);
164
165} // end of namespace parser
166
167} // end of namespace petibm
PetscErrorCode parseOneAxis(const YAML::Node &axis, PetscInt &dir, PetscReal &bg, PetscReal &ed, PetscInt &nTotal, type::RealVec1D &dL)
Parse the info of only one direction from YAML node.
Definition: parser.cpp:275
PetscErrorCode getSettings(YAML::Node &node)
Get configuration settings.
Definition: parser.cpp:175
PetscErrorCode parseSubDomains(const YAML::Node &subs, const PetscReal bg, PetscInt &nTotal, PetscReal &ed, type::RealVec1D &dL)
Parse all sub-domains in a direction.
Definition: parser.cpp:298
PetscErrorCode parseMesh(const YAML::Node &meshNode, PetscInt &dim, type::RealVec1D &bg, type::RealVec1D &ed, type::IntVec1D &nTotal, type::RealVec2D &dL)
Parse a YAML node of cartesianMesh.
Definition: parser.cpp:239
PetscErrorCode parseBCs(const YAML::Node &node, type::IntVec2D &bcTypes, type::RealVec2D &bcValues)
Parse boundary conditions from a YAML node.
Definition: parser.cpp:359
PetscErrorCode readYAMLFile(const std::string &filePath, YAML::Node &node)
Load the content of a YAML file into a YAML node.
Definition: parser.cpp:152
PetscErrorCode parseOneSubDomain(const YAML::Node &sub, const PetscReal bg, PetscInt &n, PetscReal &ed, type::RealVec1D &dL)
Parse only one sub-domain.
Definition: parser.cpp:334
PetscErrorCode parseICs(const YAML::Node &node, std::vector< SymEngine::LambdaRealDoubleVisitor > &ics)
Parse initial conditions from a YAML node.
Definition: parser.cpp:396
std::vector< IntVec1D > IntVec2D
2D std::vector holding PetscInt.
Definition: type.h:135
std::vector< RealVec1D > RealVec2D
2D std::vector holding PetscReal.
Definition: type.h:142
std::vector< PetscInt > IntVec1D
1D std::vector holding PetscInt.
Definition: type.h:133
std::vector< PetscReal > RealVec1D
1D std::vector holding PetscReal.
Definition: type.h:140
A toolbox for building flow solvers.
Definition: bodypack.h:52
Definition of user-defined types for convenience.