G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsVisitorNitsche.h
Go to the documentation of this file.
1 
14 #pragma once
15 
17 
18 namespace gismo
19 {
20 
53 template <class T>
55 {
56 public:
57 
63  : m_pde(&pde), m_dirdata_ptr( bc.function().get() ), m_penalty(-1), m_side(bc.ps)
64  { }
65 
68  {
69  gsOptionList options;
70  options.addReal ("Nitsche.Alpha",
71  "Parameter alpha for Nitsche scheme.", 1);
72  options.addReal ("Nitsche.Beta",
73  "Parameter beta for Nitsche scheme.", 1);
74  options.addReal ("Nitsche.Penalty",
75  "Penalty parameter penalty for Nitsche scheme; if negative, default 2.5(p+d)(p+1) is used.",-1);
76  options.addSwitch("Nitsche.ParameterGridSize",
77  "Use grid size on parameter domain for the penalty term.", false);
78  return options;
79  }
80 
82  void initialize(const gsBasis<T> & basis,
83  const index_t,
84  const gsOptionList & options,
85  gsQuadRule<T> & rule)
86  {
87  // Setup Quadrature (harmless slicing occurs)
88  rule = gsQuadrature::get(basis, options, m_side.direction());
89 
90  m_penalty = options.askReal("Nitsche.Penalty",-1);
91  // If not given, use default
92  if (m_penalty<0)
93  {
94  const index_t deg = basis.maxDegree();
95  m_penalty = (T)(2.5) * (T)(deg + basis.dim()) * (T)(deg + 1);
96  }
97 
98  m_alpha = options.askReal("Nitsche.Alpha", 1);
99  m_beta = options.askReal("Nitsche.Beta" , 1);
100 
101  if (options.getSwitch("Nitsche.ParameterGridSize"))
102  {
103  m_h = basis.getMinCellLength();
104  }
105  else
106  {
107  GISMO_ENSURE (m_pde, "gsVisitorNitsche::initialize: No PDE given.");
109  basis,
110  m_pde->patches()[m_side.patch],
111  m_side
112  );
113  }
114 
115  // Set Geometry evaluation flags
117 
118  }
119 
121  inline void evaluate(const gsBasis<T> & basis,
122  const gsGeometry<T> & geo,
123  // todo: add element here for efficiency
124  const gsMatrix<T> & quNodes)
125  {
126  md.points = quNodes;
127  // Compute the active basis functions
128  // Assumes actives are the same for all quadrature points on the current element
129  basis.active_into(md.points.col(0) , actives);
130  const index_t numActive = actives.rows();
131 
132  // Evaluate basis values and derivatives on element
133  basis.evalAllDers_into( md.points, 1, basisData, true);
134 
135  // Compute geometry related values
136  geo.computeMap(md);
137 
138  // Evaluate the Dirichlet data
139  m_dirdata_ptr->eval_into(md.values[0], dirData);
140 
141  // Initialize local matrix/rhs
142  localMat.setZero(numActive, numActive);
143  localRhs.setZero(numActive, m_dirdata_ptr->targetDim() );
144  }
145 
147  inline void assemble(gsDomainIterator<T> & /*element*/,
148  const gsVector<T> & quWeights)
149  {
150  gsMatrix<T> & bGrads = basisData[1];
151  const index_t numActive = actives.rows();
152 
153  for (index_t k = 0; k < quWeights.rows(); ++k) // loop over quadrature nodes
154  {
155 
156  const typename gsMatrix<T>::Block bVals =
157  basisData[0].block(0,k,numActive,1);
158 
159  // Compute the outer normal vector on the side
160  outerNormal(md, k, m_side, unormal);
161 
162  // Multiply quadrature weight by the geometry measure
163  const T weight = quWeights[k] * unormal.norm();
164 
165  // Compute the unit normal vector
166  unormal.normalize();
167 
168  // Compute physical gradients at k as a Dim x NumActive matrix
169  transformGradients(md, k, bGrads, pGrads);
170 
171  // Get penalty parameter
172  const T mu = m_penalty / m_h;
173 
174  // Sum up quadrature point evaluations
175  localRhs.noalias() -= weight * ( ( m_beta * pGrads.transpose() * unormal - mu * bVals ) * dirData.col(k).transpose() );
176 
177  localMat.noalias() -= weight * (
178  m_alpha * bVals * unormal.transpose() * pGrads
179  + m_beta * (bVals * unormal.transpose() * pGrads).transpose()
180  - mu * bVals * bVals.transpose()
181  );
182  }
183  }
184 
186  inline void localToGlobal(const index_t patchIndex,
187  const std::vector<gsMatrix<T> > & ,
188  gsSparseSystem<T> & system)
189  {
190  // Map patch-local DoFs to global DoFs
191  system.mapColIndices(actives, patchIndex, actives);
192 
193  // Add contributions to the system matrix and right-hand side
194  system.pushAllFree(localMat, localRhs, actives, 0);
195  }
196 
198  void localToGlobal(const gsDofMapper & mapper,
199  const gsMatrix<T> & eliminatedDofs,
200  const index_t patchIndex,
201  gsSparseMatrix<T> & sysMatrix,
202  gsMatrix<T> & rhsMatrix )
203  {
204  // Local DoFs to global DoFs
205  mapper.localToGlobal(actives, patchIndex, actives);
206  const index_t numActive = actives.rows();
207 
208  // Push element contribution to the global matrix and load vector
209  for (index_t j=0; j!=numActive; ++j)
210  {
211  const unsigned jj = actives(j);
212  rhsMatrix.row(jj) += localRhs.row(j);
213  for (index_t i=0; i!=numActive; ++i)
214  {
215  const unsigned ii = actives(i);
216 // if ( jj <= ii ) // assuming symmetric problem
217  sysMatrix( ii, jj ) += localMat(i,j);
218  }
219  }
220 
221  }
222 
223 private:
224 
226  const gsPde<T> * m_pde;
227 
230 
233 
236 
239 
242 
243 
244 private:
245  // Basis values
246  std::vector<gsMatrix<T> > basisData;
247  gsMatrix<T> pGrads;
248  gsMatrix<index_t> actives;
249 
250  // Normal and Neumann values
251  gsVector<T> unormal;
252  gsMatrix<T> dirData;
253 
254  // Local matrix and rhs
255  gsMatrix<T> localMat;
256  gsMatrix<T> localRhs;
257 
258  gsMapData<T> md;
259 
260  // Grid size
261  T m_h;
262 };
263 
264 
265 } // 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
const gsPde< T > * m_pde
The underlying PDE.
Definition: gsVisitorNitsche.h:226
Struct which represents a certain side of a patch.
Definition: gsBoundary.h:231
virtual void computeMap(gsMapData< T > &InOut) const
Computes map function data.
Definition: gsFunction.hpp:822
virtual void evalAllDers_into(const gsMatrix< T > &u, int n, std::vector< gsMatrix< T > > &result, bool sameElement=false) const
Evaluate the nonzero functions and their derivatives up to order n at points u into result...
Definition: gsFunctionSet.hpp:81
the gsMapData is a cache of pre-computed function (map) values.
Definition: gsFuncData.h:324
static T estimateSmallestPerpendicularCellSize(const gsBasis< T > &basis, const gsGeometry< T > &geo, patchSide side)
Estimates the gird size perpendicular to the given side on the physical domain, as required for SIPG ...
Definition: gsVisitorDg.h:232
Visitor for adding the interface conditions for the interior penalty methods of the Poisson problem...
Abstract class representing a PDE (partial differential equation).
Definition: gsPde.h:43
#define index_t
Definition: gsConfig.h:32
#define GISMO_ENSURE(cond, message)
Definition: gsDebug.h:102
void localToGlobal(const gsMatrix< index_t > &locals, index_t patchIndex, gsMatrix< index_t > &globals, index_t comp=0) const
Computes the global indices of the input local indices.
Definition: gsDofMapper.cpp:25
Maintains a mapping from patch-local dofs to global dof indices and allows the elimination of individ...
Definition: gsDofMapper.h:68
virtual T getMinCellLength() const
Get the minimum mesh size, as expected for inverse inequalities.
Definition: gsBasis.hpp:663
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
const gsFunctionSet< T > * m_dirdata_ptr
Dirichlet function.
Definition: gsVisitorNitsche.h:229
void initialize(const gsBasis< T > &basis, const index_t, const gsOptionList &options, gsQuadRule< T > &rule)
Initialize.
Definition: gsVisitorNitsche.h:82
T m_penalty
Parameter for the linear form.
Definition: gsVisitorNitsche.h:238
void localToGlobal(const gsDofMapper &mapper, const gsMatrix< T > &eliminatedDofs, const index_t patchIndex, gsSparseMatrix< T > &sysMatrix, gsMatrix< T > &rhsMatrix)
Adds the contributions to the sparse system.
Definition: gsVisitorNitsche.h:198
void pushAllFree(const gsMatrix< T > &localMat, const gsMatrix< T > &localRhs, const gsMatrix< index_t > &actives, const size_t r=0, const size_t c=0)
pushAllFree pushes the local system matrix and rhs for an element to the global system,
Definition: gsSparseSystem.h:1082
Class which enables iteration over all elements of a parameter domain.
Definition: gsDomainIterator.h:67
static gsOptionList defaultOptions()
Default options.
Definition: gsVisitorNitsche.h:67
bool getSwitch(const std::string &label) const
Reads value for option label from options.
Definition: gsOptionList.cpp:51
gsVisitorNitsche(const gsPde< T > &pde, const boundary_condition< T > &bc)
Constructor.
Definition: gsVisitorNitsche.h:62
void localToGlobal(const index_t patchIndex, const std::vector< gsMatrix< T > > &, gsSparseSystem< T > &system)
Adds the contributions to the sparse system.
Definition: gsVisitorNitsche.h:186
Interface for the set of functions defined on a domain (the total number of functions in the set equa...
Definition: gsFuncData.h:23
Visitor for adding the terms associated to weak (Nitsche-type) imposition of the Dirichlet boundary c...
Definition: gsVisitorNitsche.h:54
void assemble(gsDomainIterator< T > &, const gsVector< T > &quWeights)
Assemble on element.
Definition: gsVisitorNitsche.h:147
The density of the measure pull back.
Definition: gsForwardDeclarations.h:76
Real askReal(const std::string &label, const Real &value=0) const
Reads value for option label from options.
Definition: gsOptionList.cpp:139
Class that defines a boundary condition for a side of a patch for some unknown variable of a PDE...
Definition: gsBoundaryConditions.h:106
Value of the object.
Definition: gsForwardDeclarations.h:72
void addReal(const std::string &label, const std::string &desc, const Real &value)
Adds a option named label, with description desc and value value.
Definition: gsOptionList.cpp:211
virtual short_t maxDegree() const
If the basis is of polynomial or piecewise polynomial type, then this function returns the maximum po...
Definition: gsBasis.hpp:638
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
void evaluate(const gsBasis< T > &basis, const gsGeometry< T > &geo, const gsMatrix< T > &quNodes)
Evaluate on element.
Definition: gsVisitorNitsche.h:121
void addSwitch(const std::string &label, const std::string &desc, const bool &value)
Adds a option named label, with description desc and value value.
Definition: gsOptionList.cpp:235
T m_alpha
Parameter for the linear form.
Definition: gsVisitorNitsche.h:232
index_t patch
The index of the patch.
Definition: gsBoundary.h:234
short_t direction() const
Returns the parametric direction orthogonal to this side.
Definition: gsBoundary.h:113
patchSide m_side
Patch side.
Definition: gsVisitorNitsche.h:241
A basis represents a family of scalar basis functions defined over a common parameter domain...
Definition: gsBasis.h:78
T m_beta
Parameter for the linear form.
Definition: gsVisitorNitsche.h:235
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