v0.13.2
Loading...
Searching...
No Matches
ContactNaturalBoundaryBC.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#include <ElasticSpring.hpp>
12
14 typename OpBase>
15struct AddFluxToRhsPipelineImpl<
16
17 OpFluxRhsImpl<ContactOps::BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase>,
18 A, I, OpBase
19
20 > {
21
23
24 using T =
25 typename NaturalBC<OpBase>::template Assembly<A>::template LinearForm<I>;
26
27 using OpForce =
31 SPACE_DIM>;
32
33 static MoFEMErrorCode add(
34
35 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
36 MoFEM::Interface &m_field, std::string field_name,
37 std::vector<boost::shared_ptr<ScalingMethod>> smv, Sev sev
38
39 ) {
41 CHKERR T::template AddFluxToPipeline<OpForce>::add(
42 pipeline, m_field, field_name, smv, "FORCE", sev);
43 auto u_ptr = boost::make_shared<MatrixDouble>();
44 pipeline.push_back(
45 new OpCalculateVectorFieldValues<SPACE_DIM>(field_name, u_ptr));
46 CHKERR T::template AddFluxToPipeline<OpSpringRhs>::add(
47 pipeline, m_field, field_name, u_ptr, 1, "SPRING", sev);
49 }
50};
51
53 typename OpBase>
54struct AddFluxToLhsPipelineImpl<
55
56 OpFluxLhsImpl<ContactOps::BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase>,
57 A, I, OpBase
58
59 > {
60
62
63 using T = typename NaturalBC<OpBase>::template Assembly<
64 A>::template BiLinearForm<I>;
65
69
70 static MoFEMErrorCode add(
71
72 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
73 MoFEM::Interface &m_field, std::string field_name, Sev sev
74
75 ) {
77 CHKERR T::template AddFluxToPipeline<OpSpringLhs>::add(
78 pipeline, m_field, field_name, field_name, "SPRING", sev);
80 }
81};
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
[Define dimension]
Definition: elastic.cpp:16
IntegrationType
Form integrator integration types.
AssemblyType
[Storage and set boundary conditions]
constexpr int BASE_DIM
[Only used for dynamics]
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)
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)
Deprecated interface functions.