OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
tr_shell01.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 tr_shell01_h
36 #define tr_shell01_h
37 
38 #include "../sm/Elements/structuralelement.h"
39 #include "zznodalrecoverymodel.h"
40 #include "../sm/ErrorEstimators/zzerrorestimator.h"
41 #include "../sm/Elements/Shells/cct3d.h"
42 #include "../sm/Elements/PlaneStress/trplanrot3d.h"
43 #include "spatiallocalizer.h"
44 
45 #include <memory>
46 
47 #define _IFT_TR_SHELL01_Name "tr_shell01"
48 
49 namespace oofem {
58 {
59 protected:
61  std :: unique_ptr< CCTPlate3d > plate;
63  std :: unique_ptr< TrPlaneStrRot3d > membrane;
69  std :: unique_ptr< IntegrationRule > compositeIR;
70 
73 
74 public:
76  TR_SHELL01(int n, Domain * d);
78  virtual ~TR_SHELL01() {}
79 
80  virtual FEInterpolation *giveInterpolation() const { return plate->giveInterpolation(); }
81 
82  virtual int computeNumberOfDofs() { return 18; }
83  virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
84  { plate->giveDofManDofIDMask(inode, answer); }
85  // definition & identification
86  virtual const char *giveInputRecordName() const { return _IFT_TR_SHELL01_Name; }
87  virtual const char *giveClassName() const { return "TR_SHELL01"; }
89 
90  virtual void giveCharacteristicVector(FloatArray &answer, CharType mtrx, ValueModeType mode, TimeStep *tStep);
91  virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType mtrx, TimeStep *tStep);
92  virtual double computeVolumeAround(GaussPoint *gp);
93  virtual bool giveRotationMatrix(FloatMatrix &answer);
94 
95  virtual void updateYourself(TimeStep *tStep);
96  virtual void updateInternalState(TimeStep *tStep);
97  virtual void printOutputAt(FILE *file, TimeStep *tStep);
98  virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj = NULL);
99  virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj = NULL);
100  virtual void postInitialize();
102  void setCrossSection(int csIndx);
103 #ifdef __OOFEG
104  virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep);
105  virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type);
106  virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep);
107 #endif
108  // the membrane and plate irules are same (chacked in initializeFrom)
109  virtual int giveDefaultIntegrationRule() const { return plate->giveDefaultIntegrationRule(); }
110  virtual IntegrationRule *giveDefaultIntegrationRulePtr() { return plate->giveDefaultIntegrationRulePtr(); }
111  virtual Element_Geometry_Type giveGeometryType() const { return EGT_triangle_1; }
113  virtual MaterialMode giveMaterialMode() { return _Unknown; }
114 
116  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
117 
119  InternalStateType type, TimeStep *tStep);
120 
122  virtual void ZZErrorEstimatorI_computeLocalStress(FloatArray &answer, FloatArray &sig);
123 
124  // SpatialLocalizerI
125  virtual void SpatialLocalizerI_giveBBox(FloatArray &bb0, FloatArray &bb1);
126 
127 
128  virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords) {
129  return this->plate->computeGlobalCoordinates(answer, lcoords);
130  }
131 
132  virtual bool computeLocalCoordinates(FloatArray &answer, const FloatArray &gcoords) {
133  return this->plate->computeLocalCoordinates(answer, gcoords);
134  }
135 
136 protected:
137  virtual void computeBmatrixAt(GaussPoint *, FloatMatrix &, int = 1, int = ALL_STRAINS)
138  { OOFEM_ERROR("calling of this function is not allowed"); }
139  virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &)
140  { OOFEM_ERROR("calling of this function is not allowed"); }
141 
143 protected:
144  virtual void computeGaussPoints()
145  {
146  this->membrane->computeGaussPoints();
147  this->plate->computeGaussPoints();
148  }
149  virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
150  { OOFEM_ERROR("calling of this function is not allowed"); }
152  { OOFEM_ERROR("calling of this function is not allowed"); }
153  virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode);
154 
155 public:
156  virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
157  { OOFEM_ERROR("calling of this function is not allowed"); }
158  virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
159  { OOFEM_ERROR("calling of this function is not allowed"); }
160  virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
161  { OOFEM_ERROR("calling of this function is not allowed"); }
162 };
163 } // end namespace oofem
164 #endif
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
The element interface required by NodalAvergagingRecoveryModel.
virtual const char * giveInputRecordName() const
Definition: tr_shell01.h:86
integrationDomain
Used by integrator class to supply integration points for proper domain to be integrated (Area...
virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType mtrx, TimeStep *tStep)
Computes characteristic matrix of receiver of requested type in given time step.
Definition: tr_shell01.C:144
Class and object Domain.
Definition: domain.h:115
Element_Geometry_Type
Enumerative type used to classify element geometry Possible values are: EGT_point - point in space EG...
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: tr_shell01.C:339
virtual FEInterpolation * giveInterpolation() const
Definition: tr_shell01.h:80
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &)
Computes interpolation matrix for element unknowns.
Definition: tr_shell01.h:139
The element interface required by ZZNodalRecoveryModel.
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element&#39;s local coordinates.
Definition: tr_shell01.h:128
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: tr_shell01.C:220
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
This class implements an triangular three-node shell finite element, composed of cct3d and trplanrot3...
Definition: tr_shell01.h:57
virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes mass matrix of receiver.
Definition: tr_shell01.h:158
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual bool giveRotationMatrix(FloatMatrix &answer)
Transformation matrices updates rotation matrix between element-local and primary DOFs...
Definition: tr_shell01.C:179
#define _IFT_TR_SHELL01_Name
Definition: tr_shell01.h:47
virtual bool computeLocalCoordinates(FloatArray &answer, const FloatArray &gcoords)
Computes the element local coordinates from given global coordinates.
Definition: tr_shell01.h:132
virtual void computeGaussPoints()
Definition: tr_shell01.h:144
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
Returns equivalent nodal forces vectors.
Definition: tr_shell01.h:160
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: tr_shell01.C:252
virtual ~TR_SHELL01()
Destructor.
Definition: tr_shell01.h:78
virtual void SpatialLocalizerI_giveBBox(FloatArray &bb0, FloatArray &bb1)
Creates a bounding box of the nodes and checks if it includes the given coordinate.
Definition: tr_shell01.C:422
static IntArray loc_plate
Definition: tr_shell01.h:71
Abstract base class representing integration rule.
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
Definition: tr_shell01.h:83
virtual void updateInternalState(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: tr_shell01.C:212
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
Abstract base class for all "structural" finite elements.
std::unique_ptr< IntegrationRule > compositeIR
Element integraton rule (plate and membrane parts have their own integration rules) this one used to ...
Definition: tr_shell01.h:69
virtual int computeNumberOfDofs()
Computes or simply returns total number of element&#39;s local DOFs.
Definition: tr_shell01.h:82
virtual integrationDomain giveIntegrationDomain() const
Returns integration domain for receiver, used to initialize integration point over receiver volume...
Definition: tr_shell01.h:112
The element interface corresponding to ZZErrorEstimator.
#define OOFEM_ERROR(...)
Definition: error.h:61
virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
Computes the stress vector of receiver at given integration point, at time step tStep.
Definition: tr_shell01.h:149
static IntArray loc_membrane
Definition: tr_shell01.h:72
UnknownType
Type representing particular unknown (its physical meaning).
Definition: unknowntype.h:55
virtual void ZZErrorEstimatorI_computeLocalStress(FloatArray &answer, FloatArray &sig)
Returns stress vector in global c.s.
Definition: tr_shell01.C:365
void updateLocalNumbering(EntityRenumberingFunctor &f)
Local renumbering support.
Definition: tr_shell01.C:110
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: tr_shell01.C:64
TR_SHELL01(int n, Domain *d)
Constructor.
Definition: tr_shell01.C:56
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Definition: tr_shell01.C:529
virtual MaterialMode giveMaterialMode()
Returns material mode for receiver integration points.
Definition: tr_shell01.h:113
virtual Element_Geometry_Type giveGeometryType() const
Returns the element geometry type.
Definition: tr_shell01.h:111
#define ALL_STRAINS
virtual IntegrationRule * ZZErrorEstimatorI_giveIntegrationRule()
Returns element integration rule used to evaluate error.
Definition: tr_shell01.C:355
virtual void giveCharacteristicVector(FloatArray &answer, CharType mtrx, ValueModeType mode, TimeStep *tStep)
Computes characteristic vector of receiver of requested type in given time step.
Definition: tr_shell01.C:126
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: tr_shell01.C:461
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
CharType
Definition: chartype.h:87
virtual void computeBmatrixAt(GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Definition: tr_shell01.h:137
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes numerically stiffness matrix of receiver.
Definition: tr_shell01.h:156
Class representing the general Input Record.
Definition: inputrecord.h:101
Class Interface.
Definition: interface.h:82
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type)
Definition: tr_shell01.C:495
std::unique_ptr< CCTPlate3d > plate
Pointer to plate element.
Definition: tr_shell01.h:61
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: tr_shell01.C:230
virtual void postInitialize()
Performs post initialization steps.
Definition: tr_shell01.C:98
virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode)
Computes the load vector due to body load acting on receiver, at given time step. ...
Definition: tr_shell01.C:161
The spatial localizer element interface associated to spatial localizer.
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
Definition: tr_shell01.C:271
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
virtual const char * giveClassName() const
Definition: tr_shell01.h:87
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: tr_shell01.C:323
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
Load is base abstract class for all loads.
Definition: load.h:61
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
Definition: tr_shell01.h:151
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Access method for default integration rule.
Definition: tr_shell01.h:110
std::unique_ptr< TrPlaneStrRot3d > membrane
Pointer to membrane (plane stress) element.
Definition: tr_shell01.h:63
virtual int giveDefaultIntegrationRule() const
Returns id of default integration rule.
Definition: tr_shell01.h:109
the oofem namespace is to define a context or scope in which all oofem names are defined.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Definition: tr_shell01.C:279
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: tr_shell01.C:246
void setCrossSection(int csIndx)
Sets the cross section model of receiver.
Definition: tr_shell01.C:117

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