v0.13.2
Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Private Attributes | List of all members
Poisson2DHomogeneous Struct Reference
Collaboration diagram for Poisson2DHomogeneous:
[legend]

Public Member Functions

 Poisson2DHomogeneous (MoFEM::Interface &m_field)
 
MoFEMErrorCode runProgram ()
 [Output results] More...
 

Private Member Functions

MoFEMErrorCode readMesh ()
 [Read mesh] More...
 
MoFEMErrorCode setupProblem ()
 [Read mesh] More...
 
MoFEMErrorCode boundaryCondition ()
 [Setup problem] More...
 
MoFEMErrorCode assembleSystem ()
 [Boundary condition] More...
 
MoFEMErrorCode setIntegrationRules ()
 [Assemble system] More...
 
MoFEMErrorCode solveSystem ()
 [Set integration rules] More...
 
MoFEMErrorCode outputResults ()
 [Solve system] More...
 

Private Attributes

MoFEM::InterfacemField
 
SimplesimpleInterface
 
int oRder
 

Detailed Description

Examples
poisson_2d_homogeneous.cpp.

Definition at line 27 of file poisson_2d_homogeneous.cpp.

Constructor & Destructor Documentation

◆ Poisson2DHomogeneous()

Poisson2DHomogeneous::Poisson2DHomogeneous ( MoFEM::Interface m_field)
Examples
poisson_2d_homogeneous.cpp.

Definition at line 52 of file poisson_2d_homogeneous.cpp.

53 : mField(m_field) {}

Member Function Documentation

◆ assembleSystem()

MoFEMErrorCode Poisson2DHomogeneous::assembleSystem ( )
private

[Boundary condition]

[Assemble system]

Examples
poisson_2d_homogeneous.cpp.

Definition at line 101 of file poisson_2d_homogeneous.cpp.

101 {
103
104 auto pipeline_mng = mField.getInterface<PipelineManager>();
105
106 { // Push operators to the Pipeline that is responsible for calculating LHS
108 pipeline_mng->getOpDomainLhsPipeline(), {H1});
109 pipeline_mng->getOpDomainLhsPipeline().push_back(
111 }
112
113 { // Push operators to the Pipeline that is responsible for calculating RHS
114
115 auto set_values_to_bc_dofs = [&](auto &fe) {
116 auto get_bc_hook = [&]() {
118 return hook;
119 };
120 fe->preProcessHook = get_bc_hook();
121 };
122
123 // you can skip that if boundary condition is prescribing zero
124 auto calculate_residual_from_set_values_on_bc = [&](auto &pipeline) {
125 using DomainEle =
130
131 auto grad_u_vals_ptr = boost::make_shared<MatrixDouble>();
132 pipeline_mng->getOpDomainRhsPipeline().push_back(
134 grad_u_vals_ptr));
135 pipeline_mng->getOpDomainRhsPipeline().push_back(
136 new OpInternal(field_name, grad_u_vals_ptr,
137 [](double, double, double) constexpr { return -1; }));
138 };
139
141 pipeline_mng->getOpDomainRhsPipeline(), {H1});
142 set_values_to_bc_dofs(pipeline_mng->getDomainRhsFE());
143 calculate_residual_from_set_values_on_bc(
144 pipeline_mng->getOpDomainRhsPipeline());
145 pipeline_mng->getOpDomainRhsPipeline().push_back(
147 }
148
150}
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
constexpr auto field_name
constexpr int SPACE_DIM
Add operators pushing bases from local to physical configuration.
Class (Function) to enforce essential constrains.
Definition: Essential.hpp:39
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
PipelineManager interface.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

◆ boundaryCondition()

MoFEMErrorCode Poisson2DHomogeneous::boundaryCondition ( )
private

[Setup problem]

[Boundary condition]

Examples
poisson_2d_homogeneous.cpp.

Definition at line 85 of file poisson_2d_homogeneous.cpp.

85 {
87
88 auto bc_mng = mField.getInterface<BcManager>();
89
90 // Remove BCs from blockset name "BOUNDARY_CONDITION" or SETU, note that you
91 // can use regular expression to put list of blocksets;
92 CHKERR bc_mng->removeBlockDOFsOnEntities<BcScalarMeshsetType<BLOCKSET>>(
93 simpleInterface->getProblemName(), "(BOUNDARY_CONDITION|SETU)",
94 std::string(field_name), true);
95
97}
Simple interface for fast problem set-up.
Definition: BcManager.hpp:24
const std::string getProblemName() const
Get the Problem Name.
Definition: Simple.hpp:348

◆ outputResults()

MoFEMErrorCode Poisson2DHomogeneous::outputResults ( )
private

[Solve system]

[Output results]

Examples
poisson_2d_homogeneous.cpp.

Definition at line 196 of file poisson_2d_homogeneous.cpp.

196 {
198
199 auto pipeline_mng = mField.getInterface<PipelineManager>();
200 pipeline_mng->getDomainLhsFE().reset();
201
202 auto post_proc_fe = boost::make_shared<PostProcFaceEle>(mField);
204 post_proc_fe->getOpPtrVector(), {H1});
205
206 auto u_ptr = boost::make_shared<VectorDouble>();
207 auto grad_u_ptr = boost::make_shared<MatrixDouble>();
208 post_proc_fe->getOpPtrVector().push_back(
210
211 post_proc_fe->getOpPtrVector().push_back(
213
215
216 post_proc_fe->getOpPtrVector().push_back(
217
218 new OpPPMap(post_proc_fe->getPostProcMesh(),
219 post_proc_fe->getMapGaussPts(),
220
221 OpPPMap::DataMapVec{{"U", u_ptr}},
222
223 OpPPMap::DataMapMat{{"GRAD_U", grad_u_ptr}},
224
226
228
229 )
230
231 );
232
233 pipeline_mng->getDomainRhsFE() = post_proc_fe;
234 CHKERR pipeline_mng->loopFiniteElements();
235 CHKERR post_proc_fe->writeFile("out_result.h5m");
236
238}
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Get value at integration points for scalar field.
Post post-proc data at points from hash maps.
std::map< std::string, boost::shared_ptr< VectorDouble > > DataMapVec
std::map< std::string, boost::shared_ptr< MatrixDouble > > DataMapMat
boost::shared_ptr< FEMethod > & getDomainLhsFE()

◆ readMesh()

MoFEMErrorCode Poisson2DHomogeneous::readMesh ( )
private

[Read mesh]

Examples
poisson_2d_homogeneous.cpp.

Definition at line 56 of file poisson_2d_homogeneous.cpp.

56 {
58
62
64}
MoFEMErrorCode getOptions()
get options
Definition: Simple.cpp:180
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name)
Load mesh file.
Definition: Simple.cpp:194

◆ runProgram()

MoFEMErrorCode Poisson2DHomogeneous::runProgram ( )

[Output results]

[Run program]

Examples
poisson_2d_homogeneous.cpp.

Definition at line 242 of file poisson_2d_homogeneous.cpp.

242 {
244
252
254}
MoFEMErrorCode setupProblem()
[Read mesh]
MoFEMErrorCode boundaryCondition()
[Setup problem]
MoFEMErrorCode readMesh()
[Read mesh]
MoFEMErrorCode outputResults()
[Solve system]
MoFEMErrorCode assembleSystem()
[Boundary condition]
MoFEMErrorCode solveSystem()
[Set integration rules]
MoFEMErrorCode setIntegrationRules()
[Assemble system]

◆ setIntegrationRules()

MoFEMErrorCode Poisson2DHomogeneous::setIntegrationRules ( )
private

[Assemble system]

[Set integration rules]

Examples
poisson_2d_homogeneous.cpp.

Definition at line 154 of file poisson_2d_homogeneous.cpp.

154 {
156
157 auto rule_lhs = [](int, int, int p) -> int { return 2 * (p - 1); };
158 auto rule_rhs = [](int, int, int p) -> int { return p; };
159
160 auto pipeline_mng = mField.getInterface<PipelineManager>();
161 CHKERR pipeline_mng->setDomainLhsIntegrationRule(rule_lhs);
162 CHKERR pipeline_mng->setDomainRhsIntegrationRule(rule_rhs);
163
165}
static Index< 'p', 3 > p

◆ setupProblem()

MoFEMErrorCode Poisson2DHomogeneous::setupProblem ( )
private

[Read mesh]

[Setup problem]

Examples
poisson_2d_homogeneous.cpp.

Definition at line 68 of file poisson_2d_homogeneous.cpp.

68 {
70
73
74 int oRder = 3;
75 CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &oRder, PETSC_NULL);
77
79
81}
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
@ H1
continuous field
Definition: definitions.h:85
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
MoFEMErrorCode addDomainField(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_ZERO, int verb=-1)
Add field on domain.
Definition: Simple.cpp:264
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition: Simple.cpp:473
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition: Simple.cpp:611

◆ solveSystem()

MoFEMErrorCode Poisson2DHomogeneous::solveSystem ( )
private

[Set integration rules]

[Solve system]

Examples
poisson_2d_homogeneous.cpp.

Definition at line 169 of file poisson_2d_homogeneous.cpp.

169 {
171
172 auto pipeline_mng = mField.getInterface<PipelineManager>();
173
174 auto ksp_solver = pipeline_mng->createKSP();
175 CHKERR KSPSetFromOptions(ksp_solver);
176 CHKERR KSPSetUp(ksp_solver);
177
178 // Create RHS and solution vectors
179 auto dm = simpleInterface->getDM();
180 auto F = createDMVector(dm);
181 auto D = vectorDuplicate(F);
182
183 // Solve the system
184 CHKERR KSPSolve(ksp_solver, F, D);
185
186 // Scatter result data on the mesh
187 CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
188 CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
189 CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
190
192}
@ F
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMoFEM.cpp:509
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1003
SmartPetscObj< KSP > createKSP(SmartPetscObj< DM > dm=nullptr)
Create KSP (linear) solver.
double D
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition: Simple.cpp:667

Member Data Documentation

◆ mField

MoFEM::Interface& Poisson2DHomogeneous::mField
private
Examples
poisson_2d_homogeneous.cpp.

Definition at line 45 of file poisson_2d_homogeneous.cpp.

◆ oRder

int Poisson2DHomogeneous::oRder
private
Examples
poisson_2d_homogeneous.cpp.

Definition at line 49 of file poisson_2d_homogeneous.cpp.

◆ simpleInterface

Simple* Poisson2DHomogeneous::simpleInterface
private
Examples
poisson_2d_homogeneous.cpp.

Definition at line 46 of file poisson_2d_homogeneous.cpp.


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