OOFEM 3.0
Loading...
Searching...
No Matches
structuralelement.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 structuralelement_h
36#define structuralelement_h
37
38#include "element.h"
39#include "matresponsemode.h"
40#include "valuemodetype.h"
41#include "integrationdomain.h"
42#include "dofmantransftype.h"
43#include "floatarray.h"
44
45#include <memory>
46
47namespace oofem {
48#define ALL_STRAINS -1
49
50class TimeStep;
51class Node;
52class Material;
53class GaussPoint;
54class FloatArray;
55class IntArray;
56class FloatMatrix;
57class SparseMtrx; // required by addNonlocalStiffnessContributions declaration
58class StructuralCrossSection;
59class IDNLMaterial;
60
96{
97protected:
99 std::unique_ptr< FloatArray >initialDisplacements;
100
101public:
107 StructuralElement(int n, Domain *d);
109 virtual ~StructuralElement();
110
111 void giveCharacteristicMatrix(FloatMatrix &answer, CharType, TimeStep *tStep) override;
112 void giveCharacteristicVector(FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep) override;
113
122 virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep);
133 virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep);
146 virtual void computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity = NULL);
155 virtual void giveMassMtrxIntegrationgMask(IntArray &answer) { answer.clear(); }
183 virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep);
187 void computeStiffnessMatrix_withIRulesAsSubcells(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep);
188
197 virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
198 {
199 OOFEM_ERROR("not implemented");
200 }
201
207
208 void computeField(ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer) override;
209
210 // stress equivalent vector in nodes (vector of internal forces)
211 // - mainly for nonLinear Analysis.
226 virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord = 0);
227
231 virtual void giveInternalForcesVector_withIRulesAsSubcells(FloatArray &answer,
232 TimeStep *tStep, int useUpdatedGpRecord = 0);
233
242 virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep);
243
244 int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override;
253 virtual void computeResultingIPTemperatureAt(FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode);
262 virtual void computeResultingIPEigenstrainAt(FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode);
263
266
277 void updateBeforeNonlocalAverage(TimeStep *tStep) override;
286 virtual void giveNonlocalLocationArray(IntArray &locationArray, const UnknownNumberingScheme &us);
291 virtual void addNonlocalStiffnessContributions(SparseMtrx &dest, const UnknownNumberingScheme &s, TimeStep *tStep);
293
294 // Overloaded methods.
295 int adaptiveUpdate(TimeStep *tStep) override;
296 void updateInternalState(TimeStep *tStep) override;
297 void updateYourself(TimeStep *tStep) override;
298 int checkConsistency() override;
299 void giveInputRecord(DynamicInputRecord &input) override;
300 const char *giveClassName() const override { return "StructuralElement"; }
301
302#ifdef __OOFEG
315 int giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode,
316 int node, TimeStep *tStep) override;
318 void showSparseMtrxStructure(CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep) override;
320 void showExtendedSparseMtrxStructure(CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep) override;
321
322#endif
323
324 // Interface for body loads applied by Sets:
325 void computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep) override;
326 void computeBoundarySurfaceLoadVector(FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global = true) override;
327 void computeBoundaryEdgeLoadVector(FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global = true) override;
329 virtual void computeEdgeNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords);
341 virtual void computeSurfaceNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords);
342
343
353 MatResponseMode rMode, GaussPoint *gp,
354 TimeStep *tStep) = 0;
355
358
359 virtual void createMaterialStatus();
360
361protected:
362
363
373 virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode);
374
378
385 virtual void computePointLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode, bool global = true);
386
394 virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const { answer.clear(); }
402 virtual void giveSurfaceDofMapping(IntArray &answer, int iSurf) const { answer.clear(); }
403
409 virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge);
415 virtual double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf);
416
417 // Global to local element c.s transformation for load vector dofs
426 answer.clear();
427 return 0;
428 }
429
430 // Local edge (LE-local Edge c.s) or surface (LS-local surface c.s) c.s
431 // to element local c.s for load vector dofs
442 virtual int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp) {
443 answer.clear();
444 return 0;
445 }
446
456 virtual int computeLoadLSToLRotationMatrix(FloatMatrix &answer, int iSurf, GaussPoint *gp) {
457 answer.clear();
458 return 0;
459 }
460
461
462public:
471 virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) = 0;
472
486 virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer,
487 int lowerIndx = 1, int upperIndx = ALL_STRAINS) = 0;
488
495 virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer);
496
497protected:
504
514 void condense(FloatMatrix *stiff, FloatMatrix *mass, FloatArray *load, IntArray *what);
515
516
521
522
523
524
525 friend class IDNLMaterial;
526 friend class TrabBoneNL3D;
527 friend class MisesMatNl;
528 friend class RankineMatNl;
529 friend class GradDpElement;
530};
531} // end namespace oofem
532#endif // structuralelement_h
Element(int n, Domain *aDomain)
Definition element.C:86
*Sets size of receiver to be an empty matrix It will have zero rows and zero columns size void clear()
virtual int computeLoadLSToLRotationMatrix(FloatMatrix &answer, int iSurf, GaussPoint *gp)
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)=0
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
void giveCharacteristicVector(FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep) override
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const
virtual void computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity=NULL)
virtual void computePointLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode, bool global=true)
virtual void setupIRForMassMtrxIntegration(IntegrationRule &iRule)
virtual void giveMassMtrxIntegrationgMask(IntArray &answer)
StructuralElement(int n, Domain *d)
void giveCharacteristicMatrix(FloatMatrix &answer, CharType, TimeStep *tStep) override
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int lowerIndx=1, int upperIndx=ALL_STRAINS)=0
virtual void computeLumpedInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode)
std::unique_ptr< FloatArray > initialDisplacements
Initial displacement vector, describes the initial nodal displacements when element has been casted.
virtual int giveNumberOfIPForMassMtrxIntegration()
virtual int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp)
virtual void giveSurfaceDofMapping(IntArray &answer, int iSurf) const
virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
void condense(FloatMatrix *stiff, FloatMatrix *mass, FloatArray *load, IntArray *what)
const char * giveClassName() const override
virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)=0
#define OOFEM_ERROR(...)
Definition error.h:79
InternalStateMode
Determines the mode of internal variable.
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
#define OOFEM_EXPORT
Definition oofemcfg.h:7
#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