G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsVisitorGradGrad.h
Go to the documentation of this file.
1 
15 #pragma once
16 
17 
18 #include <gsAssembler/gsVisitorMass.h>
19 namespace gismo
20 {
21 
32 template <class T>
33 class gsVisitorGradGrad : public gsVisitorMass<T> // inherit to reuse functionality
34 {
35 public:
36  typedef gsVisitorMass<T> Base;
37 
38 public:
39 
44  { GISMO_UNUSED(pde); }
45 
47  void initialize(const gsBasis<T> & basis,
48  const index_t /*patchIndex*/,
49  const gsOptionList & options,
50  gsQuadRule<T> & rule)
51  {
52  // Setup Quadrature
53  rule = gsQuadrature::get(basis, options); // harmless slicing occurs here
54 
55  // Set Geometry evaluation flags
57  }
58 
60  inline void evaluate(const gsBasis<T> & basis,
61  const gsGeometry<T> & geo,
62  // todo: add element here for efficiency
63  gsMatrix<T> & quNodes)
64  {
65  md.points = quNodes;
66  // Compute the active basis functions
67  // Assumes actives are the same for all quadrature points on the current element
68  basis.active_into(md.points.col(0), actives);
69  const index_t numActive = actives.rows();
70 
71  // Evaluate basis functions on element
72  basis.deriv_into(md.points, basisData);
73 
74  // Compute geometry related values
75  geo.computeMap(md);
76 
77  // Initialize local matrix
78  localMat.setZero(numActive, numActive);
79  }
80 
82  inline void assemble(gsDomainIterator<T> & /*element*/,
83  gsVector<T> const & quWeights)
84  {
85  for (index_t k = 0; k < quWeights.rows(); ++k) // loop over quadrature nodes
86  {
87  // Multiply quadrature weight by the geometry measure
88  const T weight = quWeights[k] * md.measure(k);
89 
90  // Compute physical gradients at k as a Dim x NumActive matrix
91  transformGradients(md, k, basisData, basisPhGrads);
92 
93  localMat.noalias() += weight * ( basisPhGrads.transpose() * basisPhGrads );
94  }
95  }
96 
97  //Inherited from gsVisitorMass
98  //void localToGlobal( ... )
99 
100 private:
101 
102  // Gradient values
103  gsMatrix<T> basisPhGrads;
104  using Base:: basisData;
105  using Base::actives;
106 
107  // Local matrix
108  using Base::localMat;
109 
110  using Base::md;
111 };
112 
113 
114 } // namespace gismo
Abstract base class representing a geometry map.
Definition: gsGeometry.h:92
Class representing a reference quadrature rule.
Definition: gsQuadRule.h:28
The visitor computes element grad-grad integrals.
Definition: gsVisitorGradGrad.h:33
virtual void computeMap(gsMapData< T > &InOut) const
Computes map function data.
Definition: gsFunction.hpp:822
Abstract class representing a PDE (partial differential equation).
Definition: gsPde.h:43
#define index_t
Definition: gsConfig.h:32
The visitor computes element mass integrals.
Definition: gsVisitorMass.h:29
gsVisitorGradGrad(const gsPde< T > &pde)
Definition: gsVisitorGradGrad.h:43
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
Class which enables iteration over all elements of a parameter domain.
Definition: gsDomainIterator.h:67
void evaluate(const gsBasis< T > &basis, const gsGeometry< T > &geo, gsMatrix< T > &quNodes)
Evaluate on element.
Definition: gsVisitorGradGrad.h:60
void initialize(const gsBasis< T > &basis, const index_t, const gsOptionList &options, gsQuadRule< T > &rule)
Initialize.
Definition: gsVisitorGradGrad.h:47
The density of the measure pull back.
Definition: gsForwardDeclarations.h:76
#define GISMO_UNUSED(x)
Definition: gsDebug.h:112
Gradient transformation matrix.
Definition: gsForwardDeclarations.h:77
Class which holds a list of parameters/options, and provides easy access to them. ...
Definition: gsOptionList.h:32
virtual void deriv_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluates the first partial derivatives of the nonzero basis function.
Definition: gsBasis.hpp:451
void assemble(gsDomainIterator< T > &, gsVector< T > const &quWeights)
Assemble on element.
Definition: gsVisitorGradGrad.h:82
A basis represents a family of scalar basis functions defined over a common parameter domain...
Definition: gsBasis.h:78
virtual void active_into(const gsMatrix< T > &u, gsMatrix< index_t > &result) const
Returns the indices of active basis functions at points u, as a list of indices, in result...
Definition: gsBasis.hpp:293