v0.15.4
Loading...
Searching...
No Matches
Classes | Public Member Functions | Private Member Functions | Private Attributes | List of all members
Seepage Struct Reference
Collaboration diagram for Seepage:
[legend]

Classes

struct  BlockedParameters
 

Public Member Functions

 Seepage (MoFEM::Interface &m_field)
 
MoFEMErrorCode runProblem ()
 [Run problem]
 

Private Member Functions

MoFEMErrorCode setupProblem ()
 [Run problem]
 
MoFEMErrorCode createCommonData ()
 [Set up problem]
 
MoFEMErrorCode bC ()
 [Create common data]
 
MoFEMErrorCode OPs ()
 [Boundary condition]
 
MoFEMErrorCode tsSolve ()
 [Solve]
 
MoFEMErrorCode addMatBlockOps (boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, std::string block_elastic_name, std::string block_thermal_name, boost::shared_ptr< BlockedParameters > blockedParamsPtr, Sev sev)
 

Private Attributes

MoFEM::InterfacemField
 
PetscBool doEvalField
 
std::array< double, 3 > fieldEvalCoords {0.0, 0.0, 0.0}
 
boost::shared_ptr< FieldEvaluatorInterface::SetPtsDatafieldEvalData
 
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uXScatter
 
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uYScatter
 
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uZScatter
 
boost::shared_ptr< VectorDoublepressureFieldPtr
 
boost::shared_ptr< MatrixDoublefluxFieldPtr
 
boost::shared_ptr< MatrixDoubledispFieldPtr
 
boost::shared_ptr< MatrixDoubledispGradPtr
 
boost::shared_ptr< MatrixDoublestrainFieldPtr
 
boost::shared_ptr< MatrixDoublestressFieldPtr
 
boost::shared_ptr< MatrixDoublemDPtr
 

Detailed Description

Examples
mofem/tutorials/adv-5_poroelasticity/seepage.cpp.

Definition at line 194 of file seepage.cpp.

Constructor & Destructor Documentation

◆ Seepage()

Seepage::Seepage ( MoFEM::Interface m_field)
inline

Definition at line 196 of file seepage.cpp.

196: mField(m_field) {}
MoFEM::Interface & mField
Definition seepage.cpp:201

Member Function Documentation

◆ addMatBlockOps()

MoFEMErrorCode Seepage::addMatBlockOps ( boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &  pipeline,
std::string  block_elastic_name,
std::string  block_thermal_name,
boost::shared_ptr< BlockedParameters blockedParamsPtr,
Sev  sev 
)
private

[Calculate elasticity tensor]

[Calculate elasticity tensor]

Examples
mofem/tutorials/adv-5_poroelasticity/seepage.cpp.

Definition at line 265 of file seepage.cpp.

268 {
270
271 struct OpMatElasticBlocks : public DomainEleOp {
272 OpMatElasticBlocks(boost::shared_ptr<MatrixDouble> m, double bulk_modulus_K,
273 double shear_modulus_G, MoFEM::Interface &m_field,
274 Sev sev,
275 std::vector<const CubitMeshSets *> meshset_vec_ptr)
276 : DomainEleOp(NOSPACE, DomainEleOp::OPSPACE), matDPtr(m),
277 bulkModulusKDefault(bulk_modulus_K),
278 shearModulusGDefault(shear_modulus_G) {
279 CHK_THROW_MESSAGE(extractElasticBlockData(m_field, meshset_vec_ptr, sev),
280 "Can not get data from block");
281 }
282
283 MoFEMErrorCode doWork(int side, EntityType type,
286
287 for (auto &b : blockData) {
288
289 if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
290 CHKERR getMatDPtr(matDPtr, b.bulkModulusK, b.shearModulusG);
292 }
293 }
294
295 CHKERR getMatDPtr(matDPtr, bulkModulusKDefault, shearModulusGDefault);
297 }
298
299 private:
300 boost::shared_ptr<MatrixDouble> matDPtr;
301
302 struct BlockData {
303 double bulkModulusK;
304 double shearModulusG;
305 Range blockEnts;
306 };
307
308 double bulkModulusKDefault;
309 double shearModulusGDefault;
310 std::vector<BlockData> blockData;
311
313 extractElasticBlockData(MoFEM::Interface &m_field,
314 std::vector<const CubitMeshSets *> meshset_vec_ptr,
315 Sev sev) {
317
318 for (auto m : meshset_vec_ptr) {
319 MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Elastic Block") << *m;
320 std::vector<double> block_data;
321 CHKERR m->getAttributes(block_data);
322 if (block_data.size() < 2) {
323 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
324 "Expected that block has two attributes");
325 }
326 auto get_block_ents = [&]() {
327 Range ents;
328 CHKERR
329 m_field.get_moab().get_entities_by_handle(m->meshset, ents, true);
330 return ents;
331 };
332
333 double young_modulus = block_data[0];
334 double poisson_ratio = block_data[1];
335 double bulk_modulus_K = young_modulus / (3 * (1 - 2 * poisson_ratio));
336 double shear_modulus_G = young_modulus / (2 * (1 + poisson_ratio));
337
338 MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Elastic Block")
339 << m->getName() << ": E = " << young_modulus
340 << " nu = " << poisson_ratio;
341
342 blockData.push_back(
343 {bulk_modulus_K, shear_modulus_G, get_block_ents()});
344 }
345 MOFEM_LOG_CHANNEL("WORLD");
347 }
348
349 MoFEMErrorCode getMatDPtr(boost::shared_ptr<MatrixDouble> mat_D_ptr,
350 double bulk_modulus_K, double shear_modulus_G) {
352 //! [Calculate elasticity tensor]
353 auto set_material_stiffness = [&]() {
359 double A = 1.;
360 if (SPACE_DIM == 2 && !is_plane_strain) {
361 A = 2 * shear_modulus_G /
362 (bulk_modulus_K + (4. / 3.) * shear_modulus_G);
363 }
364 auto t_D = getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*mat_D_ptr);
365 t_D(i, j, k, l) =
366 2 * shear_modulus_G * ((t_kd(i, k) ^ t_kd(j, l)) / 4.) +
367 A * (bulk_modulus_K - (2. / 3.) * shear_modulus_G) * t_kd(i, j) *
368 t_kd(k, l);
369 };
370 //! [Calculate elasticity tensor]
371 constexpr auto size_symm = (SPACE_DIM * (SPACE_DIM + 1)) / 2;
372 mat_D_ptr->resize(size_symm * size_symm, 1);
373 set_material_stiffness();
375 }
376 };
377
378 double default_bulk_modulus_K =
380 double default_shear_modulus_G =
382
383 pipeline.push_back(new OpMatElasticBlocks(
384 blockedParamsPtr->getDPtr(), default_bulk_modulus_K,
385 default_shear_modulus_G, mField, sev,
386
387 // Get blockset using regular expression
388 mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
389
390 (boost::format("%s(.*)") % block_elastic_name).str()
391
392 ))
393
394 ));
395
396 struct OpMatFluidBlocks : public DomainEleOp {
397 OpMatFluidBlocks(boost::shared_ptr<double> conductivity_ptr,
398 boost::shared_ptr<double> biot_constant_ptr,
399 boost::shared_ptr<double> fluid_density_ptr,
400 boost::shared_ptr<double> mat_density_ptr,
401 boost::shared_ptr<double> storage_constant_ptr,
402 MoFEM::Interface &m_field, Sev sev,
403 std::vector<const CubitMeshSets *> meshset_vec_ptr)
404 : DomainEleOp(NOSPACE, DomainEleOp::OPSPACE),
405 conductivityPtr(conductivity_ptr), biotConstantPtr(biot_constant_ptr),
406 fluidDensityPtr(fluid_density_ptr), matDensityPtr(mat_density_ptr),
407 storageConstantPtr(storage_constant_ptr) {
408 CHK_THROW_MESSAGE(extractThermalBlockData(m_field, meshset_vec_ptr, sev),
409 "Can not get data from block");
410 }
411
412 MoFEMErrorCode doWork(int side, EntityType type,
415
416 for (auto &b : blockData) {
417
418 if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
419 *conductivityPtr = b.conductivity;
420 *biotConstantPtr = b.biot_constant;
421 *fluidDensityPtr = b.fluid_density;
422 *matDensityPtr = b.mat_density;
423 *storageConstantPtr = b.storage;
425 }
426 }
427
428 *conductivityPtr = default_conductivity;
429 *biotConstantPtr = default_biot;
430 *fluidDensityPtr = default_fluid_density;
431 *matDensityPtr = default_mat_density;
432 *storageConstantPtr = default_storage;
433
435 }
436
437 private:
438 struct BlockData {
439 double conductivity;
440 double biot_constant;
441 double storage;
442 double mat_density;
443 double fluid_density;
444 Range blockEnts;
445 };
446
447 std::vector<BlockData> blockData;
448
450 extractThermalBlockData(MoFEM::Interface &m_field,
451 std::vector<const CubitMeshSets *> meshset_vec_ptr,
452 Sev sev) {
454
455 for (auto m : meshset_vec_ptr) {
456 MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Fluid Block") << *m;
457 std::vector<double> block_data;
458 CHKERR m->getAttributes(block_data);
459 if (block_data.size() < 5) {
460 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
461 "Expected that block has five attributes");
462 }
463 auto get_block_ents = [&]() {
464 Range ents;
465 CHKERR
466 m_field.get_moab().get_entities_by_handle(m->meshset, ents, true);
467 return ents;
468 };
469
470 MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Fluid Block")
471 << m->getName() << ": conductivity = " << block_data[0]
472 << ", biot_constant = " << block_data[1]
473 << ", storage = " << block_data[2]
474 << ", mat_density = " << block_data[3]
475 << ", fluid_density = " << block_data[4];
476
477 blockData.push_back({block_data[0], block_data[1], block_data[2],
478 block_data[3], block_data[4], get_block_ents()});
479 }
480 MOFEM_LOG_CHANNEL("WORLD");
482 }
483
484 boost::shared_ptr<double> expansionPtr;
485 boost::shared_ptr<double> conductivityPtr;
486 boost::shared_ptr<double> biotConstantPtr;
487 boost::shared_ptr<double> fluidDensityPtr;
488 boost::shared_ptr<double> matDensityPtr;
489 boost::shared_ptr<double> storageConstantPtr;
490 boost::shared_ptr<double> capacityPtr;
491 };
492
493 pipeline.push_back(new OpMatFluidBlocks(
494 blockedParamsPtr->getConductivityPtr(),
495 blockedParamsPtr->getBiotConstantPtr(),
496 blockedParamsPtr->getFluidDensityPtr(),
497 blockedParamsPtr->getMatDensityPtr(),
498 blockedParamsPtr->getStorageConstantPtr(), mField, sev,
499
500 // Get blockset using regular expression
502
503 (boost::format("%s(.*)") % block_thermal_name).str()
504
505 ))
506
507 ));
508
510}
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
constexpr int SPACE_DIM
Kronecker Delta class symmetric.
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ NOSPACE
Definition definitions.h:83
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
double bulk_modulus_K
double shear_modulus_G
constexpr auto t_kd
SeverityLevel
Severity levels.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
MoFEMErrorCode getCubitMeshsetPtr(const int ms_id, const CubitBCType cubit_bc_type, const CubitMeshSets **cubit_meshset_ptr) const
get cubit meshset
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
constexpr AssemblyType A
double default_biot
Definition seepage.cpp:172
double default_conductivity
Definition seepage.cpp:171
PetscBool is_plane_strain
Definition seepage.cpp:177
double default_fluid_density
Definition seepage.cpp:174
double default_young_modulus
Definition seepage.cpp:169
double default_mat_density
Definition seepage.cpp:175
double default_storage
Definition seepage.cpp:173
double default_poisson_ratio
Definition seepage.cpp:170
FTensor::Index< 'm', 3 > m
virtual moab::Interface & get_moab()=0
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
Interface for managing meshsets containing materials and boundary conditions.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
double young_modulus
Young modulus.
Definition plastic.cpp:125
double poisson_ratio
Poisson ratio.
Definition plastic.cpp:126
constexpr auto size_symm
Definition plastic.cpp:42

◆ bC()

MoFEMErrorCode Seepage::bC ( )
private

[Create common data]

[Boundary condition]

Examples
mofem/tutorials/adv-5_poroelasticity/seepage.cpp.

Definition at line 655 of file seepage.cpp.

655 {
657
658 MOFEM_LOG("SYNC", Sev::noisy) << "bC";
660
662 auto bc_mng = mField.getInterface<BcManager>();
663
664 CHKERR bc_mng->pushMarkDOFsOnEntities<DisplacementCubitBcData>(
665 simple->getProblemName(), "U");
666 CHKERR bc_mng->pushMarkDOFsOnEntities<HeatFluxCubitBcData>(
667 simple->getProblemName(), "FLUX", false);
668
669 auto get_skin = [&]() {
670 Range body_ents;
671 CHKERR mField.get_moab().get_entities_by_dimension(0, SPACE_DIM, body_ents);
672 Skinner skin(&mField.get_moab());
673 Range skin_ents;
674 CHKERR skin.find_skin(0, body_ents, false, skin_ents);
675 return skin_ents;
676 };
677
678 auto filter_flux_blocks = [&](auto skin) {
679 auto remove_cubit_blocks = [&](auto c) {
681 for (auto m :
682
683 mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(c)
684
685 ) {
686 Range ents;
687 CHKERR mField.get_moab().get_entities_by_dimension(
688 m->getMeshset(), SPACE_DIM - 1, ents, true);
689 skin = subtract(skin, ents);
690 }
692 };
693
694 auto remove_named_blocks = [&](auto n) {
696 for (auto m : mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(
697 std::regex(
698
699 (boost::format("%s(.*)") % n).str()
700
701 ))
702
703 ) {
704 Range ents;
705 CHKERR mField.get_moab().get_entities_by_dimension(
706 m->getMeshset(), SPACE_DIM - 1, ents, true);
707 skin = subtract(skin, ents);
708 }
710 };
711
712 CHK_THROW_MESSAGE(remove_cubit_blocks(NODESET | TEMPERATURESET),
713 "remove_cubit_blocks");
714 CHK_THROW_MESSAGE(remove_cubit_blocks(SIDESET | HEATFLUXSET),
715 "remove_cubit_blocks");
716 CHK_THROW_MESSAGE(remove_named_blocks("PRESSURE"), "remove_named_blocks");
717 CHK_THROW_MESSAGE(remove_named_blocks("FLUIDFLUX"), "remove_named_blocks");
718
719 return skin;
720 };
721
722 auto filter_true_skin = [&](auto skin) {
723 Range boundary_ents;
724 ParallelComm *pcomm =
725 ParallelComm::get_pcomm(&mField.get_moab(), MYPCOMM_INDEX);
726 CHKERR pcomm->filter_pstatus(skin, PSTATUS_SHARED | PSTATUS_MULTISHARED,
727 PSTATUS_NOT, -1, &boundary_ents);
728 return boundary_ents;
729 };
730
731 auto remove_flux_ents = filter_true_skin(filter_flux_blocks(get_skin()));
732
733 CHKERR mField.getInterface<CommInterface>()->synchroniseEntities(
734 remove_flux_ents);
735
736 MOFEM_LOG("SYNC", Sev::noisy) << remove_flux_ents << endl;
738
739#ifndef NDEBUG
740
743 (boost::format("flux_remove_%d.vtk") % mField.get_comm_rank()).str(),
744 remove_flux_ents);
745
746#endif
747
748 CHKERR mField.getInterface<ProblemsManager>()->removeDofsOnEntities(
749 simple->getProblemName(), "FLUX", remove_flux_ents);
750
752}
static auto filter_true_skin(MoFEM::Interface &m_field, Range &&skin)
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
#define MYPCOMM_INDEX
default communicator number PCOMM
@ TEMPERATURESET
@ HEATFLUXSET
@ NODESET
@ SIDESET
#define MOFEM_LOG(channel, severity)
Log.
const double c
speed of light (cm/ns)
const double n
refractive index of diffusive medium
auto save_range
Definition seepage.cpp:185
constexpr int SPACE_DIM
Definition seepage.cpp:33
Boundary condition manager for finite element problem setup.
Managing BitRefLevels.
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
Definition of the displacement bc data structure.
Definition BCData.hpp:72
Definition of the heat flux bc data structure.
Definition BCData.hpp:423
Problem manager is used to build and partition problems.
Simple interface for fast problem set-up.
Definition Simple.hpp:27

◆ createCommonData()

MoFEMErrorCode Seepage::createCommonData ( )
private

[Set up problem]

[Create common data]

Examples
mofem/tutorials/adv-5_poroelasticity/seepage.cpp.

Definition at line 611 of file seepage.cpp.

611 {
613
614 auto get_command_line_parameters = [&]() {
616 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "", "-young_modulus",
617 &default_young_modulus, PETSC_NULLPTR);
618 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "", "-poisson_ratio",
619 &default_poisson_ratio, PETSC_NULLPTR);
620 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "", "-conductivity",
621 &default_conductivity, PETSC_NULLPTR);
622 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "", "-biot_constant",
623 &default_biot, PETSC_NULLPTR);
624 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "", "-mat_density",
625 &default_mat_density, PETSC_NULLPTR);
626 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "", "-fluid_density",
627 &default_fluid_density, PETSC_NULLPTR);
628 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "", "-storage",
629 &default_storage, PETSC_NULLPTR);
630
631 MOFEM_LOG("Seepage", Sev::inform)
632 << "Default Young modulus " << default_young_modulus;
633 MOFEM_LOG("Seepage", Sev::inform)
634 << "Default Poisson ratio " << default_poisson_ratio;
635 MOFEM_LOG("Seepage", Sev::inform)
636 << "Default conductivity " << default_conductivity;
637 MOFEM_LOG("Seepage", Sev::inform) << "Default biot " << default_biot;
638 MOFEM_LOG("Seepage", Sev::inform)
639 << "Default fluid density " << default_fluid_density;
640 MOFEM_LOG("Seepage", Sev::inform)
641 << "Default material density " << default_mat_density;
642 MOFEM_LOG("Seepage", Sev::inform)
643 << "Default storage constant " << default_storage;
644
646 };
647
648 CHKERR get_command_line_parameters();
649
651}
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)

◆ OPs()

MoFEMErrorCode Seepage::OPs ( )
private

[Boundary condition]

[Push operators to pipeline]

[LHS SOLID MATRIX OP]

[LHS SOLID MATRIX OP]

[RHS SOLID MATRIX OPERATORS]

[RHS SOLID MATRIX OPERATORS]

[LHS FLUID MATRIX OPERATORS]

[LHS FLUID MATRIX OPERATORS]

[RHS FLUID MATRIX OPERATORS]

[RHS FLUID MATRIX OPERATORS]

[Assembleops]

Examples
mofem/tutorials/adv-5_poroelasticity/seepage.cpp.

Definition at line 756 of file seepage.cpp.

756 {
758 auto pipeline_mng = mField.getInterface<PipelineManager>();
760 auto bc_mng = mField.getInterface<BcManager>();
761
762 auto boundary_marker =
763 bc_mng->getMergedBlocksMarker(vector<string>{"FLUIDFLUX"});
764
765 auto u_grad_ptr = boost::make_shared<MatrixDouble>();
766 auto dot_u_grad_ptr = boost::make_shared<MatrixDouble>();
767 auto trace_dot_u_grad_ptr = boost::make_shared<VectorDouble>();
768 auto p_ptr = boost::make_shared<VectorDouble>();
769 auto dot_p_ptr = boost::make_shared<VectorDouble>();
770 auto flux_ptr = boost::make_shared<MatrixDouble>();
771 auto div_flux_ptr = boost::make_shared<VectorDouble>();
772 auto strain_ptr = boost::make_shared<MatrixDouble>();
773 auto stress_ptr = boost::make_shared<MatrixDouble>();
774
775 auto time_scale = boost::make_shared<TimeScale>();
776
777 auto block_params = boost::make_shared<BlockedParameters>();
778 auto mDPtr = block_params->getDPtr();
779 auto conductivity_ptr = block_params->getConductivityPtr();
780 auto biot_constant_ptr = block_params->getBiotConstantPtr();
781 auto mat_density_ptr = block_params->getMatDensityPtr();
782 auto fluid_density_ptr = block_params->getFluidDensityPtr();
783 auto storage_constant_ptr = block_params->getStorageConstantPtr();
784
785 auto integration_rule = [](int, int, int approx_order) {
786 return 2 * approx_order;
787 };
788
789 CHKERR pipeline_mng->setDomainRhsIntegrationRule(integration_rule);
790 CHKERR pipeline_mng->setDomainLhsIntegrationRule(integration_rule);
791 CHKERR pipeline_mng->setBoundaryLhsIntegrationRule(integration_rule);
792 CHKERR pipeline_mng->setBoundaryRhsIntegrationRule(integration_rule);
793
794 auto add_domain_base_ops = [&](auto &pip) {
796
797 CHKERR addMatBlockOps(pip, "MAT_ELASTIC", "MAT_FLUID", block_params,
798 Sev::inform);
800
802 "U", u_grad_ptr));
803 pip.push_back(new OpSymmetrizeTensor<SPACE_DIM>(u_grad_ptr, strain_ptr));
805 "U", dot_u_grad_ptr));
806 pip.push_back(new OpCalculateTraceFromMat<SPACE_DIM>(dot_u_grad_ptr,
807 trace_dot_u_grad_ptr));
808
809 pip.push_back(new OpCalculateScalarFieldValues("P", p_ptr));
810 pip.push_back(new OpCalculateScalarFieldValuesDot("P", dot_p_ptr));
811 pip.push_back(new OpCalculateHVecVectorField<3>("FLUX", flux_ptr));
813 "FLUX", div_flux_ptr));
814
816 };
817
818 //! [LHS SOLID MATRIX OP]
819 auto add_domain_ops_lhs_mechanical = [&](auto &pip) {
821
822 auto biot = [biot_constant_ptr](const double, const double, const double) {
823 return -*biot_constant_ptr;
824 };
825
826 pip.push_back(new OpKCauchy("U", "U", mDPtr));
827 pip.push_back(new OpBaseDivU("P", "U", biot, true, true));
829 };
830 //! [LHS SOLID MATRIX OP]
831
832 //! [RHS SOLID MATRIX OPERATORS]
833 auto add_domain_ops_rhs_mechanical = [&](auto &pip) {
835
837 pip, mField, "U", {time_scale}, "BODY_FORCE", Sev::inform);
838 auto biot = [biot_constant_ptr](const double, const double, const double) {
839 return *biot_constant_ptr;
840 };
841 auto mat_density = [mat_density_ptr](const double, const double,
842 const double) {
843 return *mat_density_ptr;
844 };
845 // Calculate internal force
847 strain_ptr, stress_ptr, mDPtr));
848 pip.push_back(new OpInternalForceCauchy("U", stress_ptr));
850 "U", p_ptr, biot));
851
853 };
854 //! [RHS SOLID MATRIX OPERATORS]
855
856
857 //! [LHS FLUID MATRIX OPERATORS]
858 auto add_domain_ops_lhs_seepage = [&](auto &pip, auto &fe) {
860 auto resistance = [conductivity_ptr](const double, const double,
861 const double) {
862 return (1. / *(conductivity_ptr));
863 };
864 auto biot = [biot_constant_ptr](const double, const double, const double) {
865 return *biot_constant_ptr;
866 };
867 auto storage = [storage_constant_ptr](const double, const double,
868 const double) {
869 return -*storage_constant_ptr;
870 };
871 auto minus_one = []() constexpr { return -1; };
872
873 pip.push_back(new OpSetBc("FLUX", true, boundary_marker));
874
875 pip.push_back(new OpHdivHdiv("FLUX", "FLUX", resistance));
876 pip.push_back(new OpHdivQ("FLUX", "P", minus_one, true));
877
878 auto op_domain_mass = new OpDomainMass("P", "P", storage);
879 op_domain_mass->feScalingFun = [](const FEMethod *fe_ptr) {
880 return fe_ptr->ts_a;
881 };
882 pip.push_back(op_domain_mass);
883
884 auto op_base_div_u = new OpBaseDivU("P", "U", biot, false, false);
885 op_base_div_u->feScalingFun = [](const FEMethod *fe_ptr) {
886 return fe_ptr->ts_a;
887 };
888 pip.push_back(op_base_div_u);
889
890 pip.push_back(new OpUnSetBc("FLUX"));
891
893 };
894 //! [LHS FLUID MATRIX OPERATORS]
895
896 //! [RHS FLUID MATRIX OPERATORS]
897 auto add_domain_ops_rhs_seepage = [&](auto &pip) {
899 auto resistance = [conductivity_ptr](double, double, double) {
900 return (1. / (*conductivity_ptr));
901 };
902 auto minus_one = [](double, double, double) constexpr { return -1; };
903 auto biot = [biot_constant_ptr](const double, const double, const double) {
904 return *biot_constant_ptr;
905 };
906
907 auto storage = [storage_constant_ptr](const double, const double,
908 const double) {
909 return -*storage_constant_ptr;
910 };
911
912 auto fluid_density = [fluid_density_ptr](const double, const double,
913 const double) {
914 return *fluid_density_ptr;
915 };
916
917 auto dot_p_at_gauss_pts = boost::make_shared<VectorDouble>();
918 pip.push_back(new OpCalculateScalarFieldValuesDot("P", dot_p_at_gauss_pts));
919
920 pip.push_back(new OpSetBc("FLUX", true, boundary_marker));
921
922 pip.push_back(new OpHdivFlux("FLUX", flux_ptr, resistance));
923 pip.push_back(new OpHDivH("FLUX", p_ptr, minus_one));
924 pip.push_back(
925 new OpDomainTimesScalarField("P", dot_p_at_gauss_pts, storage));
926 pip.push_back(new OpBaseDotH("P", trace_dot_u_grad_ptr, biot));
927 pip.push_back(new OpBaseDivFlux("P", div_flux_ptr, minus_one));
928
929 pip.push_back(new OpUnSetBc("FLUX"));
930
932 };
933 //! [RHS FLUID MATRIX OPERATORS]
934
935 auto add_boundary_rhs_ops = [&](auto &pip) {
937
939
940 pip.push_back(new OpSetBc("FLUX", true, boundary_marker));
941
943 pipeline_mng->getOpBoundaryRhsPipeline(), mField, "U", {time_scale},
944 "FORCE", Sev::inform);
945
947 pipeline_mng->getOpBoundaryRhsPipeline(), mField, "FLUX", {time_scale},
948 "PRESSURE", Sev::inform);
949
950 pip.push_back(new OpUnSetBc("FLUX"));
951
952 auto mat_flux_ptr = boost::make_shared<MatrixDouble>();
953 pip.push_back(
954 new OpCalculateHVecVectorField<3, SPACE_DIM>("FLUX", mat_flux_ptr));
957 mField, pip, simple->getProblemName(), "FLUX", mat_flux_ptr,
958 {time_scale});
959
961 };
962
963 auto add_boundary_lhs_ops = [&](auto &pip) {
965
967
970 mField, pip, simple->getProblemName(), "FLUX");
971
973 };
974
975 //! [Assembleops]
976 // LHS
977 CHKERR add_domain_base_ops(pipeline_mng->getOpDomainLhsPipeline());
978 CHKERR add_domain_ops_lhs_mechanical(pipeline_mng->getOpDomainLhsPipeline());
979 CHKERR add_domain_ops_lhs_seepage(pipeline_mng->getOpDomainLhsPipeline(),
980 pipeline_mng->getDomainLhsFE());
981
982 // RHS
983 CHKERR add_domain_base_ops(pipeline_mng->getOpDomainRhsPipeline());
984 CHKERR add_domain_ops_rhs_mechanical(pipeline_mng->getOpDomainRhsPipeline());
985 CHKERR add_domain_ops_rhs_seepage(pipeline_mng->getOpDomainRhsPipeline());
986
987 CHKERR add_boundary_rhs_ops(pipeline_mng->getOpBoundaryRhsPipeline());
988 CHKERR add_boundary_lhs_ops(pipeline_mng->getOpBoundaryLhsPipeline());
989
991}
@ H1
continuous field
Definition definitions.h:85
@ HDIV
field with continuous normal traction
Definition definitions.h:87
auto integration_rule
OpCalculateScalarFieldValuesFromPetscVecImpl< PetscData::CTX_SET_X_T > OpCalculateScalarFieldValuesDot
static constexpr int approx_order
OpBaseDotH OpBaseDivFlux
Integrate Rhs base of temperature times divergent of flux (T)
Definition seepage.cpp:133
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, 1 > OpDomainMass
Definition seepage.cpp:139
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesScalar< 1 > OpDomainTimesScalarField
Definition seepage.cpp:141
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesSymTensor< 1, SPACE_DIM, SPACE_DIM > OpInternalForceCauchy
Definition seepage.cpp:51
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMixScalarTimesDiv< SPACE_DIM > OpBaseDivU
[Essential boundary conditions]
Definition seepage.cpp:78
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpMixDivTimesU< 3, 1, SPACE_DIM > OpHDivH
Integrate Rhs div flux base times temperature (T)
Definition seepage.cpp:119
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesVector< 3, 3, 1 > OpHdivFlux
Integrating Rhs flux base (1/k) flux (FLUX)
Definition seepage.cpp:108
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMixDivTimesScalar< SPACE_DIM > OpHdivQ
Integrate Lhs div of base of flux time base of temperature (FLUX x T) and transpose of it,...
Definition seepage.cpp:94
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 3, 3 > OpHdivHdiv
Integrate Lhs base of flux (1/k) base of flux (FLUX x FLUX)
Definition seepage.cpp:86
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradSymTensorGrad< 1, SPACE_DIM, SPACE_DIM, 0 > OpKCauchy
[Only used with Hooke equation (linear material model)]
Definition seepage.cpp:49
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesScalarField< 1 > OpBaseDotH
Integrate Rhs base of temperature time heat capacity times heat rate (T)
Definition seepage.cpp:127
Add operators pushing bases from local to physical configuration.
Essential boundary conditions.
Structure for user loop methods on finite elements.
Get vector field for H-div approximation.
Calculate divergence of vector field.
Specialization for double precision scalar field values calculation.
Operator for calculating the trace of matrices at integration points.
Get field gradients time derivative at integration pts for scalar field rank 0, i....
Get field gradients at integration pts for scalar field rank 0, i.e. vector field.
Set indices on entities on finite element.
Operator for symmetrizing tensor fields.
Unset indices on entities on finite element.
PipelineManager interface.
boost::shared_ptr< MatrixDouble > mDPtr
Definition seepage.cpp:223
MoFEMErrorCode addMatBlockOps(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, std::string block_elastic_name, std::string block_thermal_name, boost::shared_ptr< BlockedParameters > blockedParamsPtr, Sev sev)
Definition seepage.cpp:265

◆ runProblem()

MoFEMErrorCode Seepage::runProblem ( )

[Run problem]

Examples
mofem/tutorials/adv-5_poroelasticity/seepage.cpp.

Definition at line 513 of file seepage.cpp.

513 {
517 CHKERR bC();
518 CHKERR OPs();
519 CHKERR tsSolve();
521}
MoFEMErrorCode setupProblem()
[Run problem]
Definition seepage.cpp:525
MoFEMErrorCode createCommonData()
[Set up problem]
Definition seepage.cpp:611
MoFEMErrorCode bC()
[Create common data]
Definition seepage.cpp:655
MoFEMErrorCode OPs()
[Boundary condition]
Definition seepage.cpp:756
MoFEMErrorCode tsSolve()
[Solve]
Definition seepage.cpp:996

◆ setupProblem()

MoFEMErrorCode Seepage::setupProblem ( )
private

[Run problem]

[Set up problem]

Examples
mofem/tutorials/adv-5_poroelasticity/seepage.cpp.

Definition at line 525 of file seepage.cpp.

525 {
528 // Add field
530 // Mechanical fields
532 CHKERR simple->addBoundaryField("U", H1, base, SPACE_DIM);
533 // Fluid
534 const auto flux_space = (SPACE_DIM == 2) ? HCURL : HDIV;
535 CHKERR simple->addDomainField("P", L2, AINSWORTH_LEGENDRE_BASE, 1);
536 CHKERR simple->addDomainField("FLUX", flux_space, DEMKOWICZ_JACOBI_BASE, 1);
537 CHKERR simple->addBoundaryField("FLUX", flux_space, DEMKOWICZ_JACOBI_BASE, 1);
538
539 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-plane_strain",
540 &is_plane_strain, PETSC_NULLPTR);
541
542 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-save_every", &save_every,
543 PETSC_NULLPTR);
544
545 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-atom_test", &atom_test,
546 PETSC_NULLPTR);
547 int order = 2.;
548 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-order", &order, PETSC_NULLPTR);
549 CHKERR simple->setFieldOrder("U", order);
550 CHKERR simple->setFieldOrder("P", order - 1);
551 CHKERR simple->setFieldOrder("FLUX", order);
552
553 int coords_dim = 3;
554 CHKERR PetscOptionsGetRealArray(NULL, NULL, "-field_eval_coords",
555 fieldEvalCoords.data(), &coords_dim,
556 &doEvalField);
557
558 CHKERR simple->setUp();
559
560 pressureFieldPtr = boost::make_shared<VectorDouble>();
561 fluxFieldPtr = boost::make_shared<MatrixDouble>();
562 dispFieldPtr = boost::make_shared<MatrixDouble>();
563 dispGradPtr = boost::make_shared<MatrixDouble>();
564 strainFieldPtr = boost::make_shared<MatrixDouble>();
565 stressFieldPtr = boost::make_shared<MatrixDouble>();
566 mDPtr = boost::make_shared<MatrixDouble>();
567
568 if (doEvalField) {
570 mField.getInterface<FieldEvaluatorInterface>()->getData<DomainEle>();
571
572 CHKERR mField.getInterface<FieldEvaluatorInterface>()->buildTree<SPACE_DIM>(
573 fieldEvalData, simple->getDomainFEName());
574
575 fieldEvalData->setEvalPoints(fieldEvalCoords.data(), 1);
576 auto no_rule = [](int, int, int) { return -1; };
577
578 auto field_eval_fe_ptr = fieldEvalData->feMethodPtr.lock();
579 field_eval_fe_ptr->getRuleHook = no_rule;
580
581 auto block_params = boost::make_shared<BlockedParameters>();
582 mDPtr = block_params->getDPtr();
583
584 CHKERR addMatBlockOps(field_eval_fe_ptr->getOpPtrVector(), "MAT_ELASTIC",
585 "MAT_FLUID", block_params, Sev::verbose);
586
588 field_eval_fe_ptr->getOpPtrVector(), {H1, HDIV});
589
590 field_eval_fe_ptr->getOpPtrVector().push_back(
592 field_eval_fe_ptr->getOpPtrVector().push_back(
594 field_eval_fe_ptr->getOpPtrVector().push_back(
596 field_eval_fe_ptr->getOpPtrVector().push_back(
598 dispGradPtr));
599 field_eval_fe_ptr->getOpPtrVector().push_back(
601 field_eval_fe_ptr->getOpPtrVector().push_back(
604 }
605
607}
FieldApproximationBase
approximation base
Definition definitions.h:58
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ DEMKOWICZ_JACOBI_BASE
Definition definitions.h:66
@ L2
field with C-1 continuity
Definition definitions.h:88
@ HCURL
field with continuous tangents
Definition definitions.h:86
constexpr int order
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
PetscErrorCode PetscOptionsGetRealArray(PetscOptions *, const char pre[], const char name[], PetscReal dval[], PetscInt *nmax, PetscBool *set)
int save_every
Definition seepage.cpp:179
int atom_test
Definition seepage.cpp:180
Field evaluator interface.
Specialization for double precision vector field values calculation.
MoFEMErrorCode addDomainField(const std::string name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_ZERO, int verb=-1)
Add field on domain.
Definition Simple.cpp:261
std::array< double, 3 > fieldEvalCoords
Definition seepage.cpp:210
boost::shared_ptr< VectorDouble > pressureFieldPtr
Definition seepage.cpp:217
boost::shared_ptr< FieldEvaluatorInterface::SetPtsData > fieldEvalData
Definition seepage.cpp:211
boost::shared_ptr< MatrixDouble > dispFieldPtr
Definition seepage.cpp:219
boost::shared_ptr< MatrixDouble > fluxFieldPtr
Definition seepage.cpp:218
PetscBool doEvalField
Definition seepage.cpp:209
boost::shared_ptr< MatrixDouble > stressFieldPtr
Definition seepage.cpp:222
boost::shared_ptr< MatrixDouble > strainFieldPtr
Definition seepage.cpp:221
boost::shared_ptr< MatrixDouble > dispGradPtr
Definition seepage.cpp:220

◆ tsSolve()

MoFEMErrorCode Seepage::tsSolve ( )
private

[Solve]

[Assembleops] [Push operators to pipeline]

Examples
mofem/tutorials/adv-5_poroelasticity/seepage.cpp.

Definition at line 996 of file seepage.cpp.

996 {
1000
1001 auto dm = simple->getDM();
1002 auto solver = pipeline_mng->createTSIM();
1003 auto snes_ctx_ptr = getDMSnesCtx(dm);
1004
1005 auto set_section_monitor = [&](auto solver) {
1007 SNES snes;
1008 CHKERR TSGetSNES(solver, &snes);
1009 CHKERR SNESMonitorSet(snes,
1010 (MoFEMErrorCode (*)(SNES, PetscInt, PetscReal,
1011 void *))MoFEMSNESMonitorFields,
1012 (void *)(snes_ctx_ptr.get()), nullptr);
1014 };
1015
1016 auto create_post_process_element = [&]() {
1017 auto post_proc_fe = boost::make_shared<PostProcEle>(mField);
1018
1019 auto block_params = boost::make_shared<BlockedParameters>();
1020 auto mD_ptr = block_params->getDPtr();
1021 CHKERR addMatBlockOps(post_proc_fe->getOpPtrVector(), "MAT_ELASTIC",
1022 "MAT_FLUID", block_params, Sev::verbose);
1024 post_proc_fe->getOpPtrVector(), {H1, HDIV});
1025
1026 auto mat_grad_ptr = boost::make_shared<MatrixDouble>();
1027 auto mat_strain_ptr = boost::make_shared<MatrixDouble>();
1028 auto mat_stress_ptr = boost::make_shared<MatrixDouble>();
1029
1030 auto p_ptr = boost::make_shared<VectorDouble>();
1031 auto mat_flux_ptr = boost::make_shared<MatrixDouble>();
1032
1033 post_proc_fe->getOpPtrVector().push_back(
1034 new OpCalculateScalarFieldValues("P", p_ptr));
1035 post_proc_fe->getOpPtrVector().push_back(
1036 new OpCalculateHVecVectorField<3, SPACE_DIM>("FLUX", mat_flux_ptr));
1037
1038 auto u_ptr = boost::make_shared<MatrixDouble>();
1039 post_proc_fe->getOpPtrVector().push_back(
1041 post_proc_fe->getOpPtrVector().push_back(
1043 mat_grad_ptr));
1044 post_proc_fe->getOpPtrVector().push_back(
1045 new OpSymmetrizeTensor<SPACE_DIM>(mat_grad_ptr, mat_strain_ptr));
1046 post_proc_fe->getOpPtrVector().push_back(
1048 mat_strain_ptr, mat_stress_ptr, mD_ptr));
1049
1051
1052 post_proc_fe->getOpPtrVector().push_back(
1053
1054 new OpPPMap(
1055
1056 post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
1057
1058 {{"P", p_ptr}},
1059
1060 {{"U", u_ptr}, {"FLUX", mat_flux_ptr}},
1061
1062 {},
1063
1064 {{"STRAIN", mat_strain_ptr}, {"STRESS", mat_stress_ptr}}
1065
1066 )
1067
1068 );
1069
1070 return post_proc_fe;
1071 };
1072
1073 auto create_creaction_fe = [&]() {
1074 auto fe_ptr = boost::make_shared<DomainEle>(mField);
1075 fe_ptr->getRuleHook = [](int, int, int o) { return 2 * o; };
1076
1077 auto &pip = fe_ptr->getOpPtrVector();
1078
1079 auto block_params = boost::make_shared<BlockedParameters>();
1080 auto mD_ptr = block_params->getDPtr();
1081 CHKERR addMatBlockOps(pip, "MAT_ELASTIC", "MAT_FLUID", block_params,
1082 Sev::verbose);
1084
1085 auto u_grad_ptr = boost::make_shared<MatrixDouble>();
1086 auto strain_ptr = boost::make_shared<MatrixDouble>();
1087 auto stress_ptr = boost::make_shared<MatrixDouble>();
1088
1090 "U", u_grad_ptr));
1091 pip.push_back(new OpSymmetrizeTensor<SPACE_DIM>(u_grad_ptr, strain_ptr));
1092
1093 // Calculate internal force
1095 strain_ptr, stress_ptr, mD_ptr));
1096 pip.push_back(new OpInternalForceCauchy("U", stress_ptr));
1097
1098 fe_ptr->postProcessHook =
1100
1101 return fe_ptr;
1102 };
1103
1104 auto monitor_ptr = boost::make_shared<FEMethod>();
1105 auto post_proc_fe = create_post_process_element();
1106 auto res = createDMVector(dm);
1107 auto rections_fe = create_creaction_fe();
1108
1109 auto set_time_monitor = [&](auto dm, auto solver) {
1111 monitor_ptr->preProcessHook = [&]() {
1113
1114 if (save_every && (monitor_ptr->ts_step % save_every == 0)) {
1115
1116 CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(),
1117 post_proc_fe,
1118 monitor_ptr->getCacheWeakPtr());
1119 CHKERR post_proc_fe->writeFile(
1120 "out_" + boost::lexical_cast<std::string>(monitor_ptr->ts_step) +
1121 ".h5m");
1122 }
1123
1124 rections_fe->f = res;
1125 CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(),
1126 rections_fe,
1127 monitor_ptr->getCacheWeakPtr());
1128
1129 CHKERR VecAssemblyBegin(res);
1130 CHKERR VecAssemblyEnd(res);
1131 double nrm;
1132 CHKERR VecNorm(res, NORM_2, &nrm);
1133 MOFEM_LOG("Seepage", Sev::verbose)
1134 << "Residual norm " << nrm << " at time step "
1135 << monitor_ptr->ts_step;
1136
1137 struct AtomTestResult {
1138 bool fail = false;
1139 std::string msg = "";
1140 };
1141
1142 AtomTestResult atom_test_result;
1143
1144 auto fail_atom_test = [&atom_test_result](const std::string &msg) {
1145 if (!atom_test_result.fail) {
1146 atom_test_result.fail = true;
1147 atom_test_result.msg = msg;
1148 }
1149 };
1150
1151
1152 struct AtomTestData {
1153 double expected = 0.0;
1154 double tol = 0.0;
1155 };
1156
1157
1158 if (doEvalField) {
1159
1161 ->evalFEAtThePoint<SPACE_DIM>(
1162 fieldEvalCoords.data(), 1e-12, simple->getProblemName(),
1163 simple->getDomainFEName(), fieldEvalData,
1165 MF_EXIST, QUIET);
1166
1167 if (atom_test) {
1168 auto eval_num_vec =
1169 createVectorMPI(mField.get_comm(), PETSC_DECIDE, 1);
1170 CHKERR VecZeroEntries(eval_num_vec);
1171 if (pressureFieldPtr->size()) {
1172 CHKERR VecSetValue(eval_num_vec, 0, 1, ADD_VALUES);
1173 }
1174 CHKERR VecAssemblyBegin(eval_num_vec);
1175 CHKERR VecAssemblyEnd(eval_num_vec);
1176
1177 double eval_num;
1178 CHKERR VecSum(eval_num_vec, &eval_num);
1179 if (!(int)eval_num) {
1180 fail_atom_test(
1181 "did not find elements to evaluate the field, check the "
1182 "coordinates");
1183 }
1184 }
1185
1186 if (atom_test && fabs(monitor_ptr->ts_t - 0.1) < 1e-12 &&
1187 pressureFieldPtr && !pressureFieldPtr->empty()) {
1188 auto t_pressure = getFTensor0FromVec(*pressureFieldPtr);
1189 MOFEM_LOG("SeepageSync", Sev::inform)
1190 << "Rank " << mField.get_comm_rank();
1191 MOFEM_LOG("SeepageSync", Sev::inform)
1192 << "Eval point P: " << t_pressure;
1193 double pressure = 0.0;
1194
1196 auto t_flux = getFTensor1FromMat<SPACE_DIM>(*fluxFieldPtr);
1197 auto flux_mag = sqrt(t_flux(i) * t_flux(i));
1198 MOFEM_LOG("SeepageSync", Sev::inform)
1199 << "Eval point FLUX magnitude: " << flux_mag;
1200 double flux = 0.0;
1201
1202 auto t_disp = getFTensor1FromMat<SPACE_DIM>(*dispFieldPtr);
1203 auto disp_mag = sqrt(t_disp(i) * t_disp(i));
1204 MOFEM_LOG("SeepageSync", Sev::inform)
1205 << "Eval point U magnitude: " << disp_mag;
1206 double disp = 0.0;
1207
1208 auto t_strain = getFTensor2SymmetricFromMat<SPACE_DIM>(*strainFieldPtr);
1209 auto t_strain_trace = t_strain(i, i);
1210 MOFEM_LOG("SeepageSync", Sev::inform)
1211 << "Eval point trace of strain: " << t_strain_trace;
1212 double strain = 0.0;
1213
1214 double von_mises_stress;
1215 auto getVonMisesStress = [&](auto t_stress) {
1217 von_mises_stress = std::sqrt(
1218 0.5 *
1219 ((t_stress(0, 0) - t_stress(1, 1)) *
1220 (t_stress(0, 0) - t_stress(1, 1)) +
1221 (SPACE_DIM == 3 ? (t_stress(1, 1) - t_stress(2, 2)) *
1222 (t_stress(1, 1) - t_stress(2, 2))
1223 : 0) +
1224 (SPACE_DIM == 3 ? (t_stress(2, 2) - t_stress(0, 0)) *
1225 (t_stress(2, 2) - t_stress(0, 0))
1226 : 0) +
1227 6 * (t_stress(0, 1) * t_stress(0, 1) +
1228 (SPACE_DIM == 3 ? t_stress(1, 2) * t_stress(1, 2) : 0) +
1229 (SPACE_DIM == 3 ? t_stress(2, 0) * t_stress(2, 0) : 0))));
1231 };
1232 auto t_stress =
1233 getFTensor2SymmetricFromMat<SPACE_DIM>(*stressFieldPtr);
1234 CHKERR getVonMisesStress(t_stress);
1235 MOFEM_LOG("SeepageSync", Sev::inform)
1236 << "Eval point von Mises Stress: " << von_mises_stress;
1237 double stress = 0.0;
1238
1239 switch (atom_test) {
1240 case 1: // terzaghi_2d
1241 pressure = 0.0002681;
1242 flux = 0.0001734;
1243 disp = 0.0009;
1244 strain = 0.0009;
1245 stress = 53.037;
1246 break;
1247
1248 case 2: // terzaghi_3d
1249 pressure = 0.0063046;
1250 flux = 0.000093771;
1251 disp = 0.000859577;
1252 strain = 0.0009;
1253 stress = 75.00;
1254 break;
1255
1256 case 3: // mandel
1257 pressure = 0.31473;
1258 flux = 0.0009216;
1259 disp = 0.0135;
1260 strain = 0.036392;
1261 stress = 356.956;
1262 break;
1263
1264 default:
1265 fail_atom_test("unknown atom test number");
1266 }
1267 if (fabs(t_pressure + pressure) > 1e-4) {
1268 fail_atom_test("wrong pressure value");
1269 }
1270 if (fabs(flux_mag - flux) > 1e-6) {
1271 fail_atom_test("wrong flux value");
1272 }
1273 if (fabs(disp_mag - disp) > 1e-2) {
1274 fail_atom_test("wrong displacement value");
1275 }
1276 if (fabs(t_strain_trace + strain) > 1e-3) {
1277 fail_atom_test("wrong strain value");
1278 }
1279 if (fabs(von_mises_stress - stress) > 1e-2) {
1280 fail_atom_test("wrong stress value");
1281 }
1282 }
1283
1284 auto current_time_step = monitor_ptr->ts_step; // Time step
1285 PetscReal current_time;
1286 CHKERR TSGetTime(solver, &current_time);
1287
1288 if (pressureFieldPtr && !pressureFieldPtr->empty()) {
1289
1290 MOFEM_LOG("SeepageSync", Sev::inform)
1291 << "Eval point displacement: " << *dispFieldPtr;
1292 MOFEM_LOG("SeepageSync", Sev::inform)
1293 << "Eval point strain: " << *strainFieldPtr;
1294 MOFEM_LOG("SeepageSync", Sev::inform)
1295 << "Eval point symmetric stress tensor: " << *stressFieldPtr;
1296 MOFEM_LOG("SeepageSync", Sev::inform)
1297 << "Eval point pore pressure: " << *pressureFieldPtr;
1298 MOFEM_LOG("SeepageSync", Sev::inform)
1299 << "Eval point flux: " << *fluxFieldPtr;
1300 MOFEM_LOG("SeepageSync", Sev::inform)
1301 << "Eval point elasticity matrix: " << *mDPtr;
1302 }
1303 MOFEM_LOG_SEVERITY_SYNC(mField.get_comm(), Sev::inform);
1304 }
1305
1306 int fail_flag = atom_test_result.fail ? 1 : 0;
1307 MPI_Allreduce(MPI_IN_PLACE, &fail_flag, 1, MPI_INT, MPI_MAX, mField.get_comm());
1308
1309 if (fail_flag) {
1310 const int root = 0;
1311 const int MAX_MSG = 512;
1312
1313 char msg_buf[MAX_MSG];
1314 msg_buf[0] = '\0';
1315
1316 if (atom_test_result.fail) {
1317 std::snprintf(msg_buf, MAX_MSG, "%s", atom_test_result.msg.c_str());
1318 }
1319 int my_rank = mField.get_comm_rank();
1320 int send_rank = atom_test_result.fail ? my_rank : INT_MAX;
1321 MPI_Allreduce(MPI_IN_PLACE, &send_rank, 1, MPI_INT, MPI_MIN, mField.get_comm());
1322 MPI_Bcast(msg_buf, MAX_MSG, MPI_CHAR, send_rank, mField.get_comm());
1323
1324 if (my_rank == root) {
1325 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
1326 "atom test %d failed: %s", atom_test, msg_buf);
1327 }
1328
1330 "Atom test failed");
1331 }
1332
1334 };
1335
1336 auto null = boost::shared_ptr<FEMethod>();
1337 CHKERR DMMoFEMTSSetMonitor(dm, solver, simple->getDomainFEName(), null,
1338 monitor_ptr, null);
1340 };
1341
1342 auto set_fieldsplit_preconditioner = [&](auto solver) {
1344
1345 SNES snes;
1346 CHKERR TSGetSNES(solver, &snes);
1347 KSP ksp;
1348 CHKERR SNESGetKSP(snes, &ksp);
1349 PC pc;
1350 CHKERR KSPGetPC(ksp, &pc);
1351 PetscBool is_pcfs = PETSC_FALSE;
1352 PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &is_pcfs);
1353
1354 // Setup fieldsplit (block) solver - optional: yes/no
1355 if (is_pcfs == PETSC_TRUE) {
1356 auto bc_mng = mField.getInterface<BcManager>();
1357 auto is_mng = mField.getInterface<ISManager>();
1358 auto name_prb = simple->getProblemName();
1359
1360 SmartPetscObj<IS> is_u;
1361 CHKERR is_mng->isCreateProblemFieldAndRank(name_prb, ROW, "U", 0,
1362 SPACE_DIM, is_u);
1363 SmartPetscObj<IS> is_flux;
1364 CHKERR is_mng->isCreateProblemFieldAndRank(name_prb, ROW, "FLUX", 0, 0,
1365 is_flux);
1366 SmartPetscObj<IS> is_P;
1367 CHKERR is_mng->isCreateProblemFieldAndRank(name_prb, ROW, "P", 0, 0,
1368 is_P);
1369 IS is_tmp;
1370 CHKERR ISExpand(is_P, is_flux, &is_tmp);
1371 auto is_Flux = SmartPetscObj<IS>(is_tmp);
1372
1373 auto is_all_bc = bc_mng->getBlockIS(name_prb, "FLUIDFLUX", "FLUX", 0, 0);
1374 int is_all_bc_size;
1375 CHKERR ISGetSize(is_all_bc, &is_all_bc_size);
1376 MOFEM_LOG("ThermoElastic", Sev::inform)
1377 << "Field split block size " << is_all_bc_size;
1378 if (is_all_bc_size) {
1379 IS is_tmp2;
1380 CHKERR ISDifference(is_Flux, is_all_bc, &is_tmp2);
1381 is_Flux = SmartPetscObj<IS>(is_tmp2);
1382 CHKERR PCFieldSplitSetIS(pc, PETSC_NULLPTR,
1383 is_all_bc); // boundary block
1384 }
1385
1386 CHKERR ISSort(is_u);
1387 CHKERR ISSort(is_Flux);
1388 CHKERR PCFieldSplitSetIS(pc, PETSC_NULLPTR, is_Flux);
1389 CHKERR PCFieldSplitSetIS(pc, PETSC_NULLPTR, is_u);
1390 }
1391
1393 };
1394
1395 auto pre_proc_ptr = boost::make_shared<FEMethod>();
1396 auto post_proc_rhs_ptr = boost::make_shared<FEMethod>();
1397 auto post_proc_lhs_ptr = boost::make_shared<FEMethod>();
1398 auto time_scale = boost::make_shared<TimeScale>();
1399
1400 auto get_bc_hook_rhs = [this, pre_proc_ptr, time_scale]() {
1402 {time_scale}, false);
1403 return hook;
1404 };
1405
1406 auto get_post_proc_hook_rhs = [this, post_proc_rhs_ptr]() {
1409 mField, post_proc_rhs_ptr, nullptr, Sev::verbose)();
1411 mField, post_proc_rhs_ptr, 1.)();
1413 };
1414 auto get_post_proc_hook_lhs = [this, post_proc_lhs_ptr]() {
1416 post_proc_lhs_ptr, 1.);
1417 };
1418
1419 pre_proc_ptr->preProcessHook = get_bc_hook_rhs();
1420 post_proc_rhs_ptr->postProcessHook = get_post_proc_hook_rhs;
1421 post_proc_lhs_ptr->postProcessHook = get_post_proc_hook_lhs();
1422
1423 auto ts_ctx_ptr = getDMTsCtx(dm);
1424 ts_ctx_ptr->getPreProcessIFunction().push_front(pre_proc_ptr);
1425 ts_ctx_ptr->getPreProcessIJacobian().push_front(pre_proc_ptr);
1426 ts_ctx_ptr->getPostProcessIFunction().push_back(post_proc_rhs_ptr);
1427 ts_ctx_ptr->getPostProcessIJacobian().push_back(post_proc_lhs_ptr);
1428
1429 auto B = createDMMatrix(dm);
1430 CHKERR TSSetIJacobian(solver, B, B, PETSC_NULLPTR, PETSC_NULLPTR);
1431 auto D = createDMVector(dm);
1432 CHKERR TSSetSolution(solver, D);
1433 CHKERR TSSetFromOptions(solver);
1434
1435 CHKERR set_section_monitor(solver);
1436 CHKERR set_fieldsplit_preconditioner(solver);
1437 CHKERR set_time_monitor(dm, solver);
1438
1439 CHKERR TSSetUp(solver);
1440 CHKERR TSSolve(solver, NULL);
1441
1443}
@ QUIET
@ ROW
@ MF_EXIST
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
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)
Get smart vector from DM.
Definition DMMoFEM.hpp:1234
auto createDMMatrix(DM dm)
Get smart matrix from DM.
Definition DMMoFEM.hpp:1191
SmartPetscObj< TS > createTSIM(SmartPetscObj< DM > dm=nullptr)
Create TS (time) implicit solver.
double D
double tol
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:1276
MoFEMErrorCode MoFEMSNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *ctx)
Sens monitor printing residual field by field.
Definition SnesCtx.cpp:596
auto createVectorMPI(MPI_Comm comm, PetscInt n, PetscInt N)
Create MPI Vector.
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
auto getDMSnesCtx(DM dm)
Get SNES context data structure used by DM.
Definition DMMoFEM.hpp:1262
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
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
Section manager is used to create indexes and sections.
Definition ISManager.hpp:23
Post post-proc data at points from hash maps.
intrusive_ptr for managing petsc objects

Member Data Documentation

◆ dispFieldPtr

boost::shared_ptr<MatrixDouble> Seepage::dispFieldPtr
private

◆ dispGradPtr

boost::shared_ptr<MatrixDouble> Seepage::dispGradPtr
private

◆ doEvalField

PetscBool Seepage::doEvalField
private

◆ fieldEvalCoords

std::array<double, 3> Seepage::fieldEvalCoords {0.0, 0.0, 0.0}
private
Examples
mofem/tutorials/adv-5_poroelasticity/seepage.cpp.

Definition at line 210 of file seepage.cpp.

210{0.0, 0.0, 0.0};

◆ fieldEvalData

boost::shared_ptr<FieldEvaluatorInterface::SetPtsData> Seepage::fieldEvalData
private

◆ fluxFieldPtr

boost::shared_ptr<MatrixDouble> Seepage::fluxFieldPtr
private

◆ mDPtr

boost::shared_ptr<MatrixDouble> Seepage::mDPtr
private

◆ mField

MoFEM::Interface& Seepage::mField
private

◆ pressureFieldPtr

boost::shared_ptr<VectorDouble> Seepage::pressureFieldPtr
private

◆ strainFieldPtr

boost::shared_ptr<MatrixDouble> Seepage::strainFieldPtr
private

◆ stressFieldPtr

boost::shared_ptr<MatrixDouble> Seepage::stressFieldPtr
private

◆ uXScatter

std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter> > Seepage::uXScatter
private

◆ uYScatter

std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter> > Seepage::uYScatter
private

◆ uZScatter

std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter> > Seepage::uZScatter
private

The documentation for this struct was generated from the following file: