OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
cct.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 cct_h
36 #define cct_h
37 
38 #include "../sm/Elements/nlstructuralelement.h"
39 #include "../sm/CrossSections/layeredcrosssection.h"
40 #include "../sm/ErrorEstimators/zzerrorestimator.h"
41 #include "mathfem.h"
42 #include "zznodalrecoverymodel.h"
44 #include "sprnodalrecoverymodel.h"
45 
46 #define _IFT_CCTPlate_Name "cctplate"
47 
48 namespace oofem {
49 class FEI2dTrLin;
50 
62 {
63 protected:
65  //static FEI2dTrRot interp_rot;
66 
67  double area;
68 
69 public:
70  CCTPlate(int n, Domain * d);
71  virtual ~CCTPlate() { }
72 
73  virtual FEInterpolation *giveInterpolation() const;
74  virtual FEInterpolation *giveInterpolation(DofIDItem id) const;
75 
76  virtual MaterialMode giveMaterialMode() { return _2dPlate; }
77  virtual int testElementExtension(ElementExtension ext) { return ( ( ext == Element_EdgeLoadSupport ) ? 1 : 0 ); }
78  // overloaded to take into account possible element local cs (in derived cct3d)
79  virtual double computeArea();
80 
81 protected:
82  virtual void computeGaussPoints();
83  virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode);
84  virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int = 1, int = ALL_STRAINS);
85  virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer);
86 
87  virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep);
88  virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep);
89 
90 
91  virtual void giveNodeCoordinates(double &x1, double &x2, double &x3,
92  double &y1, double &y2, double &y3,
93  double &z1, double &z2, double &z3);
94 
95 
96  //virtual void computeSurfaceNMatrixAt(FloatMatrix &answer, GaussPoint *gp) { answer.clear(); }
97  virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const;
98  //virtual void giveSurfaceDofMapping(IntArray &answer, int iSurf) const { answer.clear(); }
99  //virtual IntegrationRule *GetSurfaceIntegrationRule(int i) { return NULL; }
100  virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge);
101  //virtual double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf) { return 0.; }
102  virtual int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp);
103 
104 public:
105  // definition & identification
106  virtual const char *giveClassName() const { return "CCTPlate"; }
107  virtual const char *giveInputRecordName() const { return _IFT_CCTPlate_Name; }
109 
110  virtual int computeNumberOfDofs() { return 9; }
111  virtual void giveDofManDofIDMask(int inode, IntArray &) const;
112 
113  virtual void computeMidPlaneNormal(FloatArray &answer, const GaussPoint *gp);
114 
115  virtual double giveCharacteristicLength(const FloatArray &normalToCrackPlane);
116  virtual double computeVolumeAround(GaussPoint *gp);
117 
118  virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep);
119  virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
120  { computeLumpedMassMatrix(answer, tStep); }
121 
123 
124  virtual bool computeLocalCoordinates(FloatArray &answer, const FloatArray &gcoords);
125  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
126 
128  InternalStateType type, TimeStep *tStep);
129 
131  virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap);
132  virtual int SPRNodalRecoveryMI_giveNumberOfIP() { return 1; }
134 
135  // layered cross section support functions
136  virtual void computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain,
137  GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep);
138  // boundary load support
139  virtual void computeSurfaceNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords);
140 
141 #ifdef __OOFEG
142  virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep);
143  virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type);
144  virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep);
145 #endif
146 };
147 } // end namespace oofem
148 #endif // cct_h
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
virtual void giveNodeCoordinates(double &x1, double &x2, double &x3, double &y1, double &y2, double &y3, double &z1, double &z2, double &z3)
Definition: cct.C:275
The element interface required by NodalAvergagingRecoveryModel.
virtual SPRPatchType SPRNodalRecoveryMI_givePatchType()
Definition: cct.C:528
virtual bool computeLocalCoordinates(FloatArray &answer, const FloatArray &gcoords)
Computes the element local coordinates from given global coordinates.
Definition: cct.C:411
The element interface required by ZZNodalRecoveryModel.
virtual FEInterpolation * giveInterpolation() const
Definition: cct.C:73
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
Definition: cct.C:268
Class and object Domain.
Definition: domain.h:115
The element interface required by ZZNodalRecoveryModel.
Abstract base class for "structural" finite elements with geometrical nonlinearities.
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: cct.C:389
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
virtual void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
Definition: cct.C:504
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual ~CCTPlate()
Definition: cct.h:71
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: cct.C:453
virtual int SPRNodalRecoveryMI_giveNumberOfIP()
Definition: cct.h:132
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 IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: cct.C:313
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Definition: cct.C:755
This class implements an triangular three-node plate CCT finite element.
Definition: cct.h:58
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
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...
Definition: cct.C:633
Class representing a 2d triangular linear interpolation based on area coordinates.
Definition: fei2dtrlin.h:44
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type)
Definition: cct.C:720
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
Definition: cct.C:213
The element interface corresponding to ZZErrorEstimator.
virtual double computeArea()
Computes the area (zero for all but 2d geometries).
Definition: cct.C:299
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
Definition: cct.C:321
ElementExtension
Type representing element extension.
DofIDItem
Type representing particular dof type.
Definition: dofiditem.h:86
UnknownType
Type representing particular unknown (its physical meaning).
Definition: unknowntype.h:55
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Definition: cct.C:367
virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
Computes volume related to integration point on local edge.
Definition: cct.C:625
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Definition: cct.C:145
virtual int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
Definition: cct.h:110
virtual MaterialMode giveMaterialMode()
Returns material mode for receiver integration points.
Definition: cct.h:76
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. ...
Definition: cct.C:100
#define ALL_STRAINS
virtual const char * giveClassName() const
Definition: cct.h:106
virtual void computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep)
Computes full 3D strain vector in element layer.
Definition: cct.C:538
Class representing vector of real numbers.
Definition: floatarray.h:82
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
#define _IFT_CCTPlate_Name
Definition: cct.h:46
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: cct.C:685
Class Interface.
Definition: interface.h:82
virtual const char * giveInputRecordName() const
Definition: cct.h:107
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: cct.C:88
virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
Definition: cct.C:514
CCTPlate(int n, Domain *d)
Definition: cct.C:61
double area
Definition: cct.h:67
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
Load is base abstract class for all loads.
Definition: load.h:61
The element interface required by LayeredCrossSection.
virtual void computeMidPlaneNormal(FloatArray &answer, const GaussPoint *gp)
Computes mid-plane normal of receiver at integration point.
Definition: cct.C:328
virtual int testElementExtension(ElementExtension ext)
Tests if the element implements required extension.
Definition: cct.h:77
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
Definition: cct.C:489
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: cct.C:350
virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes mass matrix of receiver.
Definition: cct.h:119
virtual double giveCharacteristicLength(const FloatArray &normalToCrackPlane)
Returns the size of element in the given direction, in some cases adjusted (e.g.
Definition: cct.C:341
static FEI2dTrLin interp_lin
Definition: cct.h:64
virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
Definition: cct.C:590
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
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: cct.C:261
Element extension for edge loads.
virtual void computeSurfaceNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords)
Computes surface interpolation matrix.
Definition: cct.C:669

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