v0.14.0
Loading...
Searching...
No Matches
schur_test_diag_mat.cpp

Test block matrix and Schur complement matrix.

Test block matrix and Schur complement matrix.

/**
* @file schur_test_diag_mat.cpp
* @example schur_test_diag_mat.cpp
* @brief Test block matrix and Schur complement matrix.
*
* @copyright Copyright (c) 2024
*
*/
#include <MoFEM.hpp>
using namespace MoFEM;
static char help[] = "...\n\n";
template <int DIM> struct ElementsAndOps {};
template <> struct ElementsAndOps<2> {
};
template <> struct ElementsAndOps<3> {
using DomainEle = VolumeElementForcesAndSourcesCore;
};
//! [Define dimension]
constexpr int SPACE_DIM = 3;
constexpr int FIELD_DIM = SPACE_DIM;
constexpr AssemblyType A =
AssemblyType::BLOCK_MAT; //< selected assembly type
constexpr IntegrationType I =
IntegrationType::GAUSS; //< selected integration type
using DomainEleOp = DomainEle::UserDataOperator;
template <>
[](ForcesAndSourcesCore::UserDataOperator *op_ptr,
const EntitiesFieldData::EntData &row_data,
const EntitiesFieldData::EntData &col_data, MatrixDouble &m) {
return MatSetValues<AssemblyTypeSelector<PETSC>>(
petsc_mat, row_data, col_data, m, ADD_VALUES);
};
template <>
[](ForcesAndSourcesCore::UserDataOperator *op_ptr,
const EntitiesFieldData::EntData &row_data,
const EntitiesFieldData::EntData &col_data, MatrixDouble &m) {
return MatSetValues<AssemblyTypeSelector<BLOCK_MAT>>(
block_mat, row_data, col_data, m, ADD_VALUES);
};
template <>
[](ForcesAndSourcesCore::UserDataOperator *op_ptr,
const EntitiesFieldData::EntData &row_data,
const EntitiesFieldData::EntData &col_data, MatrixDouble &m) {
return MatSetValues<AssemblyTypeSelector<BLOCK_SCHUR>>(
block_mat, row_data, col_data, m, ADD_VALUES);
};
template <>
DomainEleOp>::MatSetValuesHook
DomainEleOp>::matSetValuesHook =
[](ForcesAndSourcesCore::UserDataOperator *op_ptr,
const EntitiesFieldData::EntData &row_data,
const EntitiesFieldData::EntData &col_data, MatrixDouble &m) {
return MatSetValues<AssemblyTypeSelector<BLOCK_PRECONDITIONER_SCHUR>>(
block_mat, row_data, col_data, m, ADD_VALUES);
};
constexpr bool debug = false;
int main(int argc, char *argv[]) {
// initialize petsc
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
// Create MoAB database
moab::Core moab_core;
moab::Interface &moab = moab_core;
// Create MoFEM database and link it to MoAB
MoFEM::Core mofem_core(moab);
MoFEM::Interface &m_field = mofem_core;
// Register DM Manager
DMType dm_name = "DMMOFEM";
// Add logging channel for example
auto core_log = logging::core::get();
core_log->add_sink(
LogManager::createSink(LogManager::getStrmWorld(), "Timeline"));
LogManager::setLog("Timeline");
MOFEM_LOG_TAG("Timeline", "Timeline");
// Simple interface
auto simple = m_field.getInterface<Simple>();
// get options from command line
CHKERR simple->getOptions();
// load mesh file
CHKERR simple->loadFile();
// set fields order, i.e. for most first cases order is sufficient.
int order = 3;
CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
CHKERR simple->setFieldOrder("V", order);
CHKERR simple->setFieldOrder("T", order);
CHKERR simple->setFieldOrder("S", order - 1);
CHKERR simple->setFieldOrder("O", order - 2);
// setup problem
CHKERR simple->setUp();
// If block "VOL" exist, first index is removed from "T" field. That test
// blocks handling when some DOFs are removed.
auto bc_mng = m_field.getInterface<BcManager>();
CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "VOL",
"T", 0, 1);
auto schur_dm = createDM(m_field.get_comm(), "DMMOFEM");
CHKERR DMMoFEMCreateSubDM(schur_dm, simple->getDM(), "SCHUR");
CHKERR DMMoFEMSetSquareProblem(schur_dm, PETSC_TRUE);
CHKERR DMMoFEMAddElement(schur_dm, simple->getDomainFEName());
CHKERR DMMoFEMAddSubFieldRow(schur_dm, "V");
CHKERR DMMoFEMAddSubFieldCol(schur_dm, "V");
CHKERR DMSetUp(schur_dm);
auto block_dm = createDM(m_field.get_comm(), "DMMOFEM");
CHKERR DMMoFEMCreateSubDM(block_dm, simple->getDM(), "BLOCK");
CHKERR DMMoFEMSetSquareProblem(block_dm, PETSC_TRUE);
CHKERR DMMoFEMAddElement(block_dm, simple->getDomainFEName());
CHKERR DMMoFEMAddSubFieldRow(block_dm, "T");
CHKERR DMMoFEMAddSubFieldCol(block_dm, "T");
CHKERR DMMoFEMAddSubFieldRow(block_dm, "S");
CHKERR DMMoFEMAddSubFieldCol(block_dm, "S");
CHKERR DMMoFEMAddSubFieldRow(block_dm, "O");
CHKERR DMMoFEMAddSubFieldCol(block_dm, "O");
CHKERR DMSetUp(block_dm);
auto S = createDMMatrix(schur_dm);
auto shell_data =
{{simple->getDomainFEName(),
{{"V", "V"},
{"T", "T"},
{"S", "S"},
{"O", "O"},
{"V", "T"},
{"T", "V"},
{"S", "T"},
{"T", "S"},
{"T", "O"},
{"O", "T"},
{"S", "O"},
{"O", "S"}
}}}
)
);
auto [mat, block_data_ptr] = shell_data;
block_mat = mat;
std::vector<std::string> fields{"T", "S", "O"};
auto [nested_mat, nested_data_ptr] = createSchurNestedMatrix(
{schur_dm, block_dm}, block_data_ptr,
fields, {nullptr, nullptr, nullptr}, true
)
);
using OpMassPETSCAssemble = FormsIntegrators<DomainEleOp>::Assembly<
using OpMassBlockAssemble = FormsIntegrators<DomainEleOp>::Assembly<
using OpMassBlockPreconditionerAssemble =
BiLinearForm<GAUSS>::OpMass<1, FIELD_DIM>;
// get operators tester
auto pip_mng = m_field.getInterface<PipelineManager>(); // get interface to
// pipeline manager
auto close_zero = [](double, double, double) { return 1; };
auto beta = [](double, double, double) { return -1./2; };
auto gamma = [](double, double, double) { return -1. / 4; };
[](int, int, int o) { return 2 * o; });
auto &pip_lhs = pip_mng->getOpDomainLhsPipeline();
pip_lhs.push_back(new OpMassPETSCAssemble("V", "V"));
// pip_lhs.push_back(new OpMassPETSCAssemble("T", "T"));
pip_lhs.push_back(new OpMassPETSCAssemble("V", "T"));
pip_lhs.push_back(new OpMassPETSCAssemble("T", "V"));
pip_lhs.push_back(new OpMassPETSCAssemble("S", "S", close_zero));
pip_lhs.push_back(new OpMassPETSCAssemble("S", "T", beta));
pip_lhs.push_back(new OpMassPETSCAssemble("T", "S", beta));
pip_lhs.push_back(new OpMassPETSCAssemble("O", "O", close_zero));
pip_lhs.push_back(new OpMassPETSCAssemble("T", "O", beta));
pip_lhs.push_back(new OpMassPETSCAssemble("O", "T", beta));
pip_lhs.push_back(new OpMassPETSCAssemble("S", "O", gamma));
pip_lhs.push_back(new OpMassPETSCAssemble("O", "S", gamma));
pip_lhs.push_back(createOpSchurAssembleBegin());
pip_lhs.push_back(new OpMassBlockAssemble("V", "V"));
// pip_lhs.push_back(new OpMassBlockAssemble("T", "T"));
pip_lhs.push_back(new OpMassBlockAssemble("V", "T"));
pip_lhs.push_back(new OpMassBlockAssemble("T", "V"));
pip_lhs.push_back(new OpMassBlockAssemble("S", "S", close_zero));
pip_lhs.push_back(new OpMassBlockAssemble("S", "T", beta));
pip_lhs.push_back(new OpMassBlockAssemble("T", "S", beta));
pip_lhs.push_back(new OpMassBlockAssemble("O", "O", close_zero));
pip_lhs.push_back(new OpMassBlockAssemble("T", "O", beta));
pip_lhs.push_back(new OpMassBlockAssemble("O", "T", beta));
pip_lhs.push_back(new OpMassBlockAssemble("S", "O", gamma));
pip_lhs.push_back(new OpMassBlockAssemble("O", "S", gamma));
pip_lhs.push_back(new OpMassBlockPreconditionerAssemble("T", "T"));
auto schur_is = getDMSubData(schur_dm)->getSmartRowIs();
auto ao_up = createAOMappingIS(schur_is, PETSC_NULL);
pip_lhs.push_back(createOpSchurAssembleEnd(
fields, {nullptr, nullptr, nullptr}, ao_up, S, true, true)
);
{
MOFEM_LOG_CHANNEL("Timeline");
MOFEM_LOG_TAG("Timeline", "timer");
BOOST_LOG_SCOPED_THREAD_ATTR("Timeline", attrs::timer());
MOFEM_LOG("Timeline", Sev::inform) << "Assemble start";
simple->getDomainFEName(),
pip_mng->getDomainLhsFE());
MOFEM_LOG("Timeline", Sev::inform) << "Assemble end";
}
{
MOFEM_LOG_CHANNEL("Timeline");
MOFEM_LOG_TAG("Timeline", "timer");
BOOST_LOG_SCOPED_THREAD_ATTR("Timeline", attrs::timer());
MOFEM_LOG("Timeline", Sev::inform) << "Mat assemble start";
CHKERR MatAssemblyBegin(petsc_mat, MAT_FINAL_ASSEMBLY);
CHKERR MatAssemblyEnd(petsc_mat, MAT_FINAL_ASSEMBLY);
MOFEM_LOG("Timeline", Sev::inform) << "Mat assemble end";
}
auto get_random_vector = [&](auto dm) {
auto v = createDMVector(dm);
PetscRandom rctx;
PetscRandomCreate(PETSC_COMM_WORLD, &rctx);
CHK_MOAB_THROW(VecSetRandom(v, rctx), "generate rand vector");
PetscRandomDestroy(&rctx);
return v;
};
auto v = get_random_vector(simple->getDM());
auto y_petsc = createDMVector(simple->getDM());
auto y_block = createDMVector(simple->getDM());
auto test = [](auto msg, auto y, double norm0) {
double eps = 1e-10;
CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-eps", &eps, PETSC_NULL);
PetscReal norm;
CHKERR VecNorm(y, NORM_2, &norm);
norm = norm / norm0;
MOFEM_TAG_AND_LOG("WORLD", Sev::inform, "TestBlockMat")
<< msg << ": norm of difference: " << norm;
if (norm > eps || std::isnan(norm) || std::isinf(norm)) {
SETERRQ(PETSC_COMM_WORLD, 1, "norm of difference is too big");
}
};
std::vector<int> zero_rows_and_cols = {
0, 1, 10, 20,
500}; // not to remove dofs for TENSOR filed, inverse will not work
CHKERR MatZeroRowsColumns(petsc_mat, zero_rows_and_cols.size(),
&*zero_rows_and_cols.begin(), 1, PETSC_NULL,
PETSC_NULL);
CHKERR MatZeroRowsColumns(block_mat, zero_rows_and_cols.size(),
&*zero_rows_and_cols.begin(), 1, PETSC_NULL,
PETSC_NULL);
{
MOFEM_LOG_CHANNEL("Timeline");
MOFEM_LOG_TAG("Timeline", "timer");
BOOST_LOG_SCOPED_THREAD_ATTR("Timeline", attrs::timer());
MOFEM_LOG("Timeline", Sev::inform)
<< "MatMult(petsc_mat, v, y_petsc) star";
CHKERR MatMult(petsc_mat, v, y_petsc);
MOFEM_LOG("Timeline", Sev::inform)
<< "MatMult(petsc_mat, v, y_petsc) end";
}
double nrm0;
CHKERR VecNorm(y_petsc, NORM_2, &nrm0);
{
MOFEM_LOG_CHANNEL("Timeline");
MOFEM_LOG_TAG("Timeline", "timer");
BOOST_LOG_SCOPED_THREAD_ATTR("Timeline", attrs::timer());
MOFEM_LOG("Timeline", Sev::inform)
<< "MatMult(block_mat, v, y_block) star";
CHKERR MatMult(block_mat, v, y_block);
MOFEM_LOG("Timeline", Sev::inform)
<< "MatMult(block_mat, v, y_block) end";
}
CHKERR VecAXPY(y_petsc, -1.0, y_block);
CHKERR test("mult", y_petsc, nrm0);
MOFEM_LOG_CHANNEL("Timeline");
MOFEM_LOG_TAG("Timeline", "timer");
CHKERR MatMult(petsc_mat, v, y_petsc);
CHKERR MatMult(block_mat, v, y_block);
CHKERR MatMultAdd(petsc_mat, v, y_petsc, y_petsc);
CHKERR MatMultAdd(block_mat, v, y_block, y_block);
CHKERR VecAXPY(y_petsc, -1.0, y_block);
CHKERR test("mult add", y_petsc, nrm0);
CHKERR schurSwitchPreconditioner(std::get<1>(*nested_data_ptr)[3]);
auto y_nested = createDMVector(simple->getDM());
CHKERR MatMult(petsc_mat, v, y_petsc);
CHKERR MatMult(nested_mat, v, y_nested);
CHKERR VecAXPY(y_petsc, -1.0, y_nested);
CHKERR test("mult nested", y_petsc, nrm0);
CHKERR schurSwitchPreconditioner(std::get<1>(*nested_data_ptr)[3]);
auto diag_mat = std::get<0>(*nested_data_ptr)[3];
auto diag_block_x = get_random_vector(block_dm);
auto diag_block_f = createDMVector(block_dm);
auto block_solved_x = createDMVector(block_dm);
CHKERR MatMult(diag_mat, diag_block_x, diag_block_f);
// That is if one like to use MatSolve directly, not though PC, as it is
// below
// CHKERR MatSolve(diag_mat, diag_block_f, block_solved_x);
// set matrix type to shell, set data
CHKERR DMSetMatType(block_dm, MATSHELL);
CHKERR DMMoFEMSetBlocMatData(block_dm, std::get<1>(*nested_data_ptr)[3]);
// set empty operator, since block data are already calculated
CHKERR DMKSPSetComputeOperators(
block_dm,
[](KSP, Mat, Mat, void *) {
MOFEM_LOG("WORLD", Sev::inform) << "empty operator";
return 0;
},
nullptr);
auto ksp = createKSP(m_field.get_comm());
CHKERR KSPSetDM(ksp, block_dm);
CHKERR KSPSetFromOptions(ksp);
// set preconditioner to block mat
auto get_pc = [](auto ksp) {
PC pc_raw;
CHKERR KSPGetPC(ksp, &pc_raw);
return SmartPetscObj<PC>(pc_raw, true); // bump reference
};
CHKERR KSPSetUp(ksp);
CHKERR VecZeroEntries(block_solved_x);
CHKERR KSPSolve(ksp, diag_block_f, block_solved_x);
auto diag_block_f_test = createDMVector(block_dm);
CHKERR MatMult(diag_mat, block_solved_x, diag_block_f_test);
CHKERR VecAXPY(diag_block_f_test, -1.0, diag_block_f);
CHKERR test("diag solve", diag_block_f_test, nrm0);
if (m_field.get_comm_rank() == 0) {
CHKERR schurSaveBlockMesh(block_data_ptr, "block_mesh.vtk");
}
petsc_mat.reset();
block_mat.reset();
}
// finish work cleaning memory, getting statistics, etc.
return 0;
}
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
static char help[]
int main()
static const double eps
constexpr int SPACE_DIM
constexpr int FIELD_DIM
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
[Define dimension]
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base nme:nme847.
Definition definitions.h:60
@ L2
field with C-1 continuity
Definition definitions.h:88
@ 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 ...
#define CHK_MOAB_THROW(err, msg)
Check error code of MoAB function and throw MoFEM exception.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
constexpr int order
static const bool debug
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:497
PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem)
set squared problem
Definition DMMoFEM.cpp:456
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
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition DMMoFEM.cpp:586
auto createDMVector(DM dm)
Get smart vector from DM.
Definition DMMoFEM.hpp:1099
PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[])
Definition DMMoFEM.cpp:280
auto createDMMatrix(DM dm)
Get smart matrix from DM.
Definition DMMoFEM.hpp:1056
IntegrationType
Form integrator integration types.
AssemblyType
[Storage and set boundary conditions]
@ BLOCK_PRECONDITIONER_SCHUR
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
const double v
phase velocity of light in medium (cm/ns)
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
auto createKSP(MPI_Comm comm)
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:2186
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
MoFEMErrorCode schurSwitchPreconditioner(boost::shared_ptr< BlockStructure > block_mat_data)
Switch preconditioner.
Definition Schur.cpp:2365
MoFEMErrorCode setSchurA00MatSolvePC(SmartPetscObj< PC > pc)
Set PC for A00 block.
Definition Schur.cpp:2223
std::pair< SmartPetscObj< Mat >, boost::shared_ptr< NestSchurData > > createSchurNestedMatrix(boost::shared_ptr< NestSchurData > schur_net_data_ptr)
Create a Mat Diag Blocks object.
Definition Schur.cpp:2153
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)
auto getDMSubData(DM dm)
Get sub problem data structure.
Definition DMMoFEM.hpp:1157
MoFEMErrorCode DMMoFEMSetBlocMatData(DM dm, boost::shared_ptr< BlockStructure >)
Set data for block mat.
Definition DMMoFEM.cpp:1533
boost::shared_ptr< BlockStructure > createBlockMatStructure(DM dm, SchurFEOpsFEandFields schur_fe_op_vec)
Create a Mat Diag Blocks object.
Definition Schur.cpp:1009
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:1944
auto createAOMappingIS(IS isapp, IS ispetsc)
Creates an application mapping using two index sets.
MoFEMErrorCode schurSaveBlockMesh(boost::shared_ptr< BlockStructure > block_mat_data, std::string filename)
Save block matrix as a mesh.
Definition Schur.cpp:2380
auto createDM(MPI_Comm comm, const std::string dm_type_name)
Creates smart DM object.
SchurShellMatData createBlockMat(DM dm, boost::shared_ptr< BlockStructure > data)
Create a Schur Mat object.
Definition Schur.cpp:1379
OpSchurAssembleBase * createOpSchurAssembleBegin()
Definition Schur.cpp:2181
constexpr IntegrationType I
SmartPetscObj< Mat > block_mat
SmartPetscObj< Mat > petsc_mat
FTensor::Index< 'm', 3 > m
Simple interface for fast problem set-up.
Definition BcManager.hpp:25
MoFEMErrorCode removeBlockDOFsOnEntities(const std::string problem_name, const std::string block_name, const std::string field_name, int lo, int hi, bool get_low_dim_ents=true, bool is_distributed_mesh=true)
Remove DOFs from problem.
Definition BcManager.cpp:73
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() 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:112
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
boost::function< MoFEMErrorCode( ForcesAndSourcesCore::UserDataOperator *op_ptr, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, MatrixDouble &m)> MatSetValuesHook
PipelineManager interface.
MoFEMErrorCode setDomainLhsIntegrationRule(RuleHookFun rule)
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.