G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsVisitorMass.h
1 
14 #pragma once
15 
16 namespace gismo
17 {
18 
28 template <class T>
30 {
31 public:
32 
35  { }
36 
40  gsVisitorMass(const gsPde<T> & pde)
41  { GISMO_UNUSED(pde); }
42 
44  void initialize(const gsBasis<T> & basis,
45  const index_t ,
46  const gsOptionList & options,
47  gsQuadRule<T> & rule)
48  {
49  // Setup Quadrature (harmless slicing occurs)
50  rule = gsQuadrature::get(basis, options); // harmless slicing occurs here
51 
52  // Set Geometry evaluation flags
53  md.flags = NEED_MEASURE;
54  }
55 
57  inline void evaluate(const gsBasis<T> & basis, // to do: more unknowns
58  const gsGeometry<T> & geo,
59  // todo: add element here for efficiency
60  gsMatrix<T> & quNodes)
61  {
62  md.points = quNodes;
63  // Compute the active basis functions
64  // Assumes actives are the same for all quadrature points on the current element
65  basis.active_into(md.points.col(0), actives);
66  const index_t numActive = actives.rows();
67 
68  // Evaluate basis functions on element
69  basis.eval_into(md.points, basisData);
70 
71  // Compute geometry related values
72  geo.computeMap(md);
73 
74  // Initialize local matrix/rhs
75  localMat.setZero(numActive, numActive);
76  }
77 
79  inline void assemble(gsDomainIterator<T> & ,
80  gsVector<T> const & quWeights)
81  {
82  localMat.noalias() =
83  basisData * quWeights.asDiagonal() *
84  md.measures.asDiagonal() * basisData.transpose();
85  }
86 
88  inline void localToGlobal(const index_t patchIndex,
89  const std::vector<gsMatrix<T> > & eliminatedDofs,
90  gsSparseSystem<T> & system)
91  {
92  // Map patch-local DoFs to global DoFs
93  system.mapColIndices(actives, patchIndex, actives);
94 
95  // Add contributions to the system matrix
96  system.pushToMatrix(localMat, actives, eliminatedDofs.front(), 0, 0);
97  }
98 
99 /* ----------------------- to be removed later*/
100 
101  void initialize(const gsBasis<T> & basis,
102  gsQuadRule<T> & rule)
103  {
104  gsVector<short_t> numQuadNodes( basis.dim() );
105  for (short_t i = 0; i < basis.dim(); ++i)
106  numQuadNodes[i] = basis.degree(i) + 1;
107 
108  // Setup Quadrature
109  rule = gsGaussRule<T>(numQuadNodes);// harmless slicing occurs here
110 
111  // Set Geometry evaluation flags
112  md.flags = NEED_MEASURE;
113  }
114 
115 
116  void localToGlobal(const gsDofMapper & mapper,
117  const gsMatrix<T> & eliminatedDofs,
118  const index_t patchIndex,
119  gsSparseMatrix<T> & sysMatrix,
120  gsMatrix<T> & rhsMatrix )
121  {
122  mapper.localToGlobal(actives, patchIndex, actives);
123 
124  const index_t numActive = actives.rows();
125 
126  for (index_t i = 0; i < numActive; ++i)
127  {
128  const int ii = actives(i,0); // N_i
129 
130  if ( mapper.is_free_index(ii) )
131  {
132  for (index_t j = 0; j < numActive; ++j)
133  {
134  const int jj = actives(j,0); // N_j
135  if ( mapper.is_free_index(jj) )
136  //if ( jj <= ii ) // only store lower triangular part
137  sysMatrix.coeffRef(ii, jj) += localMat(i, j); // N_i*N_j
138  }
139  }
140  }
141  }
142 
143 
144 protected:
145 
146  // Basis values
147  gsMatrix<T> basisData;
148  gsMatrix<index_t> actives;
149 
150  // Local matrix
151  gsMatrix<T> localMat;
152 
153  gsMapData<T> md;
154 };
155 
156 
157 } // namespace gismo
Abstract base class representing a geometry map.
Definition: gsGeometry.h:92
Class representing a reference quadrature rule.
Definition: gsQuadRule.h:28
void mapColIndices(const gsMatrix< index_t > &actives, const index_t patchIndex, gsMatrix< index_t > &result, const index_t c=0) const
mapColIndices Maps a set of basis indices by the corresponding dofMapper.
Definition: gsSparseSystem.h:584
void pushToMatrix(const gsMatrix< T > &localMat, const gsMatrix< index_t > &actives, const gsMatrix< T > &eliminatedDofs, const size_t r=0, const size_t c=0)
pushToMatrix pushes the local matrix for an element to the global system,
Definition: gsSparseSystem.h:638
void assemble(gsDomainIterator< T > &, gsVector< T > const &quWeights)
Assemble on element.
Definition: gsVisitorMass.h:79
virtual void computeMap(gsMapData< T > &InOut) const
Computes map function data.
Definition: gsFunction.hpp:822
#define short_t
Definition: gsConfig.h:35
void evaluate(const gsBasis< T > &basis, const gsGeometry< T > &geo, gsMatrix< T > &quNodes)
Evaluate on element.
Definition: gsVisitorMass.h:57
virtual short_t degree(short_t i) const
Degree with respect to the i-th variable. If the basis is a tensor product of (piecewise) polynomial ...
Definition: gsBasis.hpp:650
void initialize(const gsBasis< T > &basis, const index_t, const gsOptionList &options, gsQuadRule< T > &rule)
Initialize.
Definition: gsVisitorMass.h:44
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
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
The density of the measure pull back.
Definition: gsForwardDeclarations.h:76
virtual void eval_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluates nonzero basis functions at point u into result.
Definition: gsBasis.hpp:443
gsVisitorMass(const gsPde< T > &pde)
Constructor.
Definition: gsVisitorMass.h:40
#define GISMO_UNUSED(x)
Definition: gsDebug.h:112
Class which holds a list of parameters/options, and provides easy access to them. ...
Definition: gsOptionList.h:32
void localToGlobal(const index_t patchIndex, const std::vector< gsMatrix< T > > &eliminatedDofs, gsSparseSystem< T > &system)
Adds the contributions to the sparse system.
Definition: gsVisitorMass.h:88
Class that represents the (tensor) Gauss-Legendre quadrature rule.
Definition: gsGaussRule.h:27
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
gsVisitorMass()
Constructor.
Definition: gsVisitorMass.h:34