G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsVisitorElasticityNeumann.h
Go to the documentation of this file.
1
16#pragma once
17
19#include <gsCore/gsFuncData.h>
20
21namespace gismo
22{
23
24template <class T>
25class gsVisitorElasticityNeumann
26{
27public:
28
29 gsVisitorElasticityNeumann(const gsPde<T> & /* pde_ */,
30 const boundary_condition<T> & s)
31 : neumannFunction_ptr(s.function().get()),
32 patchSide(s.side()) {}
33
34 void initialize(const gsBasisRefs<T> & basisRefs,
35 const index_t patchIndex,
36 const gsOptionList & options,
37 gsQuadRule<T> & rule)
38 {
39 GISMO_UNUSED(patchIndex);
40 // parametric dimension of the first displacement component
41 dim = basisRefs.front().dim();
42 // a quadrature rule is defined by the basis for the first displacement component.
43 rule = gsQuadrature::get(basisRefs.front(), options,patchSide.direction());
44 // saving necessary info
45 forceScaling = options.getReal("ForceScaling");
46 // resize containers for global indices
47 globalIndices.resize(dim);
48 blockNumbers.resize(dim);
49 }
50
51 inline void evaluate(const gsBasisRefs<T> & basisRefs,
52 const gsGeometry<T> & geo,
53 const gsMatrix<T> & quNodes)
54 {
55 // store quadrature points of the element for geometry evaluation
56 md.points = quNodes;
57 // NEED_VALUE to get points in the physical domain for evaluation of the load
58 // NEED_MEASURE to get the Jacobian determinant values for integration
59 md.flags = NEED_VALUE | NEED_MEASURE;
60 // Compute image of the quadrature points plus gradient, jacobian and other necessary data
61 geo.computeMap(md);
62 // Evaluate the Neumann functon on the images of the quadrature points
63 neumannFunction_ptr->eval_into(md.values[0], neumannValues);
64 // find local indices of the displacement basis functions active on the element
65 basisRefs.front().active_into(quNodes.col(0),localIndicesDisp);
66 N_D = localIndicesDisp.rows();
67 // Evaluate basis functions on element
68 basisRefs.front().eval_into(quNodes,basisValuesDisp);
69 }
70
71 inline void assemble(gsDomainIterator<T> & element,
72 const gsVector<T> & quWeights)
73 {
74 GISMO_UNUSED(element);
75 // Initialize local matrix/rhs
76 localRhs.setZero(dim*N_D,1);
77 // loop over the quadrature nodes
78 for (index_t q = 0; q < quWeights.rows(); ++q)
79 {
80 // Compute the outer normal vector on the side
81 // normal length equals to the local area measure
82 outerNormal(md, q, patchSide, unormal);
83 // Collect the factors here: quadrature weight and geometry measure
84 const T weight = quWeights[q] * unormal.norm() * forceScaling;
85
86 for (short_t d = 0; d < dim; ++d)
87 localRhs.middleRows(d*N_D,N_D).noalias() += weight * neumannValues(d,q) * basisValuesDisp.col(q);
88 }
89 }
90
91 inline void localToGlobal(const int patchIndex,
92 const std::vector<gsMatrix<T> > & eliminatedDofs,
93 gsSparseSystem<T> & system)
94 {
95 GISMO_UNUSED(eliminatedDofs);
96 // computes global indices for displacement components
97 for (short_t d = 0; d < dim; ++d)
98 {
99 system.mapColIndices(localIndicesDisp, patchIndex, globalIndices[d], d);
100 blockNumbers.at(d) = d;
101 }
102 // push to global system
103 system.pushToRhs(localRhs,globalIndices,blockNumbers);
104 }
105
106protected:
107 // problem info
108 short_t dim;
109 const gsFunctionSet<T> * neumannFunction_ptr;
110 T forceScaling;
111 boxSide patchSide;
112 // geometry mapping
113 gsMapData<T> md;
114 // local components of the global linear system
115 gsMatrix<T> localRhs;
116 // local indices (at the current patch) of the displacement basis functions active at the current element
117 gsMatrix<index_t> localIndicesDisp;
118 // number of displacement basis functions active at the current element
119 index_t N_D;
120 // values of displacement basis functions at quadrature points at the current element stored as a N_D x numQuadPoints matrix;
121 gsMatrix<T> basisValuesDisp;
122 // values of the boundary loading function stored as a dim x numQuadPoints matrix;
123 gsMatrix<T> neumannValues;
124
125 // all temporary matrices defined here for efficiency
126 gsVector<T> unormal;
127 // containers for global indices
128 std::vector< gsMatrix<index_t> > globalIndices;
129 gsVector<index_t> blockNumbers;
130};
131
132} // namespace gismo
#define short_t
Definition gsConfig.h:35
#define index_t
Definition gsConfig.h:32
#define GISMO_UNUSED(x)
Definition gsDebug.h:112
This object is a cache for computed values from an evaluator.
Creates a variety of quadrature rules.
The G+Smo namespace, containing all definitions for the library.
@ NEED_VALUE
Value of the object.
Definition gsForwardDeclarations.h:72
@ NEED_MEASURE
The density of the measure pull back.
Definition gsForwardDeclarations.h:76
static gsQuadRule< T > get(const gsBasis< T > &basis, const gsOptionList &options, short_t fixDir=-1)
Constructs a quadrature rule based on input options.
Definition gsQuadrature.h:51