OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
libeam3d2.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 libeam3d2_h
36 #define libeam3d2_h
37 
38 #include "../sm/Elements/nlstructuralelement.h"
39 #include "../sm/CrossSections/fiberedcs.h"
40 
42 
43 #define _IFT_LIBeam3d2_Name "libeam3d2"
44 #define _IFT_LIBeam3d2_refnode "refnode"
45 
46 
47 namespace oofem {
61 {
62 private:
63  double length;
65 
72 
73 public:
74  LIBeam3d2(int n, Domain * d);
75  virtual ~LIBeam3d2() { }
76 
78 
79  virtual void updateYourself(TimeStep *tStep);
80 
81  virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep);
82  virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
83  { computeLumpedMassMatrix(answer, tStep); }
84  virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep);
85  virtual bool computeGtoLRotationMatrix(FloatMatrix &answer);
86 
87  virtual int testElementExtension(ElementExtension ext);
88 
89  virtual int computeNumberOfDofs() { return 12; }
90  virtual void giveDofManDofIDMask(int inode, IntArray &) const;
91  virtual double computeVolumeAround(GaussPoint *gp);
92  virtual int giveLocalCoordinateSystem(FloatMatrix &answer);
93  virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep);
94 
95  virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords);
96 
97  virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep);
98  virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep);
99 
100  // Fibered cross section support functions
102  GaussPoint *slaveGp, TimeStep *tStep);
103 
105 
106  // definition & identification
107  virtual const char *giveInputRecordName() const { return _IFT_LIBeam3d2_Name; }
108  virtual const char *giveClassName() const { return "LIBeam3d2"; }
109  virtual Element_Geometry_Type giveGeometryType() const { return EGT_line_1; }
110 
111  virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj);
112  virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj);
113 
114 #ifdef __OOFEG
115  virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep);
117  virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep);
118 #endif
119 
120  virtual integrationDomain giveIntegrationDomain() const { return _Line; }
121  virtual MaterialMode giveMaterialMode() { return _3dBeam; }
122 
123 protected:
124  // edge load support
125  virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const;
126  virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge);
127  virtual int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp);
128  virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer);
129  virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode);
130 
131  virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int = 1, int = ALL_STRAINS);
132  virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer);
133  virtual void computeGaussPoints();
134  virtual double computeLength();
135 
136  // nonlinearity
137  void updateTempTriad(TimeStep *tStep);
138  void computeSMtrx(FloatMatrix &answer, FloatArray &vec);
139  void computeRotMtrx(FloatMatrix &answer, FloatArray &psi);
140  double giveCurrentLength(TimeStep *tStep);
141  void initForNewStep();
142 };
143 } // end namespace oofem
144 #endif // libeam3d2_h
virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
Definition: libeam3d2.C:327
virtual MaterialMode giveMaterialMode()
Returns material mode for receiver integration points.
Definition: libeam3d2.h:121
integrationDomain
Used by integrator class to supply integration points for proper domain to be integrated (Area...
void computeSMtrx(FloatMatrix &answer, FloatArray &vec)
Definition: libeam3d2.C:506
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: libeam3d2.C:295
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: libeam3d2.C:566
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
Definition: libeam3d2.C:254
Class and object Domain.
Definition: domain.h:115
void initForNewStep()
Initializes receivers state to new time step.
Definition: libeam3d2.C:583
Element_Geometry_Type
Enumerative type used to classify element geometry Possible values are: EGT_point - point in space EG...
Abstract base class for "structural" finite elements with geometrical nonlinearities.
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: libeam3d2.C:117
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Definition: libeam3d2.C:742
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: libeam3d2.C:661
virtual integrationDomain giveIntegrationDomain() const
Returns integration domain for receiver, used to initialize integration point over receiver volume...
Definition: libeam3d2.h:120
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
long StateCounterType
StateCounterType type used to indicate solution state.
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element's local coordinates.
Definition: libeam3d2.C:236
virtual const char * giveClassName() const
Definition: libeam3d2.h:108
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: libeam3d2.C:219
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
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: libeam3d2.C:386
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Definition: libeam3d2.C:130
virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes mass matrix of receiver.
Definition: libeam3d2.h:82
virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
Computes volume related to integration point on local edge.
Definition: libeam3d2.C:346
StateCounterType tempTcCounter
Time stamp of temporary centre triad.
Definition: libeam3d2.h:71
virtual int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
Definition: libeam3d2.h:89
virtual int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp)
Returns transformation matrix from local edge c.s to element local coordinate system of load vector c...
Definition: libeam3d2.C:395
virtual double computeLength()
Computes the length (zero for all but 1D geometries)
Definition: libeam3d2.C:275
ElementExtension
Type representing element extension.
UnknownType
Type representing particular unknown (its physical meaning).
Definition: unknowntype.h:55
virtual bool computeGtoLRotationMatrix(FloatMatrix &answer)
Returns transformation matrix from global c.s.
Definition: libeam3d2.C:184
virtual Element_Geometry_Type giveGeometryType() const
Returns the element geometry type.
Definition: libeam3d2.h:109
#define _IFT_LIBeam3d2_Name
Definition: libeam3d2.h:43
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
Definition: libeam3d2.C:358
#define ALL_STRAINS
FloatMatrix tc
Last equilibrium triad at the centre.
Definition: libeam3d2.h:67
double giveCurrentLength(TimeStep *tStep)
Definition: libeam3d2.C:547
virtual ~LIBeam3d2()
Definition: libeam3d2.h:75
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: libeam3d2.C:649
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
void computeRotMtrx(FloatMatrix &answer, FloatArray &psi)
Definition: libeam3d2.C:476
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
Definition: libeam3d2.C:145
virtual const char * giveInputRecordName() const
Definition: libeam3d2.h:107
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
Definition: libeam3d2.C:688
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: libeam3d2.C:261
virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Compute strain vector of receiver evaluated at given integration point at given time step from elemen...
Definition: libeam3d2.C:525
virtual int giveLocalCoordinateSystem(FloatMatrix &answer)
Returns local coordinate system of receiver Required by material models with ortho- and anisotrophy...
Definition: libeam3d2.C:410
Class Interface.
Definition: interface.h:82
FloatMatrix tempTc
Temporary triad at the centre.
Definition: libeam3d2.h:69
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
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 contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj)
Stores receiver state to output stream.
Definition: libeam3d2.C:593
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes the stiffness matrix of receiver.
Definition: libeam3d2.C:175
the oofem namespace is to define a context or scope in which all oofem names are defined.
This class implements a 3-dimensional Linear Isoparametric Mindlin theory beam element, with reduced integration.
Definition: libeam3d2.h:60
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Definition: libeam3d2.C:71
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
Definition: libeam3d2.C:229
LIBeam3d2(int n, Domain *d)
Definition: libeam3d2.C:60
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj)
Restores the receiver state previously written in stream.
Definition: libeam3d2.C:612
Class representing integration point in finite element program.
Definition: gausspoint.h:93
void updateTempTriad(TimeStep *tStep)
Definition: libeam3d2.C:447
void FiberedCrossSectionInterface_computeStrainVectorInFiber(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *slaveGp, TimeStep *tStep)
Computes full 3d strain vector in element fiber.
Definition: libeam3d2.C:632
Class representing solution step.
Definition: timestep.h:80
virtual int testElementExtension(ElementExtension ext)
Tests if the element implements required extension.
Definition: libeam3d2.C:268

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