30 N_M(), N_A(), localStiffening(),
44 localStiffening = options.
getReal(
"LocalStiff");
46 globalIndices.resize(2);
47 blockNumbers.resize(2);
63 basisRefs.
front().active_into(quNodes.col(0),localIndicesMain);
64 N_M = localIndicesMain.rows();
65 basisRefs.
back().active_into(quNodes.col(0), localIndicesAux);
66 N_A = localIndicesAux.rows();
68 basisRefs.
front().evalAllDers_into(quNodes,1,basisValuesMain);
69 basisRefs.
back().evalAllDers_into(quNodes,1,basisValuesAux);
71 pde_ptr->rhs()->eval_into(md.values[0],forceValues);
79 localMat.setZero(N_M + N_A, N_M + N_A);
80 localRhs.setZero(N_M + N_A,pde_ptr->numRhs());
83 for (
index_t q = 0; q < quWeights.rows(); ++q)
86 const T weightMatrix = quWeights[q] * pow(md.measure(q),1-localStiffening);
87 const T weightRHS = quWeights[q] * md.measure(q);
89 transformGradients(md, q, basisValuesAux[1], physGradAux);
90 transformGradients(md, q, basisValuesMain[1], physGradMain);
92 block = weightMatrix * basisValuesAux[0].col(q) * basisValuesAux[0].col(q).transpose();
93 localMat.block(N_M,N_M,N_A,N_A) += block.block(0,0,N_A,N_A);
95 block = weightMatrix * physGradAux.transpose()*physGradMain;
96 localMat.block(0,N_M,N_M,N_A) += block.block(0,0,N_A,N_M).transpose();
97 localMat.block(N_M,0,N_A,N_M) += block.block(0,0,N_A,N_M);
99 localRhs.middleRows(0,N_M).noalias() += weightRHS * basisValuesMain[0].col(q) * forceValues.col(q).transpose();
108 system.
mapColIndices(localIndicesMain,patchIndex,globalIndices[0],0);
109 system.
mapColIndices(localIndicesAux,patchIndex,globalIndices[1],1);
110 blockNumbers.
at(0) = 0;
111 blockNumbers.
at(1) = 1;
113 system.
pushToRhs(localRhs,globalIndices,blockNumbers);
114 system.
pushToMatrix(localMat,globalIndices,eliminatedDofs,blockNumbers,blockNumbers);
117 if (elimMat !=
nullptr)
121 for (
short_t dJ = 0; dJ < 2; ++dJ)
123 for (
short_t dI = 0; dI < 2; ++dI)
124 for (
index_t i = 0; i < N_M; ++i)
128 for (
index_t j = 0; j < N_M; ++j)
132 elimMat->coeffRef(globalI,elimSize+globalElimJ) += localMat(N_M*dI+i,N_M*dJ+j);
135 elimSize += eliminatedDofs[dJ].rows();
156 std::vector<gsMatrix<T> > basisValuesMain;
157 std::vector<gsMatrix<T> > basisValuesAux;
162 real_t localStiffening;
166 std::vector< gsMatrix<index_t> > globalIndices;
Simple class to hold a list of gsBasis which discretize a PDE system on a given patch.
Definition gsBasisRefs.h:26
const gsBasis< T > & back() const
Back.
Definition gsBasisRefs.h:54
const gsBasis< T > & front() const
Front.
Definition gsBasisRefs.h:50
bool is_free_index(index_t gl) const
Returns true if global dof gl is not eliminated.
Definition gsDofMapper.h:376
index_t global_to_bindex(index_t gl) const
Returns the boundary index of global dof gl.
Definition gsDofMapper.h:364
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
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
T at(index_t i) const
Returns the i-th element of the vectorization of the matrix.
Definition gsMatrix.h:211
Class which holds a list of parameters/options, and provides easy access to them.
Definition gsOptionList.h:33
Real getReal(const std::string &label) const
Reads value for option label from options.
Definition gsOptionList.cpp:44
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
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
A sparse linear system indexed by sets of degrees of freedom.
Definition gsSparseSystem.h:30
void pushToMatrix(const gsMatrix< T > &localMat, const gsMatrix< index_t > &actives, const gsMatrix< T > &eliminatedDofs, const size_t r=0, const size_t c=0)
pushToMatrix pushes the local matrix for an element to the global system,
Definition gsSparseSystem.h:638
const gsDofMapper & colMapper(const index_t c) const
colMapper returns the dofMapper for column block c
Definition gsSparseSystem.h:498
void mapToGlobalRowIndex(const index_t active, const index_t patchIndex, index_t &result, const index_t r=0) const
mapToGlobalRowIndex Maps a single basis index to the final position in the system,...
Definition gsSparseSystem.h:600
void pushToRhs(const gsMatrix< T > &localRhs, const gsMatrix< index_t > &actives, const size_t r=0)
pushToRhs pushes the local rhs for an element to the global system
Definition gsSparseSystem.h:884
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
T at(index_t i) const
Returns the i-th element of the vector.
Definition gsVector.h:177
Visitor for the biharmonic equation.
Definition gsVisitorBiharmonic.h:25
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
void evaluate(const gsBasis< T > &basis, const gsGeometry< T > &geo, gsMatrix< T > &quNodes)
Evaluate on element.
Definition gsVisitorBiharmonic.h:81
#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