OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
quad1mindlinshell3d.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 quad1mindlinshell3d_h
36 #define quad1mindlinshell3d_h
37 
38 #include "../sm/Elements/nlstructuralelement.h"
39 #include "zznodalrecoverymodel.h"
40 #include "sprnodalrecoverymodel.h"
41 
43 
44 #define _IFT_Quad1MindlinShell3D_Name "quad1mindlinshell3d"
45 #define _IFT_Quad1MindlinShell3D_ReducedIntegration "reducedintegration"
46 
47 
48 namespace oofem {
49 class FEI2dQuadLin;
50 
72 {
73 protected:
75  std::vector< FloatArray > lnodes;
80 
82 
87 
88 public:
89  Quad1MindlinShell3D(int n, Domain * d);
90  virtual ~Quad1MindlinShell3D();
91 
92  virtual FEInterpolation *giveInterpolation() const;
93  virtual FEInterpolation *giveInterpolation(DofIDItem id) const;
94 
95  virtual MaterialMode giveMaterialMode() { return _3dShell; }
96  virtual int testElementExtension(ElementExtension ext) { return ( ( ext == Element_EdgeLoadSupport ) ? 1 : 0 ); }
97 
98  virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep);
99  virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord = 0);
100 
101  // definition & identification
102  virtual const char *giveInputRecordName() const { return _IFT_Quad1MindlinShell3D_Name; }
103  virtual const char *giveClassName() const { return "Quad1MindlinShell3D"; }
105 
106  virtual int computeNumberOfDofs() { return 24; }
107  virtual int computeNumberOfGlobalDofs() { return 24; }
108  virtual void giveDofManDofIDMask(int inode, IntArray &) const;
109 
110  virtual void computeMidPlaneNormal(FloatArray &answer, const GaussPoint *gp);
111 
112  virtual double giveCharacteristicLength(const FloatArray &normalToCrackPlane);
113  virtual double computeVolumeAround(GaussPoint *gp);
114 
115  virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep);
116  virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
117  { computeLumpedMassMatrix(answer, tStep); }
118 
119  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
120 
121  virtual void computeLCS();
122  virtual bool computeGtoLRotationMatrix(FloatMatrix &answer);
123 
125 
127  virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap);
130 
131 
132 protected:
133  virtual void computeGaussPoints();
134  virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode);
135  virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int = 1, int = ALL_STRAINS);
136 
137  virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep);
138  virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep);
139  virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep);
140 
141  virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const;
142  virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge);
143  virtual int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp);
144  void computeVectorOfUnknowns(ValueModeType mode, TimeStep* tStep, FloatArray &shellUnknowns, FloatArray &drillUnknowns);
145  //virtual void computeSurfaceNMatrixAt(FloatMatrix &answer, GaussPoint *gp) { answer.clear(); }
146  //virtual void giveSurfaceDofMapping(IntArray &answer, int iSurf) const { answer.clear(); }
147  //virtual IntegrationRule *GetSurfaceIntegrationRule(int i) { return NULL; }
148  //virtual double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf) { return 0.; }
149 };
150 } // end namespace oofem
151 #endif // quad1mindlinshell3d_h
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. ...
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
virtual const char * giveInputRecordName() const
The element interface required by ZZNodalRecoveryModel.
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
Class and object Domain.
Definition: domain.h:115
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 FEInterpolation * giveInterpolation() const
This class implements an quadrilateral four-node shell element, using Mindlin plate theory...
virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
The element interface required by ZZNodalRecoveryModel.
Abstract base class for "structural" finite elements with geometrical nonlinearities.
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
static IntArray drillOrdering
Ordering for the drilling dofs (the out-of-plane rotations)
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual int computeNumberOfDofs()
Computes or simply returns total number of element&#39;s local DOFs.
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
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.
#define _IFT_Quad1MindlinShell3D_Name
virtual void computeMidPlaneNormal(FloatArray &answer, const GaussPoint *gp)
Computes mid-plane normal of receiver at integration point.
virtual int computeNumberOfGlobalDofs()
Computes the total number of element&#39;s global dofs.
bool reducedIntegrationFlag
Flag controlling reduced (one - point) integration for shear.
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
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.
ElementExtension
Type representing element extension.
DofIDItem
Type representing particular dof type.
Definition: dofiditem.h:86
virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes mass matrix of receiver.
virtual void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
virtual bool computeGtoLRotationMatrix(FloatMatrix &answer)
Returns transformation matrix from global c.s.
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
Evaluates nodal representation of real internal forces.
#define ALL_STRAINS
int numberOfGaussPoints
Number of integration points as specified by nip.
Definition: element.h:188
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
static IntArray shellOrdering
Ordering for the normal shell stiffness (everything but the out-of-plane rotations) ...
virtual double giveCharacteristicLength(const FloatArray &normalToCrackPlane)
Returns the size of element in the given direction, in some cases adjusted (e.g.
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual int SPRNodalRecoveryMI_giveNumberOfIP()
virtual const char * giveClassName() const
Class Interface.
Definition: interface.h:82
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
virtual MaterialMode giveMaterialMode()
Returns material mode for receiver integration points.
virtual SPRPatchType SPRNodalRecoveryMI_givePatchType()
virtual int testElementExtension(ElementExtension ext)
Tests if the element implements required extension.
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
Computes volume related to integration point on local edge.
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
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...
Load is base abstract class for all loads.
Definition: load.h:61
Class representing a 2d isoparametric linear interpolation based on natural coordinates for quadrilat...
Definition: fei2dquadlin.h:45
void computeVectorOfUnknowns(ValueModeType mode, TimeStep *tStep, FloatArray &shellUnknowns, FloatArray &drillUnknowns)
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes the stiffness matrix of receiver.
Quad1MindlinShell3D(int n, Domain *d)
Class representing integration point in finite element program.
Definition: gausspoint.h:93
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
std::vector< FloatArray > lnodes
Cached nodal coordinates in local c.s.,.
Class representing solution step.
Definition: timestep.h:80
FloatMatrix lcsMatrix
Cached coordinates in local c.s.,.
Element extension for edge loads.

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