v0.14.0
ContactNaturalBC.hpp
Go to the documentation of this file.
1 /**
2  * @file NaturalBoundaryBC.hpp
3  * @brief Implementation of natural boundary conditions
4  * @version 0.13.2
5  * @date 2022-09-22
6  *
7  * @copyright Copyright (c) 2022
8  *
9  */
10 
11 namespace ContactOps {
12 struct BoundaryBCs;
13 struct DomainBCs;
14 } // namespace ContactOps
15 
16 #include <ElasticSpring.hpp>
17 #include <FluidLevel.hpp>
18 
19 template <int BASE_DIM, int FIELD_DIM, AssemblyType A, IntegrationType I,
20  typename OpBase>
21 struct AddFluxToRhsPipelineImpl<
22 
23  OpFluxRhsImpl<ContactOps::BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase>,
24  A, I, OpBase
25 
26  > {
27 
28  AddFluxToRhsPipelineImpl() = delete;
29 
30  using T =
31  typename NaturalBC<OpBase>::template Assembly<A>::template LinearForm<I>;
32 
33  using OpForce =
35  using OpSpringRhs =
36  typename T::template OpFlux<ElasticExample::SpringBcType<BLOCKSET>, 1,
38  using OpFluidLevelRhs =
41 
43 
44  boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
45  MoFEM::Interface &m_field, std::string field_name,
46  std::vector<boost::shared_ptr<ScalingMethod>> smv, Sev sev
47 
48  ) {
50  CHKERR T::template AddFluxToPipeline<OpForce>::add(
51  pipeline, m_field, field_name, smv, "FORCE", sev);
52  CHKERR T::template AddFluxToPipeline<OpFluidLevelRhs>::add(
53  pipeline, m_field, field_name, smv, 1, "FLUID_PRESSURE", sev);
54  auto u_ptr = boost::make_shared<MatrixDouble>();
55  pipeline.push_back(
56  new OpCalculateVectorFieldValues<SPACE_DIM>(field_name, u_ptr));
57  CHKERR T::template AddFluxToPipeline<OpSpringRhs>::add(
58  pipeline, m_field, field_name, u_ptr, 1, "SPRING", sev);
60  }
61 };
62 
63 template <int BASE_DIM, int FIELD_DIM, AssemblyType A, IntegrationType I,
64  typename OpBase>
65 struct AddFluxToLhsPipelineImpl<
66 
67  OpFluxLhsImpl<ContactOps::BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase>,
68  A, I, OpBase
69 
70  > {
71 
72  AddFluxToLhsPipelineImpl() = delete;
73 
74  using T = typename NaturalBC<OpBase>::template Assembly<
75  A>::template BiLinearForm<I>;
76 
77  using OpSpringLhs =
80 
82 
83  boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
84  MoFEM::Interface &m_field, std::string field_name, Sev sev
85 
86  ) {
88  CHKERR T::template AddFluxToPipeline<OpSpringLhs>::add(
89  pipeline, m_field, field_name, field_name, "SPRING", sev);
91  }
92 };
93 
94 template <int BASE_DIM, int FIELD_DIM, AssemblyType A, IntegrationType I,
95  typename OpBase>
96 struct AddFluxToRhsPipelineImpl<
97 
98  OpFluxRhsImpl<ContactOps::DomainBCs, BASE_DIM, FIELD_DIM, A, I, OpBase>, A,
99  I, OpBase
100 
101  > {
102 
103  AddFluxToRhsPipelineImpl() = delete;
104 
105  using T =
106  typename NaturalBC<OpBase>::template Assembly<A>::template LinearForm<I>;
107  using OpBodyForce = typename T::template OpFlux<NaturalMeshsetType<BLOCKSET>,
109 
111 
112  boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
113  MoFEM::Interface &m_field, std::string field_name,
114  std::vector<boost::shared_ptr<ScalingMethod>> smv, Sev sev
115 
116  ) {
118  CHKERR T::template AddFluxToPipeline<OpBodyForce>::add(
119  pipeline, m_field, field_name, smv, "BODY_FORCE", sev);
121  }
122 };
ElasticSpring.hpp
Implementation of elastic spring bc.
AddFluxToRhsPipelineImpl< OpFluxRhsImpl< ContactOps::BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase >, A, I, OpBase >::T
typename NaturalBC< OpBase >::template Assembly< A >::template LinearForm< I > T
Definition: ContactNaturalBC.hpp:31
AddFluxToRhsPipelineImpl< OpFluxRhsImpl< ContactOps::DomainBCs, BASE_DIM, FIELD_DIM, A, I, OpBase >, A, I, OpBase >::OpBodyForce
typename T::template OpFlux< NaturalMeshsetType< BLOCKSET >, BASE_DIM, FIELD_DIM > OpBodyForce
Definition: ContactNaturalBC.hpp:108
AddFluxToRhsPipelineImpl< OpFluxRhsImpl< ContactOps::DomainBCs, BASE_DIM, FIELD_DIM, A, I, OpBase >, A, I, OpBase >::T
typename NaturalBC< OpBase >::template Assembly< A >::template LinearForm< I > T
Definition: ContactNaturalBC.hpp:106
ContactOps
Definition: contact.cpp:99
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
AddFluxToLhsPipelineImpl< OpFluxLhsImpl< ContactOps::BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase >, A, I, OpBase >::OpSpringLhs
typename T::template OpFlux< ElasticExample::SpringBcType< BLOCKSET >, BASE_DIM, FIELD_DIM > OpSpringLhs
Definition: ContactNaturalBC.hpp:79
AddFluxToRhsPipelineImpl< OpFluxRhsImpl< ContactOps::BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase >, A, I, OpBase >::OpSpringRhs
typename T::template OpFlux< ElasticExample::SpringBcType< BLOCKSET >, 1, SPACE_DIM > OpSpringRhs
Definition: ContactNaturalBC.hpp:37
BASE_DIM
constexpr int BASE_DIM
Definition: dg_projection.cpp:15
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
FIELD_DIM
constexpr int FIELD_DIM
Definition: child_and_parent.cpp:15
I
constexpr IntegrationType I
Definition: operators_tests.cpp:31
AddFluxToLhsPipelineImpl< OpFluxLhsImpl< ContactOps::BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase >, A, I, OpBase >::T
typename NaturalBC< OpBase >::template Assembly< A >::template BiLinearForm< I > T
Definition: ContactNaturalBC.hpp:75
MoFEM::Sev
MoFEM::LogManager::SeverityLevel Sev
Definition: CoreTemplates.hpp:17
MoFEM::OpBaseImpl
Definition: FormsIntegrators.hpp:178
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
SPACE_DIM
constexpr int SPACE_DIM
Definition: child_and_parent.cpp:16
OpFlux
Definition: hcurl_divergence_operator_2d.cpp:56
AddFluxToRhsPipelineImpl< OpFluxRhsImpl< ContactOps::BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase >, A, I, OpBase >::OpForce
typename T::template OpFlux< NaturalForceMeshsets, 1, SPACE_DIM > OpForce
Definition: ContactNaturalBC.hpp:34
AddFluxToRhsPipelineImpl< OpFluxRhsImpl< ContactOps::DomainBCs, BASE_DIM, FIELD_DIM, A, I, OpBase >, A, I, OpBase >::add
static MoFEMErrorCode add(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, MoFEM::Interface &m_field, std::string field_name, std::vector< boost::shared_ptr< ScalingMethod >> smv, Sev sev)
Definition: ContactNaturalBC.hpp:110
AddFluxToLhsPipelineImpl< OpFluxLhsImpl< ContactOps::BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase >, A, I, OpBase >::add
static MoFEMErrorCode add(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, MoFEM::Interface &m_field, std::string field_name, Sev sev)
Definition: ContactNaturalBC.hpp:81
BiLinearForm
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
MoFEM::IntegrationType
IntegrationType
Form integrator integration types.
Definition: FormsIntegrators.hpp:136
MoFEM::AssemblyType
AssemblyType
[Storage and set boundary conditions]
Definition: FormsIntegrators.hpp:104
AddFluxToRhsPipelineImpl< OpFluxRhsImpl< ContactOps::BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase >, A, I, OpBase >::OpFluidLevelRhs
typename T::template OpFlux< ElasticExample::FluidLevelType< BLOCKSET >, 1, SPACE_DIM > OpFluidLevelRhs
Definition: ContactNaturalBC.hpp:40
FluidLevel.hpp
Natural boundary condition applying pressure from fluid.
DomainBCs
[OpInternalForce]
Definition: elastic.cpp:45
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
BoundaryBCs
Definition: elastic.cpp:46
AddFluxToRhsPipelineImpl< OpFluxRhsImpl< ContactOps::BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase >, A, I, OpBase >::add
static MoFEMErrorCode add(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, MoFEM::Interface &m_field, std::string field_name, std::vector< boost::shared_ptr< ScalingMethod >> smv, Sev sev)
Definition: ContactNaturalBC.hpp:42