v0.14.0
ep.cpp

Implementation of mix-element for Large strains

Todo:
Implementation of plasticity
/**
* \file ep.cpp
* \example ep.cpp
*
* \brief Implementation of mix-element for Large strains
*
* \todo Implementation of plasticity
*
*/
/* 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/>. */
constexpr int adolc_tag = 1;
#define SINGULARITY
#include <MoFEM.hpp>
using namespace MoFEM;
using namespace boost::multi_index;
using namespace boost::multiprecision;
using namespace boost::numeric;
#include <cholesky.hpp>
#ifdef PYTHON_SDF
#include <boost/python.hpp>
#include <boost/python/def.hpp>
#include <boost/python/numpy.hpp>
namespace bp = boost::python;
namespace np = boost::python::numpy;
#endif
#include <TSEPAdapt.hpp>
using namespace EshelbianPlasticity;
static char help[] = "...\n\n";
int main(int argc, char *argv[]) {
// initialize petsc
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(), "EP"));
MOFEM_LOG_TAG("EP", "ep");
#ifdef PYTHON_SDF
Py_Initialize();
np::initialize();
MOFEM_LOG("EP", Sev::inform) << "Python initialised";
#else
MOFEM_LOG("EP", Sev::inform) << "Python NOT initialised";
#endif
core_log->add_sink(
LogManager::setLog("EPSYNC");
MOFEM_LOG_TAG("EPSYNC", "ep");
try {
// Get mesh file
PetscBool flg = PETSC_TRUE;
char mesh_file_name[255];
CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
255, &flg);
CHKERR PetscOptionsGetString(PETSC_NULL, "", "-file_name", mesh_file_name,
255, &flg);
char restart_file[255];
PetscBool restart_flg = PETSC_TRUE;
CHKERR PetscOptionsGetString(PETSC_NULL, "", "-restart", restart_file, 255,
&restart_flg);
double time = 0;
CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-time", &time, PETSC_NULL);
// Register DM Manager
DMType dm_name = "DMMOFEM";
DMType dm_name_mg = "DMMOFEM_MG";
// Create MoAB database
moab::Core moab_core;
moab::Interface &moab = moab_core;
ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
if (pcomm == NULL)
pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
// Read mesh to MOAB
const char *option;
if (pcomm->proc_config().proc_size() == 1)
option = "";
else
option = "PARALLEL=READ_PART;"
"PARALLEL_RESOLVE_SHARED_ENTS;"
"PARTITION=PARALLEL_PARTITION";
CHKERR moab.load_file(mesh_file_name, 0, option);
// Create MoFEM database and link it to MoAB
MoFEM::Core mofem_core(moab);
MoFEM::Interface &m_field = mofem_core;
// Register mofem DM
BitRefLevel bit_level0 = BitRefLevel().set(0);
CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
0, 3, bit_level0);
// Data stuctures
EshelbianCore ep(m_field);
CHKERR ep.addDMs();
enum MaterialModel {
StVenantKirchhoff,
MooneyRivlin,
Hencky,
Neohookean,
LastMaterial
};
const char *list_materials[LastMaterial] = {
"stvenant_kirchhoff", "mooney_rivlin", "hencky", "neo_hookean"};
PetscInt choice_material = MooneyRivlin;
CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-material", list_materials,
LastMaterial, &choice_material, PETSC_NULL);
switch (choice_material) {
case StVenantKirchhoff:
MOFEM_LOG("EP", Sev::inform) << "StVenantKirchhoff material model";
MU(1, 0.25), 0);
break;
case MooneyRivlin:
MOFEM_LOG("EP", Sev::inform) << "MooneyRivlin material model";
MOFEM_LOG("EP", Sev::inform) << "MU(1, 0)/2 = " << MU(1, 0) / 2;
MOFEM_LOG("EP", Sev::inform) << "LAMBDA(1, 0) = " << LAMBDA(1, 0);
LAMBDA(1, 0), 0);
break;
case Hencky:
MOFEM_LOG("EP", Sev::inform) << "Hencky material model";
CHKERR ep.addMaterial_Hencky(5., 0.25);
break;
case Neohookean:
MOFEM_LOG("EP", Sev::inform) << "Neohookean material model";
break;
default:
SETERRQ(m_field.get_comm(), MOFEM_DATA_INCONSISTENCY, "Unknown material");
break;
}
SmartPetscObj<Vec> x_elastic, f_elastic;
f_elastic = vectorDuplicate(x_elastic);
CHKERR VecSetOption(f_elastic, VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE);
if (restart_flg) {
PetscViewer viewer;
CHKERR PetscViewerBinaryOpen(PETSC_COMM_WORLD, restart_file,
FILE_MODE_READ, &viewer);
CHKERR VecLoad(x_elastic, viewer);
CHKERR PetscViewerDestroy(&viewer);
CHKERR DMoFEMMeshToLocalVector(ep.dmElastic, x_elastic, INSERT_VALUES,
SCATTER_REVERSE);
}
auto ts_elastic = createTS(PETSC_COMM_WORLD);
CHKERR TSSetType(ts_elastic, TSBEULER);
CHKERR TSAdaptRegister(TSADAPTEP, TSAdaptCreate_EP);
TSAdapt adapt;
CHKERR TSGetAdapt(ts_elastic, &adapt);
CHKERR TSAdaptSetType(adapt, TSADAPTEP);
// CHKERR TSSetFromOptions(ts_elastic);
CHKERR TSSetTime(ts_elastic, time);
CHKERR ep.solveElastic(ts_elastic, x_elastic);
}
// finish work cleaning memory, getting statistics, etc.
#ifdef PYTHON_SDF
if (Py_FinalizeEx() < 0) {
exit(120);
}
#endif
}
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:228
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
EshelbianPlasticity::EshelbianCore::addDMs
MoFEMErrorCode addDMs(const BitRefLevel bit=BitRefLevel().set(0), const EntityHandle meshset=0)
Definition: EshelbianPlasticity.cpp:1101
EshelbianPlasticity::EshelbianCore::projectGeometry
MoFEMErrorCode projectGeometry(const EntityHandle meshset=0)
Definition: EshelbianPlasticity.cpp:707
MoFEM::createTS
auto createTS(MPI_Comm comm)
Definition: PetscSmartObj.hpp:249
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
EshelbianPlasticity::EshelbianCore::addMaterial_Hencky
MoFEMErrorCode addMaterial_Hencky(double E, double nu)
Definition: EshelbianADOL-C.cpp:491
EshelbianPlasticity::EshelbianCore::addBoundaryFiniteElement
MoFEMErrorCode addBoundaryFiniteElement(const EntityHandle meshset=0)
Definition: EshelbianPlasticity.cpp:912
EshelbianPlasticity::EshelbianCore::setElasticElementToTs
MoFEMErrorCode setElasticElementToTs(DM dm)
Definition: EshelbianPlasticity.cpp:2192
EshelbianPlasticity::EshelbianCore::setElasticElementOps
MoFEMErrorCode setElasticElementOps(const int tag)
Definition: EshelbianPlasticity.cpp:2169
LAMBDA
#define LAMBDA(E, NU)
Definition: fem_tools.h:22
EshelbianPlasticity
Definition: CGGTonsorialBubbleBase.hpp:11
MoFEM.hpp
MoFEM::DMoFEMMeshToLocalVector
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:523
EshelbianPlasticity::EshelbianCore::getSpatialTractionFreeBc
MoFEMErrorCode getSpatialTractionFreeBc(const EntityHandle meshset=0)
Definition: EshelbianPlasticity.hpp:1061
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::LogManager::createSink
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:298
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
adolc_tag
constexpr int adolc_tag
Definition: ep.cpp:25
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MoFEM::DMCreateGlobalVector_MoFEM
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVector.
Definition: DMMoFEM.cpp:1167
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
TSEPAdapt.hpp
Step adaptation.
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
EshelbianPlasticity::EshelbianCore::addVolumeFiniteElement
MoFEMErrorCode addVolumeFiniteElement(const EntityHandle meshset=0)
Definition: EshelbianPlasticity.cpp:850
EshelbianPlasticity::EshelbianCore::dmElastic
SmartPetscObj< DM > dmElastic
Elastic problem.
Definition: EshelbianPlasticity.hpp:964
MoFEM::DMRegister_MGViaApproxOrders
MoFEMErrorCode DMRegister_MGViaApproxOrders(const char sname[])
Register DM for Multi-Grid via approximation orders.
Definition: PCMGSetUpViaApproxOrders.cpp:302
TSAdaptCreate_EP
PetscErrorCode TSAdaptCreate_EP(TSAdapt adapt)
Definition: TSEPAdapt.hpp:87
EshelbianPlasticity::EshelbianCore::addMaterial_HMHHStVenantKirchhoff
MoFEMErrorCode addMaterial_HMHHStVenantKirchhoff(const int tape, const double lambda, const double mu, const double sigma_y)
Definition: EshelbianADOL-C.cpp:464
MoFEM::LogManager::getStrmSync
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
Definition: LogManager.cpp:348
MoFEM::LogManager::getStrmWorld
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
Definition: LogManager.cpp:344
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
EshelbianPlasticity::EshelbianCore::addMaterial_HMHNeohookean
MoFEMErrorCode addMaterial_HMHNeohookean(const int tape, const double c10, const double K)
Definition: EshelbianADOL-C.cpp:482
EshelbianPlasticity::EshelbianCore::solveElastic
MoFEMErrorCode solveElastic(TS ts, Vec x)
Definition: EshelbianPlasticity.cpp:2225
EshelbianPlasticity.hpp
Eshelbian plasticity interface.
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
cholesky.hpp
cholesky decomposition
EshelbianPlasticity::EshelbianCore
Definition: EshelbianPlasticity.hpp:894
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
help
static char help[]
Definition: ep.cpp:48
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::vectorDuplicate
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
Definition: PetscSmartObj.hpp:221
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::PetscOptionsGetEList
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
Definition: DeprecatedPetsc.hpp:203
EshelbianPlasticity::EshelbianCore::addFields
MoFEMErrorCode addFields(const EntityHandle meshset=0)
Definition: EshelbianPlasticity.cpp:483
TSADAPTEP
#define TSADAPTEP
Definition: TSEPAdapt.hpp:15
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
main
int main(int argc, char *argv[])
Definition: ep.cpp:50
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MoFEM::PetscOptionsGetScalar
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:162
MoFEM::LogManager::setLog
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:389
MoFEM::SmartPetscObj< Vec >
EshelbianPlasticity::EshelbianCore::getSpatialTractionBc
MoFEMErrorCode getSpatialTractionBc()
Definition: EshelbianPlasticity.cpp:2599
EshelbianPlasticity::EshelbianCore::addMaterial_HMHMooneyRivlin
MoFEMErrorCode addMaterial_HMHMooneyRivlin(const int tape, const double alpha, const double beta, const double lambda, const double sigma_y)
Definition: EshelbianADOL-C.cpp:473
EshelbianPlasticity::EshelbianCore::getSpatialRotationBc
MoFEMErrorCode getSpatialRotationBc()
Definition: EshelbianPlasticity.hpp:1039
EshelbianPlasticity::EshelbianCore::getSpatialDispBc
MoFEMErrorCode getSpatialDispBc()
[Getting norms]
Definition: EshelbianPlasticity.cpp:2561
MU
#define MU(E, NU)
Definition: fem_tools.h:23