31using namespace boost::multi_index;
32using namespace boost::multiprecision;
33using namespace boost::numeric;
38 #include <boost/python.hpp>
39 #include <boost/python/def.hpp>
40 #include <boost/python/numpy.hpp>
41namespace bp = boost::python;
42namespace np = boost::python::numpy;
48static char help[] =
"...\n\n";
50int main(
int argc,
char *argv[]) {
53 const char param_file[] =
"param_file.petsc";
57 auto core_log = logging::core::get();
65 MOFEM_LOG(
"EP", Sev::inform) <<
"Python initialised";
67 MOFEM_LOG(
"EP", Sev::inform) <<
"Python NOT initialised";
78 PetscBool flg = PETSC_TRUE;
84 char restart_file[255];
85 PetscBool restart_flg = PETSC_TRUE;
92 DMType dm_name =
"DMMOFEM";
94 DMType dm_name_mg =
"DMMOFEM_MG";
99 moab::Interface &moab = moab_core;
101 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab,
MYPCOMM_INDEX);
103 pcomm =
new ParallelComm(&moab, PETSC_COMM_WORLD);
105 PetscBool fully_distributed = PETSC_FALSE;
107 &fully_distributed, PETSC_NULL);
108 if (fully_distributed) {
110 if (pcomm->proc_config().proc_size() == 1)
113 option =
"PARALLEL=READ_PART;"
114 "PARALLEL_RESOLVE_SHARED_ENTS;"
115 "PARTITION=PARALLEL_PARTITION";
118 CHKERR CommInterface::loadFileRootProcAllRestDistributed(
138 *meshset_ptr, CommInterface::getPartEntities(moab, pcomm->rank()));
139 auto get_adj = [&](
Range ents,
int dim) {
141 CHKERR moab.get_adjacencies(ents, dim,
false, adj,
142 moab::Interface::UNION);
147 get_adj(CommInterface::getPartEntities(moab, pcomm->rank())
152 get_adj(CommInterface::getPartEntities(moab, pcomm->rank())
157 get_adj(CommInterface::getPartEntities(moab, pcomm->rank())
185 const char *list_materials[LastMaterial] = {
186 "stvenant_kirchhoff",
"mooney_rivlin",
"hencky",
"neo_hookean"};
187 PetscInt choice_material = MooneyRivlin;
189 LastMaterial, &choice_material, PETSC_NULL);
191 switch (choice_material) {
192 case StVenantKirchhoff:
193 MOFEM_LOG(
"EP", Sev::inform) <<
"StVenantKirchhoff material model";
198 MOFEM_LOG(
"EP", Sev::inform) <<
"MooneyRivlin material model";
199 MOFEM_LOG(
"EP", Sev::inform) <<
"MU(1, 0)/2 = " <<
MU(1, 0) / 2;
205 MOFEM_LOG(
"EP", Sev::inform) <<
"Hencky material model";
209 MOFEM_LOG(
"EP", Sev::inform) <<
"Neo-Hookean material model";
223 CHKERR VecSetOption(f_elastic, VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE);
227 CHKERR PetscViewerBinaryOpen(PETSC_COMM_WORLD, restart_file,
228 FILE_MODE_READ, &viewer);
229 CHKERR VecLoad(x_elastic, viewer);
230 CHKERR PetscViewerDestroy(&viewer);
235 auto ts_elastic =
createTS(PETSC_COMM_WORLD);
236 CHKERR TSSetType(ts_elastic, TSBEULER);
239 CHKERR TSGetAdapt(ts_elastic, &adapt);
240 CHKERR TSAdaptSetType(adapt, TSADAPTNONE);
241 CHKERR TSSetTime(ts_elastic, time);
256 if (Py_FinalizeEx() < 0) {
Eshelbian plasticity interface.
#define CATCH_ERRORS
Catch errors.
#define MYPCOMM_INDEX
default communicator number PCOMM
@ MOFEM_DATA_INCONSISTENCY
#define CHKERR
Inline error check.
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
MoFEMErrorCode DMRegister_MGViaApproxOrders(const char sname[])
Register DM for Multi-Grid via approximation orders.
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVector.
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
implementation of Data Operators for Forces and Sources
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)
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
auto createTS(MPI_Comm comm)
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
auto get_temp_meshset_ptr(moab::Interface &moab)
Create smart pointer to temporary meshset.
PetscErrorCode TSAdaptCreateMoFEM(TSAdapt adapt)
Craete MOFEM adapt.
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.
virtual MPI_Comm & get_comm() const =0
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Deprecated interface functions.
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.