25class gsVisitorMassElasticity
29 gsVisitorMassElasticity(gsSparseMatrix<T> * elimMatrix =
nullptr) :
30 elimMat(elimMatrix) {}
32 void initialize(
const gsBasisRefs<T> & basisRefs,
34 const gsOptionList & options,
39 dim = basisRefs.front().dim();
43 density = options.getReal(
"Density");
45 globalIndices.resize(dim);
46 blockNumbers.resize(dim);
49 inline void evaluate(
const gsBasisRefs<T> & basisRefs,
50 const gsGeometry<T> & geo,
51 const gsMatrix<T> & quNodes)
60 basisRefs.front().eval_into(quNodes,basisValuesDisp);
62 basisRefs.front().active_into(quNodes.col(0),localIndicesDisp);
63 N_D = localIndicesDisp.rows();
67 inline void assemble(gsDomainIterator<T> & element,
68 const gsVector<T> & quWeights)
72 localMat.setZero(dim*N_D,dim*N_D);
73 block = density*basisValuesDisp * quWeights.asDiagonal() * md.measures.asDiagonal() * basisValuesDisp.transpose();
74 for (
short_t d = 0; d < dim; ++d)
75 localMat.block(d*N_D,d*N_D,N_D,N_D) = block.block(0,0,N_D,N_D);
78 inline void localToGlobal(
const int patchIndex,
79 const std::vector<gsMatrix<T> > & eliminatedDofs,
80 gsSparseSystem<T> & system)
83 for (
short_t d = 0; d < dim; ++d)
85 system.mapColIndices(localIndicesDisp, patchIndex, globalIndices[d], d);
86 blockNumbers.at(d) = d;
89 system.pushToMatrix(localMat,globalIndices,eliminatedDofs,blockNumbers,blockNumbers);
92 if (elimMat !=
nullptr)
96 for (
short_t dJ = 0; dJ < dim; ++dJ)
98 for (
short_t dI = 0; dI < dim; ++dI)
99 for (
index_t i = 0; i < N_D; ++i)
100 if (system.colMapper(dI).is_free_index(globalIndices[dI].at(i)))
102 system.mapToGlobalRowIndex(localIndicesDisp.at(i),patchIndex,globalI,dI);
103 for (
index_t j = 0; j < N_D; ++j)
104 if (!system.colMapper(dJ).is_free_index(globalIndices[dJ].at(j)))
106 globalElimJ = system.colMapper(dJ).global_to_bindex(globalIndices[dJ].at(j));
107 elimMat->coeffRef(globalI,elimSize+globalElimJ) += localMat(N_D*dI+i,N_D*dJ+j);
110 elimSize += eliminatedDofs[dJ].rows();
123 gsMatrix<T> localMat;
125 gsMatrix<index_t> localIndicesDisp;
129 gsMatrix<T> basisValuesDisp;
132 gsSparseMatrix<T> * elimMat;
137 std::vector< gsMatrix<index_t> > globalIndices;
138 gsVector<index_t> blockNumbers;
#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_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