OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
quad1platesubsoil.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 - 2014 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 quad1platesubsoil_H
36 #define quad1platesubsoil_H
37 
38 #include "../sm/Elements/structuralelement.h"
39 #include "zznodalrecoverymodel.h"
40 #include "sprnodalrecoverymodel.h"
41 
42 #define _IFT_Quad1PlateSubSoil_Name "quad1platesubsoil"
43 
44 namespace oofem {
45 class FEI2dQuadLin;
46 
60 {
61 protected:
63 
64 public:
65  Quad1PlateSubSoil(int n, Domain * d);
66  virtual ~Quad1PlateSubSoil() { }
67 
68  virtual FEInterpolation *giveInterpolation() const;
69  virtual FEInterpolation *giveInterpolation(DofIDItem id) const;
70 
71  virtual MaterialMode giveMaterialMode() { return _2dPlateSubSoil; }
72  virtual int testElementExtension(ElementExtension ext) { return ( ( ext == Element_SurfaceLoadSupport ) ? 1 : 0 ); }
73 
74  // definition & identification
75  virtual const char *giveInputRecordName() const { return _IFT_Quad1PlateSubSoil_Name; }
76  virtual const char *giveClassName() const { return "Quad1PlateSubSoil"; }
78 
79  virtual int computeNumberOfDofs() { return 4; }
80  virtual void giveDofManDofIDMask(int inode, IntArray &) const;
81 
82  virtual void computeMidPlaneNormal(FloatArray &answer, const GaussPoint *gp);
83 
84  virtual double giveCharacteristicLength(const FloatArray &normalToCrackPlane);
85  virtual double computeVolumeAround(GaussPoint *gp);
86 
87  virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep);
88  virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
89  { computeLumpedMassMatrix(answer, tStep); }
90 
91  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
93 
94 protected:
95  virtual void computeGaussPoints();
96  virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode);
97  virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int = 1, int = ALL_STRAINS);
98 
99  virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep);
100  virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep);
101 
103  virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap);
110  virtual void computeSurfaceNMatrixAt(FloatMatrix &answer, int iSurf, GaussPoint *gp);
111  virtual void giveSurfaceDofMapping(IntArray &answer, int iSurf) const;
112  virtual IntegrationRule *GetSurfaceIntegrationRule(int iSurf);
113  virtual double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf);
114  virtual int computeLoadLSToLRotationMatrix(FloatMatrix &answer, int iSurf, GaussPoint *gp);
116 };
117 } // end namespace oofem
118 #endif // quad1platesubsoil_H
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
The element interface required by ZZNodalRecoveryModel.
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.
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 computeSurfaceNMatrixAt(FloatMatrix &answer, int iSurf, GaussPoint *gp)
Class and object Domain.
Definition: domain.h:115
virtual double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf)
Computes volume related to integration point on local surface.
virtual int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
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...
The element interface required by ZZNodalRecoveryModel.
virtual IntegrationRule * GetSurfaceIntegrationRule(int iSurf)
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
virtual int testElementExtension(ElementExtension ext)
Tests if the element implements required extension.
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. ...
#define _IFT_Quad1PlateSubSoil_Name
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
virtual MaterialMode giveMaterialMode()
Returns material mode for receiver integration points.
virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes mass matrix of receiver.
This class implements an quadrilateral four-node plate subsoil element in xy plane.
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
virtual void computeMidPlaneNormal(FloatArray &answer, const GaussPoint *gp)
Computes mid-plane normal of receiver at integration point.
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 Interface * giveInterface(InterfaceType it)
Interface requesting service.
Abstract base class representing integration rule.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Element extension for surface loads.
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
virtual int SPRNodalRecoveryMI_giveNumberOfIP()
Abstract base class for all "structural" finite elements.
virtual double giveCharacteristicLength(const FloatArray &normalToCrackPlane)
Returns the size of element in the given direction, in some cases adjusted (e.g.
Quad1PlateSubSoil(int n, Domain *d)
ElementExtension
Type representing element extension.
DofIDItem
Type representing particular dof type.
Definition: dofiditem.h:86
#define ALL_STRAINS
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
int numberOfGaussPoints
Number of integration points as specified by nip.
Definition: element.h:188
virtual void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
virtual const char * giveClassName() const
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
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 FEI2dQuadLin interp_lin
Class representing the general Input Record.
Definition: inputrecord.h:101
Class Interface.
Definition: interface.h:82
virtual FEInterpolation * giveInterpolation() const
virtual SPRPatchType SPRNodalRecoveryMI_givePatchType()
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
virtual const char * giveInputRecordName() const
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
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 computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.

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