v0.14.0
Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes | List of all members
Poisson2DNonhomogeneous Struct Reference
Collaboration diagram for Poisson2DNonhomogeneous:
[legend]

Public Member Functions

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

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...
 

Static Private Member Functions

static double sourceTermFunction (const double x, const double y, const double z)
 
static double boundaryFunction (const double x, const double y, const double z)
 

Private Attributes

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

Detailed Description

Definition at line 12 of file poisson_2d_nonhomogeneous.cpp.

Constructor & Destructor Documentation

◆ Poisson2DNonhomogeneous()

Poisson2DNonhomogeneous::Poisson2DNonhomogeneous ( MoFEM::Interface m_field)

Definition at line 57 of file poisson_2d_nonhomogeneous.cpp.

58  : domainField("U"), mField(m_field) {}

Member Function Documentation

◆ assembleSystem()

MoFEMErrorCode Poisson2DNonhomogeneous::assembleSystem ( )
private

[Boundary condition]

[Assemble system]

Definition at line 108 of file poisson_2d_nonhomogeneous.cpp.

108  {
110 
111  auto pipeline_mng = mField.getInterface<PipelineManager>();
112 
113  { // Push operators to the Pipeline that is responsible for calculating LHS of
114  // domain elements
116  pipeline_mng->getOpDomainLhsPipeline(), {H1});
117  pipeline_mng->getOpDomainLhsPipeline().push_back(
118  new OpSetBc(domainField, true, boundaryMarker));
119  pipeline_mng->getOpDomainLhsPipeline().push_back(
121  pipeline_mng->getOpDomainLhsPipeline().push_back(
122  new OpUnSetBc(domainField));
123  }
124 
125  { // Push operators to the Pipeline that is responsible for calculating RHS of
126  // domain elements
127  pipeline_mng->getOpDomainRhsPipeline().push_back(
128  new OpSetBc(domainField, true, boundaryMarker));
129  pipeline_mng->getOpDomainRhsPipeline().push_back(
131  pipeline_mng->getOpDomainRhsPipeline().push_back(
132  new OpUnSetBc(domainField));
133  }
134 
135  { // Push operators to the Pipeline that is responsible for calculating LHS of
136  // boundary elements
137  pipeline_mng->getOpBoundaryLhsPipeline().push_back(
138  new OpSetBc(domainField, false, boundaryMarker));
139  pipeline_mng->getOpBoundaryLhsPipeline().push_back(
141  pipeline_mng->getOpBoundaryLhsPipeline().push_back(
142  new OpUnSetBc(domainField));
143  }
144 
145  { // Push operators to the Pipeline that is responsible for calculating RHS of
146  // boundary elements
147  pipeline_mng->getOpBoundaryRhsPipeline().push_back(
148  new OpSetBc(domainField, false, boundaryMarker));
149  pipeline_mng->getOpBoundaryRhsPipeline().push_back(
151  pipeline_mng->getOpBoundaryRhsPipeline().push_back(
152  new OpUnSetBc(domainField));
153  }
154 
156 }

◆ boundaryCondition()

MoFEMErrorCode Poisson2DNonhomogeneous::boundaryCondition ( )
private

[Setup problem]

[Boundary condition]

Definition at line 90 of file poisson_2d_nonhomogeneous.cpp.

90  {
92 
93  auto bc_mng = mField.getInterface<BcManager>();
94  CHKERR bc_mng->pushMarkDOFsOnEntities(simpleInterface->getProblemName(),
95  "BOUNDARY_CONDITION", domainField, 0, 1,
96  true);
97  // merge markers from all blocksets "BOUNDARY_CONDITION"
98  boundaryMarker = bc_mng->getMergedBlocksMarker({"BOUNDARY_CONDITION"});
99  // get entities on blocksets "BOUNDARY_CONDITION"
101  bc_mng->getMergedBlocksRange({"BOUNDARY_CONDITION"});
102 
104 }

◆ boundaryFunction()

static double Poisson2DNonhomogeneous::boundaryFunction ( const double  x,
const double  y,
const double  z 
)
inlinestaticprivate

Definition at line 36 of file poisson_2d_nonhomogeneous.cpp.

37  {
38  return sin(x * 10.) * cos(y * 10.);
39  // return 0;
40  }

◆ outputResults()

MoFEMErrorCode Poisson2DNonhomogeneous::outputResults ( )
private

[Solve system]

[Output results]

Definition at line 231 of file poisson_2d_nonhomogeneous.cpp.

231  {
233 
234  auto pipeline_mng = mField.getInterface<PipelineManager>();
235  pipeline_mng->getDomainLhsFE().reset();
236  pipeline_mng->getBoundaryLhsFE().reset();
237  pipeline_mng->getDomainRhsFE().reset();
238  pipeline_mng->getBoundaryRhsFE().reset();
239 
240  auto d_ptr = boost::make_shared<VectorDouble>();
241  auto l_ptr = boost::make_shared<VectorDouble>();
242 
244 
245  auto post_proc_fe = boost::make_shared<PostProcFaceEle>(mField);
246 
247  post_proc_fe->getOpPtrVector().push_back(
249  post_proc_fe->getOpPtrVector().push_back(new OpPPMap(
250  post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
251  {{domainField, d_ptr}}, {}, {}, {}));
252  pipeline_mng->getDomainRhsFE() = post_proc_fe;
253 
254  CHKERR pipeline_mng->loopFiniteElements();
255  CHKERR post_proc_fe->writeFile("out_result.h5m");
256 
258 }

◆ readMesh()

MoFEMErrorCode Poisson2DNonhomogeneous::readMesh ( )
private

[Read mesh]

Definition at line 60 of file poisson_2d_nonhomogeneous.cpp.

◆ runProgram()

MoFEMErrorCode Poisson2DNonhomogeneous::runProgram ( )

Definition at line 260 of file poisson_2d_nonhomogeneous.cpp.

260  {
262 
263  readMesh();
264  setupProblem();
266  assembleSystem();
268  solveSystem();
269  outputResults();
270 
272 }

◆ setIntegrationRules()

MoFEMErrorCode Poisson2DNonhomogeneous::setIntegrationRules ( )
private

[Assemble system]

[Set integration rules]

Definition at line 160 of file poisson_2d_nonhomogeneous.cpp.

160  {
162 
163  auto pipeline_mng = mField.getInterface<PipelineManager>();
164 
165  auto domain_rule_lhs = [](int, int, int p) -> int { return 2 * (p - 1); };
166  auto domain_rule_rhs = [](int, int, int p) -> int { return 2 * (p - 1); };
167  CHKERR pipeline_mng->setDomainLhsIntegrationRule(domain_rule_lhs);
168  CHKERR pipeline_mng->setDomainRhsIntegrationRule(domain_rule_rhs);
169 
170  auto boundary_rule_lhs = [](int, int, int p) -> int { return 2 * p; };
171  auto boundary_rule_rhs = [](int, int, int p) -> int { return 2 * p; };
172  CHKERR pipeline_mng->setBoundaryLhsIntegrationRule(boundary_rule_lhs);
173  CHKERR pipeline_mng->setBoundaryLhsIntegrationRule(boundary_rule_rhs);
174 
176 }

◆ setupProblem()

MoFEMErrorCode Poisson2DNonhomogeneous::setupProblem ( )
private

[Read mesh]

[Setup problem]

Definition at line 72 of file poisson_2d_nonhomogeneous.cpp.

◆ solveSystem()

MoFEMErrorCode Poisson2DNonhomogeneous::solveSystem ( )
private

[Set integration rules]

[Solve system]

Definition at line 180 of file poisson_2d_nonhomogeneous.cpp.

180  {
182 
183  auto pipeline_mng = mField.getInterface<PipelineManager>();
184 
185  auto ksp_solver = pipeline_mng->createKSP();
186  CHKERR KSPSetFromOptions(ksp_solver);
187 
188  // Create RHS and solution vectors
189  auto dm = simpleInterface->getDM();
190  auto F = createDMVector(dm);
191  auto D = vectorDuplicate(F);
192 
193  // Setup fieldsplit (block) solver - optional: yes/no
194  if (1) {
195  PC pc;
196  CHKERR KSPGetPC(ksp_solver, &pc);
197  PetscBool is_pcfs = PETSC_FALSE;
198  PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &is_pcfs);
199 
200  // Set up FIELDSPLIT, only when user set -pc_type fieldsplit
201  // Identify the index for boundary entities, remaining will be for domain
202  // Then split the fields for boundary and domain for solving
203  if (is_pcfs == PETSC_TRUE) {
204  IS is_domain, is_boundary;
205  const MoFEM::Problem *problem_ptr;
206  CHKERR DMMoFEMGetProblemPtr(dm, &problem_ptr);
207  CHKERR mField.getInterface<ISManager>()->isCreateProblemFieldAndRank(
208  problem_ptr->getName(), ROW, domainField, 0, 1, &is_boundary,
210  // CHKERR ISView(is_boundary, PETSC_VIEWER_STDOUT_SELF);
211  CHKERR PCFieldSplitSetIS(pc, NULL, is_boundary);
212  CHKERR ISDestroy(&is_boundary);
213  }
214  }
215 
216  CHKERR KSPSetUp(ksp_solver);
217 
218  // Solve the system
219  CHKERR KSPSolve(ksp_solver, F, D);
220 
221  // Scatter result data on the mesh
222  CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
223  CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
224  CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
225 
227 }

◆ sourceTermFunction()

static double Poisson2DNonhomogeneous::sourceTermFunction ( const double  x,
const double  y,
const double  z 
)
inlinestaticprivate

Definition at line 30 of file poisson_2d_nonhomogeneous.cpp.

31  {
32  return 200 * sin(x * 10.) * cos(y * 10.);
33  // return 1;
34  }

Member Data Documentation

◆ boundaryEntitiesForFieldsplit

Range Poisson2DNonhomogeneous::boundaryEntitiesForFieldsplit
private

Definition at line 54 of file poisson_2d_nonhomogeneous.cpp.

◆ boundaryMarker

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

Definition at line 51 of file poisson_2d_nonhomogeneous.cpp.

◆ domainField

std::string Poisson2DNonhomogeneous::domainField
private

Definition at line 47 of file poisson_2d_nonhomogeneous.cpp.

◆ mField

MoFEM::Interface& Poisson2DNonhomogeneous::mField
private

Definition at line 43 of file poisson_2d_nonhomogeneous.cpp.

◆ oRder

int Poisson2DNonhomogeneous::oRder
private

Definition at line 48 of file poisson_2d_nonhomogeneous.cpp.

◆ simpleInterface

Simple* Poisson2DNonhomogeneous::simpleInterface
private

Definition at line 44 of file poisson_2d_nonhomogeneous.cpp.


The documentation for this struct was generated from the following file:
Poisson2DNonhomogeneousOperators::OpBoundaryRhs
[OpBoundaryLhs]
Definition: poisson_2d_nonhomogeneous.hpp:243
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
Poisson2DNonhomogeneousOperators::OpBoundaryLhs
[OpDomainRhs]
Definition: poisson_2d_nonhomogeneous.hpp:168
H1
@ H1
continuous field
Definition: definitions.h:85
Poisson2DNonhomogeneous::oRder
int oRder
Definition: poisson_2d_nonhomogeneous.cpp:48
Poisson2DNonhomogeneousOperators::OpDomainLhs
[OpDomainLhs]
Definition: poisson_2d_nonhomogeneous.hpp:22
Poisson2DNonhomogeneousOperators::OpDomainRhs
[OpDomainLhs]
Definition: poisson_2d_nonhomogeneous.hpp:103
Poisson2DNonhomogeneous::simpleInterface
Simple * simpleInterface
Definition: poisson_2d_nonhomogeneous.cpp:44
Poisson2DNonhomogeneous::assembleSystem
MoFEMErrorCode assembleSystem()
[Boundary condition]
Definition: poisson_2d_nonhomogeneous.cpp:108
Poisson2DNonhomogeneous::mField
MoFEM::Interface & mField
Definition: poisson_2d_nonhomogeneous.cpp:43
MoFEM::Simple::loadFile
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition: Simple.cpp:194
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
MoFEM::DMoFEMMeshToLocalVector
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:527
Poisson2DNonhomogeneous::outputResults
MoFEMErrorCode outputResults()
[Solve system]
Definition: poisson_2d_nonhomogeneous.cpp:231
ROW
@ ROW
Definition: definitions.h:123
MoFEM::Simple::getOptions
MoFEMErrorCode getOptions()
get options
Definition: Simple.cpp:180
MoFEM::Simple::getDM
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition: Simple.cpp:667
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::createDMVector
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1018
Poisson2DNonhomogeneous::readMesh
MoFEMErrorCode readMesh()
[Read mesh]
Definition: poisson_2d_nonhomogeneous.cpp:60
MoFEM::ISManager
Section manager is used to create indexes and sections.
Definition: ISManager.hpp:23
MoFEM::BcManager
Simple interface for fast problem set-up.
Definition: BcManager.hpp:25
Poisson2DNonhomogeneous::domainField
std::string domainField
Definition: poisson_2d_nonhomogeneous.cpp:47
Poisson2DNonhomogeneous::boundaryFunction
static double boundaryFunction(const double x, const double y, const double z)
Definition: poisson_2d_nonhomogeneous.cpp:36
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
Poisson2DNonhomogeneous::sourceTermFunction
static double sourceTermFunction(const double x, const double y, const double z)
Definition: poisson_2d_nonhomogeneous.cpp:30
MoFEM::OpCalculateScalarFieldValues
Get value at integration points for scalar field.
Definition: UserDataOperators.hpp:82
MoFEM::Simple::addDomainField
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
MoFEM::PipelineManager::createKSP
SmartPetscObj< KSP > createKSP(SmartPetscObj< DM > dm=nullptr)
Create KSP (linear) solver.
Definition: PipelineManager.cpp:49
Poisson2DNonhomogeneous::solveSystem
MoFEMErrorCode solveSystem()
[Set integration rules]
Definition: poisson_2d_nonhomogeneous.cpp:180
Poisson2DNonhomogeneous::setIntegrationRules
MoFEMErrorCode setIntegrationRules()
[Assemble system]
Definition: poisson_2d_nonhomogeneous.cpp:160
Poisson2DNonhomogeneous::boundaryEntitiesForFieldsplit
Range boundaryEntitiesForFieldsplit
Definition: poisson_2d_nonhomogeneous.cpp:54
AINSWORTH_BERNSTEIN_BEZIER_BASE
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition: definitions.h:64
MoFEM::Simple::setFieldOrder
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition: Simple.cpp:473
MoFEM::AddHOOps
Add operators pushing bases from local to physical configuration.
Definition: HODataOperators.hpp:503
MoFEM::OpUnSetBc
Definition: FormsIntegrators.hpp:49
MoFEM::PipelineManager::getDomainLhsFE
boost::shared_ptr< FEMethod > & getDomainLhsFE()
Definition: PipelineManager.hpp:401
MoFEM::vectorDuplicate
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
Definition: PetscSmartObj.hpp:217
Poisson2DNonhomogeneous::boundaryCondition
MoFEMErrorCode boundaryCondition()
[Setup problem]
Definition: poisson_2d_nonhomogeneous.cpp:90
Poisson2DNonhomogeneous::setupProblem
MoFEMErrorCode setupProblem()
[Read mesh]
Definition: poisson_2d_nonhomogeneous.cpp:72
MoFEM::OpSetBc
Set indices on entities on finite element.
Definition: FormsIntegrators.hpp:38
MoFEM::DMMoFEMGetProblemPtr
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
Definition: DMMoFEM.cpp:430
MoFEM::Simple::addBoundaryField
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:282
Poisson2DNonhomogeneous::boundaryMarker
boost::shared_ptr< std::vector< unsigned char > > boundaryMarker
Definition: poisson_2d_nonhomogeneous.cpp:51
ReactionDiffusionEquation::D
const double D
diffusivity
Definition: reaction_diffusion.cpp:20
MoFEM::Problem
keeps basic data about problem
Definition: ProblemsMultiIndices.hpp:54
MoFEM::Simple::setUp
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition: Simple.cpp:611
MoFEM::Problem::getName
auto getName() const
Definition: ProblemsMultiIndices.hpp:372
MoFEM::Simple::getProblemName
const std::string getProblemName() const
Get the Problem Name.
Definition: Simple.hpp:362
convert.int
int
Definition: convert.py:64
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
F
@ F
Definition: free_surface.cpp:394
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698