OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
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 - 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 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 
47 namespace oofem {
48 #define ALL_STRAINS -1
49 
50 class TimeStep;
51 class Node;
52 class Material;
53 class GaussPoint;
54 class FloatArray;
55 class IntArray;
56 class FloatMatrix;
57 class SparseMtrx; // required by addNonlocalStiffnessContributions declaration
58 class StructuralCrossSection;
59 class IDNLMaterial;
60 
95 class StructuralElement : public Element
96 {
97 protected:
99  std :: unique_ptr< FloatArray >initialDisplacements;
100 
101 public:
107  StructuralElement(int n, Domain *d);
109  virtual ~StructuralElement();
110 
111  virtual void giveCharacteristicMatrix(FloatMatrix & answer, CharType, TimeStep * tStep);
112  virtual void giveCharacteristicVector(FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep);
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);
188 
197  virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
198  {
199  OOFEM_ERROR("not implemented");
200  }
201 
202  virtual void computeField(ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer);
203 
204  // stress equivalent vector in nodes (vector of internal forces)
205  // - mainly for nonLinear Analysis.
220  virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord = 0);
221 
226  TimeStep *tStep, int useUpdatedGpRecord = 0);
227 
236  virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep);
237 
238  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
247  virtual void computeResultingIPTemperatureAt(FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode);
256  virtual void computeResultingIPEigenstrainAt(FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode);
257 
271  virtual void updateBeforeNonlocalAverage(TimeStep *tStep);
280  virtual void giveNonlocalLocationArray(IntArray &locationArray, const UnknownNumberingScheme &us);
285  virtual void addNonlocalStiffnessContributions(SparseMtrx &dest, const UnknownNumberingScheme &s, TimeStep *tStep);
287 
288  // Overloaded methods.
289  virtual int adaptiveUpdate(TimeStep *tStep);
290  virtual void updateInternalState(TimeStep *tStep);
291  virtual void updateYourself(TimeStep *tStep);
292  virtual int checkConsistency();
294  virtual void giveInputRecord(DynamicInputRecord &input);
295  virtual const char *giveClassName() const { return "StructuralElement"; }
296 
297 #ifdef __OOFEG
298 
311  int node, TimeStep *tStep);
313  virtual void showSparseMtrxStructure(CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep);
316 
317 #endif
318 
319  // Interface for body loads applied by Sets:
320  virtual void computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep);
321  virtual void computeBoundarySurfaceLoadVector(FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global = true);
322  virtual void computeBoundaryEdgeLoadVector(FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global = true);
324  virtual void computeEdgeNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords);
336  virtual void computeSurfaceNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords);
337 
338 
347  virtual void computeConstitutiveMatrixAt(FloatMatrix &answer,
348  MatResponseMode rMode, GaussPoint *gp,
349  TimeStep *tStep) = 0;
352 
353  virtual void createMaterialStatus();
354 
355 protected:
356 
357 
367  virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode);
368 
379  virtual void computePointLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode, bool global = true);
380 
388  virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const { answer.clear(); }
396  virtual void giveSurfaceDofMapping(IntArray &answer, int iSurf) const { answer.clear(); }
397 
399  virtual IntegrationRule *GetSurfaceIntegrationRule(int order) { return NULL; }
405  virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge);
411  virtual double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf);
412 
413  // Global to local element c.s transformation for load vector dofs
422  answer.clear();
423  return 0;
424  }
425 
426  // Local edge (LE-local Edge c.s) or surface (LS-local surface c.s) c.s
427  // to element local c.s for load vector dofs
438  virtual int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp) {
439  answer.clear();
440  return 0;
441  }
452  virtual int computeLoadLSToLRotationMatrix(FloatMatrix &answer, int iSurf, GaussPoint *gp) {
453  answer.clear();
454  return 0;
455  }
457 
458 public:
467  virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) = 0;
468 
482  virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer,
483  int lowerIndx = 1, int upperIndx = ALL_STRAINS) = 0;
484 
491  virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer);
492 
493 protected:
500 
510  void condense(FloatMatrix *stiff, FloatMatrix *mass, FloatArray *load, IntArray *what);
511 
512 
516  virtual void setupIRForMassMtrxIntegration(IntegrationRule &iRule);
517 
518 
519 
520 
521  friend class IDNLMaterial;
522  friend class TrabBoneNL3D;
523  friend class MisesMatNl;
524  friend class RankineMatNl;
525  friend class GradDpElement;
526 };
527 } // end namespace oofem
528 #endif // structuralelement_h
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
Mises nonlocal material.
Definition: misesmatnl.h:91
virtual void computeResultingIPTemperatureAt(FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode)
Computes at given time (tStep) the the resulting temperature component array.
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Class and object Domain.
Definition: domain.h:115
Base class for all matrices stored in sparse format.
Definition: sparsemtrx.h:60
virtual int computeLoadLSToLRotationMatrix(FloatMatrix &answer, int iSurf, GaussPoint *gp)
Returns transformation matrix from local surface c.s to element local coordinate system of load vecto...
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
virtual void createMaterialStatus()
std::unique_ptr< FloatArray > initialDisplacements
Initial displacement vector, describes the initial nodal displacements when element has been casted...
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)=0
Computes constitutive matrix of receiver.
virtual void showExtendedSparseMtrxStructure(CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep)
Shows extended sparse structure (for example, due to nonlocal interactions for tangent stiffness) ...
virtual int giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *tStep)
Returns internal state variable (like stress,strain) at node of element in Reduced form...
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
Class implementing element body load, acting over whole element volume (e.g., the dead weight)...
Definition: bodyload.h:49
Trabecular bone nonlocal material model.
Definition: trabbonenl3d.h:83
virtual ~StructuralElement()
Destructor.
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual void giveSurfaceDofMapping(IntArray &answer, int iSurf) const
Assembles surface dof mapping mask, which provides mapping between surface local DOFs and "global" el...
virtual void computeEdgeNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords)
computes edge interpolation matrix
Abstract base class for all finite elements.
Definition: element.h:145
virtual void giveNonlocalLocationArray(IntArray &locationArray, const UnknownNumberingScheme &us)
Returns the "nonlocal" location array of receiver.
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
virtual int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp)
Returns transformation matrix from local edge c.s to element local coordinate system of load vector c...
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
virtual int adaptiveUpdate(TimeStep *tStep)
Updates the internal state variables stored in all IPs according to already mapped state...
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
Computes volume related to integration point on local edge.
virtual void setupIRForMassMtrxIntegration(IntegrationRule &iRule)
Setup Integration Rule Gauss Points for Mass Matrix integration.
virtual void showSparseMtrxStructure(CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep)
Shows sparse structure.
Abstract base class representing integration rule.
Abstract class for gradient formulation of coupled damage-plasticity model(GradDp).
Definition: graddpelement.h:48
void condense(FloatMatrix *stiff, FloatMatrix *mass, FloatArray *load, IntArray *what)
General service for condensation of stiffness and optionally load vector and mass or initial stress m...
Abstract base class for all "structural" finite elements.
virtual void computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep)
Computes the contribution of the given body load (volumetric).
virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode)
Computes the load vector due to body load acting on receiver, at given time step. ...
Abstract base class representing a boundary load (force, momentum, ...) that acts directly on a bound...
Definition: boundaryload.h:110
virtual void computeBoundarySurfaceLoadVector(FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true)
Computes the contribution of the given load at the given boundary surface in global coordinate system...
void computeStiffnessMatrix_withIRulesAsSubcells(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
#define OOFEM_ERROR(...)
Definition: error.h:61
void clear()
Clears the array (zero size).
Definition: intarray.h:177
virtual void giveInternalForcesVector_withIRulesAsSubcells(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
virtual void updateInternalState(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Abstract base class allowing to control the way, how equations are assigned to individual DOFs...
virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)=0
Computes the stress vector of receiver at given integration point, at time step tStep.
virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Compute strain vector of receiver evaluated at given integration point at given time step from elemen...
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes numerically stiffness matrix of receiver.
#define ALL_STRAINS
Class representing vector of real numbers.
Definition: floatarray.h:82
Rankine nonlocal material.
Definition: rankinematnl.h:90
virtual void giveCharacteristicVector(FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep)
Computes characteristic vector of receiver of requested type in given time step.
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
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
virtual void computeBoundaryEdgeLoadVector(FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true)
Computes the contribution of the given load at the given boundary edge.
CharType
Definition: chartype.h:87
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual void computeSurfaceNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords)
Computes surface interpolation matrix.
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
Returns equivalent nodal forces vectors.
Class representing the a dynamic Input Record.
StructuralElement(int n, Domain *d)
Constructor.
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int lowerIndx=1, int upperIndx=ALL_STRAINS)=0
Computes the geometrical matrix of receiver in given integration point.
virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes initial stress matrix for linear stability problem.
virtual void computeField(ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer)
Computes the unknown vector interpolated at the specified local coordinates.
virtual int checkConsistency()
Performs consistency check.
This class implements a Nonlocal Isotropic Damage Model for Concrete in Tension Model based on nonloc...
Definition: idmnl1.h:108
virtual void updateBeforeNonlocalAverage(TimeStep *tStep)
Updates internal element state (in all integration points of receiver) before nonlocal averaging take...
virtual double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf)
Computes volume related to integration point on local surface.
virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType, TimeStep *tStep)
Computes characteristic matrix of receiver of requested type in given time step.
virtual int giveNumberOfIPForMassMtrxIntegration()
Return desired number of integration points for consistent mass matrix computation, if required.
virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes mass matrix of receiver.
Load is base abstract class for all loads.
Definition: load.h:61
virtual const char * giveClassName() const
Abstract base class for all structural cross section models.
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual void computePointLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode, bool global=true)
Computes point load vector contribution of receiver for given load (should has BoundaryLoad Base)...
void clear()
Sets size of receiver to be an empty matrix. It will have zero rows and zero columns size...
Definition: floatmatrix.h:516
virtual void computeResultingIPEigenstrainAt(FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode)
Computes at given time the resulting eigenstrain component array.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual void addNonlocalStiffnessContributions(SparseMtrx &dest, const UnknownNumberingScheme &s, TimeStep *tStep)
Adds the "nonlocal" contribution to stiffness matrix, to account for nonlocality of material model...
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
InternalStateMode
Determines the mode of internal variable.
virtual IntegrationRule * GetSurfaceIntegrationRule(int order)
virtual void computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity=NULL)
Computes consistent mass matrix of receiver using numerical integration over element volume...
virtual void giveMassMtrxIntegrationgMask(IntArray &answer)
Returns mask indicating, which unknowns (their type and ordering is the same as element unknown vecto...

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