v0.14.0
Loading...
Searching...
No Matches
poisson_2d_lagrange_multiplier.cpp
Go to the documentation of this file.
1#include <stdlib.h>
2#include <MoFEM.hpp>
3
4using namespace MoFEM;
6
8
13
14static char help[] = "...\n\n";
15
17public:
19
20 // Declaration of the main function to run analysis
22
23private:
24 // Declaration of other main functions called in runProgram()
32
33 // Function to calculate the Source term
34 static double sourceTermFunction(const double x, const double y,
35 const double z) {
36 return 200 * sin(x * 10.) * cos(y * 10.);
37 // return 1;
38 }
39 // Function to calculate the Boundary term
40 static double boundaryFunction(const double x, const double y,
41 const double z) {
42 return sin(x * 10.) * cos(y * 10.);
43 // return 0;
44 }
45
46 // Main interfaces
49
50 // Field name and approximation order
51 std::string domainField; // displacement field
52 std::string boundaryField; // Lagrange multiplier field
53 int oRder;
54
55 // Object to mark boundary entities for the assembling of domain elements
56 boost::shared_ptr<std::vector<unsigned char>> boundaryMarker;
57
58 // Boundary entities marked for fieldsplit (block) solver - optional
60};
61
63 MoFEM::Interface &m_field)
64 : domainField("U"), boundaryField("L"), mField(m_field) {}
65
75
93
96
97 // Get boundary edges marked in block named "BOUNDARY_CONDITION"
98 auto get_ents_on_mesh = [&]() {
99 Range boundary_entities;
101 std::string entity_name = it->getName();
102 if (entity_name.compare(0, 18, "BOUNDARY_CONDITION") == 0) {
103 CHKERR it->getMeshsetIdEntitiesByDimension(mField.get_moab(), 1,
104 boundary_entities, true);
105 }
106 }
107 // Add vertices to boundary entities
108 Range boundary_vertices;
109 CHKERR mField.get_moab().get_connectivity(boundary_entities,
110 boundary_vertices, true);
111 boundary_entities.merge(boundary_vertices);
112
113 // Store entities for fieldsplit (block) solver
114 boundaryEntitiesForFieldsplit = boundary_entities;
115
116 return boundary_entities;
117 };
118
119 auto mark_boundary_dofs = [&](Range &&skin_edges) {
120 auto problem_manager = mField.getInterface<ProblemsManager>();
121 auto marker_ptr = boost::make_shared<std::vector<unsigned char>>();
122 problem_manager->markDofs(simpleInterface->getProblemName(), ROW,
123 skin_edges, *marker_ptr);
124 return marker_ptr;
125 };
126
127 // Get global local vector of marked DOFs. Is global, since is set for all
128 // DOFs on processor. Is local since only DOFs on processor are in the
129 // vector. To access DOFs use local indices.
130 boundaryMarker = mark_boundary_dofs(get_ents_on_mesh());
131
133}
134
137
138 auto pipeline_mng = mField.getInterface<PipelineManager>();
139 auto det_ptr = boost::make_shared<VectorDouble>();
140 auto jac_ptr = boost::make_shared<MatrixDouble>();
141 auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
142
143 { // Push operators to the Pipeline that is responsible for calculating LHS of
144 // domain elements
145 pipeline_mng->getOpDomainLhsPipeline().push_back(
146 new OpCalculateHOJac<2>(jac_ptr));
147 pipeline_mng->getOpDomainLhsPipeline().push_back(
148 new OpInvertMatrix<2>(jac_ptr, det_ptr, inv_jac_ptr));
149 pipeline_mng->getOpDomainLhsPipeline().push_back(
150 new OpSetHOInvJacToScalarBases<2>(H1, inv_jac_ptr));
151 pipeline_mng->getOpDomainLhsPipeline().push_back(
153 pipeline_mng->getOpDomainLhsPipeline().push_back(
155 }
156
157 { // Push operators to the Pipeline that is responsible for calculating RHS of
158 // domain elements
159 pipeline_mng->getOpDomainRhsPipeline().push_back(
161 }
162
163 { // Push operators to the Pipeline that is responsible for calculating LHS of
164 // boundary elements (Lagrange multiplier)
165 pipeline_mng->getOpBoundaryLhsPipeline().push_back(
167 pipeline_mng->getOpBoundaryLhsPipeline().push_back(
169 }
170
171 { // Push operators to the Pipeline that is responsible for calculating RHS of
172 // boundary elements (Lagrange multiplier)
173 pipeline_mng->getOpBoundaryRhsPipeline().push_back(
175 }
176
178}
179
182
183 auto pipeline_mng = mField.getInterface<PipelineManager>();
184
185 auto domain_rule_lhs = [](int, int, int p) -> int { return 2 * (p - 1); };
186 auto domain_rule_rhs = [](int, int, int p) -> int { return 2 * (p - 1); };
187 CHKERR pipeline_mng->setDomainLhsIntegrationRule(domain_rule_lhs);
188 CHKERR pipeline_mng->setDomainRhsIntegrationRule(domain_rule_rhs);
189
190 auto boundary_rule_lhs = [](int, int, int p) -> int { return 2 * p; };
191 auto boundary_rule_rhs = [](int, int, int p) -> int { return 2 * p; };
192 CHKERR pipeline_mng->setBoundaryLhsIntegrationRule(boundary_rule_lhs);
193 CHKERR pipeline_mng->setBoundaryLhsIntegrationRule(boundary_rule_rhs);
194
196}
197
200
201 auto pipeline_mng = mField.getInterface<PipelineManager>();
202
203 auto ksp_solver = pipeline_mng->createKSP();
204 CHKERR KSPSetFromOptions(ksp_solver);
205
206 // Create RHS and solution vectors
207 auto dm = simpleInterface->getDM();
208 auto F = createDMVector(dm);
209 auto D = vectorDuplicate(F);
210
211 // Setup fieldsplit (block) solver - optional: yes/no
212 if (1) {
213 PC pc;
214 CHKERR KSPGetPC(ksp_solver, &pc);
215 PetscBool is_pcfs = PETSC_FALSE;
216 PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &is_pcfs);
217
218 // Set up FIELDSPLIT, only when user set -pc_type fieldsplit
219 // Identify the index for boundary entities, remaining will be for domain
220 // Then split the fields for boundary and domain for solving
221 if (is_pcfs == PETSC_TRUE) {
222 IS is_domain, is_boundary;
223 cerr << "Running FIELDSPLIT..." << endl;
224 const MoFEM::Problem *problem_ptr;
225 CHKERR DMMoFEMGetProblemPtr(dm, &problem_ptr);
226 CHKERR mField.getInterface<ISManager>()->isCreateProblemFieldAndRank(
227 problem_ptr->getName(), ROW, domainField, 0, 1, &is_boundary,
229 // CHKERR ISView(is_boundary, PETSC_VIEWER_STDOUT_SELF);
230
231 CHKERR PCFieldSplitSetIS(pc, NULL, is_boundary);
232
233 CHKERR ISDestroy(&is_boundary);
234 }
235 }
236
237 CHKERR KSPSetUp(ksp_solver);
238
239 // Solve the system
240 CHKERR KSPSolve(ksp_solver, F, D);
241
242 // Scatter result data on the mesh
243 CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
244 CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
245 CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
246
248}
249
252
253 auto pipeline_mng = mField.getInterface<PipelineManager>();
254 pipeline_mng->getDomainLhsFE().reset();
255 pipeline_mng->getBoundaryLhsFE().reset();
256
257 auto d_ptr = boost::make_shared<VectorDouble>();
258 auto l_ptr = boost::make_shared<VectorDouble>();
259
261
262 auto post_proc_domain_fe = boost::make_shared<PostProcFaceEle>(mField);
263 post_proc_domain_fe->getOpPtrVector().push_back(
265 post_proc_domain_fe->getOpPtrVector().push_back(
266 new OpPPMap(post_proc_domain_fe->getPostProcMesh(),
267 post_proc_domain_fe->getMapGaussPts(), {{domainField, d_ptr}},
268 {}, {}, {}));
269 pipeline_mng->getDomainRhsFE() = post_proc_domain_fe;
270
271 auto post_proc_boundary_fe = boost::make_shared<PostProcEdgeEle>(mField);
272 post_proc_boundary_fe->getOpPtrVector().push_back(
273 new OpCalculateScalarFieldValues(boundaryField, l_ptr));
274 post_proc_boundary_fe->getOpPtrVector().push_back(
275 new OpPPMap(post_proc_boundary_fe->getPostProcMesh(),
276 post_proc_boundary_fe->getMapGaussPts(),
277 {{boundaryField, l_ptr}}, {}, {}, {}));
278 pipeline_mng->getBoundaryRhsFE() = post_proc_boundary_fe;
279
280 CHKERR pipeline_mng->loopFiniteElements();
281 CHKERR post_proc_domain_fe->writeFile("out_result_domain.h5m");
282 CHKERR post_proc_boundary_fe->writeFile("out_result_boundary.h5m");
283
285}
286
300
301int main(int argc, char *argv[]) {
302
303 // Initialisation of MoFEM/PETSc and MOAB data structures
304 const char param_file[] = "param_file.petsc";
305 MoFEM::Core::Initialize(&argc, &argv, param_file, help);
306
307 // Error handling
308 try {
309 // Register MoFEM discrete manager in PETSc
310 DMType dm_name = "DMMOFEM";
311 CHKERR DMRegister_MoFEM(dm_name);
312
313 // Create MOAB instance
314 moab::Core mb_instance; // mesh database
315 moab::Interface &moab = mb_instance; // mesh database interface
316
317 // Create MoFEM instance
318 MoFEM::Core core(moab); // finite element database
319 MoFEM::Interface &m_field = core; // finite element interface
320
321 // Run the main analysis
322 Poisson2DLagrangeMultiplier poisson_problem(m_field);
323 CHKERR poisson_problem.runProgram();
324 }
326
327 // Finish work: cleaning memory, getting statistics, etc.
329
330 return 0;
331}
int main()
@ ROW
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition definitions.h:64
@ H1
continuous field
Definition definitions.h:85
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ BLOCKSET
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
constexpr auto domainField
@ 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:523
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
Definition DMMoFEM.cpp:426
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
auto createDMVector(DM dm)
Get smart vector from DM.
Definition DMMoFEM.hpp:1099
SmartPetscObj< KSP > createKSP(SmartPetscObj< DM > dm=nullptr)
Create KSP (linear) solver.
#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.
double D
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
static char help[]
virtual moab::Interface & get_moab()=0
Core (interface) class.
Definition Core.hpp:82
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition Core.cpp:72
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:112
Deprecated interface functions.
Section manager is used to create indexes and sections.
Definition ISManager.hpp:23
Get value at integration points for scalar field.
Post post-proc data at points from hash maps.
Set inverse jacobian to base functions.
Modify integration weights on face to take in account higher-order geometry.
PipelineManager interface.
boost::shared_ptr< FEMethod > & getDomainLhsFE()
MoFEMErrorCode setDomainLhsIntegrationRule(RuleHookFun rule)
keeps basic data about problem
Problem manager is used to build and partition problems.
MoFEMErrorCode markDofs(const std::string problem_name, RowColData rc, const enum MarkOP op, const Range ents, std::vector< unsigned char > &marker) const
Create vector with marked indices.
Simple interface for fast problem set-up.
Definition Simple.hpp:27
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 loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition Simple.cpp:194
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:358
MoFEMErrorCode getOptions()
get options
Definition Simple.cpp:180
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition Simple.cpp:829
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition Simple.cpp:599
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition Simple.cpp:765
const std::string getProblemName() const
Get the Problem Name.
Definition Simple.hpp:408
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
static double sourceTermFunction(const double x, const double y, const double z)
Poisson2DLagrangeMultiplier(MoFEM::Interface &m_field)
boost::shared_ptr< std::vector< unsigned char > > boundaryMarker
static double boundaryFunction(const double x, const double y, const double z)