OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
dkt3d.C
Go to the documentation of this file.
1 /*
2  *
3  * ##### ##### ###### ###### ### ###
4  * ## ## ## ## ## ## ## ### ##
5  * ## ## ## ## #### #### ## # ##
6  * ## ## ## ## ## ## ## ##
7  * ## ## ## ## ## ## ## ##
8  * ##### ##### ## ###### ## ##
9  *
10  *
11  * OOFEM : Object Oriented Finite Element Code
12  *
13  * Copyright (C) 1993 - 2013 Borek Patzak
14  *
15  *
16  *
17  * Czech Technical University, Faculty of Civil Engineering,
18  * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19  *
20  * This library is free software; you can redistribute it and/or
21  * modify it under the terms of the GNU Lesser General Public
22  * License as published by the Free Software Foundation; either
23  * version 2.1 of the License, or (at your option) any later version.
24  *
25  * This program is distributed in the hope that it will be useful,
26  * but WITHOUT ANY WARRANTY; without even the implied warranty of
27  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28  * Lesser General Public License for more details.
29  *
30  * You should have received a copy of the GNU Lesser General Public
31  * License along with this library; if not, write to the Free Software
32  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33  */
34 
35 #include "../sm/Elements/Plates/dkt3d.h"
36 #include "../sm/Materials/structuralms.h"
37 #include "fei2dtrlin.h"
38 #include "node.h"
39 #include "load.h"
40 #include "mathfem.h"
41 #include "gaussintegrationrule.h"
42 #include "gausspoint.h"
43 #include "classfactory.h"
44 
45 #include <cstdlib>
46 
47 namespace oofem {
48 REGISTER_Element(DKTPlate3d);
49 
50 DKTPlate3d :: DKTPlate3d(int n, Domain *aDomain) : DKTPlate(n, aDomain)
51 {
52 }
53 
54 
55 void
57 // Returns global coordinates given in global vector
58 // transformed into local coordinate system of the
59 // receiver
60 {
61  FloatArray offset;
62  // test the parametr
63  if ( global.giveSize() != 3 ) {
64  OOFEM_ERROR("cannot transform coordinates - size mismatch");
65  exit(1);
66  }
67 
68  // first ensure that receiver's GtoLRotationMatrix[3,3] is defined
70 
71  offset = global;
72  offset.subtract( * this->giveNode(1)->giveCoordinates() );
73  answer.beProductOf(GtoLRotationMatrix, offset);
74 }
75 
76 
77 void
78 DKTPlate3d :: giveNodeCoordinates(double &x1, double &x2, double &x3,
79  double &y1, double &y2, double &y3,
80  double &z1, double &z2, double &z3)
81 {
82  FloatArray nc1(3), nc2(3), nc3(3);
83 
84  this->giveLocalCoordinates( nc1, * ( this->giveNode(1)->giveCoordinates() ) );
85  this->giveLocalCoordinates( nc2, * ( this->giveNode(2)->giveCoordinates() ) );
86  this->giveLocalCoordinates( nc3, * ( this->giveNode(3)->giveCoordinates() ) );
87 
88  x1 = nc1.at(1);
89  x2 = nc2.at(1);
90  x3 = nc3.at(1);
91 
92  y1 = nc1.at(2);
93  y2 = nc2.at(2);
94  y3 = nc3.at(2);
95 
96  z1 = nc1.at(3);
97  z2 = nc2.at(3);
98  z3 = nc3.at(3);
99 
100 }
101 
102 
103 void
105 {
106  answer = {D_u, D_v, D_w, R_u, R_v, R_w};
107 }
108 
109 
110 const FloatMatrix *
112 // Returns the rotation matrix of the receiver of the size [3,3]
113 // coords(local) = T * coords(global)
114 //
115 // local coordinate (described by vector triplet e1',e2',e3') is defined as follows:
116 //
117 // e1' : [N2-N1] Ni - means i - th node
118 // help : [N3-N1]
119 // e3' : e1' x help
120 // e2' : e3' x e1'
121 {
122  if ( !GtoLRotationMatrix.isNotEmpty() ) {
123  FloatArray e1, e2, e3, help;
124 
125  // compute e1' = [N2-N1] and help = [N3-N1]
126  e1.beDifferenceOf(*this->giveNode(2)->giveCoordinates(), *this->giveNode(1)->giveCoordinates());
127  help.beDifferenceOf(*this->giveNode(3)->giveCoordinates(), *this->giveNode(1)->giveCoordinates());
128 
129  // let us normalize e1'
130  e1.normalize();
131 
132  // compute e3' : vector product of e1' x help
133  e3.beVectorProductOf(e1, help);
134  // let us normalize
135  e3.normalize();
136 
137  // now from e3' x e1' compute e2'
138  e2.beVectorProductOf(e3, e1);
139 
140  //
142 
143  for ( int i = 1; i <= 3; i++ ) {
144  GtoLRotationMatrix.at(1, i) = e1.at(i);
145  GtoLRotationMatrix.at(2, i) = e2.at(i);
146  GtoLRotationMatrix.at(3, i) = e3.at(i);
147  }
148  }
149 
150  return &GtoLRotationMatrix;
151 }
152 
153 
154 bool
156 // Returns the rotation matrix of the receiver of the size [9,18]
157 // r(local) = T * r(global)
158 // for one node (r written transposed): {w,r1,r2} = T * {u,v,w,r1,r2,r3}
159 {
161 
162  answer.resize(9, 18);
163  answer.zero();
164 
165  for ( int i = 1; i <= 3; i++ ) {
166  answer.at(1, i) = answer.at(1 + 3, i + 6) = answer.at(1 + 6, i + 12) = GtoLRotationMatrix.at(3, i);
167  answer.at(2, i + 3) = answer.at(2 + 3, i + 3 + 6) = answer.at(2 + 6, i + 3 + 12) = GtoLRotationMatrix.at(1, i);
168  answer.at(3, i + 3) = answer.at(3 + 3, i + 3 + 6) = answer.at(3 + 6, i + 3 + 12) = GtoLRotationMatrix.at(2, i);
169  }
170 
171  return 1;
172 }
173 
174 void
176 // returns characteristic tensor of the receiver at given gp and tStep
177 // strain vector = (Kappa_x, Kappa_y, Kappa_xy, Gamma_zx, Gamma_zy)
178 {
179  FloatArray charVect;
181 
182  answer.resize(3, 3);
183  answer.zero();
184 
185  if ( ( type == LocalForceTensor ) || ( type == GlobalForceTensor ) ) {
186  //this->computeStressVector(charVect, gp, tStep);
187  charVect = ms->giveStressVector();
188 
189  answer.at(1, 3) = charVect.at(4);
190  answer.at(3, 1) = charVect.at(4);
191  answer.at(2, 3) = charVect.at(5);
192  answer.at(3, 2) = charVect.at(5);
193  } else if ( ( type == LocalMomentTensor ) || ( type == GlobalMomentTensor ) ) {
194  //this->computeStressVector(charVect, gp, tStep);
195  charVect = ms->giveStressVector();
196 
197  answer.at(1, 1) = charVect.at(1);
198  answer.at(2, 2) = charVect.at(2);
199  answer.at(1, 2) = charVect.at(3);
200  answer.at(2, 1) = charVect.at(3);
201  } else if ( ( type == LocalStrainTensor ) || ( type == GlobalStrainTensor ) ) {
202  //this->computeStrainVector(charVect, gp, tStep);
203  charVect = ms->giveStrainVector();
204 
205  answer.at(1, 3) = charVect.at(4) / 2.;
206  answer.at(3, 1) = charVect.at(4) / 2.;
207  answer.at(2, 3) = charVect.at(5) / 2.;
208  answer.at(3, 2) = charVect.at(5) / 2.;
209  } else if ( ( type == LocalCurvatureTensor ) || ( type == GlobalCurvatureTensor ) ) {
210  //this->computeStrainVector(charVect, gp, tStep);
211  charVect = ms->giveStrainVector();
212 
213  answer.at(1, 1) = charVect.at(1);
214  answer.at(2, 2) = charVect.at(2);
215  answer.at(1, 2) = charVect.at(3) / 2.;
216  answer.at(2, 1) = charVect.at(3) / 2.;
217  } else {
218  OOFEM_ERROR("unsupported tensor mode");
219  exit(1);
220  }
221 
222  if ( ( type == GlobalForceTensor ) || ( type == GlobalMomentTensor ) ||
223  ( type == GlobalStrainTensor ) || ( type == GlobalCurvatureTensor ) ) {
226  }
227 }
228 
229 
230 int
232 {
233  FloatMatrix globTensor;
234  CharTensor cht;
235 
236  answer.resize(6);
237 
238  if ( type == IST_CurvatureTensor || type == IST_ShellStrainTensor ) {
239  if ( type == IST_CurvatureTensor ) {
240  cht = GlobalCurvatureTensor;
241  } else {
242  cht = GlobalStrainTensor;
243  }
244 
245  this->giveCharacteristicTensor(globTensor, cht, gp, tStep);
246 
247  answer.at(1) = globTensor.at(1, 1); //xx
248  answer.at(2) = globTensor.at(2, 2); //yy
249  answer.at(3) = globTensor.at(3, 3); //zz
250  answer.at(4) = 2*globTensor.at(2, 3); //yz
251  answer.at(5) = 2*globTensor.at(1, 3); //xz
252  answer.at(6) = 2*globTensor.at(1, 2); //xy
253 
254  return 1;
255  } else if ( type == IST_ShellMomentTensor || type == IST_ShellForceTensor ) {
256  if ( type == IST_ShellMomentTensor ) {
257  cht = GlobalMomentTensor;
258  } else {
259  cht = GlobalForceTensor;
260  }
261 
262  this->giveCharacteristicTensor(globTensor, cht, gp, tStep);
263 
264  answer.at(1) = globTensor.at(1, 1); //xx
265  answer.at(2) = globTensor.at(2, 2); //yy
266  answer.at(3) = globTensor.at(3, 3); //zz
267  answer.at(4) = globTensor.at(2, 3); //yz
268  answer.at(5) = globTensor.at(1, 3); //xz
269  answer.at(6) = globTensor.at(1, 2); //xy
270 
271  return 1;
272  } else {
273  return NLStructuralElement :: giveIPValue(answer, gp, type, tStep);
274  }
275 }
276 
277 int
279 // Returns the rotation matrix of the receiver of the size [6,6]
280 // f(local) = T * f(global)
281 {
283 
284  answer.resize(6, 6);
285  answer.zero();
286 
287  for ( int i = 1; i <= 3; i++ ) {
288  answer.at(1, i) = answer.at(4, i + 3) = GtoLRotationMatrix.at(1, i);
289  answer.at(2, i) = answer.at(5, i + 3) = GtoLRotationMatrix.at(2, i);
290  answer.at(3, i) = answer.at(6, i + 3) = GtoLRotationMatrix.at(3, i);
291  }
292 
293  return 1;
294 }
295 
296 void
298 {
299  FloatMatrix ne;
300  this->computeNmatrixAt(sgp->giveNaturalCoordinates(), ne);
301 
302  answer.resize(6, 18);
303  answer.zero();
304  int ri[] = {
305  2, 3, 4
306  };
307  int ci[] = {
308  2, 3, 4, 8, 9, 10, 14, 15, 16
309  };
310 
311  for ( int i = 0; i < 3; i++ ) {
312  for ( int j = 0; j < 9; j++ ) {
313  answer(ri [ i ], ci [ j ]) = ne(i, j);
314  }
315  }
316 }
317 
318 void
320 {
321  answer.resize(18);
322  answer.zero();
323  if ( iSurf == 1 ) {
324  answer.at(3) = 1; // node 1
325  answer.at(4) = 2;
326  answer.at(5) = 3;
327 
328  answer.at(9) = 4; // node 2
329  answer.at(10) = 5;
330  answer.at(11) = 6;
331 
332  answer.at(15) = 7; // node 3
333  answer.at(16) = 8;
334  answer.at(17) = 9;
335  } else {
336  OOFEM_ERROR("wrong surface number");
337  }
338 }
339 
342 {
343  IntegrationRule *iRule = new GaussIntegrationRule(1, this, 1, 1);
344  int npoints = iRule->getRequiredNumberOfIntegrationPoints(_Triangle, approxOrder);
345  iRule->SetUpPointsOnTriangle(npoints, _Unknown);
346  return iRule;
347 }
348 
349 double
351 {
352  return this->computeVolumeAround(gp);
353 }
354 
355 
356 int
358 {
359  return 0;
360 }
361 
362 void
364 // Performs end-of-step operations.
365 {
366  FloatArray v;
367 
368  fprintf( file, "element %d (%8d) :\n", this->giveLabel(), this->giveNumber() );
369 
370  for ( int i = 0; i < (int)integrationRulesArray.size(); i++ ) {
371  for ( GaussPoint *gp: *integrationRulesArray [ i ] ) {
372 
373  fprintf( file, " GP %2d.%-2d :", i + 1, gp->giveNumber() );
374 
375  this->giveIPValue(v, gp, IST_ShellStrainTensor, tStep);
376  fprintf(file, " strains ");
377  // eps_x, eps_y, eps_z, eps_yz, eps_xz, eps_xy (global)
378  for ( auto &val : v ) fprintf(file, " %.4e", val);
379 
380  this->giveIPValue(v, gp, IST_CurvatureTensor, tStep);
381  fprintf(file, "\n curvatures ");
382  for ( auto &val : v ) fprintf(file, " %.4e", val);
383 
384  // Forces - Moments
385  this->giveIPValue(v, gp, IST_ShellForceTensor, tStep);
386  fprintf(file, "\n stresses ");
387  for ( auto &val : v ) fprintf(file, " %.4e", val);
388 
389  this->giveIPValue(v, gp, IST_ShellMomentTensor, tStep);
390  fprintf(file, "\n moments ");
391  for ( auto &val : v ) fprintf(file, " %.4e", val);
392 
393  fprintf(file, "\n");
394  }
395  }
396 }
397 
398 
399 void
401 {
402  /*
403  * provides dof mapping of local edge dofs (only nonzero are taken into account)
404  * to global element dofs
405  */
406 
407  answer.resize(12);
408  answer.zero();
409  if ( iEdge == 1 ) { // edge between nodes 1,2
410  answer.at(3) = 3;
411  answer.at(4) = 4;
412  answer.at(5) = 5;
413  answer.at(9) = 9;
414  answer.at(10) = 10;
415  answer.at(11) = 11;
416  } else if ( iEdge == 2 ) { // edge between nodes 2 3
417  answer.at(3) = 9;
418  answer.at(4) = 10;
419  answer.at(5) = 11;
420  answer.at(9) = 15;
421  answer.at(10) = 16;
422  answer.at(11) = 17;
423  } else if ( iEdge == 3 ) { // edge between nodes 3 1
424  answer.at(3) = 15;
425  answer.at(4) = 16;
426  answer.at(5) = 17;
427  answer.at(9) = 4;
428  answer.at(10) = 5;
429  answer.at(11) = 6;
430  } else {
431  OOFEM_ERROR("wrong edge number");
432  }
433 }
434 
435 double
437 {
439  return detJ *gp->giveWeight();
440 }
441 
442 
443 int
445 {
446  // returns transformation matrix from
447  // edge local coordinate system
448  // to element local c.s
449  // (same as global c.s in this case)
450  //
451  // i.e. f(element local) = T * f(edge local)
452  //
453  double dx, dy, length;
454  IntArray edgeNodes;
455  Node *nodeA, *nodeB;
456 
457  answer.resize(3, 3);
458  answer.zero();
459 
460  this->interp_lin.computeLocalEdgeMapping(edgeNodes, iEdge);
461 
462  nodeA = this->giveNode( edgeNodes.at(1) );
463  nodeB = this->giveNode( edgeNodes.at(2) );
464 
465  FloatArray cb(3), ca(3);
466  this->giveLocalCoordinates (ca, * (nodeA->giveCoordinates() ) );
467  this->giveLocalCoordinates (cb, * (nodeB->giveCoordinates() ) );
468 
469  dx = cb.at(1) - ca.at(1);
470  dy = cb.at(2) - ca.at(2);
471  length = sqrt(dx * dx + dy * dy);
472 
473  answer.at(1, 1) = 1.0;
474  answer.at(2, 2) = dx / length;
475  answer.at(2, 3) = -dy / length;
476  answer.at(3, 2) = dy / length;
477  answer.at(3, 3) = dx / length;
478 
479  return 1;
480 }
481 
482 bool
484 //converts global coordinates to local planar area coordinates,
485 //does not return a coordinate in the thickness direction, but
486 //does check that the point is in the element thickness
487 {
488  // rotate the input point Coordinate System into the element CS
489  FloatArray inputCoords_ElCS;
490  std::vector< FloatArray > lc(3);
491  FloatArray llc;
492  this->giveLocalCoordinates( inputCoords_ElCS, const_cast< FloatArray & >(coords) );
493  for ( int _i = 0; _i < 3; _i++ ) {
494  this->giveLocalCoordinates( lc [ _i ], * this->giveNode(_i + 1)->giveCoordinates() );
495  }
496  FEI2dTrLin _interp(1, 2);
497  bool inplane = _interp.global2local(llc, inputCoords_ElCS, FEIVertexListGeometryWrapper(lc)) > 0;
498  answer.resize(2);
499  answer.at(1) = inputCoords_ElCS.at(1);
500  answer.at(2) = inputCoords_ElCS.at(2);
501  GaussPoint _gp(NULL, 1, answer, 2.0, _2dPlate);
502  // now check if the third local coordinate is within the thickness of element
503  bool outofplane = ( fabs( inputCoords_ElCS.at(3) ) <= this->giveCrossSection()->give(CS_Thickness, & _gp) / 2. );
504 
505  return inplane && outofplane;
506 }
507 
508 
509 int
511 {
512  double l1 = lcoords.at(1);
513  double l2 = lcoords.at(2);
514  double l3 = 1. - l2 - l1;
515 
516  answer.resize(3);
517  for ( int _i = 1; _i <= 3; _i++ ) {
518  answer.at(_i) = l1 * this->giveNode(1)->giveCoordinate(_i) + l2 *this->giveNode(2)->giveCoordinate(_i) + l3 *this->giveNode(3)->giveCoordinate(_i);
519  }
520  return true;
521 }
522 
523 void
525 // Computes numerically the load vector of the receiver due to the body loads, at tStep.
526 // load is assumed to be in global cs.
527 // load vector is then transformed to coordinate system in each node.
528 // (should be global coordinate system, but there may be defined
529 // different coordinate system in each node)
530 {
531  double dens, dV, load;
532  FloatArray force;
533  FloatMatrix T;
534 
535  if ( ( forLoad->giveBCGeoType() != BodyLoadBGT ) || ( forLoad->giveBCValType() != ForceLoadBVT ) ) {
536  OOFEM_ERROR("unknown load type");
537  }
538 
539  GaussIntegrationRule irule(1, this, 1, 5);
540  irule.SetUpPointsOnTriangle(1, _2dPlate);
541 
542  // note: force is assumed to be in global coordinate system.
543  forLoad->computeComponentArrayAt(force, tStep, mode);
544 
545  if ( force.giveSize() ) {
546  GaussPoint *gp = irule.getIntegrationPoint(0);
547 
548  dens = this->giveStructuralCrossSection()->give('d', gp); // constant density assumed
549  dV = this->computeVolumeAround(gp) * this->giveCrossSection()->give(CS_Thickness, gp); // constant thickness assumed
550 
551  answer.resize(18);
552  answer.zero();
553 
554  load = force.at(1) * dens * dV / 3.0;
555  answer.at(1) = load;
556  answer.at(7) = load;
557  answer.at(13) = load;
558 
559  load = force.at(2) * dens * dV / 3.0;
560  answer.at(2) = load;
561  answer.at(8) = load;
562  answer.at(14) = load;
563 
564  load = force.at(3) * dens * dV / 3.0;
565  answer.at(3) = load;
566  answer.at(9) = load;
567  answer.at(15) = load;
568 
569  // transform result from global cs to local element cs.
570  if ( this->computeGtoLRotationMatrix(T) ) {
571  answer.rotatedWith(T, 'n');
572  }
573  } else {
574  answer.clear(); // nil resultant
575  }
576 }
577 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
This class implements an triangular Discrete Kirchhoff Theory (DKT) element.
Definition: dkt.h:67
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
void giveLocalCoordinates(FloatArray &answer, FloatArray &global)
Definition: dkt3d.C:56
void subtract(const FloatArray &src)
Subtracts array src to receiver.
Definition: floatarray.C:258
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
void beVectorProductOf(const FloatArray &v1, const FloatArray &v2)
Computes vector product (or cross product) of vectors given as parameters, , and stores the result in...
Definition: floatarray.C:415
Class and object Domain.
Definition: domain.h:115
Wrapper around cell with vertex coordinates stored in FloatArray**.
Definition: feinterpol.h:115
virtual void computeLocalEdgeMapping(IntArray &edgeNodes, int iedge)
Definition: fei2dtrlin.C:230
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: dkt.C:389
void zero()
Sets all component to zero.
Definition: intarray.C:52
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
This class implements a structural material status information.
Definition: structuralms.h:65
void giveCharacteristicTensor(FloatMatrix &answer, CharTensor type, GaussPoint *gp, TimeStep *tStep)
Definition: dkt3d.C:175
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: dkt3d.C:231
virtual double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf)
Computes volume related to integration point on local surface.
Definition: dkt3d.C:350
virtual double edgeGiveTransformationJacobian(int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the edge Jacobian of transformation between local and global coordinates.
Definition: feinterpol2d.C:175
virtual IntegrationRule * GetSurfaceIntegrationRule(int iSurf)
Definition: dkt3d.C:341
virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
Definition: dkt3d.C:400
virtual void computeSurfaceNMatrixAt(FloatMatrix &answer, int iSurf, GaussPoint *gp)
Definition: dkt3d.C:297
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Definition: dkt3d.C:363
virtual void computeComponentArrayAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Computes boundary condition value - its components values at given time.
Definition: load.C:82
virtual double giveCoordinate(int i)
Definition: node.C:82
virtual void giveNodeCoordinates(double &x1, double &x2, double &x3, double &y1, double &y2, double &y3, double &z1, double &z2, double &z3)
Definition: dkt3d.C:78
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
virtual int computeLoadLSToLRotationMatrix(FloatMatrix &answer, int iSurf, GaussPoint *gp)
Returns transformation matrix from local surface c.s to element local coordinate system of load vecto...
Definition: dkt3d.C:357
void rotatedWith(FloatMatrix &r, char mode)
Returns the receiver a rotated according the change-of-base matrix r.
Definition: floatarray.C:799
Abstract base class representing integration rule.
Distributed body load.
Definition: bcgeomtype.h:43
void beDifferenceOf(const FloatArray &a, const FloatArray &b)
Sets receiver to be a - b.
Definition: floatarray.C:341
virtual void giveSurfaceDofMapping(IntArray &answer, int iSurf) const
Assembles surface dof mapping mask, which provides mapping between surface local DOFs and "global" el...
Definition: dkt3d.C:319
bool isNotEmpty() const
Tests for empty matrix.
Definition: floatmatrix.h:162
Class representing a 2d triangular linear interpolation based on area coordinates.
Definition: fei2dtrlin.h:44
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
FloatMatrix GtoLRotationMatrix
Transformation Matrix form GtoL(3,3) is stored at the element level for computation efficiency...
Definition: dkt3d.h:77
virtual double giveWeight()
Returns integration weight of receiver.
Definition: gausspoint.h:181
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
virtual int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp)
Returns transformation matrix from local edge c.s to element local coordinate system of load vector c...
Definition: dkt3d.C:444
void rotatedWith(const FloatMatrix &r, char mode= 'n')
Returns the receiver &#39;a&#39; transformed using give transformation matrix r.
Definition: floatmatrix.C:1557
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
Definition: dkt3d.C:278
virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode)
Computes the load vector due to body load acting on receiver, at given time step. ...
Definition: dkt3d.C:524
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
Definition: floatarray.C:676
virtual int global2local(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Default implementation using Newton&#39;s method to find the local coordinates.
Definition: fei2dtrlin.C:103
virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
Computes volume related to integration point on local edge.
Definition: dkt3d.C:436
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
void resize(int n)
Checks size of receiver towards requested bounds.
Definition: intarray.C:124
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
Definition: dkt.C:293
Class representing vector of real numbers.
Definition: floatarray.h:82
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
virtual double give(CrossSectionProperty a, GaussPoint *gp)
Returns the value of cross section property at given point.
Definition: crosssection.C:151
GaussPoint * getIntegrationPoint(int n)
Access particular integration point of receiver.
virtual const FloatMatrix * computeGtoLRotationMatrix()
Definition: dkt3d.C:111
const FloatArray & giveStressVector() const
Returns the const pointer to receiver&#39;s stress vector.
Definition: structuralms.h:107
IntegrationPointStatus * giveMaterialStatus()
Returns reference to associated material status (NULL if not defined).
Definition: gausspoint.h:205
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1358
virtual int getRequiredNumberOfIntegrationPoints(integrationDomain dType, int approxOrder)
Abstract service.
static FEI2dTrLin interp_lin
Element geometry approximation.
Definition: dkt.h:74
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
virtual bcGeomType giveBCGeoType() const
Returns geometry character of boundary condition.
std::vector< std::unique_ptr< IntegrationRule > > integrationRulesArray
List of integration rules of receiver (each integration rule contains associated integration points a...
Definition: element.h:170
virtual int SetUpPointsOnTriangle(int, MaterialMode mode)
Sets up receiver&#39;s integration points on triangular (area coords) integration domain.
virtual FloatArray * giveCoordinates()
Definition: node.h:114
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
virtual bool computeLocalCoordinates(FloatArray &answer, const FloatArray &gcoords)
Computes the element local coordinates from given global coordinates.
Definition: dkt3d.C:483
Load is base abstract class for all loads.
Definition: load.h:61
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
int giveLabel() const
Definition: element.h:1055
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual bcValType giveBCValType() const
Returns receiver load type.
Class implementing node in finite element mesh.
Definition: node.h:87
int giveNumber() const
Definition: femcmpnn.h:107
double normalize()
Normalizes receiver.
Definition: floatarray.C:828
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element&#39;s local coordinates.
Definition: dkt3d.C:510
const FloatArray & giveStrainVector() const
Returns the const pointer to receiver&#39;s strain vector.
Definition: structuralms.h:105
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
Definition: dkt3d.C:104
DKTPlate3d(int n, Domain *d)
Definition: dkt3d.C:50
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
virtual int SetUpPointsOnTriangle(int nPoints, MaterialMode mode)
Sets up receiver&#39;s integration points on triangular (area coords) integration domain.
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:28 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011