OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
libeam3dnl.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 libeam3dnl_h
36 #define libeam3dnl_h
37 
38 #include "../sm/Elements/nlstructuralelement.h"
39 
41 
42 #define _IFT_LIBeam3dNL_Name "libeam3dnl"
43 #define _IFT_LIBeam3dNL_refnode "refnode"
44 
45 
46 namespace oofem {
55 {
56 private:
58  double l0;
61  // curvature at the centre
62  // FloatArray kappa;
69 
70 public:
71  LIBeam3dNL(int n, Domain * d);
72  virtual ~LIBeam3dNL() { }
73 
74  virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep);
75  virtual void computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity = NULL)
76  { computeLumpedMassMatrix(answer, tStep); }
77  virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep);
78  //int computeGtoLRotationMatrix(FloatMatrix &answer);
79  //void computeInitialStressMatrix(FloatMatrix& answer, TimeStep* tStep);
80 
81  virtual int computeNumberOfDofs() { return 12; }
82  virtual void giveDofManDofIDMask(int inode, IntArray &) const;
83  virtual double computeVolumeAround(GaussPoint *gp);
84  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
85  virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords);
86 
87  virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep);
88  virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep);
89 
90  // definition & identification
91  virtual const char *giveInputRecordName() const { return _IFT_LIBeam3dNL_Name; }
92  virtual const char *giveClassName() const { return "LIBeam3dNL"; }
94  virtual Element_Geometry_Type giveGeometryType() const { return EGT_line_1; }
95 
96 #ifdef __OOFEG
97  virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep);
99 #endif
100 
101  virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep);
102  virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord = 0);
103 
104  virtual integrationDomain giveIntegrationDomain() const { return _Line; }
105  virtual MaterialMode giveMaterialMode() { return _3dBeam; }
106 
107 
108 protected:
109  // edge load support
110  virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const;
111  virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge);
112  virtual int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp);
113  virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer);
114  virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode);
115 
116  virtual void updateYourself(TimeStep *tStep);
117  virtual void initForNewStep();
118  virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int lowerIndx, int upperIndx)
119  { OOFEM_ERROR("not implemented"); }
120  //int computeGtoLRotationMatrix(FloatMatrix& answer);
121 
122  virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer);
123  virtual void computeGaussPoints();
124  virtual double computeLength();
125  //double givePitch();
126  virtual int giveLocalCoordinateSystem(FloatMatrix &answer);
133  void updateTempTriad(TimeStep *tStep);
139  void computeTempCurv(FloatArray &answer, TimeStep *tStep);
140 
146  void computeSMtrx(FloatMatrix &answer, FloatArray &vec);
153  void computeRotMtrx(FloatMatrix &answer, FloatArray &psi);
159  void computeXMtrx(FloatMatrix &answer, TimeStep *tStep);
165  void computeXdVector(FloatArray &answer, TimeStep *tStep);
166 };
167 } // end namespace oofem
168 #endif // libeam3dnl_h
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 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: libeam3dnl.C:652
Class and object Domain.
Definition: domain.h:115
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
Evaluates nodal representation of real internal forces.
Definition: libeam3dnl.C:199
FloatMatrix tc
Last equilibrium triad at the centre.
Definition: libeam3dnl.h:60
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
Definition: libeam3dnl.C:623
Element_Geometry_Type
Enumerative type used to classify element geometry Possible values are: EGT_point - point in space EG...
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: libeam3dnl.C:337
Abstract base class for "structural" finite elements with geometrical nonlinearities.
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Definition: libeam3dnl.C:405
double l0
Initial length.
Definition: libeam3dnl.h:58
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
virtual void initForNewStep()
Initializes receivers state to new time step.
Definition: libeam3dnl.C:692
long StateCounterType
StateCounterType type used to indicate solution state.
virtual const char * giveClassName() const
Definition: libeam3dnl.h:92
void computeRotMtrx(FloatMatrix &answer, FloatArray &psi)
Evaluates the rotation matrix for large rotations according to Rodrigues formula for given pseudovect...
Definition: libeam3dnl.C:87
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
Definition: libeam3dnl.C:516
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: libeam3dnl.C:560
FloatMatrix tempTc
Temporary triad at the centre.
Definition: libeam3dnl.h:64
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
Definition: libeam3dnl.C:429
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
Definition: libeam3dnl.C:774
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: libeam3dnl.C:676
#define _IFT_LIBeam3dNL_Name
Definition: libeam3dnl.h:42
virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
Definition: libeam3dnl.C:530
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element's local coordinates.
Definition: libeam3dnl.C:498
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes the stiffness matrix of receiver.
Definition: libeam3dnl.C:251
void updateTempTriad(TimeStep *tStep)
Updates the temporary triad at the centre to the state identified by given solution step...
Definition: libeam3dnl.C:117
virtual const char * giveInputRecordName() const
Definition: libeam3dnl.h:91
#define OOFEM_ERROR(...)
Definition: error.h:61
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: libeam3dnl.C:474
virtual double computeLength()
Computes the length (zero for all but 1D geometries)
Definition: libeam3dnl.C:385
UnknownType
Type representing particular unknown (its physical meaning).
Definition: unknowntype.h:55
void computeTempCurv(FloatArray &answer, TimeStep *tStep)
Compute the temporary curvature at the centre to the state identified by given solution step...
Definition: libeam3dnl.C:702
void computeSMtrx(FloatMatrix &answer, FloatArray &vec)
Evaluates the S matrix from given vector vec.
Definition: libeam3dnl.C:68
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: libeam3dnl.C:349
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: libeam3dnl.C:748
This class implements a 3-dimensional Linear Isoparametric Mindlin theory beam element, with reduced integration.
Definition: libeam3dnl.h:54
virtual ~LIBeam3dNL()
Definition: libeam3dnl.h:72
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
Definition: libeam3dnl.h:81
virtual integrationDomain giveIntegrationDomain() const
Returns integration domain for receiver, used to initialize integration point over receiver volume...
Definition: libeam3dnl.h:104
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
LIBeam3dNL(int n, Domain *d)
Definition: libeam3dnl.C:56
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual MaterialMode giveMaterialMode()
Returns material mode for receiver integration points.
Definition: libeam3dnl.h:105
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: libeam3dnl.C:464
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: libeam3dnl.C:523
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: libeam3dnl.C:146
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
Definition: libeam3dnl.C:491
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: libeam3dnl.C:667
void computeXdVector(FloatArray &answer, TimeStep *tStep)
Computes x_21' vector for given solution state.
Definition: libeam3dnl.C:233
Load is base abstract class for all loads.
Definition: load.h:61
int referenceNode
Reference node.
Definition: libeam3dnl.h:68
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
Computes volume related to integration point on local edge.
Definition: libeam3dnl.C:548
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: libeam3dnl.h:75
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int lowerIndx, int upperIndx)
Computes the geometrical matrix of receiver in given integration point.
Definition: libeam3dnl.h:118
virtual Element_Geometry_Type giveGeometryType() const
Returns the element geometry type.
Definition: libeam3dnl.h:94
Class representing integration point in finite element program.
Definition: gausspoint.h:93
StateCounterType tempTcCounter
Time stamp of temporary centre triad.
Definition: libeam3dnl.h:66
Class representing solution step.
Definition: timestep.h:80
void computeXMtrx(FloatMatrix &answer, TimeStep *tStep)
Computes X mtrx at given solution state.
Definition: libeam3dnl.C:174

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