OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
truss3d.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/Bars/truss3d.h"
36 #include "../sm/CrossSections/structuralcrosssection.h"
37 #include "fei3dlinelin.h"
38 #include "node.h"
39 #include "material.h"
40 #include "gausspoint.h"
41 #include "gaussintegrationrule.h"
42 #include "floatmatrix.h"
43 #include "floatarray.h"
44 #include "intarray.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(Truss3d);
54 
55 FEI3dLineLin Truss3d :: interp;
56 
57 Truss3d :: Truss3d(int n, Domain *aDomain) :
59 {
60  numberOfDofMans = 2;
61 }
62 
63 
65 
66 
67 Interface *
69 {
70  if ( interface == ZZNodalRecoveryModelInterfaceType ) {
71  return static_cast< ZZNodalRecoveryModelInterface * >(this);
72  } else if ( interface == NodalAveragingRecoveryModelInterfaceType ) {
73  return static_cast< NodalAveragingRecoveryModelInterface * >(this);
74  }
75 
76  //OOFEM_LOG_INFO("Interface on Truss3d element not supported");
77  return NULL;
78 }
79 
80 
81 void
83 {
84  answer.clear();
85  OOFEM_WARNING("IP values will not be transferred to nodes. Use ZZNodalRecovery instead (parameter stype 1)");
86 }
87 
88 
89 void
90 Truss3d :: computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int li, int ui)
91 //
92 // Returns linear part of geometrical equations of the receiver at gp.
93 // Returns the linear part of the B matrix
94 //
95 {
96  FloatMatrix dN;
98 
99  answer.resize(1, 6);
100  answer.at(1, 1) = dN.at(1, 1);
101  answer.at(1, 2) = dN.at(1, 2);
102  answer.at(1, 3) = dN.at(1, 3);
103  answer.at(1, 4) = dN.at(2, 1);
104  answer.at(1, 5) = dN.at(2, 2);
105  answer.at(1, 6) = dN.at(2, 3);
106 }
107 
108 
109 void
111 // Sets up the array of Gauss Points of the receiver.
112 {
113  if ( integrationRulesArray.size() == 0 ) {
114  integrationRulesArray.resize( 1 );
115  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 2) );
117  }
118 }
119 
120 
121 double
123 {
124  return this->interp.giveLength( FEIElementGeometryWrapper(this) );
125 }
126 
127 
128 void
130 // Returns the lumped mass matrix of the receiver. This expression is
131 // valid in both local and global axes.
132 {
133  answer.resize(6, 6);
134  answer.zero();
135  if ( !this->isActivated(tStep) ) {
136  return;
137  }
138 
139  GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
140  double density = this->giveStructuralCrossSection()->give('d', gp);
141  double halfMass = density * this->giveCrossSection()->give(CS_Area, gp) * this->computeLength() * 0.5;
142  answer.at(1, 1) = halfMass;
143  answer.at(2, 2) = halfMass;
144  answer.at(3, 3) = halfMass;
145  answer.at(4, 4) = halfMass;
146  answer.at(5, 5) = halfMass;
147  answer.at(6, 6) = halfMass;
148 }
149 
150 
151 void
153 // Returns the displacement interpolation matrix {N} of the receiver, eva-
154 // luated at gp.
155 {
156  FloatArray n;
157  this->interp.evalN( n, iLocCoord, FEIElementGeometryWrapper(this) );
158  answer.beNMatrixOf(n, 3);
159 }
160 
161 
162 double
164 // Returns the length of the receiver. This method is valid only if 1
165 // Gauss point is used.
166 {
168  double weight = gp->giveWeight();
169  return detJ *weight *this->giveCrossSection()->give(CS_Area, gp);
170 }
171 
172 
173 int
175 //
176 // returns a unit vectors of local coordinate system at element
177 // stored rowwise (mainly used by some materials with ortho and anisotrophy)
178 //
179 {
180  FloatArray lx, ly(3), lz;
181 
182  lx.beDifferenceOf( * this->giveNode(2)->giveCoordinates(), * this->giveNode(1)->giveCoordinates() );
183  lx.normalize();
184 
185  ly(0) = lx(1);
186  ly(1) = -lx(2);
187  ly(2) = lx(0);
188 
189  // Construct orthogonal vector
190  double npn = ly.dotProduct(lx);
191  ly.add(-npn, lx);
192  ly.normalize();
193  lz.beVectorProductOf(ly, lx);
194 
195  answer.resize(3, 3);
196  for ( int i = 1; i <= 3; i++ ) {
197  answer.at(1, i) = lx.at(i);
198  answer.at(2, i) = ly.at(i);
199  answer.at(3, i) = lz.at(i);
200  }
201 
202  return 1;
203 }
204 
205 
208 {
210 }
211 
212 
213 void
215 {
216  this->giveStructuralCrossSection()->giveRealStress_1d(answer, gp, strain, tStep);
217 }
218 
219 void
221 {
222  this->giveStructuralCrossSection()->giveStiffnessMatrix_1d(answer, rMode, gp, tStep);
223 }
224 
225 void
226 Truss3d :: giveDofManDofIDMask(int inode, IntArray &answer) const
227 {
228  answer = {D_u, D_v, D_w};
229 }
230 
231 
232 void
233 Truss3d :: giveEdgeDofMapping(IntArray &answer, int iEdge) const
234 {
235  /*
236  * provides dof mapping of local edge dofs (only nonzero are taken into account)
237  * to global element dofs
238  */
239 
240  if ( iEdge != 1 ) {
241  OOFEM_ERROR("wrong edge number");
242  }
243 
244  answer.resize(6);
245  answer.at(1) = 1;
246  answer.at(2) = 2;
247  answer.at(3) = 3;
248  answer.at(4) = 4;
249  answer.at(5) = 5;
250  answer.at(6) = 6;
251 }
252 
253 
254 double
256 {
257  if ( iEdge != 1 ) { // edge between nodes 1 2
258  OOFEM_ERROR("wrong edge number");
259  }
260 
261  double weight = gp->giveWeight();
263 }
264 
265 
266 int
268 {
269  // returns transformation matrix from
270  // edge local coordinate system
271  // to element local c.s
272  // (same as global c.s in this case)
273  //
274  // i.e. f(element local) = T * f(edge local)
275  //
276  FloatMatrix lcs;
277  this->giveLocalCoordinateSystem(lcs);
278  answer.beTranspositionOf(lcs);
279 
280  return 1;
281 }
282 
283 
284 #ifdef __OOFEG
286 {
287  GraphicObj *go;
288  // if (!go) { // create new one
289  WCRec p [ 2 ]; /* point */
290  if ( !gc.testElementGraphicActivity(this) ) {
291  return;
292  }
293 
294  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
295  EASValsSetColor( gc.getElementColor() );
296  EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
297  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
298  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
299  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
300  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
301  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
302  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
303  go = CreateLine3D(p);
304  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
305  EGAttachObject(go, ( EObjectP ) this);
306  EMAddGraphicsToModel(ESIModel(), go);
307 }
308 
309 
311 {
312  GraphicObj *go;
313  double defScale = gc.getDefScale();
314  // if (!go) { // create new one
315  WCRec p [ 2 ]; /* point */
316  if ( !gc.testElementGraphicActivity(this) ) {
317  return;
318  }
319 
320  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
321  EASValsSetColor( gc.getDeformedElementColor() );
322  EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
323  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
324  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
325  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
326 
327  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
328  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
329  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
330  go = CreateLine3D(p);
331  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
332  EMAddGraphicsToModel(ESIModel(), go);
333 }
334 #endif
335 } // end namespace oofem
virtual double computeEdgeVolumeAround(GaussPoint *gp, int)
Computes volume related to integration point on local edge.
Definition: truss3d.C:255
CrossSection * giveCrossSection()
Definition: element.C:495
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Definition: truss3d.C:129
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.
The element interface required by NodalAvergagingRecoveryModel.
virtual bool isActivated(TimeStep *tStep)
Definition: element.C:793
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: truss3d.C:207
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
The element interface required by ZZNodalRecoveryModel.
Abstract base class for "structural" finite elements with geometrical nonlinearities.
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
Definition: truss3d.C:220
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
#define OOFEG_RAW_GEOMETRY_LAYER
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual FEInterpolation * giveInterpolation() const
Definition: truss3d.C:64
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
Definition: truss3d.C:310
virtual double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
Definition: fei3dlinelin.C:112
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
MatResponseMode
Describes the character of characteristic material matrix.
virtual double evaldNdx(FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the matrix of derivatives of interpolation functions (shape functions) at given point...
Definition: fei3dlinelin.C:59
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: truss3d.C:267
static FEI3dLineLin interp
Definition: truss3d.h:57
#define OOFEG_DEFORMED_GEOMETRY_LAYER
void beDifferenceOf(const FloatArray &a, const FloatArray &b)
Sets receiver to be a - b.
Definition: floatarray.C:341
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
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: truss3d.C:214
virtual double computeLength()
Computes the length (zero for all but 1D geometries)
Definition: truss3d.C:122
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
Definition: truss3d.C:152
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
Definition: fei3dlinelin.C:49
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: truss3d.C:68
#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.
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
virtual double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Returns updated ic-th coordinate of receiver.
Definition: node.C:245
virtual int giveLocalCoordinateSystem(FloatMatrix &answer)
Returns local coordinate system of receiver Required by material models with ortho- and anisotrophy...
Definition: truss3d.C:174
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
Truss3d(int n, Domain *d)
Definition: truss3d.C:57
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
Definition: crosssection.C:54
void beNMatrixOf(const FloatArray &n, int nsd)
Assigns the receiver to be a repeated diagonal matrix.
Definition: floatmatrix.C:505
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
Definition: truss3d.C:82
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: truss3d.C:163
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 void giveEdgeDofMapping(IntArray &answer, int iEdge) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
Definition: truss3d.C:233
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 giveRealStress_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)=0
virtual double giveLength(const FEICellGeometry &cellgeo) const
Definition: fei3dlinelin.C:43
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: truss3d.C:110
void beTranspositionOf(const FloatMatrix &src)
Assigns to the receiver the transposition of parameter.
Definition: floatmatrix.C:323
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 IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
virtual void giveStiffnessMatrix_1d(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: truss3d.C:285
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Definition: truss3d.C:90
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
Definition: truss3d.C:226
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
Class representing integration point in finite element program.
Definition: gausspoint.h:93
#define OOFEM_WARNING(...)
Definition: error.h:62
Class representing solution step.
Definition: timestep.h:80
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
Class representing Gaussian-quadrature integration rule.

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