OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
libeam2d.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/libeam2d.h"
36 #include "../sm/Materials/structuralms.h"
37 #include "fei2dlinelin.h"
38 #include "node.h"
39 #include "material.h"
40 #include "crosssection.h"
41 #include "gausspoint.h"
42 #include "gaussintegrationrule.h"
43 #include "floatmatrix.h"
44 #include "intarray.h"
45 #include "floatarray.h"
46 #include "mathfem.h"
47 #include "classfactory.h"
48 
49 namespace oofem {
50 REGISTER_Element(LIBeam2d);
51 
52 // Set up interpolation coordinates
53 FEI2dLineLin LIBeam2d :: interpolation(1, 3);
54 
56 {
57  numberOfDofMans = 2;
58  length = 0.;
59  pitch = 10.; // a dummy value
60 }
61 
62 
64 
65 
66 Interface *
68 {
69  if ( interface == LayeredCrossSectionInterfaceType ) {
70  return static_cast< LayeredCrossSectionInterface * >(this);
71  }
72 
73  return NULL;
74 }
75 
76 
77 void
79 // Returns the strain matrix of the receiver.
80 {
81  double l, ksi, n1x, n4x, n3xx, n6xx, n2xxx, n3xxx, n5xxx, n6xxx;
82 
83  l = this->computeLength();
84  ksi = gp->giveNaturalCoordinate(1);
85  n1x = -1.0 / l;
86  n4x = 1.0 / l;
87  n3xx = -1.0 / l;
88  n6xx = 1.0 / l;
89  n2xxx = -1.0 / l;
90  n3xxx = 0.5 * ( 1 - ksi );
91  n5xxx = 1.0 / l;
92  n6xxx = 0.5 * ( 1. + ksi );
93 
94  answer.resize(3, 6);
95  answer.zero();
96 
97  answer.at(1, 1) = n1x;
98  answer.at(1, 4) = n4x;
99  answer.at(2, 3) = n3xx;
100  answer.at(2, 6) = n6xx;
101  answer.at(3, 2) = n2xxx;
102  answer.at(3, 3) = n3xxx;
103  answer.at(3, 5) = n5xxx;
104  answer.at(3, 6) = n6xxx;
105 }
106 
107 
108 void
110 // Sets up the array of Gauss Points of the receiver.
111 {
112  if ( integrationRulesArray.size() == 0 ) {
113  integrationRulesArray.resize( 1 );
114  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 2) );
116  }
117 }
118 
119 
120 void
122 {
123  this->giveStructuralCrossSection()->give2dBeamStiffMtrx(answer, rMode, gp, tStep);
124 }
125 
126 
127 void
129 {
130  this->giveStructuralCrossSection()->giveGeneralizedStress_Beam2d(answer, gp, strain, tStep);
131 }
132 
133 
134 int
136 {
137  double ksi, n1, n2;
138 
139  ksi = lcoords.at(1);
140  n1 = ( 1. - ksi ) * 0.5;
141  n2 = ( 1. + ksi ) * 0.5;
142 
143  answer.resize(3);
144  answer.at(1) = n1 * this->giveNode(1)->giveCoordinate(1) + n2 *this->giveNode(2)->giveCoordinate(1);
145  answer.at(3) = n1 * this->giveNode(1)->giveCoordinate(3) + n2 *this->giveNode(2)->giveCoordinate(3);
146 
147  return 1;
148 }
149 
150 
151 void
153 // Returns the lumped mass matrix of the receiver. This expression is
154 // valid in both local and global axes.
155 {
156  GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
157  double density = this->giveStructuralCrossSection()->give('d', gp);
158  double halfMass = density * this->giveCrossSection()->give(CS_Area, gp) * this->computeLength() / 2.;
159  answer.resize(6, 6);
160  answer.zero();
161  answer.at(1, 1) = halfMass;
162  answer.at(2, 2) = halfMass;
163  answer.at(4, 4) = halfMass;
164  answer.at(5, 5) = halfMass;
165 }
166 
167 
168 void
170 // Returns the displacement interpolation matrix {N} of the receiver, eva-
171 // luated at gp.
172 {
173  double ksi, n1, n2;
174 
175  ksi = iLocCoord.at(1);
176  n1 = ( 1. - ksi ) * 0.5;
177  n2 = ( 1. + ksi ) * 0.5;
178 
179  answer.resize(3, 6);
180  answer.zero();
181 
182  answer.at(1, 1) = n1;
183  answer.at(1, 4) = n2;
184  answer.at(2, 2) = n1;
185  answer.at(2, 5) = n2;
186  answer.at(3, 3) = n1;
187  answer.at(3, 6) = n2;
188 }
189 
190 
191 void
193 // Returns the stiffness matrix of the receiver, expressed in the global
194 // axes.
195 {
196  StructuralElement :: computeStiffnessMatrix(answer, rMode, tStep);
197 }
198 
199 
200 bool
202 {
203  double sine, cosine;
204 
205  sine = sin( this->givePitch() );
206  cosine = cos(pitch);
207 
208  answer.resize(6, 6);
209  answer.zero();
210 
211  answer.at(1, 1) = cosine;
212  answer.at(1, 2) = sine;
213  answer.at(2, 1) = -sine;
214  answer.at(2, 2) = cosine;
215  answer.at(3, 3) = 1.;
216  answer.at(4, 4) = cosine;
217  answer.at(4, 5) = sine;
218  answer.at(5, 4) = -sine;
219  answer.at(5, 5) = cosine;
220  answer.at(6, 6) = 1.;
221 
222  return true;
223 }
224 
225 
226 double
228 // Returns the length of the receiver. This method is valid only if 1
229 // Gauss point is used.
230 {
231  double weight = gp->giveWeight();
232  return weight * 0.5 * this->computeLength();
233 }
234 
235 
236 void
237 LIBeam2d :: computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep)
238 //
239 // returns full 3d strain vector of given layer (whose z-coordinate from center-line is
240 // stored in slaveGp) for given tStep
241 //
242 {
243  double layerZeta, layerZCoord, top, bottom;
244 
245  top = this->giveCrossSection()->give(CS_TopZCoord, masterGp);
246  bottom = this->giveCrossSection()->give(CS_BottomZCoord, masterGp);
247  layerZeta = slaveGp->giveNaturalCoordinate(3);
248  layerZCoord = 0.5 * ( ( 1. - layerZeta ) * bottom + ( 1. + layerZeta ) * top );
249 
250  answer.resize(2); // {Exx,GMzx}
251 
252  answer.at(1) = masterGpStrain.at(1) + masterGpStrain.at(2) * layerZCoord;
253  answer.at(2) = masterGpStrain.at(3);
254 }
255 
256 
257 int
259 {
260  if ( type == IST_BeamForceMomentTensor ) {
261  answer = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
262  return 1;
263  } else if ( type == IST_BeamStrainCurvatureTensor ) {
264  answer = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
265  return 1;
266  } else {
267  return StructuralElement :: giveIPValue(answer, gp, type, tStep);
268  }
269 }
270 
271 
272 void
273 LIBeam2d :: giveDofManDofIDMask(int inode, IntArray &answer) const
274 {
275  answer = {D_u, D_w, R_v};
276 }
277 
278 
279 double
281 // Returns the length of the receiver.
282 {
283  double dx, dy;
284  Node *nodeA, *nodeB;
285 
286  if ( length == 0. ) {
287  nodeA = this->giveNode(1);
288  nodeB = this->giveNode(2);
289  dx = nodeB->giveCoordinate(1) - nodeA->giveCoordinate(1);
290  dy = nodeB->giveCoordinate(3) - nodeA->giveCoordinate(3);
291  length = sqrt(dx * dx + dy * dy);
292  }
293 
294  return length;
295 }
296 
297 
298 double
300 // Returns the pitch of the receiver.
301 {
302  double xA, xB, yA, yB;
303  Node *nodeA, *nodeB;
304 
305  if ( pitch == 10. ) { // 10. : dummy initialization value
306  nodeA = this->giveNode(1);
307  nodeB = this->giveNode(2);
308  xA = nodeA->giveCoordinate(1);
309  xB = nodeB->giveCoordinate(1);
310  yA = nodeA->giveCoordinate(3);
311  yB = nodeB->giveCoordinate(3);
312  pitch = atan2(yB - yA, xB - xA);
313  }
314 
315  return pitch;
316 }
317 
318 
321 {
323 }
324 
325 
326 
327 void
328 LIBeam2d :: giveEdgeDofMapping(IntArray &answer, int iEdge) const
329 {
330  /*
331  * provides dof mapping of local edge dofs (only nonzero are taken into account)
332  * to global element dofs
333  */
334 
335  if ( iEdge != 1 ) {
336  OOFEM_ERROR("wrong edge number");
337  }
338 
339  answer.resize(6);
340  answer.at(1) = 1;
341  answer.at(2) = 2;
342  answer.at(3) = 3;
343  answer.at(4) = 4;
344  answer.at(5) = 5;
345  answer.at(6) = 6;
346 }
347 
348 
349 double
351 {
352  if ( iEdge != 1 ) { // edge between nodes 1 2
353  OOFEM_ERROR("wrong egde number");
354  }
355 
356  double weight = gp->giveWeight();
357  return 0.5 * this->computeLength() * weight;
358 }
359 
360 
361 void
363 {
364  FloatArray lc(1);
365  StructuralElement :: computeBodyLoadVectorAt(answer, load, tStep, mode);
366  answer.times( this->giveCrossSection()->give(CS_Area, lc, this) );
367 }
368 
369 
370 int
372 {
373  /*
374  * Returns transformation matrix from global coordinate system to local
375  * element coordinate system for element load vector components.
376  * If no transformation is necessary, answer is empty matrix (default);
377  */
378 
379  // f(elemLocal) = T * f (global)
380 
381  double sine, cosine;
382 
383  answer.resize(3, 3);
384  answer.zero();
385 
386  sine = sin( this->givePitch() );
387  cosine = cos(pitch);
388 
389  answer.at(1, 1) = cosine;
390  answer.at(1, 2) = -sine;
391  answer.at(2, 1) = sine;
392  answer.at(2, 2) = cosine;
393  answer.at(3, 3) = 1.0;
394 
395  return 1;
396 }
397 
398 
399 int
401 {
402  // returns transformation matrix from
403  // edge local coordinate system
404  // to element local c.s
405  // (same as global c.s in this case)
406  //
407  // i.e. f(element local) = T * f(edge local)
408  //
409  answer.clear();
410  return 0;
411 }
412 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
virtual void computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep)
Computes full 3D strain vector in element layer.
Definition: libeam2d.C:237
double givePitch()
Definition: libeam2d.C:299
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: libeam2d.C:67
double length
Definition: libeam2d.h:57
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
virtual int computeLoadLEToLRotationMatrix(FloatMatrix &, int, GaussPoint *)
Returns transformation matrix from local edge c.s to element local coordinate system of load vector c...
Definition: libeam2d.C:400
virtual void giveGeneralizedStress_Beam2d(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)=0
Computes the generalized stress vector for given strain and integration point.
Class and object Domain.
Definition: domain.h:115
Bottom z coordinate.
Definition: crosssection.h:72
LIBeam2d(int n, Domain *aDomain)
Definition: libeam2d.C:55
virtual void computeBmatrixAt(GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Definition: libeam2d.C:78
virtual double computeLength()
Computes the length (zero for all but 1D geometries)
Definition: libeam2d.C:280
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: libeam2d.C:227
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes numerically stiffness matrix of receiver.
Definition: libeam2d.C:192
This class implements a structural material status information.
Definition: structuralms.h:65
virtual FEInterpolation * giveInterpolation() const
Definition: libeam2d.C:63
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
Definition: libeam2d.C:121
Top z coordinate.
Definition: crosssection.h:71
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 computeEdgeVolumeAround(GaussPoint *, int)
Computes volume related to integration point on local edge.
Definition: libeam2d.C:350
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &)
Computes interpolation matrix for element unknowns.
Definition: libeam2d.C:169
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
Definition: libeam2d.C:371
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 bool computeGtoLRotationMatrix(FloatMatrix &answer)
Returns transformation matrix from global c.s.
Definition: libeam2d.C:201
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 void giveEdgeDofMapping(IntArray &answer, int) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
Definition: libeam2d.C:328
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
virtual double giveWeight()
Returns integration weight of receiver.
Definition: gausspoint.h:181
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: libeam2d.C:258
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
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 computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: libeam2d.C:109
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
Definition: crosssection.C:54
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
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 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: libeam2d.C:128
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
static FEI2dLineLin interpolation
Interpolation.
Definition: libeam2d.h:54
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
double pitch
Definition: libeam2d.h:57
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: libeam2d.C:362
Load is base abstract class for all loads.
Definition: load.h:61
The element interface required by LayeredCrossSection.
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
Definition: libeam2d.C:273
the oofem namespace is to define a context or scope in which all oofem names are defined.
void clear()
Sets size of receiver to be an empty matrix. It will have zero rows and zero columns size...
Definition: floatmatrix.h:516
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: libeam2d.C:320
Class implementing node in finite element mesh.
Definition: node.h:87
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: libeam2d.C:135
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual void give2dBeamStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
Computes the stiffness matrix for 2d beams.
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 void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Definition: libeam2d.C:152
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:29 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011