OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
libeam3d.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/Beams/libeam3d.h"
36 #include "../sm/Materials/structuralms.h"
37 #include "node.h"
38 #include "material.h"
39 #include "crosssection.h"
40 #include "gausspoint.h"
41 #include "gaussintegrationrule.h"
42 #include "floatmatrix.h"
43 #include "intarray.h"
44 #include "floatarray.h"
45 #include "mathfem.h"
46 #include "classfactory.h"
47 
48 #ifdef __OOFEG
49  #include "oofeggraphiccontext.h"
50 #endif
51 
52 namespace oofem {
53 REGISTER_Element(LIBeam3d);
54 
55 LIBeam3d :: LIBeam3d(int n, Domain *aDomain) : StructuralElement(n, aDomain)
56  // Constructor.
57 {
58  numberOfDofMans = 2;
59  referenceNode = 0;
60  length = 0.;
61 }
62 
63 
64 void
66 // Returns the strain matrix of the receiver.
67 // eeps = {\eps_x, \gamma_xz, \gamma_xy, \der{phi_x}{x}, \kappa_y, \kappa_z}^T
68 {
69  double l, ksi, n1, n2, n1x, n2x;
70 
71  l = this->computeLength();
72  ksi = gp->giveNaturalCoordinate(1);
73 
74  answer.resize(6, 12);
75  answer.zero();
76 
77  n1 = 0.5 * ( 1 - ksi );
78  n2 = 0.5 * ( 1. + ksi );
79  n1x = -1.0 / l;
80  n2x = 1.0 / l;
81 
82  answer.at(1, 1) = -1. / l;
83  answer.at(1, 7) = 1. / l;
84 
85  answer.at(2, 3) = n1x;
86  answer.at(2, 5) = n1;
87  answer.at(2, 9) = n2x;
88  answer.at(2, 11) = n2;
89 
90  answer.at(3, 2) = n1x;
91  answer.at(3, 6) = -n1;
92  answer.at(3, 8) = n2x;
93  answer.at(3, 12) = -n2;
94 
95  answer.at(4, 4) = -1. / l;
96  answer.at(4, 10) = 1. / l;
97 
98  answer.at(5, 5) = n1x;
99  answer.at(5, 11) = n2x;
100 
101  answer.at(6, 6) = n1x;
102  answer.at(6, 12) = n2x;
103 }
104 
105 
106 void
108 // Sets up the array of Gauss Points of the receiver.
109 {
110  if ( integrationRulesArray.size() == 0 ) {
111  integrationRulesArray.resize( 1 );
112  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 2) );
114  }
115 }
116 
117 
118 void
120 {
121  this->giveStructuralCrossSection()->give3dBeamStiffMtrx(answer, rMode, gp, tStep);
122 }
123 
124 
125 void
127 {
128  this->giveStructuralCrossSection()->giveGeneralizedStress_Beam3d(answer, gp, strain, tStep);
129 }
130 
131 
132 void
134 // Returns the lumped mass matrix of the receiver. This expression is
135 // valid in both local and global axes.
136 {
137  GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
138  double density = this->giveStructuralCrossSection()->give('d', gp);
139  double halfMass = density * this->giveCrossSection()->give(CS_Area, gp) * this->computeLength() / 2.;
140  answer.resize(12, 12);
141  answer.zero();
142  answer.at(1, 1) = answer.at(2, 2) = answer.at(3, 3) = halfMass;
143  answer.at(7, 7) = answer.at(8, 8) = answer.at(9, 9) = halfMass;
144 }
145 
146 
147 void
149 // Returns the displacement interpolation matrix {N} of the receiver, eva-
150 // luated at gp.
151 {
152  double ksi, n1, n2;
153 
154  ksi = iLocCoord.at(1);
155  n1 = ( 1. - ksi ) * 0.5;
156  n2 = ( 1. + ksi ) * 0.5;
157 
158  answer.resize(6, 12);
159  answer.zero();
160 
161  answer.at(1, 1) = n1;
162  answer.at(1, 7) = n2;
163  answer.at(2, 2) = n1;
164  answer.at(2, 8) = n2;
165  answer.at(3, 3) = n1;
166  answer.at(3, 9) = n2;
167 
168  answer.at(4, 4) = n1;
169  answer.at(4, 10) = n2;
170  answer.at(5, 5) = n1;
171  answer.at(5, 11) = n2;
172  answer.at(6, 6) = n1;
173  answer.at(6, 12) = n2;
174 }
175 
176 void
178 // Returns the stiffness matrix of the receiver, expressed in the global
179 // axes.
180 {
181  StructuralElement :: computeStiffnessMatrix(answer, rMode, tStep);
182 }
183 
184 
185 bool
187 {
188  FloatMatrix lcs;
189 
190  answer.resize(12, 12);
191  answer.zero();
192 
193  this->giveLocalCoordinateSystem(lcs);
194  for ( int i = 1; i <= 3; i++ ) {
195  for ( int j = 1; j <= 3; j++ ) {
196  answer.at(i, j) = lcs.at(i, j);
197  answer.at(i + 3, j + 3) = lcs.at(i, j);
198  answer.at(i + 6, j + 6) = lcs.at(i, j);
199  answer.at(i + 9, j + 9) = lcs.at(i, j);
200  }
201  }
202  return true;
203 }
204 
205 
206 double
208 // Returns the length of the receiver. This method is valid only if 1
209 // Gauss point is used.
210 {
211  double weight = gp->giveWeight();
212  return weight * 0.5 * this->computeLength();
213 }
214 
215 
216 int
218 {
221  if ( type == IST_BeamForceMomentTensor ) {
222  answer = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
223  return 1;
224  } else if ( type == IST_BeamStrainCurvatureTensor ) {
225  answer = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
226  return 1;
227  } else {
228  return StructuralElement :: giveIPValue(answer, gp, type, tStep);
229  }
230 }
231 
232 
233 void
234 LIBeam3d :: giveDofManDofIDMask(int inode, IntArray &answer) const
235 {
236  answer = {D_u, D_v, D_w, R_u, R_v, R_w};
237 }
238 
239 
240 int
242 {
243  double ksi, n1, n2;
244 
245  ksi = lcoords.at(1);
246  n1 = ( 1. - ksi ) * 0.5;
247  n2 = ( 1. + ksi ) * 0.5;
248 
249  answer.resize(3);
250  answer.at(1) = n1 * this->giveNode(1)->giveCoordinate(1) + n2 *this->giveNode(2)->giveCoordinate(1);
251  answer.at(2) = n1 * this->giveNode(1)->giveCoordinate(2) + n2 *this->giveNode(2)->giveCoordinate(2);
252  answer.at(3) = n1 * this->giveNode(1)->giveCoordinate(3) + n2 *this->giveNode(2)->giveCoordinate(3);
253 
254  return 1;
255 }
256 
257 
258 int
260 {
261  return ( ( ext == Element_EdgeLoadSupport ) ? 1 : 0 );
262 }
263 
264 
265 double
267 // Returns the length of the receiver.
268 {
269  double dx, dy, dz;
270  Node *nodeA, *nodeB;
271 
272  if ( length == 0. ) {
273  nodeA = this->giveNode(1);
274  nodeB = this->giveNode(2);
275  dx = nodeB->giveCoordinate(1) - nodeA->giveCoordinate(1);
276  dy = nodeB->giveCoordinate(2) - nodeA->giveCoordinate(2);
277  dz = nodeB->giveCoordinate(3) - nodeA->giveCoordinate(3);
278  length = sqrt(dx * dx + dy * dy + dz * dz);
279  }
280 
281  return length;
282 }
283 
284 
287 {
288  IRResultType result; // Required by IR_GIVE_FIELD macro
289 
291  if ( referenceNode == 0 ) {
292  OOFEM_ERROR("wrong reference node specified");
293  }
294 
295  // if (this->hasString (initString, "dofstocondense")) {
296  // dofsToCondense = this->ReadIntArray (initString, "dofstocondense");
297  // if (dofsToCondense->giveSize() >= 12)
298  // OOFEM_ERROR("wrong input data for condensed dofs");
299  // } else {
300  // dofsToCondense = NULL;
301  // }
302 
304 }
305 
306 
307 void
308 LIBeam3d :: giveEdgeDofMapping(IntArray &answer, int iEdge) const
309 {
310  /*
311  * provides dof mapping of local edge dofs (only nonzero are taken into account)
312  * to global element dofs
313  */
314  if ( iEdge != 1 ) {
315  OOFEM_ERROR("wrong edge number");
316  }
317 
318 
319  answer.resize(12);
320  for ( int i = 1; i <= 12; i++ ) {
321  answer.at(i) = i;
322  }
323 }
324 
325 double
327 {
328  if ( iEdge != 1 ) { // edge between nodes 1 2
329  OOFEM_ERROR("wrong egde number");
330  }
331 
332  double weight = gp->giveWeight();
333  return 0.5 * this->computeLength() * weight;
334 }
335 
336 int
338 {
339  /*
340  * Returns transformation matrix from global coordinate system to local
341  * element coordinate system for element load vector components.
342  * If no transformation is necessary, answer is empty matrix (default);
343  */
344 
345  // f(elemLocal) = T * f (global)
346 
347  FloatMatrix lcs;
348 
349  answer.resize(6, 6);
350  answer.zero();
351 
352  this->giveLocalCoordinateSystem(lcs);
353  for ( int i = 1; i <= 3; i++ ) {
354  for ( int j = 1; j <= 3; j++ ) {
355  answer.at(i, j) = lcs.at(i, j);
356  answer.at(3 + i, 3 + j) = lcs.at(i, j);
357  }
358  }
359  return 1;
360 }
361 
362 
363 int
365 {
366  // returns transformation matrix from
367  // edge local coordinate system
368  // to element local c.s
369  // (same as global c.s in this case)
370  //
371  // i.e. f(element local) = T * f(edge local)
372  //
373  answer.clear();
374  return 0;
375 }
376 
377 void
379 {
380  FloatArray lc(1);
381  StructuralElement :: computeBodyLoadVectorAt(answer, load, tStep, mode);
382  answer.times( this->giveCrossSection()->give(CS_Area, lc, this) );
383 }
384 
385 int
387 //
388 // returns a unit vectors of local coordinate system at element
389 // stored rowwise (mainly used by some materials with ortho and anisotrophy)
390 //
391 {
392  FloatArray lx(3), ly(3), lz(3), help(3);
393  double length = this->computeLength();
394  Node *nodeA, *nodeB, *refNode;
395 
396  answer.resize(3, 3);
397  answer.zero();
398  nodeA = this->giveNode(1);
399  nodeB = this->giveNode(2);
400  refNode = this->giveDomain()->giveNode(this->referenceNode);
401 
402  for ( int i = 1; i <= 3; i++ ) {
403  lx.at(i) = ( nodeB->giveCoordinate(i) - nodeA->giveCoordinate(i) ) / length;
404  help.at(i) = ( refNode->giveCoordinate(i) - nodeA->giveCoordinate(i) );
405  }
406 
407  lz.beVectorProductOf(lx, help);
408  lz.normalize();
409  ly.beVectorProductOf(lz, lx);
410  ly.normalize();
411 
412  for ( int i = 1; i <= 3; i++ ) {
413  answer.at(1, i) = lx.at(i);
414  answer.at(2, i) = ly.at(i);
415  answer.at(3, i) = lz.at(i);
416  }
417 
418  return 1;
419 }
420 
421 void
423  GaussPoint *slaveGp, TimeStep *tStep)
424 {
425  double layerYCoord, layerZCoord;
426 
427  layerZCoord = slaveGp->giveNaturalCoordinate(2);
428  layerYCoord = slaveGp->giveNaturalCoordinate(1);
429 
430  answer.resize(3); // {Exx,GMzx,GMxy}
431 
432  answer.at(1) = masterGpStrain.at(1) + masterGpStrain.at(5) * layerZCoord - masterGpStrain.at(6) * layerYCoord;
433  answer.at(2) = masterGpStrain.at(2);
434  answer.at(3) = masterGpStrain.at(3);
435 }
436 
437 Interface *
439 {
440  if ( interface == FiberedCrossSectionInterfaceType ) {
441  return static_cast< FiberedCrossSectionInterface * >(this);
442  }
443 
444  return NULL;
445 }
446 
447 
448 #ifdef __OOFEG
449 void
451 {
452  GraphicObj *go;
453 
454  if ( !gc.testElementGraphicActivity(this) ) {
455  return;
456  }
457 
458  // if (!go) { // create new one
459  WCRec p [ 2 ]; /* poin */
460  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
461  EASValsSetColor( gc.getElementColor() );
462  EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
463  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
464  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
465  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
466  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
467  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
468  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
469  go = CreateLine3D(p);
470  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
471  EGAttachObject(go, ( EObjectP ) this);
472  EMAddGraphicsToModel(ESIModel(), go);
473 }
474 
475 
476 void
478 {
479  GraphicObj *go;
480 
481  if ( !gc.testElementGraphicActivity(this) ) {
482  return;
483  }
484 
485  double defScale = gc.getDefScale();
486  // if (!go) { // create new one
487  WCRec p [ 2 ]; /* poin */
488  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
489  EASValsSetColor( gc.getDeformedElementColor() );
490  EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
491  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
492  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
493  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
494 
495  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
496  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
497  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
498  go = CreateLine3D(p);
499  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
500  EMAddGraphicsToModel(ESIModel(), go);
501 }
502 #endif
503 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
int testElementGraphicActivity(Element *)
Test if particular element passed fulfills various filtering criteria for its graphics output...
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Class and object Domain.
Definition: domain.h:115
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Definition: libeam3d.C:133
int referenceNode
Definition: libeam3d.h:56
virtual int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int, GaussPoint *gp)
Returns transformation matrix from local edge c.s to element local coordinate system of load vector c...
Definition: libeam3d.C:364
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: libeam3d.C:107
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 Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: libeam3d.C:438
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: libeam3d.C:217
#define OOFEG_RAW_GEOMETRY_LAYER
This class implements a structural material status information.
Definition: structuralms.h:65
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual void giveGeneralizedStress_Beam3d(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)=0
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
Definition: libeam3d.C:119
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: libeam3d.C:450
virtual double giveCoordinate(int i)
Definition: node.C:82
LIBeam3d(int n, Domain *d)
Definition: libeam3d.C:55
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
MatResponseMode
Describes the character of characteristic material matrix.
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
Definition: libeam3d.C:148
#define OOFEG_DEFORMED_GEOMETRY_LAYER
#define _IFT_LIBeam3d_refnode
Definition: libeam3d.h:44
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
Definition: libeam3d.C:477
Abstract base class for all "structural" finite elements.
double giveNaturalCoordinate(int i) const
Returns i-th natural element coordinate of receiver.
Definition: gausspoint.h:136
virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode)
Computes the load vector due to body load acting on receiver, at given time step. ...
virtual double computeLength()
Computes the length (zero for all but 1D geometries)
Definition: libeam3d.C:266
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: libeam3d.C:286
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: libeam3d.C:207
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes numerically stiffness matrix of receiver.
Definition: libeam3d.C:177
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
ElementExtension
Type representing element extension.
#define OOFEG_RAW_GEOMETRY_WIDTH
virtual double giveWeight()
Returns integration weight of receiver.
Definition: gausspoint.h:181
UnknownType
Type representing particular unknown (its physical meaning).
Definition: unknowntype.h:55
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
Definition: libeam3d.C:308
double length
Definition: libeam3d.h:55
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element&#39;s local coordinates.
Definition: libeam3d.C:241
virtual double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Returns updated ic-th coordinate of receiver.
Definition: node.C:245
virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
Computes the stress vector of receiver at given integration point, at time step tStep.
Definition: libeam3d.C:126
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes numerically stiffness matrix of receiver.
void resize(int n)
Checks size of receiver towards requested bounds.
Definition: intarray.C:124
virtual void give3dBeamStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
Computes the stiffness matrix for 2d beams.
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
Definition: libeam3d.C:337
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
Definition: crosssection.C:54
virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode)
Computes the load vector due to body load acting on receiver, at given time step. ...
Definition: libeam3d.C:378
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Definition: libeam3d.C:65
Class representing vector of real numbers.
Definition: floatarray.h:82
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
virtual double give(CrossSectionProperty a, GaussPoint *gp)
Returns the value of cross section property at given point.
Definition: crosssection.C:151
virtual int giveLocalCoordinateSystem(FloatMatrix &answer)
Returns local coordinate system of receiver Required by material models with ortho- and anisotrophy...
Definition: libeam3d.C:386
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
Class representing the general Input Record.
Definition: inputrecord.h:101
Class Interface.
Definition: interface.h:82
virtual void FiberedCrossSectionInterface_computeStrainVectorInFiber(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *slaveGp, TimeStep *tStep)
Computes full 3d strain vector in element fiber.
Definition: libeam3d.C:422
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
void times(double s)
Multiplies receiver with scalar.
Definition: floatarray.C:818
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
The element interface required by FiberedCrossSection.
Definition: fiberedcs.h:220
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
Domain * giveDomain() const
Definition: femcmpnn.h:100
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
Load is base abstract class for all loads.
Definition: load.h:61
virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
Computes volume related to integration point on local edge.
Definition: libeam3d.C:326
Node * giveNode(int n)
Service for accessing particular domain node.
Definition: domain.h:371
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual int testElementExtension(ElementExtension ext)
Tests if the element implements required extension.
Definition: libeam3d.C:259
void clear()
Sets size of receiver to be an empty matrix. It will have zero rows and zero columns size...
Definition: floatmatrix.h:516
Class implementing node in finite element mesh.
Definition: node.h:87
#define IR_GIVE_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:69
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
virtual bool computeGtoLRotationMatrix(FloatMatrix &answer)
Returns transformation matrix from global c.s.
Definition: libeam3d.C:186
Element extension for edge loads.
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
Definition: libeam3d.C:234

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:29 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011