33template <
class T,
bool paramCoef = false>
52 rhs_ptr = &pde_ptr->rhs()->piece(patchIndex);
69 basis.active_into(md.points.col(0), actives);
70 numActive = actives.rows();
73 basis.evalAllDers_into( md.points, 1, basisData,
true);
81 rhs_ptr->eval_into( (paramCoef ? md.points : md.values[0] ), rhsVals );
84 localMat.setZero(numActive, numActive );
85 localRhs.setZero(numActive, rhsVals.rows() );
95 for (
index_t k = 0; k < quWeights.rows(); ++k)
98 const T weight = quWeights[k] * md.measure(k);
101 transformGradients(md, k, bGrads, physGrad);
103 localRhs.noalias() += weight * ( bVals.col(k) * rhsVals.col(k).transpose() ) ;
104 localMat.noalias() += weight * (physGrad.transpose() * physGrad);
117 system.
push(localMat, localRhs, actives, eliminatedDofs.front(), 0, 0);
126 std::vector<gsMatrix<T> > basisData;
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
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
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
A Poisson PDE.
Definition gsPoissonPde.h:35
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
Visitor for the Poisson equation.
Definition gsVisitorPoisson.h:35
void localToGlobal(const index_t patchIndex, const std::vector< gsMatrix< T > > &eliminatedDofs, gsSparseSystem< T > &system)
Adds the contributions to the sparse system.
Definition gsVisitorPoisson.h:109
void assemble(gsDomainIterator< T > &, gsVector< T > const &quWeights)
Assemble on element.
Definition gsVisitorPoisson.h:89
gsVisitorPoisson(const gsPde< T > &pde)
Constructor.
Definition gsVisitorPoisson.h:41
void initialize(const gsBasis< T > &basis, const index_t patchIndex, const gsOptionList &options, gsQuadRule< T > &rule)
Initialize.
Definition gsVisitorPoisson.h:46
void evaluate(const gsBasis< T > &basis, const gsGeometry< T > &geo, const gsMatrix< T > &quNodes)
Evaluate on element.
Definition gsVisitorPoisson.h:62
#define index_t
Definition gsConfig.h:32
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