OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
nodalspringelement.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 nodalspringelement_h
36 #define nodalspringelement_h
37 
38 #include "../sm/Elements/structuralelement.h"
39 
41 
42 #define _IFT_NodalSpringElement_Name "nodalspring"
43 #define _IFT_NodalSpringElement_dofmask "dofmask"
44 #define _IFT_NodalSpringElement_springConstants "k"
45 #define _IFT_NodalSpringElement_masses "m"
46 
48 
49 namespace oofem {
62 {
63  protected:
70 
71 public:
72  NodalSpringElement(int n, Domain * d);
73  virtual ~NodalSpringElement() { }
74 
75  virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep);
76  virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
77  { computeLumpedMassMatrix(answer, tStep); }
78  virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep);
79  virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
80  { answer.clear(); }
81 
82  virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord = 0);
83 
84  virtual int computeNumberOfDofs() { return dofMask.giveSize(); }
85  virtual int computeNumberOfGlobalDofs();
86 
87  virtual void giveDofManDofIDMask(int inode, IntArray &answer) const;
88 
89  virtual void updateInternalState(TimeStep *tStep) { }
90  virtual void updateYourself(TimeStep *tStep) { }
91  virtual int checkConsistency() { return 1; }
92  virtual void printOutputAt(FILE *file, TimeStep *tStep);
93 
94 #ifdef __OOFEG
95  //void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep);
96  //void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType);
97  //void drawScalar(oofegGraphicContext &gc, TimeStep *tStep);
98 #endif
99 
100  // definition & identification
101  virtual const char *giveInputRecordName() const { return _IFT_NodalSpringElement_Name; }
102  virtual const char *giveClassName() const { return "NodalSpringElement"; }
104  virtual Element_Geometry_Type giveGeometryType() const { return EGT_point; }
105  virtual bool isCast(TimeStep *tStep) {return true;}
106 
107 protected:
108  virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
109  { answer.clear(); }
111  { answer.clear(); }
112 
113  virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer,
114  int lowerIndx = 1, int upperIndx = ALL_STRAINS)
115  { answer.clear(); }
116  virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) { answer.clear(); }
117  virtual bool computeGtoLRotationMatrix(FloatMatrix &answer);
118 };
119 } // end namespace oofem
120 #endif // nodalspringelement_h
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
IntArray dofMask
Dof mask.
Class and object Domain.
Definition: domain.h:115
Element_Geometry_Type
Enumerative type used to classify element geometry Possible values are: EGT_point - point in space EG...
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
virtual bool computeGtoLRotationMatrix(FloatMatrix &answer)
Returns transformation matrix from global c.s.
FloatArray springConstants
Spring constants.
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
virtual int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes numerically stiffness matrix of receiver.
Abstract base class for all "structural" finite elements.
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
virtual void updateInternalState(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
virtual const char * giveInputRecordName() const
virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes initial stress matrix for linear stability problem.
#define ALL_STRAINS
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual int checkConsistency()
Performs consistency check.
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
NodalSpringElement(int n, Domain *d)
virtual bool isCast(TimeStep *tStep)
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
This class implements a simple linear spring element connecting the given node and the ground...
Class representing the general Input Record.
Definition: inputrecord.h:101
FloatArray masses
total mass of the spring; to be distributed to nodes
virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes mass matrix of receiver.
virtual const char * giveClassName() const
int giveSize() const
Definition: intarray.h:203
virtual Element_Geometry_Type giveGeometryType() const
Returns the element geometry type.
the oofem namespace is to define a context or scope in which all oofem names are defined.
void clear()
Sets size of receiver to be an empty matrix. It will have zero rows and zero columns size...
Definition: floatmatrix.h:516
#define _IFT_NodalSpringElement_Name
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
Returns equivalent nodal forces vectors.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
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 solution step.
Definition: timestep.h:80
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.
virtual int computeNumberOfGlobalDofs()
Computes the total number of element's global dofs.

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