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

Public Member Functions

 Poisson (moab::Core &mb_instance, MoFEM::Core &core, int order)
 
MoFEMErrorCode runAnalysis ()
 

Private Member Functions

MoFEMErrorCode readMesh ()
 
MoFEMErrorCode setupSystem ()
 
MoFEMErrorCode setIntegrationRule ()
 
MoFEMErrorCode assembleSystem ()
 
MoFEMErrorCode applyBc ()
 
MoFEMErrorCode solveSystem ()
 
MoFEMErrorCode setOutput ()
 

Private Attributes

boost::shared_ptr< VolElevolPipelineLhs
 
boost::shared_ptr< VolElevolPipelineRhs
 
boost::shared_ptr< VolElepostFaceVolume
 
MoFEM::InterfacemField
 
moab::Interface & moab
 
SimplesimpleInterface
 
DM dm
 
KSP ksp
 
MatrixDouble invJac
 
MPI_Comm mpiComm
 
const int mpiRank
 
int order
 
std::string domainField
 

Detailed Description

Definition at line 10 of file poisson_prob.cpp.

Constructor & Destructor Documentation

◆ Poisson()

Poisson::Poisson ( moab::Core &  mb_instance,
MoFEM::Core core,
int  order 
)

Definition at line 46 of file poisson_prob.cpp.

47 : domainField("U"), moab(mb_instance), mField(core),
49 order(order) {
50 volPipelineLhs = boost::shared_ptr<VolEle>(new VolEle(mField));
51 volPipelineRhs = boost::shared_ptr<VolEle>(new VolEle(mField));
52}
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
const int mpiRank
std::string domainField
MPI_Comm mpiComm
boost::shared_ptr< VolEle > volPipelineLhs
MoFEM::Interface & mField
boost::shared_ptr< VolEle > volPipelineRhs
moab::Interface & moab
VolumeElementForcesAndSourcesCore VolEle

Member Function Documentation

◆ applyBc()

MoFEMErrorCode Poisson::applyBc ( )
private

Definition at line 171 of file poisson_prob.cpp.

171 {
173 Range bdry_entities;
175 string name = it->getName();
176 if (name.compare(0, 14, "ESSENTIAL") == 0) {
177 CHKERR it->getMeshsetIdEntitiesByDimension(mField.get_moab(), 1,
178 bdry_entities, true);
179 }
180 }
181
182 cerr << bdry_entities;
183 Range bdr_verts;
184 CHKERR mField.get_moab().get_connectivity(bdry_entities, bdr_verts, true);
185 bdry_entities.merge(bdr_verts);
186
187 CHKERR mField.getInterface<ProblemsManager>()->removeDofsOnEntities(
188 simpleInterface->getProblemName(), domainField, bdry_entities);
189
190
192}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ BLOCKSET
Definition: definitions.h:148
#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
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
virtual moab::Interface & get_moab()=0
Problem manager is used to build and partition problems.
const std::string getProblemName() const
Get the Problem Name.
Definition: Simple.hpp:348
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Simple * simpleInterface

◆ assembleSystem()

MoFEMErrorCode Poisson::assembleSystem ( )
private

Definition at line 84 of file poisson_prob.cpp.

84 {
86 auto det_ptr = boost::make_shared<VectorDouble>();
87 auto jac_ptr = boost::make_shared<MatrixDouble>();
88 auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
89 volPipelineLhs->getOpPtrVector().push_back(
90 new OpCalculateHOJacForFace(jac_ptr));
91 volPipelineLhs->getOpPtrVector().push_back(
92 new OpInvertMatrix<2>(jac_ptr, det_ptr, inv_jac_ptr));
93 volPipelineLhs->getOpPtrVector().push_back(
95 volPipelineRhs->getOpPtrVector().push_back(new OpRhsU(domainField));
96
97 // dm = simpleInterface->getDM();
99 boost::shared_ptr<VolEle> null;
100
102
106}
PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set KSP right hand side evaluation function
Definition: DMMoFEM.cpp:625
PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
Set KSP operators and push mofem finite element methods.
Definition: DMMoFEM.cpp:666
OpCalculateHOJacForFaceImpl< 2 > OpCalculateHOJacForFace
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition: Simple.cpp:671
const std::string getDomainFEName() const
Get the Domain FE Name.
Definition: Simple.hpp:327
MoFEMErrorCode setIntegrationRule()

◆ readMesh()

MoFEMErrorCode Poisson::readMesh ( )
private

Definition at line 54 of file poisson_prob.cpp.

54 {
56
60
62}
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

◆ runAnalysis()

MoFEMErrorCode Poisson::runAnalysis ( )

Definition at line 194 of file poisson_prob.cpp.

194 {
196
199 CHKERR applyBc();
200
202
204
206
208}
MoFEMErrorCode applyBc()
MoFEMErrorCode assembleSystem()
MoFEMErrorCode setupSystem()
MoFEMErrorCode setOutput()
MoFEMErrorCode readMesh()
MoFEMErrorCode solveSystem()

◆ setIntegrationRule()

MoFEMErrorCode Poisson::setIntegrationRule ( )
private

Definition at line 75 of file poisson_prob.cpp.

75 {
77 auto vol_ruleLhs = [](int, int, int p) -> int { return 2 * (p - 1); };
78 auto vol_ruleRhs = [](int, int, int p) -> int { return 2 * (p - 1); };
79 volPipelineLhs->getRuleHook = vol_ruleLhs;
80 volPipelineRhs->getRuleHook = vol_ruleRhs;
82}
static Index< 'p', 3 > p

◆ setOutput()

MoFEMErrorCode Poisson::setOutput ( )
private

Definition at line 152 of file poisson_prob.cpp.

152 {
154
156 boost::shared_ptr<VolEle>(new PostProcFaceOnRefinedMesh(mField));
157
158 CHKERR
159 boost::static_pointer_cast<PostProcFaceOnRefinedMesh>(postFaceVolume)
160 ->generateReferenceElementMesh();
161 CHKERR
162 boost::static_pointer_cast<PostProcFaceOnRefinedMesh>(postFaceVolume)
163 ->addFieldValuesPostProc(domainField);
164 // CHKERR
165 // boost::static_pointer_cast<PostProcFaceOnRefinedMesh>(postFaceVolume)
166 // ->addFieldValuesGradientPostProc(domainField);
167
169}
boost::shared_ptr< VolEle > postFaceVolume
Postprocess on face.

◆ setupSystem()

MoFEMErrorCode Poisson::setupSystem ( )
private

Definition at line 64 of file poisson_prob.cpp.

64 {
66
71
73}
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
@ H1
continuous field
Definition: definitions.h:85
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:476
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition: Simple.cpp:614

◆ solveSystem()

MoFEMErrorCode Poisson::solveSystem ( )
private

Definition at line 108 of file poisson_prob.cpp.

108 {
110
111 Vec global_rhs, global_solution;
112 CHKERR DMCreateGlobalVector(dm, &global_rhs);
113 CHKERR VecDuplicate(global_rhs, &global_solution);
114
115
116 // ksp = createKSP(mField.get_comm());
117
118 CHKERR KSPCreate(PETSC_COMM_WORLD, &ksp);
119
120 // CHKERR KSPSetOperators(ksp, volPipelineLhs->getFEMethod()->ksp_B,
121 // volPipelineLhs->getFEMethod()->ksp_B);
122
123 CHKERR KSPSetFromOptions(ksp);
124
125 CHKERR KSPSetDM(ksp, dm);
126
127 CHKERR KSPSetUp(ksp);
128
129 CHKERR KSPSolve(ksp, global_rhs, global_solution);
130
131 CHKERR DMoFEMMeshToGlobalVector(dm, global_solution, INSERT_VALUES,
132 SCATTER_REVERSE);
133
134
135
136 CHKERR KSPDestroy(&ksp);
137 CHKERR VecDestroy(&global_solution);
138 CHKERR VecDestroy(&global_rhs);
139
142
143 CHKERR
144 boost::static_pointer_cast<PostProcFaceOnRefinedMesh>(postFaceVolume)
145 ->writeFile("out_vol.h5m");
146
147 CHKERR DMDestroy(&dm);
148
150}
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMoFEM.cpp:574
PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode)
set ghosted vector values on all existing mesh entities
Definition: DMMoFEM.cpp:523
const FTensor::Tensor2< T, Dim, Dim > Vec

Member Data Documentation

◆ dm

DM Poisson::dm
private

Definition at line 33 of file poisson_prob.cpp.

◆ domainField

std::string Poisson::domainField
private

Definition at line 43 of file poisson_prob.cpp.

◆ invJac

MatrixDouble Poisson::invJac
private

Definition at line 36 of file poisson_prob.cpp.

◆ ksp

KSP Poisson::ksp
private

Definition at line 34 of file poisson_prob.cpp.

◆ mField

MoFEM::Interface& Poisson::mField
private

Definition at line 29 of file poisson_prob.cpp.

◆ moab

moab::Interface& Poisson::moab
private

Definition at line 30 of file poisson_prob.cpp.

◆ mpiComm

MPI_Comm Poisson::mpiComm
private

Definition at line 38 of file poisson_prob.cpp.

◆ mpiRank

const int Poisson::mpiRank
private

Definition at line 39 of file poisson_prob.cpp.

◆ order

int Poisson::order
private

Definition at line 41 of file poisson_prob.cpp.

◆ postFaceVolume

boost::shared_ptr<VolEle> Poisson::postFaceVolume
private

Definition at line 27 of file poisson_prob.cpp.

◆ simpleInterface

Simple* Poisson::simpleInterface
private

Definition at line 31 of file poisson_prob.cpp.

◆ volPipelineLhs

boost::shared_ptr<VolEle> Poisson::volPipelineLhs
private

Definition at line 25 of file poisson_prob.cpp.

◆ volPipelineRhs

boost::shared_ptr<VolEle> Poisson::volPipelineRhs
private

Definition at line 26 of file poisson_prob.cpp.


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