G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsVisitorBiharmonicMixed.h
Go to the documentation of this file.
1
15#pragma once
16
18#include <gsCore/gsFuncData.h>
19
20namespace gismo
21{
22
23template <class T>
24class gsVisitorBiharmonicMixed
25{
26public:
27
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(),
31 elimMat(elimMatrix)
32 {}
33
34 void initialize(const gsBasisRefs<T> & basisRefs,
35 const index_t patchIndex,
36 const gsOptionList & options,
37 gsQuadRule<T> & rule)
38 {
39 GISMO_UNUSED(patchIndex);
40 // a quadrature rule is defined by the basis for the auxiliary variable.
41 // the same rule is used for the main variable
42 rule = gsQuadrature::get(basisRefs.back(), options);
43 // saving necessary info
44 localStiffening = options.getReal("LocalStiff");
45 // resize containers for global indices
46 globalIndices.resize(2);
47 blockNumbers.resize(2);
48 }
49
50 inline void evaluate(const gsBasisRefs<T> & basisRefs,
51 const gsGeometry<T> & geo,
52 const gsMatrix<T> & quNodes)
53 {
54 // store quadrature points of the element for geometry evaluation
55 md.points = quNodes;
56 // NEED_VALUE to get points in the physical domain for evaluation of the RHS
57 // NEED_MEASURE to get the Jacobian determinant values for integration
58 // NEED_GRAD_TRANSFORM to get the Jacobian matrix to transform gradient from the parametric to physical domain
60 // Compute image of the quadrature points plus gradient, jacobian and other necessary data
61 geo.computeMap(md);
62 // find local indices of the main and auxiliary basis functions active on the element
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();
67 // Evaluate basis functions and their derivatives on the element
68 basisRefs.front().evalAllDers_into(quNodes,1,basisValuesMain);
69 basisRefs.back().evalAllDers_into(quNodes,1,basisValuesAux);
70 // Evaluate right-hand side at the image of the quadrature points
71 pde_ptr->rhs()->eval_into(md.values[0],forceValues);
72 }
73
74 inline void assemble(gsDomainIterator<T> & element,
75 const gsVector<T> & quWeights)
76 {
77 GISMO_UNUSED(element);
78 // Initialize local matrix/rhs // 0 | B^T = L
79 localMat.setZero(N_M + N_A, N_M + N_A); // --|-- matrix structure
80 localRhs.setZero(N_M + N_A,pde_ptr->numRhs()); // B | A = 0
81
82 // Loop over the quadrature nodes
83 for (index_t q = 0; q < quWeights.rows(); ++q)
84 {
85 // Multiply quadrature weight by the geometry measure
86 const T weightMatrix = quWeights[q] * pow(md.measure(q),1-localStiffening);
87 const T weightRHS = quWeights[q] * md.measure(q);
88 // Compute physical gradients of the basis functions at q as a 1 x numActiveFunction matrix
89 transformGradients(md, q, basisValuesAux[1], physGradAux);
90 transformGradients(md, q, basisValuesMain[1], physGradMain);
91 // matrix A
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);
94 // matrix B
95 block = weightMatrix * physGradAux.transpose()*physGradMain; // N_A x N_M
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);
98 // rhs contribution
99 localRhs.middleRows(0,N_M).noalias() += weightRHS * basisValuesMain[0].col(q) * forceValues.col(q).transpose();
100 }
101 }
102
103 inline void localToGlobal(const int patchIndex,
104 const std::vector<gsMatrix<T> > & eliminatedDofs,
105 gsSparseSystem<T> & system)
106 {
107 // compute global indices
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;
112 // push to global system
113 system.pushToRhs(localRhs,globalIndices,blockNumbers);
114 system.pushToMatrix(localMat,globalIndices,eliminatedDofs,blockNumbers,blockNumbers);
115
116 // push to the elimination system
117 if (elimMat != nullptr)
118 {
119 index_t globalI,globalElimJ;
120 index_t elimSize = 0;
121 for (short_t dJ = 0; dJ < 2; ++dJ)
122 {
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)))
126 {
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)))
130 {
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);
133 }
134 }
135 elimSize += eliminatedDofs[dJ].rows();
136 }
137 }
138 }
139
140protected:
141 // problem info
142 const gsPoissonPde<T> * pde_ptr;
143 // geometry mapping
144 gsMapData<T> md;
145 // local components of the global linear system
146 gsMatrix<T> localMat;
147 gsMatrix<T> localRhs;
148 // local indices (at the current patch) of basis functions active at the current element
149 gsMatrix<index_t> localIndicesMain;
150 gsMatrix<index_t> localIndicesAux;
151 // number of main and auxiliary basis functions active at the current element
152 index_t N_M, N_A;
153 // values and derivatives of main basis functions at quadrature points at the current element
154 // values are stored as a N_M x numQuadPoints matrix; not sure about derivatives, must be smth like N_M x numQuadPoints
155 // same for the auxiliary basis functions
156 std::vector<gsMatrix<T> > basisValuesMain;
157 std::vector<gsMatrix<T> > basisValuesAux;
158 // RHS values at quadrature points at the current element; stored as a 1 x numQuadPoints matrix
159 gsMatrix<T> forceValues;
160 // all temporary matrices defined here for efficiency
161 gsMatrix<T> block, physGradMain, physGradAux;
162 real_t localStiffening;
163 // elimination matrix to efficiently change Dirichlet degrees of freedom
164 gsSparseMatrix<T> * elimMat;
165 // containers for global indices
166 std::vector< gsMatrix<index_t> > globalIndices;
167 gsVector<index_t> blockNumbers;
168};
169
170} // namespace gismo
#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