20 ? AssemblyType::BLOCK_SCHUR
24 IntegrationType::GAUSS;
89 static boost::shared_ptr<SetUpSchur>
97static char help[] =
"...\n\n";
99int main(
int argc,
char *argv[]) {
102 const char param_file[] =
"param_file.petsc";
105 auto core_log = logging::core::get();
117 DMType dm_name =
"DMMOFEM";
119 DMType dm_name_mg =
"DMMOFEM_MG";
124 moab::Core mb_instance;
125 moab::Interface &moab = mb_instance;
149 if (
A == AssemblyType::BLOCK_SCHUR ||
A == AssemblyType::SCHUR) {
151 CHKERR schur_ptr->setUp(solver);
155 BOOST_LOG_SCOPED_THREAD_ATTR(
"Timeline", attrs::timer());
156 MOFEM_LOG(
"TIMER", Sev::inform) <<
"KSPSetUp";
158 MOFEM_LOG(
"TIMER", Sev::inform) <<
"KSPSetUp <= Done";
161 CHKERR KSPGetDM(solver, &dm);
165 MOFEM_LOG(
"TIMER", Sev::inform) <<
"KSPSolve";
167 MOFEM_LOG(
"TIMER", Sev::inform) <<
"KSPSolve <= Done";
169 CHKERR VecGhostUpdateBegin(
D, INSERT_VALUES, SCATTER_FORWARD);
170 CHKERR VecGhostUpdateEnd(
D, INSERT_VALUES, SCATTER_FORWARD);
182 "Is expected that schur matrix is not allocated. This is "
183 "possible only is PC is set up twice");
211 CHKERR KSPGetPC(solver, &pc);
212 PetscBool is_pcfs = PETSC_FALSE;
213 PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &is_pcfs);
218 "Is expected that schur matrix is not allocated. This is "
219 "possible only is PC is set up twice");
226 CHKERR MatSetDM(
S, PETSC_NULLPTR);
228 CHKERR MatSetOption(
S, MAT_SYMMETRIC, PETSC_TRUE);
233 if constexpr (
A == AssemblyType::BLOCK_SCHUR) {
236 CHKERR KSPGetDM(solver, &solver_dm);
237 CHKERR DMSetMatType(solver_dm, MATSHELL);
266 auto create_dm = [&](
const char *name,
auto &ents,
auto dm_type) {
268 auto create_dm_imp = [&]() {
273 auto sub_ents_ptr = boost::make_shared<Range>(ents);
281 "Error in creating schurDM. It is possible that schurDM is "
289 if constexpr (
A == AssemblyType::BLOCK_SCHUR) {
291 auto get_nested_mat_data = [&]() -> boost::shared_ptr<NestSchurData> {
292 auto block_mat_data =
297 simple->getDomainFEName(),
309 {
"U"}, {boost::make_shared<Range>(
volEnts)}
314 auto nested_mat_data = get_nested_mat_data();
332 {
"U"}, {boost::make_shared<Range>(
volEnts)}, ao_up,
S,
true,
true));
336 {
"U"}, {boost::make_shared<Range>(
volEnts)}, ao_up,
S,
true,
true));
338 auto pre_proc_schur_lhs_ptr = boost::make_shared<FEMethod>();
339 auto post_proc_schur_lhs_ptr = boost::make_shared<FEMethod>();
341 pre_proc_schur_lhs_ptr->preProcessHook = [
this]() {
346 MOFEM_LOG(
"TIMER", Sev::inform) <<
"Lhs Assemble Begin";
350 post_proc_schur_lhs_ptr->postProcessHook = [
this, post_proc_schur_lhs_ptr,
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";
362 ksp_ctx_ptr->getPreProcSetOperators().push_front(pre_proc_schur_lhs_ptr);
363 ksp_ctx_ptr->getPostProcSetOperators().push_back(post_proc_schur_lhs_ptr);
374 CHKERR PCFieldSplitSetIS(pc, NULL, vol_is);
375 CHKERR PCFieldSplitSetSchurPre(pc, PC_FIELDSPLIT_SCHUR_PRE_USER,
S);
382 auto get_pc = [](
auto ksp) {
384 CHKERR KSPGetPC(ksp, &pc_raw);
388 if constexpr (
A == AssemblyType::BLOCK_SCHUR) {
392 CHKERR PCSetOperators(pc, A, P);
395 CHKERR PCFieldSplitSchurGetSubKSP(pc, PETSC_NULLPTR, &subksp);
398 auto set_pc_p_mg = [](
auto dm,
auto pc,
auto S) {
401 PetscBool same = PETSC_FALSE;
402 PetscObjectTypeCompare((PetscObject)pc, PCMG, &same);
404 MOFEM_LOG(
"TIMER", Sev::inform) <<
"Set up MG";
407 CHKERR PCSetFromOptions(pc);
412 auto set_pc_ksp = [&](
auto dm,
auto pc,
auto S) {
414 PetscBool same = PETSC_FALSE;
415 PetscObjectTypeCompare((PetscObject)pc, PCKSP, &same);
417 MOFEM_LOG(
"TIMER", Sev::inform) <<
"Set up inner KSP for PCKSP";
418 CHKERR PCSetFromOptions(pc);
420 CHKERR PCKSPGetKSP(pc, &inner_ksp);
421 CHKERR KSPSetFromOptions(inner_ksp);
423 CHKERR KSPGetPC(inner_ksp, &inner_pc);
424 CHKERR PCSetFromOptions(inner_pc);
425 CHKERR set_pc_p_mg(dm, inner_pc,
S);
438boost::shared_ptr<SetUpSchur>
440 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)
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
#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
#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
PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[])
Must be called by user to set Sub DM MoFEM data structures.
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem)
set squared problem
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[])
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
MoFEMErrorCode DMRegister_MGViaApproxOrders(const char sname[])
Register DM for Multi-Grid via approximation orders.
auto createDMVector(DM dm)
Get smart vector from DM.
PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[])
auto createDMMatrix(DM dm)
Get smart matrix from DM.
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.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
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.
MoFEMErrorCode setSchurA00MatSolvePC(SmartPetscObj< PC > pc)
Set PC for A00 block.
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.
auto getDMSubData(DM dm)
Get sub problem data structure.
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.
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.
auto createAOMappingIS(IS isapp, IS ispetsc)
Creates an application mapping using two index sets.
MoFEMErrorCode DMMoFEMSetNestSchurData(DM dm, boost::shared_ptr< NestSchurData >)
auto createDMNestSchurMat(DM dm)
auto createDM(MPI_Comm comm, const std::string dm_type_name)
Creates smart DM object.
OpSchurAssembleBase * createOpSchurAssembleBegin()
auto createDMBlockMat(DM dm)
#define EXECUTABLE_DIMENSION
ElementsAndOps< SPACE_DIM >::SideEle SideEle
constexpr int SPACE_DIM
[Define dimension]
constexpr IntegrationType I
constexpr AssemblyType A
[Define dimension]
Boundary conditions marker.
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
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
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.
Section manager is used to create indexes and sections.
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.
Template struct for dimension-specific finite element types.
PipelineManager interface.
Simple interface for fast problem set-up.
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
MoFEMErrorCode createSubDM()
MoFEMErrorCode setDiagonalPC(PC pc)
SetUpSchurImpl(MoFEM::Interface &m_field)
SmartPetscObj< DM > schurDM
MoFEMErrorCode setEntities()
MoFEMErrorCode setUp(SmartPetscObj< KSP >)
SmartPetscObj< DM > blockDM
MoFEMErrorCode setPC(PC pc)
MoFEMErrorCode setOperator()
MoFEM::Interface & mField
[Push operators to pipeline]
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)