OOFEM 3.0
Loading...
Searching...
No Matches
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 - 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 nodalspringelement_h
36#define nodalspringelement_h
37
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
49namespace oofem {
50class ParamKey;
63{
64 protected:
71
75public:
76 NodalSpringElement(int n, Domain * d);
77 virtual ~NodalSpringElement() { }
78
79 void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep) override;
80 void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep) override
81 { computeLumpedMassMatrix(answer, tStep); }
82 void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) override;
83 void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep) override
84 { answer.clear(); }
85
86 void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord = 0) override;
87
88 int computeNumberOfDofs() override { return dofMask.giveSize(); }
89 int computeNumberOfGlobalDofs() override;
90
91 void giveDofManDofIDMask(int inode, IntArray &answer) const override;
92
93 void updateInternalState(TimeStep *tStep) override { }
94 void updateYourself(TimeStep *tStep) override { }
95 int checkConsistency() override { return 1; }
96 void printOutputAt(FILE *file, TimeStep *tStep) override;
97
98#ifdef __OOFEG
99 //void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep);
100 //void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType);
101 //void drawScalar(oofegGraphicContext &gc, TimeStep *tStep);
102#endif
103
104 // definition & identification
105 const char *giveInputRecordName() const override { return _IFT_NodalSpringElement_Name; }
106 const char *giveClassName() const override { return "NodalSpringElement"; }
107 void initializeFrom(InputRecord &ir, int priority) override;
108 void initializeFinish() override;
109 void postInitialize() override;
110 Element_Geometry_Type giveGeometryType() const override { return EGT_point; }
111 bool isCast(TimeStep *tStep) override { return true; }
112
113protected:
114 void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
115 { answer.clear(); }
116 void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
117 { answer.clear(); }
118
120 int lowerIndx = 1, int upperIndx = ALL_STRAINS) override
121 { answer.clear(); }
122 void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override { answer.clear(); }
123 bool computeGtoLRotationMatrix(FloatMatrix &answer) override;
124};
125} // end namespace oofem
126#endif // nodalspringelement_h
*Sets size of receiver to be an empty matrix It will have zero rows and zero columns size void clear()
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
void initializeFrom(InputRecord &ir, int priority) override
void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) override
void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int lowerIndx=1, int upperIndx=ALL_STRAINS) override
void printOutputAt(FILE *file, TimeStep *tStep) override
FloatArray masses
total mass of the spring; to be distributed to nodes
const char * giveInputRecordName() const override
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
static ParamKey IPK_NodalSpringElement_springConstants
void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep) override
void updateYourself(TimeStep *tStep) override
NodalSpringElement(int n, Domain *d)
void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep) override
void updateInternalState(TimeStep *tStep) override
bool computeGtoLRotationMatrix(FloatMatrix &answer) override
static ParamKey IPK_NodalSpringElement_masses
void giveDofManDofIDMask(int inode, IntArray &answer) const override
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
void postInitialize() override
Performs post initialization steps.
static ParamKey IPK_NodalSpringElement_dofmask
int computeNumberOfGlobalDofs() override
FloatArray springConstants
Spring constants.
Element_Geometry_Type giveGeometryType() const override
const char * giveClassName() const override
bool isCast(TimeStep *tStep) override
void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0) override
void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep) override
StructuralElement(int n, Domain *d)
#define _IFT_NodalSpringElement_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