44 : dirdata_ptr(&dirdata), m_penalty(penalty), side(s)
53 for (
int i = 0; i < basis.dim(); ++i)
54 numQuadNodes[i] = basis.degree(i) + 1;
55 numQuadNodes[dir] = 1;
73 basis.active_into(md.points.col(0), actives);
74 const index_t numActive = actives.rows();
77 basis.evalAllDers_into(md.points, 2, basisData,
true);
83 dirdata_ptr->eval_into(md.values[0], dirData);
86 localMat.setZero(numActive, numActive);
87 localRhs.setZero(numActive, dirdata_ptr->targetDim());
94 const unsigned d = element.
dim();
96 const index_t numActive = actives.rows();
98 for (
index_t k = 0; k < quWeights.rows(); ++k)
101 const typename gsMatrix<T>::Block basisGrads =
102 basisData.middleRows(numActive, numActive * d);
103 const typename gsMatrix<T>::Block basis2ndDerivs =
104 basisData.bottomRows(numActive * (d * (d + 1)) / 2);
107 outerNormal(md, k, side, unormal);
110 const T weight = quWeights[k] * unormal.norm();
116 transformGradients(md, k, basisGrads, physBasisGrads);
118 transformLaplaceHgrad(md, k, basisGrads, basis2ndDerivs, physBasisLaplace);
122 const T mu = m_penalty / (0 != h ? h : 1);
125 localRhs.noalias() += weight * ((physBasisLaplace.transpose() + mu * physBasisGrads.transpose() * unormal)
126 * dirData.col(k).transpose());
128 localMat.noalias() += weight * (physBasisGrads.transpose() * unormal * physBasisLaplace
129 + (physBasisGrads.transpose() * unormal * physBasisLaplace).transpose()
130 - mu * physBasisGrads.transpose() * physBasisGrads);
143 const index_t numActive = actives.rows();
146 for (
index_t j=0; j!=numActive; ++j)
148 const unsigned jj = actives(j);
149 rhsMatrix.row(jj) -= localRhs.row(j);
150 for (
index_t i=0; i!=numActive; ++i)
152 const unsigned ii = actives(i);
154 sysMatrix( ii, jj ) -= localMat(i,j);
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
Class which enables iteration over all elements of a parameter domain.
Definition gsDomainIterator.h:68
short_t dim() const
Return dimension of the elements.
Definition gsDomainIterator.h:115
T getCellSize() const
Return the diagonal of the element.
Definition gsDomainIterator.h:171
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 representing a reference quadrature rule.
Definition gsQuadRule.h:29
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
Visitor for the weak imposition of the first-type dirichlet boundary condition.
Definition gsVisitorNitscheBiharmonic.h:35
void initialize(const gsBasis< T > &basis, gsQuadRule< T > &rule)
Initialize.
Definition gsVisitorNitscheBiharmonic.h:48
gsVisitorNitscheBiharmonic(const gsFunction< T > &dirdata, T penalty, boxSide s)
Constructor.
Definition gsVisitorNitscheBiharmonic.h:43
void assemble(gsDomainIterator< T > &element, const gsVector< T > &quWeights)
Assemble on element.
Definition gsVisitorNitscheBiharmonic.h:91
void evaluate(const gsBasis< T > &basis, const gsGeometry< T > &geo, gsMatrix< T > &quNodes)
Evaluate on element.
Definition gsVisitorNitscheBiharmonic.h:65
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 gsVisitorNitscheBiharmonic.h:135
#define index_t
Definition gsConfig.h:32
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