v0.14.0
seepage.cpp

Hydromechanical problem

/**
* \file seepage.cpp
* \example seepage.cpp
*
* Hydromechanical problem
*
*/
/* This file is part of MoFEM.
* MoFEM is free software: you can redistribute it and/or modify it under
* the terms of the GNU Lesser General Public License as published by the
* Free Software Foundation, either version 3 of the License, or (at your
* option) any later version.
*
* MoFEM is distributed in the hope that it will be useful, but WITHOUT
* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
* FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
* License for more details.
*
* You should have received a copy of the GNU Lesser General Public
* License along with MoFEM. If not, see <http://www.gnu.org/licenses/>. */
#ifndef EXECUTABLE_DIMENSION
#define EXECUTABLE_DIMENSION 2
#endif
#include <MoFEM.hpp>
using namespace MoFEM;
constexpr int SPACE_DIM =
EXECUTABLE_DIMENSION; //< Space dimension of problem, mesh
using BoundaryEle =
//! [Only used with Hooke equation (linear material model)]
GAUSS>::OpGradSymTensorGrad<1, SPACE_DIM, SPACE_DIM, 0>;
PETSC>::LinearForm<GAUSS>::OpGradTimesSymTensor<1, SPACE_DIM, SPACE_DIM>;
//! [Only used with Hooke equation (linear material model)]
//! [Only used for dynamics]
PETSC>::LinearForm<GAUSS>::OpBaseTimesVector<1, SPACE_DIM, 1>;
//! [Only used for dynamics]
//! [Only used with Hencky/nonlinear material]
GAUSS>::OpGradTensorGrad<1, SPACE_DIM, SPACE_DIM, 1>;
//! [Only used with Hencky/nonlinear material]
//! [Essential boundary conditions]
PETSC>::LinearForm<GAUSS>::OpBaseTimesVector<1, SPACE_DIM, 0>;
PETSC>::LinearForm<GAUSS>::OpBaseTimesVector<1, SPACE_DIM, 1>;
//! [Essential boundary conditions]
// Thermal operators
/**
* @brief Integrate Lhs base of flux (1/k) base of flux (FLUX x FLUX)
*
*/
/**
* @brief Integrate Lhs div of base of flux time base of temperature (FLUX x T)
* and transpose of it, i.e. (T x FLAX)
*
*/
GAUSS>::OpMixDivTimesScalar<SPACE_DIM>;
/**
* @brief Integrate Lhs base of temperature times (heat capacity) times base of
* temperature (T x T)
*
*/
/**
* @brief Integrating Rhs flux base (1/k) flux (FLUX)
*/
GAUSS>::OpBaseTimesVector<3, 3, 1>;
/**
* @brief Integrate Rhs div flux base times temperature (T)
*
*/
GAUSS>::OpMixDivTimesU<3, 1, 2>;
/**
* @brief Integrate Rhs base of temperature time heat capacity times heat rate
* (T)
*
*/
GAUSS>::OpBaseTimesScalarField<1>;
/**
* @brief Integrate Rhs base of temperature times divergent of flux (T)
*
*/
//! [Body and heat source]
using OpBodyForce =
using OpHeatSource =
//! [Body and heat source]
//! [Natural boundary conditions]
//! [Natural boundary conditions]
//! [Essential boundary conditions (Least square approach)]
GAUSS>::OpEssentialRhs<HeatFluxCubitBcData, 3, SPACE_DIM>;
GAUSS>::OpEssentialLhs<HeatFluxCubitBcData, 3, SPACE_DIM>;
//! [Essential boundary conditions (Least square approach)]
double scale = 1.;
double default_poisson_ratio = 0.25; // zero for verification
double fluid_density = 10;
// #include <OpPostProcElastic.hpp>
#include <SeepageOps.hpp>
auto save_range = [](moab::Interface &moab, const std::string name,
const Range r) {
auto out_meshset = get_temp_meshset_ptr(moab);
CHKERR moab.add_entities(*out_meshset, r);
CHKERR moab.write_file(name.c_str(), "VTK", "", out_meshset->get_ptr(), 1);
};
struct Seepage {
Seepage(MoFEM::Interface &m_field) : mField(m_field) {}
MoFEMErrorCode runProblem();
private:
MoFEMErrorCode setupProblem();
MoFEMErrorCode createCommonData();
MoFEMErrorCode tsSolve();
PetscBool doEvalField;
std::array<double, SPACE_DIM> fieldEvalCoords;
std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uXScatter;
std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uYScatter;
std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uZScatter;
struct BlockedParameters
: public boost::enable_shared_from_this<BlockedParameters> {
double fluidConductivity;
inline auto getDPtr() {
return boost::shared_ptr<MatrixDouble>(shared_from_this(), &mD);
}
inline auto getConductivityPtr() {
return boost::shared_ptr<double>(shared_from_this(), &fluidConductivity);
}
};
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
std::string block_elastic_name, std::string block_thermal_name,
boost::shared_ptr<BlockedParameters> blockedParamsPtr, Sev sev);
};
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
std::string block_elastic_name, std::string block_thermal_name,
boost::shared_ptr<BlockedParameters> blockedParamsPtr, Sev sev) {
struct OpMatElasticBlocks : public DomainEleOp {
OpMatElasticBlocks(boost::shared_ptr<MatrixDouble> m, double bulk_modulus_K,
Sev sev,
std::vector<const CubitMeshSets *> meshset_vec_ptr)
: DomainEleOp(NOSPACE, DomainEleOp::OPSPACE), matDPtr(m),
bulkModulusKDefault(bulk_modulus_K),
shearModulusGDefault(shear_modulus_G) {
CHK_THROW_MESSAGE(extractElasticBlockData(m_field, meshset_vec_ptr, sev),
"Can not get data from block");
}
MoFEMErrorCode doWork(int side, EntityType type,
for (auto &b : blockData) {
if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
CHKERR getMatDPtr(matDPtr, b.bulkModulusK, b.shearModulusG);
}
}
CHKERR getMatDPtr(matDPtr, bulkModulusKDefault, shearModulusGDefault);
}
private:
boost::shared_ptr<MatrixDouble> matDPtr;
struct BlockData {
double bulkModulusK;
double shearModulusG;
Range blockEnts;
};
double bulkModulusKDefault;
double shearModulusGDefault;
std::vector<BlockData> blockData;
extractElasticBlockData(MoFEM::Interface &m_field,
std::vector<const CubitMeshSets *> meshset_vec_ptr,
Sev sev) {
for (auto m : meshset_vec_ptr) {
MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Elastic Block") << *m;
std::vector<double> block_data;
CHKERR m->getAttributes(block_data);
if (block_data.size() < 2) {
SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
"Expected that block has two attributes");
}
auto get_block_ents = [&]() {
Range ents;
m_field.get_moab().get_entities_by_handle(m->meshset, ents, true);
return ents;
};
double young_modulus = block_data[0];
double poisson_ratio = block_data[1];
double bulk_modulus_K = young_modulus / (3 * (1 - 2 * poisson_ratio));
MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Elastic Block")
<< m->getName() << ": E = " << young_modulus
<< " nu = " << poisson_ratio;
blockData.push_back(
{bulk_modulus_K, shear_modulus_G, get_block_ents()});
}
}
MoFEMErrorCode getMatDPtr(boost::shared_ptr<MatrixDouble> mat_D_ptr,
double bulk_modulus_K, double shear_modulus_G) {
//! [Calculate elasticity tensor]
auto set_material_stiffness = [&]() {
double A = (SPACE_DIM == 2)
: 1;
auto t_D = getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*mat_D_ptr);
t_D(i, j, k, l) =
2 * shear_modulus_G * ((t_kd(i, k) ^ t_kd(j, l)) / 4.) +
A * (bulk_modulus_K - (2. / 3.) * shear_modulus_G) * t_kd(i, j) *
t_kd(k, l);
};
//! [Calculate elasticity tensor]
constexpr auto size_symm = (SPACE_DIM * (SPACE_DIM + 1)) / 2;
mat_D_ptr->resize(size_symm * size_symm, 1);
set_material_stiffness();
}
};
double default_bulk_modulus_K =
double default_shear_modulus_G =
pipeline.push_back(new OpMatElasticBlocks(
blockedParamsPtr->getDPtr(), default_bulk_modulus_K,
default_bulk_modulus_K, mField, sev,
// Get blockset using regular expression
mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
(boost::format("%s(.*)") % block_elastic_name).str()
))
));
struct OpMatFluidBlocks : public DomainEleOp {
OpMatFluidBlocks(boost::shared_ptr<double> conductivity_ptr,
MoFEM::Interface &m_field, Sev sev,
std::vector<const CubitMeshSets *> meshset_vec_ptr)
conductivityPtr(conductivity_ptr) {
CHK_THROW_MESSAGE(extractThermalBlockData(m_field, meshset_vec_ptr, sev),
"Can not get data from block");
}
MoFEMErrorCode doWork(int side, EntityType type,
for (auto &b : blockData) {
if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
*conductivityPtr = b.conductivity;
}
}
*conductivityPtr = default_conductivity;
}
private:
struct BlockData {
double conductivity;
Range blockEnts;
};
std::vector<BlockData> blockData;
extractThermalBlockData(MoFEM::Interface &m_field,
std::vector<const CubitMeshSets *> meshset_vec_ptr,
Sev sev) {
for (auto m : meshset_vec_ptr) {
MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Thermal Block") << *m;
std::vector<double> block_data;
CHKERR m->getAttributes(block_data);
if (block_data.size() < 1) {
SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
"Expected that block has two attributes");
}
auto get_block_ents = [&]() {
Range ents;
m_field.get_moab().get_entities_by_handle(m->meshset, ents, true);
return ents;
};
MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Thermal Block")
<< m->getName() << ": conductivity = " << block_data[0];
blockData.push_back({block_data[0], get_block_ents()});
}
}
boost::shared_ptr<double> expansionPtr;
boost::shared_ptr<double> conductivityPtr;
boost::shared_ptr<double> capacityPtr;
};
pipeline.push_back(new OpMatFluidBlocks(
blockedParamsPtr->getConductivityPtr(), mField, sev,
// Get blockset using regular expression
mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
(boost::format("%s(.*)") % block_thermal_name).str()
))
));
}
//! [Run problem]
CHKERR setupProblem();
CHKERR createCommonData();
CHKERR bC();
CHKERR OPs();
CHKERR tsSolve();
}
//! [Run problem]
//! [Set up problem]
Simple *simple = mField.getInterface<Simple>();
// Add field
// Mechanical fields
CHKERR simple->addDomainField("U", H1, base, SPACE_DIM);
CHKERR simple->addBoundaryField("U", H1, base, SPACE_DIM);
// Temerature
const auto flux_space = (SPACE_DIM == 2) ? HCURL : HDIV;
CHKERR simple->addDomainField("H", L2, AINSWORTH_LEGENDRE_BASE, 1);
CHKERR simple->addDomainField("FLUX", flux_space, DEMKOWICZ_JACOBI_BASE, 1);
CHKERR simple->addBoundaryField("FLUX", flux_space, DEMKOWICZ_JACOBI_BASE, 1);
int order = 2.;
CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
CHKERR simple->setFieldOrder("U", order);
CHKERR simple->setFieldOrder("H", order - 1);
CHKERR simple->setFieldOrder("FLUX", order);
CHKERR simple->setUp();
}
//! [Set up problem]
//! [Create common data]
auto get_command_line_parameters = [&]() {
CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-young_modulus",
&default_young_modulus, PETSC_NULL);
CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-poisson_ratio",
&default_poisson_ratio, PETSC_NULL);
CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-conductivity",
&default_conductivity, PETSC_NULL);
CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-fluid_density",
&fluid_density, PETSC_NULL);
MOFEM_LOG("Seepage", Sev::inform)
<< "Default Young modulus " << default_young_modulus;
MOFEM_LOG("Seepage", Sev::inform)
<< "Default Poisson ratio " << default_poisson_ratio;
MOFEM_LOG("Seepage", Sev::inform)
<< "Default conductivity " << default_conductivity;
MOFEM_LOG("Seepage", Sev::inform) << "Fluid denisty " << fluid_density;
int coords_dim = SPACE_DIM;
CHKERR PetscOptionsGetRealArray(NULL, NULL, "-field_eval_coords",
fieldEvalCoords.data(), &coords_dim,
};
CHKERR get_command_line_parameters();
}
//! [Create common data]
//! [Boundary condition]
MOFEM_LOG("SYNC", Sev::noisy) << "bC";
MOFEM_LOG_SEVERITY_SYNC(mField.get_comm(), Sev::noisy);
auto simple = mField.getInterface<Simple>();
auto bc_mng = mField.getInterface<BcManager>();
CHKERR bc_mng->pushMarkDOFsOnEntities<DisplacementCubitBcData>(
simple->getProblemName(), "U");
CHKERR bc_mng->pushMarkDOFsOnEntities<HeatFluxCubitBcData>(
simple->getProblemName(), "FLUX", false);
auto get_skin = [&]() {
Range body_ents;
CHKERR mField.get_moab().get_entities_by_dimension(0, SPACE_DIM, body_ents);
Skinner skin(&mField.get_moab());
Range skin_ents;
CHKERR skin.find_skin(0, body_ents, false, skin_ents);
return skin_ents;
};
auto filter_flux_blocks = [&](auto skin) {
auto remove_cubit_blocks = [&](auto c) {
for (auto m :
mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(c)
) {
Range ents;
CHKERR mField.get_moab().get_entities_by_dimension(
m->getMeshset(), SPACE_DIM - 1, ents, true);
skin = subtract(skin, ents);
}
};
auto remove_named_blocks = [&](auto n) {
for (auto m : mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(
std::regex(
(boost::format("%s(.*)") % n).str()
))
) {
Range ents;
CHKERR mField.get_moab().get_entities_by_dimension(
m->getMeshset(), SPACE_DIM - 1, ents, true);
skin = subtract(skin, ents);
}
};
CHK_THROW_MESSAGE(remove_cubit_blocks(NODESET | TEMPERATURESET),
"remove_cubit_blocks");
CHK_THROW_MESSAGE(remove_cubit_blocks(SIDESET | HEATFLUXSET),
"remove_cubit_blocks");
CHK_THROW_MESSAGE(remove_named_blocks("PRESSURE"), "remove_named_blocks");
CHK_THROW_MESSAGE(remove_named_blocks("FLUIDFLUX"), "remove_named_blocks");
return skin;
};
auto filter_true_skin = [&](auto skin) {
Range boundary_ents;
ParallelComm *pcomm =
ParallelComm::get_pcomm(&mField.get_moab(), MYPCOMM_INDEX);
CHKERR pcomm->filter_pstatus(skin, PSTATUS_SHARED | PSTATUS_MULTISHARED,
PSTATUS_NOT, -1, &boundary_ents);
return boundary_ents;
};
auto remove_flux_ents = filter_true_skin(filter_flux_blocks(get_skin()));
CHKERR mField.getInterface<CommInterface>()->synchroniseEntities(
remove_flux_ents);
MOFEM_LOG("SYNC", Sev::noisy) << remove_flux_ents << endl;
MOFEM_LOG_SEVERITY_SYNC(mField.get_comm(), Sev::noisy);
#ifndef NDEBUG
mField.get_moab(),
(boost::format("flux_remove_%d.vtk") % mField.get_comm_rank()).str(),
remove_flux_ents);
#endif
CHKERR mField.getInterface<ProblemsManager>()->removeDofsOnEntities(
simple->getProblemName(), "FLUX", remove_flux_ents);
}
//! [Boundary condition]
//! [Push operators to pipeline]
auto pipeline_mng = mField.getInterface<PipelineManager>();
auto simple = mField.getInterface<Simple>();
auto bc_mng = mField.getInterface<BcManager>();
auto boundary_marker =
bc_mng->getMergedBlocksMarker(vector<string>{"FLUIDFLUX"});
auto u_grad_ptr = boost::make_shared<MatrixDouble>();
auto dot_u_grad_ptr = boost::make_shared<MatrixDouble>();
auto trace_dot_u_grad_ptr = boost::make_shared<VectorDouble>();
auto h_ptr = boost::make_shared<VectorDouble>();
auto dot_h_ptr = boost::make_shared<VectorDouble>();
auto flux_ptr = boost::make_shared<MatrixDouble>();
auto div_flux_ptr = boost::make_shared<VectorDouble>();
auto strain_ptr = boost::make_shared<MatrixDouble>();
auto stress_ptr = boost::make_shared<MatrixDouble>();
auto time_scale = boost::make_shared<TimeScale>();
auto block_params = boost::make_shared<BlockedParameters>();
auto mDPtr = block_params->getDPtr();
auto conductivity_ptr = block_params->getConductivityPtr();
auto integration_rule = [](int, int, int approx_order) {
return 2 * approx_order;
};
CHKERR pipeline_mng->setDomainRhsIntegrationRule(integration_rule);
CHKERR pipeline_mng->setDomainLhsIntegrationRule(integration_rule);
CHKERR pipeline_mng->setBoundaryLhsIntegrationRule(integration_rule);
CHKERR pipeline_mng->setBoundaryRhsIntegrationRule(integration_rule);
auto add_domain_base_ops = [&](auto &pip) {
CHKERR addMatBlockOps(pip, "MAT_ELASTIC", "MAT_FLUID", block_params,
Sev::inform);
"U", u_grad_ptr));
pip.push_back(new OpSymmetrizeTensor<SPACE_DIM>(u_grad_ptr, strain_ptr));
"U", dot_u_grad_ptr));
pip.push_back(new OpCalculateTraceFromMat<SPACE_DIM>(dot_u_grad_ptr,
trace_dot_u_grad_ptr));
pip.push_back(new OpCalculateScalarFieldValues("H", h_ptr));
pip.push_back(new OpCalculateScalarFieldValuesDot("H", dot_h_ptr));
pip.push_back(new OpCalculateHVecVectorField<3>("FLUX", flux_ptr));
"FLUX", div_flux_ptr));
};
auto add_domain_ops_lhs_mechanical = [&](auto &pip) {
pip.push_back(new OpKCauchy("U", "U", mDPtr));
pip.push_back(new OpBaseDivU(
"H", "U",
[](const double, const double, const double) { return -9.81; }, true,
true));
};
auto add_domain_ops_rhs_mechanical = [&](auto &pip) {
pip, mField, "U", {time_scale}, "BODY_FORCE", Sev::inform);
// Calculate internal forece
"U", strain_ptr, stress_ptr, mDPtr));
pip.push_back(new OpInternalForceCauchy("U", stress_ptr));
pip.push_back(
};
auto add_domain_ops_lhs_seepage = [&](auto &pip, auto &fe) {
auto resistance = [conductivity_ptr](const double, const double,
const double) {
return (1. / *(conductivity_ptr));
};
auto unity = []() constexpr { return -1; };
pip.push_back(new OpHdivHdiv("FLUX", "FLUX", resistance));
pip.push_back(new OpHdivQ("FLUX", "H", unity, true));
auto op_base_div_u = new OpBaseDivU(
"H", "U", [](double, double, double) constexpr { return -1; }, false,
false);
op_base_div_u->feScalingFun = [](const FEMethod *fe_ptr) {
return fe_ptr->ts_a;
};
pip.push_back(op_base_div_u);
};
auto add_domain_ops_rhs_seepage = [&](auto &pip) {
auto resistance = [conductivity_ptr](double, double, double) {
return (1. / (*conductivity_ptr));
};
auto minus_one = [](double, double, double) constexpr { return -1; };
pip.push_back(new OpHdivFlux("FLUX", flux_ptr, resistance));
pip.push_back(new OpHDivH("FLUX", h_ptr, minus_one));
pip.push_back(new OpBaseDotH("H", trace_dot_u_grad_ptr, minus_one));
pip.push_back(new OpBaseDivFlux("H", div_flux_ptr, minus_one));
};
auto add_boundary_rhs_ops = [&](auto &pip) {
pip.push_back(new OpSetBc("FLUX", true, boundary_marker));
pipeline_mng->getOpBoundaryRhsPipeline(), mField, "U", {time_scale},
"FORCE", Sev::inform);
pipeline_mng->getOpBoundaryRhsPipeline(), mField, "FLUX", {time_scale},
"PRESSURE", Sev::inform);
pip.push_back(new OpUnSetBc("FLUX"));
auto mat_flux_ptr = boost::make_shared<MatrixDouble>();
pip.push_back(
new OpCalculateHVecVectorField<3, SPACE_DIM>("FLUX", mat_flux_ptr));
mField, pip, simple->getProblemName(), "FLUX", mat_flux_ptr,
{time_scale});
};
auto add_boundary_lhs_ops = [&](auto &pip) {
mField, pip, simple->getProblemName(), "FLUX");
};
// LHS
CHKERR add_domain_base_ops(pipeline_mng->getOpDomainLhsPipeline());
CHKERR add_domain_ops_lhs_mechanical(pipeline_mng->getOpDomainLhsPipeline());
CHKERR add_domain_ops_lhs_seepage(pipeline_mng->getOpDomainLhsPipeline(),
pipeline_mng->getDomainLhsFE());
// RHS
CHKERR add_domain_base_ops(pipeline_mng->getOpDomainRhsPipeline());
CHKERR add_domain_ops_rhs_mechanical(pipeline_mng->getOpDomainRhsPipeline());
CHKERR add_domain_ops_rhs_seepage(pipeline_mng->getOpDomainRhsPipeline());
CHKERR add_boundary_rhs_ops(pipeline_mng->getOpBoundaryRhsPipeline());
CHKERR add_boundary_lhs_ops(pipeline_mng->getOpBoundaryLhsPipeline());
}
//! [Push operators to pipeline]
//! [Solve]
auto simple = mField.getInterface<Simple>();
PipelineManager *pipeline_mng = mField.getInterface<PipelineManager>();
ISManager *is_manager = mField.getInterface<ISManager>();
auto dm = simple->getDM();
auto solver = pipeline_mng->createTSIM();
auto snes_ctx_ptr = getDMSnesCtx(dm);
auto set_section_monitor = [&](auto solver) {
SNES snes;
CHKERR TSGetSNES(solver, &snes);
CHKERR SNESMonitorSet(snes,
(MoFEMErrorCode(*)(SNES, PetscInt, PetscReal,
(void *)(snes_ctx_ptr.get()), nullptr);
};
auto create_post_process_element = [&]() {
auto post_proc_fe = boost::make_shared<PostProcEle>(mField);
auto block_params = boost::make_shared<BlockedParameters>();
auto mDPtr = block_params->getDPtr();
CHKERR addMatBlockOps(post_proc_fe->getOpPtrVector(), "MAT_ELASTIC",
"MAT_FLUID", block_params, Sev::verbose);
post_proc_fe->getOpPtrVector(), {H1, HDIV});
auto mat_grad_ptr = boost::make_shared<MatrixDouble>();
auto mat_strain_ptr = boost::make_shared<MatrixDouble>();
auto mat_stress_ptr = boost::make_shared<MatrixDouble>();
auto h_ptr = boost::make_shared<VectorDouble>();
auto mat_flux_ptr = boost::make_shared<MatrixDouble>();
post_proc_fe->getOpPtrVector().push_back(
new OpCalculateScalarFieldValues("H", h_ptr));
post_proc_fe->getOpPtrVector().push_back(
new OpCalculateHVecVectorField<3, SPACE_DIM>("FLUX", mat_flux_ptr));
auto u_ptr = boost::make_shared<MatrixDouble>();
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getOpPtrVector().push_back(
mat_grad_ptr));
post_proc_fe->getOpPtrVector().push_back(
new OpSymmetrizeTensor<SPACE_DIM>(mat_grad_ptr, mat_strain_ptr));
post_proc_fe->getOpPtrVector().push_back(
"U", mat_strain_ptr, mat_stress_ptr, mDPtr));
post_proc_fe->getOpPtrVector().push_back(
new OpPPMap(
post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
{{"H", h_ptr}},
{{"U", u_ptr}, {"FLUX", mat_flux_ptr}},
{},
{{"STRAIN", mat_strain_ptr}, {"STRESS", mat_stress_ptr}}
)
);
return post_proc_fe;
};
auto create_creaction_fe = [&]() {
auto fe_ptr = boost::make_shared<DomainEle>(mField);
fe_ptr->getRuleHook = [](int, int, int o) { return 2 * o; };
auto &pip = fe_ptr->getOpPtrVector();
auto block_params = boost::make_shared<BlockedParameters>();
auto mDPtr = block_params->getDPtr();
CHKERR addMatBlockOps(pip, "MAT_ELASTIC", "MAT_FLUID", block_params,
Sev::verbose);
auto u_grad_ptr = boost::make_shared<MatrixDouble>();
auto strain_ptr = boost::make_shared<MatrixDouble>();
auto stress_ptr = boost::make_shared<MatrixDouble>();
"U", u_grad_ptr));
pip.push_back(new OpSymmetrizeTensor<SPACE_DIM>(u_grad_ptr, strain_ptr));
// Calculate internal forece
"U", strain_ptr, stress_ptr, mDPtr));
pip.push_back(new OpInternalForceCauchy("U", stress_ptr));
fe_ptr->postProcessHook =
return fe_ptr;
};
auto monitor_ptr = boost::make_shared<FEMethod>();
auto post_proc_fe = create_post_process_element();
auto res = createDMVector(dm);
auto rections_fe = create_creaction_fe();
auto set_time_monitor = [&](auto dm, auto solver) {
monitor_ptr->preProcessHook = [&]() {
CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(),
post_proc_fe,
monitor_ptr->getCacheWeakPtr());
CHKERR post_proc_fe->writeFile(
"out_" + boost::lexical_cast<std::string>(monitor_ptr->ts_step) +
".h5m");
rections_fe->f = res;
CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(),
rections_fe,
monitor_ptr->getCacheWeakPtr());
double nrm;
CHKERR VecNorm(res, NORM_2, &nrm);
MOFEM_LOG("Seepage", Sev::verbose)
<< "Residual norm " << nrm << " at time step "
<< monitor_ptr->ts_step;
if (doEvalField) {
auto scalar_field_ptr = boost::make_shared<VectorDouble>();
auto vector_field_ptr = boost::make_shared<MatrixDouble>();
auto tensor_field_ptr = boost::make_shared<MatrixDouble>();
auto field_eval_data = mField.getInterface<FieldEvaluatorInterface>()
->getData<DomainEle>();
if constexpr (SPACE_DIM == 3) {
CHKERR mField.getInterface<FieldEvaluatorInterface>()->buildTree3D(
field_eval_data, simple->getDomainFEName());
} else {
CHKERR mField.getInterface<FieldEvaluatorInterface>()->buildTree2D(
field_eval_data, simple->getDomainFEName());
}
field_eval_data->setEvalPoints(fieldEvalCoords.data(), 1);
auto no_rule = [](int, int, int) { return -1; };
auto field_eval_ptr = field_eval_data->feMethodPtr.lock();
field_eval_ptr->getRuleHook = no_rule;
auto u_grad_ptr = boost::make_shared<MatrixDouble>();
auto strain_ptr = boost::make_shared<MatrixDouble>();
auto stress_ptr = boost::make_shared<MatrixDouble>();
auto h_ptr = boost::make_shared<VectorDouble>();
auto block_params = boost::make_shared<BlockedParameters>();
auto mDPtr = block_params->getDPtr();
CHKERR addMatBlockOps(field_eval_ptr->getOpPtrVector(), "MAT_ELASTIC",
"MAT_FLUID", block_params, Sev::noisy);
field_eval_ptr->getOpPtrVector(), {H1, HDIV});
field_eval_ptr->getOpPtrVector().push_back(
"U", u_grad_ptr));
field_eval_ptr->getOpPtrVector().push_back(
new OpSymmetrizeTensor<SPACE_DIM>(u_grad_ptr, strain_ptr));
field_eval_ptr->getOpPtrVector().push_back(
new OpCalculateScalarFieldValues("H", h_ptr));
field_eval_ptr->getOpPtrVector().push_back(
"U", strain_ptr, stress_ptr, mDPtr));
if constexpr (SPACE_DIM == 3) {
CHKERR mField.getInterface<FieldEvaluatorInterface>()
->evalFEAtThePoint3D(
fieldEvalCoords.data(), 1e-12, simple->getProblemName(),
simple->getDomainFEName(), field_eval_data,
mField.get_comm_rank(), mField.get_comm_rank(), nullptr,
} else {
CHKERR mField.getInterface<FieldEvaluatorInterface>()
->evalFEAtThePoint2D(
fieldEvalCoords.data(), 1e-12, simple->getProblemName(),
simple->getDomainFEName(), field_eval_data,
mField.get_comm_rank(), mField.get_comm_rank(), nullptr,
}
MOFEM_LOG("SeepageSync", Sev::inform)
<< "Eval point hydrostatic hight: " << *h_ptr;
MOFEM_LOG("SeepageSync", Sev::inform)
<< "Eval point skeleton stress pressure: " << *stress_ptr;
MOFEM_LOG_SEVERITY_SYNC(mField.get_comm(), Sev::inform);
}
};
auto null = boost::shared_ptr<FEMethod>();
CHKERR DMMoFEMTSSetMonitor(dm, solver, simple->getDomainFEName(), null,
monitor_ptr, null);
};
auto set_fieldsplit_preconditioner = [&](auto solver) {
SNES snes;
CHKERR TSGetSNES(solver, &snes);
KSP ksp;
CHKERR SNESGetKSP(snes, &ksp);
PC pc;
CHKERR KSPGetPC(ksp, &pc);
PetscBool is_pcfs = PETSC_FALSE;
PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &is_pcfs);
// Setup fieldsplit (block) solver - optional: yes/no
if (is_pcfs == PETSC_TRUE) {
auto bc_mng = mField.getInterface<BcManager>();
auto is_mng = mField.getInterface<ISManager>();
auto name_prb = simple->getProblemName();
CHKERR is_mng->isCreateProblemFieldAndRank(name_prb, ROW, "U", 0,
SPACE_DIM, is_u);
CHKERR is_mng->isCreateProblemFieldAndRank(name_prb, ROW, "FLUX", 0, 0,
is_flux);
CHKERR is_mng->isCreateProblemFieldAndRank(name_prb, ROW, "H", 0, 0,
is_H);
IS is_tmp;
CHKERR ISExpand(is_H, is_flux, &is_tmp);
auto is_Flux = SmartPetscObj<IS>(is_tmp);
auto is_all_bc = bc_mng->getBlockIS(name_prb, "FLUIDFLUX", "FLUX", 0, 0);
int is_all_bc_size;
CHKERR ISGetSize(is_all_bc, &is_all_bc_size);
MOFEM_LOG("ThermoElastic", Sev::inform)
<< "Field split block size " << is_all_bc_size;
if (is_all_bc_size) {
IS is_tmp2;
CHKERR ISDifference(is_Flux, is_all_bc, &is_tmp2);
is_Flux = SmartPetscObj<IS>(is_tmp2);
CHKERR PCFieldSplitSetIS(pc, PETSC_NULL,
is_all_bc); // boundary block
}
CHKERR ISSort(is_u);
CHKERR ISSort(is_Flux);
CHKERR PCFieldSplitSetIS(pc, PETSC_NULL, is_Flux);
CHKERR PCFieldSplitSetIS(pc, PETSC_NULL, is_u);
}
};
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<TimeScale>();
auto get_bc_hook_rhs = [this, pre_proc_ptr, time_scale]() {
EssentialPreProc<DisplacementCubitBcData> hook(mField, pre_proc_ptr,
{time_scale}, false);
return hook;
};
auto get_post_proc_hook_rhs = [this, post_proc_rhs_ptr]() {
mField, post_proc_rhs_ptr, nullptr, Sev::verbose)();
mField, post_proc_rhs_ptr, 1.)();
};
auto get_post_proc_hook_lhs = [this, post_proc_lhs_ptr]() {
post_proc_lhs_ptr, 1.);
};
pre_proc_ptr->preProcessHook = get_bc_hook_rhs();
post_proc_rhs_ptr->postProcessHook = get_post_proc_hook_rhs;
post_proc_lhs_ptr->postProcessHook = get_post_proc_hook_lhs();
auto ts_ctx_ptr = getDMTsCtx(dm);
ts_ctx_ptr->getPreProcessIFunction().push_front(pre_proc_ptr);
ts_ctx_ptr->getPreProcessIJacobian().push_front(pre_proc_ptr);
ts_ctx_ptr->getPostProcessIFunction().push_back(post_proc_rhs_ptr);
ts_ctx_ptr->getPostProcessIJacobian().push_back(post_proc_lhs_ptr);
auto D = createDMVector(dm);
CHKERR TSSetSolution(solver, D);
CHKERR TSSetFromOptions(solver);
CHKERR set_section_monitor(solver);
CHKERR set_fieldsplit_preconditioner(solver);
CHKERR set_time_monitor(dm, solver);
CHKERR TSSetUp(solver);
CHKERR TSSolve(solver, NULL);
}
//! [Solve]
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);
// Add logging channel for example
auto core_log = logging::core::get();
core_log->add_sink(
LogManager::createSink(LogManager::getStrmWorld(), "Seepage"));
LogManager::setLog("Seepage");
MOFEM_LOG_TAG("Seepage", "seepage");
core_log->add_sink(
LogManager::createSink(LogManager::getStrmSync(), "SeepageSync"));
LogManager::setLog("SeepageSync");
MOFEM_LOG_TAG("SeepageSync", "seepage");
try {
//! [Register MoFEM discrete manager in PETSc]
DMType dm_name = "DMMOFEM";
//! [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 insterface
//! [Create MoFEM]
//! [Load mesh]
CHKERR simple->getOptions();
CHKERR simple->loadFile();
//! [Load mesh]
//! [Seepage]
Seepage ex(m_field);
CHKERR ex.runProblem();
//! [Seepage]
}
}
NOSPACE
@ NOSPACE
Definition: definitions.h:83
MoFEM::NaturalBC::Assembly::LinearForm
Definition: Natural.hpp:67
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
SIDESET
@ SIDESET
Definition: definitions.h:160
BoundaryEle
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
Definition: child_and_parent.cpp:39
save_range
auto save_range
Definition: seepage.cpp:170
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:228
MoFEM::EssentialPreProcReaction< DisplacementCubitBcData >
Specialization for DisplacementCubitBcData.
Definition: EssentialDisplacementCubitBcData.hpp:157
Seepage::createCommonData
MoFEMErrorCode createCommonData()
[Set up problem]
Definition: seepage.cpp:476
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
MOFEM_LOG_SEVERITY_SYNC
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
Definition: LogManager.hpp:352
SeepageOps.hpp
OpHdivHdiv
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 3, SPACE_DIM > OpHdivHdiv
[Linear elastic problem]
Definition: thermo_elastic.cpp:49
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
help
static char help[]
[Solve]
Definition: seepage.cpp:1108
young_modulus
double young_modulus
Young modulus.
Definition: plastic.cpp:121
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
MoFEM::FEMethod
structure for User Loop Methods on finite elements
Definition: LoopMethods.hpp:369
MoFEM::OpFluxRhsImpl
Definition: Natural.hpp:39
default_conductivity
double default_conductivity
Definition: seepage.cpp:164
SPACE_DIM
constexpr int SPACE_DIM
Definition: seepage.cpp:34
MoFEM::PipelineManager::ElementsAndOpsByDim
Definition: PipelineManager.hpp:38
EntData
EntitiesFieldData::EntData EntData
Definition: seepage.cpp:37
MoFEM::EssentialPostProcLhs< DisplacementCubitBcData >
Specialization for DisplacementCubitBcData.
Definition: EssentialDisplacementCubitBcData.hpp:134
MoFEM::NaturalBC::Assembly::BiLinearForm
Definition: Natural.hpp:74
Seepage::OPs
MoFEMErrorCode OPs()
[Boundary condition]
Definition: seepage.cpp:614
OpKCauchy
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradSymTensorGrad< 1, SPACE_DIM, SPACE_DIM, 0 > OpKCauchy
[Only used with Hooke equation (linear material model)]
Definition: seepage.cpp:50
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:609
L2
@ L2
field with C-1 continuity
Definition: definitions.h:88
OpHdivFlux
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesVector< 3, 3, 1 > OpHdivFlux
Integrating Rhs flux base (1/k) flux (FLUX)
Definition: seepage.cpp:109
MoFEM::EssentialBC
Essential boundary conditions.
Definition: Essential.hpp:111
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::OpCalculateVectorFieldValues
Get values at integration pts for tensor filed rank 1, i.e. vector field.
Definition: UserDataOperators.hpp:466
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
MoFEM::OpCalculateScalarFieldValuesDot
OpCalculateScalarFieldValuesFromPetscVecImpl< PetscData::CTX_SET_X_T > OpCalculateScalarFieldValuesDot
Definition: UserDataOperators.hpp:273
MoFEM::EssentialPostProcRhs< DisplacementCubitBcData >
Specialization for DisplacementCubitBcData.
Definition: EssentialDisplacementCubitBcData.hpp:113
OpKPiola
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradTensorGrad< 1, SPACE_DIM, SPACE_DIM, 1 > OpKPiola
[Only used for dynamics]
Definition: seepage.cpp:64
MoFEM::PETSC
@ PETSC
Definition: FormsIntegrators.hpp:105
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
MoFEM.hpp
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
bulk_modulus_K
double bulk_modulus_K
Definition: dynamic_first_order_con_law.cpp:96
MoFEM::DisplacementCubitBcData
Definition of the displacement bc data structure.
Definition: BCData.hpp:76
OpHdivHdiv
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 3, 3 > OpHdivHdiv
Integrate Lhs base of flux (1/k) base of flux (FLUX x FLUX)
Definition: seepage.cpp:87
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
Seepage::addMatBlockOps
MoFEMErrorCode addMatBlockOps(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, std::string block_elastic_name, std::string block_thermal_name, boost::shared_ptr< BlockedParameters > blockedParamsPtr, Sev sev)
Definition: seepage.cpp:221
OpBaseDotH
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesScalarField< 1 > OpBaseDotH
Integrate Rhs base of temperature time heat capacity times heat rate (T)
Definition: seepage.cpp:124
OpBase
OpBaseImpl< PETSC, EdgeEleOp > OpBase
Definition: radiation.cpp:29
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
OpInternalForceCauchy
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesSymTensor< 1, SPACE_DIM, SPACE_DIM > OpInternalForceCauchy
Definition: seepage.cpp:52
sdf.r
int r
Definition: sdf.py:8
Seepage::setupProblem
MoFEMErrorCode setupProblem()
[Run problem]
Definition: seepage.cpp:449
Seepage::bC
MoFEMErrorCode bC()
[Create common data]
Definition: seepage.cpp:513
order
constexpr int order
Definition: dg_projection.cpp:18
OpMixScalarTimesDiv
FormsIntegrators< DomainEleOp >::Assembly< A >::BiLinearForm< I >::OpMixScalarTimesDiv< SPACE_DIM, COORD_TYPE > OpMixScalarTimesDiv
Definition: free_surface.cpp:125
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
OpBoundaryMass
FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpBoundaryMass
[Only used with Hencky/nonlinear material]
Definition: seepage.cpp:71
ROW
@ ROW
Definition: definitions.h:136
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
c
const double c
speed of light (cm/ns)
Definition: initial_diffusion.cpp:39
scale
double scale
[Essential boundary conditions (Least square approach)]
Definition: seepage.cpp:160
OpInertiaForce
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesVector< 1, SPACE_DIM, 1 > OpInertiaForce
Definition: dynamic_first_order_con_law.cpp:63
MoFEM::PostProcBrokenMeshInMoab
Definition: PostProcBrokenMeshInMoabBase.hpp:667
NODESET
@ NODESET
Definition: definitions.h:159
MoFEM::OpBaseImpl
Definition: FormsIntegrators.hpp:178
OpHDivH
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpMixDivTimesU< 3, 1, 2 > OpHDivH
Integrate Rhs div flux base times temperature (T)
Definition: seepage.cpp:116
default_young_modulus
double default_young_modulus
Definition: seepage.cpp:162
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::createDMVector
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1099
MoFEM::EssentialPreProc< DisplacementCubitBcData >
Specialization for DisplacementCubitBcData.
Definition: EssentialDisplacementCubitBcData.hpp:91
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
OpInternalForcePiola
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, SPACE_DIM, SPACE_DIM > OpInternalForcePiola
Definition: seepage.cpp:66
MoFEM::ISManager
Section manager is used to create indexes and sections.
Definition: ISManager.hpp:23
MoFEM::BcManager
Simple interface for fast problem set-up.
Definition: BcManager.hpp:25
OpMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMass
[Only used with Hooke equation (linear material model)]
Definition: seepage.cpp:57
BoundaryEleOp
MoFEM::OpCalculateHVecVectorField
Get vector field for H-div approximation.
Definition: UserDataOperators.hpp:2115
MoFEM::get_temp_meshset_ptr
auto get_temp_meshset_ptr(moab::Interface &moab)
Create smart pointer to temporary meshset.
Definition: Templates.hpp:1886
MoFEM::FieldEvaluatorInterface
Field evaluator interface.
Definition: FieldEvaluator.hpp:21
simple
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
double
MoFEM::OpEssentialLhsImpl
Enforce essential constrains on lhs.
Definition: Essential.hpp:81
convert.type
type
Definition: convert.py:64
Seepage
Definition: seepage.cpp:179
MoFEM::FormsIntegrators::Assembly
Assembly methods.
Definition: FormsIntegrators.hpp:317
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
OpBaseDivFlux
OpBaseDotT OpBaseDivFlux
Integrate Rhs base of temperature times divergence of flux (T)
Definition: thermo_elastic.cpp:92
MoFEM::OpCalculateScalarFieldValues
Get value at integration points for scalar field.
Definition: UserDataOperators.hpp:82
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
poisson_ratio
double poisson_ratio
Poisson ratio.
Definition: plastic.cpp:122
MoFEM::GAUSS
@ GAUSS
Definition: FormsIntegrators.hpp:136
EXECUTABLE_DIMENSION
#define EXECUTABLE_DIMENSION
Definition: seepage.cpp:24
MatrixFunction.hpp
OpKCauchy
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradSymTensorGrad< 1, SPACE_DIM, SPACE_DIM, 0 > OpKCauchy
[Linear elastic problem]
Definition: thermo_elastic.cpp:36
doEvalField
PetscBool doEvalField
Definition: incompressible_elasticity.cpp:41
Seepage::tsSolve
MoFEMErrorCode tsSolve()
[Push operators to pipeline]
Definition: seepage.cpp:793
MoFEM::LogManager::SeverityLevel
SeverityLevel
Severity levels.
Definition: LogManager.hpp:33
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
MoFEM::AddFluxToRhsPipelineImpl
Definition: Natural.hpp:46
size_symm
constexpr auto size_symm
Definition: plastic.cpp:42
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MoFEM::OpTensorTimesSymmetricTensor
Calculate .
Definition: UserDataOperators.hpp:1883
t_kd
constexpr auto t_kd
Definition: free_surface.cpp:137
BiLinearForm
OpBoundaryInternal
FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesVector< 1, SPACE_DIM, 1 > OpBoundaryInternal
Definition: seepage.cpp:75
default_poisson_ratio
double default_poisson_ratio
Definition: seepage.cpp:163
OpInternalForceCauchy
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesSymTensor< 1, SPACE_DIM, SPACE_DIM > OpInternalForceCauchy
Definition: thermo_elastic.cpp:40
OpGradTimesTensor
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor
Definition: operators_tests.cpp:48
FTensor::Index< 'i', SPACE_DIM >
MoFEM::AddHOOps
Add operators pushing bases from local to physical configuration.
Definition: HODataOperators.hpp:413
MoFEM::OpCalculateVectorFieldGradient
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Definition: UserDataOperators.hpp:1535
convert.n
n
Definition: convert.py:82
fluid_density
double fluid_density
Definition: seepage.cpp:165
integration_rule
auto integration_rule
Definition: free_surface.cpp:185
Seepage::runProblem
MoFEMErrorCode runProblem()
[Run problem]
Definition: seepage.cpp:437
MoFEM::OpUnSetBc
Definition: FormsIntegrators.hpp:49
Range
DomainEleOp
OpBaseDivFlux
OpBaseDotH OpBaseDivFlux
Integrate Rhs base of temperature times divergent of flux (T)
Definition: seepage.cpp:130
MoFEM::PetscOptionsGetRealArray
PetscErrorCode PetscOptionsGetRealArray(PetscOptions *, const char pre[], const char name[], PetscReal dval[], PetscInt *nmax, PetscBool *set)
Definition: DeprecatedPetsc.hpp:192
MoFEM::CoreTmp< 0 >::Initialize
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
MOFEM_TAG_AND_LOG
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
Definition: LogManager.hpp:362
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::NaturalBC::Assembly
Assembly methods.
Definition: Natural.hpp:65
MoFEM::OpEssentialRhsImpl
Enforce essential constrains on rhs.
Definition: Essential.hpp:65
OpHdivQ
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMixDivTimesScalar< SPACE_DIM > OpHdivQ
Integrate Lhs div of base of flux time base of temperature (FLUX x T) and transpose of it,...
Definition: seepage.cpp:95
shear_modulus_G
double shear_modulus_G
Definition: dynamic_first_order_con_law.cpp:97
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
MoFEM::OpSymmetrizeTensor
Definition: UserDataOperators.hpp:1948
MoFEM::CommInterface
Managing BitRefLevels.
Definition: CommInterface.hpp:21
IntegrationRules.hpp
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
approx_order
int approx_order
Definition: test_broken_space.cpp:50
MoFEM::OpCalculateTraceFromMat
Calculates the trace of an input matrix.
Definition: UserDataOperators.hpp:3319
MoFEM::DMMoFEMTSSetMonitor
PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
Set Monitor To TS solver.
Definition: DMMoFEM.cpp:1056
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::OpSetBc
Set indices on entities on finite element.
Definition: FormsIntegrators.hpp:38
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
TEMPERATURESET
@ TEMPERATURESET
Definition: definitions.h:168
MoFEM::MeshsetsManager
Interface for managing meshsets containing materials and boundary conditions.
Definition: MeshsetsManager.hpp:104
MoFEM::BlockData
Definition: MeshsetsManager.cpp:755
MoFEM::OpCalculateHdivVectorDivergence
Calculate divergence of vector field.
Definition: UserDataOperators.hpp:2212
MoFEM::MoFEMSNESMonitorFields
MoFEMErrorCode MoFEMSNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *snes_ctx)
Sens monitor printing residual field by field.
Definition: SnesCtx.cpp:232
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
MoFEM::EssentialBC::Assembly
Assembly methods.
Definition: Essential.hpp:119
ReactionDiffusionEquation::D
const double D
diffusivity
Definition: reaction_diffusion.cpp:20
OpCapacity
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, 1 > OpCapacity
Integrate Lhs base of temperature times (heat capacity) times base of temperature (T x T)
Definition: thermo_elastic.cpp:65
MoFEM::HeatFluxCubitBcData
Definition of the heat flux bc data structure.
Definition: BCData.hpp:427
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
main
int main(int argc, char *argv[])
Definition: seepage.cpp:1110
FTensor::Kronecker_Delta_symmetric
Kronecker Delta class symmetric.
Definition: Kronecker_Delta.hpp:49
MoFEM::getDMSnesCtx
auto getDMSnesCtx(DM dm)
Get SNES context data structure used by DM.
Definition: DMMoFEM.hpp:1127
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
HEATFLUXSET
@ HEATFLUXSET
Definition: definitions.h:169
QUIET
@ QUIET
Definition: definitions.h:221
MoFEM::PetscOptionsGetScalar
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:162
DomainEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
Definition: child_and_parent.cpp:34
MoFEM::MeshsetsManager::getCubitMeshsetPtr
MoFEMErrorCode getCubitMeshsetPtr(const int ms_id, const CubitBCType cubit_bc_type, const CubitMeshSets **cubit_meshset_ptr) const
get cubit meshset
Definition: MeshsetsManager.cpp:578
MoFEM::SmartPetscObj< VecScatter >
OpBoundaryVec
FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesVector< 1, SPACE_DIM, 0 > OpBoundaryVec
Definition: seepage.cpp:73
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
MoFEM::DMoFEMLoopFiniteElements
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
MF_EXIST
@ MF_EXIST
Definition: definitions.h:113
convert.int
int
Definition: convert.py:64
MoFEM::OpCalculateVectorFieldGradientDot
Get field gradients time derivative at integration pts for scalar filed rank 0, i....
Definition: UserDataOperators.hpp:1550
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
HDIV
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
MoFEM::getDMTsCtx
auto getDMTsCtx(DM dm)
Get TS context data structure used by DM.
Definition: DMMoFEM.hpp:1141
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21
OpBaseDivU
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMixScalarTimesDiv< SPACE_DIM > OpBaseDivU
[Essential boundary conditions]
Definition: seepage.cpp:79
SeepageOps::OpDomainRhsHydrostaticStress
Definition: SeepageOps.hpp:22
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698
OpHdivFlux
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesVector< 3, SPACE_DIM, 1 > OpHdivFlux
Integrating Rhs flux base (1/k) flux (FLUX)
Definition: thermo_elastic.cpp:71
PlasticOps::addMatBlockOps
MoFEMErrorCode addMatBlockOps(MoFEM::Interface &m_field, std::string block_name, Pip &pip, boost::shared_ptr< MatrixDouble > mat_D_Ptr, boost::shared_ptr< CommonData::BlockParams > mat_params_ptr, double scale_value, Sev sev)
Definition: PlasticOps.hpp:248