v0.15.0
Loading...
Searching...
No Matches
ExampleArcLength.hpp
Date
2025-10-5
/**
* @file ExampleArcLength.hpp
* @example ExampleArcLength.hpp
* @date 2025-10-5
*
* @copyright Anonymous authors (c) 2025 under the MIT license (see LICENSE for
* details)
*
*/
MoFEMErrorCode runProblem();
protected:
MoFEMErrorCode ArcLengthSolve();
MoFEMErrorCode checkResults();
};
//! [Run problem]
MoFEMErrorCode ExampleArcLength::runProblem() {
ExampleTimeScale::snesPtr.reset();
}
//! [Run problem]
//! [Arc Length Solve]
auto *simple = mField.getInterface<Simple>();
auto *pipeline_mng = mField.getInterface<PipelineManager>();
auto dm = simple->getDM();
auto snes = pipeline_mng->createSNES();
ExampleTimeScale::snesPtr = snes;
PetscBool post_proc_vol;
PetscBool post_proc_skin;
if constexpr (SPACE_DIM == 2) {
post_proc_vol = PETSC_TRUE;
post_proc_skin = PETSC_FALSE;
} else {
post_proc_vol = PETSC_FALSE;
post_proc_skin = PETSC_TRUE;
}
CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-post_proc_vol",
&post_proc_vol, PETSC_NULLPTR);
CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-post_proc_skin",
&post_proc_skin, PETSC_NULLPTR);
// Setup postprocessing
auto create_post_proc_fe = [dm, this, simple, post_proc_vol,
post_proc_skin]() {
auto post_proc_ele_domain = [dm, this](auto &pip_domain) {
"GEOMETRY");
auto common_ptr = commonDataFactory<SPACE_DIM, GAUSS, DomainEleOp>(
mField, pip_domain, "U", "MAT_ELASTIC", Sev::inform);
return common_ptr;
};
auto post_proc_map = [this](auto &pip, auto u_ptr, auto common_ptr) {
using OpPPMap = OpPostProcMapInMoab<SPACE_DIM, SPACE_DIM>;
auto X_ptr = boost::make_shared<MatrixDouble>();
pip->getOpPtrVector().push_back(
new OpCalculateVectorFieldValues<SPACE_DIM>("GEOMETRY", X_ptr));
pip->getOpPtrVector().push_back(
new OpPPMap(pip->getPostProcMesh(), pip->getMapGaussPts(), {},
{{"U", u_ptr}, {"GEOMETRY", X_ptr}},
{{"GRAD", common_ptr->matGradPtr},
{"FIRST_PIOLA", common_ptr->getMatFirstPiolaStress()}},
{}));
};
auto push_post_proc_bdy = [dm, this, simple, post_proc_skin,
&post_proc_ele_domain,
&post_proc_map](auto &pip_bdy) {
if (post_proc_skin == PETSC_FALSE)
return boost::shared_ptr<PostProcEleBdy>();
auto u_ptr = boost::make_shared<MatrixDouble>();
pip_bdy->getOpPtrVector().push_back(
new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_ptr));
auto op_loop_side =
new OpLoopSide<SideEle>(mField, simple->getDomainFEName(), SPACE_DIM);
auto common_ptr = post_proc_ele_domain(op_loop_side->getOpPtrVector());
pip_bdy->getOpPtrVector().push_back(op_loop_side);
CHKERR post_proc_map(pip_bdy, u_ptr, common_ptr);
return pip_bdy;
};
auto push_post_proc_domain = [dm, this, simple, post_proc_vol,
&post_proc_ele_domain,
&post_proc_map](auto &pip_domain) {
if (post_proc_vol == PETSC_FALSE)
return boost::shared_ptr<PostProcEleDomain>();
auto u_ptr = boost::make_shared<MatrixDouble>();
pip_domain->getOpPtrVector().push_back(
new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_ptr));
auto common_ptr = post_proc_ele_domain(pip_domain->getOpPtrVector());
CHKERR post_proc_map(pip_domain, u_ptr, common_ptr);
return pip_domain;
};
auto post_proc_fe_domain = boost::make_shared<PostProcEleDomain>(mField);
auto post_proc_fe_bdy = boost::make_shared<PostProcEleBdy>(mField);
return std::make_pair(push_post_proc_domain(post_proc_fe_domain),
push_post_proc_bdy(post_proc_fe_bdy));
};
auto add_arc_length_load_tangent_operator = [&](auto snes) {
auto bdy_pipeline = boost::make_shared<BoundaryEle>(mField);
auto domain_pipeline = boost::make_shared<DomainEle>(mField);
auto integration_rule = [](int, int, int approx_order) {
return 2 * (approx_order - 1);
};
bdy_pipeline->getRuleHook = integration_rule;
domain_pipeline->getRuleHook = integration_rule;
struct TimeOne : TimeScale {
double getScale(const double) override { return -1.; }
};
auto time_one = boost::make_shared<TimeOne>();
// The trick is here, that is assumed that load is applied linearly
// proportional to the load factor (time in this case). So derivative of the
// load with respect to the load factor is the same as the load at time
// = 1.0
CHKERR BoundaryNaturalBC::AddFluxToPipeline<OpForce>::add(
bdy_pipeline->getOpPtrVector(), mField, "U", {time_one}, "FORCE",
Sev::inform);
//! [Define gravity vector]
CHKERR DomainNaturalBC::AddFluxToPipeline<OpBodyForce>::add(
domain_pipeline->getOpPtrVector(), mField, "U", {time_one},
"BODY_FORCE", Sev::inform);
auto snes_ctx_ptr = getDMSnesCtx(simple->getDM());
snes_ctx_ptr->getLoadTangent().push_back(
{simple->getBoundaryFEName(), bdy_pipeline});
snes_ctx_ptr->getLoadTangent().push_back(
{simple->getDomainFEName(), domain_pipeline});
auto time_scale = boost::make_shared<ExampleTimeScale>();
auto pre_proc_ptr = boost::make_shared<FEMethod>();
auto get_bc_hook_rhs = [this, pre_proc_ptr, time_scale]() {
CHKERR EssentialPreProc<DisplacementCubitBcData>(mField, pre_proc_ptr,
{time_scale}, false)();
};
pre_proc_ptr->preProcessHook = get_bc_hook_rhs;
snes_ctx_ptr->getPreProcLoadTangent().push_back(pre_proc_ptr);
auto post_proc_rhs_ptr = boost::make_shared<FEMethod>();
auto get_post_proc_hook_rhs = [this, post_proc_rhs_ptr]() {
CHKERR EssentialPostProcRhs<DisplacementCubitBcData>(
mField, post_proc_rhs_ptr, 0.)();
};
post_proc_rhs_ptr->postProcessHook = get_post_proc_hook_rhs;
snes_ctx_ptr->getPostProcLoadTangent().push_back(post_proc_rhs_ptr);
CHKERR SNESNewtonALSetFunction(snes, SnesLoadTangent, snes_ctx_ptr.get());
};
auto add_extra_finite_elements_to_solver_pipelines = [&]() {
auto pre_proc_ptr = boost::make_shared<FEMethod>();
auto post_proc_rhs_ptr = boost::make_shared<FEMethod>();
auto post_proc_lhs_ptr = boost::make_shared<FEMethod>();
auto time_scale = boost::make_shared<ExampleTimeScale>();
auto get_bc_hook_rhs = [this, pre_proc_ptr, time_scale]() {
CHKERR EssentialPreProc<DisplacementCubitBcData>(mField, pre_proc_ptr,
{time_scale}, false)();
};
pre_proc_ptr->preProcessHook = get_bc_hook_rhs;
auto get_post_proc_hook_rhs = [this, post_proc_rhs_ptr]() {
CHKERR EssentialPreProcReaction<DisplacementCubitBcData>(
mField, post_proc_rhs_ptr, nullptr, Sev::verbose)();
CHKERR EssentialPostProcRhs<DisplacementCubitBcData>(
mField, post_proc_rhs_ptr, 1.)();
};
auto get_post_proc_hook_lhs = [this, post_proc_lhs_ptr]() {
CHKERR EssentialPostProcLhs<DisplacementCubitBcData>(
mField, post_proc_lhs_ptr, 1.)();
};
post_proc_rhs_ptr->postProcessHook = get_post_proc_hook_rhs;
post_proc_lhs_ptr->postProcessHook = get_post_proc_hook_lhs;
// This is low level pushing finite elements (pipelines) to solver
auto snes_ctx_ptr = getDMSnesCtx(simple->getDM());
snes_ctx_ptr->getPreProcComputeRhs().push_front(pre_proc_ptr);
snes_ctx_ptr->getPreProcSetOperators().push_front(pre_proc_ptr);
snes_ctx_ptr->getPostProcComputeRhs().push_back(post_proc_rhs_ptr);
snes_ctx_ptr->getPostProcSetOperators().push_back(post_proc_lhs_ptr);
};
auto set_log_and_post_hook = [this](auto post_proc_fe, auto post_proc_bdy) {
auto make_vtk = [this, post_proc_fe, post_proc_bdy](auto iter,
auto cache_ptr) {
auto simple = mField.getInterface<Simple>();
auto dm = simple->getDM();
if (post_proc_fe) {
CHKERR DMoFEMLoopFiniteElements(dm, "dFE", post_proc_fe, cache_ptr);
CHKERR post_proc_fe->writeFile(
"out_step_" + boost::lexical_cast<std::string>(iter) + ".h5m");
}
if (post_proc_bdy) {
CHKERR DMoFEMLoopFiniteElements(dm, "bFE", post_proc_bdy, cache_ptr);
CHKERR post_proc_bdy->writeFile(
"out_skin_" + boost::lexical_cast<std::string>(iter) + ".h5m");
}
};
auto hook = [this, make_vtk](SNES snes, Vec x, Vec F,
boost::shared_ptr<CacheTuple> cache_ptr,
void *ctx) -> MoFEMErrorCode {
double atol, rtol, stol;
int maxit, maxf;
CHKERR SNESGetTolerances(snes, &atol, &rtol, &stol, &maxit, &maxf);
double lambda;
CHKERR SNESNewtonALGetLoadParameter(snes, &lambda);
double nrm;
CHKERR VecNorm(F, NORM_2, &nrm);
int iter;
CHKERR SNESGetIterationNumber(snes, &iter);
SNESConvergedReason reason;
CHKERR SNESGetConvergedReason(snes, &reason);
MOFEM_LOG_C("EXAMPLE", Sev::inform,
"SNES iteration %d, reason %d nrm2 %g lambda %g", iter,
reason, nrm, lambda);
CHKERR make_vtk(iter, cache_ptr);
};
auto simple = mField.getInterface<Simple>();
auto snes_ctx_ptr = getDMSnesCtx(simple->getDM());
snes_ctx_ptr->getRhsHook() = hook;
};
// Set arc length solver
auto D = createDMVector(simple->getDM());
CHKERR SNESSetType(snes, SNESNEWTONLS);
CHKERR SNESSetFromOptions(snes);
auto B = createDMMatrix(dm);
CHKERR SNESSetJacobian(snes, B, B, PETSC_NULLPTR, PETSC_NULLPTR);
// Add extra finite elements to SNES solver pipelines to resolve essential
// boundary conditions
CHKERR add_extra_finite_elements_to_solver_pipelines();
CHKERR add_arc_length_load_tangent_operator(snes);
auto [post_proc_fe, post_proc_bdy] = create_post_proc_fe();
CHKERR set_log_and_post_hook(post_proc_fe, post_proc_bdy);
CHKERR SNESSolve(snes, PETSC_NULLPTR, D);
}
//! [Arc Length Solve]
//! [Check]
PetscInt test_nb = 0;
PetscBool test_flg = PETSC_FALSE;
CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-test", &test_nb, &test_flg);
if (test_flg) {
auto simple = mField.getInterface<Simple>();
auto T = createDMVector(simple->getDM());
CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
SCATTER_FORWARD);
double nrm2;
CHKERR VecNorm(T, NORM_2, &nrm2);
MOFEM_LOG("EXAMPLE", Sev::verbose) << "Regression norm " << nrm2;
double regression_value = 0;
switch (test_nb) {
case 1:
regression_value = 1.215134e-04;
break;
default:
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Wrong test number.");
break;
}
if (fabs(nrm2 - regression_value) > 1e-2)
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"Regression test field; wrong norm value. %6.4e != %6.4e", nrm2,
regression_value);
}
}
//! [Check]
#define MOFEM_LOG_C(channel, severity, format,...)
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
constexpr int SPACE_DIM
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ 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 ...
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ F
auto integration_rule
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 DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition DMMoFEM.cpp:576
auto createDMVector(DM dm)
Get smart vector from DM.
Definition DMMoFEM.hpp:1102
auto createDMMatrix(DM dm)
Get smart matrix from DM.
Definition DMMoFEM.hpp:1059
#define MOFEM_LOG(channel, severity)
Log.
static double lambda
double D
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
auto getDMSnesCtx(DM dm)
Get SNES context data structure used by DM.
Definition DMMoFEM.hpp:1130
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
static constexpr int approx_order
MoFEMErrorCode checkResults()
[Arc Length Solve]
MoFEMErrorCode ArcLengthSolve()
[Run problem]
ExampleArcLength(MoFEM::Interface &m_field)
MoFEMErrorCode runProblem()
[Run problem]
MoFEMErrorCode readMesh()
[Run problem]
MoFEMErrorCode gettingNorms()
[TS Solve]
MoFEMErrorCode outputResults()
[Getting norms]
MoFEMErrorCode assembleSystem()
[Boundary condition]
MoFEMErrorCode setupProblem()
[Read mesh]
ExampleNonlinearElastic(MoFEM::Interface &m_field)
MoFEMErrorCode boundaryCondition()
[Set up problem]
Deprecated interface functions.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.