OOFEM 3.0
Loading...
Searching...
No Matches
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 - 2025 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
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
50namespace oofem {
51class ParamKey;
58{
59public:
69
70protected:
74 double mass=0.0;
82
87
88public:
89 SpringElement(int n, Domain * d);
90 virtual ~SpringElement() { }
91
92 void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep) override;
93 void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep) override
94 { computeLumpedMassMatrix(answer, tStep); }
95 void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) override;
96 void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep) override
97 { answer.clear(); }
98
99 void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord = 0) override;
100
101 int computeNumberOfDofs() override { return 2; }
102 int computeNumberOfGlobalDofs() override;
103
104 void giveDofManDofIDMask(int inode, IntArray &answer) const override;
105
106 void updateInternalState(TimeStep *tStep) override { }
107 void updateYourself(TimeStep *tStep) override { }
108 int checkConsistency() override { return 1; }
109 void printOutputAt(FILE *file, TimeStep *tStep) override;
110 bool isCast(TimeStep *tStep) override { return true; }
111
112#ifdef __OOFEG
113 //void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep);
114 //void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType);
115 //void drawScalar(oofegGraphicContext &gc, TimeStep *tStep);
116#endif
117
118 // definition & identification
119 const char *giveInputRecordName() const override { return _IFT_SpringElement_Name; }
120 const char *giveClassName() const override { return "SpringElement"; }
121 void initializeFrom(InputRecord &ir, int priority) override;
122 Element_Geometry_Type giveGeometryType() const override { return EGT_point; }
123
124protected:
125 void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
126 { answer.clear(); }
127 void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
128 { answer.clear(); }
129
131 int lowerIndx = 1, int upperIndx = ALL_STRAINS) override
132 { answer.clear(); }
133 void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override { answer.clear(); }
134 bool computeGtoLRotationMatrix(FloatMatrix &answer) override;
136};
137} // end namespace oofem
138#endif // springelement_h
*Sets size of receiver to be an empty matrix It will have zero rows and zero columns size void clear()
void printOutputAt(FILE *file, TimeStep *tStep) override
int computeNumberOfDofs() override
double springConstant
The longitudinal spring constant [Force/Length], torsional spring constant [Force*Length/Radians].
static ParamKey IPK_SpringElement_mass
SpringElement(int n, Domain *d)
static ParamKey IPK_SpringElement_orientation
double computeSpringInternalForce(TimeStep *tStep)
int computeNumberOfGlobalDofs() override
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
SpringElementType mode
Mode.
void updateInternalState(TimeStep *tStep) override
int checkConsistency() override
void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int lowerIndx=1, int upperIndx=ALL_STRAINS) override
bool isCast(TimeStep *tStep) override
void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep) override
static ParamKey IPK_SpringElement_springConstant
void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0) override
void giveDofManDofIDMask(int inode, IntArray &answer) const override
double mass
total mass of the spring; to be distributed to nodes
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
SpringElementType
Defines type of spring element (longitudinal/rotational) spring.
@ SE_1D_SPRING
1D spring element along x-axis.
@ SE_3D_TORSIONALSPRING
3D torsional spring in space, requires R_u, R_v, and R_w DOFs in each node.
@ SE_2D_TORSIONALSPRING_XZ
< 2D spring element in xz plane, requires D_u and D_w DOFs in each node (orientation vector should be...
@ SE_3D_SPRING
3D spring element in space, requires D_u, D_v, and D_w DOFs in each node.
@ SE_2D_SPRING_XY
2D spring element in xy plane, requires D_u and D_v DOFs in each node (orientation vector should be i...
void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep) override
const char * giveClassName() const override
bool computeGtoLRotationMatrix(FloatMatrix &answer) override
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
void initializeFrom(InputRecord &ir, int priority) override
const char * giveInputRecordName() const override
void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) override
void updateYourself(TimeStep *tStep) override
static ParamKey IPK_SpringElement_mode
void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep) override
Element_Geometry_Type giveGeometryType() const override
StructuralElement(int n, Domain *d)
#define _IFT_SpringElement_Name
#define ALL_STRAINS

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011