OOFEM 3.0
Loading...
Searching...
No Matches
truss1d.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
37#include "fei1dlin.h"
38#include "node.h"
39#include "material.h"
40#include "gausspoint.h"
42#include "floatmatrix.h"
43#include "floatarray.h"
44#include "intarray.h"
45#include "mathfem.h"
46#include "classfactory.h"
47
48#ifdef __OOFEG
49 #include "oofeggraphiccontext.h"
50#endif
51
52namespace oofem {
54
55FEI1dLin Truss1d::interp(1); // Initiates the static interpolator
56
57
67
68
70// Sets up the array of Gauss Points of the receiver.
71{
72 if ( integrationRulesArray.size() == 0 ) {
73 integrationRulesArray.resize(1);
74 integrationRulesArray [ 0 ] = std::make_unique< GaussIntegrationRule >(1, this, 1, 2);
76 }
77}
78
81
82void
84// Returns the lumped mass matrix of the receiver. This expression is
85// valid in both local and global axes.
86{
87 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
88 double density = this->giveStructuralCrossSection()->give('d', gp);
89 double halfMass = density * this->giveCrossSection()->give(CS_Area, gp) * this->computeLength() * 0.5;
90 answer.resize(2, 2);
91 answer.zero();
92 answer.at(1, 1) = halfMass;
93 answer.at(2, 2) = halfMass;
94}
95
96
97void
99//
100// Returns linear part of geometrical equations of the receiver at gp.
101// Returns the linear part of the B matrix
102//
103{
104 FloatMatrix dN;
105 this->interp.evaldNdx(dN, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
106 answer.beTranspositionOf(dN);
107}
108
109
110void
112//
113// Returns the [1x2] displacement gradient matrix {BH} of the receiver,
114// evaluated at gp.
115// @todo not checked if correct
116{
117 this->computeBmatrixAt(gp, answer);
118}
119
120
121
122void
124// Returns the displacement interpolation matrix {N} of the receiver,
125// evaluated at gp.
126{
127 FloatArray n;
128 this->interp.evalN(n, iLocCoord, FEIElementGeometryWrapper(this) );
129 // Reshape
130 answer.resize(1, 2);
131 answer.at(1, 1) = n.at(1);
132 answer.at(1, 2) = n.at(2);
133}
134
135
136double
138// Returns the length of the receiver. This method is valid only if 1
139// Gauss point is used.
140{
141 double detJ = fabs(this->interp.giveTransformationJacobian(gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) ) );
142 return detJ * gp->giveWeight() * this->giveCrossSection()->give(CS_Area, gp);
143}
144
145
146void
148{
149 answer = this->giveStructuralCrossSection()->giveRealStress_1d(strain, gp, tStep);
150}
151
152void
153Truss1d::computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
154{
155 answer = this->giveStructuralCrossSection()->giveStiffnessMatrix_1d(rMode, gp, tStep);
156}
157
158void
160{
161 answer = this->giveStructuralCrossSection()->giveStiffnessMatrix_dPdF_1d(rMode, gp, tStep);
162}
163
164void
165Truss1d::giveDofManDofIDMask(int inode, IntArray &answer) const
166{
167 answer = {
168 D_u
169 };
170}
171
172
173void
175 IntArray &localNodeIdArray, IntArray &globalNodeIdArray,
177 int &localNodeId, int &localElemId, int &localBcId,
178 IntArray &controlNode, IntArray &controlDof,
180{
181 int nodes = 2;
182
183 FloatArray corner[ 2 ], midNode, cor[ 2 ];
186 double x = 0.0;
187 for ( int inode = 0; inode < nodes; inode++ ) {
188 corner [ inode ] = this->giveNode(inode + 1)->giveCoordinates();
189 if ( corner [ inode ].giveSize() != 3 ) {
190 cor [ inode ].resize(3);
191 cor [ inode ].at(1) = corner [ inode ].at(1);
192 cor [ inode ].at(2) = 0.0;
193 cor [ inode ].at(3) = 0.0;
194
195 corner [ inode ] = cor [ inode ];
196 }
197
198 x += corner [ inode ].at(1);
199 }
200
201 midNode.resize(3);
202
203 midNode.at(1) = x / nodes;
204 midNode.at(2) = 0.0;
205 midNode.at(3) = 0.0;
206 }
207
208 this->setupRefinedElementProblem1D(this, refinedElement, level, nodeId, localNodeIdArray, globalNodeIdArray,
209 sMode, tStep, nodes, corner, midNode,
210 localNodeId, localElemId, localBcId,
211 controlNode, controlDof, aMode, "Truss1d");
212}
213
214
219
220
221#ifdef __OOFEG
223{
224 GraphicObj *go;
225 // if (!go) { // create new one
226 WCRec p[ 2 ]; /* point */
227 if ( !gc.testElementGraphicActivity(this) ) {
228 return;
229 }
230
231 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
232 EASValsSetColor(gc.getElementColor() );
233 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
234 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
235 p [ 0 ].y = 0.;
236 p [ 0 ].z = 0.0;
237 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
238 p [ 1 ].y = 0.;
239 p [ 1 ].z = 0.0;
240 go = CreateLine3D(p);
241 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
242 EGAttachObject(go, ( EObjectP ) this);
243 EMAddGraphicsToModel(ESIModel(), go);
244}
245
246
248{
249 GraphicObj *go;
250 double defScale = gc.getDefScale();
251 // if (!go) { // create new one
252 WCRec p[ 2 ]; /* point */
253 if ( !gc.testElementGraphicActivity(this) ) {
254 return;
255 }
256
257 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
258 EASValsSetColor(gc.getDeformedElementColor() );
259 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
260 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
261 p [ 0 ].y = 0.;
262 p [ 0 ].z = 0.;
263
264 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
265 p [ 1 ].y = 0.;
266 p [ 1 ].z = 0.;
267 go = CreateLine3D(p);
268 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
269 EMAddGraphicsToModel(ESIModel(), go);
270}
271
272
274{
275 int i, indx, result = 0;
276 WCRec p[ 2 ];
277 GraphicObj *tr;
278 FloatArray v1, v2;
279 double s[ 2 ], defScale;
280
281 if ( !gc.testElementGraphicActivity(this) ) {
282 return;
283 }
284
285 if ( gc.giveIntVarMode() == ISM_recovered ) {
286 result += this->giveInternalStateAtNode(v1, gc.giveIntVarType(), gc.giveIntVarMode(), 1, tStep);
287 result += this->giveInternalStateAtNode(v2, gc.giveIntVarType(), gc.giveIntVarMode(), 2, tStep);
288 } else if ( gc.giveIntVarMode() == ISM_local ) {
289 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
290 result += giveIPValue(v1, gp, gc.giveIntVarType(), tStep);
291 v2 = v1;
292 result *= 2;
293 }
294
295 if ( result != 2 ) {
296 return;
297 }
298
299 indx = gc.giveIntVarIndx();
300
301 s [ 0 ] = v1.at(indx);
302 s [ 1 ] = v2.at(indx);
303
304 EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
305
306 if ( ( gc.getScalarAlgo() == SA_ISO_SURF ) || ( gc.getScalarAlgo() == SA_ISO_LINE ) ) {
307 for ( i = 0; i < 2; i++ ) {
308 if ( gc.getInternalVarsDefGeoFlag() ) {
309 // use deformed geometry
310 defScale = gc.getDefScale();
311 p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(1, tStep, defScale);
312 p [ i ].y = 0.;
313 p [ i ].z = 0.;
314 } else {
315 p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
316 p [ i ].y = 0.;
317 p [ i ].z = 0.;
318 }
319 }
320
321 //EASValsSetColor(gc.getYieldPlotColor(ratio));
322 tr = CreateLine3D(p);
323 EGWithMaskChangeAttributes(LAYER_MASK, tr);
324 EMAddGraphicsToModel(ESIModel(), tr);
325 } else if ( ( gc.getScalarAlgo() == SA_ZPROFILE ) || ( gc.getScalarAlgo() == SA_COLORZPROFILE ) ) {
326 double landScale = gc.getLandScale();
327
328 for ( i = 0; i < 2; i++ ) {
329 if ( gc.getInternalVarsDefGeoFlag() ) {
330 // use deformed geometry
331 defScale = gc.getDefScale();
332 p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(1, tStep, defScale);
333 p [ i ].y = 0.0;
334 p [ i ].z = s [ i ] * landScale;
335 } else {
336 p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
337 p [ i ].y = 0.0;
338 p [ i ].z = s [ i ] * landScale;
339 }
340 }
341
342 if ( gc.getScalarAlgo() == SA_ZPROFILE ) {
343 /*
344 * EASValsSetColor(gc.getDeformedElementColor());
345 * EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
346 * tr = CreateLine3D(p);
347 * EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, tr);
348 */
349 WCRec pp[ 4 ];
350 pp [ 0 ].x = p [ 0 ].x;
351 pp [ 0 ].y = 0.0;
352 pp [ 0 ].z = 0.0;
353 pp [ 1 ].x = p [ 0 ].x;
354 pp [ 1 ].y = 0.0;
355 pp [ 1 ].z = p [ 0 ].z;
356 pp [ 2 ].x = p [ 1 ].x;
357 pp [ 2 ].y = 0.0;
358 pp [ 2 ].z = p [ 1 ].z;
359 pp [ 3 ].x = p [ 1 ].x;
360 pp [ 3 ].y = 0.0;
361 pp [ 3 ].z = 0.0;
362 tr = CreateQuad3D(pp);
363 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
364 EASValsSetColor(gc.getDeformedElementColor() );
365 //EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
366 EASValsSetFillStyle(FILL_HOLLOW);
367 EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | LAYER_MASK, tr);
368 EMAddGraphicsToModel(ESIModel(), tr);
369 } else {
370 //tr = CreateTriangleWD3D(p, s[0], s[1], s[2]);
371 EASValsSetColor(gc.getDeformedElementColor() );
372 tr = CreateLine3D(p);
373 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, tr);
374 EMAddGraphicsToModel(ESIModel(), tr);
375 }
376 }
377}
378
379#endif
380
381
382Interface *
384{
385 if ( interface == ZZNodalRecoveryModelInterfaceType ) {
386 return static_cast< ZZNodalRecoveryModelInterface * >( this );
387 } else if ( interface == NodalAveragingRecoveryModelInterfaceType ) {
388 return static_cast< NodalAveragingRecoveryModelInterface * >( this );
389 } else if ( interface == SpatialLocalizerInterfaceType ) {
390 return static_cast< SpatialLocalizerInterface * >( this );
391 } else if ( interface == ZZErrorEstimatorInterfaceType ) {
392 return static_cast< ZZErrorEstimatorInterface * >( this );
393 } else if ( interface == HuertaErrorEstimatorInterfaceType ) {
394 return static_cast< HuertaErrorEstimatorInterface * >( this );
395 }
396
397 return NULL;
398}
399
400
401void
403 InternalStateType type, TimeStep *tStep)
404{
405 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
406 this->giveIPValue(answer, gp, type, tStep);
407}
408} // end namespace oofem
#define REGISTER_Element(class)
virtual double give(CrossSectionProperty a, GaussPoint *gp) const
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
double giveCoordinate(int i) const
Definition dofmanager.h:383
const FloatArray & giveCoordinates() const
Definition dofmanager.h:390
Node * giveNode(int i) const
Definition element.h:629
int numberOfDofMans
Number of dofmanagers.
Definition element.h:136
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
CrossSection * giveCrossSection()
Definition element.C:534
virtual double computeLength()
Definition element.C:1149
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
void beTranspositionOf(const FloatMatrix &src)
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
const FloatArray & giveSubPatchCoordinates() const
Returns local sub-patch coordinates of the receiver.
Definition gausspoint.h:142
double giveWeight()
Returns integration weight of receiver.
Definition gausspoint.h:180
void setupRefinedElementProblem1D(Element *element, RefinedElement *refinedElement, int level, int nodeId, IntArray &localNodeIdArray, IntArray &globalNodeIdArray, HuertaErrorEstimatorInterface ::SetupMode mode, TimeStep *tStep, int nodes, FloatArray *corner, FloatArray &midNode, int &localNodeId, int &localElemId, int &localBcId, IntArray &controlNode, IntArray &controlDof, HuertaErrorEstimator ::AnalysisMode aMode, const char *edgetype)
NLStructuralElement(int n, Domain *d)
virtual double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Definition node.C:234
SpatialLocalizerInterface(Element *element)
virtual FloatMatrixF< 1, 1 > giveStiffnessMatrix_dPdF_1d(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0
virtual FloatArrayF< 1 > giveRealStress_1d(const FloatArrayF< 1 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const =0
virtual FloatMatrixF< 1, 1 > giveStiffnessMatrix_1d(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
int giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *tStep) override
void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType) override
Definition truss1d.C:247
void computeGaussPoints() override
Definition truss1d.C:69
void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep) override
Definition truss1d.C:83
Truss1d(int n, Domain *d)
Definition truss1d.C:58
void HuertaErrorEstimatorI_setupRefinedElementProblem(RefinedElement *refinedElement, int level, int nodeId, IntArray &localNodeIdArray, IntArray &globalNodeIdArray, HuertaErrorEstimatorInterface::SetupMode sMode, TimeStep *tStep, int &localNodeId, int &localElemId, int &localBcId, IntArray &controlNode, IntArray &controlDof, HuertaErrorEstimator::AnalysisMode aMode) override
Definition truss1d.C:174
void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep) override
Definition truss1d.C:402
void computeBHmatrixAt(GaussPoint *gp, FloatMatrix &answer) override
Definition truss1d.C:111
void giveDofManDofIDMask(int inode, IntArray &) const override
Definition truss1d.C:165
void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep) override
Definition truss1d.C:222
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
Definition truss1d.C:123
void computeConstitutiveMatrix_dPdF_At(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
Definition truss1d.C:159
FEInterpolation * giveInterpolation() const override
Definition truss1d.C:80
static FEI1dLin interp
Definition truss1d.h:62
double computeVolumeAround(GaussPoint *gp) override
Definition truss1d.C:137
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
Definition truss1d.C:153
Interface * giveInterface(InterfaceType it) override
Definition truss1d.C:383
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
Definition truss1d.C:147
void HuertaErrorEstimatorI_computeNmatrixAt(GaussPoint *gp, FloatMatrix &answer) override
Definition truss1d.C:215
void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS) override
Definition truss1d.C:98
void drawScalar(oofegGraphicContext &gc, TimeStep *tStep) override
Definition truss1d.C:273
ZZErrorEstimatorInterface(Element *element)
Constructor.
ZZNodalRecoveryModelInterface(Element *element)
Constructor.
@ CS_Area
Area.
@ HuertaErrorEstimatorInterfaceType
@ ZZNodalRecoveryModelInterfaceType
@ ZZErrorEstimatorInterfaceType
@ SpatialLocalizerInterfaceType
@ NodalAveragingRecoveryModelInterfaceType
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
#define OOFEG_VARPLOT_PATTERN_LAYER
#define OOFEG_DEFORMED_GEOMETRY_LAYER
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
#define OOFEG_RAW_GEOMETRY_WIDTH
#define OOFEG_RAW_GEOMETRY_LAYER

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