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

Implementation of mix-element for Large strains.

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 <cmath>
#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
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"));
LogManager::setLog("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::createSink(LogManager::getStrmSync(), "EPSYNC"));
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";
CHKERR DMRegister_MoFEM(dm_name);
DMType dm_name_mg = "DMMOFEM_MG";
CHKERR DMRegister_MGViaApproxOrders(dm_name_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
PetscBool fully_distributed = PETSC_FALSE;
CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-fully_distributed",
&fully_distributed, PETSC_NULL);
if (fully_distributed) {
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);
} else {
CHKERR CommInterface::loadFileRootProcAllRestDistributed(
}
// Create MoFEM database and link it to MoAB
MoFEM::Core mofem_core(moab);
MoFEM::Interface &m_field = mofem_core;
// Register mofem DM
CHKERR DMRegister_MoFEM("DMMOFEM");
BitRefLevel bit_level0 = BitRefLevel().set(0);
CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
0, 3, bit_level0);
// Data stuctures
EshelbianCore ep(m_field);
auto meshset_ptr = get_temp_meshset_ptr(moab);
CHKERR moab.add_entities(
*meshset_ptr, CommInterface::getPartEntities(moab, pcomm->rank()));
auto get_adj = [&](Range ents, int dim) {
Range adj;
CHKERR moab.get_adjacencies(ents, dim, false, adj,
moab::Interface::UNION);
return adj;
};
CHKERR moab.add_entities(
*meshset_ptr,
get_adj(CommInterface::getPartEntities(moab, pcomm->rank())
.subset_by_dimension(SPACE_DIM),
2));
CHKERR moab.add_entities(
*meshset_ptr,
get_adj(CommInterface::getPartEntities(moab, pcomm->rank())
.subset_by_dimension(SPACE_DIM),
1));
CHKERR moab.add_entities(
*meshset_ptr,
get_adj(CommInterface::getPartEntities(moab, pcomm->rank())
.subset_by_dimension(SPACE_DIM),
0));
// create growing crack surface meshset
CHKERR ep.createExchangeVectors(Sev::inform);
CHKERR ep.addFields(*meshset_ptr);
CHKERR ep.projectGeometry(*meshset_ptr);
CHKERR ep.addVolumeFiniteElement(*meshset_ptr);
CHKERR ep.addBoundaryFiniteElement(*meshset_ptr);
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) << "Neo-Hookean material model";
break;
default:
SETERRQ(m_field.get_comm(), MOFEM_DATA_INCONSISTENCY, "Unknown material");
break;
}
SmartPetscObj<Vec> x_elastic, f_elastic;
CHKERR DMCreateGlobalVector_MoFEM(ep.dmElastic, x_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(TSADAPTMOFEM, TSAdaptCreateMoFEM);
TSAdapt adapt;
CHKERR TSGetAdapt(ts_elastic, &adapt);
CHKERR TSAdaptSetType(adapt, TSADAPTNONE);
CHKERR TSSetTime(ts_elastic, time);
CHKERR ep.solveDynamicRelaxation(ts_elastic, x_elastic);
} else {
CHKERR ep.solveElastic(ts_elastic, x_elastic);
}
}
// finish work cleaning memory, getting statistics, etc.
#ifdef PYTHON_SDF
if (Py_FinalizeEx() < 0) {
exit(120);
}
#endif
}
Eshelbian plasticity interface.
#define TSADAPTMOFEM
Definition TsCtx.hpp:10
static char help[]
int main()
constexpr int SPACE_DIM
cholesky decomposition
#define CATCH_ERRORS
Catch errors.
#define MYPCOMM_INDEX
default communicator number PCOMM
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define CHKERR
Inline error check.
constexpr int adolc_tag
Definition ep.cpp:25
#define MU(E, NU)
Definition fem_tools.h:23
#define LAMBDA(E, NU)
Definition fem_tools.h:22
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
char mesh_file_name[255]
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
MoFEMErrorCode setElasticElementOps(const int tag)
MoFEMErrorCode addVolumeFiniteElement(const EntityHandle meshset=0)
MoFEMErrorCode addFields(const EntityHandle meshset=0)
MoFEMErrorCode createCrackSurfaceMeshset()
MoFEMErrorCode addBoundaryFiniteElement(const EntityHandle meshset=0)
MoFEMErrorCode getSpatialRotationBc()
MoFEMErrorCode addMaterial_Hencky(double E, double nu)
MoFEMErrorCode projectGeometry(const EntityHandle meshset=0)
MoFEMErrorCode setBlockTagsOnSkin()
MoFEMErrorCode addMaterial_HMHHStVenantKirchhoff(const int tape, const double lambda, const double mu, const double sigma_y)
MoFEMErrorCode addMaterial_HMHMooneyRivlin(const int tape, const double alpha, const double beta, const double lambda, const double sigma_y)
static PetscBool dynamicRelaxation
MoFEMErrorCode solveElastic(TS ts, Vec x)
MoFEMErrorCode setElasticElementToTs(DM dm)
MoFEMErrorCode solveDynamicRelaxation(TS ts, Vec x)
MoFEMErrorCode getSpatialTractionFreeBc(const EntityHandle meshset=0)
MoFEMErrorCode getSpatialTractionBc()
MoFEMErrorCode addDMs(const BitRefLevel bit=BitRefLevel().set(0), const EntityHandle meshset=0)
MoFEMErrorCode getSpatialDispBc()
MoFEMErrorCode createExchangeVectors(Sev sev)
MoFEMErrorCode addMaterial_HMHNeohookean(const int tape, const double c10, const double K)
SmartPetscObj< DM > dmElastic
Elastic problem.
Managing BitRefLevels.
virtual MPI_Comm & get_comm() const =0
Core (interface) class.
Definition Core.hpp:82
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition Core.cpp:72
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:112
Deprecated interface functions.
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.