v0.15.0
Loading...
Searching...
No Matches
StandardPoisson Struct Reference
Collaboration diagram for StandardPoisson:
[legend]

Classes

struct  CommonData
 [Source function] More...
 
struct  OpError
 

Public Member Functions

 StandardPoisson (MoFEM::Interface &m_field)
 
MoFEMErrorCode runProgram ()
 

Private Member Functions

MoFEMErrorCode readMesh ()
 
MoFEMErrorCode setupProblem ()
 
MoFEMErrorCode boundaryCondition ()
 
MoFEMErrorCode createCommonData ()
 [Create common data]
 
MoFEMErrorCode assembleSystem ()
 [Create common data]
 
MoFEMErrorCode setIntegrationRules ()
 
MoFEMErrorCode solveSystem ()
 
MoFEMErrorCode checkError ()
 [Check error]
 
MoFEMErrorCode outputResults ()
 [Check error]
 

Static Private Member Functions

static double analyticalFunction (const double x, const double y, const double z)
 [Analytical function]
 
static VectorDouble analyticalFunctionGrad (const double x, const double y, const double z)
 [Analytical function]
 
static double sourceFunction (const double x, const double y, const double z)
 [Analytical function gradient]
 

Private Attributes

boost::shared_ptr< CommonDatacommonDataPtr
 
MoFEM::InterfacemField
 
SimplesimpleInterface
 
std::string domainField
 
int oRder
 
boost::shared_ptr< std::vector< unsigned char > > boundaryMarker
 
Range boundaryEntitiesForFieldsplit
 

Detailed Description

Definition at line 24 of file standard_poisson.cpp.

Constructor & Destructor Documentation

◆ StandardPoisson()

StandardPoisson::StandardPoisson ( MoFEM::Interface & m_field)

Definition at line 114 of file standard_poisson.cpp.

115 : domainField("U"), mField(m_field) {}
MoFEM::Interface & mField

Member Function Documentation

◆ analyticalFunction()

static double StandardPoisson::analyticalFunction ( const double x,
const double y,
const double z )
inlinestaticprivate

[Analytical function]

Definition at line 44 of file standard_poisson.cpp.

45 {
46 return exp(-100. * (sqr(x) + sqr(y))) * cos(M_PI * x) * cos(M_PI * y);
47 }
double sqr(double x)

◆ analyticalFunctionGrad()

static VectorDouble StandardPoisson::analyticalFunctionGrad ( const double x,
const double y,
const double z )
inlinestaticprivate

[Analytical function]

[Analytical function gradient]

Definition at line 51 of file standard_poisson.cpp.

52 {
53 VectorDouble res;
54 res.resize(2);
55 res[0] = -exp(-100. * (sqr(x) + sqr(y))) *
56 (200. * x * cos(M_PI * x) + M_PI * sin(M_PI * x)) * cos(M_PI * y);
57 res[1] = -exp(-100. * (sqr(x) + sqr(y))) *
58 (200. * y * cos(M_PI * y) + M_PI * sin(M_PI * y)) * cos(M_PI * x);
59 return res;
60 }

◆ assembleSystem()

MoFEMErrorCode StandardPoisson::assembleSystem ( )
private

[Create common data]

Definition at line 175 of file standard_poisson.cpp.

175 {
177
178 auto pipeline_mng = mField.getInterface<PipelineManager>();
179
180 { // Push operators to the Pipeline that is responsible for calculating LHS of
181 // domain elements
182 CHKERR AddHOOps<2, 2, 2>::add(pipeline_mng->getOpDomainLhsPipeline(), {H1});
183 auto unity = [](const double, const double, const double) { return 1; };
184 pipeline_mng->getOpDomainLhsPipeline().push_back(
186 }
187
188 { // Push operators to the Pipeline that is responsible for calculating RHS of
189 // domain elements
190 pipeline_mng->getOpDomainRhsPipeline().push_back(
192 }
193
195}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, 1 > OpDomainSource
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradGrad< 1, 1, SPACE_DIM > OpDomainGradGrad
Add operators pushing bases from local to physical configuration.
PipelineManager interface.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
static double sourceFunction(const double x, const double y, const double z)
[Analytical function gradient]

◆ boundaryCondition()

MoFEMErrorCode StandardPoisson::boundaryCondition ( )
private

Definition at line 144 of file standard_poisson.cpp.

144 {
146
147 auto bc_mng = mField.getInterface<BcManager>();
148
149 // Remove BCs from blockset name "BOUNDARY_CONDITION" or SETU, note that you
150 // can use regular expression to put list of blocksets;
152 simpleInterface->getProblemName(), "BOUNDARY", std::string(domainField),
153 true);
154
156}
Simple interface for fast problem set-up.
Definition BcManager.hpp:29
MoFEMErrorCode removeBlockDOFsOnEntities(const std::string problem_name, const std::string block_name, const std::string field_name, int lo, int hi, bool get_low_dim_ents=true, bool is_distributed_mesh=true)
Remove DOFs from problem.
Definition BcManager.cpp:72
const std::string getProblemName() const
Get the Problem Name.
Definition Simple.hpp:410

◆ checkError()

MoFEMErrorCode StandardPoisson::checkError ( )
private

[Check error]

Definition at line 242 of file standard_poisson.cpp.

242 {
245 pipeline_mng->getDomainLhsFE().reset();
246 pipeline_mng->getDomainRhsFE().reset();
247 pipeline_mng->getOpDomainRhsPipeline().clear();
248
250
251 pipeline_mng->getOpDomainRhsPipeline().push_back(
253 pipeline_mng->getOpDomainRhsPipeline().push_back(
255 commonDataPtr->approxValsGrad));
256 pipeline_mng->getOpDomainRhsPipeline().push_back(
258
259 CHKERR VecZeroEntries(commonDataPtr->petscVec);
260
261 CHKERR pipeline_mng->loopFiniteElements();
262
263 CHKERR VecAssemblyBegin(commonDataPtr->petscVec);
264 CHKERR VecAssemblyEnd(commonDataPtr->petscVec);
265 const double *array;
266 CHKERR VecGetArrayRead(commonDataPtr->petscVec, &array);
267 MOFEM_LOG("EXAMPLE", Sev::inform)
268 << "Global error L2 norm: " << std::sqrt(array[0]);
269 MOFEM_LOG("EXAMPLE", Sev::inform)
270 << "Global error H1 seminorm: " << std::sqrt(array[1]);
271 CHKERR VecRestoreArrayRead(commonDataPtr->petscVec, &array);
272
274}
MoFEMErrorCode loopFiniteElements(SmartPetscObj< DM > dm=nullptr)
Iterate finite elements.
boost::ptr_deque< UserDataOperator > & getOpDomainRhsPipeline()
Get the Op Domain Rhs Pipeline object.
#define MOFEM_LOG(channel, severity)
Log.
Get field gradients at integration pts for scalar field rank 0, i.e. vector field.
Get value at integration points for scalar field.
boost::shared_ptr< FEMethod > & getDomainRhsFE()
boost::shared_ptr< FEMethod > & getDomainLhsFE()
boost::shared_ptr< CommonData > commonDataPtr

◆ createCommonData()

MoFEMErrorCode StandardPoisson::createCommonData ( )
private

[Create common data]

Definition at line 159 of file standard_poisson.cpp.

159 {
161 commonDataPtr = boost::make_shared<CommonData>();
162 PetscInt ghosts[2] = {0, 1};
163 if (!mField.get_comm_rank())
164 commonDataPtr->petscVec =
165 createGhostVector(mField.get_comm(), 2, 2, 0, ghosts);
166 else
167 commonDataPtr->petscVec =
168 createGhostVector(mField.get_comm(), 0, 2, 2, ghosts);
169 commonDataPtr->approxVals = boost::make_shared<VectorDouble>();
170 commonDataPtr->approxValsGrad = boost::make_shared<MatrixDouble>();
172}
auto createGhostVector(MPI_Comm comm, PetscInt n, PetscInt N, PetscInt nghost, const PetscInt ghosts[])
Create smart ghost vector.
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0

◆ outputResults()

MoFEMErrorCode StandardPoisson::outputResults ( )
private

[Check error]

Definition at line 277 of file standard_poisson.cpp.

277 {
279
280 auto pipeline_mng = mField.getInterface<PipelineManager>();
281 pipeline_mng->getDomainLhsFE().reset();
282 pipeline_mng->getBoundaryLhsFE().reset();
283 pipeline_mng->getDomainRhsFE().reset();
284 pipeline_mng->getBoundaryRhsFE().reset();
285
287
288 auto post_proc_fe = boost::make_shared<PostProcFaceEle>(mField);
289
290 CHKERR AddHOOps<2, 2, 2>::add(post_proc_fe->getOpPtrVector(), {H1});
291
292 post_proc_fe->getOpPtrVector().push_back(
294 post_proc_fe->getOpPtrVector().push_back(
296 commonDataPtr->approxValsGrad));
297
298 post_proc_fe->getOpPtrVector().push_back(new OpPPMap(
299 post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
300 {{domainField, commonDataPtr->approxVals}},
301 {{domainField + "_GRAD", commonDataPtr->approxValsGrad}}, {}, {}));
302 pipeline_mng->getDomainRhsFE() = post_proc_fe;
303
304 CHKERR pipeline_mng->loopFiniteElements();
305 CHKERR post_proc_fe->writeFile("out_result.h5m");
306
308}
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Post post-proc data at points from hash maps.

◆ readMesh()

MoFEMErrorCode StandardPoisson::readMesh ( )
private

Definition at line 117 of file standard_poisson.cpp.

117 {
119
123
125}
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition Simple.cpp:191
MoFEMErrorCode getOptions()
get options
Definition Simple.cpp:180

◆ runProgram()

MoFEMErrorCode StandardPoisson::runProgram ( )

Definition at line 310 of file standard_poisson.cpp.

310 {
312
322
324}
MoFEMErrorCode outputResults()
[Check error]
MoFEMErrorCode setIntegrationRules()
MoFEMErrorCode createCommonData()
[Create common data]
MoFEMErrorCode checkError()
[Check error]
MoFEMErrorCode solveSystem()
MoFEMErrorCode assembleSystem()
[Create common data]
MoFEMErrorCode setupProblem()
MoFEMErrorCode boundaryCondition()
MoFEMErrorCode readMesh()

◆ setIntegrationRules()

MoFEMErrorCode StandardPoisson::setIntegrationRules ( )
private

Definition at line 197 of file standard_poisson.cpp.

197 {
199
200 auto pipeline_mng = mField.getInterface<PipelineManager>();
201
202 auto domain_rule_lhs = [](int, int, int p) -> int { return 2 * (p - 1); };
203 auto domain_rule_rhs = [](int, int, int p) -> int { return 2 * (p - 1); };
204 CHKERR pipeline_mng->setDomainLhsIntegrationRule(domain_rule_lhs);
205 CHKERR pipeline_mng->setDomainRhsIntegrationRule(domain_rule_rhs);
206
207 auto boundary_rule_lhs = [](int, int, int p) -> int { return 2 * p; };
208 auto boundary_rule_rhs = [](int, int, int p) -> int { return 2 * p; };
209 CHKERR pipeline_mng->setBoundaryLhsIntegrationRule(boundary_rule_lhs);
210 CHKERR pipeline_mng->setBoundaryLhsIntegrationRule(boundary_rule_rhs);
211
213}
MoFEMErrorCode setDomainLhsIntegrationRule(RuleHookFun rule)

◆ setupProblem()

MoFEMErrorCode StandardPoisson::setupProblem ( )
private

Definition at line 127 of file standard_poisson.cpp.

127 {
129
134
135 int oRder = 2;
136 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-order", &oRder, PETSC_NULLPTR);
138
140
142}
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base nme:nme847.
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:261
MoFEMErrorCode addBoundaryField(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 boundary.
Definition Simple.cpp:355
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition Simple.cpp:575
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition Simple.cpp:736

◆ solveSystem()

MoFEMErrorCode StandardPoisson::solveSystem ( )
private

Definition at line 215 of file standard_poisson.cpp.

215 {
217
218 auto pipeline_mng = mField.getInterface<PipelineManager>();
219
220 auto ksp_solver = pipeline_mng->createKSP();
221 CHKERR KSPSetFromOptions(ksp_solver);
222
223 // Create RHS and solution vectors
224 auto dm = simpleInterface->getDM();
225 auto F = createDMVector(dm);
226 auto D = vectorDuplicate(F);
227
228 CHKERR KSPSetUp(ksp_solver);
229
230 // Solve the system
231 CHKERR KSPSolve(ksp_solver, F, D);
232
233 // Scatter result data on the mesh
234 CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
235 CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
236 CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
237
239}
@ 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:514
auto createDMVector(DM dm)
Get smart vector from DM.
Definition DMMoFEM.hpp:1102
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:800

◆ sourceFunction()

static double StandardPoisson::sourceFunction ( const double x,
const double y,
const double z )
inlinestaticprivate

[Analytical function gradient]

[Source function]

Definition at line 64 of file standard_poisson.cpp.

64 {
65 return -exp(-100. * (sqr(x) + sqr(y))) *
66 (400. * M_PI *
67 (x * cos(M_PI * y) * sin(M_PI * x) +
68 y * cos(M_PI * x) * sin(M_PI * y)) +
69 2. * (20000. * (sqr(x) + sqr(y)) - 200. - sqr(M_PI)) *
70 cos(M_PI * x) * cos(M_PI * y));
71 }

Member Data Documentation

◆ boundaryEntitiesForFieldsplit

Range StandardPoisson::boundaryEntitiesForFieldsplit
private

Definition at line 111 of file standard_poisson.cpp.

◆ boundaryMarker

boost::shared_ptr<std::vector<unsigned char> > StandardPoisson::boundaryMarker
private

Definition at line 108 of file standard_poisson.cpp.

◆ commonDataPtr

boost::shared_ptr<CommonData> StandardPoisson::commonDataPtr
private

Definition at line 82 of file standard_poisson.cpp.

◆ domainField

std::string StandardPoisson::domainField
private

Definition at line 104 of file standard_poisson.cpp.

◆ mField

MoFEM::Interface& StandardPoisson::mField
private

Definition at line 100 of file standard_poisson.cpp.

◆ oRder

int StandardPoisson::oRder
private

Definition at line 105 of file standard_poisson.cpp.

◆ simpleInterface

Simple* StandardPoisson::simpleInterface
private

Definition at line 101 of file standard_poisson.cpp.


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