G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsVisitorNeumann.h
Go to the documentation of this file.
1
14#pragma once
15
16#include <gsCore/gsFuncData.h>
17
18namespace gismo
19{
20
32template <class T>
34{
35public:
36
42 : neudata_ptr( bc.function().get() ), side(bc.side())
43 { GISMO_UNUSED(pde); }
44
50 : neudata_ptr(&neudata), side(s)
51 { }
52
54 void initialize(const gsBasis<T> & basis,
55 gsQuadRule<T> & rule)
56 {
57 const int dir = side.direction();
58 gsVector<int> numQuadNodes ( basis.dim() );
59 for (int i = 0; i < basis.dim(); ++i)
60 numQuadNodes[i] = basis.degree(i) + 1;
61 numQuadNodes[dir] = 1;
62
63 // Setup Quadrature
64 rule = gsGaussRule<T>(numQuadNodes);// harmless slicing occurs here
65
66 // Set Geometry evaluation flags
67 md.flags = NEED_VALUE|NEED_JACOBIAN;
68 }
69
71 void initialize(const gsBasis<T> & basis,
72 const index_t ,
73 const gsOptionList & options,
74 gsQuadRule<T> & rule)
75 {
76 // Setup Quadrature (harmless slicing occurs)
77 rule = gsQuadrature::get(basis, options, side.direction());
78
79 // Set Geometry evaluation flags
81 }
82
84 inline void evaluate(const gsBasis<T> & basis, // to do: more unknowns
85 const gsGeometry<T> & geo,
86 // todo: add element here for efficiency
87 const gsMatrix<T> & quNodes)
88 {
89 md.points = quNodes;
90 // Compute the active basis functions
91 // Assumes actives are the same for all quadrature points on the current element
92 basis.active_into(md.points.col(0) , actives);
93 const index_t numActive = actives.rows();
94
95 // Evaluate basis functions on element
96 basis.eval_into(md.points, basisData);
97
98 // Compute geometry related values
99 geo.computeMap(md);
100
101 // Evaluate the Neumann data
102 neudata_ptr->eval_into(md.values[0], neuData);
103
104 // Initialize local matrix/rhs
105 localRhs.setZero(numActive, neudata_ptr->targetDim() );
106 }
107
110 const gsVector<T> & quWeights)
111 {
112 for (index_t k = 0; k < quWeights.rows(); ++k) // loop over quadrature nodes
113 {
114 // Compute the outer normal vector on the side
115 outerNormal(md, k, side, unormal);
116
117 // Multiply quadrature weight by the measure of normal
118 const T weight = quWeights[k] * unormal.norm();
119
120 localRhs.noalias() += weight * basisData.col(k) * neuData.col(k).transpose() ;
121 }
122 }
123
125 inline void localToGlobal(const index_t patchIndex,
126 const std::vector<gsMatrix<T> > & ,
127 gsSparseSystem<T> & system)
128 {
129 // Map patch-local DoFs to global DoFs
130 system.mapColIndices(actives, patchIndex, actives);
131
132 // Add contributions to the system matrix and right-hand side
133 system.pushToRhs(localRhs, actives, 0);
134 }
135
137 void localToGlobal(const gsDofMapper & mapper,
138 const gsMatrix<T> & eliminatedDofs,
139 const index_t patchIndex,
140 gsSparseMatrix<T> & sysMatrix,
141 gsMatrix<T> & rhsMatrix )
142 {
143 // Local DoFs to global DoFs
144 mapper.localToGlobal(actives, patchIndex, actives);
145 const index_t numActive = actives.rows();
146
147 // Push element contribution to the global load vector
148 for (index_t j=0; j!=numActive; ++j)
149 {
150 // convert local DoF index to global DoF index
151 const unsigned jj = actives(j);
152 if (mapper.is_free_index(jj))
153 rhsMatrix.row(jj) += localRhs.row(j);
154 }
155 }
156
157protected:
158
159
160 // Neumann function
161 const gsFunctionSet<T> * neudata_ptr;
162 boxSide side;
163
164 // Basis values
165 gsMatrix<T> basisData;
166 gsMatrix<index_t> actives;
167
168 // Normal and Neumann values
169 gsVector<T> unormal;
170 gsMatrix<T> neuData;
171
172 // Local matrix and rhs
173 gsMatrix<T> localMat;
174 gsMatrix<T> localRhs;
175
176 gsMapData<T> md;
177};
178
179
180} // namespace gismo
Struct which represents a certain side of a box.
Definition gsBoundary.h:85
short_t direction() const
Returns the parametric direction orthogonal to this side.
Definition gsBoundary.h:113
A basis represents a family of scalar basis functions defined over a common parameter domain.
Definition gsBasis.h:79
Maintains a mapping from patch-local dofs to global dof indices and allows the elimination of individ...
Definition gsDofMapper.h:69
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
bool is_free_index(index_t gl) const
Returns true if global dof gl is not eliminated.
Definition gsDofMapper.h:376
Class which enables iteration over all elements of a parameter domain.
Definition gsDomainIterator.h:68
Interface for the set of functions defined on a domain (the total number of functions in the set equa...
Definition gsFunctionSet.h:219
A function from a n-dimensional domain to an m-dimensional image.
Definition gsFunction.h:60
virtual void computeMap(gsMapData< T > &InOut) const
Computes map function data.
Definition gsFunction.hpp:817
Class that represents the (tensor) Gauss-Legendre quadrature rule.
Definition gsGaussRule.h:28
Abstract base class representing a geometry map.
Definition gsGeometry.h:93
the gsMapData is a cache of pre-computed function (map) values.
Definition gsFuncData.h:349
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Class which holds a list of parameters/options, and provides easy access to them.
Definition gsOptionList.h:33
Abstract class representing a PDE (partial differential equation).
Definition gsPde.h:44
Class representing a reference quadrature rule.
Definition gsQuadRule.h:29
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
A sparse linear system indexed by sets of degrees of freedom.
Definition gsSparseSystem.h:30
void pushToRhs(const gsMatrix< T > &localRhs, const gsMatrix< index_t > &actives, const size_t r=0)
pushToRhs pushes the local rhs for an element to the global system
Definition gsSparseSystem.h:884
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
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
Implementation of a Neumann BC for elliptic assemblers.
Definition gsVisitorNeumann.h:34
void initialize(const gsBasis< T > &basis, const index_t, const gsOptionList &options, gsQuadRule< T > &rule)
Initialize.
Definition gsVisitorNeumann.h:71
void assemble(gsDomainIterator< T > &, const gsVector< T > &quWeights)
Assemble on element.
Definition gsVisitorNeumann.h:109
void initialize(const gsBasis< T > &basis, gsQuadRule< T > &rule)
Initialize.
Definition gsVisitorNeumann.h:54
gsVisitorNeumann(const gsFunction< T > &neudata, boxSide s)
Constructor.
Definition gsVisitorNeumann.h:49
gsVisitorNeumann(const gsPde< T > &pde, const boundary_condition< T > &bc)
Constructor.
Definition gsVisitorNeumann.h:41
void localToGlobal(const index_t patchIndex, const std::vector< gsMatrix< T > > &, gsSparseSystem< T > &system)
Adds the contributions to the sparse system.
Definition gsVisitorNeumann.h:125
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 gsVisitorNeumann.h:137
void evaluate(const gsBasis< T > &basis, const gsGeometry< T > &geo, const gsMatrix< T > &quNodes)
Evaluate on element.
Definition gsVisitorNeumann.h:84
#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.
The G+Smo namespace, containing all definitions for the library.
@ NEED_VALUE
Value of the object.
Definition gsForwardDeclarations.h:72
@ NEED_JACOBIAN
Jacobian of the object.
Definition gsForwardDeclarations.h:75
@ NEED_GRAD_TRANSFORM
Gradient transformation matrix.
Definition gsForwardDeclarations.h:77
@ NEED_MEASURE
The density of the measure pull back.
Definition gsForwardDeclarations.h:76
Class that defines a boundary condition for a side of a patch for some unknown variable of a PDE.
Definition gsBoundaryConditions.h:107
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