28class gsVisitorMixedNonLinearElasticity
31 gsVisitorMixedNonLinearElasticity(
const gsPde<T> & pde_,
const gsMultiPatch<T> & displacement_,
32 const gsMultiPatch<T> & pressure_)
33 : pde_ptr(static_cast<const gsBasePde<T>*>(&pde_)),
34 displacement(displacement_),
37 void initialize(
const gsBasisRefs<T> & basisRefs,
39 const gsOptionList & options,
44 dim = basisRefs.front().dim();
50 T YM = options.getReal(
"YoungsModulus");
51 T PR = options.getReal(
"PoissonsRatio");
52 lambda_inv = ( 1. + PR ) * ( 1. - 2. * PR ) / YM / PR ;
53 mu = YM / ( 2. * ( 1. + PR ) );
54 forceScaling = options.getReal(
"ForceScaling");
55 I = gsMatrix<T>::Identity(dim,dim);
57 globalIndices.resize(dim+1);
58 blockNumbers.resize(dim+1);
61 inline void evaluate(
const gsBasisRefs<T> & basisRefs,
62 const gsGeometry<T> & geo,
63 const gsMatrix<T> & quNodes)
74 basisRefs.front().active_into(quNodes.col(0),localIndicesDisp);
75 N_D = localIndicesDisp.rows();
76 basisRefs.back().active_into(quNodes.col(0), localIndicesPres);
77 N_P = localIndicesPres.rows();
79 basisRefs.front().evalAllDers_into(quNodes,1,basisValuesDisp);
81 basisRefs.back().eval_into(quNodes,basisValuesPres);
83 pde_ptr->rhs()->eval_into(md.values[0],forceValues);
85 mdDisplacement.points = quNodes;
89 displacement.patch(patch).computeMap(mdDisplacement);
92 pressure.patch(patch).eval_into(quNodes,pressureValues);
95 inline void assemble(gsDomainIterator<T> & element,
96 const gsVector<T> & quWeights)
100 localMat.setZero(dim*N_D + N_P, dim*N_D + N_P);
101 localRhs.setZero(dim*N_D + N_P,1);
103 for (
index_t q = 0; q < quWeights.rows(); ++q)
106 const T weight = quWeights[q] * md.measure(q);
108 transformGradients(md,q,basisValuesDisp[1],physGradDisp);
110 physDispJac = mdDisplacement.jacobian(q)*(md.jacobian(q).cramerInverse());
114 T J = F.determinant();
116 RCG = F.transpose() * F;
119 RCGinv = RCG.cramerInverse();
121 S = (pressureValues.at(q)-mu)*RCGinv + mu*I;
123 symmetricIdentityTensor<T>(C,RCGinv);
124 C *= mu-pressureValues.at(q);
126 for (
index_t i = 0; i < N_D; i++)
128 setB<T>(B_i,F,physGradDisp.col(i));
129 materialTangentTemp = B_i.transpose() * C;
131 geometricTangentTemp = S * physGradDisp.col(i);
133 for (
index_t j = 0; j < N_D; j++)
135 setB<T>(B_j,F,physGradDisp.col(j));
136 materialTangent = materialTangentTemp * B_j;
137 T geometricTangent = geometricTangentTemp.transpose() * physGradDisp.col(j);
139 for (
short_t d = 0; d < dim; ++d)
140 materialTangent(d,d) += geometricTangent;
142 for (
short_t di = 0; di < dim; ++di)
143 for (
short_t dj = 0; dj < dim; ++dj)
144 localMat(di*N_D+i, dj*N_D+j) += weight * materialTangent(di,dj);
148 voigtStress<T>(Svec,S);
150 localResidual = B_i.transpose() * Svec;
151 for (
short_t d = 0; d < dim; d++)
152 localRhs(d*N_D+i) -= weight * localResidual(d);
155 divV = F.cramerInverse().transpose() * physGradDisp;
156 for (
short_t d = 0; d < dim; ++d)
158 block = weight*basisValuesPres.col(q)*divV.row(d);
159 localMat.block(dim*N_D,d*N_D,N_P,N_D) += block.block(0,0,N_P,N_D);
160 localMat.block(d*N_D,dim*N_D,N_D,N_P) += block.transpose().block(0,0,N_D,N_P);
163 if (
abs(lambda_inv) > 0)
164 localMat.block(dim*N_D,dim*N_D,N_P,N_P) -=
165 (weight*lambda_inv*basisValuesPres.col(q)*basisValuesPres.col(q).transpose()).block(0,0,N_P,N_P);
167 localRhs.middleRows(dim*N_D,N_P) += weight*basisValuesPres.col(q)*(lambda_inv*pressureValues.at(q)-log(J));
169 for (
short_t d = 0; d < dim; ++d)
170 localRhs.middleRows(d*N_D,N_D).noalias() += weight * forceScaling * forceValues(d,q) * basisValuesDisp[0].col(q) ;
174 inline void localToGlobal(
const int patchIndex,
175 const std::vector<gsMatrix<T> > & eliminatedDofs,
176 gsSparseSystem<T> & system)
179 for (
short_t d = 0; d < dim; ++d)
181 system.mapColIndices(localIndicesDisp,patchIndex,globalIndices[d],d);
182 blockNumbers.at(d) = d;
185 system.mapColIndices(localIndicesPres, patchIndex, globalIndices[dim], dim);
186 blockNumbers.at(dim) = dim;
188 system.pushToRhs(localRhs,globalIndices,blockNumbers);
189 system.pushToMatrix(localMat,globalIndices,eliminatedDofs,blockNumbers,blockNumbers);
195 const gsBasePde<T> * pde_ptr;
198 T lambda_inv, mu, forceScaling;
202 gsMatrix<T> localMat;
203 gsMatrix<T> localRhs;
205 gsMatrix<index_t> localIndicesDisp;
206 gsMatrix<index_t> localIndicesPres;
211 std::vector<gsMatrix<T> > basisValuesDisp;
214 gsMatrix<T> basisValuesPres;
216 gsMatrix<T> forceValues;
218 const gsMultiPatch<T> & displacement;
220 gsMapData<T> mdDisplacement;
222 const gsMultiPatch<T> & pressure;
224 gsMatrix<T> pressureValues;
227 gsMatrix<T> C, Ctemp, physGradDisp, physDispJac, F, RCG, E, S, RCGinv, B_i, materialTangentTemp, B_j, materialTangent, divV, block, I;
228 gsVector<T> geometricTangentTemp, Svec, localResidual;
230 std::vector< gsMatrix<index_t> > globalIndices;
231 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
#define GISMO_ENSURE(cond, message)
Definition gsDebug.h:102
This object is a cache for computed values from an evaluator.
Creates a variety of quadrature rules.
Tensor operations for elasticity.
EIGEN_STRONG_INLINE abs_expr< E > abs(const E &u)
Absolute value.
Definition gsExpressions.h:4486
The G+Smo namespace, containing all definitions for the library.
@ NEED_VALUE
Value of the object.
Definition gsForwardDeclarations.h:72
@ NEED_DERIV
Gradient of the object.
Definition gsForwardDeclarations.h:73
@ 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