OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
tria1platesubsoil.C
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 #include "../sm/Elements/tria1platesubsoil.h"
36 #include "../sm/Materials/structuralms.h"
37 #include "../sm/CrossSections/structuralcrosssection.h"
38 #include "fei2dtrlin.h"
39 #include "node.h"
40 #include "material.h"
41 #include "crosssection.h"
42 #include "gausspoint.h"
43 #include "gaussintegrationrule.h"
44 #include "floatmatrix.h"
45 #include "floatarray.h"
46 #include "intarray.h"
47 #include "load.h"
48 #include "mathfem.h"
49 #include "classfactory.h"
50 
51 namespace oofem {
52 REGISTER_Element(Tria1PlateSubSoil);
53 
54 FEI2dTrLin Tria1PlateSubSoil :: interp_lin(1, 2);
55 
59 {
61  numberOfDofMans = 3;
62 }
63 
64 
67 {
68  return & interp_lin;
69 }
70 
71 
74 
75 
76 void
78 // Sets up the array containing the four Gauss points of the receiver.
79 {
80  if ( integrationRulesArray.size() == 0 ) {
81  integrationRulesArray.resize( 1 );
82  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 5) );
84  }
85 }
86 
87 
88 void
90 {
91  OOFEM_ERROR("Body load not supported, use surface load instead");
92 }
93 
94 
95 void
97 // Returns the [3x3] strain-displacement matrix {B} of the receiver,
98 // evaluated at gp.
99 {
100  FloatArray n;
101  FloatMatrix dn;
102 
105 
106  answer.resize(3, 3);
107  answer.zero();
108 
110  for ( int i = 0; i < 3; ++i ) {
111  answer(0, i) = n(i); // eps_z
112  answer(1, i) = dn(i, 0); // gamma_xz
113  answer(2, i) = dn(i, 1); // gamma_yz
114  }
115 }
116 
117 
118 void
120 {
121  this->giveStructuralCrossSection()->giveGeneralizedStress_PlateSubSoil(answer, gp, strain, tStep);
122 }
123 
124 
125 void
127 {
128  this->giveStructuralCrossSection()->give2dPlateSubSoilStiffMtrx(answer, rMode, gp, tStep);
129 }
130 
131 
134 {
135  this->numberOfGaussPoints = 1;
137 }
138 
139 
140 void
142 {
143  answer = {D_w};
144 }
145 
146 
147 void
149 {
150  FloatArray u, v;
151  u.beDifferenceOf( * this->giveNode(2)->giveCoordinates(), * this->giveNode(1)->giveCoordinates() );
152  v.beDifferenceOf( * this->giveNode(3)->giveCoordinates(), * this->giveNode(1)->giveCoordinates() );
153 
154  answer.beVectorProductOf(u, v);
155  answer.normalize();
156 }
157 
158 
159 double
161 //
162 // returns receiver's characteristic length for crack band models
163 // for a crack formed in the plane with normal normalToCrackPlane.
164 //
165 {
166  return this->giveCharacteristicLengthForPlaneElements(normalToCrackPlane);
167 }
168 
169 
170 double
172 {
173  double detJ, weight;
174 
175  weight = gp->giveWeight();
177  return detJ * weight;
178 }
179 
180 
181 void
183 // Returns the lumped mass matrix of the receiver.
184 {
185  OOFEM_ERROR("Mass matrix not provided");
186 }
187 
188 
189 int
191 {
192  return StructuralElement :: giveIPValue(answer, gp, type, tStep);
193 }
194 
195 Interface *
197 {
198  if ( interface == ZZNodalRecoveryModelInterfaceType ) {
199  return static_cast< ZZNodalRecoveryModelInterface * >(this);
200  } else if ( interface == SPRNodalRecoveryModelInterfaceType ) {
201  return static_cast< SPRNodalRecoveryModelInterface * >(this);
202  }
203 
204  return NULL;
205 }
206 
207 void
209 {
210  pap.resize(3);
211  for ( int i = 1; i < 4; i++ ) {
212  pap.at(i) = this->giveNode(i)->giveNumber();
213  }
214 }
215 
216 void
218 {
219  int found = 0;
220  answer.resize(1);
221 
222  for ( int i = 1; i < 4; i++ ) {
223  if ( pap == this->giveNode(i)->giveNumber() ) {
224  found = 1;
225  }
226  }
227 
228  if ( found ) {
229  answer.at(1) = pap;
230  } else {
231  OOFEM_ERROR("node unknown");
232  }
233 }
234 
235 
236 void
238 {
239  this->computeNmatrixAt(sgp->giveNaturalCoordinates(), answer);
240 }
241 
242 void
244 {
245  answer.resize(3);
246  answer.zero();
247  if ( iSurf == 1 ) {
248  for (int i = 1; i<=3; i++) {
249  answer.at(i) = i;
250  }
251  } else {
252  OOFEM_ERROR("wrong surface number");
253  }
254 }
255 
258 {
259  IntegrationRule *iRule = new GaussIntegrationRule(1, this, 1, 1);
260  int npoints = iRule->getRequiredNumberOfIntegrationPoints(_Triangle, approxOrder);
261  iRule->SetUpPointsOnTriangle(npoints, _Unknown);
262  return iRule;
263 }
264 
265 double
267 {
268  return this->computeVolumeAround(gp);
269 }
270 
271 
272 int
274 {
275  return 0;
276 }
277 
278 
279 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
The element interface required by ZZNodalRecoveryModel.
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
virtual double evaldNdx(FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the matrix of derivatives of interpolation functions (shape functions) at given point...
Definition: fei2dtrlin.C:53
void beVectorProductOf(const FloatArray &v1, const FloatArray &v2)
Computes vector product (or cross product) of vectors given as parameters, , and stores the result in...
Definition: floatarray.C:415
Class and object Domain.
Definition: domain.h:115
virtual FEInterpolation * giveInterpolation() const
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
Definition: fei2dtrlin.C:43
virtual void giveGeneralizedStress_PlateSubSoil(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)=0
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
The element interface required by ZZNodalRecoveryModel.
void zero()
Sets all component to zero.
Definition: intarray.C:52
static FEI2dTrLin interp_lin
Tria1PlateSubSoil(int n, Domain *d)
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
virtual IntegrationRule * GetSurfaceIntegrationRule(int iSurf)
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
virtual void computeMidPlaneNormal(FloatArray &answer, const GaussPoint *gp)
Computes mid-plane normal of receiver at integration point.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
MatResponseMode
Describes the character of characteristic material matrix.
Abstract base class representing integration rule.
void beDifferenceOf(const FloatArray &a, const FloatArray &b)
Sets receiver to be a - b.
Definition: floatarray.C:341
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
virtual void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
Abstract base class for all "structural" finite elements.
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
#define OOFEM_ERROR(...)
Definition: error.h:61
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...
REGISTER_Element(LSpace)
DofIDItem
Type representing particular dof type.
Definition: dofiditem.h:86
virtual void give2dPlateSubSoilStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
Method for computing subsoil stiffness matrix for plates.
virtual double giveWeight()
Returns integration weight of receiver.
Definition: gausspoint.h:181
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. ...
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
void resize(int n)
Checks size of receiver towards requested bounds.
Definition: intarray.C:124
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
Definition: crosssection.C:54
int numberOfGaussPoints
Number of integration points as specified by nip.
Definition: element.h:188
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
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
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 double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
Definition: fei2dtrlin.C:147
virtual double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf)
Computes volume related to integration point on local surface.
double giveCharacteristicLengthForPlaneElements(const FloatArray &normalToCrackPlane)
Returns the size of element in the given direction if the direction is in the XY plane, otherwise gives the mean size defined as the square root of the element area.
Definition: element.C:1170
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1358
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual int getRequiredNumberOfIntegrationPoints(integrationDomain dType, int approxOrder)
Abstract service.
Class Interface.
Definition: interface.h:82
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
std::vector< std::unique_ptr< IntegrationRule > > integrationRulesArray
List of integration rules of receiver (each integration rule contains associated integration points a...
Definition: element.h:170
virtual int SetUpPointsOnTriangle(int, MaterialMode mode)
Sets up receiver&#39;s integration points on triangular (area coords) integration domain.
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
virtual void giveSurfaceDofMapping(IntArray &answer, int iSurf) const
Assembles surface dof mapping mask, which provides mapping between surface local DOFs and "global" el...
Load is base abstract class for all loads.
Definition: load.h:61
the oofem namespace is to define a context or scope in which all oofem names are defined.
int giveNumber() const
Definition: femcmpnn.h:107
double normalize()
Normalizes receiver.
Definition: floatarray.C:828
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
virtual double giveCharacteristicLength(const FloatArray &normalToCrackPlane)
Returns the size of element in the given direction, in some cases adjusted (e.g.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
Class representing Gaussian-quadrature integration rule.
virtual void computeSurfaceNMatrixAt(FloatMatrix &answer, int iSurf, GaussPoint *gp)

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