OOFEM 3.0
Loading...
Searching...
No Matches
linedistributedspring.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
38#include "fei3dlinelin.h"
39#include "node.h"
40#include "material.h"
41#include "crosssection.h"
42#include "gausspoint.h"
44#include "floatmatrix.h"
45#include "floatarray.h"
46#include "intarray.h"
47#include "load.h"
48#include "mathfem.h"
49#include "classfactory.h"
50#include "parametermanager.h"
51#include "paramkey.h"
52
53namespace oofem {
55
58
59FEI3dLineLin LineDistributedSpring :: interp_lin;
60
61LineDistributedSpring :: LineDistributedSpring(int n, Domain *aDomain) :
64{
67}
68
69
71LineDistributedSpring :: giveInterpolation(DofIDItem id) const
72{
73 return & interp_lin;
74}
75
76
78LineDistributedSpring :: giveInterpolation() const { return & interp_lin; }
79
80
81void
82LineDistributedSpring :: computeGaussPoints()
83// Sets up the array containing the four Gauss points of the receiver.
84{
85 if ( integrationRulesArray.size() == 0 ) {
86 integrationRulesArray.resize( 1 );
87 integrationRulesArray [ 0 ] = std::make_unique<GaussIntegrationRule>(1, this, 1, 5);
88 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 0 ], numberOfGaussPoints, this);
89 }
90}
91
92
93void
94LineDistributedSpring :: computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode)
95{
96 OOFEM_ERROR("Body load not supported, use surface load instead");
97}
98
99
100void
101LineDistributedSpring :: computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int li, int ui)
102// Returns the [3x3] strain-displacement matrix {B} of the receiver,
103// evaluated at gp.
104{
105 FloatArray n;
106 int ndofs = this->dofs.giveSize();
107
108 this->interp_lin.evalN( n, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
109
110 answer.resize(ndofs, ndofs*2);
111 answer.zero();
112
113 for (int idof=1; idof<=ndofs; idof++) {
114 answer.at(idof, idof) = n.at(1);
115 answer.at(idof, ndofs+idof) = n.at(2);
116 }
117}
118
119
120void
121LineDistributedSpring :: computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
122{
123 int ndofs = this->dofs.giveSize();
124 answer.resize(ndofs);
125
126 for (int idof=1; idof<=ndofs; idof++) {
127 answer.at(idof) = strain.at(idof)*this->springStiffnesses.at(idof);
128 }
129 //this->giveStructuralCrossSection()->giveGeneralizedStress_DistributedSpring(answer, gp, strain, tStep);
130}
131
132
133void
134LineDistributedSpring :: computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
135{
136 //this->giveStructuralCrossSection()->give2dPlateSubSoilStiffMtrx(answer, rMode, gp, tStep);
137 int ndofs = this->dofs.giveSize();
138 answer.resize(ndofs, ndofs);
139 answer.beDiagonal(this->springStiffnesses);
140}
141
142
143void
145 TimeStep *tStep, int useUpdatedGpRecord)
146{
147 FloatArray u;
148 FloatMatrix k;
149
150 this->computeVectorOf(VM_Total, tStep, u);
151 this->computeStiffnessMatrix(k, TangentStiffness, tStep);
152 answer.beProductOf(k, u);
153
154}
155
156
157
158void
159LineDistributedSpring :: initializeFrom(InputRecord &ir, int priority)
160{
161 ParameterManager &ppm = this->giveDomain()->elementPPM;
165}
166
167void
168LineDistributedSpring :: postInitialize()
169{
171 ParameterManager &ppm = this->giveDomain()->elementPPM;
174
175 if (dofs.giveSize() != springStiffnesses.giveSize()) {
176 OOFEM_ERROR ("dofs and k params size mismatch");
177 }
178}
179
180
181int
183{
184 // skip StructuralElement consistency as there is checjk for structurak cross section capability
186}
187
188void LineDistributedSpring :: printOutputAt(FILE *File, TimeStep *tStep)
189{
190 FloatArray stress, strain;
191
192 fprintf(File, "element %d (%8d) :\n", this->giveLabel(), this->giveNumber() );
193
194 for ( GaussPoint *gp: *this->giveDefaultIntegrationRulePtr() ) {
195 fprintf(File, "\tGP 1.%d :\n", gp->giveNumber());
196 this->computeStrainVector(strain, gp, tStep);
197 this->computeStressVector(stress, strain, gp, tStep);
198 fprintf(File, "\t\tStrain");
199 for (int i=1; i<=strain.giveSize(); i++) fprintf (File, " %e", strain.at(i));
200 fprintf(File, "\n\t\tStress");
201 for (int i=1; i<=stress.giveSize(); i++) fprintf (File, " %e", stress.at(i));
202 fprintf(File, "\n");
203 }
204
205}
206
207
208
209void
210LineDistributedSpring :: giveDofManDofIDMask(int inode, IntArray &answer) const
211{
212 answer = this->dofs;
213}
214
215
216double
217LineDistributedSpring :: computeVolumeAround(GaussPoint *gp)
218{
219 double detJ, weight;
220
221 weight = gp->giveWeight();
222 detJ = fabs( this->interp_lin.giveTransformationJacobian( gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) ) );
223 return detJ * weight;
224}
225
226
227void
228LineDistributedSpring :: computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
229// Returns the lumped mass matrix of the receiver.
230{
231 OOFEM_ERROR("Mass matrix not provided");
232}
233
234
235int
236LineDistributedSpring :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
237{
238 return StructuralElement :: giveIPValue(answer, gp, type, tStep);
239}
240
241Interface *
242LineDistributedSpring :: giveInterface(InterfaceType interface)
243{
244 if ( interface == ZZNodalRecoveryModelInterfaceType ) {
245 return static_cast< ZZNodalRecoveryModelInterface * >(this);
246 } else if ( interface == SPRNodalRecoveryModelInterfaceType ) {
247 return static_cast< SPRNodalRecoveryModelInterface * >(this);
248 }
249
250 return NULL;
251}
252
253void
254LineDistributedSpring :: SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
255{
256 pap.resize(2);
257 for ( int i = 1; i < 3; i++ ) {
258 pap.at(i) = this->giveNode(i)->giveNumber();
259 }
260}
261
262void
263LineDistributedSpring :: SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
264{
265 int found = 0;
266 answer.resize(1);
267
268 for ( int i = 1; i < 3; i++ ) {
269 if ( pap == this->giveNode(i)->giveNumber() ) {
270 found = 1;
271 }
272 }
273
274 if ( found ) {
275 answer.at(1) = pap;
276 } else {
277 OOFEM_ERROR("node unknown");
278 }
279}
280
281
282} // end namespace oofem
#define REGISTER_Element(class)
Node * giveNode(int i) const
Definition element.h:629
void initializeFrom(InputRecord &ir, int priority) override
Definition element.C:687
int numberOfDofMans
Number of dofmanagers.
Definition element.h:136
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Definition element.C:103
void postInitialize() override
Performs post initialization steps.
Definition element.C:797
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
int numberOfGaussPoints
Definition element.h:175
CrossSection * giveCrossSection()
Definition element.C:534
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Definition element.h:886
int checkConsistency() override
Definition element.h:807
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
Index giveSize() const
Returns the size of receiver.
Definition floatarray.h:261
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
void beDiagonal(const FloatArray &diag)
void zero()
Zeroes all coefficient of receiver.
double at(std::size_t i, std::size_t j) const
const FloatArray & giveNaturalCoordinates() const
Returns coordinate array of receiver.
Definition gausspoint.h:138
double giveWeight()
Returns integration weight of receiver.
Definition gausspoint.h:180
void resize(int n)
Definition intarray.C:73
int & at(std::size_t i)
Definition intarray.h:104
void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord) override
static ParamKey IPK_LineDistributedSpring_Stiffnesses
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
StructuralElement(int n, Domain *d)
virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
ZZNodalRecoveryModelInterface(Element *element)
Constructor.
#define OOFEM_ERROR(...)
Definition error.h:79
@ SPRNodalRecoveryModelInterfaceType
@ ZZNodalRecoveryModelInterfaceType
#define PM_ELEMENT_ERROR_IFNOTSET(_pm, _componentnum, _paramkey)
#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