OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
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 - 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 #include "../sm/Elements/springelement.h"
36 #include "floatmatrix.h"
37 #include "intarray.h"
38 #include "floatarray.h"
39 #include "classfactory.h"
40 
41 #ifdef __OOFEG
42  #include "oofeggraphiccontext.h"
43 #endif
44 
45 namespace oofem {
46 REGISTER_Element(SpringElement);
47 
49 {
50  numberOfDofMans = 2;
51  springConstant = 0.0;
52 }
53 
54 void
56 {
57  /* spring stiffness matrix in local coordinate system (along orientation axis) */
58  answer.resize(2, 2);
59  answer.at(1, 1) = answer.at(2, 2) = this->springConstant;
60  answer.at(1, 2) = answer.at(2, 1) = -this->springConstant;
61 }
62 
63 
64 void
65 SpringElement :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
66 {
67  double f = this->computeSpringInternalForce(tStep);
68  answer.resize(2);
69  answer.at(1) = -f;
70  answer.at(2) = f;
71 }
72 
73 
74 bool
76 {
77  /*
78  * Spring is defined as 1D element along orientation axis (or around orientation axis for torsional springs)
79  * The transformation from local (1d) to global system typically expand dimensions
80  */
81  if ( ( this->mode == SE_1D_SPRING ) || ( this->mode == SE_2D_TORSIONALSPRING_XZ ) ) {
82  answer.resize(2, 2);
83  answer.at(1, 1) = answer.at(2, 2) = 1.0;
84  } else if ( this->mode == SE_2D_SPRING_XY ) {
85  answer.resize(2, 4);
86  answer.at(1, 1) = this->dir.at(1);
87  answer.at(1, 2) = this->dir.at(2);
88  answer.at(2, 3) = this->dir.at(1);
89  answer.at(2, 4) = this->dir.at(2);
90  } else if ( this->mode == SE_2D_SPRING_XZ ) {
91  answer.resize(2, 4);
92  answer.at(1, 1) = this->dir.at(1);
93  answer.at(1, 2) = this->dir.at(3);
94  answer.at(2, 3) = this->dir.at(1);
95  answer.at(2, 4) = this->dir.at(3);
96  } else if ( ( this->mode == SE_3D_SPRING ) || ( this->mode == SE_3D_TORSIONALSPRING ) ) {
97  answer.resize(2, 6);
98  answer.at(1, 1) = this->dir.at(1);
99  answer.at(1, 2) = this->dir.at(2);
100  answer.at(1, 3) = this->dir.at(3);
101  answer.at(2, 4) = this->dir.at(1);
102  answer.at(2, 5) = this->dir.at(2);
103  answer.at(2, 6) = this->dir.at(3);
104  }
105  return 1;
106 }
107 
108 
109 void
111 {
112  if ( this->mode == SE_1D_SPRING ) {
113  answer = {D_u};
114  } else if ( this->mode == SE_2D_SPRING_XY ) {
115  answer = {D_u, D_v};
116  } else if ( this->mode == SE_2D_SPRING_XZ ) {
117  answer = {D_u, D_w};
118  } else if ( this->mode == SE_2D_TORSIONALSPRING_XZ ) {
119  answer = {R_v};
120  } else if ( this->mode == SE_3D_SPRING ) {
121  answer = {D_u, D_v, D_w};
122  } else if ( this->mode == SE_3D_TORSIONALSPRING ) {
123  answer = {R_u, R_v, R_w};
124  }
125 }
126 
127 double
129 {
130  FloatArray u;
131  this->computeVectorOf(VM_Total, tStep, u);
132  return ( this->springConstant * ( u.at(2) - u.at(1) ) );
133 }
134 
135 void
137 {
138  answer.resize(2, 2);
139  answer.at(1,1)=answer.at(2,2) = this->mass/2.0;
140  answer.at(1,2)=answer.at(2,1) = 0.0;
141 }
142 
143 int
145 {
146  if ( ( this->mode == SE_1D_SPRING ) || ( this->mode == SE_2D_TORSIONALSPRING_XZ ) ) {
147  return 2;
148  } else if ( ( this->mode == SE_2D_SPRING_XY ) || ( this->mode == SE_2D_SPRING_XZ ) ) {
149  return 4;
150  } else if ( ( this->mode == SE_3D_SPRING ) || ( this->mode == SE_3D_TORSIONALSPRING ) ) {
151  return 6;
152  }
153  return 0;
154 }
155 
156 
159 {
160  IRResultType result; // Required by IR_GIVE_FIELD macro
161 
162  int _mode;
165  this->mass = 0.0;
167 
168  this->mode = ( SpringElementType ) _mode;
169  if ( mode != SE_1D_SPRING ) {
171  this->dir.normalize();
172  }
174 }
175 
176 void SpringElement :: printOutputAt(FILE *File, TimeStep *tStep)
177 {
178  fprintf(File, "spring element %d (%8d) :\n", this->giveLabel(), this->giveNumber() );
179  fprintf(File, " spring force or moment %.4e", this->computeSpringInternalForce(tStep) );
180  fprintf(File, "\n");
181 }
182 } // end namespace oofem
< 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
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Returns local vector of unknowns.
Definition: element.C:86
virtual int computeNumberOfGlobalDofs()
Computes the total number of element&#39;s global dofs.
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
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
#define _IFT_SpringElement_mass
Definition: springelement.h:46
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 void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Abstract base class for all "structural" finite elements.
REGISTER_Element(LSpace)
#define _IFT_SpringElement_orientation
Definition: springelement.h:44
double computeSpringInternalForce(TimeStep *tStep)
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
1D spring element along x-axis.
Definition: springelement.h:61
SpringElementType mode
Mode.
Definition: springelement.h:80
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
Class representing vector of real numbers.
Definition: floatarray.h:82
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
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
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1358
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
#define IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:78
int giveLabel() const
Definition: element.h:1055
#define _IFT_SpringElement_mode
Definition: springelement.h:43
the oofem namespace is to define a context or scope in which all oofem names are defined.
#define IR_GIVE_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:69
int giveNumber() const
Definition: femcmpnn.h:107
double normalize()
Normalizes receiver.
Definition: floatarray.C:828
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Class representing solution step.
Definition: timestep.h:80
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
#define _IFT_SpringElement_springConstant
Definition: springelement.h:45
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631

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