OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
trplanestressrotallman3d.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/PlaneStress/trplanestressrotallman3d.h"
36 #include "../sm/CrossSections/structuralcrosssection.h"
37 #include "../sm/Materials/structuralms.h"
38 #include "material.h"
39 #include "node.h"
40 #include "load.h"
41 #include "mathfem.h"
42 #include "classfactory.h"
43 #include "gaussintegrationrule.h"
44 #include "gausspoint.h"
45 #include "fei2dtrlin.h"
46 
47 namespace oofem {
48 REGISTER_Element(TrPlanestressRotAllman3d);
49 
51 {
52  GtoLRotationMatrix = NULL;
53 }
54 
55 
56 void
58 // Returns global coordinates given in global vector
59 // transformed into local coordinate system of the
60 // receiver
61 {
62  // test if pereviously computed
63  if ( GtoLRotationMatrix == NULL ) {
65  }
66 
67 
68  lxy.resize(6);
69  const FloatArray *nc;
70  for ( int i = 0; i < 3; i++ ) {
71  nc = this->giveNode(i + 1)->giveCoordinates();
72  lxy [ i ].beProductOf(* GtoLRotationMatrix, * nc);
73  }
74  lxy [ 3 ].resize(3);
75  lxy [ 4 ].resize(3);
76  lxy [ 5 ].resize(3);
77  for ( int i = 1; i <= 3; i++ ) {
78  lxy [ 3 ].at(i) = 0.5 * ( lxy [ 0 ].at(i) + lxy [ 1 ].at(i) );
79  lxy [ 4 ].at(i) = 0.5 * ( lxy [ 1 ].at(i) + lxy [ 2 ].at(i) );
80  lxy [ 5 ].at(i) = 0.5 * ( lxy [ 2 ].at(i) + lxy [ 0 ].at(i) );
81  }
82 }
83 
84 
85 double
87 {
88  double detJ, weight;
89 
90  std :: vector< FloatArray >lc;
92 
93  weight = gp->giveWeight();
95  return detJ * weight * this->giveStructuralCrossSection()->give(CS_Thickness, gp);
96 }
97 
98 
99 void
101 {
102  answer = {
103  D_u, D_v, D_w, R_u, R_v, R_w
104  };
105 }
106 
107 
108 const FloatMatrix *
110 // Returns the rotation matrix of the receiver of the size [3,3]
111 // coords(local) = T * coords(global)
112 //
113 // local coordinate (described by vector triplet e1',e2',e3') is defined as follows:
114 //
115 // e1' : [N2-N1] Ni - means i - th node
116 // help : [N3-N1]
117 // e3' : e1' x help
118 // e2' : e3' x e1'
119 {
120  if ( GtoLRotationMatrix == NULL ) {
121  FloatArray e1(3), e2(3), e3(3), help(3);
122 
123  // compute e1' = [N2-N1] and help = [N3-N1]
124  e1.beDifferenceOf( * this->giveNode(2)->giveCoordinates(), * this->giveNode(1)->giveCoordinates() );
125  help.beDifferenceOf( * this->giveNode(3)->giveCoordinates(), * this->giveNode(1)->giveCoordinates() );
126 
127  // let us normalize e1'
128  e1.normalize();
129 
130  // compute e3' : vector product of e1' x help
131  e3.beVectorProductOf(e1, help);
132  // let us normalize
133  e3.normalize();
134 
135  // now from e3' x e1' compute e2'
136  e2.beVectorProductOf(e3, e1);
137 
138  //
139  GtoLRotationMatrix = new FloatMatrix(3, 3);
140 
141  for ( int i = 1; i <= 3; i++ ) {
142  GtoLRotationMatrix->at(1, i) = e1.at(i);
143  GtoLRotationMatrix->at(2, i) = e2.at(i);
144  GtoLRotationMatrix->at(3, i) = e3.at(i);
145  }
146  }
147 
148  return GtoLRotationMatrix;
149 }
150 
151 
152 bool
154 // Returns the rotation matrix of the receiver of the size [9,18]
155 // r(local) = T * r(global)
156 // for one node (r written transposed): {u,v,r3} = T * {u,v,w,r1,r2,r3}
157 {
158  // test if pereviously computed
159  if ( GtoLRotationMatrix == NULL ) {
161  }
162 
163  answer.resize(9, 18);
164  answer.zero();
165 
166  for ( int i = 1; i <= 3; i++ ) {
167  answer.at(1, i) = answer.at(1 + 3, i + 6) = answer.at(1 + 6, i + 12) = GtoLRotationMatrix->at(1, i);
168  answer.at(2, i) = answer.at(2 + 3, i + 6) = answer.at(2 + 6, i + 12) = GtoLRotationMatrix->at(2, i);
169  answer.at(3, i + 3) = answer.at(3 + 3, i + 3 + 6) = answer.at(3 + 6, i + 3 + 12) = GtoLRotationMatrix->at(3, i);
170  }
171 
172  return 1;
173 }
174 
175 
176 void
178 // returns characteristic tensor of the receiver at given gp and tStep
179 // strain vector = (Eps_X, Eps_y, Gamma_xy, Kappa_z)
180 {
181  FloatArray charVect;
183 
184  answer.resize(3, 3);
185  answer.zero();
186 
187  if ( ( type == LocalForceTensor ) || ( type == GlobalForceTensor ) ) {
188  //this->computeStressVector(charVect, gp, tStep);
189  charVect = ms->giveStressVector();
190 
191  answer.at(1, 1) = charVect.at(1);
192  answer.at(2, 2) = charVect.at(2);
193  answer.at(1, 2) = charVect.at(3);
194  answer.at(2, 1) = charVect.at(3);
195  } else if ( ( type == LocalMomentTensor ) || ( type == GlobalMomentTensor ) ) {} else if ( ( type == LocalStrainTensor ) || ( type == GlobalStrainTensor ) ) {
196  //this->computeStrainVector(charVect, gp, tStep);
197  charVect = ms->giveStrainVector();
198 
199  answer.at(1, 1) = charVect.at(1);
200  answer.at(2, 2) = charVect.at(2);
201  answer.at(1, 2) = charVect.at(3) / 2.;
202  answer.at(2, 1) = charVect.at(3) / 2.;
203  } else if ( ( type == LocalCurvatureTensor ) || ( type == GlobalCurvatureTensor ) ) {} else {
204  OOFEM_ERROR("unsupported tensor mode");
205  exit(1);
206  }
207 
208  if ( ( type == GlobalForceTensor ) || ( type == GlobalMomentTensor ) ||
209  ( type == GlobalStrainTensor ) || ( type == GlobalCurvatureTensor ) ) {
212  }
213 }
214 
215 
216 int
218 {
219  FloatMatrix globTensor;
220  CharTensor cht;
221 
222  if ( type == IST_ShellForceTensor || type == IST_ShellStrainTensor ) {
223  double c = 1.0;
224  if ( type == IST_ShellForceTensor ) {
225  cht = GlobalForceTensor;
226  } else {
227  cht = GlobalStrainTensor;
228  c = 2.0; // tensor components reported
229  }
230 
231  this->giveCharacteristicTensor(globTensor, cht, gp, tStep);
232 
233  answer.resize(6);
234  answer.at(1) = globTensor.at(1, 1); //xx
235  answer.at(2) = globTensor.at(2, 2); //yy
236  answer.at(3) = globTensor.at(3, 3); //zz
237  answer.at(4) = c * globTensor.at(2, 3); //yz
238  answer.at(5) = c * globTensor.at(1, 3); //xz
239  answer.at(6) = c * globTensor.at(1, 2); //xy
240  // mutiply stresses by thickness to get forces
241 
242  return 1;
243  } else if ( type == IST_ShellMomentTensor || type == IST_CurvatureTensor ) {
244  answer.clear();
245  return 1;
246  } else {
247  return TrPlanestressRotAllman :: giveIPValue(answer, gp, type, tStep);
248  }
249 }
250 
251 int
253 // Returns the rotation matrix of the receiver of the size [6,6]
254 // f(local) = T * f(global)
255 {
256  // test if previously computed
257  if ( GtoLRotationMatrix == NULL ) {
259  }
260 
261  answer.resize(6, 6);
262  answer.zero();
263 
264  for ( int i = 1; i <= 3; i++ ) {
265  answer.at(1, i) = answer.at(4, i + 3) = GtoLRotationMatrix->at(1, i);
266  answer.at(2, i) = answer.at(5, i + 3) = GtoLRotationMatrix->at(2, i);
267  answer.at(3, i) = answer.at(6, i + 3) = GtoLRotationMatrix->at(3, i);
268  }
269 
270  return 1;
271 }
272 
273 
274 void
276 // Computes numerically the load vector of the receiver due to the body loads, at tStep.
277 // load is assumed to be in global cs.
278 // load vector is then transformed to coordinate system in each node.
279 // (should be global coordinate system, but there may be defined
280 // different coordinate system in each node)
281 {
282  double dens, dV, load;
283  GaussPoint *gp = NULL;
284  FloatArray force;
285  FloatMatrix T;
286 
287  if ( ( forLoad->giveBCGeoType() != BodyLoadBGT ) || ( forLoad->giveBCValType() != ForceLoadBVT ) ) {
288  OOFEM_ERROR("unknown load type");
289  }
290 
291  GaussIntegrationRule irule(0, this);
292  irule.SetUpPointsOnTriangle( 1, this->giveMaterialMode() );
293 
294  // note: force is assumed to be in global coordinate system.
295  forLoad->computeComponentArrayAt(force, tStep, mode);
296 
297  if ( force.giveSize() ) {
298  gp = irule.getIntegrationPoint(0);
299 
300  dens = this->giveStructuralCrossSection()->give('d', gp);
301  dV = this->computeVolumeAround(gp);
302 
303  answer.resize(18);
304  answer.zero();
305 
306  load = force.at(1) * dens * dV / 3.0;
307  answer.at(1) = load;
308  answer.at(7) = load;
309  answer.at(13) = load;
310 
311  load = force.at(2) * dens * dV / 3.0;
312  answer.at(2) = load;
313  answer.at(8) = load;
314  answer.at(14) = load;
315 
316  load = force.at(3) * dens * dV / 3.0;
317  answer.at(3) = load;
318  answer.at(9) = load;
319  answer.at(15) = load;
320 
321  // transform result from global cs to local element cs.
322  if ( this->computeGtoLRotationMatrix(T) ) {
323  answer.rotatedWith(T, 'n');
324  }
325  } else {
326  answer.clear(); // nil resultant
327  }
328 }
329 
330 void
332 // Performs end-of-step operations.
333 {
334  FloatArray v;
335 
336  fprintf( file, "element %d (%8d) :\n", this->giveLabel(), this->giveNumber() );
337 
338  for ( int i = 0; i < ( int ) integrationRulesArray.size(); i++ ) {
339  for ( GaussPoint *gp : *integrationRulesArray [ i ] ) {
340  fprintf( file, " GP %2d.%-2d :", i + 1, gp->giveNumber() );
341 
342  this->giveIPValue(v, gp, IST_ShellStrainTensor, tStep);
343  fprintf(file, " strains ");
344  for ( auto &val : v ) {
345  fprintf(file, " %.4e", val);
346  }
347 
348  this->giveIPValue(v, gp, IST_CurvatureTensor, tStep);
349  fprintf(file, "\n curvatures ");
350  for ( auto &val : v ) {
351  fprintf(file, " %.4e", val);
352  }
353 
354  // Forces - Moments
355  this->giveIPValue(v, gp, IST_ShellForceTensor, tStep);
356  fprintf(file, "\n stresses ");
357  for ( auto &val : v ) {
358  fprintf(file, " %.4e", val);
359  }
360 
361  this->giveIPValue(v, gp, IST_ShellMomentTensor, tStep);
362  fprintf(file, "\n moments ");
363  for ( auto &val : v ) {
364  fprintf(file, " %.4e", val);
365  }
366 
367  fprintf(file, "\n");
368  }
369  }
370 }
371 } // end namespace oofem
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Class and object Domain.
Definition: domain.h:115
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
void giveCharacteristicTensor(FloatMatrix &answer, CharTensor type, GaussPoint *gp, TimeStep *tStep)
Wrapper around cell with vertex coordinates stored in FloatArray**.
Definition: feinterpol.h:115
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 void computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode)
Computes the load vector due to body load acting on receiver, at given time step. ...
virtual void computeComponentArrayAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Computes boundary condition value - its components values at given time.
Definition: load.C:82
virtual MaterialMode giveMaterialMode()
Returns material mode for receiver integration points.
FloatMatrix * GtoLRotationMatrix
Transformation Matrix form GtoL(3,3) is stored at the element level for computation efficiency...
Class implementing an array of integers.
Definition: intarray.h:61
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
void rotatedWith(FloatMatrix &r, char mode)
Returns the receiver a rotated according the change-of-base matrix r.
Definition: floatarray.C:799
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 int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
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.
void rotatedWith(const FloatMatrix &r, char mode= 'n')
Returns the receiver &#39;a&#39; transformed using give transformation matrix r.
Definition: floatmatrix.C:1557
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
Definition: floatarray.C:676
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
const FloatMatrix * computeGtoLRotationMatrix()
Class representing vector of real numbers.
Definition: floatarray.h:82
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
virtual double give(CrossSectionProperty a, GaussPoint *gp)
Returns the value of cross section property at given point.
Definition: crosssection.C:151
GaussPoint * getIntegrationPoint(int n)
Access particular integration point of receiver.
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
virtual double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
Definition: fei2dtrlin.C:147
const FloatArray & giveStressVector() const
Returns the const pointer to receiver&#39;s stress vector.
Definition: structuralms.h:107
IntegrationPointStatus * giveMaterialStatus()
Returns reference to associated material status (NULL if not defined).
Definition: gausspoint.h:205
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1358
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
virtual bcGeomType giveBCGeoType() const
Returns geometry character of boundary condition.
void computeLocalNodalCoordinates(std::vector< FloatArray > &lxy)
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 FloatArray * giveCoordinates()
Definition: node.h:114
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
Load is base abstract class for all loads.
Definition: load.h:61
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
int giveLabel() const
Definition: element.h:1055
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.
int giveNumber() const
Definition: femcmpnn.h:107
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
const FloatArray & giveStrainVector() const
Returns the const pointer to receiver&#39;s strain vector.
Definition: structuralms.h:105
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
Class implements an triangular three-node plane- stress elasticity finite element with independentver...
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
virtual int SetUpPointsOnTriangle(int nPoints, MaterialMode mode)
Sets up receiver&#39;s integration points on triangular (area coords) integration domain.
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631
static FEI2dTrLin interp
Definition: trplanstrss.h:67

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