OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
truss2d.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 truss2d_h
36 #define truss2d_h
37 
39 
41 
42 #define _IFT_Truss2d_Name "truss2d"
43 #define _IFT_Truss2d_cs "cs"
44 
45 
46 namespace oofem {
47 class FEI2dLineLin;
61 {
62 protected:
63  double length;
64  double pitch;
65  int cs_mode;
66  // array for diffrent cs_modes
67  static FEI2dLineLin interp[3]; // only defined it so far...
68 
69 public:
70  Truss2d(int n, Domain * d);
71  virtual ~Truss2d() { }
72 
73  virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep);
74  virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
75  { computeLumpedMassMatrix(answer, tStep); }
76  virtual int giveLocalCoordinateSystem(FloatMatrix &answer);
77 
78  virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords);
79 
80  virtual int computeNumberOfDofs() { return 4; }
81  virtual void giveDofManDofIDMask(int inode, IntArray &answer) const;
82 
83  virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep);
84  virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep);
85 
86  // characteristic length (for crack band approach)
87  virtual double giveCharacteristicLength(const FloatArray &)
88  { return this->computeLength(); }
89 
90  virtual double computeVolumeAround(GaussPoint *gp);
91 
92  virtual int testElementExtension(ElementExtension ext) { return ( ( ext == Element_EdgeLoadSupport ) ? 1 : 0 ); }
93 
94 #ifdef __OOFEG
95  virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep);
97 #endif
98 
99  // definition & identification
100  virtual const char *giveInputRecordName() const { return _IFT_Truss2d_Name; }
101  virtual const char *giveClassName() const { return "Truss2d"; }
104  virtual Element_Geometry_Type giveGeometryType() const { return EGT_line_1; }
105  virtual integrationDomain giveIntegrationDomain() const { return _Line; }
106  virtual MaterialMode giveMaterialMode() { return _1dMat; }
107  virtual FEInterpolation *giveInterpolation() const;
108 
109 protected:
110  // edge load support
111  void resolveCoordIndices(int &c1, int &c2);
112  virtual void giveEdgeDofMapping(IntArray &answer, int) const;
113  virtual double computeEdgeVolumeAround(GaussPoint *gp, int);
114  virtual int computeLoadLEToLRotationMatrix(FloatMatrix &, int, GaussPoint *gp);
115  virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &, int = 1, int = ALL_STRAINS);
116  virtual void computeBHmatrixAt(GaussPoint *gp, FloatMatrix &);
117  virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &);
118  virtual void computeGaussPoints();
119 
120  virtual double computeLength();
121  double givePitch();
122 };
123 } // end namespace oofem
124 #endif // truss2d_h
Truss2d(int n, Domain *d)
Definition: truss2d.C:59
integrationDomain
Used by integrator class to supply integration points for proper domain to be integrated (Area...
double givePitch()
Definition: truss2d.C:216
Class and object Domain.
Definition: domain.h:115
virtual ~Truss2d()
Definition: truss2d.h:71
Element_Geometry_Type
Enumerative type used to classify element geometry Possible values are: EGT_point - point in space EG...
virtual FEInterpolation * giveInterpolation() const
Definition: truss2d.C:387
Abstract base class for "structural" finite elements with geometrical nonlinearities.
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Definition: truss2d.C:70
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
Definition: truss2d.C:322
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: truss2d.C:285
virtual int computeLoadLEToLRotationMatrix(FloatMatrix &, int, GaussPoint *gp)
Returns transformation matrix from local edge c.s to element local coordinate system of load vector c...
Definition: truss2d.C:360
virtual const char * giveClassName() const
Definition: truss2d.h:101
void resolveCoordIndices(int &c1, int &c2)
Definition: truss2d.C:265
virtual double computeEdgeVolumeAround(GaussPoint *gp, int)
Computes volume related to integration point on local edge.
Definition: truss2d.C:349
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
virtual Element_Geometry_Type giveGeometryType() const
Definition: truss2d.h:104
virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes mass matrix of receiver.
Definition: truss2d.h:74
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
static FEI2dLineLin interp[3]
Definition: truss2d.h:67
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
virtual int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
Definition: truss2d.h:80
virtual double computeLength()
Computes the length (zero for all but 1D geometries)
Definition: truss2d.C:194
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: truss2d.C:395
virtual int testElementExtension(ElementExtension ext)
Tests if the element implements required extension.
Definition: truss2d.h:92
ElementExtension
Type representing element extension.
UnknownType
Type representing particular unknown (its physical meaning).
Definition: unknowntype.h:55
virtual void computeBHmatrixAt(GaussPoint *gp, FloatMatrix &)
Computes a matrix which, multiplied by the column matrix of nodal displacements, gives the displaceme...
Definition: truss2d.C:103
virtual void giveEdgeDofMapping(IntArray &answer, int) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
Definition: truss2d.C:329
Class representing a 2d line with linear interpolation.
Definition: fei2dlinelin.h:46
#define ALL_STRAINS
#define _IFT_Truss2d_Name
Definition: truss2d.h:42
double length
Definition: truss2d.h:63
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: truss2d.C:111
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual const char * giveInputRecordName() const
Definition: truss2d.h:100
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 MaterialMode giveMaterialMode()
Returns material mode for receiver integration points.
Definition: truss2d.h:106
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual integrationDomain giveIntegrationDomain() const
Returns integration domain for receiver, used to initialize integration point over receiver volume...
Definition: truss2d.h:105
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: truss2d.C:185
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element's local coordinates.
Definition: truss2d.C:164
virtual double giveCharacteristicLength(const FloatArray &)
Returns the size of element in the given direction, in some cases adjusted (e.g.
Definition: truss2d.h:87
This class implements a two-node truss bar element for two-dimensional analysis.
Definition: truss2d.h:60
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
Definition: truss2d.C:302
double pitch
Definition: truss2d.h:64
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
Definition: truss2d.C:440
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: truss2d.C:315
Class representing integration point in finite element program.
Definition: gausspoint.h:93
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &)
Computes interpolation matrix for element unknowns.
Definition: truss2d.C:144
Class representing solution step.
Definition: timestep.h:80
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Definition: truss2d.C:123
Element extension for edge loads.
virtual int giveLocalCoordinateSystem(FloatMatrix &answer)
Returns local coordinate system of receiver Required by material models with ortho- and anisotrophy...
Definition: truss2d.C:240

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