OOFEM 3.0
Loading...
Searching...
No Matches
springelement.C
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
36#include "floatmatrix.h"
37#include "intarray.h"
38#include "floatarray.h"
39#include "classfactory.h"
40#include "parametermanager.h"
41#include "paramkey.h"
42
43#ifdef __OOFEG
44 #include "oofeggraphiccontext.h"
45#endif
46
47namespace oofem {
49
54
55SpringElement :: SpringElement(int n, Domain *aDomain) : StructuralElement(n, aDomain)
56{
58 springConstant = 0.0;
59}
60
61void
62SpringElement :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
63{
64 /* spring stiffness matrix in local coordinate system (along orientation axis) */
65 answer.resize(2, 2);
66 answer.at(1, 1) = answer.at(2, 2) = this->springConstant;
67 answer.at(1, 2) = answer.at(2, 1) = -this->springConstant;
68}
69
70
71void
72SpringElement :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
73{
74 double f = this->computeSpringInternalForce(tStep);
75 answer.resize(2);
76 answer.at(1) = -f;
77 answer.at(2) = f;
78}
79
80
81bool
82SpringElement :: computeGtoLRotationMatrix(FloatMatrix &answer)
83{
84 /*
85 * Spring is defined as 1D element along orientation axis (or around orientation axis for torsional springs)
86 * The transformation from local (1d) to global system typically expand dimensions
87 */
88 if ( ( this->mode == SE_1D_SPRING ) || ( this->mode == SE_2D_TORSIONALSPRING_XZ ) ) {
89 answer.resize(2, 2);
90 answer.at(1, 1) = answer.at(2, 2) = 1.0;
91 } else if ( this->mode == SE_2D_SPRING_XY ) {
92 answer.resize(2, 4);
93 answer.at(1, 1) = this->dir.at(1);
94 answer.at(1, 2) = this->dir.at(2);
95 answer.at(2, 3) = this->dir.at(1);
96 answer.at(2, 4) = this->dir.at(2);
97 } else if ( this->mode == SE_2D_SPRING_XZ ) {
98 answer.resize(2, 4);
99 answer.at(1, 1) = this->dir.at(1);
100 answer.at(1, 2) = this->dir.at(3);
101 answer.at(2, 3) = this->dir.at(1);
102 answer.at(2, 4) = this->dir.at(3);
103 } else if ( ( this->mode == SE_3D_SPRING ) || ( this->mode == SE_3D_TORSIONALSPRING ) ) {
104 answer.resize(2, 6);
105 answer.at(1, 1) = this->dir.at(1);
106 answer.at(1, 2) = this->dir.at(2);
107 answer.at(1, 3) = this->dir.at(3);
108 answer.at(2, 4) = this->dir.at(1);
109 answer.at(2, 5) = this->dir.at(2);
110 answer.at(2, 6) = this->dir.at(3);
111 }
112 return 1;
113}
114
115
116void
117SpringElement :: giveDofManDofIDMask(int inode, IntArray &answer) const
118{
119 if ( this->mode == SE_1D_SPRING ) {
120 answer = {D_u};
121 } else if ( this->mode == SE_2D_SPRING_XY ) {
122 answer = {D_u, D_v};
123 } else if ( this->mode == SE_2D_SPRING_XZ ) {
124 answer = {D_u, D_w};
125 } else if ( this->mode == SE_2D_TORSIONALSPRING_XZ ) {
126 answer = {R_v};
127 } else if ( this->mode == SE_3D_SPRING ) {
128 answer = {D_u, D_v, D_w};
129 } else if ( this->mode == SE_3D_TORSIONALSPRING ) {
130 answer = {R_u, R_v, R_w};
131 }
132}
133
134double
135SpringElement :: computeSpringInternalForce(TimeStep *tStep)
136{
137 FloatArray u;
138 this->computeVectorOf(VM_Total, tStep, u);
139 return ( this->springConstant * ( u.at(2) - u.at(1) ) );
140}
141
142void
144{
145 answer.resize(2, 2);
146 answer.at(1,1)=answer.at(2,2) = this->mass/2.0;
147 answer.at(1,2)=answer.at(2,1) = 0.0;
148}
149
150int
151SpringElement :: computeNumberOfGlobalDofs()
152{
153 if ( ( this->mode == SE_1D_SPRING ) || ( this->mode == SE_2D_TORSIONALSPRING_XZ ) ) {
154 return 2;
155 } else if ( ( this->mode == SE_2D_SPRING_XY ) || ( this->mode == SE_2D_SPRING_XZ ) ) {
156 return 4;
157 } else if ( ( this->mode == SE_3D_SPRING ) || ( this->mode == SE_3D_TORSIONALSPRING ) ) {
158 return 6;
159 }
160 return 0;
161}
162
163
164void
165SpringElement :: initializeFrom(InputRecord &ir, int priority)
166{
167 ParameterManager &ppm = this->giveDomain()->elementPPM;
168
169 StructuralElement :: initializeFrom(ir, priority);
170 int _mode;
171 bool modeFlag, dirFlag;
172 PM_UPDATE_PARAMETER_AND_REPORT(_mode, ppm, ir, this->number, IPK_SpringElement_mode, priority, modeFlag) ;
173 if (modeFlag) {
174 this->mode = ( SpringElementType ) _mode;
175 }
177 PM_UPDATE_PARAMETER(mass, ppm, ir, this->number, IPK_SpringElement_mass, priority) ;
178 PM_UPDATE_PARAMETER_AND_REPORT(dir, ppm, ir, this->number, IPK_SpringElement_orientation, priority, dirFlag) ;
179 if (dirFlag) {
180 this->dir.normalize();
181 }
182}
183
184void SpringElement :: printOutputAt(FILE *File, TimeStep *tStep)
185{
186 fprintf(File, "spring element %d (%8d) :\n", this->giveLabel(), this->giveNumber() );
187 fprintf(File, " spring force or moment %.4e", this->computeSpringInternalForce(tStep) );
188 fprintf(File, "\n");
189}
190} // end namespace oofem
#define REGISTER_Element(class)
int numberOfDofMans
Number of dofmanagers.
Definition element.h:136
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Definition element.C:103
int giveLabel() const
Definition element.h:1125
Domain * giveDomain() const
Definition femcmpnn.h:97
int number
Component number.
Definition femcmpnn.h:77
int giveNumber() const
Definition femcmpnn.h:104
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
double at(std::size_t i, std::size_t j) const
double springConstant
The longitudinal spring constant [Force/Length], torsional spring constant [Force*Length/Radians].
static ParamKey IPK_SpringElement_mass
static ParamKey IPK_SpringElement_orientation
double computeSpringInternalForce(TimeStep *tStep)
SpringElementType mode
Mode.
static ParamKey IPK_SpringElement_springConstant
double mass
total mass of the spring; to be distributed to nodes
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...
static ParamKey IPK_SpringElement_mode
void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep) override
StructuralElement(int n, Domain *d)
#define PM_UPDATE_PARAMETER_AND_REPORT(_val, _pm, _ir, _componentnum, _paramkey, _prio, _flag)
#define PM_UPDATE_PARAMETER(_val, _pm, _ir, _componentnum, _paramkey, _prio)

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