v0.14.0
Loading...
Searching...
No Matches
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
13namespace PlasticOps {
14
15struct DomainBCs;
16struct BoundaryBCs;
17
18} // namespace PlasticOps
19
20
21template <int BASE_DIM, int FIELD_DIM, AssemblyType A, IntegrationType I,
22 typename OpBase>
23struct AddFluxToRhsPipelineImpl<
24
25 OpFluxRhsImpl<PlasticOps::DomainBCs, BASE_DIM, FIELD_DIM, A, I, OpBase>, A,
26 I, OpBase
27
28 > {
29
31
32 using T =
33 typename NaturalBC<OpBase>::template Assembly<A>::template LinearForm<I>;
36
37 static MoFEMErrorCode add(
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
52 typename OpBase>
53struct AddFluxToRhsPipelineImpl<
54
55 OpFluxRhsImpl<PlasticOps::BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase>,
56 A, I, OpBase
57
58 > {
59
61
62 using T =
63 typename NaturalBC<OpBase>::template Assembly<A>::template LinearForm<I>;
64
65 using OpForce =
69 SPACE_DIM>;
70
71 static MoFEMErrorCode add(
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
91 typename OpBase>
92struct AddFluxToLhsPipelineImpl<
93
94 OpFluxLhsImpl<PlasticOps::BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase>,
95 A, I, OpBase
96
97 > {
98
100
101 using T = typename NaturalBC<OpBase>::template Assembly<
102 A>::template BiLinearForm<I>;
103
107
108 static MoFEMErrorCode add(
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};
Implementation of elastic spring bc.
constexpr int SPACE_DIM
constexpr int FIELD_DIM
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
constexpr int BASE_DIM
constexpr int BASE_DIM
Definition: elastic.cpp:15
IntegrationType
Form integrator integration types.
AssemblyType
[Storage and set boundary conditions]
constexpr IntegrationType I
constexpr AssemblyType A
constexpr auto field_name
typename T::template OpFlux< ElasticExample::SpringBcType< BLOCKSET >, BASE_DIM, FIELD_DIM > OpSpringLhs
static MoFEMErrorCode add(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, MoFEM::Interface &m_field, std::string field_name, Sev sev)
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)
typename T::template OpFlux< ElasticExample::SpringBcType< BLOCKSET >, 1, SPACE_DIM > OpSpringRhs
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)
typename T::template OpFlux< NaturalMeshsetType< BLOCKSET >, BASE_DIM, FIELD_DIM > OpBodyForce
[OpInternalForce]
Definition: elastic.cpp:45
Deprecated interface functions.