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