OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
libeam2dnl.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 libeam2dnl_h
36 #define libeam2dnl_h
37 
38 #include "../sm/Elements/nlstructuralelement.h"
39 #include "../sm/CrossSections/layeredcrosssection.h"
40 
41 #define _IFT_LIBeam2dNL_Name "libeam2dNL"
42 
43 namespace oofem {
50 {
51 protected:
52  double pitch, length;
53 
54 public:
55  LIBeam2dNL(int n, Domain * d);
56  virtual ~LIBeam2dNL() { }
57 
58  virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep);
59  virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
60  { computeLumpedMassMatrix(answer, tStep); }
61  virtual bool computeGtoLRotationMatrix(FloatMatrix &answer);
62  virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep);
63 
64  virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords);
65 
66  virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep);
67  virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep);
68 
69  // layered cross section support functions
70  virtual void computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain,
71  GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep);
72 
74 
75  virtual int computeNumberOfDofs() { return 6; }
76  virtual void giveDofManDofIDMask(int inode, IntArray &) const;
77  virtual double computeVolumeAround(GaussPoint *gp);
78 
79  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
80 
81  // definition & identification
82  virtual const char *giveInputRecordName() const { return _IFT_LIBeam2dNL_Name; }
83  virtual const char *giveClassName() const { return "LIBeam2dNL"; }
85 
86 #ifdef __OOFEG
87  virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep);
89 #endif
90 
91  virtual integrationDomain giveIntegrationDomain() const { return _Line; }
92  virtual MaterialMode giveMaterialMode() { return _2dBeam; }
93  virtual Element_Geometry_Type giveGeometryType() const { return EGT_line_1; }
94 
95 protected:
96  // edge load support
97  virtual void giveEdgeDofMapping(IntArray &answer, int) const;
98  virtual double computeEdgeVolumeAround(GaussPoint *, int);
100  virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer);
101  virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode);
102 
103  virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int = 1, int = ALL_STRAINS);
104  // nonlinear part of geometrical eqs. for i-th component of strain vector.
105  void computeNLBMatrixAt(FloatMatrix &answer, GaussPoint *gp, int);
106  virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer);
107  virtual void computeGaussPoints();
108  virtual double computeLength();
109  double givePitch();
110 };
111 } // end namespace oofem
112 #endif // libeam2dnl_h
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Definition: libeam2dnl.C:295
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
integrationDomain
Used by integrator class to supply integration points for proper domain to be integrated (Area...
virtual bool computeGtoLRotationMatrix(FloatMatrix &answer)
Returns transformation matrix from global c.s.
Definition: libeam2dnl.C:342
Class and object Domain.
Definition: domain.h:115
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.
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual MaterialMode giveMaterialMode()
Returns material mode for receiver integration points.
Definition: libeam2dnl.h:92
virtual int computeLoadLEToLRotationMatrix(FloatMatrix &, int, GaussPoint *)
Returns transformation matrix from local edge c.s to element local coordinate system of load vector c...
Definition: libeam2dnl.C:559
#define _IFT_LIBeam2dNL_Name
Definition: libeam2dnl.h:41
virtual const char * giveInputRecordName() const
Definition: libeam2dnl.h:82
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: libeam2dnl.C:491
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
virtual ~LIBeam2dNL()
Definition: libeam2dnl.h:56
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: libeam2dnl.C:583
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element's local coordinates.
Definition: libeam2dnl.C:436
double givePitch()
Definition: libeam2dnl.C:470
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
Definition: libeam2dnl.C:401
virtual void giveEdgeDofMapping(IntArray &answer, int) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
Definition: libeam2dnl.C:498
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: libeam2dnl.C:283
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: libeam2dnl.C:574
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
Definition: libeam2dnl.C:312
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: libeam2dnl.C:377
virtual Element_Geometry_Type giveGeometryType() const
Returns the element geometry type.
Definition: libeam2dnl.h:93
UnknownType
Type representing particular unknown (its physical meaning).
Definition: unknowntype.h:55
virtual void computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep)
Computes full 3D strain vector in element layer.
Definition: libeam2dnl.C:408
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: libeam2dnl.C:63
This class implements a 2-dimensional Linear Isoparametric Mindlin theory beam element, with reduced integration.
Definition: libeam2dnl.h:49
#define ALL_STRAINS
void computeNLBMatrixAt(FloatMatrix &answer, GaussPoint *gp, int)
Definition: libeam2dnl.C:106
LIBeam2dNL(int n, Domain *d)
Definition: libeam2dnl.C:55
Class representing vector of real numbers.
Definition: floatarray.h:82
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Definition: libeam2dnl.C:74
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
virtual double computeLength()
Computes the length (zero for all but 1D geometries)
Definition: libeam2dnl.C:452
virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes mass matrix of receiver.
Definition: libeam2dnl.h:59
virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes the initial stiffness matrix of receiver.
Definition: libeam2dnl.C:256
virtual const char * giveClassName() const
Definition: libeam2dnl.h:83
virtual int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
Definition: libeam2dnl.h:75
Class representing the general Input Record.
Definition: inputrecord.h:101
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: libeam2dnl.C:394
Class Interface.
Definition: interface.h:82
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
Definition: libeam2dnl.C:609
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: libeam2dnl.C:367
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
virtual double computeEdgeVolumeAround(GaussPoint *, int)
Computes volume related to integration point on local edge.
Definition: libeam2dnl.C:520
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: libeam2dnl.C:429
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
Definition: libeam2dnl.C:532
the oofem namespace is to define a context or scope in which all oofem names are defined.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
virtual integrationDomain giveIntegrationDomain() const
Returns integration domain for receiver, used to initialize integration point over receiver volume...
Definition: libeam2dnl.h:91

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