v0.14.0
Loading...
Searching...
No Matches
ep.cpp File Reference
#include <MoFEM.hpp>
#include <cmath>
#include <cholesky.hpp>
#include <EshelbianPlasticity.hpp>

Go to the source code of this file.

Macros

#define SINGULARITY
 

Functions

int main (int argc, char *argv[])
 

Variables

constexpr int adolc_tag = 1
 
static char help [] = "...\n\n"
 

Macro Definition Documentation

◆ SINGULARITY

#define SINGULARITY

Definition at line 27 of file ep.cpp.

Function Documentation

◆ main()

int main ( int argc,
char * argv[] )

Definition at line 50 of file ep.cpp.

50 {
51
52 // initialize petsc
53 const char param_file[] = "param_file.petsc";
54 MoFEM::Core::Initialize(&argc, &argv, param_file, help);
55
56 // Add logging channel for example
57 auto core_log = logging::core::get();
58 core_log->add_sink(LogManager::createSink(LogManager::getStrmWorld(), "EP"));
60 MOFEM_LOG_TAG("EP", "ep");
61
62#ifdef PYTHON_SDF
63 Py_Initialize();
64 np::initialize();
65 MOFEM_LOG("EP", Sev::inform) << "Python initialised";
66#else
67 MOFEM_LOG("EP", Sev::inform) << "Python NOT initialised";
68#endif
69
70 core_log->add_sink(
72 LogManager::setLog("EPSYNC");
73 MOFEM_LOG_TAG("EPSYNC", "ep");
74
75 try {
76
77 // Get mesh file
78 PetscBool flg = PETSC_TRUE;
79 char mesh_file_name[255];
80 CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
81 255, &flg);
82 CHKERR PetscOptionsGetString(PETSC_NULL, "", "-file_name", mesh_file_name,
83 255, &flg);
84 char restart_file[255];
85 PetscBool restart_flg = PETSC_TRUE;
86 CHKERR PetscOptionsGetString(PETSC_NULL, "", "-restart", restart_file, 255,
87 &restart_flg);
88 double time = 0;
89 CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-time", &time, PETSC_NULL);
90
91 // Register DM Manager
92 DMType dm_name = "DMMOFEM";
93 CHKERR DMRegister_MoFEM(dm_name);
94 DMType dm_name_mg = "DMMOFEM_MG";
96
97 // Create MoAB database
98 moab::Core moab_core;
99 moab::Interface &moab = moab_core;
100
101 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
102 if (pcomm == NULL)
103 pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
104 // Read mesh to MOAB
105 PetscBool fully_distributed = PETSC_FALSE;
106 CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-fully_distributed",
107 &fully_distributed, PETSC_NULL);
108 if (fully_distributed) {
109 const char *option;
110 if (pcomm->proc_config().proc_size() == 1)
111 option = "";
112 else
113 option = "PARALLEL=READ_PART;"
114 "PARALLEL_RESOLVE_SHARED_ENTS;"
115 "PARTITION=PARALLEL_PARTITION";
116 CHKERR moab.load_file(mesh_file_name, 0, option);
117 } else {
118 CHKERR CommInterface::loadFileRootProcAllRestDistributed(
120 }
121
122 // Create MoFEM database and link it to MoAB
123 MoFEM::Core mofem_core(moab);
124 MoFEM::Interface &m_field = mofem_core;
125
126 // Register mofem DM
127 CHKERR DMRegister_MoFEM("DMMOFEM");
128
129 BitRefLevel bit_level0 = BitRefLevel().set(0);
130 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
131 0, 3, bit_level0);
132
133 // Data stuctures
134 EshelbianCore ep(m_field);
135
136 auto meshset_ptr = get_temp_meshset_ptr(moab);
137 CHKERR moab.add_entities(
138 *meshset_ptr, CommInterface::getPartEntities(moab, pcomm->rank()));
139 auto get_adj = [&](Range ents, int dim) {
140 Range adj;
141 CHKERR moab.get_adjacencies(ents, dim, false, adj,
142 moab::Interface::UNION);
143 return adj;
144 };
145 CHKERR moab.add_entities(
146 *meshset_ptr,
147 get_adj(CommInterface::getPartEntities(moab, pcomm->rank())
148 .subset_by_dimension(SPACE_DIM),
149 2));
150 CHKERR moab.add_entities(
151 *meshset_ptr,
152 get_adj(CommInterface::getPartEntities(moab, pcomm->rank())
153 .subset_by_dimension(SPACE_DIM),
154 1));
155 CHKERR moab.add_entities(
156 *meshset_ptr,
157 get_adj(CommInterface::getPartEntities(moab, pcomm->rank())
158 .subset_by_dimension(SPACE_DIM),
159 0));
160
161 // create growing crack surface meshset
162 CHKERR ep.createCrackSurfaceMeshset();
163
164 CHKERR ep.getSpatialDispBc();
165 CHKERR ep.getSpatialRotationBc();
166 CHKERR ep.getSpatialTractionBc();
167 CHKERR ep.getSpatialTractionFreeBc();
168 CHKERR ep.setBlockTagsOnSkin();
169
170 CHKERR ep.createExchangeVectors(Sev::inform);
171
172 CHKERR ep.addFields(*meshset_ptr);
173 CHKERR ep.projectGeometry(*meshset_ptr);
174 CHKERR ep.addVolumeFiniteElement(*meshset_ptr);
175 CHKERR ep.addBoundaryFiniteElement(*meshset_ptr);
176 CHKERR ep.addDMs();
177
178 enum MaterialModel {
179 StVenantKirchhoff,
180 MooneyRivlin,
181 Hencky,
182 Neohookean,
183 LastMaterial
184 };
185 const char *list_materials[LastMaterial] = {
186 "stvenant_kirchhoff", "mooney_rivlin", "hencky", "neo_hookean"};
187 PetscInt choice_material = MooneyRivlin;
188 CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-material", list_materials,
189 LastMaterial, &choice_material, PETSC_NULL);
190
191 switch (choice_material) {
192 case StVenantKirchhoff:
193 MOFEM_LOG("EP", Sev::inform) << "StVenantKirchhoff material model";
194 CHKERR ep.addMaterial_HMHHStVenantKirchhoff(adolc_tag, LAMBDA(1, 0.25),
195 MU(1, 0.25), 0);
196 break;
197 case MooneyRivlin:
198 MOFEM_LOG("EP", Sev::inform) << "MooneyRivlin material model";
199 MOFEM_LOG("EP", Sev::inform) << "MU(1, 0)/2 = " << MU(1, 0) / 2;
200 MOFEM_LOG("EP", Sev::inform) << "LAMBDA(1, 0) = " << LAMBDA(1, 0);
201 CHKERR ep.addMaterial_HMHMooneyRivlin(adolc_tag, MU(1, 0) / 2., 0,
202 LAMBDA(1, 0), 0);
203 break;
204 case Hencky:
205 MOFEM_LOG("EP", Sev::inform) << "Hencky material model";
206 CHKERR ep.addMaterial_Hencky(5., 0.25);
207 break;
208 case Neohookean:
209 MOFEM_LOG("EP", Sev::inform) << "Neo-Hookean material model";
210 CHKERR ep.addMaterial_HMHNeohookean(adolc_tag, 1.0, 1.0);
211 break;
212 default:
213 SETERRQ(m_field.get_comm(), MOFEM_DATA_INCONSISTENCY, "Unknown material");
214 break;
215 }
216
217 CHKERR ep.setElasticElementOps(adolc_tag);
218 CHKERR ep.setElasticElementToTs(ep.dmElastic);
219
220 SmartPetscObj<Vec> x_elastic, f_elastic;
221 CHKERR DMCreateGlobalVector_MoFEM(ep.dmElastic, x_elastic);
222 f_elastic = vectorDuplicate(x_elastic);
223 CHKERR VecSetOption(f_elastic, VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE);
224
225 if (restart_flg) {
226 PetscViewer viewer;
227 CHKERR PetscViewerBinaryOpen(PETSC_COMM_WORLD, restart_file,
228 FILE_MODE_READ, &viewer);
229 CHKERR VecLoad(x_elastic, viewer);
230 CHKERR PetscViewerDestroy(&viewer);
231 CHKERR DMoFEMMeshToLocalVector(ep.dmElastic, x_elastic, INSERT_VALUES,
232 SCATTER_REVERSE);
233 }
234
235 auto ts_elastic = createTS(PETSC_COMM_WORLD);
236 CHKERR TSSetType(ts_elastic, TSBEULER);
237 CHKERR TSAdaptRegister(TSADAPTMOFEM, TSAdaptCreateMoFEM);
238 TSAdapt adapt;
239 CHKERR TSGetAdapt(ts_elastic, &adapt);
240 CHKERR TSAdaptSetType(adapt, TSADAPTNONE);
241 CHKERR TSSetTime(ts_elastic, time);
242
243 if(ep.dynamicRelaxation) {
244 CHKERR ep.solveDynamicRelaxation(ts_elastic, x_elastic);
245 } else {
246 CHKERR ep.solveElastic(ts_elastic, x_elastic);
247 }
248
249 }
251
252 // finish work cleaning memory, getting statistics, etc.
254
255#ifdef PYTHON_SDF
256 if (Py_FinalizeEx() < 0) {
257 exit(120);
258 }
259#endif
260}
#define TSADAPTMOFEM
Definition TsCtx.hpp:10
constexpr int SPACE_DIM
#define CATCH_ERRORS
Catch errors.
#define MYPCOMM_INDEX
default communicator number PCOMM
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define CHKERR
Inline error check.
static char help[]
Definition ep.cpp:48
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
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
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
MoFEMErrorCode DMRegister_MGViaApproxOrders(const char sname[])
Register DM for Multi-Grid via approximation orders.
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVector.
Definition DMMoFEM.cpp:1167
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.
char mesh_file_name[255]
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
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.
Definition TsCtx.cpp:797
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.
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.

Variable Documentation

◆ adolc_tag

int adolc_tag = 1
constexpr
Examples
ep.cpp.

Definition at line 25 of file ep.cpp.

◆ help

char help[] = "...\n\n"
static

Definition at line 48 of file ep.cpp.