G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsVisitorNeumannBiharmonic.h
Go to the documentation of this file.
1
14#pragma once
15
16namespace gismo
17{
18
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
68
69 }
70
72 void initialize(const gsBasis<T> & basis,
73 const index_t ,
74 const gsOptionList & options,
75 gsQuadRule<T> & rule)
76 {
77 // Setup Quadrature (harmless slicing occurs)
78 rule = gsGaussRule<T>(basis, options.getReal("quA"),
79 options.getInt("quB"),
80 side.direction() );
81
82 // Set Geometry evaluation flags
84 }
85
87 inline void evaluate(const gsBasis<T> & basis, // to do: more unknowns
88 const gsGeometry<T> & geo,
89 // todo: add element here for efficiency
90 gsMatrix<T> & quNodes)
91 {
92 md.points = quNodes;
93 // Compute the active basis functions
94 // Assumes actives are the same for all quadrature points on the elements
95 basis.active_into(md.points.col(0), actives);
96 numActive = actives.rows();
97
98 // Evaluate basis gradients on element
99 basis.deriv_into( md.points, basisGrads);
100
101 // Compute geometry related values
102 geo.computeMap(md);
103
104 // Evaluate the Neumann data
105 neudata_ptr->eval_into(md.values[0], neuData);
106
107 // Initialize local matrix/rhs
108 localRhs.setZero(numActive, neudata_ptr->targetDim() );
109 }
110
113 gsVector<T> const & quWeights)
114 {
115 for (index_t k = 0; k < quWeights.rows(); ++k) // loop over quadrature nodes
116 {
117 // Compute the outer normal vector on the side
118 outerNormal(md, k, side, unormal);
119
120 // Multiply quadrature weight by the measure of normal
121 const T weight = quWeights[k] * unormal.norm();
122 unormal.normalize();
123 //Get gradients of the physical space
124 transformGradients(md, k, basisGrads, physBasisGrad);
125
126 localRhs.noalias() += weight *(( physBasisGrad.transpose() * unormal )* neuData.col(k).transpose());
127 }
128 }
129
131 inline void localToGlobal(const index_t patchIndex,
132 const std::vector<gsMatrix<T> > & ,
133 gsSparseSystem<T> & system)
134 {
135 // Map patch-local DoFs to global DoFs
136 system.mapColIndices(actives, patchIndex, actives);
137
138 // Add contributions to the system matrix and right-hand side
139 system.pushToRhs(localRhs, actives, 0);
140 }
141
142 /*
143 void localToGlobal(const gsDofMapper & mapper,
144 const gsMatrix<T> & eliminatedDofs,
145 const index_t patchIndex,
146 gsSparseMatrix<T> & sysMatrix,
147 gsMatrix<T> & rhsMatrix )
148 {
149 // Local DoFs to global DoFs
150 mapper.localToGlobal(actives, patchIndex, actives);
151
152 // Push element contribution to the global load vector
153 for (index_t j=0; j!=numActive; ++j)
154 {
155 // convert local DoF index to global DoF index
156 const unsigned jj = actives(j);
157 if (mapper.is_free_index(jj))
158 {
159 rhsMatrix.row(jj) += localRhs.row(j);
160 }
161 }
162 }
163 */
164
165protected:
166
167
168 // Neumann function
169 const gsFunctionSet<T> * neudata_ptr;
170 boxSide side;
171
172 // Basis values
173 gsMatrix<T> basisGrads;
174 gsMatrix<index_t> actives;
175
176 // Normal and Neumann values
177 gsMatrix<T> physBasisGrad;
178
179 gsVector<T> unormal;
180 gsMatrix<T> neuData;
181 index_t numActive;
182
183
184 // Local matrix and rhs
185 gsMatrix<T> localRhs;
186
187 gsMapData<T> md;
188};
189
190
191} // 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
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
const index_t & getInt(const std::string &label) const
Reads value for option label from options.
Definition gsOptionList.cpp:37
Real getReal(const std::string &label) const
Reads value for option label from options.
Definition gsOptionList.cpp:44
Abstract class representing a PDE (partial differential equation).
Definition gsPde.h:44
Class representing a reference quadrature rule.
Definition gsQuadRule.h:29
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
Visitor for Neumann boundary condition for the biharmonic equation.
Definition gsVisitorNeumannBiharmonic.h:34
void initialize(const gsBasis< T > &basis, const index_t, const gsOptionList &options, gsQuadRule< T > &rule)
Initialize.
Definition gsVisitorNeumannBiharmonic.h:72
void initialize(const gsBasis< T > &basis, gsQuadRule< T > &rule)
Initialize.
Definition gsVisitorNeumannBiharmonic.h:54
void localToGlobal(const index_t patchIndex, const std::vector< gsMatrix< T > > &, gsSparseSystem< T > &system)
Adds the contributions to the sparse system.
Definition gsVisitorNeumannBiharmonic.h:131
void evaluate(const gsBasis< T > &basis, const gsGeometry< T > &geo, gsMatrix< T > &quNodes)
Evaluate on element.
Definition gsVisitorNeumannBiharmonic.h:87
void assemble(gsDomainIterator< T > &, gsVector< T > const &quWeights)
Assemble on element.
Definition gsVisitorNeumannBiharmonic.h:112
gsVisitorNeumannBiharmonic(const gsFunction< T > &neudata, boxSide s)
Constructor.
Definition gsVisitorNeumannBiharmonic.h:49
gsVisitorNeumannBiharmonic(const gsPde< T > &pde, const boundary_condition< T > &bc)
Constructor.
Definition gsVisitorNeumannBiharmonic.h:41
#define index_t
Definition gsConfig.h:32
#define GISMO_UNUSED(x)
Definition gsDebug.h:112
The G+Smo namespace, containing all definitions for the library.
@ NEED_VALUE
Value of the object.
Definition gsForwardDeclarations.h:72
@ 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