OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
springelement.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 springelement_h
36 #define springelement_h
37 
38 #include "../sm/Elements/structuralelement.h"
39 
41 
42 #define _IFT_SpringElement_Name "spring"
43 #define _IFT_SpringElement_mode "mode"
44 #define _IFT_SpringElement_orientation "orientation"
45 #define _IFT_SpringElement_springConstant "k"
46 #define _IFT_SpringElement_mass "m"
47 
49 
50 namespace oofem {
57 {
58 public:
67  };
68 
69 protected:
73  double mass;
81 
82 public:
83  SpringElement(int n, Domain * d);
84  virtual ~SpringElement() { }
85 
86  virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep);
87  virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
88  { computeLumpedMassMatrix(answer, tStep); }
89  virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep);
90  virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
91  { answer.clear(); }
92 
93  virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord = 0);
94 
95  virtual int computeNumberOfDofs() { return 2; }
96  virtual int computeNumberOfGlobalDofs();
97 
98  virtual void giveDofManDofIDMask(int inode, IntArray &answer) const;
99 
100  virtual void updateInternalState(TimeStep *tStep) { }
101  virtual void updateYourself(TimeStep *tStep) { }
102  virtual int checkConsistency() { return 1; }
103  virtual void printOutputAt(FILE *file, TimeStep *tStep);
104  virtual bool isCast(TimeStep *tStep) {return true;}
105 
106 #ifdef __OOFEG
107  //void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep);
108  //void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType);
109  //void drawScalar(oofegGraphicContext &gc, TimeStep *tStep);
110 #endif
111 
112  // definition & identification
113  virtual const char *giveInputRecordName() const { return _IFT_SpringElement_Name; }
114  virtual const char *giveClassName() const { return "SpringElement"; }
116  virtual Element_Geometry_Type giveGeometryType() const { return EGT_point; }
117 
118 protected:
119  virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
120  { answer.clear(); }
122  { answer.clear(); }
123 
124  virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer,
125  int lowerIndx = 1, int upperIndx = ALL_STRAINS)
126  { answer.clear(); }
127  virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) { answer.clear(); }
128  virtual bool computeGtoLRotationMatrix(FloatMatrix &answer);
129  double computeSpringInternalForce(TimeStep *tStep);
130 };
131 } // end namespace oofem
132 #endif // springelement_h
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
< 2D spring element in xz plane, requires D_u and D_w DOFs in each node (orientation vector should be...
Definition: springelement.h:64
Class and object Domain.
Definition: domain.h:115
virtual bool isCast(TimeStep *tStep)
Element_Geometry_Type
Enumerative type used to classify element geometry Possible values are: EGT_point - point in space EG...
virtual int computeNumberOfGlobalDofs()
Computes the total number of element&#39;s global dofs.
This class implements a simple spring element.
Definition: springelement.h:56
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.
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes mass matrix of receiver.
Definition: springelement.h:87
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
Returns equivalent nodal forces vectors.
Definition: springelement.C:65
virtual bool computeGtoLRotationMatrix(FloatMatrix &answer)
Returns transformation matrix from global c.s.
Definition: springelement.C:75
SpringElement(int n, Domain *d)
Definition: springelement.C:48
Class implementing an array of integers.
Definition: intarray.h:61
2D spring element in xy plane, requires D_u and D_v DOFs in each node (orientation vector should be i...
Definition: springelement.h:62
MatResponseMode
Describes the character of characteristic material matrix.
double mass
total mass of the spring; to be distributed to nodes
Definition: springelement.h:73
virtual ~SpringElement()
Definition: springelement.h:84
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Abstract base class for all "structural" finite elements.
virtual Element_Geometry_Type giveGeometryType() const
Returns the element geometry type.
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
double computeSpringInternalForce(TimeStep *tStep)
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
1D spring element along x-axis.
Definition: springelement.h:61
#define _IFT_SpringElement_Name
Definition: springelement.h:42
SpringElementType mode
Mode.
Definition: springelement.h:80
#define ALL_STRAINS
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual void updateInternalState(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes numerically stiffness matrix of receiver.
Definition: springelement.C:55
virtual int computeNumberOfDofs()
Computes or simply returns total number of element&#39;s local DOFs.
Definition: springelement.h:95
virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes initial stress matrix for linear stability problem.
Definition: springelement.h:90
FloatArray dir
Orientation vector.
Definition: springelement.h:78
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
SpringElementType
Defines type of spring element (longitudinal/rotational) spring.
Definition: springelement.h:60
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
virtual const char * giveInputRecordName() const
Class representing the general Input Record.
Definition: inputrecord.h:101
double springConstant
The longitudinal spring constant [Force/Length], torsional spring constant [Force*Length/Radians].
Definition: springelement.h:71
3D spring element in space, requires D_u, D_v, and D_w DOFs in each node.
Definition: springelement.h:65
3D torsional spring in space, requires R_u, R_v, and R_w DOFs in each node.
Definition: springelement.h:66
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual const char * giveClassName() const
void clear()
Sets size of receiver to be an empty matrix. It will have zero rows and zero columns size...
Definition: floatmatrix.h:516
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int lowerIndx=1, int upperIndx=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
virtual int checkConsistency()
Performs consistency check.

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:31 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011