v0.15.5
Loading...
Searching...
No Matches
BasicBoundaryConditionsInterface.hpp
Go to the documentation of this file.
1/** \file BasicBoundaryConditionsInterface.hpp
2 \brief Header file for BasicBoundaryConditionsInterface element implementation
3*/
4
5
6
7// TODO: This is still work in progress !!!
8
9// FIX and ROTATE Boundary conditions
10// include SurfacePressureComplexForLazy for following load
11
12
13#ifndef __BASICBOUNDARYCONDIONSINTERFACE_HPP__
14#define __BASICBOUNDARYCONDIONSINTERFACE_HPP__
15
16
17/** \brief Set of functions declaring elements and setting operators
18 * for basic boundary conditions interface
19 */
21
22 using DomainEle = VolumeElementForcesAndSourcesCore;
23 using DomainEleOp = DomainEle::UserDataOperator;
25 using BoundaryEleOp = BoundaryEle::UserDataOperator;
26
27 using OpBodyForce = FormsIntegrators<DomainEleOp>::Assembly<
30 NaturalBC<DomainEleOp>::Assembly<PETSC>::LinearForm<GAUSS>;
33
34 using OpMass = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::BiLinearForm<
35 GAUSS>::OpMass<1, 3>;
36 using OpInertiaForce = FormsIntegrators<DomainEleOp>::Assembly<
38
42
44 FTensor::Index<'i', 3> i;
45 tForce(i) *= sCale;
46 return tForce;
47 }
48
49 private:
50 double sCale;
52 };
53
55 double sCale;
56 BasicBCVectorScale(double scale, std::string file_name)
57 : sCale(scale), TimeScaleVector3(file_name, false) {}
58
61 auto vec2 = MoFEM::TimeScaleVector3::getVector(time);
62 FTensor::Index<'i', 3> i;
63 vec(i) = sCale * vec2(i);
64 return vec;
65 }
66 };
67
69 SmartPetscObj<DM> dM;
70
72
74 double *lAmbda;
75 LoadScale(double *my_lambda) : lAmbda(my_lambda){};
76 MoFEMErrorCode scaleNf(const FEMethod *fe, VectorDouble &nf) {
78 nf *= *lAmbda;
80 }
81 };
82
87
88 BitRefLevel bIt;
89
92
93 boost::ptr_map<std::string, NeumannForcesSurface> neumann_forces;
94 boost::ptr_map<std::string, NodalForce> nodal_forces;
95 boost::ptr_map<std::string, EdgeForce> edge_forces;
96
97 boost::shared_ptr<FluidPressure> fluidPressureElementPtr;
98
99 boost::shared_ptr<FaceElementForcesAndSourcesCore> springRhsPtr;
100 boost::shared_ptr<FaceElementForcesAndSourcesCore> springLhsPtr;
101
102 boost::shared_ptr<VolumeElementForcesAndSourcesCore> bodyForceRhsPtr;
103 boost::shared_ptr<VolumeElementForcesAndSourcesCore> bodyForceLhsPtr;
104 // boost::shared_ptr<FEMethod> dirichletBcPtr;
105
106 boost::shared_ptr<DirichletDisplacementBc> dirichletBcPtr;
107 boost::shared_ptr<KelvinVoigtDamper> damperElementPtr;
108
109 const string domainProblemName;
110 const string domainElementName;
111
112 // boost::shared_ptr<PostProcVolumeOnRefinedMesh> postProc;
113 // boost::shared_ptr<PostProcFaceOnRefinedMesh> postProcSkin;
114
116 MoFEM::Interface &m_field, string postion_field,
117 string mesh_pos_field_name = "MESH_NODE_POSITIONS",
118 string problem_name = "ELASTIC",
119 string domain_element_name = "ELASTIC_FE",
120 bool is_displacement_field = true, bool is_quasi_static = true,
121 double *snes_load_factor = nullptr, bool is_partitioned = true)
122 : mField(m_field), positionField(postion_field),
123 meshNodeField(mesh_pos_field_name), domainProblemName(problem_name),
124 domainElementName(domain_element_name),
125 isDisplacementField(is_displacement_field),
127 snesLambdaLoadFactorPtr(snes_load_factor),
128 isPartitioned(is_partitioned), isLinear(PETSC_FALSE) {}
129
131
132 MoFEMErrorCode getCommandLineParameters() override {
134
135 PetscBool quasi_static = PETSC_FALSE;
136 PetscBool is_linear = PETSC_FALSE;
137 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "-is_quasi_static", &quasi_static,
138 PETSC_NULLPTR);
139 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "-is_linear", &is_linear,
140 PETSC_NULLPTR);
141
142 isQuasiStatic = quasi_static;
143 isLinear = is_linear;
144
146 };
147
148 MoFEMErrorCode addElementFields() override { return 0;};
149
150 MoFEMErrorCode createElements() override {
152
158
160 dirichletBcPtr = boost::make_shared<DirichletSpatialRemoveDofsBc>(
162 "DISPLACEMENT", isPartitioned);
163 else
164 dirichletBcPtr = boost::make_shared<DirichletDisplacementRemoveDofsBc>(
165 mField, positionField, domainProblemName, "DISPLACEMENT",
167 // CHKERR dynamic_cast<DirichletDisplacementRemoveDofsBc &>(
168 // *dirichletBcPtr).iNitialize();
169
179
180 // CHKERR mField.add_finite_element(domainElementName, "FLUID_PRESSURE_FE");
181
182 fluidPressureElementPtr = boost::make_shared<FluidPressure>(mField);
183 fluidPressureElementPtr->addNeumannFluidPressureBCElements(positionField);
184 CHKERR addHOOpsFace3D(meshNodeField, fluidPressureElementPtr->getLoopFe(),
185 false, false);
186
187 damperElementPtr = boost::make_shared<KelvinVoigtDamper>(mField);
188 damperElementPtr->commonData.meshNodePositionName = meshNodeField;
189
190 auto &common_data = damperElementPtr->commonData;
191
192 common_data.spatialPositionName = positionField;
193 common_data.spatialPositionNameDot = "DOT_" + positionField;
194 damperElementPtr->setBlockDataMap(); FIXME:
195
196 for (auto &[id, data] : damperElementPtr->blockMaterialDataMap) {
197 data.lInear = isLinear;
198 int cid = id;
199 damperElementPtr->constitutiveEquationMap.insert(
202 "DAMPER_FE");
203 }
204
206 };
207
208 MoFEMErrorCode setOperators() override {
214
216 mField, neumann_forces, PETSC_NULLPTR, positionField);
217 springLhsPtr = boost::make_shared<FaceElementForcesAndSourcesCore>(mField);
218 springRhsPtr = boost::make_shared<FaceElementForcesAndSourcesCore>(mField);
220 boost::make_shared<VolumeElementForcesAndSourcesCore>(mField);
222 boost::make_shared<VolumeElementForcesAndSourcesCore>(mField);
223
227
228
229 fluidPressureElementPtr->setNeumannFluidPressureFiniteElementOperators(
230 positionField, PETSC_NULLPTR, true, true);
231
232 // KelvinVoigtDamper::CommonData &common_data =
233 // damperElementPtr->commonData;
234 CHKERR damperElementPtr->setOperators(3);
235
236 // auto dm = mField.getInterface<Simple>()->getDM();
237
238 auto get_id_block_param = [&](string base_name, int id) {
239 char load_hist_file[255] = "hist.in";
240 PetscBool ctg_flag = PETSC_FALSE;
241 string param_name_with_id = "-" + base_name + "_" + to_string(id);
242 CHKERR PetscOptionsGetString(PETSC_NULLPTR, PETSC_NULLPTR,
243 param_name_with_id.c_str(), load_hist_file,
244 255, &ctg_flag);
245 if (ctg_flag)
246 return param_name_with_id;
247
248 param_name_with_id = "-" + base_name;
249 CHKERR PetscOptionsGetString(PETSC_NULLPTR, PETSC_NULLPTR,
250 param_name_with_id.c_str(), load_hist_file,
251 255, &ctg_flag);
252 if (ctg_flag) {
253 MOFEM_LOG("WORLD", Sev::verbose)
254 << "Setting one accelerogram for all blocks!";
255 return param_name_with_id;
256 }
257
258 return string("");
259 };
260
261 auto get_adj_ents = [&](const Range &ents) {
262 Range verts;
263 CHKERR mField.get_moab().get_connectivity(ents, verts, true);
264 for (size_t d = 1; d < 3; ++d)
265 CHKERR mField.get_moab().get_adjacencies(ents, d, false, verts,
266 moab::Interface::UNION);
267 verts.merge(ents);
268 CHKERR mField.getInterface<CommInterface>()->synchroniseEntities(verts);
269 return verts;
270 };
271
273 const std::string block_name = "BODY_FORCE";
274 if (it->getName().compare(0, block_name.size(), block_name) == 0) {
275 std::vector<double> attr;
276 CHKERR it->getAttributes(attr);
277 if (attr.size() > 3) {
278
279 const int id = it->getMeshsetId();
280
281 Range bc_ents;
282 CHKERR it->getMeshsetIdEntitiesByDimension(mField.get_moab(), 3,
283 bc_ents, true);
284 auto bc_ents_ptr = boost::make_shared<Range>(get_adj_ents(bc_ents));
285
286 // the first parameter is density!!! FIXME:
287
288 VectorDouble acc({attr[1], attr[2], attr[3]});
289 double density = attr[0];
290 bool inertia_flag =
291 attr.size() > 4 ? bool(std::floor(attr[4])) : true;
292 // if accelerogram is provided then change the acceleration,
293 // otherwise use whatever is in that block TODO:
294 std::vector<boost::shared_ptr<TimeScaleVector3>> methods_for_scaling;
295 std::string param_name_for_scaling =
296 get_id_block_param("accelerogram", id);
297
298 if (!param_name_for_scaling.empty())
299 methods_for_scaling.push_back(
300 boost::make_shared<BasicBCVectorScale>(density,
301 param_name_for_scaling));
302 else
303 methods_for_scaling.push_back(
304 boost::make_shared<BasicBCVectorConst>(
305 density,
306 FTensor::Tensor1<double, 3>({acc(0), acc(1), acc(2)})));
307
308 // FIXME: this require correction for large strains, multiply by det_F
309 auto get_rho = [&](double, double, double) {
310 auto &fe_domain_lhs = bodyForceLhsPtr;
311 return density * fe_domain_lhs->ts_aa;
312 };
313
314 auto &pipeline_rhs = bodyForceRhsPtr->getOpPtrVector();
315 auto &pipeline_lhs = bodyForceLhsPtr->getOpPtrVector();
316
317 CHKERR addHOOpsVol(meshNodeField, *bodyForceRhsPtr, true, false,
318 false, false);
319 CHKERR addHOOpsVol(meshNodeField, *bodyForceLhsPtr, true, false,
320 false, false);
321
322 //FIXME: fix for large strains
323 pipeline_rhs.push_back(
324 new OpSetBc(positionField, true, mBoundaryMarker));
325
326 // using new way of adding BCs
327 CHKERR
329 pipeline_rhs, mField, "U", {}, methods_for_scaling,
330 it->getName(), Sev::inform);
331
332 if (!isQuasiStatic && inertia_flag) {
333 pipeline_lhs.push_back(
334 new OpSetBc(positionField, true, mBoundaryMarker));
335 pipeline_lhs.push_back(
336 new OpMass(positionField, positionField, get_rho, bc_ents_ptr));
337 auto mat_acceleration = boost::make_shared<MatrixDouble>();
338 pipeline_rhs.push_back(new OpCalculateVectorFieldValuesDotDot<3>(
339 positionField, mat_acceleration));
340 pipeline_rhs.push_back(new OpInertiaForce(
341 positionField, mat_acceleration,
342 [&](double, double, double) { return density; }, bc_ents_ptr));
343 pipeline_lhs.push_back(new OpUnSetBc(positionField));
344 }
345 pipeline_rhs.push_back(new OpUnSetBc(positionField));
346
347 } else {
348 SETERRQ(
349 PETSC_COMM_SELF, MOFEM_INVALID_DATA,
350 "There should be (1 density + 3 accelerations ) attributes in "
351 "BODY_FORCE blockset, but is %ld. Optionally, you can set 5th "
352 "parameter to inertia flag.",
353 attr.size());
354 }
355 }
356 }
357
358 auto integration_rule_vol = [](int, int, int approx_order) {
359 return 2 * approx_order + 1;
360 };
361 auto integration_rule_boundary = [](int, int, int approx_order) {
362 return 2 * approx_order + 1;
363 };
364
365 springLhsPtr->getRuleHook = integration_rule_boundary;
366 springRhsPtr->getRuleHook = integration_rule_boundary;
367 bodyForceLhsPtr->getRuleHook = integration_rule_vol;
368 bodyForceRhsPtr->getRuleHook = integration_rule_vol;
369
370 // set other boundary conditions
371 auto bc_mng = mField.getInterface<BcManager>();
372 auto *pipeline_mng = mField.getInterface<PipelineManager>();
373
374 CHKERR bc_mng->removeBlockDOFsOnEntities(domainProblemName, "REMOVE_X",
375 positionField, 0, 0, true,
377 CHKERR bc_mng->removeBlockDOFsOnEntities(domainProblemName, "REMOVE_Y",
378 positionField, 1, 1, true,
380 CHKERR bc_mng->removeBlockDOFsOnEntities(domainProblemName, "REMOVE_Z",
381 positionField, 2, 2, true,
383 CHKERR bc_mng->removeBlockDOFsOnEntities(domainProblemName, "REMOVE_ALL",
384 positionField, 0, 2, true,
386
387 CHKERR bc_mng->pushMarkDOFsOnEntities(domainProblemName, "FIX_X",
388 positionField, 0, 0);
389 CHKERR bc_mng->pushMarkDOFsOnEntities(domainProblemName, "FIX_Y",
390 positionField, 1, 1);
391 CHKERR bc_mng->pushMarkDOFsOnEntities(domainProblemName, "FIX_Z",
392 positionField, 2, 2);
393 CHKERR bc_mng->pushMarkDOFsOnEntities(domainProblemName, "FIX_ALL",
394 positionField, 0, 2);
395 CHKERR bc_mng->pushMarkDOFsOnEntities(domainProblemName, "ROTATE",
396 positionField, 0, 2);
397
399 bc_mng->getMergedBlocksMarker(vector<string>{"FIX_", "ROTATE"});
400
402 };
403
404 MoFEMErrorCode addElementsToDM(SmartPetscObj<DM> dm) override {
406 this->dM = dm;
407 auto simple = mField.getInterface<Simple>();
408 vector<const char *> element_list{"FORCE_FE", "PRESSURE_FE",
409 "FLUID_PRESSURE_FE",
410 "SPRING",
411 "DAMPER_FE"};
412 for (auto &el : element_list) {
413 CHKERR DMMoFEMAddElement(dM, el);
414 simple->getOtherFiniteElements().push_back(el);
415 }
416 // if (!fluidPressureElementPtr->setOfFluids.empty())
417 // FIXME:
418 // CHKERR mField.modify_problem_add_finite_element(domainProblemName,
419 // "FLUID_PRESSURE_FE");
420 // CHKERR dynamic_cast<DirichletDisplacementRemoveDofsBc &>(
421 // *dirichletBcPtr).iNitialize();
423 };
424
425 MoFEMErrorCode updateElementVariables() override { return 0; };
426 MoFEMErrorCode postProcessElement(int step) override { return 0; };
427
428 string getHistoryParam(string prefix) {
429 char load_hist_file[255] = "hist.in";
430 PetscBool ctg_flag = PETSC_FALSE;
431 string new_param_file = string("-") + prefix + string("_history");
432 CHKERR PetscOptionsGetString(PETSC_NULLPTR, PETSC_NULLPTR, new_param_file.c_str(),
433 load_hist_file, 255, &ctg_flag);
434 if (ctg_flag)
435 return new_param_file;
436 return string("-load_history");
437 };
438
439 template <typename T>
440 MoFEMErrorCode setupSolverFunction(const TSType type = IM) {
441 CHKERR setupSolverImpl<T, true>(type);
443 }
444
445 template <typename T>
446 MoFEMErrorCode setupSolverJacobian(const TSType type = IM) {
447 CHKERR setupSolverImpl<T, false>(type);
449 }
450
451 template <typename T, bool RHS>
452 MoFEMErrorCode setupSolverImpl(const TSType type = IM) {
454 // auto dm = dM;
455 // boost::shared_ptr<FEMethod> null;
456
457 auto set_solver_pipelines =
458 [&](PetscErrorCode (*function)(DM, const char fe_name[], MoFEM::FEMethod *,
461 PetscErrorCode (*jacobian)(DM, const char fe_name[], MoFEM::FEMethod *,
464
466 if (RHS) {
467 if (std::is_same_v<T, TS>)
468 dirichletBcPtr->methodsOp.push_back(
469 new TimeForceScale(getHistoryParam("dirichlet"), false));
470
471 CHKERR DMoFEMPreProcessFiniteElements(dM, dirichletBcPtr.get());
472
473 // auto push_fmethods = [&](auto method, string element_name) {
474 // CHKERR function(dm, element_name.c_str(), method, method, method);
475 // };
476
477 auto set_neumann_methods = [&](auto &neumann_el, string hist_name,
478 int dim) {
480 for (auto &&mit : neumann_el) {
481 if constexpr (std::is_same_v<T, SNES>)
482 mit->second->methodsOp.push_back(
484 if constexpr (std::is_same_v<T, TS>)
485 mit->second->methodsOp.push_back(
486 new TimeForceScale(getHistoryParam(hist_name), false));
487 string element_name = mit->first;
488 switch (dim) {
489 case 2:
491 mit->second->getLoopFe().getOpPtrVector(), {},
493 break;
494 case 1:
496 mit->second->getLoopFe().getOpPtrVector(), {},
498 break;
499 case 0:
500 break;
501 default:
502 break;
503 }
504 // CHKERR push_fmethods(&mit->second->getLoopFe(),
505 // element_name);
506 CHKERR function(dM, element_name.c_str(),
507 &mit->second->getLoopFe(), NULL, NULL);
508 }
510 };
511
512 CHKERR set_neumann_methods(neumann_forces, "force", 2);
513 CHKERR set_neumann_methods(nodal_forces, "force", 0);
514 CHKERR set_neumann_methods(edge_forces, "force", 1);
515
516 CHKERR function(dM, domainElementName.c_str(), dirichletBcPtr.get(),
517 dirichletBcPtr.get(), dirichletBcPtr.get());
518 CHKERR function(dM, domainElementName.c_str(),
519 bodyForceRhsPtr.get(), NULL, NULL);
520 CHKERR function(dM, "SPRING", springRhsPtr.get(), NULL, NULL);
521 CHKERR function(dM, "DAMPER_FE", &damperElementPtr->feRhs, NULL,
522 NULL);
523 CHKERR function(dM, "FLUID_PRESSURE_FE",
524 &fluidPressureElementPtr->getLoopFe(), NULL, NULL);
525 } else {
526
527 CHKERR jacobian(dM, domainElementName.c_str(), dirichletBcPtr.get(),
528 dirichletBcPtr.get(), dirichletBcPtr.get());
529 CHKERR jacobian(dM, domainElementName.c_str(),
530 bodyForceLhsPtr.get(), NULL, NULL);
531 CHKERR jacobian(dM, "SPRING", springLhsPtr.get(), NULL, NULL);
532
533 CHKERR jacobian(dM, "DAMPER_FE", &damperElementPtr->feLhs, NULL,
534 NULL);
535 }
536
538 };
539
540 if constexpr (std::is_same_v<T, SNES>) {
541
543 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
544 "SNES lambda factor pointer not set in the module constructor");
545
546 CHKERR set_solver_pipelines(&DMMoFEMSNESSetFunction,
547 &DMMoFEMSNESSetJacobian);
548
549 } else if constexpr (std::is_same_v<T, TS>) {
550
551 switch (type) {
552 case IM:
553 CHKERR set_solver_pipelines(&DMMoFEMTSSetIFunction,
554 &DMMoFEMTSSetIJacobian);
555 break;
556 case IM2:
557 CHKERR set_solver_pipelines(&DMMoFEMTSSetI2Function,
558 &DMMoFEMTSSetI2Jacobian);
559 break;
560 case EX:
561 CHKERR set_solver_pipelines(&DMMoFEMTSSetRHSFunction,
562 &DMMoFEMTSSetRHSJacobian);
563 break;
564 default:
565 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
566 "This TS is not yet implemented for basic BCs");
567 break;
568 }
569
570 } else
571 static_assert(!std::is_same_v<T, KSP>,
572 "this solver has not been implemented for basic BCs yet");
573
574
576 }
577
578 MoFEMErrorCode setupSolverFunctionSNES() override {
580 CHKERR this->setupSolverFunction<SNES>();
582 }
583 MoFEMErrorCode setupSolverJacobianSNES() override {
585 CHKERR this->setupSolverJacobian<SNES>();
587 }
588 MoFEMErrorCode
589 setupSolverFunctionTS(const TSType type) override {
591 CHKERR this->setupSolverFunction<TS>(type);
593 }
594 MoFEMErrorCode
595 setupSolverJacobianTS(const TSType type) override {
597 CHKERR this->setupSolverJacobian<TS>(type);
599 }
600
601};
602
603#endif //__BASICBOUNDARYCONDIONSINTERFACE_HPP__
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
@ MF_ZERO
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ BLOCKSET
@ MOFEM_INVALID_DATA
Definition definitions.h:36
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
#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 ...
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string name, const bool recursive=true)=0
add entities to finite element
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
#define MOFEM_LOG(channel, severity)
Log.
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
FTensor::Index< 'i', SPACE_DIM > i
static constexpr int approx_order
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMass
[Only used with Hooke equation (linear material model)]
Definition seepage.cpp:56
FTensor::Tensor1< double, 3 > getVector(const double time)
BasicBCVectorConst(double scale, FTensor::Tensor1< double, 3 > t_vec)
FTensor::Tensor1< double, 3 > getVector(const double time)
MoFEMErrorCode scaleNf(const FEMethod *fe, VectorDouble &nf)
Set of functions declaring elements and setting operators for basic boundary conditions interface.
BasicBoundaryConditionsInterface(MoFEM::Interface &m_field, string postion_field, string mesh_pos_field_name="MESH_NODE_POSITIONS", string problem_name="ELASTIC", string domain_element_name="ELASTIC_FE", bool is_displacement_field=true, bool is_quasi_static=true, double *snes_load_factor=nullptr, bool is_partitioned=true)
boost::shared_ptr< FaceElementForcesAndSourcesCore > springLhsPtr
MoFEMErrorCode setupSolverFunction(const TSType type=IM)
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, 3 > OpBodyForce
MoFEMErrorCode setupSolverJacobianTS(const TSType type) override
boost::shared_ptr< VolumeElementForcesAndSourcesCore > bodyForceLhsPtr
boost::shared_ptr< FluidPressure > fluidPressureElementPtr
boost::shared_ptr< FaceElementForcesAndSourcesCore > springRhsPtr
boost::ptr_map< std::string, EdgeForce > edge_forces
MoFEMErrorCode setupSolverImpl(const TSType type=IM)
boost::ptr_map< std::string, NodalForce > nodal_forces
boost::ptr_map< std::string, NeumannForcesSurface > neumann_forces
MoFEMErrorCode setupSolverFunctionTS(const TSType type) override
NaturalBC< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS > DomainNaturalBC
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, 3 > OpMass
boost::shared_ptr< VolumeElementForcesAndSourcesCore > bodyForceRhsPtr
MoFEMErrorCode addElementsToDM(SmartPetscObj< DM > dm) override
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesVector< 1, 3, 1 > OpInertiaForce
MoFEMErrorCode setupSolverJacobian(const TSType type=IM)
boost::shared_ptr< KelvinVoigtDamper > damperElementPtr
boost::shared_ptr< DirichletDisplacementBc > dirichletBcPtr
MoFEMErrorCode postProcessElement(int step) override
Set of functions declaring elements and setting operators for generic element interface.
static MoFEMErrorCode setOperators(MoFEM::Interface &m_field, boost::ptr_map< std::string, EdgeForce > &edge_forces, Vec F, const std::string field_name, std::string mesh_node_positions="MESH_NODE_POSITIONS")
Set integration point operators.
Definition EdgeForce.hpp:97
static MoFEMErrorCode addElement(MoFEM::Interface &m_field, const std::string field_name, Range *intersect_ptr=NULL)
Add element taking information from NODESET.
Definition EdgeForce.hpp:62
static MoFEMErrorCode addNeumannBCElements(MoFEM::Interface &m_field, const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS", Range *intersect_ptr=NULL)
Declare finite element.
static MoFEMErrorCode setMomentumFluxOperators(MoFEM::Interface &m_field, boost::ptr_map< std::string, NeumannForcesSurface > &neumann_forces, Vec F, const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
Set operators to finite elements calculating right hand side vector.
static MoFEMErrorCode setOperators(MoFEM::Interface &m_field, boost::ptr_map< std::string, NodalForce > &nodal_forces, Vec F, const std::string field_name)
Set integration point operators.
static MoFEMErrorCode addElement(MoFEM::Interface &m_field, const std::string field_name, Range *intersect_ptr=NULL)
Add element taking information from NODESET.
static MoFEMErrorCode setSpringOperators(MoFEM::Interface &m_field, boost::shared_ptr< FaceElementForcesAndSourcesCore > fe_spring_lhs_ptr, boost::shared_ptr< FaceElementForcesAndSourcesCore > fe_spring_rhs_ptr, const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS", double stiffness_scale=1.)
Implementation of spring element. Set operators to calculate LHS and RHS.
static MoFEMErrorCode addSpringElements(MoFEM::Interface &m_field, const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
Declare spring element.
Class used to scale loads, f.e. in arc-length control.
Data structure to exchange data between MoFEM and User Loop Methods.
virtual moab::Interface & get_moab()=0
Deprecated interface functions.
Structure for user loop methods on finite elements.
Force scale operator for reading four columns (time and vector)
virtual FTensor::Tensor1< double, SPACE_DIM > getVector(const double time)
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Force scale operator for reading two columns.
PetscBool is_quasi_static
Definition plastic.cpp:143
double scale
Definition plastic.cpp:123