24class gsVisitorBiharmonicMixed
28 gsVisitorBiharmonicMixed(
const gsPde<T> & pde_, gsSparseMatrix<T> * elimMatrix =
nullptr)
29 : pde_ptr(static_cast<const gsPoissonPde<T>*>(&pde_)),
30 N_M(), N_A(), localStiffening(),
34 void initialize(
const gsBasisRefs<T> & basisRefs,
36 const gsOptionList & options,
44 localStiffening = options.getReal(
"LocalStiff");
46 globalIndices.resize(2);
47 blockNumbers.resize(2);
50 inline void evaluate(
const gsBasisRefs<T> & basisRefs,
51 const gsGeometry<T> & geo,
52 const gsMatrix<T> & quNodes)
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);
74 inline void assemble(gsDomainIterator<T> & element,
75 const gsVector<T> & quWeights)
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();
103 inline void localToGlobal(
const int patchIndex,
104 const std::vector<gsMatrix<T> > & eliminatedDofs,
105 gsSparseSystem<T> & system)
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)
125 if (system.colMapper(dI).is_free_index(globalIndices[dI].at(i)))
127 system.mapToGlobalRowIndex(localIndicesMain.at(i),patchIndex,globalI,dI);
128 for (
index_t j = 0; j < N_M; ++j)
129 if (!system.colMapper(dJ).is_free_index(globalIndices[dJ].at(j)))
131 globalElimJ = system.colMapper(dJ).global_to_bindex(globalIndices[dJ].at(j));
132 elimMat->coeffRef(globalI,elimSize+globalElimJ) += localMat(N_M*dI+i,N_M*dJ+j);
135 elimSize += eliminatedDofs[dJ].rows();
142 const gsPoissonPde<T> * pde_ptr;
146 gsMatrix<T> localMat;
147 gsMatrix<T> localRhs;
149 gsMatrix<index_t> localIndicesMain;
150 gsMatrix<index_t> localIndicesAux;
156 std::vector<gsMatrix<T> > basisValuesMain;
157 std::vector<gsMatrix<T> > basisValuesAux;
159 gsMatrix<T> forceValues;
161 gsMatrix<T> block, physGradMain, physGradAux;
162 real_t localStiffening;
164 gsSparseMatrix<T> * elimMat;
166 std::vector< gsMatrix<index_t> > globalIndices;
167 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_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