OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
cct3d.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/Shells/cct3d.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(CCTPlate3d);
49 
50 CCTPlate3d :: CCTPlate3d(int n, Domain *aDomain) : CCTPlate(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 CCTPlate3d :: 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 bool
112 //converts global coordinates to local planar area coordinates,
113 //does not return a coordinate in the thickness direction, but
114 //does check that the point is in the element thickness
115 {
116  // rotate the input point Coordinate System into the element CS
117  FloatArray inputCoords_ElCS;
118  std::vector< FloatArray > lc(3);
119  FloatArray llc;
120  this->giveLocalCoordinates( inputCoords_ElCS, const_cast< FloatArray & >(coords) );
121  for ( int _i = 0; _i < 3; _i++ ) {
122  this->giveLocalCoordinates( lc [ _i ], * this->giveNode(_i + 1)->giveCoordinates() );
123  }
124  FEI2dTrLin _interp(1, 2);
125  bool inplane = _interp.global2local(llc, inputCoords_ElCS, FEIVertexListGeometryWrapper(lc)) > 0;
126  answer.resize(2);
127  answer.at(1) = inputCoords_ElCS.at(1);
128  answer.at(2) = inputCoords_ElCS.at(2);
129  GaussPoint _gp(NULL, 1, answer, 2.0, _2dPlate);
130  // now check if the third local coordinate is within the thickness of element
131  bool outofplane = ( fabs( inputCoords_ElCS.at(3) ) <= this->giveCrossSection()->give(CS_Thickness, & _gp) / 2. );
132 
133  return inplane && outofplane;
134 }
135 
136 
137 int
139 {
140  double l1 = lcoords.at(1);
141  double l2 = lcoords.at(2);
142  double l3 = 1. - l2 - l1;
143 
144  answer.resize(3);
145  for ( int _i = 1; _i <= 3; _i++ ) {
146  answer.at(_i) = l1 * this->giveNode(1)->giveCoordinate(_i) + l2 *this->giveNode(2)->giveCoordinate(_i) + l3 *this->giveNode(3)->giveCoordinate(_i);
147  }
148  return true;
149 }
150 
151 
152 const FloatMatrix *
154 // Returns the rotation matrix of the receiver of the size [3,3]
155 // coords(local) = T * coords(global)
156 //
157 // local coordinate (described by vector triplet e1',e2',e3') is defined as follows:
158 //
159 // e1' : [N2-N1] Ni - means i - th node
160 // help : [N3-N1]
161 // e3' : e1' x help
162 // e2' : e3' x e1'
163 {
164  if ( !GtoLRotationMatrix.isNotEmpty() ) {
165  FloatArray e1, e2, e3, help;
166 
167  // compute e1' = [N2-N1] and help = [N3-N1]
168  e1.beDifferenceOf(*this->giveNode(2)->giveCoordinates(), *this->giveNode(1)->giveCoordinates());
169  help.beDifferenceOf(*this->giveNode(3)->giveCoordinates(), *this->giveNode(1)->giveCoordinates());
170 
171  // let us normalize e1'
172  e1.normalize();
173 
174  // compute e3' : vector product of e1' x help
175  e3.beVectorProductOf(e1, help);
176  // let us normalize
177  e3.normalize();
178 
179  // now from e3' x e1' compute e2'
180  e2.beVectorProductOf(e3, e1);
181 
182  //
184 
185  for ( int i = 1; i <= 3; i++ ) {
186  GtoLRotationMatrix.at(1, i) = e1.at(i);
187  GtoLRotationMatrix.at(2, i) = e2.at(i);
188  GtoLRotationMatrix.at(3, i) = e3.at(i);
189  }
190  }
191 
192  return &GtoLRotationMatrix;
193 }
194 
195 
196 bool
198 // Returns the rotation matrix of the receiver of the size [9,18]
199 // r(local) = T * r(global)
200 // for one node (r written transposed): {w,r1,r2} = T * {u,v,w,r1,r2,r3}
201 {
203 
204  answer.resize(9, 18);
205  answer.zero();
206 
207  for ( int i = 1; i <= 3; i++ ) {
208  answer.at(1, i) = answer.at(1 + 3, i + 6) = answer.at(1 + 6, i + 12) = GtoLRotationMatrix.at(3, i);
209  answer.at(2, i + 3) = answer.at(2 + 3, i + 3 + 6) = answer.at(2 + 6, i + 3 + 12) = GtoLRotationMatrix.at(1, i);
210  answer.at(3, i + 3) = answer.at(3 + 3, i + 3 + 6) = answer.at(3 + 6, i + 3 + 12) = GtoLRotationMatrix.at(2, i);
211  }
212 
213  return 1;
214 }
215 
216 void
218 // returns characteristic tensor of the receiver at given gp and tStep
219 // strain vector = (Kappa_x, Kappa_y, Kappa_xy, Gamma_zx, Gamma_zy)
220 {
221  FloatArray charVect;
223 
224  answer.resize(3, 3);
225  answer.zero();
226 
227  if ( ( type == LocalForceTensor ) || ( type == GlobalForceTensor ) ) {
228  //this->computeStressVector(charVect, gp, tStep);
229  charVect = ms->giveStressVector();
230 
231  answer.at(1, 3) = charVect.at(4);
232  answer.at(3, 1) = charVect.at(4);
233  answer.at(2, 3) = charVect.at(5);
234  answer.at(3, 2) = charVect.at(5);
235  } else if ( ( type == LocalMomentTensor ) || ( type == GlobalMomentTensor ) ) {
236  //this->computeStressVector(charVect, gp, tStep);
237  charVect = ms->giveStressVector();
238 
239  answer.at(1, 1) = charVect.at(1);
240  answer.at(2, 2) = charVect.at(2);
241  answer.at(1, 2) = charVect.at(3);
242  answer.at(2, 1) = charVect.at(3);
243  } else if ( ( type == LocalStrainTensor ) || ( type == GlobalStrainTensor ) ) {
244  //this->computeStrainVector(charVect, gp, tStep);
245  charVect = ms->giveStrainVector();
246 
247  answer.at(1, 3) = charVect.at(4) / 2.;
248  answer.at(3, 1) = charVect.at(4) / 2.;
249  answer.at(2, 3) = charVect.at(5) / 2.;
250  answer.at(3, 2) = charVect.at(5) / 2.;
251  } else if ( ( type == LocalCurvatureTensor ) || ( type == GlobalCurvatureTensor ) ) {
252  //this->computeStrainVector(charVect, gp, tStep);
253  charVect = ms->giveStrainVector();
254 
255  answer.at(1, 1) = charVect.at(1);
256  answer.at(2, 2) = charVect.at(2);
257  answer.at(1, 2) = charVect.at(3) / 2.;
258  answer.at(2, 1) = charVect.at(3) / 2.;
259  } else {
260  OOFEM_ERROR("unsupported tensor mode");
261  exit(1);
262  }
263 
264  if ( ( type == GlobalForceTensor ) || ( type == GlobalMomentTensor ) ||
265  ( type == GlobalStrainTensor ) || ( type == GlobalCurvatureTensor ) ) {
268  }
269 }
270 
271 
272 int
274 {
275  FloatMatrix globTensor;
276  CharTensor cht;
277 
278  answer.resize(6);
279 
280  if ( type == IST_CurvatureTensor || type == IST_ShellStrainTensor ) {
281  if ( type == IST_CurvatureTensor ) {
282  cht = GlobalCurvatureTensor;
283  } else {
284  cht = GlobalStrainTensor;
285  }
286 
287  this->giveCharacteristicTensor(globTensor, cht, gp, tStep);
288 
289  answer.at(1) = globTensor.at(1, 1); //xx
290  answer.at(2) = globTensor.at(2, 2); //yy
291  answer.at(3) = globTensor.at(3, 3); //zz
292  answer.at(4) = 2 * globTensor.at(2, 3); //yz
293  answer.at(5) = 2 * globTensor.at(1, 3); //xz
294  answer.at(6) = 2 * globTensor.at(1, 2); //xy
295 
296  return 1;
297  } else if ( type == IST_ShellMomentTensor || type == IST_ShellForceTensor ) {
298  if ( type == IST_ShellMomentTensor ) {
299  cht = GlobalMomentTensor;
300  } else {
301  cht = GlobalForceTensor;
302  }
303 
304  this->giveCharacteristicTensor(globTensor, cht, gp, tStep);
305 
306  answer.at(1) = globTensor.at(1, 1); //xx
307  answer.at(2) = globTensor.at(2, 2); //yy
308  answer.at(3) = globTensor.at(3, 3); //zz
309  answer.at(4) = globTensor.at(2, 3); //yz
310  answer.at(5) = globTensor.at(1, 3); //xz
311  answer.at(6) = globTensor.at(1, 2); //xy
312 
313  return 1;
314  } else {
315  return NLStructuralElement :: giveIPValue(answer, gp, type, tStep);
316  }
317 }
318 
319 int
321 // Returns the rotation matrix of the receiver of the size [6,6]
322 // f(local) = T * f(global)
323 {
325 
326  answer.resize(6, 6);
327  answer.zero();
328 
329  for ( int i = 1; i <= 3; i++ ) {
330  answer.at(1, i) = answer.at(4, i + 3) = GtoLRotationMatrix.at(1, i);
331  answer.at(2, i) = answer.at(5, i + 3) = GtoLRotationMatrix.at(2, i);
332  answer.at(3, i) = answer.at(6, i + 3) = GtoLRotationMatrix.at(3, i);
333  }
334 
335  return 1;
336 }
337 
338 void
340 // Computes numerically the load vector of the receiver due to the body loads, at tStep.
341 // load is assumed to be in global cs.
342 // load vector is then transformed to coordinate system in each node.
343 // (should be global coordinate system, but there may be defined
344 // different coordinate system in each node)
345 {
346  double dens, dV, load;
347  FloatArray force;
348  FloatMatrix T;
349 
350  if ( ( forLoad->giveBCGeoType() != BodyLoadBGT ) || ( forLoad->giveBCValType() != ForceLoadBVT ) ) {
351  OOFEM_ERROR("unknown load type");
352  }
353 
354  GaussIntegrationRule irule(1, this, 1, 5);
355  irule.SetUpPointsOnTriangle(1, _2dPlate);
356 
357  // note: force is assumed to be in global coordinate system, 6 components
358  forLoad->computeComponentArrayAt(force, tStep, mode);
359  // get it in local c.s.
361  force.rotatedWith(T, 'n');
362 
363  if ( force.giveSize() ) {
364  GaussPoint *gp = irule.getIntegrationPoint(0);
365 
366  dens = this->giveStructuralCrossSection()->give('d', gp); // constant density assumed
367  dV = this->computeVolumeAround(gp) * this->giveCrossSection()->give(CS_Thickness, gp); // constant thickness assumed
368 
369  answer.resize(9);
370  answer.zero();
371 
372  load = force.at(3) * dens * dV / 3.0;
373  answer.at(1) = load;
374  answer.at(4) = load;
375  answer.at(7) = load;
376 
377  load = force.at(4) * dens * dV / 3.0;
378  answer.at(2) = load;
379  answer.at(5) = load;
380  answer.at(8) = load;
381 
382  load = force.at(5) * dens * dV / 3.0;
383  answer.at(3) = load;
384  answer.at(6) = load;
385  answer.at(9) = load;
386 
387  // transform result from global cs to local element cs.
388  //if ( this->computeGtoLRotationMatrix(T) ) {
389  // answer.rotatedWith(T, 'n');
390  //}
391  } else {
392  answer.clear(); // nil resultant
393  }
394 }
395 
396 void
398 {
399  FloatMatrix ne;
400  this->computeNmatrixAt(sgp->giveNaturalCoordinates(), ne);
401 
402  answer.resize(6, 18);
403  answer.zero();
404  int ri[] = {
405  2, 3, 4
406  };
407  int ci[] = {
408  2, 3, 4, 8, 9, 10, 14, 15, 16
409  };
410 
411  for ( int i = 0; i < 3; i++ ) {
412  for ( int j = 0; j < 9; j++ ) {
413  answer(ri [ i ], ci [ j ]) = ne(i, j);
414  }
415  }
416 }
417 
418 void
420 {
421  answer.resize(18);
422  answer.zero();
423  if ( iSurf == 1 ) {
424  answer.at(3) = 1; // node 1
425  answer.at(4) = 2;
426  answer.at(5) = 3;
427 
428  answer.at(9) = 4; // node 2
429  answer.at(10) = 5;
430  answer.at(11) = 6;
431 
432  answer.at(15) = 7; // node 3
433  answer.at(16) = 8;
434  answer.at(17) = 9;
435  } else {
436  OOFEM_ERROR("wrong surface number");
437  }
438 }
439 
442 {
443  IntegrationRule *iRule = new GaussIntegrationRule(1, this, 1, 1);
444  int npoints = iRule->getRequiredNumberOfIntegrationPoints(_Triangle, approxOrder);
445  iRule->SetUpPointsOnTriangle(npoints, _Unknown);
446  return iRule;
447 }
448 
449 double
451 {
452  return this->computeVolumeAround(gp);
453 }
454 
455 
456 int
458 {
459  return 0;
460 }
461 
462 void
464 // Performs end-of-step operations.
465 {
466  FloatArray v;
467 
468  fprintf( file, "element %d (%8d) :\n", this->giveLabel(), this->giveNumber() );
469 
470  for ( int i = 0; i < (int)integrationRulesArray.size(); i++ ) {
471  for ( GaussPoint *gp: *integrationRulesArray [ i ] ) {
472 
473  fprintf( file, " GP %2d.%-2d :", i + 1, gp->giveNumber() );
474 
475  this->giveIPValue(v, gp, IST_ShellStrainTensor, tStep);
476  fprintf(file, " strains ");
477  for ( auto &val : v ) fprintf(file, " %.4e", val);
478 
479  this->giveIPValue(v, gp, IST_CurvatureTensor, tStep);
480  fprintf(file, "\n curvatures ");
481  for ( auto &val : v ) fprintf(file, " %.4e", val);
482 
483  // Forces - Moments
484  this->giveIPValue(v, gp, IST_ShellForceTensor, tStep);
485  fprintf(file, "\n stresses ");
486  for ( auto &val : v ) fprintf(file, " %.4e", val);
487 
488  this->giveIPValue(v, gp, IST_ShellMomentTensor, tStep);
489  fprintf(file, "\n moments ");
490  for ( auto &val : v ) fprintf(file, " %.4e", val);
491 
492  fprintf(file, "\n");
493  }
494  }
495 }
496 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
void subtract(const FloatArray &src)
Subtracts array src to receiver.
Definition: floatarray.C:258
void giveLocalCoordinates(FloatArray &answer, FloatArray &global)
Definition: cct3d.C:56
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
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
Definition: cct3d.C:104
Wrapper around cell with vertex coordinates stored in FloatArray**.
Definition: feinterpol.h:115
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
virtual double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf)
Computes volume related to integration point on local surface.
Definition: cct3d.C:450
This class implements a structural material status information.
Definition: structuralms.h:65
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element&#39;s local coordinates.
Definition: cct3d.C:138
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
Definition: cct3d.C:320
virtual void computeComponentArrayAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Computes boundary condition value - its components values at given time.
Definition: load.C:82
void giveCharacteristicTensor(FloatMatrix &answer, CharTensor type, GaussPoint *gp, TimeStep *tStep)
Definition: cct3d.C:217
virtual double giveCoordinate(int i)
Definition: node.C:82
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
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
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: cct3d.C:457
void beDifferenceOf(const FloatArray &a, const FloatArray &b)
Sets receiver to be a - b.
Definition: floatarray.C:341
CCTPlate3d(int n, Domain *d)
Definition: cct3d.C:50
This class implements an triangular three-node plate CCT finite element.
Definition: cct.h:58
bool isNotEmpty() const
Tests for empty matrix.
Definition: floatmatrix.h:162
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Definition: cct3d.C:463
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: cct3d.C:273
Class representing a 2d triangular linear interpolation based on area coordinates.
Definition: fei2dtrlin.h:44
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
Definition: cct.C:213
FloatMatrix GtoLRotationMatrix
Transformation Matrix form GtoL(3,3) is stored at the element level for computation efficiency...
Definition: cct3d.h:79
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
void rotatedWith(const FloatMatrix &r, char mode= 'n')
Returns the receiver &#39;a&#39; transformed using give transformation matrix r.
Definition: floatmatrix.C:1557
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
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
Class representing vector of real numbers.
Definition: floatarray.h:82
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: cct3d.C:339
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.
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.
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
virtual bcGeomType giveBCGeoType() const
Returns geometry character of boundary condition.
virtual void giveNodeCoordinates(double &x1, double &x2, double &x3, double &y1, double &y2, double &y3, double &z1, double &z2, double &z3)
Definition: cct3d.C:78
virtual const FloatMatrix * computeGtoLRotationMatrix()
Definition: cct3d.C:153
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.
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
virtual void giveSurfaceDofMapping(IntArray &answer, int iSurf) const
Assembles surface dof mapping mask, which provides mapping between surface local DOFs and "global" el...
Definition: cct3d.C:419
Load is base abstract class for all loads.
Definition: load.h:61
virtual void computeSurfaceNMatrixAt(FloatMatrix &answer, int iSurf, GaussPoint *gp)
Definition: cct3d.C:397
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 double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: cct.C:350
virtual bcValType giveBCValType() const
Returns receiver load type.
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
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 bool computeLocalCoordinates(FloatArray &answer, const FloatArray &gcoords)
Computes the element local coordinates from given global coordinates.
Definition: cct3d.C:111
virtual IntegrationRule * GetSurfaceIntegrationRule(int iSurf)
Definition: cct3d.C:441
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:27 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011