29    gsVisitorStokes(
const gsPde<T> & pde_)
 
   30    : dim(), pde_ptr(static_cast<const gsBasePde<T>*>(&pde_)),
 
   31      viscosity(), density(), N_V(), N_P()
 
   34    void initialize(
const gsBasisRefs<T> & basisRefs,
 
   36                    const gsOptionList & options,
 
   41        dim = basisRefs.front().dim();
 
   46        viscosity = options.getReal(
"Viscosity");
 
   47        density = options.getReal(
"Density");
 
   49        globalIndices.resize(dim+1);
 
   50        blockNumbers.resize(dim+1);
 
   53    inline void evaluate(
const gsBasisRefs<T> & basisRefs,
 
   54                         const gsGeometry<T> & geo,
 
   55                         const gsMatrix<T> & quNodes)
 
   66        basisRefs.front().active_into(quNodes.col(0),localIndicesVel);
 
   67        N_V = localIndicesVel.rows();
 
   68        basisRefs.back().active_into(quNodes.col(0), localIndicesPres);
 
   69        N_P = localIndicesPres.rows();
 
   71        basisRefs.front().evalAllDers_into(quNodes,1,basisValuesVel);
 
   73        basisRefs.back().eval_into(quNodes,basisValuesPres);
 
   75        pde_ptr->rhs()->eval_into(md.values[0],forceValues);
 
   78    inline void assemble(gsDomainIterator<T> & element,
 
   79                         const gsVector<T> & quWeights)
 
   83        localMat.setZero(dim*N_V + N_P, dim*N_V + N_P);         
 
   84        localRhs.setZero(dim*N_V + N_P,1);                      
 
   87        for (
index_t q = 0; q < quWeights.rows(); ++q)
 
   90            const T weight = quWeights[q] * md.measure(q);
 
   92            transformGradients(md, q, basisValuesVel[1], physGradVel);
 
   94            block = weight*viscosity*density * physGradVel.transpose()*physGradVel;
 
   95            for (
short_t d = 0; d < dim; ++d)
 
   96                localMat.block(d*N_V,d*N_V,N_V,N_V) += block.block(0,0,N_V,N_V);
 
   98            for (
short_t d = 0; d < dim; ++d)
 
  100                block = weight*basisValuesPres.col(q)*physGradVel.row(d);
 
  101                localMat.block(dim*N_V,d*N_V,N_P,N_V) -= block.block(0,0,N_P,N_V);
 
  102                localMat.block(d*N_V,dim*N_V,N_V,N_P) -= block.transpose().block(0,0,N_V,N_P);
 
  105            for (
short_t d = 0; d < dim; ++d)
 
  106                localRhs.middleRows(d*N_V,N_V).noalias() += weight *density * forceValues(d,q) * basisValuesVel[0].col(q) ;
 
  110    inline void localToGlobal(
const int patchIndex,
 
  111                              const std::vector<gsMatrix<T> > & eliminatedDofs,
 
  112                              gsSparseSystem<T> & system)
 
  115        for (
short_t d = 0; d < dim; ++d)
 
  117            system.mapColIndices(localIndicesVel, patchIndex, globalIndices[d], d);
 
  118            blockNumbers.at(d) = d;
 
  121        system.mapColIndices(localIndicesPres, patchIndex, globalIndices[dim], dim);
 
  122        blockNumbers.at(dim) = dim;
 
  124        system.pushToRhs(localRhs,globalIndices,blockNumbers);
 
  125        system.pushToMatrix(localMat,globalIndices,eliminatedDofs,blockNumbers,blockNumbers);
 
  131    const gsBasePde<T> * pde_ptr;
 
  132    T viscosity, density;
 
  136    gsMatrix<T> localMat;
 
  137    gsMatrix<T> localRhs;
 
  139    gsMatrix<index_t> localIndicesVel;
 
  140    gsMatrix<index_t> localIndicesPres;
 
  145    std::vector<gsMatrix<T> > basisValuesVel;
 
  148    gsMatrix<T> basisValuesPres;
 
  150    gsMatrix<T> forceValues;
 
  153    gsMatrix<T> block, physGradVel;
 
  155    std::vector< gsMatrix<index_t> > globalIndices;
 
  156    gsVector<index_t> blockNumbers;
 
IMHO, a useless class, but it is necessary to use the gsAssembler class. Contains proper information ...
 
#define short_t
Definition gsConfig.h:35
 
#define index_t
Definition gsConfig.h:32
 
#define GISMO_UNUSED(x)
Definition gsDebug.h:112
 
This object is a cache for computed values from an evaluator.
 
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