![]() |
MechSys
1.0
Computing library for simulations in continuum and discrete mechanics
|
00001 /************************************************************************ 00002 * MechSys - Open Library for Mechanical Systems * 00003 * Copyright (C) 2005 Dorival M. Pedroso, Raul Durand * 00004 * Copyright (C) 2009 Sergio Galindo * 00005 * * 00006 * This program is free software: you can redistribute it and/or modify * 00007 * it under the terms of the GNU General Public License as published by * 00008 * the Free Software Foundation, either version 3 of the License, or * 00009 * any later version. * 00010 * * 00011 * This program is distributed in the hope that it will be useful, * 00012 * but WITHOUT ANY WARRANTY; without even the implied warranty of * 00013 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * 00014 * GNU General Public License for more details. * 00015 * * 00016 * You should have received a copy of the GNU General Public License * 00017 * along with this program. If not, see <http://www.gnu.org/licenses/> * 00018 ************************************************************************/ 00019 00020 #ifndef MECHSYS_FEM_USIGEPSELEM_H 00021 #define MECHSYS_FEM_USIGEPSELEM_H 00022 00023 // MechSys 00024 #include <mechsys/fem/element.h> 00025 #include <mechsys/models/equilibstate.h> 00026 #include <mechsys/models/stressupdate.h> 00027 00028 namespace FEM 00029 { 00030 00031 class USigEpsElem : public Element 00032 { 00033 public: 00034 // Static 00035 static size_t NCo; 00036 static size_t NDs; 00037 static size_t NDu; 00038 00039 // Constructor 00040 USigEpsElem (int NDim, 00041 Mesh::Cell const & Cell, 00042 Model const * Mdl, 00043 Model const * XMdl, 00044 SDPair const & Prp, 00045 SDPair const & Ini, 00046 Array<Node*> const & Nodes); 00047 00048 // Destructor 00049 ~USigEpsElem () { if (GEs!=NULL) delete GEs; } 00050 00051 // Methods 00052 void IncNLocDOF (size_t & NEq) const { FirstEQ = NEq; NEq += 2*NDs; } 00053 void GetLoc (Array<size_t> & Loc) const; 00054 void SetBCs (size_t IdxEdgeOrFace, SDPair const & BCs, BCFuncs * BCF); 00055 void CalcK (Mat_t & K) const; 00056 void UpdateState (Vec_t const & dU, Vec_t * F_int=NULL) const; 00057 void StateKeys (Array<String> & Keys) const; 00058 void StateAtIP (SDPair & KeysVals, int IdxIP) const; 00059 00060 // Internal methods 00061 void Matrices (Mat_t & A, Mat_t & E, Mat_t & Q) const; 00062 void Interp (Mat_t const & C, IntegPoint const & IP, Mat_t & B, Mat_t & N, Mat_t & Ns, double & detJ, double & Coef) const; 00063 00064 // Constants 00065 GeomElem * GEs; 00066 mutable long FirstEQ; 00067 double rho; 00068 }; 00069 00070 size_t USigEpsElem::NCo = 0; 00071 size_t USigEpsElem::NDs = 0; 00072 size_t USigEpsElem::NDu = 0; 00073 00074 00076 00077 00078 inline USigEpsElem::USigEpsElem (int NDim, Mesh::Cell const & Cell, Model const * Mdl, Model const * XMdl, SDPair const & Prp, SDPair const & Ini, Array<Node*> const & Nodes) 00079 : Element(NDim,Cell,Mdl,XMdl,Prp,Ini,Nodes) 00080 { 00081 // check GE 00082 if (GE==NULL) throw new Fatal("USigEpsElem::USigEpsElem: GE (geometry element) must be defined"); 00083 if (GTy==pse_t) throw new Fatal("USigEpsElem::USigEpsElem: This element does not work for plane-stress (pse)"); 00084 if (Mdl==NULL) throw new Fatal("USigEpsElem::USigEpsElem: Model must be defined"); 00085 00086 // local nodes 00087 if (Prp.HasKey("geom_sig")) 00088 { 00089 String geom_sig; 00090 GEOM.Val2Key (Prp("geom_sig"), geom_sig); 00091 GEs = AllocGeomElem (geom_sig, NDim); 00092 } 00093 else GEs = AllocGeomElem (GE->Name, NDim); 00094 00095 // set constants of this class (just once) 00096 if (NDs==0) 00097 { 00098 NCo = 2*NDim; 00099 NDs = NCo*GEs->NN; 00100 NDu = NDim*GE->NN; 00101 } 00102 00103 // parameters/properties 00104 rho = (Prp.HasKey("rho") ? Prp("rho") : 1.0); 00105 00106 // allocate and initialize state at each IP 00107 for (size_t i=0; i<GE->NIP; ++i) 00108 { 00109 Sta.Push (new EquilibState(NDim)); 00110 Mdl->InitIvs (Ini, Sta[i]); 00111 } 00112 00113 // set initial values 00114 if (Ini.HasKey("geostatic")) 00115 { 00116 if (!Ini.HasKey("K0")) throw new Fatal("USigEpsElem::USigEpsElem: For geostatic stresses, 'K0' must be provided in 'Ini' dictionary"); 00117 if (!Ini.HasKey("gam")) throw new Fatal("USigEpsElem::USigEpsElem: For geostatic stresses, 'gam' must be provided in 'Ini' dictionary"); 00118 if (NDim==2 && !Ini.HasKey("y_surf")) throw new Fatal("USigEpsElem::USigEpsElem: For geostatic stresses in 2D, 'y_surf' must be provided in 'Ini' dictionary"); 00119 if (NDim==3 && !Ini.HasKey("z_surf")) throw new Fatal("USigEpsElem::USigEpsElem: For geostatic stresses in 3D, 'z_surf' must be provided in 'Ini' dictionary"); 00120 double K0 = Ini("K0"); 00121 double gam = Ini("gam"); 00122 double surf = (NDim==2 ? Ini("y_surf") : Ini("z_surf")); 00123 Vec_t X; // coords of IP 00124 for (size_t i=0; i<GE->NIP; ++i) 00125 { 00126 CoordsOfIP (i, X); 00127 Vec_t & sig = static_cast<EquilibState *>(Sta[i])->Sig; 00128 if (NDim==2) 00129 { 00130 double hei = fabs(surf-X(1)); 00131 sig(1) = -gam*hei; // sy 00132 sig(0) = K0*sig(1); // sx 00133 sig(2) = K0*sig(1); // sz 00134 sig(3) = 0.0; // sxy*sq2 00135 } 00136 else // 3D 00137 { 00138 double hei = fabs(surf-X(2)); 00139 sig(2) = -gam*hei; // sz 00140 sig(0) = K0*sig(2); // sx 00141 sig(1) = K0*sig(2); // sy 00142 sig(3) = 0.0; // sxy*sq2 00143 sig(4) = 0.0; // syz*sq2 00144 sig(5) = 0.0; // szx*sq2 00145 } 00146 } 00147 } 00148 00149 // set F in nodes due to initial stresses 00150 double detJ, coef; 00151 Mat_t C, B, N, Ns; 00152 Vec_t Fe(NDu); 00153 CoordMatrix (C); 00154 set_to_zero (Fe); 00155 for (size_t i=0; i<GE->NIP; ++i) 00156 { 00157 Interp (C, GE->IPs[i], B, N, Ns, detJ, coef); 00158 Vec_t const & sig = static_cast<EquilibState const *>(Sta[i])->Sig; 00159 Fe += (coef) * trans(B)*sig; 00160 } 00161 for (size_t i=0; i<GE->NN; ++i) 00162 { 00163 Con[i]->F("fx") += Fe(0+i*NDim); 00164 Con[i]->F("fy") += Fe(1+i*NDim); if (NDim==3) 00165 Con[i]->F("fz") += Fe(2+i*NDim); 00166 } 00167 } 00168 00169 inline void USigEpsElem::SetBCs (size_t IdxEdgeOrFace, SDPair const & BCs, BCFuncs * BCF) 00170 { 00171 // deactivate/activate commands 00172 if (BCs.HasKey("deactivate")) throw new Fatal("USigEpsElem::SetBCs: 'deactivate' command failed: method not implemented yet"); 00173 if (BCs.HasKey("activate")) throw new Fatal("USigEpsElem::SetBCs: 'activate' command failed: method not implemented yet"); 00174 00175 // check 00176 if (!Active) throw new Fatal("USigEpsElem::SetBCs: Element %d is inactive",Cell.ID); 00177 00178 bool has_bx = BCs.HasKey("bx"); // x component of body force 00179 bool has_by = BCs.HasKey("by"); // y component of body force 00180 bool has_bz = BCs.HasKey("bz"); // z component of body force 00181 bool has_cbx = BCs.HasKey("cbx"); // centrifugal body force along x (in axisymmetric problems) 00182 bool has_qx = BCs.HasKey("qx"); // x component of distributed loading 00183 bool has_qy = BCs.HasKey("qy"); // y component of distributed loading 00184 bool has_qz = BCs.HasKey("qz"); // z component of distributed loading 00185 bool has_qn = BCs.HasKey("qn"); // normal distributed loading 00186 bool has_qt = BCs.HasKey("qt"); // tangential distributed loading (2D only) 00187 bool has_ux = BCs.HasKey("ux"); // x displacement 00188 bool has_uy = BCs.HasKey("uy"); // y displacement 00189 bool has_uz = BCs.HasKey("uz"); // z displacement 00190 bool has_sup = BCs.HasKey("incsup"); // inclined support 00191 bool has_gra = BCs.HasKey("gravity"); // gravity 00192 00193 // force components specified 00194 if (has_bx || has_by || has_bz || has_cbx || has_gra || 00195 has_qx || has_qy || has_qz || has_qn || has_qt) 00196 { 00197 // body forces 00198 if (has_bx || has_by || has_bz || has_cbx || has_gra) // prescribed body forces 00199 { 00200 // matrix of coordinates of nodes 00201 Mat_t C; 00202 CoordMatrix (C); 00203 00204 // loading 00205 double bx = (has_bx ? BCs("bx") : 0.0); 00206 double by = (has_by ? BCs("by") : 0.0); 00207 double bz = (has_bz ? BCs("bz") : 0.0); 00208 bx = (has_cbx ? BCs("cbx") : bx ); 00209 if (has_gra) 00210 { 00211 if (NDim==2) by = -rho*BCs("gravity"); 00212 else bz = -rho*BCs("gravity"); 00213 } 00214 00215 // set 00216 for (size_t i=0; i<GE->NIP; ++i) 00217 { 00218 // geometric data 00219 GE->Shape (GE->IPs[i].r, GE->IPs[i].s, GE->IPs[i].t); 00220 GE->Derivs (GE->IPs[i].r, GE->IPs[i].s, GE->IPs[i].t); 00221 00222 // Jacobian and its determinant 00223 Mat_t J(GE->dNdR * C); // J = dNdR * C 00224 double detJ = Det(J); 00225 00226 // coefficient used during integration 00227 double coef = detJ*GE->IPs[i].w; 00228 if (GTy==axs_t) 00229 { 00230 // calculate radius=x at this IP 00231 double radius = 0.0; 00232 for (size_t j=0; j<GE->NN; ++j) radius += GE->N(j)*Con[j]->Vert.C[0]; 00233 00234 // correct coef 00235 if (has_cbx) coef *= radius*radius; 00236 else coef *= radius; 00237 } 00238 00239 // set boundary conditions 00240 for (size_t j=0; j<GE->NN; ++j) 00241 { 00242 Con[j]->AddToPF("fx", coef*GE->N(j)*bx, BCF); 00243 Con[j]->AddToPF("fy", coef*GE->N(j)*by, BCF); if (NDim==3) 00244 Con[j]->AddToPF("fz", coef*GE->N(j)*bz, BCF); 00245 } 00246 } 00247 } 00248 00249 // surface loading 00250 if (has_qx || has_qy || has_qz || has_qn || has_qt) 00251 { 00252 // matrix of coordinates of edge/face 00253 Mat_t Cf; 00254 FCoordMatrix (IdxEdgeOrFace, Cf); 00255 00256 // loading 00257 double qx = (has_qx ? BCs("qx") : 0.0); 00258 double qy = (has_qy ? BCs("qy") : 0.0); 00259 double qz = (has_qz ? BCs("qz") : 0.0); 00260 double qn = (has_qn ? BCs("qn") : 0.0); 00261 double qt = (has_qt ? BCs("qt") : 0.0); 00262 00263 // set 00264 for (size_t i=0; i<GE->NFIP; ++i) 00265 { 00266 // geometric data 00267 GE->FaceShape (GE->FIPs[i].r, GE->FIPs[i].s); 00268 GE->FaceDerivs (GE->FIPs[i].r, GE->FIPs[i].s); 00269 00270 // face/edge Jacobian and its determinant 00271 Mat_t J(GE->FdNdR * Cf); 00272 00273 // coefficient used during integration 00274 double coef = GE->FIPs[i].w; // *detJ is not necessary since qx,qy,qz are already multiplied by detJ (due to normal) 00275 00276 if (GTy==axs_t) 00277 { 00278 // calculate radius=x at this FIP 00279 double radius = 0.0; 00280 for (size_t j=0; j<GE->NFN; ++j) radius += GE->FN(j)*Con[GE->FNode(IdxEdgeOrFace,j)]->Vert.C[0]; 00281 coef *= radius; // correct coef 00282 } 00283 00284 // calculate qx, qy and qz from qn and qt 00285 if (has_qn || has_qt) 00286 { 00287 // normal to edge/face 00288 Vec_t n(NDim); // normal multiplied by detJ 00289 if (NDim==2) n = J(0,1), -J(0,0); 00290 else 00291 { 00292 // vectorial product 00293 Vec_t a(3); a = J(0,0), J(0,1), J(0,2); 00294 Vec_t b(3); b = J(1,0), J(1,1), J(1,2); 00295 n = a(1)*b(2) - a(2)*b(1), 00296 a(2)*b(0) - a(0)*b(2), 00297 a(0)*b(1) - a(1)*b(0); 00298 } 00299 00300 // loading 00301 if (NDim==2) 00302 { 00303 qx = n(0)*qn - n(1)*qt; 00304 qy = n(1)*qn + n(0)*qt; 00305 } 00306 else 00307 { 00308 qx = n(0)*qn; 00309 qy = n(1)*qn; 00310 qz = n(2)*qn; 00311 } 00312 } 00313 00314 // set boundary conditions 00315 for (size_t j=0; j<GE->NFN; ++j) 00316 { 00317 size_t k = GE->FNode(IdxEdgeOrFace,j); 00318 Con[k]->AddToPF("fx", coef*GE->FN(j)*qx, BCF); 00319 Con[k]->AddToPF("fy", coef*GE->FN(j)*qy, BCF); if (NDim==3) 00320 Con[k]->AddToPF("fz", coef*GE->FN(j)*qz, BCF); 00321 } 00322 } 00323 } 00324 } 00325 00326 // prescribed displacements 00327 else if (has_ux || has_uy || has_uz || has_sup) 00328 { 00329 double ux = (has_ux ? BCs("ux") : 0.0); 00330 double uy = (has_uy ? BCs("uy") : 0.0); 00331 double uz = (has_uz ? BCs("uz") : 0.0); 00332 double alpha = (has_sup ? BCs("alpha") : 0.0); 00333 for (size_t j=0; j<GE->NFN; ++j) 00334 { 00335 size_t k = GE->FNode(IdxEdgeOrFace,j); 00336 if (has_ux) Con[k]->SetPU("ux", ux, BCF); 00337 if (has_uy) Con[k]->SetPU("uy", uy, BCF); 00338 if (has_uz) Con[k]->SetPU("uz", uz, BCF); 00339 if (has_sup) Con[k]->SetIncSup (alpha); 00340 } 00341 } 00342 } 00343 00344 inline void USigEpsElem::GetLoc (Array<size_t> & Loc) const 00345 { 00346 // Sig and Eps DOFs 00347 Loc.Resize (2*NDs + NDu); 00348 for (size_t i=0; i<2*NDs; ++i) Loc[i] = FirstEQ + i; 00349 00350 // U DOFs 00351 for (size_t i=0; i<GE->NN; ++i) 00352 { 00353 Loc[2*NDs+i*NDim+0] = Con[i]->Eq("ux"); 00354 Loc[2*NDs+i*NDim+1] = Con[i]->Eq("uy"); if (NDim==3) 00355 Loc[2*NDs+i*NDim+2] = Con[i]->Eq("uz"); 00356 } 00357 } 00358 00359 inline void USigEpsElem::Matrices (Mat_t & A, Mat_t & E, Mat_t & Q) const 00360 { 00361 // submatrices 00362 A.change_dim (NDs,NDs); 00363 E.change_dim (NDs,NDu); 00364 Q.change_dim (NDs,NDs); 00365 set_to_zero (A); 00366 set_to_zero (E); 00367 set_to_zero (Q); 00368 00369 // auxiliar matrices 00370 double detJ, coef; 00371 Mat_t N, Ns, B, C, D, Di; 00372 CoordMatrix (C); 00373 00374 for (size_t i=0; i<GE->NIP; ++i) 00375 { 00376 // interpolation 00377 Interp (C, GE->IPs[i], B, N, Ns, detJ, coef); 00378 00379 // calc D 00380 Mdl->Stiffness (Sta[i], D); 00381 00382 // matrices 00383 A += (coef) * trans(Ns)*D*Ns; 00384 E += (coef) * trans(Ns)*B; 00385 Q += (-coef) * trans(Ns)*Ns; 00386 } 00387 } 00388 00389 inline void USigEpsElem::CalcK (Mat_t & K) const 00390 { 00391 Mat_t A, E, Q; 00392 Matrices (A, E, Q); 00393 K.change_dim (2*NDs+NDu, 2*NDs+NDu); 00394 set_to_zero (K); 00395 for (size_t i=0; i<NDs; ++i) 00396 { 00397 for (size_t j=0; j<NDs; ++j) 00398 { 00399 K(i, j) = A(i,j); 00400 K(i,NDs+j) = Q(i,j); 00401 K(NDs+i,j) = Q(j,i); 00402 } 00403 for (size_t j=0; j<NDu; ++j) 00404 { 00405 K( NDs+i, 2*NDs+j) = E(i,j); 00406 K(2*NDs+j, NDs+i) = E(i,j); 00407 } 00408 } 00409 } 00410 00411 inline void USigEpsElem::Interp (Mat_t const & C, IntegPoint const & IP, Mat_t & B, Mat_t & N, Mat_t & Ns, double & detJ, double & Coef) const 00412 { 00413 // deriv of shape func w.r.t natural coordinates 00414 GE ->Shape (IP.r, IP.s, IP.t); 00415 GE ->Derivs (IP.r, IP.s, IP.t); 00416 GEs->Shape (IP.r, IP.s, IP.t); 00417 00418 // Jacobian and its determinant 00419 Mat_t J(GE->dNdR * C); // J = dNdR * C 00420 detJ = Det(J); 00421 00422 // deriv of shape func w.r.t real coordinates 00423 Mat_t Ji; 00424 Inv (J, Ji); 00425 Mat_t dNdX(Ji * GE->dNdR); // dNdX = Inv(J) * dNdR 00426 00427 // coefficient used during integration 00428 Coef = detJ*IP.w; 00429 00430 // B matrix 00431 B.change_dim (NCo, NDu); 00432 set_to_zero (B); 00433 if (NDim==2) 00434 { 00435 if (GTy==axs_t) 00436 { 00437 // correct Coef 00438 double radius = 0.0; // radius=x at this IP 00439 for (size_t j=0; j<GE->NN; ++j) radius += GE->N(j)*Con[j]->Vert.C[0]; 00440 Coef *= radius; 00441 00442 // B matrix 00443 for (size_t i=0; i<GE->NN; ++i) 00444 { 00445 B(0,0+i*NDim) = dNdX(0,i); 00446 B(1,1+i*NDim) = dNdX(1,i); 00447 B(2,0+i*NDim) = GE->N(i)/radius; 00448 B(3,0+i*NDim) = dNdX(1,i)/sqrt(2.0); 00449 B(3,1+i*NDim) = dNdX(0,i)/sqrt(2.0); 00450 } 00451 } 00452 else // pse_t, psa_t 00453 { 00454 for (size_t i=0; i<GE->NN; ++i) 00455 { 00456 B(0,0+i*NDim) = dNdX(0,i); 00457 B(1,1+i*NDim) = dNdX(1,i); 00458 B(3,0+i*NDim) = dNdX(1,i)/sqrt(2.0); 00459 B(3,1+i*NDim) = dNdX(0,i)/sqrt(2.0); 00460 } 00461 } 00462 } 00463 else // 3D 00464 { 00465 for (size_t i=0; i<GE->NN; ++i) 00466 { 00467 B(0,0+i*NDim) = dNdX(0,i); 00468 B(1,1+i*NDim) = dNdX(1,i); 00469 B(2,2+i*NDim) = dNdX(2,i); 00470 B(3,0+i*NDim) = dNdX(1,i)/sqrt(2.0); B(3,1+i*NDim) = dNdX(0,i)/sqrt(2.0); 00471 B(4,1+i*NDim) = dNdX(2,i)/sqrt(2.0); B(4,2+i*NDim) = dNdX(1,i)/sqrt(2.0); 00472 B(5,2+i*NDim) = dNdX(0,i)/sqrt(2.0); B(5,0+i*NDim) = dNdX(2,i)/sqrt(2.0); 00473 } 00474 } 00475 00476 // N matrix 00477 N.change_dim (NDim,NDu); 00478 set_to_zero (N); 00479 for (int i=0; i<NDim; ++i) 00480 for (size_t j=0; j<GE->NN; ++j) 00481 N(i,i+j*NDim) = GE->N(j); 00482 00483 // Ns matrix 00484 Ns.change_dim (NCo,NDs); 00485 set_to_zero (Ns); 00486 for (size_t i=0; i<NCo; ++i) 00487 for (size_t j=0; j<GEs->NN; ++j) 00488 Ns(i,i+j*NCo) = GEs->N(j); 00489 } 00490 00491 inline void USigEpsElem::UpdateState (Vec_t const & dU, Vec_t * F_int) const 00492 { 00493 // stress and strain increments 00494 Vec_t dSe(NDs); 00495 Vec_t dEe(NDs); 00496 for (size_t i=0; i<NDs; ++i) 00497 { 00498 dEe(i) = dU(FirstEQ+i); 00499 dSe(i) = dU(FirstEQ+NDs+i); 00500 } 00501 00502 // displacement increment 00503 Vec_t dUe(NDu); 00504 for (size_t i=0; i<GE->NN; ++i) 00505 { 00506 dUe(i*NDim+0) = dU(Con[i]->Eq("ux")); 00507 dUe(i*NDim+1) = dU(Con[i]->Eq("uy")); if (NDim==3) 00508 dUe(i*NDim+2) = dU(Con[i]->Eq("uz")); 00509 } 00510 00511 // update F_int 00512 if (F_int!=NULL) 00513 { 00514 double detJ, coef; 00515 Mat_t C, B, N, Ns; 00516 Vec_t dsig(NCo), deps(NCo), deps_u(NCo), dF1(NDs), dF2(NDs), dF3(NDu); 00517 CoordMatrix (C); 00518 set_to_zero (dF1); 00519 set_to_zero (dF2); 00520 set_to_zero (dF3); 00521 for (size_t i=0; i<GE->NIP; ++i) 00522 { 00523 // interpolation 00524 Interp (C, GE->IPs[i], B, N, Ns, detJ, coef); 00525 00526 // (external) stress and strain increment 00527 dsig = Ns * dSe; 00528 deps = Ns * dEe; 00529 deps_u = B * dUe; 00530 00531 // (internal) strain increment 00532 Vec_t dsig_i(NCo); 00533 Mdl->SUp.Update (deps, Sta[i], dsig_i); 00534 00535 // (internal) forces 00536 dF1 += (coef) * trans(Ns)*(dsig_i - dsig); 00537 dF2 += (coef) * trans(Ns)*(deps_u - deps); 00538 dF3 += (coef) * trans(B)*dsig; 00539 } 00540 00541 // add to F_int 00542 for (size_t i=0; i<NDs; ++i) 00543 { 00544 (*F_int)(FirstEQ+i) += dF1(i); 00545 (*F_int)(FirstEQ+NDs+i) += dF2(i); 00546 } 00547 for (size_t i=0; i<GE->NN; ++i) 00548 { 00549 (*F_int)(Con[i]->Eq("ux")) += dF3(i*NDim+0); 00550 (*F_int)(Con[i]->Eq("uy")) += dF3(i*NDim+1); if (NDim==3) 00551 (*F_int)(Con[i]->Eq("uz")) += dF3(i*NDim+2); 00552 } 00553 } 00554 } 00555 00556 inline void USigEpsElem::StateKeys (Array<String> & Keys) const 00557 { 00558 Keys = EquilibState::Keys; 00559 for (size_t i=0; i<Mdl->NIvs; ++i) Keys.Push (Mdl->IvNames[i]); 00560 } 00561 00562 inline void USigEpsElem::StateAtIP (SDPair & KeysVals, int IdxIP) const 00563 { 00564 Vec_t const & sig = static_cast<EquilibState const *>(Sta[IdxIP])->Sig; 00565 Vec_t const & eps = static_cast<EquilibState const *>(Sta[IdxIP])->Eps; 00566 Vec_t const & ivs = static_cast<EquilibState const *>(Sta[IdxIP])->Ivs; 00567 00568 if (NDim==2) 00569 { 00570 KeysVals.Set("sx sy sz sxy ex ey ez exy pcam qcam ev ed", 00571 sig(0), sig(1), sig(2), sig(3)/Util::SQ2, 00572 eps(0), eps(1), eps(2), eps(3)/Util::SQ2, 00573 Calc_pcam(sig), Calc_qcam(sig), Calc_ev(eps), Calc_ed(eps)); 00574 } 00575 else 00576 { 00577 KeysVals.Set("sx sy sz sxy syz szx ex ey ez exy eyz ezx pcam qcam ev ed", 00578 sig(0), sig(1), sig(2), sig(3)/Util::SQ2, sig(4)/Util::SQ2, sig(5)/Util::SQ2, 00579 eps(0), eps(1), eps(2), eps(3)/Util::SQ2, eps(4)/Util::SQ2, eps(5)/Util::SQ2, 00580 Calc_pcam(sig), Calc_qcam(sig), Calc_ev(eps), Calc_ed(eps)); 00581 } 00582 for (size_t k=0; k<Mdl->NIvs; ++k) KeysVals.Set(Mdl->IvNames[k].CStr(), ivs(k)); 00583 } 00584 00585 00587 00588 00589 // Allocate a new element 00590 Element * USigEpsElemMaker(int NDim, Mesh::Cell const & Cell, Model const * Mdl, Model const * XMdl, SDPair const & Prp, SDPair const & Ini, Array<Node*> const & Nodes) { return new USigEpsElem(NDim,Cell,Mdl,XMdl,Prp,Ini,Nodes); } 00591 00592 // Register element 00593 int USigEpsElemRegister() 00594 { 00595 ElementFactory["USigEps"] = USigEpsElemMaker; 00596 ElementVarKeys["USigEps2D"] = std::make_pair ("ux uy", "fx fy"); 00597 ElementVarKeys["USigEps3D"] = std::make_pair ("ux uy uz", "fx fy fz"); 00598 PROB.Set ("USigEps", (double)PROB.Keys.Size()); 00599 return 0; 00600 } 00601 00602 // Call register 00603 int __USigEpsElem_dummy_int = USigEpsElemRegister(); 00604 00605 }; // namespace FEM 00606 00607 #endif // MECHSYS_FEM_USIGEPSELEM_H