v0.15.0
Loading...
Searching...
No Matches
mofem/tutorials/adv-2/thermo_elastic.cpp

Thermo plasticity

/**
* \file thermo_elastic.cpp
* \example mofem/tutorials/adv-2/thermo_elastic.cpp
*
* Thermo plasticity
*
*/
#ifndef EXECUTABLE_DIMENSION
#define EXECUTABLE_DIMENSION 3
#endif
#ifndef FINITE_DEFORMATION_FLAG
#define FINITE_DEFORMATION_FLAG true
#endif
#include <MoFEM.hpp>
using namespace MoFEM;
constexpr int SPACE_DIM =
EXECUTABLE_DIMENSION; //< Space dimension of problem, mesh
constexpr bool IS_LARGE_STRAINS =
FINITE_DEFORMATION_FLAG; //< Flag to turn off/on geometric nonlinearities
using DomainEleOp = DomainEle::UserDataOperator;
using BoundaryEle =
using BoundaryEleOp = BoundaryEle::UserDataOperator;
constexpr AssemblyType AT = AssemblyType::PETSC; //< selected assembly type
constexpr IntegrationType IT =
IntegrationType::GAUSS; //< selected integration type
//! [Linear elastic problem]
#include <HenckyOps.hpp> // Include Hencky operators
using namespace HenckyOps;
//! [Linear elastic problem]
// Include finite deformation operators
#include <FiniteThermalOps.hpp> // Include operators for finite strain diffusion problem
using namespace FiniteThermalOps;
//! [Thermal problem]
#include <ThermalOps.hpp> // Include operators for thermal problem which are agnostic to small/large strains
using namespace ThermalOps;
//! [Thermal problem]
//! [Body and heat source]
using OpBodyForce =
DomainNaturalBCRhs::OpFlux<NaturalMeshsetType<BLOCKSET>, 1, SPACE_DIM>;
using OpHeatSource =
DomainNaturalBCRhs::OpFlux<NaturalMeshsetType<BLOCKSET>, 1, 1>;
//! [Body and heat source]
//! [Natural boundary conditions]
using OpForce = BoundaryNaturalBC::OpFlux<NaturalForceMeshsets, 1, SPACE_DIM>;
//! [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)]
//! [Default input parameters]
double default_poisson_ratio = 0.25;
double default_ref_temp = 0.0;
double default_init_temp = 0.0;
1; // Force / (time temperature ) or Power /
// (length temperature) // Time unit is hour. force unit kN
double default_heat_capacity = 1; // length^2/(time^2 temperature) // length is
// millimeter time is hour
int order_temp = 2; //< default approximation order for the temperature field
int order_flux = 3; //< default approximation order for heat flux field
int order_disp = 3; //< default approximation order for the displacement field
int atom_test = 0;
int save_every = 1; //< Save every n-th step
PetscBool do_output_domain;
PetscBool do_output_skin;
//! [Default input parameters]
#include <ThermoElasticOps.hpp> //< additional coupling operators
using namespace ThermoElasticOps; //< name space of coupling operators
DomainNaturalBCRhs::OpFlux<SetTargetTemperature, 1, 1>;
DomainNaturalBCLhs::OpFlux<SetTargetTemperature, 1, 1>;
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);
};
private:
PetscBool doEvalField;
std::array<double, 3> fieldEvalCoords = {0.0, 0.0, 0.0};
boost::shared_ptr<FieldEvaluatorInterface::SetPtsData> fieldEvalData;
boost::shared_ptr<VectorDouble> tempFieldPtr;
boost::shared_ptr<MatrixDouble> fluxFieldPtr;
boost::shared_ptr<MatrixDouble> dispFieldPtr;
boost::shared_ptr<MatrixDouble> dispGradPtr;
boost::shared_ptr<MatrixDouble> strainFieldPtr;
boost::shared_ptr<MatrixDouble> stressFieldPtr;
MoFEMErrorCode setupProblem(); ///< add fields
getCommandLineParameters(); //< read parameters from command line
MoFEMErrorCode bC(); //< read boundary conditions
MoFEMErrorCode OPs(); //< add operators to pipeline
MoFEMErrorCode tsSolve(); //< time solver
struct BlockedThermalParameters
: public boost::enable_shared_from_this<BlockedThermalParameters> {
double heatCapacity;
inline auto getHeatConductivityPtr() {
return boost::shared_ptr<double>(shared_from_this(), &heatConductivity);
}
inline auto getHeatCapacityPtr() {
return boost::shared_ptr<double>(shared_from_this(), &heatCapacity);
}
};
struct BlockedThermoElasticParameters
: public boost::enable_shared_from_this<BlockedThermoElasticParameters> {
double refTemp;
inline auto getCoeffExpansionPtr() {
return boost::shared_ptr<VectorDouble>(shared_from_this(),
}
inline auto getRefTempPtr() {
return boost::shared_ptr<double>(shared_from_this(), &refTemp);
}
};
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
std::string block_name,
boost::shared_ptr<BlockedThermalParameters> blockedParamsPtr, Sev sev);
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
std::string block_name,
boost::shared_ptr<BlockedThermoElasticParameters> blockedParamsPtr,
Sev sev);
template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
MoFEM::Interface &m_field,
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
std::string field_name,
boost::shared_ptr<HenckyOps::CommonData> elastic_common_ptr,
boost::shared_ptr<ThermoElasticProblem::BlockedThermalParameters>
thermal_common_ptr,
boost::shared_ptr<ThermoElasticProblem::BlockedThermoElasticParameters>
thermoelastic_common_ptr,
Sev sev) {
using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
A>::template LinearForm<I>;
auto vec_temp_ptr = boost::make_shared<VectorDouble>();
pip.push_back(new OpCalculateScalarFieldValues("T", vec_temp_ptr));
auto coeff_expansion_ptr = thermoelastic_common_ptr->getCoeffExpansionPtr();
auto ref_temp_ptr = thermoelastic_common_ptr->getRefTempPtr();
pip.push_back(
new typename H::template OpCalculateHenckyThermalStress<DIM, I, 0>(
"U", vec_temp_ptr, elastic_common_ptr, coeff_expansion_ptr,
ref_temp_ptr));
typename B::template OpGradTimesTensor<1, DIM, DIM>;
pip.push_back(new OpInternalForcePiola(
"U", elastic_common_ptr->getMatFirstPiolaStress()));
}
template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
MoFEM::Interface &m_field,
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
std::string field_name, std::string elastic_block_name,
std::string thermal_block_name, std::string thermoelastic_block_name,
Sev sev, double scale = 1) {
auto elastic_common_ptr = commonDataFactory<DIM, I, DomainEleOp>(
m_field, pip, field_name, elastic_block_name, sev, scale);
auto thermal_common_ptr = boost::make_shared<BlockedThermalParameters>();
CHKERR addMatThermalBlockOps(pip, thermal_block_name, thermal_common_ptr,
Sev::inform);
auto thermoelastic_common_ptr =
boost::make_shared<BlockedThermoElasticParameters>();
CHKERR addMatThermoElasticBlockOps(pip, thermoelastic_block_name,
thermoelastic_common_ptr, Sev::inform);
CHKERR opThermoElasticFactoryDomainRhs<DIM, A, I, DomainEleOp>(
m_field, pip, field_name, elastic_common_ptr, thermal_common_ptr,
thermoelastic_common_ptr, sev);
}
template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
MoFEM::Interface &m_field,
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
std::string field_name, std::string coupled_field_name,
boost::shared_ptr<HenckyOps::CommonData> elastic_common_ptr,
boost::shared_ptr<ThermoElasticProblem::BlockedThermalParameters>
thermal_common_ptr,
boost::shared_ptr<ThermoElasticProblem::BlockedThermoElasticParameters>
thermoelastic_common_ptr,
Sev sev) {
using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
A>::template BiLinearForm<I>;
using OpKPiola = typename B::template OpGradTensorGrad<1, DIM, DIM, 1>;
auto vec_temp_ptr = boost::make_shared<VectorDouble>();
pip.push_back(new OpCalculateScalarFieldValues("T", vec_temp_ptr));
auto coeff_expansion_ptr = thermoelastic_common_ptr->getCoeffExpansionPtr();
auto ref_temp_ptr = thermoelastic_common_ptr->getRefTempPtr();
pip.push_back(
new typename H::template OpCalculateHenckyThermalStress<DIM, I, 0>(
"U", vec_temp_ptr, elastic_common_ptr, coeff_expansion_ptr,
ref_temp_ptr));
pip.push_back(new typename H::template OpHenckyTangent<DIM, I, 0>(
field_name, elastic_common_ptr));
pip.push_back(new OpKPiola(field_name, field_name,
elastic_common_ptr->getMatTangent()));
pip.push_back(new typename H::template OpCalculateHenckyThermalStressdT<
field_name, coupled_field_name, elastic_common_ptr,
coeff_expansion_ptr));
}
template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
MoFEM::Interface &m_field,
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
std::string field_name, std::string coupled_field_name,
std::string elastic_block_name, std::string thermal_block_name,
std::string thermoelastic_block_name, Sev sev, double scale = 1) {
auto elastic_common_ptr = commonDataFactory<DIM, I, DomainEleOp>(
m_field, pip, field_name, elastic_block_name, sev, scale);
auto thermal_common_ptr = boost::make_shared<BlockedThermalParameters>();
CHKERR addMatThermalBlockOps(pip, thermal_block_name, thermal_common_ptr,
Sev::inform);
auto thermoelastic_common_ptr =
boost::make_shared<BlockedThermoElasticParameters>();
CHKERR addMatThermoElasticBlockOps(pip, thermoelastic_block_name,
thermoelastic_common_ptr, Sev::inform);
CHKERR opThermoElasticFactoryDomainLhs<DIM, A, I, DomainEleOp>(
m_field, pip, field_name, coupled_field_name, elastic_common_ptr,
thermal_common_ptr, thermoelastic_common_ptr, sev);
}
};
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
std::string block_name,
boost::shared_ptr<BlockedThermalParameters> blockedParamsPtr, Sev sev) {
struct OpMatThermalBlocks : public DomainEleOp {
OpMatThermalBlocks(boost::shared_ptr<double> conductivity_ptr,
boost::shared_ptr<double> capacity_ptr,
MoFEM::Interface &m_field, Sev sev,
std::vector<const CubitMeshSets *> meshset_vec_ptr)
conductivityPtr(conductivity_ptr), capacityPtr(capacity_ptr) {
CHK_THROW_MESSAGE(extractThermalBlockData(m_field, meshset_vec_ptr, sev),
"Cannot get data from thermal block");
}
MoFEMErrorCode doWork(int side, EntityType type,
for (auto &b : blockData) {
if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
*conductivityPtr = b.conductivity;
*capacityPtr = b.capacity;
}
}
*conductivityPtr = default_heat_conductivity;
*capacityPtr = default_heat_capacity;
}
private:
struct BlockData {
double conductivity;
double capacity;
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() < 2) {
SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
"Expected that block has at least 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]
<< " capacity = " << block_data[1];
blockData.push_back({block_data[0], block_data[1], get_block_ents()});
}
}
boost::shared_ptr<double> conductivityPtr;
boost::shared_ptr<double> capacityPtr;
};
pipeline.push_back(new OpMatThermalBlocks(
blockedParamsPtr->getHeatConductivityPtr(),
blockedParamsPtr->getHeatCapacityPtr(), mField, sev,
// Get blockset using regular expression
(boost::format("%s(.*)") % block_name).str()
))
));
}
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
std::string block_name,
boost::shared_ptr<BlockedThermoElasticParameters> blockedParamsPtr,
Sev sev) {
struct OpMatThermoElasticBlocks : public DomainEleOp {
OpMatThermoElasticBlocks(boost::shared_ptr<VectorDouble> expansion_ptr,
boost::shared_ptr<double> ref_temp_ptr,
MoFEM::Interface &m_field, Sev sev,
std::vector<const CubitMeshSets *> meshset_vec_ptr)
expansionPtr(expansion_ptr), refTempPtr(ref_temp_ptr) {
extractThermoElasticBlockData(m_field, meshset_vec_ptr, sev),
"Cannot get data from thermoelastic block");
}
MoFEMErrorCode doWork(int side, EntityType type,
for (auto &b : blockData) {
if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
*expansionPtr = b.expansion;
*refTempPtr = b.ref_temp;
}
}
*refTempPtr = default_ref_temp;
}
private:
struct BlockData {
double ref_temp;
VectorDouble expansion;
Range blockEnts;
};
std::vector<BlockData> blockData;
MoFEMErrorCode extractThermoElasticBlockData(
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 Thermoelastic 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 at least two attributes");
}
auto get_block_ents = [&]() {
Range ents;
m_field.get_moab().get_entities_by_handle(m->meshset, ents, true);
return ents;
};
auto get_expansion = [&]() {
VectorDouble expansion(SPACE_DIM, block_data[1]);
if (block_data.size() > 2) {
expansion[1] = block_data[2];
}
if (SPACE_DIM == 3 && block_data.size() > 3) {
expansion[2] = block_data[3];
}
return expansion;
};
auto coeff_exp_vec = get_expansion();
MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Thermoelastic Block")
<< " ref_temp = " << block_data[0]
<< " expansion = " << coeff_exp_vec;
blockData.push_back({block_data[0], coeff_exp_vec, get_block_ents()});
}
}
boost::shared_ptr<VectorDouble> expansionPtr;
boost::shared_ptr<double> refTempPtr;
};
pipeline.push_back(new OpMatThermoElasticBlocks(
blockedParamsPtr->getCoeffExpansionPtr(),
blockedParamsPtr->getRefTempPtr(), mField, sev,
// Get blockset using regular expression
(boost::format("%s(.*)") % block_name).str()
))
));
}
//! [Run problem]
}
//! [Run problem]
//! [Get command line parameters]
auto get_command_line_parameters = [&]() {
CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-order", &order_temp,
PETSC_NULLPTR);
CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-order_temp", &order_temp,
PETSC_NULLPTR);
CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-order_flux", &order_flux,
PETSC_NULLPTR);
CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-order_disp", &order_disp,
PETSC_NULLPTR);
CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-atom_test", &atom_test,
PETSC_NULLPTR);
CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-save_every", &save_every,
PETSC_NULLPTR);
CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "", "-young_modulus",
&default_young_modulus, PETSC_NULLPTR);
CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "", "-poisson_ratio",
&default_poisson_ratio, PETSC_NULLPTR);
CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "", "-coeff_expansion",
&default_coeff_expansion, PETSC_NULLPTR);
CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "", "-ref_temp", &default_ref_temp,
PETSC_NULLPTR);
CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "", "-init_temp",
&default_init_temp, PETSC_NULLPTR);
CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "", "-capacity",
&default_heat_capacity, PETSC_NULLPTR);
CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "", "-conductivity",
&default_heat_conductivity, PETSC_NULLPTR);
if constexpr (SPACE_DIM == 2) {
do_output_domain = PETSC_TRUE;
do_output_skin = PETSC_FALSE;
} else {
do_output_domain = PETSC_FALSE;
do_output_skin = PETSC_TRUE;
}
CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-output_domain",
&do_output_domain, PETSC_NULLPTR);
CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-output_skin", &do_output_skin,
PETSC_NULLPTR);
MOFEM_LOG("ThermoElastic", Sev::inform)
<< "Young's modulus " << default_young_modulus;
MOFEM_LOG("ThermoElastic", Sev::inform)
<< "Poisson's ratio " << default_poisson_ratio;
MOFEM_LOG("ThermoElastic", Sev::inform)
<< "Coeff of expansion " << default_coeff_expansion;
MOFEM_LOG("ThermoElastic", Sev::inform)
<< "Default heat capacity " << default_heat_capacity;
MOFEM_LOG("ThermoElastic", Sev::inform)
<< "Heat conductivity " << default_heat_conductivity;
MOFEM_LOG("ThermoElastic", Sev::inform)
<< "Reference temperature " << default_ref_temp;
MOFEM_LOG("ThermoElastic", Sev::inform)
<< "Initial temperature " << default_init_temp;
};
CHKERR get_command_line_parameters();
}
//! [Get command line parameters]
//! [Set up problem]
// Add field
// Mechanical fields
CHKERR simple->addBoundaryField("U", H1, AINSWORTH_LEGENDRE_BASE, SPACE_DIM);
// Temperature
constexpr auto flux_space = (SPACE_DIM == 2) ? HCURL : HDIV;
CHKERR simple->addDomainField("T", L2, base, 1);
CHKERR simple->addDomainField("FLUX", flux_space, base, 1);
CHKERR simple->addBoundaryField("FLUX", flux_space, base, 1);
CHKERR simple->addBoundaryField("TBC", L2, base, 1);
CHKERR simple->setFieldOrder("U", order_disp);
CHKERR simple->setFieldOrder("FLUX", order_flux);
CHKERR simple->setFieldOrder("T", order_temp);
CHKERR simple->setFieldOrder("TBC", order_temp);
CHKERR simple->setUp();
int coords_dim = 3;
CHKERR PetscOptionsGetRealArray(NULL, NULL, "-field_eval_coords",
fieldEvalCoords.data(), &coords_dim,
tempFieldPtr = boost::make_shared<VectorDouble>();
fluxFieldPtr = boost::make_shared<MatrixDouble>();
dispFieldPtr = boost::make_shared<MatrixDouble>();
dispGradPtr = boost::make_shared<MatrixDouble>();
strainFieldPtr = boost::make_shared<MatrixDouble>();
stressFieldPtr = boost::make_shared<MatrixDouble>();
if (doEvalField) {
mField.getInterface<FieldEvaluatorInterface>()->getData<DomainEle>();
fieldEvalData, simple->getDomainFEName());
fieldEvalData->setEvalPoints(fieldEvalCoords.data(), 1);
auto no_rule = [](int, int, int) { return -1; };
auto field_eval_fe_ptr = fieldEvalData->feMethodPtr.lock();
field_eval_fe_ptr->getRuleHook = no_rule;
auto block_thermal_params = boost::make_shared<BlockedThermalParameters>();
auto block_thermoelastic_params =
boost::make_shared<BlockedThermoElasticParameters>();
auto coeff_expansion_ptr =
block_thermoelastic_params->getCoeffExpansionPtr();
auto ref_temp_ptr = block_thermoelastic_params->getRefTempPtr();
CHKERR addMatThermalBlockOps(field_eval_fe_ptr->getOpPtrVector(),
"MAT_THERMAL", block_thermal_params,
Sev::verbose);
field_eval_fe_ptr->getOpPtrVector(), "MAT_THERMOELASTIC",
block_thermoelastic_params, Sev::verbose);
field_eval_fe_ptr->getOpPtrVector(), {H1, HDIV});
auto hencky_common_data_ptr =
HenckyOps::commonDataFactory<SPACE_DIM, IT, DomainEleOp>(
mField, field_eval_fe_ptr->getOpPtrVector(), "U", "MAT_ELASTIC",
Sev::inform);
auto mat_D_ptr = hencky_common_data_ptr->matDPtr;
auto dispGradPtr = hencky_common_data_ptr->matGradPtr;
auto mat_stress_ptr = boost::make_shared<MatrixDouble>();
field_eval_fe_ptr->getOpPtrVector().push_back(
field_eval_fe_ptr->getOpPtrVector().push_back(
field_eval_fe_ptr->getOpPtrVector().push_back(
field_eval_fe_ptr->getOpPtrVector().push_back(
field_eval_fe_ptr->getOpPtrVector().push_back(
new
typename H::template OpCalculateHenckyThermalStress<SPACE_DIM, IT, 0>(
"U", tempFieldPtr, hencky_common_data_ptr, coeff_expansion_ptr,
ref_temp_ptr));
field_eval_fe_ptr->getOpPtrVector().push_back(
hencky_common_data_ptr->getMatFirstPiolaStress(),
field_eval_fe_ptr->getOpPtrVector().push_back(
} else {
field_eval_fe_ptr->getOpPtrVector().push_back(
new typename H::template OpCalculateLogC<SPACE_DIM, IT>(
"U", hencky_common_data_ptr));
stressFieldPtr = hencky_common_data_ptr->getMatFirstPiolaStress();
strainFieldPtr = hencky_common_data_ptr->getMatLogC();
};
}
}
//! [Set up problem]
//! [Boundary condition]
MOFEM_LOG("SYNC", Sev::noisy) << "bC";
auto bc_mng = mField.getInterface<BcManager>();
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, bool temp_bc = false) {
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);
}
};
if (!temp_bc) {
CHK_THROW_MESSAGE(remove_cubit_blocks(NODESET | TEMPERATURESET),
"remove_cubit_blocks");
CHK_THROW_MESSAGE(remove_named_blocks("TEMPERATURE"),
"remove_named_blocks");
}
CHK_THROW_MESSAGE(remove_cubit_blocks(SIDESET | HEATFLUXSET),
"remove_cubit_blocks");
CHK_THROW_MESSAGE(remove_named_blocks("HEATFLUX"), "remove_named_blocks");
CHK_THROW_MESSAGE(remove_named_blocks("CONVECTION"), "remove_named_blocks");
CHK_THROW_MESSAGE(remove_named_blocks("RADIATION"), "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()));
auto remove_temp_bc_ents =
filter_true_skin(filter_flux_blocks(get_skin(), true));
CHKERR mField.getInterface<CommInterface>()->synchroniseEntities(
remove_flux_ents);
CHKERR mField.getInterface<CommInterface>()->synchroniseEntities(
remove_temp_bc_ents);
MOFEM_LOG("SYNC", Sev::noisy) << remove_flux_ents << endl;
MOFEM_LOG("SYNC", Sev::noisy) << remove_temp_bc_ents << endl;
#ifndef NDEBUG
(boost::format("flux_remove_%d.vtk") % mField.get_comm_rank()).str(),
remove_flux_ents);
(boost::format("temp_bc_remove_%d.vtk") % mField.get_comm_rank()).str(),
remove_temp_bc_ents);
#endif
CHKERR mField.getInterface<ProblemsManager>()->removeDofsOnEntities(
simple->getProblemName(), "FLUX", remove_flux_ents);
CHKERR mField.getInterface<ProblemsManager>()->removeDofsOnEntities(
simple->getProblemName(), "TBC", remove_temp_bc_ents);
auto set_init_temp = [](boost::shared_ptr<FieldEntity> field_entity_ptr) {
field_entity_ptr->getEntFieldData()[0] = default_init_temp;
return 0;
};
CHKERR mField.getInterface<FieldBlas>()->fieldLambdaOnEntities(set_init_temp,
"T");
CHKERR bc_mng->removeBlockDOFsOnEntities<DisplacementCubitBcData>(
simple->getProblemName(), "U");
CHKERR bc_mng->pushMarkDOFsOnEntities<HeatFluxCubitBcData>(
simple->getProblemName(), "FLUX", false);
}
//! [Boundary condition]
//! [Push operators to pipeline]
MOFEM_LOG("SYNC", Sev::noisy) << "OPs";
auto pipeline_mng = mField.getInterface<PipelineManager>();
auto bc_mng = mField.getInterface<BcManager>();
auto boundary_marker =
bc_mng->getMergedBlocksMarker(vector<string>{"HEATFLUX"});
auto integration_rule = [](int, int, int approx_order) {
return 2 * approx_order;
};
CHKERR pipeline_mng->setDomainLhsIntegrationRule(integration_rule);
CHKERR pipeline_mng->setBoundaryRhsIntegrationRule(integration_rule);
CHKERR pipeline_mng->setBoundaryLhsIntegrationRule(integration_rule);
auto block_thermal_params = boost::make_shared<BlockedThermalParameters>();
auto heat_conductivity_ptr = block_thermal_params->getHeatConductivityPtr();
auto heat_capacity_ptr = block_thermal_params->getHeatCapacityPtr();
auto block_thermoelastic_params =
boost::make_shared<BlockedThermoElasticParameters>();
auto coeff_expansion_ptr = block_thermoelastic_params->getCoeffExpansionPtr();
auto ref_temp_ptr = block_thermoelastic_params->getRefTempPtr();
// Default time scaling of BCs and sources, from command line arguments
auto time_scale =
boost::make_shared<TimeScale>("", false, [](const double) { return 1; });
auto def_time_scale = [time_scale](const double time) {
return (!time_scale->argFileScale) ? time : 1;
};
auto def_file_name = [time_scale](const std::string &&name) {
return (!time_scale->argFileScale) ? name : "";
};
// Files which define scaling for separate variables, if provided
auto time_bodyforce_scaling = boost::make_shared<TimeScale>(
def_file_name("bodyforce_scale.txt"), false, def_time_scale);
auto time_heatsource_scaling = boost::make_shared<TimeScale>(
def_file_name("heatsource_scale.txt"), false, def_time_scale);
auto time_temperature_scaling = boost::make_shared<TimeScale>(
def_file_name("temperature_bc_scale.txt"), false, def_time_scale);
auto time_displacement_scaling = boost::make_shared<TimeScale>(
def_file_name("displacement_bc_scale.txt"), false, def_time_scale);
auto time_flux_scaling = boost::make_shared<TimeScale>(
def_file_name("flux_bc_scale.txt"), false, def_time_scale);
auto time_force_scaling = boost::make_shared<TimeScale>(
def_file_name("force_bc_scale.txt"), false, def_time_scale);
auto add_domain_rhs_ops = [&](auto &pipeline) {
CHKERR addMatThermalBlockOps(pipeline, "MAT_THERMAL", block_thermal_params,
Sev::inform);
CHKERR addMatThermoElasticBlockOps(pipeline, "MAT_THERMOELASTIC",
block_thermoelastic_params, Sev::inform);
auto hencky_common_data_ptr =
HenckyOps::commonDataFactory<SPACE_DIM, IT, DomainEleOp>(
mField, pipeline, "U", "MAT_ELASTIC", Sev::inform);
auto mat_D_ptr = hencky_common_data_ptr->matDPtr;
auto mat_grad_ptr = hencky_common_data_ptr->matGradPtr;
auto mat_strain_ptr = boost::make_shared<MatrixDouble>();
auto mat_stress_ptr = boost::make_shared<MatrixDouble>();
auto vec_temp_ptr = boost::make_shared<VectorDouble>();
auto vec_temp_dot_ptr = boost::make_shared<VectorDouble>();
auto mat_flux_ptr = boost::make_shared<MatrixDouble>();
auto vec_temp_div_ptr = boost::make_shared<VectorDouble>();
pipeline.push_back(new OpCalculateScalarFieldValues("T", vec_temp_ptr));
pipeline.push_back(
new OpCalculateScalarFieldValuesDot("T", vec_temp_dot_ptr));
"FLUX", vec_temp_div_ptr));
pipeline.push_back(
new OpCalculateHVecVectorField<3, SPACE_DIM>("FLUX", mat_flux_ptr));
CHKERR opThermoElasticFactoryDomainRhs<SPACE_DIM, AT, IT, DomainEleOp>(
mField, pipeline, "U", "MAT_ELASTIC", "MAT_THERMAL",
"MAT_THERMOELASTIC", Sev::inform);
auto resistance = [heat_conductivity_ptr](const double, const double,
const double) {
return (1. / (*heat_conductivity_ptr));
};
// negative value is consequence of symmetric system
auto capacity = [heat_capacity_ptr](const double, const double,
const double) {
return -(*heat_capacity_ptr);
};
auto unity = [](const double, const double, const double) constexpr {
return -1.;
};
pipeline.push_back(
new OpHdivFlux("FLUX", mat_flux_ptr, resistance, mat_grad_ptr));
pipeline.push_back(new OpHDivTemp("FLUX", vec_temp_ptr, unity));
pipeline.push_back(new OpBaseDivFlux("T", vec_temp_div_ptr, unity));
pipeline.push_back(new OpBaseDotT("T", vec_temp_dot_ptr, capacity));
// Set source terms
CHKERR DomainNaturalBCRhs::AddFluxToPipeline<OpHeatSource>::add(
pipeline, mField, "T", {time_scale, time_heatsource_scaling},
"HEAT_SOURCE", Sev::inform);
CHKERR DomainNaturalBCRhs::AddFluxToPipeline<OpBodyForce>::add(
pipeline, mField, "U", {time_scale, time_bodyforce_scaling},
"BODY_FORCE", Sev::inform);
CHKERR DomainNaturalBCRhs::AddFluxToPipeline<OpSetTemperatureRhs>::add(
pipeline, mField, "T", vec_temp_ptr, "SETTEMP", Sev::inform);
};
auto add_domain_lhs_ops = [&](auto &pipeline) {
CHKERR addMatThermalBlockOps(pipeline, "MAT_THERMAL", block_thermal_params,
Sev::verbose);
CHKERR addMatThermoElasticBlockOps(pipeline, "MAT_THERMOELASTIC",
block_thermoelastic_params,
Sev::verbose);
auto hencky_common_data_ptr =
HenckyOps::commonDataFactory<SPACE_DIM, IT, DomainEleOp>(
mField, pipeline, "U", "MAT_ELASTIC", Sev::inform, 1);
auto mat_D_ptr = hencky_common_data_ptr->matDPtr;
auto mat_grad_ptr = hencky_common_data_ptr->matGradPtr;
auto resistance = [heat_conductivity_ptr](const double, const double,
const double) {
return (1. / (*heat_conductivity_ptr));
};
auto capacity = [heat_capacity_ptr](const double, const double,
const double) {
return -(*heat_capacity_ptr);
};
pipeline.push_back(
new OpHdivHdiv("FLUX", "FLUX", resistance, mat_grad_ptr));
pipeline.push_back(
new OpHdivT("FLUX", "T", []() constexpr { return -1; }, true));
auto mat_flux_ptr = boost::make_shared<MatrixDouble>();
pipeline.push_back(
new OpCalculateHVecVectorField<3, SPACE_DIM>("FLUX", mat_flux_ptr));
pipeline.push_back(
new OpHdivU("FLUX", "U", mat_flux_ptr, resistance, mat_grad_ptr));
CHKERR opThermoElasticFactoryDomainLhs<SPACE_DIM, AT, IT, DomainEleOp>(
mField, pipeline, "U", "T", "MAT_ELASTIC", "MAT_THERMAL",
"MAT_THERMOELASTIC", Sev::inform);
auto op_capacity = new OpCapacity("T", "T", capacity);
op_capacity->feScalingFun = [](const FEMethod *fe_ptr) {
return fe_ptr->ts_a;
};
pipeline.push_back(op_capacity);
auto vec_temp_ptr = boost::make_shared<VectorDouble>();
pipeline.push_back(new OpCalculateScalarFieldValues("T", vec_temp_ptr));
CHKERR DomainNaturalBCLhs::AddFluxToPipeline<OpSetTemperatureLhs>::add(
pipeline, mField, "T", vec_temp_ptr, "SETTEMP", Sev::verbose);
};
auto add_boundary_rhs_ops = [&](auto &pipeline) {
CHKERR BoundaryNaturalBC::AddFluxToPipeline<OpForce>::add(
pipeline, mField, "U", {time_scale, time_force_scaling}, "FORCE",
Sev::inform);
// Temperature BC
BoundaryNaturalBC::OpFlux<NaturalTemperatureMeshsets, 3, SPACE_DIM>;
CHKERR BoundaryNaturalBC::AddFluxToPipeline<OpTemperatureBC>::add(
pipeline, mField, "FLUX", {time_scale, time_temperature_scaling},
"TEMPERATURE", Sev::inform);
// Note: fluxes for temperature should be aggregated. Here separate is
// NaturalMeshsetType<HEATFLUXSET>, we should add version with BLOCKSET,
// convection, see example tutorials/vec-0/src/NaturalBoundaryBC.hpp.
// and vec-0/elastic.cpp
using OpFluxBC =
BoundaryNaturalBC::OpFlux<NaturalMeshsetType<HEATFLUXSET>, 1, 1>;
CHKERR BoundaryNaturalBC::AddFluxToPipeline<OpFluxBC>::add(
pipeline, mField, "TBC", {time_scale, time_flux_scaling}, "FLUX",
Sev::inform);
using OpConvectionRhsBC =
T::OpFlux<ThermoElasticOps::ConvectionBcType<BLOCKSET>, 1, 1>;
using OpRadiationRhsBC =
T::OpFlux<ThermoElasticOps::RadiationBcType<BLOCKSET>, 1, 1>;
auto temp_bc_ptr = boost::make_shared<VectorDouble>();
pipeline.push_back(new OpCalculateScalarFieldValues("TBC", temp_bc_ptr));
T::AddFluxToPipeline<OpConvectionRhsBC>::add(
pipeline, mField, "TBC", temp_bc_ptr, "CONVECTION", Sev::inform);
T::AddFluxToPipeline<OpRadiationRhsBC>::add(
pipeline, mField, "TBC", temp_bc_ptr, "RADIATION", Sev::inform);
auto mat_flux_ptr = boost::make_shared<MatrixDouble>();
pipeline.push_back(
new OpCalculateHVecVectorField<3, SPACE_DIM>("FLUX", mat_flux_ptr));
// This is temporary implementation. It be moved to LinearFormsIntegrators,
// however for hybridised case, what is goal of this changes, such function
// is implemented for fluxes in broken space. Thus ultimately this operator
// would be not needed.
struct OpTBCTimesNormalFlux
: public FormsIntegrators<BoundaryEleOp>::Assembly<PETSC>::OpBase {
OpTBCTimesNormalFlux(const std::string field_name,
boost::shared_ptr<MatrixDouble> vec,
boost::shared_ptr<Range> ents_ptr = nullptr)
: OP(field_name, field_name, OP::OPROW, ents_ptr), sourceVec(vec) {}
// get integration weights
auto t_w = OP::getFTensor0IntegrationWeight();
// get base function values on rows
auto t_row_base = row_data.getFTensor0N();
// get normal vector
auto t_normal = OP::getFTensor1NormalsAtGaussPts();
// get vector values
auto t_vec = getFTensor1FromMat<SPACE_DIM, 1>(*sourceVec);
// loop over integration points
for (int gg = 0; gg != OP::nbIntegrationPts; gg++) {
// take into account Jacobian
const double alpha = t_w * (t_vec(i) * t_normal(i));
// loop over rows base functions
int rr = 0;
for (; rr != OP::nbRows; ++rr) {
OP::locF[rr] += alpha * t_row_base;
++t_row_base;
}
for (; rr < OP::nbRowBaseFunctions; ++rr)
++t_row_base;
++t_w; // move to another integration weight
++t_vec;
++t_normal;
}
EntityType fe_type = OP::getNumeredEntFiniteElementPtr()->getEntType();
if (fe_type == MBTRI) {
OP::locF /= 2;
}
}
protected:
boost::shared_ptr<MatrixDouble> sourceVec;
};
pipeline.push_back(new OpTBCTimesNormalFlux("TBC", mat_flux_ptr));
struct OpBaseNormalTimesTBC
: public FormsIntegrators<BoundaryEleOp>::Assembly<PETSC>::OpBase {
OpBaseNormalTimesTBC(const std::string field_name,
boost::shared_ptr<VectorDouble> vec,
boost::shared_ptr<Range> ents_ptr = nullptr)
: OP(field_name, field_name, OP::OPROW, ents_ptr), sourceVec(vec) {}
// get integration weights
auto t_w = OP::getFTensor0IntegrationWeight();
// get base function values on rows
auto t_row_base = row_data.getFTensor1N<3>();
// get normal vector
auto t_normal = OP::getFTensor1NormalsAtGaussPts();
// get vector values
auto t_vec = getFTensor0FromVec(*sourceVec);
// loop over integration points
for (int gg = 0; gg != OP::nbIntegrationPts; gg++) {
// take into account Jacobian
const double alpha = t_w * t_vec;
// loop over rows base functions
int rr = 0;
for (; rr != OP::nbRows; ++rr) {
OP::locF[rr] += alpha * (t_row_base(i) * t_normal(i));
++t_row_base;
}
for (; rr < OP::nbRowBaseFunctions; ++rr)
++t_row_base;
++t_w; // move to another integration weight
++t_vec;
++t_normal;
}
EntityType fe_type = OP::getNumeredEntFiniteElementPtr()->getEntType();
if (fe_type == MBTRI) {
OP::locF /= 2;
}
}
protected:
boost::shared_ptr<VectorDouble> sourceVec;
};
pipeline.push_back(new OpBaseNormalTimesTBC("FLUX", temp_bc_ptr));
};
auto add_boundary_lhs_ops = [&](auto &pipeline) {
using T =
NaturalBC<BoundaryEleOp>::template Assembly<PETSC>::BiLinearForm<GAUSS>;
using OpConvectionLhsBC =
T::OpFlux<ThermoElasticOps::ConvectionBcType<BLOCKSET>, 1, 1>;
using OpRadiationLhsBC =
T::OpFlux<ThermoElasticOps::RadiationBcType<BLOCKSET>, 1, 1>;
auto temp_bc_ptr = boost::make_shared<VectorDouble>();
pipeline.push_back(new OpCalculateScalarFieldValues("TBC", temp_bc_ptr));
T::AddFluxToPipeline<OpConvectionLhsBC>::add(pipeline, mField, "TBC", "TBC",
"CONVECTION", Sev::verbose);
T::AddFluxToPipeline<OpRadiationLhsBC>::add(
pipeline, mField, "TBC", "TBC", temp_bc_ptr, "RADIATION", Sev::verbose);
// This is temporary implementation. It be moved to LinearFormsIntegrators,
// however for hybridised case, what is goal of this changes, such function
// is implemented for fluxes in broken space. Thus ultimately this operator
// would be not needed.
struct OpTBCTimesNormalFlux
: public FormsIntegrators<BoundaryEleOp>::Assembly<PETSC>::OpBase {
OpTBCTimesNormalFlux(const std::string row_field_name,
const std::string col_field_name,
boost::shared_ptr<Range> ents_ptr = nullptr)
: OP(row_field_name, col_field_name, OP::OPROWCOL, ents_ptr) {
this->sYmm = false;
this->assembleTranspose = true;
this->onlyTranspose = false;
}
// get integration weights
auto t_w = OP::getFTensor0IntegrationWeight();
// get base function values on rows
auto t_row_base = row_data.getFTensor0N();
// get normal vector
auto t_normal = OP::getFTensor1NormalsAtGaussPts();
// loop over integration points
for (int gg = 0; gg != OP::nbIntegrationPts; gg++) {
// loop over rows base functions
auto a_mat_ptr = &*OP::locMat.data().begin();
int rr = 0;
for (; rr != OP::nbRows; rr++) {
// take into account Jacobian
const double alpha = t_w * t_row_base;
// get column base functions values at gauss point gg
auto t_col_base = col_data.getFTensor1N<3>(gg, 0);
// loop over columns
for (int cc = 0; cc != OP::nbCols; cc++) {
// calculate element of local matrix
// using L2 norm of normal vector for convenient area calculation
// for quads, tris etc.
*a_mat_ptr += alpha * (t_col_base(i) * t_normal(i));
++t_col_base;
++a_mat_ptr;
}
++t_row_base;
}
for (; rr < OP::nbRowBaseFunctions; ++rr)
++t_row_base;
++t_normal;
++t_w; // move to another integration weight
}
EntityType fe_type = OP::getNumeredEntFiniteElementPtr()->getEntType();
if (fe_type == MBTRI) {
OP::locMat /= 2;
}
}
};
pipeline.push_back(new OpTBCTimesNormalFlux("TBC", "FLUX"));
};
// Set BCs by eliminating degrees of freedom
auto get_bc_hook_rhs = [&]() {
mField, pipeline_mng->getDomainRhsFE(),
{time_scale, time_displacement_scaling});
return hook;
};
auto get_bc_hook_lhs = [&]() {
mField, pipeline_mng->getDomainLhsFE(),
{time_scale, time_displacement_scaling});
return hook;
};
pipeline_mng->getDomainRhsFE()->preProcessHook = get_bc_hook_rhs();
pipeline_mng->getDomainLhsFE()->preProcessHook = get_bc_hook_lhs();
CHKERR add_domain_rhs_ops(pipeline_mng->getOpDomainRhsPipeline());
CHKERR add_domain_lhs_ops(pipeline_mng->getOpDomainLhsPipeline());
CHKERR add_boundary_rhs_ops(pipeline_mng->getOpBoundaryRhsPipeline());
CHKERR add_boundary_lhs_ops(pipeline_mng->getOpBoundaryLhsPipeline());
}
//! [Push operators to pipeline]
//! [Solve]
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_elements = [&]() {
auto block_thermal_params = boost::make_shared<BlockedThermalParameters>();
auto block_thermoelastic_params =
boost::make_shared<BlockedThermoElasticParameters>();
auto coeff_expansion_ptr =
block_thermoelastic_params->getCoeffExpansionPtr();
auto ref_temp_ptr = block_thermoelastic_params->getRefTempPtr();
auto u_ptr = boost::make_shared<MatrixDouble>();
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 vec_temp_ptr = boost::make_shared<VectorDouble>();
auto mat_flux_ptr = boost::make_shared<MatrixDouble>();
auto push_domain_ops = [&](auto &pp_fe) {
auto &pip = pp_fe->getOpPtrVector();
CHKERR addMatThermalBlockOps(pip, "MAT_THERMAL", block_thermal_params,
Sev::verbose);
pip, "MAT_THERMOELASTIC", block_thermoelastic_params, Sev::verbose);
pip.push_back(new OpCalculateScalarFieldValues("T", vec_temp_ptr));
pip.push_back(
new OpCalculateHVecVectorField<3, SPACE_DIM>("FLUX", mat_flux_ptr));
pip.push_back(new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_ptr));
"U", mat_grad_ptr));
auto elastic_common_ptr = commonDataFactory<SPACE_DIM, IT, DomainEleOp>(
mField, pip, "U", "MAT_ELASTIC", Sev::inform);
pip.push_back(
new
typename H::template OpCalculateHenckyThermalStress<SPACE_DIM, IT, 0>(
"U", vec_temp_ptr, elastic_common_ptr, coeff_expansion_ptr,
ref_temp_ptr));
pip.push_back(new OpSymmetrizeTensor<SPACE_DIM>(
elastic_common_ptr->getMatFirstPiolaStress(), mat_stress_ptr));
pip.push_back(
new OpSymmetrizeTensor<SPACE_DIM>(mat_grad_ptr, mat_strain_ptr));
} else {
mat_stress_ptr = elastic_common_ptr->getMatFirstPiolaStress();
mat_strain_ptr = elastic_common_ptr->getMatLogC();
}
};
auto push_post_proc_ops = [&](auto &pp_fe) {
auto &pip = pp_fe->getOpPtrVector();
pip.push_back(
new OpPPMap(
pp_fe->getPostProcMesh(), pp_fe->getMapGaussPts(),
{{"T", vec_temp_ptr}},
{{"U", u_ptr}, {"FLUX", mat_flux_ptr}},
{},
{{"CAUCHY", mat_stress_ptr}, {"STRAIN", mat_strain_ptr}}
)
);
} else {
pip.push_back(
new OpPPMap(
pp_fe->getPostProcMesh(), pp_fe->getMapGaussPts(),
{{"T", vec_temp_ptr}},
{{"U", u_ptr}, {"FLUX", mat_flux_ptr}},
{{"PIOLA", mat_stress_ptr}},
{{"HENCKY_STRAIN", mat_strain_ptr}}
)
);
}
};
auto domain_post_proc = [&]() {
if (do_output_domain == PETSC_FALSE)
return boost::shared_ptr<PostProcEle>();
auto pp_fe = boost::make_shared<PostProcEle>(mField);
CHK_MOAB_THROW(push_domain_ops(pp_fe),
"push domain ops to domain element");
CHK_MOAB_THROW(push_post_proc_ops(pp_fe),
"push post proc ops to domain element");
return pp_fe;
};
auto skin_post_proc = [&]() {
if (do_output_skin == PETSC_FALSE)
return boost::shared_ptr<SkinPostProcEle>();
auto pp_fe = boost::make_shared<SkinPostProcEle>(mField);
auto simple = mField.getInterface<Simple>();
auto op_side = new OpLoopSide<SideEle>(mField, simple->getDomainFEName(),
SPACE_DIM, Sev::verbose);
CHK_MOAB_THROW(push_domain_ops(op_side),
"push domain ops to side element");
pp_fe->getOpPtrVector().push_back(op_side);
CHK_MOAB_THROW(push_post_proc_ops(pp_fe),
"push post proc ops to skin element");
return pp_fe;
};
return std::make_pair(domain_post_proc(), skin_post_proc());
};
auto monitor_ptr = boost::make_shared<FEMethod>();
auto set_time_monitor = [&](auto dm, auto solver, auto domain_post_proc_fe,
auto skin_post_proc_fe) {
monitor_ptr->preProcessHook = [&]() {
if (save_every && (monitor_ptr->ts_step % save_every == 0)) {
CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(),
domain_post_proc_fe,
monitor_ptr->getCacheWeakPtr());
CHKERR domain_post_proc_fe->writeFile(
"out_" + boost::lexical_cast<std::string>(monitor_ptr->ts_step) +
".h5m");
}
CHKERR DMoFEMLoopFiniteElements(dm, simple->getBoundaryFEName(),
skin_post_proc_fe,
monitor_ptr->getCacheWeakPtr());
CHKERR skin_post_proc_fe->writeFile(
"out_skin_" +
boost::lexical_cast<std::string>(monitor_ptr->ts_step) + ".h5m");
}
}
if (doEvalField) {
->evalFEAtThePoint<SPACE_DIM>(
fieldEvalCoords.data(), 1e-12, simple->getProblemName(),
simple->getDomainFEName(), fieldEvalData,
mField.get_comm_rank(), mField.get_comm_rank(), nullptr,
if (atom_test) {
auto eval_num_vec =
createVectorMPI(mField.get_comm(), PETSC_DECIDE, 1);
CHKERR VecZeroEntries(eval_num_vec);
if (tempFieldPtr->size()) {
CHKERR VecSetValue(eval_num_vec, 0, 1, ADD_VALUES);
}
CHKERR VecAssemblyBegin(eval_num_vec);
CHKERR VecAssemblyEnd(eval_num_vec);
double eval_num;
CHKERR VecSum(eval_num_vec, &eval_num);
if (!(int)eval_num) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: did not find elements to evaluate "
"the field, check the coordinates",
}
}
if (tempFieldPtr->size()) {
auto t_temp = getFTensor0FromVec(*tempFieldPtr);
MOFEM_LOG("ThermoElasticSync", Sev::inform)
<< "Eval point T: " << t_temp;
if (atom_test && fabs(monitor_ptr->ts_t - 10) < 1e-12) {
if (atom_test <= 3 && fabs(t_temp - 554.48) > 1e-2) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: wrong temperature value",
}
if (atom_test == 4 && fabs(t_temp - 250) > 1e-2) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: wrong temperature value",
}
if (atom_test == 5 && fabs(t_temp - 1) > 1e-2) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: wrong temperature value",
}
}
}
if (fluxFieldPtr->size1()) {
auto t_flux = getFTensor1FromMat<SPACE_DIM>(*fluxFieldPtr);
auto flux_mag = sqrt(t_flux(i) * t_flux(i));
MOFEM_LOG("ThermoElasticSync", Sev::inform)
<< "Eval point FLUX magnitude: " << flux_mag;
if (atom_test && fabs(monitor_ptr->ts_t - 10) < 1e-12) {
if (atom_test <= 3 && fabs(flux_mag - 27008.0) > 2e1) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: wrong flux value", atom_test);
}
if (atom_test == 4 && fabs(flux_mag - 150e3) > 1e-1) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: wrong flux value", atom_test);
}
if (atom_test == 5 && fabs(flux_mag) > 1e-6) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: wrong flux value", atom_test);
}
}
}
if (dispFieldPtr->size1()) {
auto t_disp = getFTensor1FromMat<SPACE_DIM>(*dispFieldPtr);
auto disp_mag = sqrt(t_disp(i) * t_disp(i));
MOFEM_LOG("ThermoElasticSync", Sev::inform)
<< "Eval point U magnitude: " << disp_mag;
if (atom_test && fabs(monitor_ptr->ts_t - 10) < 1e-12) {
if (atom_test == 1 && fabs(disp_mag - 0.00345) > 1e-5) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: wrong displacement value",
}
if ((atom_test == 2 || atom_test == 3) &&
fabs(disp_mag - 0.00265) > 1e-5) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: wrong displacement value",
}
if ((atom_test == 5) &&
fabs(t_disp(0) - (std::sqrt(std::exp(0.2)) - 1)) > 1e-5 &&
fabs(t_disp(1) - (std::sqrt(std::exp(0.2)) - 1)) > 1e-5) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: wrong displacement value",
}
}
}
if (strainFieldPtr->size1()) {
auto t_strain =
getFTensor2SymmetricFromMat<SPACE_DIM>(*strainFieldPtr);
auto t_strain_trace = t_strain(i, i);
if (atom_test && fabs(monitor_ptr->ts_t - 10) < 1e-12) {
if (atom_test == 1 && fabs(t_strain_trace - 0.00679) > 1e-5) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: wrong strain value", atom_test);
}
if ((atom_test == 2 || atom_test == 3) &&
fabs(t_strain_trace - 0.00522) > 1e-5) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: wrong strain value", atom_test);
}
if ((atom_test == 5) && fabs(t_strain_trace - 0.2) > 1e-5) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: wrong strain value", atom_test);
}
}
}
if (stressFieldPtr->size1()) {
double von_mises_stress;
auto getVonMisesStress = [&](auto t_stress) {
von_mises_stress = std::sqrt(
0.5 *
((t_stress(0, 0) - t_stress(1, 1)) *
(t_stress(0, 0) - t_stress(1, 1)) +
(SPACE_DIM == 3 ? (t_stress(1, 1) - t_stress(2, 2)) *
(t_stress(1, 1) - t_stress(2, 2))
: 0) +
(SPACE_DIM == 3 ? (t_stress(2, 2) - t_stress(0, 0)) *
(t_stress(2, 2) - t_stress(0, 0))
: 0) +
6 * (t_stress(0, 1) * t_stress(0, 1) +
(SPACE_DIM == 3 ? t_stress(1, 2) * t_stress(1, 2) : 0) +
(SPACE_DIM == 3 ? t_stress(2, 0) * t_stress(2, 0) : 0))));
};
auto t_stress =
getFTensor2SymmetricFromMat<SPACE_DIM>(*stressFieldPtr);
CHKERR getVonMisesStress(t_stress);
} else {
auto t_stress =
getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*stressFieldPtr);
CHKERR getVonMisesStress(t_stress);
}
MOFEM_LOG("ThermoElasticSync", Sev::inform)
<< "Eval point von Mises Stress: " << von_mises_stress;
if (atom_test && fabs(monitor_ptr->ts_t - 10) < 1e-12) {
if (atom_test == 1 && fabs(von_mises_stress - 523.0) > 5e-1) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: wrong von Mises stress value",
}
if (atom_test == 2 && fabs(von_mises_stress - 16.3) > 5e-2) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: wrong von Mises stress value",
}
if (atom_test == 3 && fabs(von_mises_stress - 14.9) > 5e-2) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: wrong von Mises stress value",
}
if (atom_test == 5 && fabs(von_mises_stress) > 5e-2) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"atom test %d failed: wrong von Mises stress value",
}
}
}
}
};
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 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, "T", 0, 0,
is_T);
CHKERR is_mng->isCreateProblemFieldAndRank(name_prb, ROW, "TBC", 0, 0,
is_TBC);
IS is_tmp, is_tmp2;
CHKERR ISExpand(is_T, is_flux, &is_tmp);
CHKERR ISExpand(is_TBC, is_tmp, &is_tmp2);
CHKERR ISDestroy(&is_tmp);
auto is_TFlux = SmartPetscObj<IS>(is_tmp2);
CHKERR ISSort(is_u);
CHKERR ISSort(is_TFlux);
CHKERR PCFieldSplitSetIS(pc, PETSC_NULLPTR, is_TFlux);
CHKERR PCFieldSplitSetIS(pc, PETSC_NULLPTR, is_u);
}
};
auto B = createDMMatrix(dm);
CHKERR TSSetIJacobian(solver, B, B, PETSC_NULLPTR, PETSC_NULLPTR);
auto D = createDMVector(dm);
CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_FORWARD);
CHKERR TSSetSolution(solver, D);
CHKERR TSSetFromOptions(solver);
CHKERR set_section_monitor(solver);
CHKERR set_fieldsplit_preconditioner(solver);
auto [domain_post_proc_fe, skin_post_proc_fe] =
create_post_process_elements();
CHKERR set_time_monitor(dm, solver, domain_post_proc_fe, skin_post_proc_fe);
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(), "ThermoElastic"));
LogManager::setLog("ThermoElastic");
MOFEM_LOG_TAG("ThermoElastic", "ThermoElastic");
core_log->add_sink(
LogManager::createSink(LogManager::getStrmSync(), "ThermoElasticSync"));
LogManager::setLog("ThermoElasticSync");
MOFEM_LOG_TAG("ThermoElasticSync", "ThermoElasticSync");
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 interface
//! [Create MoFEM]
//! [Load mesh]
CHKERR simple->loadFile();
//! [Load mesh]
//! [ThermoElasticProblem]
ThermoElasticProblem ex(m_field);
CHKERR ex.runProblem();
//! [ThermoElasticProblem]
}
}
static auto filter_true_skin(MoFEM::Interface &m_field, Range &&skin)
Operators for Piola transformation of the thermal conductivity.
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
#define MOFEM_LOG_SYNCHRONISE(comm)
Synchronise "SYNC" channel.
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
#define FTENSOR_INDEX(DIM, I)
Implementation of thermal convection bc.
Thermal operators agnostic to small/large deformations.
Implementation of thermal radiation bc.
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
@ QUIET
@ ROW
#define CATCH_ERRORS
Catch errors.
@ MF_EXIST
FieldApproximationBase
approximation base
Definition definitions.h:58
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ DEMKOWICZ_JACOBI_BASE
Definition definitions.h:66
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ L2
field with C-1 continuity
Definition definitions.h:88
@ H1
continuous field
Definition definitions.h:85
@ NOSPACE
Definition definitions.h:83
@ HCURL
field with continuous tangents
Definition definitions.h:86
@ HDIV
field with continuous normal traction
Definition definitions.h:87
#define MYPCOMM_INDEX
default communicator number PCOMM
#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.
@ TEMPERATURESET
@ HEATFLUXSET
@ NODESET
@ SIDESET
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
@ 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.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
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 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:576
auto createDMVector(DM dm)
Get smart vector from DM.
Definition DMMoFEM.hpp:1234
auto createDMMatrix(DM dm)
Get smart matrix from DM.
Definition DMMoFEM.hpp:1191
IntegrationType
Form integrator integration types.
AssemblyType
[Storage and set boundary conditions]
@ GAUSS
Gaussian quadrature integration.
@ PETSC
Standard PETSc assembly.
#define MOFEM_LOG(channel, severity)
Log.
SeverityLevel
Severity levels.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
MoFEMErrorCode getCubitMeshsetPtr(const int ms_id, const CubitBCType cubit_bc_type, const CubitMeshSets **cubit_meshset_ptr) const
get cubit meshset
FTensor::Index< 'i', SPACE_DIM > i
PetscBool doEvalField
const double c
speed of light (cm/ns)
double D
const double n
refractive index of diffusive medium
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasVector< double > VectorDouble
Definition Types.hpp:68
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
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:1046
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
MoFEMErrorCode MoFEMSNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *ctx)
Sens monitor printing residual field by field.
Definition SnesCtx.cpp:592
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)
PetscErrorCode PetscOptionsGetRealArray(PetscOptions *, const char pre[], const char name[], PetscReal dval[], PetscInt *nmax, PetscBool *set)
auto createVectorMPI(MPI_Comm comm, PetscInt n, PetscInt N)
Create MPI Vector.
OpCalculateScalarFieldValuesFromPetscVecImpl< PetscData::CTX_SET_X_T > OpCalculateScalarFieldValuesDot
auto getDMSnesCtx(DM dm)
Get SNES context data structure used by DM.
Definition DMMoFEM.hpp:1262
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)
OpCalculateQdotQLhs_dU< SPACE_DIM, GAUSS, AssemblyDomainEleOp, IS_LARGE_STRAINS > OpHdivU
Integrate Lhs of flux term coupled to displacement field.
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMixDivTimesScalar< SPACE_DIM > OpHdivT
Integrate Lhs div of base of flux times base of temperature (FLUX x T) and transpose of it,...
OpHdivFluxImpl< SPACE_DIM, IS_LARGE_STRAINS > OpHdivFlux
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpMixDivTimesU< 3, 1, SPACE_DIM > OpHDivTemp
Integrate Rhs div flux base times temperature (T)
OpBaseDotT OpBaseDivFlux
Integrate Rhs base of temperature times divergence of flux (T)
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesScalar< 1 > OpBaseDotT
Integrate Rhs base of temperature time heat capacity times heat rate (T)
OpHdivHdivImpl< SPACE_DIM, IS_LARGE_STRAINS > OpHdivHdiv
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor
constexpr IntegrationType I
constexpr AssemblyType A
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
PostProcBrokenMeshInMoab< BoundaryEle > SkinPostProcEle
Definition plastic.cpp:60
int atom_test
Atom test.
Definition plastic.cpp:121
#define EXECUTABLE_DIMENSION
Definition plastic.cpp:13
double scale
Definition plastic.cpp:123
double H
Hardening.
Definition plastic.cpp:128
ElementsAndOps< SPACE_DIM >::SideEle SideEle
Definition plastic.cpp:61
constexpr auto field_name
static constexpr int approx_order
OpBaseImpl< PETSC, EdgeEleOp > OpBase
Definition radiation.cpp:29
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, SPACE_DIM, SPACE_DIM > OpInternalForcePiola
Definition seepage.cpp:65
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradTensorGrad< 1, SPACE_DIM, SPACE_DIM, 1 > OpKPiola
[Only used for dynamics]
Definition seepage.cpp:63
BoundaryNaturalBC::OpFlux< NaturalTemperatureMeshsets, 3, SPACE_DIM > OpTemperatureBC
Definition seepage.cpp:147
FTensor::Index< 'm', 3 > m
Boundary condition manager for finite element problem setup.
Managing BitRefLevels.
virtual moab::Interface & get_moab()=0
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:118
Deprecated interface functions.
Definition of the displacement bc data structure.
Definition BCData.hpp:72
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
Class (Function) to enforce essential constrains.
Definition Essential.hpp:25
Basic algebra on fields.
Definition FieldBlas.hpp:21
Field evaluator interface.
Definition of the heat flux bc data structure.
Definition BCData.hpp:423
Section manager is used to create indexes and sections.
Definition ISManager.hpp:23
Interface for managing meshsets containing materials and boundary conditions.
Assembly methods.
Definition Natural.hpp:65
Natural boundary conditions.
Definition Natural.hpp:57
Get vector field for H-div approximation.
Calculate divergence of vector field.
Specialization for double precision scalar field values calculation.
Get field gradients at integration pts for scalar field rank 0, i.e. vector field.
Specialization for double precision vector field values calculation.
Element used to execute operators on side of the element.
Post post-proc data at points from hash maps.
Operator for symmetrizing tensor fields.
Template struct for dimension-specific finite element types.
PipelineManager interface.
MoFEMErrorCode setDomainRhsIntegrationRule(RuleHookFun rule)
Set integration rule for domain right-hand side finite element.
Problem manager is used to build and partition problems.
Simple interface for fast problem set-up.
Definition Simple.hpp:27
MoFEMErrorCode addDomainField(const std::string name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_ZERO, int verb=-1)
Add field on domain.
Definition Simple.cpp:261
MoFEMErrorCode getOptions()
get options
Definition Simple.cpp:180
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
boost::shared_ptr< FieldEvaluatorInterface::SetPtsData > fieldEvalData
MoFEMErrorCode getCommandLineParameters()
[Run problem]
MoFEMErrorCode runProblem()
[Run problem]
boost::shared_ptr< MatrixDouble > fluxFieldPtr
boost::shared_ptr< VectorDouble > tempFieldPtr
boost::shared_ptr< MatrixDouble > stressFieldPtr
MoFEM::Interface & mField
std::array< double, 3 > fieldEvalCoords
boost::shared_ptr< MatrixDouble > dispGradPtr
MoFEMErrorCode opThermoElasticFactoryDomainRhs(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string field_name, boost::shared_ptr< HenckyOps::CommonData > elastic_common_ptr, boost::shared_ptr< ThermoElasticProblem::BlockedThermalParameters > thermal_common_ptr, boost::shared_ptr< ThermoElasticProblem::BlockedThermoElasticParameters > thermoelastic_common_ptr, Sev sev)
MoFEMErrorCode setupProblem()
add fields
boost::shared_ptr< MatrixDouble > dispFieldPtr
MoFEMErrorCode bC()
[Set up problem]
MoFEMErrorCode OPs()
[Boundary condition]
boost::shared_ptr< MatrixDouble > strainFieldPtr
MoFEMErrorCode tsSolve()
[Push operators to pipeline]
MoFEMErrorCode opThermoElasticFactoryDomainLhs(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string field_name, std::string coupled_field_name, boost::shared_ptr< HenckyOps::CommonData > elastic_common_ptr, boost::shared_ptr< ThermoElasticProblem::BlockedThermalParameters > thermal_common_ptr, boost::shared_ptr< ThermoElasticProblem::BlockedThermoElasticParameters > thermoelastic_common_ptr, Sev sev)
MoFEMErrorCode addMatThermalBlockOps(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, std::string block_name, boost::shared_ptr< BlockedThermalParameters > blockedParamsPtr, Sev sev)
MoFEMErrorCode addMatThermoElasticBlockOps(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, std::string block_name, boost::shared_ptr< BlockedThermoElasticParameters > blockedParamsPtr, Sev sev)
constexpr AssemblyType AT
constexpr IntegrationType IT
EssentialBC< BoundaryEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpEssentialLhs< HeatFluxCubitBcData, 3, SPACE_DIM > OpEssentialFluxLhs
#define FINITE_DEFORMATION_FLAG
DomainNaturalBCLhs::OpFlux< SetTargetTemperature, 1, 1 > OpSetTemperatureLhs
NaturalBC< BoundaryEleOp >::Assembly< PETSC >::LinearForm< GAUSS > BoundaryNaturalBC
[Body and heat source]
EssentialBC< BoundaryEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpEssentialRhs< HeatFluxCubitBcData, 3, SPACE_DIM > OpEssentialFluxRhs
[Natural boundary conditions]
int save_every
int order_flux
constexpr bool IS_LARGE_STRAINS
int atom_test
BoundaryNaturalBC::OpFlux< NaturalForceMeshsets, 1, SPACE_DIM > OpForce
constexpr int SPACE_DIM
NaturalBC< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS > DomainNaturalBCRhs
[Thermal problem]
double default_ref_temp
double default_heat_capacity
DomainNaturalBCRhs::OpFlux< NaturalMeshsetType< BLOCKSET >, 1, 1 > OpHeatSource
int order_disp
double default_young_modulus
[Essential boundary conditions (Least square approach)]
PetscBool do_output_skin
double default_coeff_expansion
int order_temp
DomainNaturalBCRhs::OpFlux< SetTargetTemperature, 1, 1 > OpSetTemperatureRhs
auto save_range
PetscBool do_output_domain
double default_heat_conductivity
NaturalBC< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS > DomainNaturalBCLhs
DomainNaturalBCRhs::OpFlux< NaturalMeshsetType< BLOCKSET >, 1, SPACE_DIM > OpBodyForce
double default_init_temp
double default_poisson_ratio