v0.15.0
Loading...
Searching...
No Matches
schur_elastic.cpp
Go to the documentation of this file.
1/**
2 * @file schur_elastic.cpp
3 * @example mofem/tutorials/vec-10/schur_elastic.cpp
4 *
5 * @copyright Anonymous authors (c) 2025 under the MIT license (see LICENSE for
6 * details)
7 */
8
9#include <MoFEM.hpp>
10
11using namespace MoFEM;
12
13constexpr int BASE_DIM = 1; //< Dimension of the base functions
14
15//! [Define dimension]
16constexpr int SPACE_DIM =
17 EXECUTABLE_DIMENSION; //< Space dimension of problem, mesh
18//! [Define dimension]
20 ? AssemblyType::BLOCK_SCHUR
21 : AssemblyType::PETSC; //< selected assembly type
22
23constexpr IntegrationType I =
24 IntegrationType::GAUSS; //< selected integration type
25
26//! [Define entities]
31using DomainEleOp = DomainEle::UserDataOperator;
32using BoundaryEleOp = BoundaryEle::UserDataOperator;
33//! [Define entities]
34
35struct DomainBCs {};
36struct BoundaryBCs {};
37
44
45template <int DIM> struct PostProcEleByDim;
46
47template <> struct PostProcEleByDim<2> {
51};
52
53template <> struct PostProcEleByDim<3> {
57};
58
62
63#include <ElasticExample.hpp>
64
73
74//! [Run problem]
86//! [Run problem]
87
88struct SetUpSchur {
89 static boost::shared_ptr<SetUpSchur>
92
93protected:
94 SetUpSchur() = default;
95};
96
97static char help[] = "...\n\n";
98
99int main(int argc, char *argv[]) {
100
101 // Initialisation of MoFEM/PETSc and MOAB data structures
102 const char param_file[] = "param_file.petsc";
103 MoFEM::Core::Initialize(&argc, &argv, param_file, help);
104
105 auto core_log = logging::core::get();
106 core_log->add_sink(
108
109 core_log->add_sink(
111 LogManager::setLog("FieldEvaluator");
112 MOFEM_LOG_TAG("FieldEvaluator", "field_eval");
113
114 try {
115
116 //! [Register MoFEM discrete manager in PETSc]
117 DMType dm_name = "DMMOFEM";
118 CHKERR DMRegister_MoFEM(dm_name);
119 DMType dm_name_mg = "DMMOFEM_MG";
121 //! [Register MoFEM discrete manager in PETSc
122
123 //! [Create MoAB]
124 moab::Core mb_instance; ///< mesh database
125 moab::Interface &moab = mb_instance; ///< mesh database interface
126 //! [Create MoAB]
127
128 //! [Create MoFEM]
129 MoFEM::Core core(moab); ///< finite element database
130 MoFEM::Interface &m_field = core; ///< finite element database interface
131 //! [Create MoFEM]
132
133 //! [Example]
134 ElasticSchurExample ex(m_field);
135 CHKERR ex.runProblem();
136 //! [Example]
137 }
139
141}
142
145
146 MOFEM_LOG_CHANNEL("TIMER");
147 MOFEM_LOG_TAG("TIMER", "timer");
148
149 if (A == AssemblyType::BLOCK_SCHUR || A == AssemblyType::SCHUR) {
150 auto schur_ptr = SetUpSchur::createSetUpSchur(mField);
151 CHKERR schur_ptr->setUp(solver);
152 }
153
154
155 BOOST_LOG_SCOPED_THREAD_ATTR("Timeline", attrs::timer());
156 MOFEM_LOG("TIMER", Sev::inform) << "KSPSetUp";
157 CHKERR KSPSetUp(solver);
158 MOFEM_LOG("TIMER", Sev::inform) << "KSPSetUp <= Done";
159
160 DM dm;
161 CHKERR KSPGetDM(solver, &dm);
162 auto D = createDMVector(dm);
163 auto F = vectorDuplicate(D);
164
165 MOFEM_LOG("TIMER", Sev::inform) << "KSPSolve";
166 CHKERR KSPSolve(solver, F, D);
167 MOFEM_LOG("TIMER", Sev::inform) << "KSPSolve <= Done";
168
169 CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
170 CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
171 CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
172
174};
175
176struct SetUpSchurImpl : public SetUpSchur {
177
179 if (S) {
182 "Is expected that schur matrix is not allocated. This is "
183 "possible only is PC is set up twice");
184 }
185 }
186 virtual ~SetUpSchurImpl() = default;
187
189
190private:
196
198
200
205};
206
209 auto pip = mField.getInterface<PipelineManager>();
210 PC pc;
211 CHKERR KSPGetPC(solver, &pc);
212 PetscBool is_pcfs = PETSC_FALSE;
213 PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &is_pcfs);
214 if (is_pcfs) {
215 if (S) {
218 "Is expected that schur matrix is not allocated. This is "
219 "possible only is PC is set up twice");
220 }
223
224 // Add data to DM storage
226 CHKERR MatSetDM(S, PETSC_NULLPTR);
227 CHKERR MatSetBlockSize(S, SPACE_DIM);
228 CHKERR MatSetOption(S, MAT_SYMMETRIC, PETSC_TRUE);
229
231 CHKERR setPC(pc);
232
233 if constexpr (A == AssemblyType::BLOCK_SCHUR) {
234 // Set DM to use shell block matrix
235 DM solver_dm;
236 CHKERR KSPGetDM(solver, &solver_dm);
237 CHKERR DMSetMatType(solver_dm, MATSHELL);
238 }
239 CHKERR KSPSetUp(solver);
241
242 } else {
243 pip->getOpBoundaryLhsPipeline().push_front(createOpSchurAssembleBegin());
244 pip->getOpBoundaryLhsPipeline().push_back(createOpSchurAssembleEnd({}, {}));
245 pip->getOpDomainLhsPipeline().push_front(createOpSchurAssembleBegin());
246 pip->getOpDomainLhsPipeline().push_back(createOpSchurAssembleEnd({}, {}));
247 }
249}
250
254 CHKERR mField.get_moab().get_entities_by_dimension(simple->getMeshset(),
256 CHKERR mField.get_moab().get_entities_by_handle(simple->getMeshset(),
257 subEnts);
258 subEnts = subtract(subEnts, volEnts);
260};
261
265
266 auto create_dm = [&](const char *name, auto &ents, auto dm_type) {
267 auto dm = createDM(mField.get_comm(), dm_type);
268 auto create_dm_imp = [&]() {
270 CHKERR DMMoFEMCreateSubDM(dm, simple->getDM(), name);
271 CHKERR DMMoFEMSetSquareProblem(dm, PETSC_TRUE);
272 CHKERR DMMoFEMAddElement(dm, simple->getDomainFEName());
273 auto sub_ents_ptr = boost::make_shared<Range>(ents);
274 CHKERR DMMoFEMAddSubFieldRow(dm, "U", sub_ents_ptr);
275 CHKERR DMMoFEMAddSubFieldCol(dm, "U", sub_ents_ptr);
276 CHKERR DMSetUp(dm);
278 };
280 create_dm_imp(),
281 "Error in creating schurDM. It is possible that schurDM is "
282 "already created");
283 return dm;
284 };
285
286 schurDM = create_dm("SCHUR", subEnts, "DMMOFEM_MG");
287 blockDM = create_dm("BLOCK", volEnts, "DMMOFEM");
288
289 if constexpr (A == AssemblyType::BLOCK_SCHUR) {
290
291 auto get_nested_mat_data = [&]() -> boost::shared_ptr<NestSchurData> {
292 auto block_mat_data =
294
295 {{
296
297 simple->getDomainFEName(),
298
299 {{"U", "U"}
300
301 }}}
302
303 );
304
306
307 {schurDM, blockDM}, block_mat_data,
308
309 {"U"}, {boost::make_shared<Range>(volEnts)}
310
311 );
312 };
313
314 auto nested_mat_data = get_nested_mat_data();
315
316 CHKERR DMMoFEMSetNestSchurData(simple->getDM(), nested_mat_data);
317 }
318
320}
321
325 auto pip = mField.getInterface<PipelineManager>();
326
327 // Boundary
328 auto dm_is = getDMSubData(schurDM)->getSmartRowIs();
329 auto ao_up = createAOMappingIS(dm_is, PETSC_NULLPTR);
330 pip->getOpBoundaryLhsPipeline().push_front(createOpSchurAssembleBegin());
331 pip->getOpBoundaryLhsPipeline().push_back(createOpSchurAssembleEnd(
332 {"U"}, {boost::make_shared<Range>(volEnts)}, ao_up, S, true, true));
333 // Domain
334 pip->getOpDomainLhsPipeline().push_front(createOpSchurAssembleBegin());
335 pip->getOpDomainLhsPipeline().push_back(createOpSchurAssembleEnd(
336 {"U"}, {boost::make_shared<Range>(volEnts)}, ao_up, S, true, true));
337
338 auto pre_proc_schur_lhs_ptr = boost::make_shared<FEMethod>();
339 auto post_proc_schur_lhs_ptr = boost::make_shared<FEMethod>();
340
341 pre_proc_schur_lhs_ptr->preProcessHook = [this]() {
343 if (S) {
344 CHKERR MatZeroEntries(S);
345 }
346 MOFEM_LOG("TIMER", Sev::inform) << "Lhs Assemble Begin";
348 };
349
350 post_proc_schur_lhs_ptr->postProcessHook = [this, post_proc_schur_lhs_ptr,
351 ao_up]() {
353 CHKERR MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY);
354 CHKERR MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY);
356 mField, post_proc_schur_lhs_ptr, 1, S, ao_up)();
357 MOFEM_LOG("TIMER", Sev::inform) << "Lhs Assemble End";
359 };
360
361 auto ksp_ctx_ptr = getDMKspCtx(simple->getDM());
362 ksp_ctx_ptr->getPreProcSetOperators().push_front(pre_proc_schur_lhs_ptr);
363 ksp_ctx_ptr->getPostProcSetOperators().push_back(post_proc_schur_lhs_ptr);
364
366}
367
371 SmartPetscObj<IS> vol_is;
372 mField.getInterface<ISManager>()->isCreateProblemFieldAndRank(
373 simple->getProblemName(), ROW, "U", 0, SPACE_DIM, vol_is, &volEnts);
374 CHKERR PCFieldSplitSetIS(pc, NULL, vol_is);
375 CHKERR PCFieldSplitSetSchurPre(pc, PC_FIELDSPLIT_SCHUR_PRE_USER, S);
377}
378
381
382 auto get_pc = [](auto ksp) {
383 PC pc_raw;
384 CHKERR KSPGetPC(ksp, &pc_raw);
385 return SmartPetscObj<PC>(pc_raw, true); // bump reference
386 };
387
388 if constexpr (A == AssemblyType::BLOCK_SCHUR) {
390 auto A = createDMBlockMat(simple->getDM());
391 auto P = createDMNestSchurMat(simple->getDM());
392 CHKERR PCSetOperators(pc, A, P);
393
394 KSP *subksp;
395 CHKERR PCFieldSplitSchurGetSubKSP(pc, PETSC_NULLPTR, &subksp);
396 CHKERR setSchurA00MatSolvePC(get_pc(subksp[0]));
397
398 auto set_pc_p_mg = [](auto dm, auto pc, auto S) {
400 CHKERR PCSetDM(pc, dm);
401 PetscBool same = PETSC_FALSE;
402 PetscObjectTypeCompare((PetscObject)pc, PCMG, &same);
403 if (same) {
404 MOFEM_LOG("TIMER", Sev::inform) << "Set up MG";
406 pc, createPCMGSetUpViaApproxOrdersCtx(dm, S, true));
407 CHKERR PCSetFromOptions(pc);
408 }
410 };
411
412 auto set_pc_ksp = [&](auto dm, auto pc, auto S) {
414 PetscBool same = PETSC_FALSE;
415 PetscObjectTypeCompare((PetscObject)pc, PCKSP, &same);
416 if (same) {
417 MOFEM_LOG("TIMER", Sev::inform) << "Set up inner KSP for PCKSP";
418 CHKERR PCSetFromOptions(pc);
419 KSP inner_ksp;
420 CHKERR PCKSPGetKSP(pc, &inner_ksp);
421 CHKERR KSPSetFromOptions(inner_ksp);
422 PC inner_pc;
423 CHKERR KSPGetPC(inner_ksp, &inner_pc);
424 CHKERR PCSetFromOptions(inner_pc);
425 CHKERR set_pc_p_mg(dm, inner_pc, S);
426 }
428 };
429
430 CHKERR set_pc_ksp(schurDM, get_pc(subksp[1]), S);
431 CHKERR set_pc_p_mg(schurDM, get_pc(subksp[1]), S);
432
433 CHKERR PetscFree(subksp);
434 }
436}
437
438boost::shared_ptr<SetUpSchur>
440 return boost::shared_ptr<SetUpSchur>(new SetUpSchurImpl(m_field));
441}
Implementation of elastic example class.
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
int main()
constexpr int SPACE_DIM
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
#define SCHUR_ASSEMBLE
Definition contact.cpp:18
@ ROW
#define CATCH_ERRORS
Catch errors.
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
PostProcEleByDim< SPACE_DIM >::PostProcEleDomain PostProcEleDomain
PostProcEleByDim< SPACE_DIM >::PostProcEleBdy PostProcEleBdy
@ F
PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[])
Must be called by user to set Sub DM MoFEM data structures.
Definition DMMoFEM.cpp:215
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition DMMoFEM.cpp:488
PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem)
set squared problem
Definition DMMoFEM.cpp:450
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
PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[])
Definition DMMoFEM.cpp:238
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
MoFEMErrorCode DMRegister_MGViaApproxOrders(const char sname[])
Register DM for Multi-Grid via approximation orders.
auto createDMVector(DM dm)
Get smart vector from DM.
Definition DMMoFEM.hpp:1234
PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[])
Definition DMMoFEM.cpp:280
auto createDMMatrix(DM dm)
Get smart matrix from DM.
Definition DMMoFEM.hpp:1191
IntegrationType
Form integrator integration types.
AssemblyType
[Storage and set boundary conditions]
@ PETSC
Standard PETSc assembly.
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
double D
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
OpSchurAssembleBase * createOpSchurAssembleEnd(std::vector< std::string > fields_name, std::vector< boost::shared_ptr< Range > > field_ents, SmartPetscObj< AO > ao, SmartPetscObj< Mat > schur, bool sym_schur, bool symm_op)
Construct a new Op Schur Assemble End object.
Definition Schur.cpp:2585
MoFEMErrorCode setSchurA00MatSolvePC(SmartPetscObj< PC > pc)
Set PC for A00 block.
Definition Schur.cpp:2627
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
boost::shared_ptr< PCMGSetUpViaApproxOrdersCtx > createPCMGSetUpViaApproxOrdersCtx(DM dm, Mat A, bool use_shell_mat)
createPCMGSetUpViaApproxOrdersCtx
auto getDMKspCtx(DM dm)
Get KSP context data structure used by DM.
Definition DMMoFEM.hpp:1248
auto getDMSubData(DM dm)
Get sub problem data structure.
Definition DMMoFEM.hpp:1292
MoFEMErrorCode PCMGSetUpViaApproxOrders(PC pc, boost::shared_ptr< PCMGSetUpViaApproxOrdersCtx > ctx, int verb)
Function build MG structure.
boost::shared_ptr< BlockStructure > createBlockMatStructure(DM dm, SchurFEOpsFEandFields schur_fe_op_vec)
Create a Mat Diag Blocks object.
Definition Schur.cpp:1082
boost::shared_ptr< NestSchurData > createSchurNestedMatrixStruture(std::pair< SmartPetscObj< DM >, SmartPetscObj< DM > > dms, boost::shared_ptr< BlockStructure > block_mat_data_ptr, std::vector< std::string > fields_names, std::vector< boost::shared_ptr< Range > > field_ents, bool add_preconditioner_block)
Get the Schur Nest Mat Array object.
Definition Schur.cpp:2343
auto createAOMappingIS(IS isapp, IS ispetsc)
Creates an application mapping using two index sets.
MoFEMErrorCode DMMoFEMSetNestSchurData(DM dm, boost::shared_ptr< NestSchurData >)
Definition DMMoFEM.cpp:1554
auto createDMNestSchurMat(DM dm)
Definition DMMoFEM.hpp:1218
auto createDM(MPI_Comm comm, const std::string dm_type_name)
Creates smart DM object.
OpSchurAssembleBase * createOpSchurAssembleBegin()
Definition Schur.cpp:2580
auto createDMBlockMat(DM dm)
Definition DMMoFEM.hpp:1211
constexpr AssemblyType A
#define EXECUTABLE_DIMENSION
Definition plastic.cpp:13
ElementsAndOps< SPACE_DIM >::SideEle SideEle
Definition plastic.cpp:61
static char help[]
constexpr int SPACE_DIM
[Define dimension]
constexpr int BASE_DIM
constexpr IntegrationType I
constexpr AssemblyType A
[Define dimension]
Boundary conditions marker.
Definition elastic.cpp:39
[Define entities]
Definition elastic.cpp:38
MoFEMErrorCode solveSystem()
[Solve]
MoFEMErrorCode outputResults()
[Solve]
MoFEMErrorCode assembleSystem()
[Boundary condition]
MoFEMErrorCode setupProblem()
[Read mesh]
MoFEM::Interface & mField
MoFEMErrorCode boundaryCondition()
[Set up problem]
MoFEMErrorCode checkResults()
[Postprocess results]
MoFEMErrorCode readMesh()
[Run problem]
MoFEMErrorCode kspSetUpAndSolve(SmartPetscObj< KSP > solver)
[Push operators to pipeline]
MoFEMErrorCode runProblem()
[Run problem]
ElasticSchurExample(MoFEM::Interface &field)
virtual moab::Interface & get_moab()=0
virtual MPI_Comm & get_comm() const =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:118
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
Class (Function) to enforce essential constrains on the left hand side diagonal.
Definition Essential.hpp:33
Section manager is used to create indexes and sections.
Definition ISManager.hpp:23
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
Assembly methods.
Definition Natural.hpp:65
Template struct for dimension-specific finite element types.
PipelineManager interface.
Simple interface for fast problem set-up.
Definition Simple.hpp:27
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
PipelineManager::ElementsAndOpsByDim< 2 >::FaceSideEle SideEle
PipelineManager::ElementsAndOpsByDim< 3 >::FaceSideEle SideEle
MoFEMErrorCode setUp(SmartPetscObj< KSP > solver)
virtual ~SetUpSchurImpl()=default
SmartPetscObj< Mat > S
MoFEMErrorCode createSubDM()
MoFEMErrorCode setDiagonalPC(PC pc)
SetUpSchurImpl(MoFEM::Interface &m_field)
SmartPetscObj< DM > schurDM
Definition contact.cpp:1025
MoFEMErrorCode setEntities()
MoFEMErrorCode setUp(SmartPetscObj< KSP >)
SmartPetscObj< DM > blockDM
Definition contact.cpp:1026
MoFEMErrorCode setPC(PC pc)
MoFEMErrorCode setOperator()
MoFEM::Interface & mField
[Push operators to pipeline]
SetUpSchur()=default
virtual MoFEMErrorCode setUp(SmartPetscObj< KSP > solver)=0
static boost::shared_ptr< SetUpSchur > createSetUpSchur(MoFEM::Interface &m_field)
static boost::shared_ptr< SetUpSchur > createSetUpSchur(MoFEM::Interface &m_field)