v0.14.0
NaturalBoundaryBC.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 template <int BASE_DIM, int FIELD_DIM, AssemblyType A, IntegrationType I,
12  typename OpBase>
13 struct AddFluxToRhsPipelineImpl<
14 
15  OpFluxRhsImpl<BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase>, A, I, OpBase
16 
17  > {
18 
19  AddFluxToRhsPipelineImpl() = delete;
20 
21  using T =
22  typename NaturalBC<OpBase>::template Assembly<A>::template LinearForm<I>;
23 
24  using OpForce =
26  using OpSpringRhs =
27  typename T::template OpFlux<ElasticExample::SpringBcType<BLOCKSET>, 1,
29  using OpFluidLevelRhs =
32 
34 
35  boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
36  MoFEM::Interface &m_field, std::string field_name, double scale, Sev sev
37 
38  ) {
40  CHKERR T::template AddFluxToPipeline<OpForce>::add(
41  pipeline, m_field, field_name, {}, "FORCE", sev);
42  CHKERR T::template AddFluxToPipeline<OpFluidLevelRhs>::add(
43  pipeline, m_field, field_name, {}, scale, "FLUID_PRESSURE", sev);
44  auto u_ptr = boost::make_shared<MatrixDouble>();
45  pipeline.push_back(
46  new OpCalculateVectorFieldValues<SPACE_DIM>(field_name, u_ptr));
47  CHKERR T::template AddFluxToPipeline<OpSpringRhs>::add(
48  pipeline, m_field, field_name, u_ptr, scale, "SPRING", sev);
49 
51  }
52 };
53 
54 template <int BASE_DIM, int FIELD_DIM, AssemblyType A, IntegrationType I,
55  typename OpBase>
56 struct AddFluxToLhsPipelineImpl<
57 
58  OpFluxLhsImpl<BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase>, A, I, OpBase
59 
60  > {
61 
62  AddFluxToLhsPipelineImpl() = delete;
63 
64  using T = typename NaturalBC<OpBase>::template Assembly<
65  A>::template BiLinearForm<I>;
66 
67  using OpSpringLhs =
70 
72 
73  boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
74  MoFEM::Interface &m_field, std::string field_name, Sev sev
75 
76  ) {
78  CHKERR T::template AddFluxToPipeline<OpSpringLhs>::add(
79  pipeline, m_field, field_name, field_name, "SPRING", sev);
81  }
82 };
AddFluxToRhsPipelineImpl< OpFluxRhsImpl< BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase >, A, I, OpBase >::OpFluidLevelRhs
typename T::template OpFlux< ElasticExample::FluidLevelType< BLOCKSET >, 1, SPACE_DIM > OpFluidLevelRhs
Definition: NaturalBoundaryBC.hpp:31
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
scale
double scale
Definition: plastic.cpp:119
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
AddFluxToRhsPipelineImpl< OpFluxRhsImpl< 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, double scale, Sev sev)
Definition: NaturalBoundaryBC.hpp:33
AddFluxToRhsPipelineImpl< OpFluxRhsImpl< BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase >, A, I, OpBase >::T
typename NaturalBC< OpBase >::template Assembly< A >::template LinearForm< I > T
Definition: NaturalBoundaryBC.hpp:22
AddFluxToLhsPipelineImpl< OpFluxLhsImpl< 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: NaturalBoundaryBC.hpp:69
OpFlux
Definition: hcurl_divergence_operator_2d.cpp:56
AddFluxToRhsPipelineImpl< OpFluxRhsImpl< BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase >, A, I, OpBase >::OpForce
typename T::template OpFlux< NaturalForceMeshsets, 1, SPACE_DIM > OpForce
Definition: NaturalBoundaryBC.hpp:25
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< BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase >, A, I, OpBase >::OpSpringRhs
typename T::template OpFlux< ElasticExample::SpringBcType< BLOCKSET >, 1, SPACE_DIM > OpSpringRhs
Definition: NaturalBoundaryBC.hpp:28
AddFluxToLhsPipelineImpl< OpFluxLhsImpl< BoundaryBCs, BASE_DIM, FIELD_DIM, A, I, OpBase >, A, I, OpBase >::T
typename NaturalBC< OpBase >::template Assembly< A >::template BiLinearForm< I > T
Definition: NaturalBoundaryBC.hpp:65
AddFluxToLhsPipelineImpl< OpFluxLhsImpl< 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: NaturalBoundaryBC.hpp:71
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