v0.15.5
Loading...
Searching...
No Matches
nonlinear_elastic.cpp
Go to the documentation of this file.
1/**
2 * \file nonlinear_elastic.cpp
3 * \example nonlinear_elastic.cpp
4 *
5 * Plane stress elastic dynamic problem
6 *
7 */
8
9#include <MoFEM.hpp>
10#include <MatrixFunction.hpp>
11
12using namespace MoFEM;
13
14constexpr int SPACE_DIM =
15 EXECUTABLE_DIMENSION; //< Space dimension of problem, mesh
16
21using DomainEleOp = DomainEle::UserDataOperator;
22using BoundaryEleOp = BoundaryEle::UserDataOperator;
23
28
32
33template <int DIM> struct PostProcEleByDim;
34
35template <> struct PostProcEleByDim<2> {
39};
40
41template <> struct PostProcEleByDim<3> {
45};
46
51
52#include <HenckyOps.hpp>
53#include <AdolCOps.hpp>
54using namespace HenckyOps;
55
56struct Example {
57
58 Example(MoFEM::Interface &m_field) : mField(m_field) {}
59
61
62private:
64
73};
74
75//! [Run problem]
87}
88//! [Run problem]
89
90//! [Read mesh]
95 char meshFileName[255];
96 CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-file_name",
97 meshFileName, 255, PETSC_NULL);
98 CHKERR simple->loadFile("", meshFileName,
101}
102//! [Read mesh]
103
104//! [Set up problem]
108
109 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOPT };
110 const char *list_bases[LASBASETOPT] = {"ainsworth", "demkowicz"};
111 PetscInt choice_base_value = AINSWORTH;
112 CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
113 LASBASETOPT, &choice_base_value, PETSC_NULL);
114
116 switch (choice_base_value) {
117 case AINSWORTH:
119 MOFEM_LOG("WORLD", Sev::inform)
120 << "Set AINSWORTH_LEGENDRE_BASE for displacements";
121 break;
122 case DEMKOWICZ:
124 MOFEM_LOG("WORLD", Sev::inform)
125 << "Set DEMKOWICZ_JACOBI_BASE for displacements";
126 break;
127 default:
128 base = LASTBASE;
129 break;
130 }
131
132 // Add field
133 CHKERR simple->addDomainField("U", H1, base, SPACE_DIM);
134 CHKERR simple->addBoundaryField("U", H1, base, SPACE_DIM);
135 int order = 2;
136 CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
137 CHKERR simple->setFieldOrder("U", order);
138 CHKERR simple->setUp();
140}
141//! [Set up problem]
142
143//! [Boundary condition]
146 auto *pipeline_mng = mField.getInterface<PipelineManager>();
148 auto bc_mng = mField.getInterface<BcManager>();
149 auto time_scale = boost::make_shared<TimeScale>();
150
151 auto integration_rule = [](int, int, int approx_order) {
152 return 2 * (approx_order - 1);
153 };
154
155 CHKERR pipeline_mng->setDomainRhsIntegrationRule(integration_rule);
156 CHKERR pipeline_mng->setDomainLhsIntegrationRule(integration_rule);
157 CHKERR pipeline_mng->setBoundaryRhsIntegrationRule(integration_rule);
158
160 pipeline_mng->getOpBoundaryRhsPipeline(), mField, "U", {time_scale},
161 "FORCE", Sev::inform);
162
163 //! [Define gravity vector]
165 pipeline_mng->getOpDomainRhsPipeline(), mField, "U", {time_scale},
166 "BODY_FORCE", Sev::inform);
167
168 // Essential BC
169 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_X",
170 "U", 0, 0);
171 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_Y",
172 "U", 1, 1);
173 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_Z",
174 "U", 2, 2);
175 CHKERR bc_mng->pushMarkDOFsOnEntities<DisplacementCubitBcData>(
176 simple->getProblemName(), "U");
177
178 // Essential MPCs BC
179 CHKERR bc_mng->addBlockDOFsToMPCs(simple->getProblemName(), "U");
180
182}
183//! [Boundary condition]
184
185//! [Push operators to pipeline]
188 auto *simple = mField.getInterface<Simple>();
189 auto *pipeline_mng = mField.getInterface<PipelineManager>();
190
191 auto add_domain_ops_lhs = [&](auto &pip) {
194 CHKERR HenckyOps::opFactoryDomainLhs<SPACE_DIM, PETSC, GAUSS, DomainEleOp>(
195 mField, pip, "U", "MAT_ELASTIC", Sev::inform);
197 };
198
199 auto add_domain_ops_rhs = [&](auto &pip) {
202 CHKERR HenckyOps::opFactoryDomainRhs<SPACE_DIM, PETSC, GAUSS, DomainEleOp>(
203 mField, pip, "U", "MAT_ELASTIC", Sev::inform);
205 };
206
207 CHKERR add_domain_ops_lhs(pipeline_mng->getOpDomainLhsPipeline());
208 CHKERR add_domain_ops_rhs(pipeline_mng->getOpDomainRhsPipeline());
209
211}
212//! [Push operators to pipeline]
213
214/**
215 * @brief Monitor solution
216 *
217 * This functions is called by TS solver at the end of each step. It is used
218 * to output results to the hard drive.
219 */
220struct Monitor : public FEMethod {
221 Monitor(SmartPetscObj<DM> dm, boost::shared_ptr<PostProcEleBdy> post_proc)
222 : dM(dm), postProc(post_proc){};
225 constexpr int save_every_nth_step = 1;
226 if (ts_step % save_every_nth_step == 0) {
229 postProc->mField, postProc->getPostProcMesh(), {"U"});
230
231 CHKERR postProc->writeFile(
232 "out_step_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
233 }
235 }
236
237private:
239 boost::shared_ptr<PostProcEleBdy> postProc;
240};
241
242//! [Solve]
245 auto *simple = mField.getInterface<Simple>();
246 auto *pipeline_mng = mField.getInterface<PipelineManager>();
247
248 auto dm = simple->getDM();
249 auto ts = pipeline_mng->createTSIM();
250
251 // Setup postprocessing
252 auto create_post_proc_fe = [dm, this, simple]() {
253 auto post_proc_ele_domain = [dm, this](auto &pip_domain) {
255 auto common_ptr = commonDataFactory<SPACE_DIM, GAUSS, DomainEleOp>(
256 mField, pip_domain, "U", "MAT_ELASTIC", Sev::inform);
257 return common_ptr;
258 };
259
260 auto post_proc_fe_bdy = boost::make_shared<PostProcEleBdy>(mField);
261 auto u_ptr = boost::make_shared<MatrixDouble>();
262 post_proc_fe_bdy->getOpPtrVector().push_back(
264 auto op_loop_side =
265 new OpLoopSide<SideEle>(mField, simple->getDomainFEName(), SPACE_DIM);
266 auto common_ptr = post_proc_ele_domain(op_loop_side->getOpPtrVector());
267 post_proc_fe_bdy->getOpPtrVector().push_back(op_loop_side);
268
270
271 post_proc_fe_bdy->getOpPtrVector().push_back(
272
273 new OpPPMap(
274
275 post_proc_fe_bdy->getPostProcMesh(),
276 post_proc_fe_bdy->getMapGaussPts(),
277
278 {},
279
280 {{"U", u_ptr}},
281
282 {{"GRAD", common_ptr->matGradPtr},
283 {"FIRST_PIOLA", common_ptr->getMatFirstPiolaStress()}},
284
285 {}
286
287 )
288
289 );
290 return post_proc_fe_bdy;
291 };
292
293 auto add_extra_finite_elements_to_ksp_solver_pipelines = [&]() {
295
296 auto pre_proc_ptr = boost::make_shared<FEMethod>();
297 auto post_proc_rhs_ptr = boost::make_shared<FEMethod>();
298 auto post_proc_lhs_ptr = boost::make_shared<FEMethod>();
299
300 auto time_scale = boost::make_shared<TimeScale>();
301
302 auto get_bc_hook_rhs = [this, pre_proc_ptr, time_scale]() {
305 {time_scale}, false)();
307 };
308
309 pre_proc_ptr->preProcessHook = get_bc_hook_rhs;
310
311 auto get_post_proc_hook_rhs = [this, post_proc_rhs_ptr]() {
314 mField, post_proc_rhs_ptr, nullptr, Sev::verbose)();
316 mField, post_proc_rhs_ptr, 1.)();
317 CHKERR EssentialPostProcRhs<MPCsType>(mField, post_proc_rhs_ptr)();
319 };
320 auto get_post_proc_hook_lhs = [this, post_proc_lhs_ptr]() {
323 mField, post_proc_lhs_ptr, 1.)();
324 CHKERR EssentialPostProcLhs<MPCsType>(mField, post_proc_lhs_ptr)();
326 };
327 post_proc_rhs_ptr->postProcessHook = get_post_proc_hook_rhs;
328 post_proc_lhs_ptr->postProcessHook = get_post_proc_hook_lhs;
329
330 // This is low level pushing finite elements (pipelines) to solver
331 auto ts_ctx_ptr = getDMTsCtx(simple->getDM());
332 ts_ctx_ptr->getPreProcessIFunction().push_front(pre_proc_ptr);
333 ts_ctx_ptr->getPreProcessIJacobian().push_front(pre_proc_ptr);
334 ts_ctx_ptr->getPostProcessIFunction().push_back(post_proc_rhs_ptr);
335 ts_ctx_ptr->getPostProcessIJacobian().push_back(post_proc_lhs_ptr);
337 };
338
339 // Add extra finite elements to SNES solver pipelines to resolve essential
340 // boundary conditions
341 CHKERR add_extra_finite_elements_to_ksp_solver_pipelines();
342
343 auto create_monitor_fe = [dm](auto &&post_proc_fe) {
344 return boost::make_shared<Monitor>(dm, post_proc_fe);
345 };
346
347 // Set monitor which postprocessing results and saves them to the hard drive
348 boost::shared_ptr<FEMethod> null_fe;
349 auto monitor_ptr = create_monitor_fe(create_post_proc_fe());
350 CHKERR DMMoFEMTSSetMonitor(dm, ts, simple->getDomainFEName(), null_fe,
351 null_fe, monitor_ptr);
352
353 // Set time solver
354 double ftime = 1;
355 CHKERR TSSetDuration(ts, PETSC_DEFAULT, ftime);
356 CHKERR TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);
357
358 auto D = createDMVector(simple->getDM());
359
360 CHKERR TSSetSolution(ts, D);
361 CHKERR TSSetFromOptions(ts);
362
363 CHKERR TSSolve(ts, NULL);
364 CHKERR TSGetTime(ts, &ftime);
365
366 PetscInt steps, snesfails, rejects, nonlinits, linits;
367 CHKERR TSGetStepNumber(ts, &steps);
368 CHKERR TSGetSNESFailures(ts, &snesfails);
369 CHKERR TSGetStepRejections(ts, &rejects);
370 CHKERR TSGetSNESIterations(ts, &nonlinits);
371 CHKERR TSGetKSPIterations(ts, &linits);
372 MOFEM_LOG_C("EXAMPLE", Sev::inform,
373 "steps %d (%d rejected, %d SNES fails), ftime %g, nonlinits "
374 "%d, linits %d",
375 steps, rejects, snesfails, ftime, nonlinits, linits);
376
378}
379//! [Solve]
380
381//! [Getting norms]
384
386 auto dm = simple->getDM();
387
388 auto T = createDMVector(simple->getDM());
389 CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
390 SCATTER_FORWARD);
391 double nrm2;
392 CHKERR VecNorm(T, NORM_2, &nrm2);
393 MOFEM_LOG("EXAMPLE", Sev::inform) << "Solution norm " << nrm2;
394
395 auto post_proc_norm_fe = boost::make_shared<DomainEle>(mField);
396
397 auto post_proc_norm_rule_hook = [](int, int, int p) -> int { return 2 * p; };
398 post_proc_norm_fe->getRuleHook = post_proc_norm_rule_hook;
399
401 post_proc_norm_fe->getOpPtrVector(), {H1});
402
403 enum NORMS { U_NORM_L2 = 0, PIOLA_NORM, LAST_NORM };
404 auto norms_vec =
406 (mField.get_comm_rank() == 0) ? LAST_NORM : 0, LAST_NORM);
407 CHKERR VecZeroEntries(norms_vec);
408
409 auto u_ptr = boost::make_shared<MatrixDouble>();
410 post_proc_norm_fe->getOpPtrVector().push_back(
412
413 post_proc_norm_fe->getOpPtrVector().push_back(
414 new OpCalcNormL2Tensor1<SPACE_DIM>(u_ptr, norms_vec, U_NORM_L2));
415
416 auto common_ptr = commonDataFactory<SPACE_DIM, GAUSS, DomainEleOp>(
417 mField, post_proc_norm_fe->getOpPtrVector(), "U", "MAT_ELASTIC",
418 Sev::inform);
419
420 post_proc_norm_fe->getOpPtrVector().push_back(
422 common_ptr->getMatFirstPiolaStress(), norms_vec, PIOLA_NORM));
423
425 post_proc_norm_fe);
426
427 CHKERR VecAssemblyBegin(norms_vec);
428 CHKERR VecAssemblyEnd(norms_vec);
429
430 MOFEM_LOG_CHANNEL("SELF"); // Clear channel from old tags
431 if (mField.get_comm_rank() == 0) {
432 const double *norms;
433 CHKERR VecGetArrayRead(norms_vec, &norms);
434 MOFEM_TAG_AND_LOG("SELF", Sev::inform, "example")
435 << "norm_u: " << std::scientific << std::sqrt(norms[U_NORM_L2]);
436 MOFEM_TAG_AND_LOG("SELF", Sev::inform, "example")
437 << "norm_piola: " << std::scientific << std::sqrt(norms[PIOLA_NORM]);
438 CHKERR VecRestoreArrayRead(norms_vec, &norms);
439 }
440
442}
443//! [Getting norms]
444
445//! [Postprocessing results]
448 PetscInt test_nb = 0;
449 PetscBool test_flg = PETSC_FALSE;
450 CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-test", &test_nb, &test_flg);
451
452 if (test_flg) {
454 auto T = createDMVector(simple->getDM());
455 CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
456 SCATTER_FORWARD);
457 double nrm2;
458 CHKERR VecNorm(T, NORM_2, &nrm2);
459 MOFEM_LOG("EXAMPLE", Sev::verbose) << "Regression norm " << nrm2;
460 double regression_value = 0;
461 switch (test_nb) {
462 case 1:
463 regression_value = 1.02789;
464 break;
465 case 2:
466 regression_value = 1.8841e+00;
467 break;
468 case 3:
469 regression_value = 1.8841e+00;
470 break;
471 case 4: // just links
472 regression_value = 4.9625e+00;
473 break;
474 case 5: // link master
475 regression_value = 6.6394e+00;
476 break;
477 case 6: // link master swap
478 regression_value = 4.98764e+00;
479 break;
480 case 7: // link Y
481 regression_value = 4.9473e+00;
482 break;
483 case 8: // link_3D_repr
484 regression_value = 2.5749e-01;
485 break;
486
487 default:
488 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Wrong test number.");
489 break;
490 }
491 if (fabs(nrm2 - regression_value) > 1e-2)
492 SETERRQ2(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
493 "Regression test field; wrong norm value. %6.4e != %6.4e", nrm2,
494 regression_value);
495 }
497}
498//! [Postprocessing results]
499
500//! [Check]
504}
505//! [Check]
506
507static char help[] = "...\n\n";
508
509int main(int argc, char *argv[]) {
510
511 // Initialisation of MoFEM/PETSc and MOAB data structures
512 const char param_file[] = "param_file.petsc";
513 MoFEM::Core::Initialize(&argc, &argv, param_file, help);
514
515 // Add logging channel for example
516 auto core_log = logging::core::get();
517 core_log->add_sink(
519 LogManager::setLog("EXAMPLE");
520 MOFEM_LOG_TAG("EXAMPLE", "example");
521
522 try {
523
524 //! [Register MoFEM discrete manager in PETSc]
525 DMType dm_name = "DMMOFEM";
526 CHKERR DMRegister_MoFEM(dm_name);
527 //! [Register MoFEM discrete manager in PETSc
528
529 //! [Create MoAB]
530 moab::Core mb_instance; ///< mesh database
531 moab::Interface &moab = mb_instance; ///< mesh database interface
532 //! [Create MoAB]
533
534 //! [Create MoFEM]
535 MoFEM::Core core(moab); ///< finite element database
536 MoFEM::Interface &m_field = core; ///< finite element database interface
537 //! [Create MoFEM]
538
539 //! [Example]
540 Example ex(m_field);
541 CHKERR ex.runProblem();
542 //! [Example]
543 }
545
547}
#define MOFEM_LOG_C(channel, severity, format,...)
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
int main()
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
#define CATCH_ERRORS
Catch errors.
FieldApproximationBase
approximation base
Definition definitions.h:58
@ LASTBASE
Definition definitions.h:69
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ DEMKOWICZ_JACOBI_BASE
Definition definitions.h:66
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ H1
continuous field
Definition definitions.h:85
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
constexpr int order
PostProcEleByDim< SPACE_DIM >::PostProcEleBdy PostProcEleBdy
auto integration_rule
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode, RowColData rc=RowColData::COL)
set local (or ghosted) vector values on mesh for partition only
Definition DMMoFEM.cpp:514
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition DMMoFEM.cpp:576
auto createDMVector(DM dm, RowColData rc=RowColData::COL)
Get smart vector from DM.
Definition DMMoFEM.hpp:1237
@ PETSC
Standard PETSc assembly.
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.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
double D
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
Set Monitor To TS solver.
Definition DMMoFEM.cpp:1046
auto getDMTsCtx(DM dm)
Get TS context data structure used by DM.
Definition DMMoFEM.hpp:1279
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
auto createVectorMPI(MPI_Comm comm, PetscInt n, PetscInt N)
Create MPI Vector.
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)
int save_every_nth_step
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
static constexpr int approx_order
[Example]
Definition plastic.cpp:216
MoFEMErrorCode boundaryCondition()
MoFEMErrorCode assembleSystem()
MoFEMErrorCode readMesh()
FieldApproximationBase base
Choice of finite element basis functions
Definition plot_base.cpp:68
MoFEMErrorCode checkResults()
MoFEMErrorCode solveSystem()
Example(MoFEM::Interface &m_field)
MoFEMErrorCode runProblem()
MoFEMErrorCode gettingNorms()
[Solve]
MoFEM::Interface & mField
Reference to MoFEM interface.
Definition plastic.cpp:226
MoFEMErrorCode setupProblem()
MoFEMErrorCode outputResults()
Add operators pushing bases from local to physical configuration.
Boundary condition manager for finite element problem setup.
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() 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:118
Deprecated interface functions.
Definition of the displacement bc data structure.
Definition BCData.hpp:72
Data on single entity (This is passed as argument to DataOperator::doWork)
Class (Function) to enforce essential constrains on the left hand side diagonal.
Definition Essential.hpp:33
Class (Function) to enforce essential constrains on the right hand side diagonal.
Definition Essential.hpp:41
Class (Function) to calculate residual side diagonal.
Definition Essential.hpp:49
Class (Function) to enforce essential constrains.
Definition Essential.hpp:25
Structure for user loop methods on finite elements.
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.
Assembly methods.
Definition Natural.hpp:65
Get norm of input MatrixDouble for Tensor1.
Get norm of input MatrixDouble for Tensor2.
Specialization for double precision vector field values calculation.
Element used to execute operators on side of the element.
Post post-proc data at points from hash maps.
Template struct for dimension-specific finite element types.
PipelineManager interface.
Simple interface for fast problem set-up.
Definition Simple.hpp:27
MoFEMErrorCode getOptions()
get options
Definition Simple.cpp:180
const std::string getDomainFEName() const
Get the Domain FE Name.
Definition Simple.hpp:429
intrusive_ptr for managing petsc objects
PetscInt ts_step
Current time step number.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
[Push operators to pipeline]
boost::shared_ptr< PostProcEleBdy > postProc
SmartPetscObj< DM > dM
MoFEMErrorCode postProcess()
Post-processing function executed at loop completion.
boost::shared_ptr< PostProcEle > postProc
Monitor(SmartPetscObj< DM > dm, boost::shared_ptr< PostProcEleBdy > post_proc)
PipelineManager::ElementsAndOpsByDim< 2 >::FaceSideEle SideEle
PipelineManager::ElementsAndOpsByDim< 3 >::FaceSideEle SideEle
#define EXECUTABLE_DIMENSION
Definition plastic.cpp:13
ElementsAndOps< SPACE_DIM >::SideEle SideEle
Definition plastic.cpp:61
static char help[]
[Check]
constexpr int SPACE_DIM