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