16#include <gsUnstructuredSplines/src/gsC1SurfGluingData.h> 
   21    class gsC1SurfVisitorBasisEdge
 
   25        gsC1SurfVisitorBasisEdge()
 
   33            for (index_t i = 0; i < basis.dim(); ++i) 
 
   34                numQuadNodes[i] = basis.degree(i) + 1;
 
   44        inline void evaluate(
const index_t bfID,
 
   53                             gsC1SurfGluingData<T>  & gluingData,
 
   60            basis.active_into(md.points.col(0), actives);
 
   63            basis.eval_into(md.points, basisData);
 
   68            numActive = actives.rows();
 
   74            T tau_1 = bsp_temp.
knots().at(p + 2);
 
   85                alpha = gluingData.evalAlpha_R(md.points.bottomRows(1));
 
   86                beta = gluingData.evalBeta_R(md.points.bottomRows(1));
 
   98                    beta = isBoundary ? beta.setZero() : beta; 
 
  101                    rhsVals = N_i_plus.cwiseProduct(N_0 + N_1) - temp.cwiseProduct(der_N_i_plus) * tau_1 / p;
 
  103                    localMat.setZero(numActive, numActive);
 
  104                    localRhs.setZero(numActive, rhsVals.rows());
 
  107                else if (typeBf == 
"minus")
 
  112                    alpha = isBoundary ? alpha.setOnes() : alpha; 
 
  114                    rhsVals = alpha.cwiseProduct(N_j_minus.cwiseProduct(N_1));
 
  116                    localMat.setZero(numActive, numActive);
 
  117                    localRhs.setZero(numActive, rhsVals.rows());
 
  124                alpha = gluingData.evalAlpha_L(md.points.topRows(1));
 
  125                beta = gluingData.evalBeta_L(md.points.topRows(1));
 
  131                if (typeBf == 
"plus")
 
  136                    beta = isBoundary ? beta.setZero() : beta; 
 
  139                    rhsVals = N_i_plus.cwiseProduct(N_0 + N_1) - temp.cwiseProduct(der_N_i_plus) * tau_1 / p;
 
  141                    localMat.setZero(numActive, numActive);
 
  142                    localRhs.setZero(numActive, rhsVals.rows());
 
  145                else if (typeBf == 
"minus")
 
  150                    alpha = isBoundary ? alpha.setOnes() : alpha; 
 
  152                    rhsVals = - alpha.cwiseProduct(N_j_minus.cwiseProduct(N_1));
 
  154                    localMat.setZero(numActive, numActive);
 
  155                    localRhs.setZero(numActive, rhsVals.rows());
 
  170                    basisData * quWeights.asDiagonal() *
 
  171                    md.measures.asDiagonal() * basisData.transpose();
 
  173            for (index_t k = 0; k < quWeights.rows(); ++k) 
 
  176                const T weight = quWeights[k] * md.measure(k);
 
  177                localRhs.noalias() += weight * (basisVals.col(k) * rhsVals.col(k).transpose());
 
  182        inline void localToGlobal(
const index_t patchIndex,
 
  189            system.
push(localMat, localRhs, actives, eliminatedDofs[0], 0, 0);
 
A univariate B-spline basis.
Definition gsBSplineBasis.h:700
A basis represents a family of scalar basis functions defined over a common parameter domain.
Definition gsBasis.h:79
virtual short_t maxDegree() const
If the basis is of polynomial or piecewise polynomial type, then this function returns the maximum po...
Definition gsBasis.hpp:687
virtual void evalSingle_into(index_t i, const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluate the i-th basis function at points u into result.
Definition gsBasis.hpp:470
virtual void derivSingle_into(index_t i, const gsMatrix< T > &u, gsMatrix< T > &result) const
Evaluates the (partial) derivatives of the i-th basis function at points u into result.
Definition gsBasis.hpp:478
Class which enables iteration over all elements of a parameter domain.
Definition gsDomainIterator.h:68
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
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
const KnotVectorType & knots(int i=0) const
Returns the knot vector of the basis.
Definition gsBSplineBasis.h:371
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
#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_MEASURE
The density of the measure pull back.
Definition gsForwardDeclarations.h:76