28class gsVisitorNonLinearElasticity
31 gsVisitorNonLinearElasticity(
const gsPde<T> & pde_,
const gsMultiPatch<T> & displacement_)
32 : pde_ptr(static_cast<const gsBasePde<T>*>(&pde_)),
33 displacement(displacement_) { }
35 void initialize(
const gsBasisRefs<T> & basisRefs,
37 const gsOptionList & options,
42 dim = basisRefs.front().dim();
47 materialLaw = options.getInt(
"MaterialLaw");
48 T YM = options.getReal(
"YoungsModulus");
49 T PR = options.getReal(
"PoissonsRatio");
50 lambda = YM * PR / ( ( 1. + PR ) * ( 1. - 2. * PR ) );
51 mu = YM / ( 2. * ( 1. + PR ) );
52 forceScaling = options.getReal(
"ForceScaling");
53 localStiffening = options.getReal(
"LocalStiff");
55 I = gsMatrix<T>::Identity(dim,dim);
58 matrixTraceTensor<T>(C,I,I);
60 symmetricIdentityTensor<T>(Ctemp,I);
64 globalIndices.resize(dim);
65 blockNumbers.resize(dim);
68 inline void evaluate(
const gsBasisRefs<T> & basisRefs,
69 const gsGeometry<T> & geo,
70 const gsMatrix<T> & quNodes)
81 basisRefs.front().active_into(quNodes.col(0),localIndicesDisp);
82 N_D = localIndicesDisp.rows();
84 basisRefs.front().evalAllDers_into(quNodes,1,basisValuesDisp);
86 pde_ptr->rhs()->eval_into(md.values[0],forceValues);
88 mdDisplacement.points = quNodes;
92 displacement.patch(patch).computeMap(mdDisplacement);
95 inline void assemble(gsDomainIterator<T> & element,
96 const gsVector<T> & quWeights)
100 localMat.setZero(dim*N_D,dim*N_D);
101 localRhs.setZero(dim*N_D,1);
103 for (
index_t q = 0; q < quWeights.rows(); ++q)
105 const T weightForce = quWeights[q] * md.measure(q);
107 transformGradients(md,q,basisValuesDisp[1],physGrad);
109 physDispJac = mdDisplacement.jacobian(q)*(md.jacobian(q).cramerInverse());
113 T J = F.determinant();
115 RCG = F.transpose() * F;
118 const T weightBody = quWeights[q] * pow(md.measure(q),-1.*localStiffening) * md.measure(q);
120 if (materialLaw == 0)
121 S = lambda*E.trace()*I + 2*mu*E;
122 if (materialLaw == 1)
125 RCGinv = RCG.cramerInverse();
126 S = (lambda*log(J)-mu)*RCGinv + mu*I;
128 matrixTraceTensor<T>(C,RCGinv,RCGinv);
130 symmetricIdentityTensor<T>(Ctemp,RCGinv);
131 C += (mu-lambda*log(J))*Ctemp;
133 if (materialLaw == 2)
135 RCGinv = RCG.cramerInverse();
136 S = (lambda*(J*J-1)/2-mu)*RCGinv + mu*I;
138 matrixTraceTensor<T>(C,RCGinv,RCGinv);
140 symmetricIdentityTensor<T>(Ctemp,RCGinv);
141 C += (mu-lambda*(J*J-1)/2)*Ctemp;
144 for (
index_t i = 0; i < N_D; i++)
147 setB<T>(B_i,F,physGrad.col(i));
148 materialTangentTemp = B_i.transpose() * C;
150 geometricTangentTemp = S * physGrad.col(i);
152 for (
index_t j = 0; j < N_D; j++)
154 setB<T>(B_j,F,physGrad.col(j));
156 materialTangent = materialTangentTemp * B_j;
157 T geometricTangent = geometricTangentTemp.transpose() * physGrad.col(j);
159 for (
short_t d = 0; d < dim; ++d)
160 materialTangent(d,d) += geometricTangent;
162 for (
short_t di = 0; di < dim; ++di)
163 for (
short_t dj = 0; dj < dim; ++dj)
164 localMat(di*N_D+i, dj*N_D+j) += weightBody * materialTangent(di,dj);
167 voigtStress<T>(Svec,S);
169 localResidual = B_i.transpose() * Svec;
170 for (
short_t d = 0; d < dim; d++)
171 localRhs(d*N_D+i) -= weightBody * localResidual(d);
174 for (
short_t d = 0; d < dim; ++d)
175 localRhs.middleRows(d*N_D,N_D).noalias() += weightForce * forceScaling * forceValues(d,q) * basisValuesDisp[0].col(q);
179 inline void localToGlobal(
const int patchIndex,
180 const std::vector<gsMatrix<T> > & eliminatedDofs,
181 gsSparseSystem<T> & system)
184 for (
short_t d = 0; d < dim; ++d)
186 system.mapColIndices(localIndicesDisp, patchIndex, globalIndices[d], d);
187 blockNumbers.at(d) = d;
190 system.pushToRhs(localRhs,globalIndices,blockNumbers);
191 system.pushToMatrix(localMat,globalIndices,eliminatedDofs,blockNumbers,blockNumbers);
198 const gsBasePde<T> * pde_ptr;
201 T lambda, mu, forceScaling;
205 gsMatrix<T> localMat;
206 gsMatrix<T> localRhs;
208 gsMatrix<index_t> localIndicesDisp;
213 std::vector<gsMatrix<T> > basisValuesDisp;
215 gsMatrix<T> forceValues;
217 const gsMultiPatch<T> & displacement;
219 gsMapData<T> mdDisplacement;
222 gsMatrix<T> C, Ctemp, physGrad, physDispJac, F, RCG, E, S, RCGinv, B_i, materialTangentTemp, B_j, materialTangent, I;
223 gsVector<T> geometricTangentTemp, Svec, localResidual;
226 std::vector< gsMatrix<index_t> > globalIndices;
227 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.
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