32class gsVisitorBiharmonic
57 for (
int i = 0; i < basis.dim(); ++i)
58 numQuadNodes[i] = basis.degree(i) + 1;
88 basis.active_into(md.points.col(0), actives);
89 numActive = actives.rows();
95 basis.evalAllDers_into(md.points, 2, basisData,
true);
101 rhs_ptr->eval_into(md.values[0], rhsVals);
104 localMat.setZero(numActive, numActive);
105 localRhs.setZero(numActive, rhsVals.rows());
116 for (
index_t k = 0; k < quWeights.rows(); ++k)
119 const T weight = quWeights[k] * md.measure(k);
122 transformLaplaceHgrad(md, k, basisGrads, basis2ndDerivs, physBasisLaplace);
125 localMat.noalias() += weight * (physBasisLaplace.transpose() * physBasisLaplace);
127 localRhs.noalias() += weight * ( basisVals.col(k) * rhsVals.col(k).transpose() ) ;
140 system.
push(localMat, localRhs, actives, eliminatedDofs[0], 0, 0);
186 std::vector<gsMatrix<T> > basisData;
A basis represents a family of scalar basis functions defined over a common parameter domain.
Definition gsBasis.h:79
A Biharmonic PDE.
Definition gsBiharmonicPde.h:48
Class which enables iteration over all elements of a parameter domain.
Definition gsDomainIterator.h:68
A function from a n-dimensional domain to an m-dimensional image.
Definition gsFunction.h:60
virtual short_t targetDim() const
Dimension of the target space.
Definition gsFunctionSet.h:595
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
A sparse linear system indexed by sets of degrees of freedom.
Definition gsSparseSystem.h:30
void push(const gsMatrix< T > &localMat, const gsMatrix< T > &localRhs, const gsMatrix< index_t > &actives, const gsMatrix< T > &eliminatedDofs, const size_t r=0, const size_t c=0)
push pushes the local system matrix and rhs for an element to the global system,
Definition gsSparseSystem.h:972
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
void initialize(const gsBasis< T > &basis, const index_t, const gsOptionList &options, gsQuadRule< T > &rule)
Initialize.
Definition gsVisitorBiharmonic.h:68
void assemble(gsDomainIterator< T > &, const gsVector< T > &quWeights)
Assemble on element.
Definition gsVisitorBiharmonic.h:109
void initialize(const gsBasis< T > &basis, gsQuadRule< T > &rule)
Initialize.
Definition gsVisitorBiharmonic.h:53
void localToGlobal(const index_t patchIndex, const std::vector< gsMatrix< T > > &eliminatedDofs, gsSparseSystem< T > &system)
Adds the contributions to the sparse system.
Definition gsVisitorBiharmonic.h:132
gsVisitorBiharmonic(const gsFunction< T > &rhs)
Constructor.
Definition gsVisitorBiharmonic.h:46
void evaluate(const gsBasis< T > &basis, const gsGeometry< T > &geo, gsMatrix< T > &quNodes)
Evaluate on element.
Definition gsVisitorBiharmonic.h:81
gsVisitorBiharmonic(const gsPde< T > &pde)
Constructor.
Definition gsVisitorBiharmonic.h:39
#define index_t
Definition gsConfig.h:32
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
Creates a variety of quadrature rules.
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
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