OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
beam3d.h
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 #ifndef beam3d_h
36 #define beam3d_h
37 
38 #include "../sm/Elements/Beams/beambaseelement.h"
39 #include "../sm/CrossSections/fiberedcs.h"
40 #include "../sm/Materials/winklermodel.h"
41 #include "dofmanager.h"
42 #include "vtkxmlexportmodule.h"
43 
45 
46 #define _IFT_Beam3d_Name "beam3d"
47 #define _IFT_Beam3d_dofstocondense "dofstocondense"
48 #define _IFT_Beam3d_refnode "refnode"
49 #define _IFT_Beam3d_refangle "refangle"
50 #define _IFT_Beam3d_zaxis "zaxis"
51 #define _IFT_Beam3d_subsoilmat "subsoilmat"
52 
53 
54 #define Beam3d_nSubBeams 10
55 
56 namespace oofem {
57 
58 class FEI3dLineLin;
59 
75 {
76 protected:
79 
80  double kappay, kappaz, length;
83  double referenceAngle = 0;
84  //IntArray *dofsToCondense;
85  /*
86  * Ghost nodes are used to introduce additional DOFs at element.
87  * These are needed as we actually do not want to condense selected DOFs, but rather
88  * allocate an extra equation to these. This allows to get cooresponding DOFs directly from
89  * the global system, avoiding the need to postprocess local displacements at element.
90  */
94 
97 
98 public:
99  Beam3d(int n, Domain *d);
100  virtual ~Beam3d();
101 
102  virtual FEInterpolation *giveInterpolation() const;
103 
104  virtual void computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity = NULL);
105  virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep);
106  virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep);
107  virtual int giveLocalCoordinateSystem(FloatMatrix &answer);
108  virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *, int useUpdatedGpRecord = 0);
109  void giveEndForcesVector(FloatArray &answer, TimeStep *tStep);
110 
111  virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords);
112 
113 
115  return ( ( ext == Element_EdgeLoadSupport ) ? 1 : 0 );
116  }
117  //int hasLayeredSupport () {return 1;}
118 
119  virtual int computeNumberOfDofs() { return 12; }
120  virtual int computeNumberOfGlobalDofs() { return 12 + this->numberOfCondensedDofs; }
121  virtual void giveDofManDofIDMask(int inode, IntArray &) const;
122  virtual int giveNumberOfInternalDofManagers() const { return ( ghostNodes [ 0 ] != NULL ) + ( ghostNodes [ 1 ] != NULL ); }
123  virtual DofManager *giveInternalDofManager(int i) const {
124  if ( i == 1 ) {
125  if ( ghostNodes [ 0 ] ) { return ghostNodes [ 0 ]; } else { return ghostNodes [ 1 ]; }
126  } else if ( i == 2 ) { // i==2
127  return ghostNodes [ 1 ];
128  } else {
129  OOFEM_ERROR("No such DOF available on Element %d", number);
130  return NULL;
131  }
132  }
133  virtual void giveInternalDofManDofIDMask(int i, IntArray &answer) const {
134  if ( i == 1 ) {
135  if ( ghostNodes [ 0 ] ) {
136  ghostNodes [ 0 ]->giveCompleteMasterDofIDArray(answer);
137  } else {
138  ghostNodes [ 1 ]->giveCompleteMasterDofIDArray(answer);
139  }
140  } else if ( i == 2 ) { // i==2
141  ghostNodes [ 1 ]->giveCompleteMasterDofIDArray(answer);
142  } else {
143  OOFEM_ERROR("No such DOF available on Element %d", number);
144  }
145  }
146  virtual void giveBoundaryLocationArray(IntArray &locationArray, const IntArray &bNodes, const UnknownNumberingScheme &s, IntArray *dofIds = NULL) {
147  giveLocationArray (locationArray, s, dofIds);
148  }
149 
150  virtual void giveBoundaryLocationArray(IntArray &locationArray, const IntArray &bNodes, const IntArray &dofIDMask, const UnknownNumberingScheme &s, IntArray *dofIds = NULL) {
151  giveLocationArray (locationArray, dofIDMask, s, dofIds);
152  }
153 
154  virtual double computeVolumeAround(GaussPoint *gp);
155 
156  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
157 
158  virtual void printOutputAt(FILE *file, TimeStep *tStep);
159 
160  //
161  // fibered cross section support functions
162  //
163  virtual void FiberedCrossSectionInterface_computeStrainVectorInFiber(FloatArray &answer, const FloatArray &masterGpStrain,
164  GaussPoint *slaveGp, TimeStep *tStep);
165 
167 
168  // definition & identification
169  virtual const char *giveClassName() const { return "Beam3d"; }
170  virtual const char *giveInputRecordName() const { return _IFT_Beam3d_Name; }
173  virtual integrationDomain giveIntegrationDomain() const { return _Line; }
174  //virtual Element_Geometry_Type giveGeometryType() const { return EGT_line_1; }
175  virtual Element_Geometry_Type giveGeometryType() const { return EGT_Composite; }
177 
178  /*
180 
181  virtual void B3SSI_getNMatrix (FloatMatrix &answer, GaussPoint *gp) {
182  this->computeNmatrixAt(gp->giveNaturalCoordinates(), answer);
183  }
184  virtual void B3SSI_getGtoLRotationMatrix (FloatMatrix &answer) {
185  this->computeGtoLRotationMatrix(answer);
186  }
187  virtual double B3SSI_computeVolumeAround (GaussPoint* gp) {
188  return this->computeVolumeAround(gp);
189  }
190  */
192 
193  virtual void giveCompositeExportData(std::vector< VTKPiece > &vtkPieces, IntArray &primaryVarsToExport, IntArray &internalVarsToExport, IntArray cellVarsToExport, TimeStep *tStep );
194 
195 
196 #ifdef __OOFEG
197  virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep);
198  virtual void drawDeformedGeometry(oofegGraphicContext & gc, TimeStep * tStep, UnknownType);
199 #endif
200 
201 protected:
202  virtual void computeBoundaryEdgeLoadVector(FloatArray &answer, BoundaryLoad *load, int edge, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true);
203  virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer);
204  virtual void computeBmatrixAt(GaussPoint *, FloatMatrix &, int = 1, int = ALL_STRAINS);
205  virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &);
206  virtual bool computeGtoLRotationMatrix(FloatMatrix &answer);
207  virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode);
208 
209  double giveKappayCoeff(TimeStep *tStep);
210  double giveKappazCoeff(TimeStep *tStep);
211  void computeKappaCoeffs(TimeStep *tStep);
212  virtual double computeLength();
213  virtual void computeGaussPoints();
214  virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep);
215  virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep);
216 
217  virtual MaterialMode giveMaterialMode() { return _3dBeam; }
218  virtual int giveNumberOfIPForMassMtrxIntegration() { return 4; }
219 
220  bool hasDofs2Condense() { return ( ghostNodes [ 0 ] || ghostNodes [ 1 ] ); }
221 
222 
225  MatResponseMode rMode, TimeStep *tStep);
226 
227  void giveInternalForcesVectorAtPoint(FloatArray &answer, TimeStep *tStep, FloatArray &coords);
229  TimeStep *tStep, FloatArray &pointCoords, double ds, bool global);
230  void computeInternalForcesFromBodyLoadVectorAtPoint(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode, FloatArray &pointCoords, double ds);
231  virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords, const FloatArray &point);
232 
233 };
234 } // end namespace oofem
235 #endif // beam3d_h
void giveInternalForcesVectorAtPoint(FloatArray &answer, TimeStep *tStep, FloatArray &coords)
Definition: beam3d.C:1002
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
Definition: beam3d.C:931
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
int number
Component number.
Definition: femcmpnn.h:80
integrationDomain
Used by integrator class to supply integration points for proper domain to be integrated (Area...
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes numerically stiffness matrix of receiver.
Definition: beam3d.C:194
void computeSubSoilStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Definition: beam3d.C:967
virtual Element_Geometry_Type giveGeometryType() const
Returns the element geometry type.
Definition: beam3d.h:175
bool hasDofs2Condense()
Definition: beam3d.h:220
Class and object Domain.
Definition: domain.h:115
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: beam3d.C:361
Element_Geometry_Type
Enumerative type used to classify element geometry Possible values are: EGT_point - point in space EG...
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: beam3d.C:638
virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes initial stress matrix for linear stability problem.
Definition: beam3d.C:791
void computeKappaCoeffs(TimeStep *tStep)
Definition: beam3d.C:431
void computeSubSoilNMatrixAt(GaussPoint *gp, FloatMatrix &answer)
Definition: beam3d.C:959
double giveKappazCoeff(TimeStep *tStep)
Definition: beam3d.C:467
Beam3d(int n, Domain *d)
Definition: beam3d.C:66
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
virtual void giveBoundaryLocationArray(IntArray &locationArray, const IntArray &bNodes, const IntArray &dofIDMask, const UnknownNumberingScheme &s, IntArray *dofIds=NULL)
Definition: beam3d.h:150
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
Definition: beam3d.C:262
double referenceAngle
Definition: beam3d.h:83
Elements with geometry defined as EGT_Composite are exported using individual pieces.
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual int giveLocalCoordinateSystem(FloatMatrix &answer)
Returns local coordinate system of receiver Required by material models with ortho- and anisotrophy...
Definition: beam3d.C:478
virtual void updateLocalNumbering(EntityRenumberingFunctor &f)
Local renumbering support.
Definition: beam3d.C:895
virtual FEInterpolation * giveInterpolation() const
Definition: beam3d.C:85
Base class for dof managers.
Definition: dofmanager.h:113
virtual integrationDomain giveIntegrationDomain() const
Definition: beam3d.h:173
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Definition: beam3d.C:672
double kappaz
Definition: beam3d.h:80
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *, int useUpdatedGpRecord=0)
Returns equivalent nodal forces vectors.
Definition: beam3d.C:605
virtual void B3SSMI_getUnknownsGtoLRotationMatrix(FloatMatrix &answer)
Evaluate transformation matrix for reciver unknowns.
Definition: beam3d.C:335
virtual MaterialMode giveMaterialMode()
Returns material mode for receiver integration points.
Definition: beam3d.h:217
virtual int computeNumberOfGlobalDofs()
Computes the total number of element&#39;s global dofs.
Definition: beam3d.h:120
double giveKappayCoeff(TimeStep *tStep)
Definition: beam3d.C:456
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: beam3d.C:545
virtual const char * giveClassName() const
Definition: beam3d.h:169
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
virtual ~Beam3d()
Definition: beam3d.C:79
Abstract base class representing a boundary load (force, momentum, ...) that acts directly on a bound...
Definition: boundaryload.h:110
void computeInternalForcesFromBodyLoadVectorAtPoint(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode, FloatArray &pointCoords, double ds)
Definition: beam3d.C:1157
int numberOfCondensedDofs
number of condensed DOFs
Definition: beam3d.h:93
virtual int computeNumberOfDofs()
Computes or simply returns total number of element&#39;s local DOFs.
Definition: beam3d.h:119
virtual void giveCompositeExportData(std::vector< VTKPiece > &vtkPieces, IntArray &primaryVarsToExport, IntArray &internalVarsToExport, IntArray cellVarsToExport, TimeStep *tStep)
Definition: beam3d.C:1206
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element&#39;s local coordinates.
Definition: beam3d.C:773
#define OOFEM_ERROR(...)
Definition: error.h:61
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: beam3d.C:905
ElementExtension
Type representing element extension.
void computeInternalForcesFromBoundaryEdgeLoadVectorAtPoint(FloatArray &answer, BoundaryLoad *load, int edge, CharType type, ValueModeType mode, TimeStep *tStep, FloatArray &pointCoords, double ds, bool global)
Definition: beam3d.C:1109
DofManager * ghostNodes[2]
Definition: beam3d.h:91
virtual bool computeGtoLRotationMatrix(FloatMatrix &answer)
Returns transformation matrix from global c.s.
Definition: beam3d.C:281
UnknownType
Type representing particular unknown (its physical meaning).
Definition: unknowntype.h:55
void giveLocationArray(IntArray &locationArray, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const
Returns the location array (array of code numbers) of receiver for given numbering scheme...
Definition: element.C:390
virtual void FiberedCrossSectionInterface_computeStrainVectorInFiber(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *slaveGp, TimeStep *tStep)
Computes full 3d strain vector in element fiber.
Definition: beam3d.C:863
Abstract base class allowing to control the way, how equations are assigned to individual DOFs...
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: beam3d.C:702
void giveCompleteMasterDofIDArray(IntArray &dofIDArray) const
Returns the full dof ID array of receiver.
Definition: dofmanager.C:249
static FEI3dLineLin interp
Geometry interpolator only.
Definition: beam3d.h:78
Class representing a linear line interpolation in 3D.
Definition: fei3dlinelin.h:44
virtual const char * giveInputRecordName() const
Definition: beam3d.h:170
#define ALL_STRAINS
Class representing vector of real numbers.
Definition: floatarray.h: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: beam3d.C:354
This class implements a 2-dimensional beam element with cubic lateral displacement interpolation (rot...
Definition: beam3d.h:74
int subsoilMat
Subsoil material.
Definition: beam3d.h:96
virtual DofManager * giveInternalDofManager(int i) const
Returns i-th internal element dof manager of the receiver.
Definition: beam3d.h:123
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
FloatArray zaxis
Definition: beam3d.h:82
CharType
Definition: chartype.h:87
double kappay
Definition: beam3d.h:80
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual void computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity=NULL)
Computes consistent mass matrix of receiver using numerical integration over element volume...
Definition: beam3d.C:711
int referenceNode
Definition: beam3d.h:81
Class Interface.
Definition: interface.h:82
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: beam3d.C:880
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: beam3d.C:133
virtual void computeBmatrixAt(GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Definition: beam3d.C:88
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &)
Computes interpolation matrix for element unknowns.
Definition: beam3d.C:146
The element interface required by FiberedCrossSection.
Definition: fiberedcs.h:220
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 int testElementExtension(ElementExtension ext)
Tests if the element implements required extension.
Definition: beam3d.h:114
void giveEndForcesVector(FloatArray &answer, TimeStep *tStep)
Definition: beam3d.C:645
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
Definition: beam3d.C:631
virtual void giveBoundaryLocationArray(IntArray &locationArray, const IntArray &bNodes, const UnknownNumberingScheme &s, IntArray *dofIds=NULL)
Returns the location array for the boundary of the element.
Definition: beam3d.h:146
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual void giveInternalDofManDofIDMask(int i, IntArray &answer) const
Returns internal dofmanager dof mask for node.
Definition: beam3d.h:133
#define _IFT_Beam3d_Name
Definition: beam3d.h:46
Interface defining required functionality from associated element.
Definition: winklermodel.h:103
double length
Definition: beam3d.h:80
virtual int giveNumberOfIPForMassMtrxIntegration()
Return desired number of integration points for consistent mass matrix computation, if required.
Definition: beam3d.h:218
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
Definition: beam3d.C:403
Class representing integration point in finite element program.
Definition: gausspoint.h:93
virtual void computeBoundaryEdgeLoadVector(FloatArray &answer, BoundaryLoad *load, int edge, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true)
Computes the contribution of the given load at the given boundary edge.
Definition: beam3d.C:217
Class representing solution step.
Definition: timestep.h:80
This class implements a base beam intented to be a base class for beams based on lagrangian interpola...
Element extension for edge loads.
virtual int giveNumberOfInternalDofManagers() const
Definition: beam3d.h:122
virtual double computeLength()
Computes the length (zero for all but 1D geometries)
Definition: beam3d.C:412

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