v0.15.0
Loading...
Searching...
No Matches
mofem/tutorials/vec-10/schur_elastic.cpp
/**
* @file schur_elastic.cpp
* @example mofem/tutorials/vec-10/schur_elastic.cpp
*
* @copyright Anonymous authors (c) 2025 under the MIT license (see LICENSE for
* details)
*/
#include <MoFEM.hpp>
using namespace MoFEM;
constexpr int BASE_DIM = 1; //< Dimension of the base functions
//! [Define dimension]
constexpr int SPACE_DIM =
EXECUTABLE_DIMENSION; //< Space dimension of problem, mesh
//! [Define dimension]
? AssemblyType::BLOCK_SCHUR
: AssemblyType::PETSC; //< selected assembly type
constexpr IntegrationType I =
IntegrationType::GAUSS; //< selected integration type
//! [Define entities]
using BoundaryEle =
using DomainEleOp = DomainEle::UserDataOperator;
using BoundaryEleOp = BoundaryEle::UserDataOperator;
//! [Define entities]
struct DomainBCs {};
struct BoundaryBCs {};
using OpDomainRhsBCs = DomainRhsBCs::OpFlux<DomainBCs, 1, SPACE_DIM>;
using OpBoundaryRhsBCs = BoundaryRhsBCs::OpFlux<BoundaryBCs, 1, SPACE_DIM>;
using OpBoundaryLhsBCs = BoundaryLhsBCs::OpFlux<BoundaryBCs, 1, SPACE_DIM>;
template <int DIM> struct PostProcEleByDim;
template <> struct PostProcEleByDim<2> {
};
template <> struct PostProcEleByDim<3> {
};
protected:
};
//! [Run problem]
}
//! [Run problem]
struct SetUpSchur {
static boost::shared_ptr<SetUpSchur>
protected:
SetUpSchur() = default;
};
static char help[] = "...\n\n";
int main(int argc, char *argv[]) {
// Initialisation of MoFEM/PETSc and MOAB data structures
const char param_file[] = "param_file.petsc";
MoFEM::Core::Initialize(&argc, &argv, param_file, help);
auto core_log = logging::core::get();
core_log->add_sink(
LogManager::createSink(LogManager::getStrmWorld(), "TIMER"));
core_log->add_sink(
LogManager::createSink(LogManager::getStrmSync(), "FieldEvaluator"));
LogManager::setLog("FieldEvaluator");
MOFEM_LOG_TAG("FieldEvaluator", "field_eval");
try {
//! [Register MoFEM discrete manager in PETSc]
DMType dm_name = "DMMOFEM";
DMType dm_name_mg = "DMMOFEM_MG";
//! [Register MoFEM discrete manager in PETSc
//! [Create MoAB]
moab::Core mb_instance; ///< mesh database
moab::Interface &moab = mb_instance; ///< mesh database interface
//! [Create MoAB]
//! [Create MoFEM]
MoFEM::Core core(moab); ///< finite element database
MoFEM::Interface &m_field = core; ///< finite element database interface
//! [Create MoFEM]
//! [Example]
ElasticSchurExample ex(m_field);
CHKERR ex.runProblem();
//! [Example]
}
}
MOFEM_LOG_TAG("TIMER", "timer");
if (A == AssemblyType::BLOCK_SCHUR || A == AssemblyType::SCHUR) {
CHKERR schur_ptr->setUp(solver);
}
BOOST_LOG_SCOPED_THREAD_ATTR("Timeline", attrs::timer());
MOFEM_LOG("TIMER", Sev::inform) << "KSPSetUp";
CHKERR KSPSetUp(solver);
MOFEM_LOG("TIMER", Sev::inform) << "KSPSetUp <= Done";
DM dm;
CHKERR KSPGetDM(solver, &dm);
auto D = createDMVector(dm);
auto F = vectorDuplicate(D);
MOFEM_LOG("TIMER", Sev::inform) << "KSPSolve";
CHKERR KSPSolve(solver, F, D);
MOFEM_LOG("TIMER", Sev::inform) << "KSPSolve <= Done";
CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
};
struct SetUpSchurImpl : public SetUpSchur {
if (S) {
"Is expected that schur matrix is not allocated. This is "
"possible only is PC is set up twice");
}
}
virtual ~SetUpSchurImpl() = default;
private:
};
PC pc;
CHKERR KSPGetPC(solver, &pc);
PetscBool is_pcfs = PETSC_FALSE;
PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &is_pcfs);
if (is_pcfs) {
if (S) {
"Is expected that schur matrix is not allocated. This is "
"possible only is PC is set up twice");
}
// Add data to DM storage
CHKERR MatSetDM(S, PETSC_NULLPTR);
CHKERR MatSetBlockSize(S, SPACE_DIM);
CHKERR MatSetOption(S, MAT_SYMMETRIC, PETSC_TRUE);
if constexpr (A == AssemblyType::BLOCK_SCHUR) {
// Set DM to use shell block matrix
DM solver_dm;
CHKERR KSPGetDM(solver, &solver_dm);
CHKERR DMSetMatType(solver_dm, MATSHELL);
}
CHKERR KSPSetUp(solver);
} else {
pip->getOpBoundaryLhsPipeline().push_front(createOpSchurAssembleBegin());
pip->getOpBoundaryLhsPipeline().push_back(createOpSchurAssembleEnd({}, {}));
pip->getOpDomainLhsPipeline().push_front(createOpSchurAssembleBegin());
pip->getOpDomainLhsPipeline().push_back(createOpSchurAssembleEnd({}, {}));
}
}
CHKERR mField.get_moab().get_entities_by_dimension(simple->getMeshset(),
CHKERR mField.get_moab().get_entities_by_handle(simple->getMeshset(),
subEnts = subtract(subEnts, volEnts);
};
auto create_dm = [&](const char *name, auto &ents, auto dm_type) {
auto dm = createDM(mField.get_comm(), dm_type);
auto create_dm_imp = [&]() {
CHKERR DMMoFEMCreateSubDM(dm, simple->getDM(), name);
CHKERR DMMoFEMAddElement(dm, simple->getDomainFEName());
auto sub_ents_ptr = boost::make_shared<Range>(ents);
CHKERR DMMoFEMAddSubFieldRow(dm, "U", sub_ents_ptr);
CHKERR DMMoFEMAddSubFieldCol(dm, "U", sub_ents_ptr);
CHKERR DMSetUp(dm);
};
create_dm_imp(),
"Error in creating schurDM. It is possible that schurDM is "
"already created");
return dm;
};
schurDM = create_dm("SCHUR", subEnts, "DMMOFEM_MG");
blockDM = create_dm("BLOCK", volEnts, "DMMOFEM");
if constexpr (A == AssemblyType::BLOCK_SCHUR) {
auto get_nested_mat_data = [&]() -> boost::shared_ptr<NestSchurData> {
auto block_mat_data =
{{
simple->getDomainFEName(),
{{"U", "U"}
}}}
);
{schurDM, blockDM}, block_mat_data,
{"U"}, {boost::make_shared<Range>(volEnts)}
);
};
auto nested_mat_data = get_nested_mat_data();
CHKERR DMMoFEMSetNestSchurData(simple->getDM(), nested_mat_data);
}
}
// Boundary
auto dm_is = getDMSubData(schurDM)->getSmartRowIs();
auto ao_up = createAOMappingIS(dm_is, PETSC_NULLPTR);
pip->getOpBoundaryLhsPipeline().push_front(createOpSchurAssembleBegin());
pip->getOpBoundaryLhsPipeline().push_back(createOpSchurAssembleEnd(
{"U"}, {boost::make_shared<Range>(volEnts)}, ao_up, S, true, true));
// Domain
pip->getOpDomainLhsPipeline().push_front(createOpSchurAssembleBegin());
pip->getOpDomainLhsPipeline().push_back(createOpSchurAssembleEnd(
{"U"}, {boost::make_shared<Range>(volEnts)}, ao_up, S, true, true));
auto pre_proc_schur_lhs_ptr = boost::make_shared<FEMethod>();
auto post_proc_schur_lhs_ptr = boost::make_shared<FEMethod>();
pre_proc_schur_lhs_ptr->preProcessHook = [this]() {
if (S) {
CHKERR MatZeroEntries(S);
}
MOFEM_LOG("TIMER", Sev::inform) << "Lhs Assemble Begin";
};
post_proc_schur_lhs_ptr->postProcessHook = [this, post_proc_schur_lhs_ptr,
ao_up]() {
CHKERR MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY);
CHKERR MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY);
mField, post_proc_schur_lhs_ptr, 1, S, ao_up)();
MOFEM_LOG("TIMER", Sev::inform) << "Lhs Assemble End";
};
auto ksp_ctx_ptr = getDMKspCtx(simple->getDM());
ksp_ctx_ptr->getPreProcSetOperators().push_front(pre_proc_schur_lhs_ptr);
ksp_ctx_ptr->getPostProcSetOperators().push_back(post_proc_schur_lhs_ptr);
}
mField.getInterface<ISManager>()->isCreateProblemFieldAndRank(
simple->getProblemName(), ROW, "U", 0, SPACE_DIM, vol_is, &volEnts);
CHKERR PCFieldSplitSetIS(pc, NULL, vol_is);
CHKERR PCFieldSplitSetSchurPre(pc, PC_FIELDSPLIT_SCHUR_PRE_USER, S);
}
auto get_pc = [](auto ksp) {
PC pc_raw;
CHKERR KSPGetPC(ksp, &pc_raw);
return SmartPetscObj<PC>(pc_raw, true); // bump reference
};
if constexpr (A == AssemblyType::BLOCK_SCHUR) {
auto A = createDMBlockMat(simple->getDM());
auto P = createDMNestSchurMat(simple->getDM());
CHKERR PCSetOperators(pc, A, P);
KSP *subksp;
CHKERR PCFieldSplitSchurGetSubKSP(pc, PETSC_NULLPTR, &subksp);
CHKERR setSchurA00MatSolvePC(get_pc(subksp[0]));
auto set_pc_p_mg = [](auto dm, auto pc, auto S) {
CHKERR PCSetDM(pc, dm);
PetscBool same = PETSC_FALSE;
PetscObjectTypeCompare((PetscObject)pc, PCMG, &same);
if (same) {
MOFEM_LOG("TIMER", Sev::inform) << "Set up MG";
CHKERR PCSetFromOptions(pc);
}
};
auto set_pc_ksp = [&](auto dm, auto pc, auto S) {
PetscBool same = PETSC_FALSE;
PetscObjectTypeCompare((PetscObject)pc, PCKSP, &same);
if (same) {
MOFEM_LOG("TIMER", Sev::inform) << "Set up inner KSP for PCKSP";
CHKERR PCSetFromOptions(pc);
KSP inner_ksp;
CHKERR PCKSPGetKSP(pc, &inner_ksp);
CHKERR KSPSetFromOptions(inner_ksp);
PC inner_pc;
CHKERR KSPGetPC(inner_ksp, &inner_pc);
CHKERR PCSetFromOptions(inner_pc);
CHKERR set_pc_p_mg(dm, inner_pc, S);
}
};
CHKERR set_pc_ksp(schurDM, get_pc(subksp[1]), S);
CHKERR set_pc_p_mg(schurDM, get_pc(subksp[1]), S);
CHKERR PetscFree(subksp);
}
}
boost::shared_ptr<SetUpSchur>
return boost::shared_ptr<SetUpSchur>(new SetUpSchurImpl(m_field));
}
Implementation of elastic example class.
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
static char help[]
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.
constexpr int BASE_DIM
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]
#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 IntegrationType I
constexpr AssemblyType A
NaturalBC< BoundaryEleOp >::Assembly< AT >::BiLinearForm< IT > BoundaryLhsBCs
Definition plastic.cpp:175
NaturalBC< BoundaryEleOp >::Assembly< AT >::LinearForm< IT > BoundaryRhsBCs
Definition plastic.cpp:172
#define EXECUTABLE_DIMENSION
Definition plastic.cpp:13
BoundaryLhsBCs::OpFlux< PlasticOps::BoundaryBCs, 1, SPACE_DIM > OpBoundaryLhsBCs
Definition plastic.cpp:177
NaturalBC< DomainEleOp >::Assembly< AT >::LinearForm< IT > DomainRhsBCs
Definition plastic.cpp:169
DomainRhsBCs::OpFlux< PlasticOps::DomainBCs, 1, SPACE_DIM > OpDomainRhsBCs
Definition plastic.cpp:171
BoundaryRhsBCs::OpFlux< PlasticOps::BoundaryBCs, 1, SPACE_DIM > OpBoundaryRhsBCs
Definition plastic.cpp:174
ElementsAndOps< SPACE_DIM >::SideEle SideEle
Definition plastic.cpp:61
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]
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
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.
virtual ~SetUpSchurImpl()=default
SmartPetscObj< Mat > S
MoFEMErrorCode createSubDM()
Definition contact.cpp:1104
MoFEMErrorCode setDiagonalPC(PC pc)
Definition contact.cpp:1254
SmartPetscObj< DM > schurDM
Definition contact.cpp:1025
MoFEMErrorCode setEntities()
MoFEMErrorCode setUp(SmartPetscObj< KSP >)
SmartPetscObj< DM > blockDM
Definition contact.cpp:1026
MoFEMErrorCode setPC(PC pc)
Definition contact.cpp:1246
MoFEMErrorCode setOperator()
Definition contact.cpp:1171
MoFEM::Interface & mField
[Push operators to pipeline]
SetUpSchur()=default
virtual MoFEMErrorCode setUp(SmartPetscObj< KSP >)=0
static boost::shared_ptr< SetUpSchur > createSetUpSchur(MoFEM::Interface &m_field)