OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
quad1mindlin.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 - 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 #include "../sm/Elements/Plates/quad1mindlin.h"
36 #include "../sm/Materials/structuralms.h"
37 #include "../sm/CrossSections/structuralcrosssection.h"
38 #include "node.h"
39 #include "material.h"
40 #include "crosssection.h"
41 #include "gausspoint.h"
42 #include "gaussintegrationrule.h"
43 #include "floatmatrix.h"
44 #include "floatarray.h"
45 #include "intarray.h"
46 #include "load.h"
47 #include "mathfem.h"
48 #include "fei2dquadlin.h"
49 #include "classfactory.h"
50 
51 namespace oofem {
52 REGISTER_Element(Quad1Mindlin);
53 
54 FEI2dQuadLin Quad1Mindlin :: interp_lin(1, 2);
55 
59 {
61  numberOfDofMans = 4;
62  this->reducedIntegrationFlag = false;
63 }
64 
67 {
68  return & interp_lin;
69 }
70 
73 {
74  return & interp_lin;
75 }
76 
77 void
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  // Only gravity load
92  double dV, load;
93  FloatArray force, gravity, n;
94 
95  if ( ( forLoad->giveBCGeoType() != BodyLoadBGT ) || ( forLoad->giveBCValType() != ForceLoadBVT ) ) {
96  OOFEM_ERROR("unknown load type");
97  }
98 
99  // note: force is assumed to be in global coordinate system.
100  forLoad->computeComponentArrayAt(gravity, tStep, mode);
101 
102  force.clear();
103  if ( gravity.giveSize() ) {
105  for ( GaussPoint *gp: *integrationRulesArray [ 0 ] ) {
106 
107  this->interp_lin.evalN( n, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
108  dV = this->computeVolumeAround(gp) * this->giveCrossSection()->give(CS_Thickness, gp);
109  load = this->giveStructuralCrossSection()->give('d', gp) * gravity.at(3) * dV;
110 
111  force.add(load, n);
112  }
113 
114  answer.resize(12);
115  answer.zero();
116 
117  answer.at(1) = force.at(1);
118  answer.at(4) = force.at(2);
119  answer.at(7) = force.at(3);
120  answer.at(10) = force.at(4);
121  } else {
122  answer.clear();
123  }
124 }
125 
126 
127 void
129 // Returns the [5x9] strain-displacement matrix {B} of the receiver,
130 // evaluated at gp.
131 {
132  FloatArray n, ns;
133  FloatMatrix dn, dns;
134 
137 
138  answer.resize(5, 12);
139  answer.zero();
140 
141  // enforce one-point reduced integration if requested
142  if ( this->reducedIntegrationFlag ) {
143  FloatArray lc(2);
144  lc.zero(); // set to element center coordinates
145 
146  this->interp_lin.evaldNdx( dns, lc, FEIElementGeometryWrapper(this));
147  this->interp_lin.evalN( ns, lc, FEIElementGeometryWrapper(this));
148  } else {
149  dns = dn;
150  ns = n;
151  }
152 
154  for ( int i = 0; i < 4; ++i ) {
155  answer(0, 2 + i * 3) = dn(i, 0); // kappa_x = d(fi_y)/dx
156  answer(1, 1 + i * 3) = -dn(i, 1); // kappa_y = -d(fi_x)/dy
157  answer(2, 2 + i * 3) = dn(i, 1); // kappa_xy=d(fi_y)/dy-d(fi_x)/dx
158  answer(2, 1 + i * 3) = -dn(i, 0);
159 
160  answer(3, 0 + i * 3) = dns(i, 0); // gamma_xz = fi_y+dw/dx
161  answer(3, 2 + i * 3) = ns(i);
162  answer(4, 0 + i * 3) = dns(i, 1);// gamma_yz = -fi_x+dw/dy
163  answer(4, 1 + i * 3) = -ns(i);
164  }
165 }
166 
167 
168 void
170 {
171  this->giveStructuralCrossSection()->giveGeneralizedStress_Plate(answer, gp, strain, tStep);
172 }
173 
174 
175 void
177 {
178  this->giveStructuralCrossSection()->give2dPlateStiffMtrx(answer, rMode, gp, tStep);
179 }
180 
181 
184 {
185  this->numberOfGaussPoints = 4;
188 }
189 
190 
191 void
193 {
194  answer = {D_w, R_u, R_v};
195 }
196 
197 
198 void
200 {
201  FloatArray u, v;
202  u.beDifferenceOf( * this->giveNode(2)->giveCoordinates(), * this->giveNode(1)->giveCoordinates() );
203  v.beDifferenceOf( * this->giveNode(3)->giveCoordinates(), * this->giveNode(1)->giveCoordinates() );
204 
205  answer.beVectorProductOf(u, v);
206  answer.normalize();
207 }
208 
209 
210 double
212 {
213  return this->giveCharacteristicLengthForPlaneElements(normalToCrackPlane);
214 }
215 
216 
217 double
219 {
220  double detJ, weight;
221 
222  weight = gp->giveWeight();
224  return detJ * weight;
225 }
226 
227 
228 void
230 // Returns the lumped mass matrix of the receiver.
231 {
232  double dV, mass = 0.;
233 
235  for ( GaussPoint *gp: *integrationRulesArray [ 0 ] ) {
236 
237  dV = this->computeVolumeAround(gp);
238  mass += dV * this->giveStructuralCrossSection()->give('d', gp);
239  }
240 
241  answer.resize(12, 12);
242  answer.zero();
243  answer.at(1, 1) = mass * 0.25;
244  answer.at(4, 4) = mass * 0.25;
245  answer.at(7, 7) = mass * 0.25;
246  answer.at(10, 10) = mass * 0.25;
247 }
248 
249 
250 int
252 {
253  FloatArray help;
254  answer.resize(6);
255  if ( type == IST_ShellForceTensor || type == IST_ShellStrainTensor ) {
256  if ( type == IST_ShellForceTensor ) {
257  help = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
258  } else {
259  help = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
260  }
261  answer.at(1) = 0.0; // nx
262  answer.at(2) = 0.0; // ny
263  answer.at(3) = 0.0; // nz
264  answer.at(4) = help.at(5); // vyz
265  answer.at(5) = help.at(4); // vxz
266  answer.at(6) = 0.0; // vxy
267  return 1;
268  } else if ( type == IST_ShellMomentTensor || type == IST_CurvatureTensor ) {
269  if ( type == IST_ShellMomentTensor ) {
270  help = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
271  } else {
272  help = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
273  }
274  answer.at(1) = help.at(1); // mx
275  answer.at(2) = help.at(2); // my
276  answer.at(3) = 0.0; // mz
277  answer.at(4) = 0.0; // myz
278  answer.at(5) = 0.0; // mxz
279  answer.at(6) = help.at(3); // mxy
280  return 1;
281  } else {
282  return NLStructuralElement :: giveIPValue(answer, gp, type, tStep);
283  }
284 }
285 
286 void
288 {
289  if ( iEdge == 1 ) { // edge between nodes 1 2
290  answer = {1, 2, 3, 4, 5, 6};
291  } else if ( iEdge == 2 ) { // edge between nodes 2 3
292  answer = {4, 5, 6, 7, 8, 9};
293  } else if ( iEdge == 3 ) { // edge between nodes 3 4
294  answer = {7, 8, 9, 10, 11, 12};
295  } else if ( iEdge == 4 ) { // edge between nodes 4 1
296  answer = {10, 11, 12, 1, 2, 3};
297  } else {
298  OOFEM_ERROR("wrong edge number");
299  }
300 }
301 
302 
303 double
305 {
307  return detJ *gp->giveWeight();
308 }
309 
310 
311 int
313 {
314  double dx, dy, length;
315  IntArray edgeNodes;
316  Node *nodeA, *nodeB;
317 
318  answer.resize(3, 3);
319  answer.zero();
320 
321  this->interp_lin.computeLocalEdgeMapping(edgeNodes, iEdge);
322 
323  nodeA = this->giveNode( edgeNodes.at(1) );
324  nodeB = this->giveNode( edgeNodes.at(2) );
325 
326  dx = nodeB->giveCoordinate(1) - nodeA->giveCoordinate(1);
327  dy = nodeB->giveCoordinate(2) - nodeA->giveCoordinate(2);
328  length = sqrt(dx * dx + dy * dy);
329 
330  answer.at(1, 1) = 1.0;
331  answer.at(2, 2) = dx / length;
332  answer.at(2, 3) = -dy / length;
333  answer.at(3, 2) = -answer.at(2, 3);
334  answer.at(3, 3) = answer.at(2, 2);
335 
336  return 1;
337 }
338 Interface *
340 {
341  if ( interface == ZZNodalRecoveryModelInterfaceType ) {
342  return static_cast< ZZNodalRecoveryModelInterface * >(this);
343  } else if ( interface == SPRNodalRecoveryModelInterfaceType ) {
344  return static_cast< SPRNodalRecoveryModelInterface * >(this);
345  }
346 
347  return NULL;
348 }
349 
350 
351 
352 void
354 {
355  pap.resize(4);
356  for ( int i = 1; i < 5; i++ ) {
357  pap.at(i) = this->giveNode(i)->giveNumber();
358  }
359 }
360 
361 void
363 {
364  int found = 0;
365  answer.resize(1);
366 
367  for ( int i = 1; i < 5; i++ ) {
368  if ( pap == this->giveNode(i)->giveNumber() ) {
369  found = 1;
370  }
371  }
372 
373  if ( found ) {
374  answer.at(1) = pap;
375  } else {
376  OOFEM_ERROR("node unknown");
377  }
378 }
379 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
Definition: quad1mindlin.C:192
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.
Definition: quad1mindlin.C:176
The element interface required by ZZNodalRecoveryModel.
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: quad1mindlin.C:78
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
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 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: quad1mindlin.C:312
virtual void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
Definition: quad1mindlin.C:353
virtual double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
Definition: feinterpol.C:43
virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
Definition: quad1mindlin.C:287
The element interface required by ZZNodalRecoveryModel.
Abstract base class for "structural" finite elements with geometrical nonlinearities.
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
This class implements a structural material status information.
Definition: structuralms.h:65
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
virtual bool hasField(InputFieldType id)=0
Returns true if record contains field identified by idString keyword.
virtual double edgeGiveTransformationJacobian(int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the edge Jacobian of transformation between local and global coordinates.
Definition: feinterpol2d.C:175
virtual void computeComponentArrayAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Computes boundary condition value - its components values at given time.
Definition: load.C:82
virtual double giveCoordinate(int i)
Definition: node.C:82
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.
#define _IFT_Quad1Mindlin_ReducedIntegration
Definition: quad1mindlin.h:43
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: quad1mindlin.C:89
Distributed body load.
Definition: bcgeomtype.h:43
void beDifferenceOf(const FloatArray &a, const FloatArray &b)
Sets receiver to be a - b.
Definition: floatarray.C:341
virtual void give2dPlateStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
Method for computing 2d plate stiffness matrix.
virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
Computes volume related to integration point on local edge.
Definition: quad1mindlin.C:304
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
virtual void computeLocalEdgeMapping(IntArray &edgeNodes, int iedge)
Definition: fei2dquadlin.C:295
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
DofIDItem
Type representing particular dof type.
Definition: dofiditem.h:86
virtual double giveWeight()
Returns integration weight of receiver.
Definition: gausspoint.h:181
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
virtual void giveGeneralizedStress_Plate(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)=0
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: quad1mindlin.C:251
Quad1Mindlin(int n, Domain *d)
Definition: quad1mindlin.C:56
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
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
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 double give(CrossSectionProperty a, GaussPoint *gp)
Returns the value of cross section property at given point.
Definition: crosssection.C:151
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
IntegrationPointStatus * giveMaterialStatus()
Returns reference to associated material status (NULL if not defined).
Definition: gausspoint.h:205
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Definition: quad1mindlin.C:128
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Definition: quad1mindlin.C:229
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 double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: quad1mindlin.C:218
virtual double giveCharacteristicLength(const FloatArray &normalToCrackPlane)
Returns the size of element in the given direction, in some cases adjusted (e.g.
Definition: quad1mindlin.C:211
static FEI2dQuadLin interp_lin
Definition: quad1mindlin.h:69
Class Interface.
Definition: interface.h:82
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
virtual bcGeomType giveBCGeoType() const
Returns geometry character of boundary condition.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: quad1mindlin.C:183
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 IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
Definition: quad1mindlin.C:362
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 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: quad1mindlin.C:169
Load is base abstract class for all loads.
Definition: load.h:61
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual bcValType giveBCValType() const
Returns receiver load type.
Class implementing node in finite element mesh.
Definition: node.h:87
virtual void computeMidPlaneNormal(FloatArray &answer, const GaussPoint *gp)
Computes mid-plane normal of receiver at integration point.
Definition: quad1mindlin.C:199
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 evaldNdx(FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the matrix of derivatives of interpolation functions (shape functions) at given point...
Definition: fei2dquadlin.C:75
Class representing integration point in finite element program.
Definition: gausspoint.h:93
virtual FEInterpolation * giveInterpolation() const
Definition: quad1mindlin.C:66
Class representing solution step.
Definition: timestep.h:80
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: quad1mindlin.C:339
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
void add(const FloatArray &src)
Adds array src to receiver.
Definition: floatarray.C:156
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
Definition: fei2dquadlin.C:59
bool reducedIntegrationFlag
Flag controlling reduced (one - point) integration for shear.
Definition: quad1mindlin.h:71
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631

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