32 gsVisitorMuscle(
const gsPde<T> & pde_,
33 const gsPiecewiseFunction<T> & muscleTendon_,
34 const gsVector<T> & fiberDir_,
35 const gsMultiPatch<T> & displacement_,
36 const gsMultiPatch<T> & pressure_)
37 : pde_ptr(static_cast<const gsBasePde<T>*>(&pde_)),
38 muscleTendon(muscleTendon_),
40 displacement(displacement_),
43 void initialize(
const gsBasisRefs<T> & basisRefs,
45 const gsOptionList & options,
49 dim = basisRefs.front().dim();
55 T YM_muscle = options.getReal(
"MuscleYoungsModulus");
56 T PR_muscle = options.getReal(
"MusclePoissonsRatio");
57 T YM_tendon = options.getReal(
"TendonYoungsModulus");
58 T PR_tendon = options.getReal(
"TendonPoissonsRatio");
59 lambda_inv_muscle = ( 1. + PR_muscle ) * ( 1. - 2. * PR_muscle ) / YM_muscle / PR_muscle ;
60 mu_muscle = YM_muscle / ( 2. * ( 1. + PR_muscle ) );
61 lambda_inv_tendon = ( 1. + PR_tendon) * ( 1. - 2. * PR_tendon) / YM_tendon / PR_tendon;
62 mu_tendon = YM_tendon / ( 2. * ( 1. + PR_tendon ) );
63 forceScaling = options.getReal(
"ForceScaling");
64 maxMuscleStress = options.getReal(
"MaxMuscleStress");
65 optFiberStretch = options.getReal(
"OptFiberStretch");
66 deltaW = options.getReal(
"DeltaW");
67 powerNu = options.getReal(
"PowerNu");
68 alpha = options.getReal(
"Alpha");
69 I = gsMatrix<T>::Identity(dim,dim);
71 globalIndices.resize(dim+1);
72 blockNumbers.resize(dim+1);
75 inline void evaluate(
const gsBasisRefs<T> & basisRefs,
76 const gsGeometry<T> & geo,
77 const gsMatrix<T> & quNodes)
88 basisRefs.front().active_into(quNodes.col(0),localIndicesDisp);
89 N_D = localIndicesDisp.rows();
90 basisRefs.back().active_into(quNodes.col(0), localIndicesPres);
91 N_P = localIndicesPres.rows();
93 basisRefs.front().evalAllDers_into(quNodes,1,basisValuesDisp);
95 basisRefs.back().eval_into(quNodes,basisValuesPres);
97 pde_ptr->rhs()->eval_into(md.values[0],forceValues);
99 mdDisplacement.points = quNodes;
103 displacement.patch(patch).computeMap(mdDisplacement);
106 pressure.patch(patch).eval_into(quNodes,pressureValues);
108 muscleTendon.piece(patch).eval_into(quNodes,muscleTendonValues);
111 inline void assemble(gsDomainIterator<T> & element,
112 const gsVector<T> & quWeights)
115 localMat.setZero(dim*N_D + N_P, dim*N_D + N_P);
116 localRhs.setZero(dim*N_D + N_P,1);
118 for (
index_t q = 0; q < quWeights.rows(); ++q)
121 const T mu = muscleTendonValues.at(q) * mu_muscle + (1-muscleTendonValues.at(q))*mu_tendon;
122 const T lambda_inv = muscleTendonValues.at(q) * lambda_inv_muscle + (1-muscleTendonValues.at(q))*lambda_inv_tendon;
124 const T weight = quWeights[q] * md.measure(q);
126 transformGradients(md,q,basisValuesDisp[1],physGradDisp);
128 physDispJac = mdDisplacement.jacobian(q)*(md.jacobian(q).cramerInverse());
132 T J = F.determinant();
134 RCG = F.transpose() * F;
137 RCGinv = RCG.cramerInverse();
139 S = (pressureValues.at(q)-mu)*RCGinv + mu*I;
142 fiberDirPhys = md.jacobian(q)*fiberDir;
143 fiberDirPhys /= fiberDirPhys.norm();
145 M = fiberDirPhys * fiberDirPhys.transpose();
147 T fiberStretch = sqrt((M*RCG).trace());
148 T ratioInExp = (fiberStretch/optFiberStretch-1)/deltaW;
149 T megaExp = exp(-1*pow(
abs(ratioInExp),powerNu));
150 S += M * maxMuscleStress * alpha * muscleTendonValues.at(q)/ pow(fiberStretch,2) * megaExp;
153 symmetricIdentityTensor<T>(C,RCGinv);
154 C *= mu-pressureValues.at(q);
157 matrixTraceTensor<T>(Ctemp,M,M);
158 C += -1*Ctemp*alpha*maxMuscleStress*megaExp/pow(fiberStretch,3)* muscleTendonValues.at(q)*
159 (2 + powerNu*pow(ratioInExp,powerNu-1)/deltaW/optFiberStretch);
162 for (
index_t i = 0; i < N_D; i++)
164 setB<T>(B_i,F,physGradDisp.col(i));
165 materialTangentTemp = B_i.transpose() * C;
167 geometricTangentTemp = S * physGradDisp.col(i);
169 for (
index_t j = 0; j < N_D; j++)
171 setB<T>(B_j,F,physGradDisp.col(j));
172 materialTangent = materialTangentTemp * B_j;
173 T geometricTangent = geometricTangentTemp.transpose() * physGradDisp.col(j);
175 for (
short_t d = 0; d < dim; ++d)
176 materialTangent(d,d) += geometricTangent;
178 for (
short_t di = 0; di < dim; ++di)
179 for (
short_t dj = 0; dj < dim; ++dj)
180 localMat(di*N_D+i, dj*N_D+j) += weight * materialTangent(di,dj);
184 voigtStress<T>(Svec,S);
186 localResidual = B_i.transpose() * Svec;
187 for (
short_t d = 0; d < dim; d++)
188 localRhs(d*N_D+i) -= weight * localResidual(d);
191 divV = F.cramerInverse().transpose() * physGradDisp;
192 for (
short_t d = 0; d < dim; ++d)
194 block = weight*basisValuesPres.col(q)*divV.row(d);
195 localMat.block(dim*N_D,d*N_D,N_P,N_D) += block.block(0,0,N_P,N_D);
196 localMat.block(d*N_D,dim*N_D,N_D,N_P) += block.transpose().block(0,0,N_D,N_P);
199 if (
abs(lambda_inv) > 0)
200 localMat.block(dim*N_D,dim*N_D,N_P,N_P) -=
201 (weight*lambda_inv*basisValuesPres.col(q)*basisValuesPres.col(q).transpose()).block(0,0,N_P,N_P);
203 localRhs.middleRows(dim*N_D,N_P) += weight*basisValuesPres.col(q)*(lambda_inv*pressureValues.at(q)-log(J));
205 for (
short_t d = 0; d < dim; ++d)
206 localRhs.middleRows(d*N_D,N_D).noalias() += weight * forceScaling * forceValues(d,q) * basisValuesDisp[0].col(q) ;
210 inline void localToGlobal(
const int patchIndex,
211 const std::vector<gsMatrix<T> > & eliminatedDofs,
212 gsSparseSystem<T> & system)
215 for (
short_t d = 0; d < dim; ++d)
217 system.mapColIndices(localIndicesDisp,patchIndex,globalIndices[d],d);
218 blockNumbers.at(d) = d;
221 system.mapColIndices(localIndicesPres, patchIndex, globalIndices[dim], dim);
222 blockNumbers.at(dim) = dim;
224 system.pushToRhs(localRhs,globalIndices,blockNumbers);
225 system.pushToMatrix(localMat,globalIndices,eliminatedDofs,blockNumbers,blockNumbers);
231 const gsBasePde<T> * pde_ptr;
233 const gsPiecewiseFunction<T> & muscleTendon;
235 const gsVector<T> & fiberDir;
238 T lambda_inv_muscle, mu_muscle, lambda_inv_tendon, mu_tendon, forceScaling;
240 T maxMuscleStress, optFiberStretch, deltaW, powerNu, alpha;
244 gsMatrix<T> localMat;
245 gsMatrix<T> localRhs;
247 gsMatrix<index_t> localIndicesDisp;
248 gsMatrix<index_t> localIndicesPres;
253 std::vector<gsMatrix<T> > basisValuesDisp;
256 gsMatrix<T> basisValuesPres;
258 gsMatrix<T> forceValues;
260 const gsMultiPatch<T> & displacement;
262 gsMapData<T> mdDisplacement;
264 const gsMultiPatch<T> & pressure;
266 gsMatrix<T> pressureValues;
268 gsMatrix<T> muscleTendonValues;
271 gsMatrix<T> C, Ctemp, physGradDisp, physDispJac, F, RCG, E, S, RCGinv, B_i, materialTangentTemp, B_j, materialTangent, divV, block, I, M;
272 gsVector<T> geometricTangentTemp, Svec, localResidual, fiberDirPhys;
274 std::vector< gsMatrix<index_t> > globalIndices;
275 gsVector<index_t> blockNumbers;
Gradient of the object.
Definition: gsForwardDeclarations.h:73
#define short_t
Definition: gsConfig.h:35
#define index_t
Definition: gsConfig.h:32
#define GISMO_ENSURE(cond, message)
Definition: gsDebug.h:102
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:48
Creates a variety of quadrature rules.
IMHO, a useless class, but it is necessary to use the gsAssembler class. Contains proper information ...
The density of the measure pull back.
Definition: gsForwardDeclarations.h:76
Tensor operations for elasticity.
Value of the object.
Definition: gsForwardDeclarations.h:72
EIGEN_STRONG_INLINE abs_expr< E > abs(const E &u)
Absolute value.
Definition: gsExpressions.h:4486
Gradient transformation matrix.
Definition: gsForwardDeclarations.h:77
This object is a cache for computed values from an evaluator.