v0.14.0
ep.cpp
Go to the documentation of this file.
1 /**
2  * \file ep.cpp
3  * \example ep.cpp
4  *
5  * \brief Implementation of mix-element for Large strains
6  *
7  * \todo Implementation of plasticity
8  *
9  */
10 
11 /* This file is part of MoFEM.
12  * MoFEM is free software: you can redistribute it and/or modify it under
13  * the terms of the GNU Lesser General Public License as published by the
14  * Free Software Foundation, either version 3 of the License, or (at your
15  * option) any later version.
16  *
17  * MoFEM is distributed in the hope that it will be useful, but WITHOUT
18  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
19  * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
20  * License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with MoFEM. If not, see <http://www.gnu.org/licenses/>. */
24 
25 
26 constexpr int adolc_tag = 1;
27 
28 #include <cholesky.hpp>
29 #ifdef PYTHON_SDF
30 #include <boost/python.hpp>
31 #include <boost/python/def.hpp>
32 #include <boost/python/def.hpp>
33 #include <boost/python/numpy.hpp>
34 namespace bp = boost::python;
35 namespace np = boost::python::numpy;
36 #endif
37 
38 #include <BasicFiniteElements.hpp>
39 #include <EshelbianPlasticity.hpp>
40 #include <TSEPAdapt.hpp>
41 
42 using namespace MoFEM;
43 using namespace EshelbianPlasticity;
44 
45 static char help[] = "...\n\n";
46 
47 int main(int argc, char *argv[]) {
48 
49  // initialize petsc
50  const char param_file[] = "param_file.petsc";
51  MoFEM::Core::Initialize(&argc, &argv, param_file, help);
52 
53  // Add logging channel for example
54  auto core_log = logging::core::get();
55  core_log->add_sink(LogManager::createSink(LogManager::getStrmWorld(), "EP"));
56  LogManager::setLog("EP");
57  MOFEM_LOG_TAG("EP", "ep");
58 
59 #ifdef PYTHON_SDF
60  Py_Initialize();
61  np::initialize();
62  MOFEM_LOG("EP", Sev::inform) << "Python initialised";
63 #else
64  MOFEM_LOG("EP", Sev::inform) << "Python NOT initialised";
65 #endif
66 
67 
68 
69  core_log->add_sink(
71  LogManager::setLog("EPSYNC");
72  MOFEM_LOG_TAG("EPSYNC", "ep");
73 
74  try {
75 
76  // Get mesh file
77  PetscBool flg = PETSC_TRUE;
78  char mesh_file_name[255];
79  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
80  255, &flg);
81  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-file_name", mesh_file_name,
82  255, &flg);
83  char restart_file[255];
84  PetscBool restart_flg = PETSC_TRUE;
85  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-restart", restart_file, 255,
86  &restart_flg);
87  double time = 0;
88  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-time", &time, PETSC_NULL);
89 
90  // Register DM Manager
91  DMType dm_name = "DMMOFEM";
92  CHKERR DMRegister_MoFEM(dm_name);
93 
94  // Create MoAB database
95  moab::Core moab_core;
96  moab::Interface &moab = moab_core;
97 
98  ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
99  if (pcomm == NULL)
100  pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
101  // Read mesh to MOAB
102  const char *option;
103  if(pcomm->proc_config().proc_size() == 1)
104  option = "";
105  else
106  option = "PARALLEL=READ_PART;"
107  "PARALLEL_RESOLVE_SHARED_ENTS;"
108  "PARTITION=PARALLEL_PARTITION";
109  CHKERR moab.load_file(mesh_file_name, 0, option);
110 
111  // Create MoFEM database and link it to MoAB
112  MoFEM::Core mofem_core(moab);
113  MoFEM::Interface &m_field = mofem_core;
114 
115  // Register mofem DM
116  CHKERR DMRegister_MoFEM("DMMOFEM");
117 
118  BitRefLevel bit_level0 = BitRefLevel().set(0);
119  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
120  0, 3, bit_level0);
121 
122  // Data stuctures
123  EshelbianCore ep(m_field);
124 
129 
130  CHKERR ep.addFields();
133  CHKERR ep.addDMs();
134 
135  enum MaterialModel { StVenantKirchhoff, MooneyRivlin, Hencky, LastMaterial };
136  const char *list_materials[LastMaterial] = {"stvenant_kirchhoff",
137  "mooney_rivlin", "hencky"};
138  PetscInt choice_material = MooneyRivlin;
139  CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-material", list_materials,
140  LastMaterial, &choice_material, PETSC_NULL);
141 
142  switch (choice_material) {
143  case StVenantKirchhoff:
144  MOFEM_LOG("EP", Sev::inform)
145  << "StVenantKirchhoff material model";
147  MU(1, 0.25), 0);
148  break;
149  case MooneyRivlin:
150  MOFEM_LOG("EP", Sev::inform) << "MooneyRivlin material model";
152  LAMBDA(1, 0), 0);
153  break;
154  case Hencky:
155  MOFEM_LOG("EP", Sev::inform) << "Hencky material model";
156  CHKERR ep.addMaterial_Hencky(5., 0.25);
157  break;
158  default:
159  SETERRQ(m_field.get_comm(), MOFEM_DATA_INCONSISTENCY, "Unknown material");
160  break;
161  }
162 
165 
166  SmartPetscObj<Vec> x_elastic, f_elastic;
168  f_elastic = vectorDuplicate(x_elastic);
169  CHKERR VecSetOption(f_elastic, VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE);
170  SmartPetscObj<Mat> m_elastic;
171  CHKERR DMCreateMatrix_MoFEM(ep.dmElastic, m_elastic);
172 
173  if (restart_flg) {
174  PetscViewer viewer;
175  CHKERR PetscViewerBinaryOpen(PETSC_COMM_WORLD, restart_file,
176  FILE_MODE_READ, &viewer);
177  CHKERR VecLoad(x_elastic, viewer);
178  CHKERR PetscViewerDestroy(&viewer);
179  CHKERR DMoFEMMeshToLocalVector(ep.dmElastic, x_elastic, INSERT_VALUES,
180  SCATTER_REVERSE);
181  }
182 
183  auto ts_elastic = createTS(PETSC_COMM_WORLD);
184  CHKERR TSSetType(ts_elastic, TSBEULER);
185  CHKERR TSAdaptRegister(TSADAPTEP, TSAdaptCreate_EP);
186  TSAdapt adapt;
187  CHKERR TSGetAdapt(ts_elastic, &adapt);
188  CHKERR TSAdaptSetType(adapt, TSADAPTEP);
189  // CHKERR TSSetFromOptions(ts_elastic);
190 
191  CHKERR TSSetTime(ts_elastic, time);
192  CHKERR ep.solveElastic(ts_elastic, m_elastic, f_elastic, x_elastic);
193  }
194  CATCH_ERRORS;
195 
196  // finish work cleaning memory, getting statistics, etc.
198 
199 #ifdef PYTHON_SDF
200  if (Py_FinalizeEx() < 0) {
201  exit(120);
202  }
203 #endif
204 }
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
EshelbianPlasticity::EshelbianCore::addDMs
MoFEMErrorCode addDMs(const BitRefLevel bit=BitRefLevel().set(0))
Definition: EshelbianPlasticity.cpp:523
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:215
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
MoFEM::createTS
auto createTS(MPI_Comm comm)
Definition: PetscSmartObj.hpp:245
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:643
EshelbianPlasticity::EshelbianCore::addBoundaryFiniteElement
MoFEMErrorCode addBoundaryFiniteElement(const EntityHandle meshset=0)
Definition: EshelbianPlasticity.cpp:436
EshelbianPlasticity::EshelbianCore::setElasticElementToTs
MoFEMErrorCode setElasticElementToTs(DM dm)
Definition: EshelbianPlasticity.cpp:1265
EshelbianPlasticity::EshelbianCore::setElasticElementOps
MoFEMErrorCode setElasticElementOps(const int tag)
Definition: EshelbianPlasticity.cpp:1239
LAMBDA
#define LAMBDA(E, NU)
Definition: fem_tools.h:22
EshelbianPlasticity
Definition: CGGTonsorialBubbleBase.hpp:11
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:527
BasicFiniteElements.hpp
EshelbianPlasticity::EshelbianCore::getSpatialTractionFreeBc
MoFEMErrorCode getSpatialTractionFreeBc(const EntityHandle meshset=0)
Definition: EshelbianPlasticity.hpp:1021
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:26
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
MoFEM::DMCreateGlobalVector_MoFEM
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVector.
Definition: DMMoFEM.cpp:1171
MoFEM::DMCreateMatrix_MoFEM
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
Definition: DMMoFEM.cpp:1201
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
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:396
EshelbianPlasticity::EshelbianCore::dmElastic
SmartPetscObj< DM > dmElastic
Elastic problem.
Definition: EshelbianPlasticity.hpp:927
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:446
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:47
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:23
cholesky.hpp
cholesky decomposition
EshelbianPlasticity::EshelbianCore
Definition: EshelbianPlasticity.hpp:862
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:45
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:217
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
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:256
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:47
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:1648
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:455
EshelbianPlasticity::EshelbianCore::solveElastic
MoFEMErrorCode solveElastic(TS ts, Mat m, Vec f, Vec x)
Definition: EshelbianPlasticity.cpp:1309
EshelbianPlasticity::EshelbianCore::getSpatialRotationBc
MoFEMErrorCode getSpatialRotationBc()
Definition: EshelbianPlasticity.hpp:999
EshelbianPlasticity::EshelbianCore::getSpatialDispBc
MoFEMErrorCode getSpatialDispBc()
[Getting norms]
Definition: EshelbianPlasticity.cpp:1610
MU
#define MU(E, NU)
Definition: fem_tools.h:23