OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
trplanestressrotallman.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/trplanestressrotallman.h"
36 #include "../sm/CrossSections/structuralcrosssection.h"
37 #include "fei2dtrquad.h"
38 #include "fei2dtrlin.h"
39 #include "node.h"
40 #include "gausspoint.h"
41 #include "gaussintegrationrule.h"
42 #include "floatmatrix.h"
43 #include "floatarray.h"
44 #include "intarray.h"
45 #include "mathfem.h"
46 #include "boundaryload.h"
47 #include "classfactory.h"
48 
49 #ifdef __OOFEG
50  #include "oofeggraphiccontext.h"
51  #include "oofegutils.h"
52 #endif
53 
54 namespace oofem {
55 REGISTER_Element(TrPlanestressRotAllman);
56 
58 
61 {
62  numberOfDofMans = 3;
64 }
65 
66 Interface *
68 {
69  if ( interface == LayeredCrossSectionInterfaceType ) {
70  return static_cast< LayeredCrossSectionInterface * >(this);
71  } else if ( interface == ZZNodalRecoveryModelInterfaceType ) {
72  return static_cast< ZZNodalRecoveryModelInterface * >(this);
73  } else if ( interface == SPRNodalRecoveryModelInterfaceType ) {
74  return static_cast< SPRNodalRecoveryModelInterface * >(this);
75  } else if ( interface == SpatialLocalizerInterfaceType ) {
76  return static_cast< SpatialLocalizerInterface * >(this);
77 
78  }
79  return NULL;
80 }
81 
82 void
84 {
85  lxy.resize(6);
86  for ( int i = 0; i < 3; i++ ) {
87  lxy [ i ] = * this->giveNode(i + 1)->giveCoordinates();
88  }
89  lxy [ 3 ].resize(2);
90  lxy [ 4 ].resize(2);
91  lxy [ 5 ].resize(2);
92  for ( int i = 1; i <= 2; i++ ) {
93  lxy [ 3 ].at(i) = 0.5 * ( lxy [ 0 ].at(i) + lxy [ 1 ].at(i) );
94  lxy [ 4 ].at(i) = 0.5 * ( lxy [ 1 ].at(i) + lxy [ 2 ].at(i) );
95  lxy [ 5 ].at(i) = 0.5 * ( lxy [ 2 ].at(i) + lxy [ 0 ].at(i) );
96  }
97 }
98 
99 
100 
101 void
103 // Returns the displacement interpolation matrix {N} of the receiver, eva-
104 // luated at gp.
105 {
106  FloatArray L, n;
107  std::vector< FloatArray > lxy;
108 
109  answer.resize(3, 9);
110  answer.zero();
111 
112  this->computeLocalNodalCoordinates(lxy); // get ready for tranformation into 3d
113  this->qinterpolation.evalN( n, iLocCoord, FEIVertexListGeometryWrapper(lxy) );
114  this->interp.evalN( L, iLocCoord, FEIVertexListGeometryWrapper(lxy));
115 
116  answer.at(1, 1) = answer.at(2, 2) = n.at(1) + n.at(4) / 2. + n.at(6) / 2.;
117  answer.at(1, 4) = answer.at(2, 5) = n.at(2) + n.at(4) / 2. + n.at(5) / 2.;
118  answer.at(1, 7) = answer.at(2, 8) = n.at(3) + n.at(5) / 2. + n.at(6) / 2.;
119  answer.at(1, 3) = n.at(6) * ( lxy [ 0 ].at(2) - lxy [ 2 ].at(2) ) / 8.0 - n.at(4) * ( lxy [ 1 ].at(2) - lxy [ 0 ].at(2) ) / 8.0;
120  answer.at(1, 6) = n.at(4) * ( lxy [ 1 ].at(2) - lxy [ 0 ].at(2) ) / 8.0 - n.at(5) * ( lxy [ 2 ].at(2) - lxy [ 1 ].at(2) ) / 8.0;
121  answer.at(1, 9) = n.at(5) * ( lxy [ 2 ].at(2) - lxy [ 1 ].at(2) ) / 8.0 - n.at(6) * ( lxy [ 0 ].at(2) - lxy [ 2 ].at(2) ) / 8.0;
122  answer.at(2, 3) = -n.at(6) * ( lxy [ 0 ].at(1) - lxy [ 2 ].at(1) ) / 8.0 + n.at(4) * ( lxy [ 1 ].at(1) - lxy [ 0 ].at(1) ) / 8.0;
123  answer.at(2, 6) = -n.at(4) * ( lxy [ 1 ].at(1) - lxy [ 0 ].at(1) ) / 8.0 + n.at(5) * ( lxy [ 2 ].at(1) - lxy [ 1 ].at(1) ) / 8.0;
124  answer.at(2, 9) = -n.at(5) * ( lxy [ 2 ].at(1) - lxy [ 1 ].at(1) ) / 8.0 + n.at(6) * ( lxy [ 0 ].at(1) - lxy [ 2 ].at(1) ) / 8.0;
125  // linear approx for rotations
126  answer.at(3, 3) = L.at(1);
127  answer.at(3, 6) = L.at(2);
128  answer.at(3, 9) = L.at(3);
129 }
130 
131 void
133 // Returns the [3x12] strain-displacement matrix {B} of the receiver, eva-
134 // luated at gp.
135 {
136  FloatMatrix dnx;
137  std::vector< FloatArray > lxy;
138 
139  this->computeLocalNodalCoordinates(lxy); // get ready for tranformation into 3d
141 
142  answer.resize(3, 9);
143  answer.zero();
144 
145  // epsilon_x
146  answer.at(1, 1) = dnx.at(1, 1) + 0.5 * dnx.at(4, 1) + 0.5 * dnx.at(6, 1);
147  answer.at(1, 4) = dnx.at(2, 1) + 0.5 * dnx.at(4, 1) + 0.5 * dnx.at(5, 1);
148  answer.at(1, 7) = dnx.at(3, 1) + 0.5 * dnx.at(5, 1) + 0.5 * dnx.at(6, 1);
149  answer.at(1, 3) =+dnx.at(6, 1) * ( lxy [ 0 ].at(2) - lxy [ 2 ].at(2) ) / 8.0 - dnx.at(4, 1) * ( lxy [ 1 ].at(2) - lxy [ 0 ].at(2) ) / 8.0;
150  answer.at(1, 6) =+dnx.at(4, 1) * ( lxy [ 1 ].at(2) - lxy [ 0 ].at(2) ) / 8.0 - dnx.at(5, 1) * ( lxy [ 2 ].at(2) - lxy [ 1 ].at(2) ) / 8.0;
151  answer.at(1, 9) =+dnx.at(5, 1) * ( lxy [ 2 ].at(2) - lxy [ 1 ].at(2) ) / 8.0 - dnx.at(6, 1) * ( lxy [ 0 ].at(2) - lxy [ 2 ].at(2) ) / 8.0;
152 
153  // epsilon_y
154  answer.at(2, 2) = dnx.at(1, 2) + 0.5 * dnx.at(4, 2) + 0.5 * dnx.at(6, 2);
155  answer.at(2, 5) = dnx.at(2, 2) + 0.5 * dnx.at(4, 2) + 0.5 * dnx.at(5, 2);
156  answer.at(2, 8) = dnx.at(3, 2) + 0.5 * dnx.at(5, 2) + 0.5 * dnx.at(6, 2);
157  answer.at(2, 3) =-dnx.at(6, 2) * ( lxy [ 0 ].at(1) - lxy [ 2 ].at(1) ) / 8.0 + dnx.at(4, 2) * ( lxy [ 1 ].at(1) - lxy [ 0 ].at(1) ) / 8.0;
158  answer.at(2, 6) =-dnx.at(4, 2) * ( lxy [ 1 ].at(1) - lxy [ 0 ].at(1) ) / 8.0 + dnx.at(5, 2) * ( lxy [ 2 ].at(1) - lxy [ 1 ].at(1) ) / 8.0;
159  answer.at(2, 9) =-dnx.at(5, 2) * ( lxy [ 2 ].at(1) - lxy [ 1 ].at(1) ) / 8.0 + dnx.at(6, 2) * ( lxy [ 0 ].at(1) - lxy [ 2 ].at(1) ) / 8.0;
160 
161  // gamma_xy (shear)
162  answer.at(3, 1) = dnx.at(1, 2) + 0.5 * dnx.at(4, 2) + 0.5 * dnx.at(6, 2);
163  answer.at(3, 2) = dnx.at(1, 1) + 0.5 * dnx.at(4, 1) + 0.5 * dnx.at(6, 1);
164  answer.at(3, 4) = dnx.at(2, 2) + 0.5 * dnx.at(4, 2) + 0.5 * dnx.at(5, 2);
165  answer.at(3, 5) = dnx.at(2, 1) + 0.5 * dnx.at(4, 1) + 0.5 * dnx.at(5, 1);
166  answer.at(3, 7) = dnx.at(3, 2) + 0.5 * dnx.at(5, 2) + 0.5 * dnx.at(6, 2);
167  answer.at(3, 8) = dnx.at(3, 1) + 0.5 * dnx.at(5, 1) + 0.5 * dnx.at(6, 1);
168 
169  answer.at(3, 3) = dnx.at(6, 2) * ( lxy [ 0 ].at(2) - lxy [ 2 ].at(2) ) / 8.0 - dnx.at(4, 2) * ( lxy [ 1 ].at(2) - lxy [ 0 ].at(2) ) / 8.0;
170  answer.at(3, 3) += -dnx.at(6, 1) * ( lxy [ 0 ].at(1) - lxy [ 2 ].at(1) ) / 8.0 + dnx.at(4, 1) * ( lxy [ 1 ].at(1) - lxy [ 0 ].at(1) ) / 8.0;
171  answer.at(3, 6) = dnx.at(4, 2) * ( lxy [ 1 ].at(2) - lxy [ 0 ].at(2) ) / 8.0 - dnx.at(5, 2) * ( lxy [ 2 ].at(2) - lxy [ 1 ].at(2) ) / 8.0;
172  answer.at(3, 6) += -dnx.at(4, 1) * ( lxy [ 1 ].at(1) - lxy [ 0 ].at(1) ) / 8.0 + dnx.at(5, 1) * ( lxy [ 2 ].at(1) - lxy [ 1 ].at(1) ) / 8.0;
173  answer.at(3, 9) = dnx.at(5, 2) * ( lxy [ 2 ].at(2) - lxy [ 1 ].at(2) ) / 8.0 - dnx.at(6, 2) * ( lxy [ 0 ].at(2) - lxy [ 2 ].at(2) ) / 8.0;
174  answer.at(3, 9) += -dnx.at(5, 1) * ( lxy [ 2 ].at(1) - lxy [ 1 ].at(1) ) / 8.0 + dnx.at(6, 1) * ( lxy [ 0 ].at(1) - lxy [ 2 ].at(1) ) / 8.0;
175 }
176 
177 
178 void
180 {
181  // compute standard stiffness matrix
182  TrPlaneStress2d :: computeStiffnessMatrix(answer, rMode, tStep);
183  // add zero energy mode stabilization
184  FloatMatrix ks;
185  this->computeStiffnessMatrixZeroEnergyStabilization(ks, rMode, tStep);
186  answer.add(ks);
187 }
188 
189 void
191 {
192  FloatMatrix b(1, 9);
193  FloatMatrix dnx;
194  FloatArray lec = {0.333333333333, 0.333333333333, 0.333333333333}; // element center in local coordinates
195  std::vector< FloatArray > lxy;
196 
197  this->computeLocalNodalCoordinates(lxy); // get ready for tranformation into 3d
198  this->qinterpolation.evaldNdx( dnx, lec, FEIVertexListGeometryWrapper(lxy) );
199 
200  // evaluate (dv/dx-du/dy)/2. at element center
201  b.at(1, 1) = -1.0 * ( dnx.at(1, 2) + 0.5 * dnx.at(4, 2) + 0.5 * dnx.at(6, 2) );
202  b.at(1, 2) = dnx.at(1, 1) + 0.5 * dnx.at(4, 1) + 0.5 * dnx.at(6, 1);
203  b.at(1, 4) = -1.0 * ( dnx.at(2, 2) + 0.5 * dnx.at(4, 2) + 0.5 * dnx.at(5, 2) );
204  b.at(1, 5) = dnx.at(2, 1) + 0.5 * dnx.at(4, 1) + 0.5 * dnx.at(5, 1);
205  b.at(1, 7) = -1.0 * ( dnx.at(3, 2) + 0.5 * dnx.at(5, 2) + 0.5 * dnx.at(6, 2) );
206  b.at(1, 8) = dnx.at(3, 1) + 0.5 * dnx.at(5, 1) + 0.5 * dnx.at(6, 1);
207 
208  b.at(1, 3) = -dnx.at(6, 2) * ( lxy [ 0 ].at(2) - lxy [ 2 ].at(2) ) / 8.0 + dnx.at(4, 2) * ( lxy [ 1 ].at(2) - lxy [ 0 ].at(2) ) / 8.0;
209  b.at(1, 3) += -dnx.at(6, 1) * ( lxy [ 0 ].at(1) - lxy [ 2 ].at(1) ) / 8.0 + dnx.at(4, 1) * ( lxy [ 1 ].at(1) - lxy [ 0 ].at(1) ) / 8.0;
210  b.at(1, 6) = -dnx.at(4, 2) * ( lxy [ 1 ].at(2) - lxy [ 0 ].at(2) ) / 8.0 + dnx.at(5, 2) * ( lxy [ 2 ].at(2) - lxy [ 1 ].at(2) ) / 8.0;
211  b.at(1, 6) += -dnx.at(4, 1) * ( lxy [ 1 ].at(1) - lxy [ 0 ].at(1) ) / 8.0 + dnx.at(5, 1) * ( lxy [ 2 ].at(1) - lxy [ 1 ].at(1) ) / 8.0;
212  b.at(1, 9) = -dnx.at(5, 2) * ( lxy [ 2 ].at(2) - lxy [ 1 ].at(2) ) / 8.0 + dnx.at(6, 2) * ( lxy [ 0 ].at(2) - lxy [ 2 ].at(2) ) / 8.0;
213  b.at(1, 9) += -dnx.at(5, 1) * ( lxy [ 2 ].at(1) - lxy [ 1 ].at(1) ) / 8.0 + dnx.at(6, 1) * ( lxy [ 0 ].at(1) - lxy [ 2 ].at(1) ) / 8.0;
214  b.times(0.5);
215  // add -1.0*sum(r_w)/3.0
216  b.at(1, 3) -= 1.0 / 3.0;
217  b.at(1, 6) -= 1.0 / 3.0;
218  b.at(1, 9) -= 1.0 / 3.0;
219  // add alpha*Volume*B^T[G]B to element stiffness matrix
220  double G = this->giveStructuralCrossSection()->give( Gxy, this->giveDefaultIntegrationRulePtr()->getIntegrationPoint(0) );
221  double coeff = G * this->giveArea() * this->giveCrossSection()->give( CS_Thickness, this->giveDefaultIntegrationRulePtr()->getIntegrationPoint(0) ) * 1.e-6;
222  answer.beTProductOf(b, b);
223  answer.times(coeff);
224 }
225 
226 
227 
228 void
230 {
231  answer = {D_u, D_v, R_w};
232 }
233 
234 
235 
236 double
238 // returns the area occupied by the receiver
239 {
240  if ( area > 0 ) {
241  return area; // check if previously computed
242  }
243 
244  return ( area = fabs( this->interp.giveArea( FEIElementGeometryWrapper(this) ) ) );
245 }
246 
248 // Sets up the array containing the four Gauss points of the receiver.
249 {
250  if ( integrationRulesArray.size() == 0 ) {
251  integrationRulesArray.resize( 1 );
252  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 3) );
254  }
255 }
256 
257 void
259 {
260  std::vector< FloatArray > lxy;
261  FloatArray l, n;
262  IntArray en;
263  FEI2dTrQuad qi(1, 2);
264 
265  this->computeLocalNodalCoordinates(lxy); // get ready for tranformation into 3d
267  qi.computeLocalEdgeMapping(en, iedge); // get edge mapping
268  this->interp.edgeEvalN( l, iedge, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
269  answer.resize(3, 6);
270 
271  answer.at(1, 1) = answer.at(2, 2) = n.at(1) + n.at(3) / 2.0;
272  answer.at(1, 4) = answer.at(2, 5) = n.at(2) + n.at(3) / 2.0;
273  answer.at(1, 3) = n.at(3) * ( lxy [ en.at(2) - 1 ].at(2) - lxy [ en.at(1) - 1 ].at(2) ) / 8.0;
274  answer.at(1, 6) = -answer.at(1, 3);
275  answer.at(2, 3) = n.at(3) * ( lxy [ en.at(2) - 1 ].at(1) - lxy [ en.at(1) - 1 ].at(1) ) / 8.0;
276  answer.at(2, 6) = -answer.at(2, 3);
277  answer.at(3, 3) = l.at(1);
278  answer.at(3, 6) = l.at(2);
279 }
280 
281 void
283 {
284  /*
285  * provides dof mapping of local edge dofs (only nonzero are taken into account)
286  * to global element dofs
287  */
288 
289  answer.resize(6);
290  if ( iEdge == 1 ) { // edge between nodes 1,2
291  answer.at(1) = 1;
292  answer.at(2) = 2;
293  answer.at(3) = 3;
294  answer.at(4) = 4;
295  answer.at(5) = 5;
296  answer.at(6) = 6;
297  } else if ( iEdge == 2 ) { // edge between nodes 2 3
298  answer.at(1) = 4;
299  answer.at(2) = 5;
300  answer.at(3) = 6;
301  answer.at(4) = 7;
302  answer.at(5) = 8;
303  answer.at(6) = 9;
304  } else if ( iEdge == 3 ) { // edge between nodes 3 1
305  answer.at(1) = 7;
306  answer.at(2) = 8;
307  answer.at(3) = 9;
308  answer.at(4) = 1;
309  answer.at(5) = 2;
310  answer.at(6) = 3;
311  } else {
312  OOFEM_ERROR("wrong edge number");
313  }
314 }
315 
316 //
317 // layered cross section support functions
318 //
319 void
321 // returns full 3d strain vector of given layer (whose z-coordinate from center-line is
322 // stored in slaveGp) for given tStep
323 {
324  double layerZeta, layerZCoord, top, bottom;
325 
326  top = this->giveCrossSection()->give(CS_TopZCoord, masterGp);
327  bottom = this->giveCrossSection()->give(CS_BottomZCoord, masterGp);
328  layerZeta = slaveGp->giveNaturalCoordinate(3);
329  layerZCoord = 0.5 * ( ( 1. - layerZeta ) * bottom + ( 1. + layerZeta ) * top );
330  answer.resize(3); // {eps_xx,eps_yy,gamma_yz}
331 
332  answer.at(1) = masterGpStrain.at(1) * layerZCoord;
333  answer.at(2) = masterGpStrain.at(2) * layerZCoord;
334  answer.at(3) = masterGpStrain.at(3);
335 }
336 
337 void TrPlanestressRotAllman :: computeBoundaryEdgeLoadVector(FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global)
338 {
339  answer.clear();
340  if ( type != ExternalForcesVector ) {
341  return;
342  }
343 
344  double dV;
345  FloatMatrix T;
346  FloatArray globalIPcoords;
347 
348  EdgeLoad *edgeLoad = dynamic_cast< EdgeLoad * >(load);
349  if ( edgeLoad ) {
350  int approxOrder = edgeLoad->giveApproxOrder() + this->giveInterpolation()->giveInterpolationOrder();
351  int numberOfGaussPoints = ( int ) ceil( ( approxOrder + 1. ) / 2. );
352  GaussIntegrationRule iRule(1, this, 1, 1);
353  iRule.SetUpPointsOnLine(numberOfGaussPoints, _Unknown);
354  FloatArray reducedAnswer, force, ntf;
355  IntArray mask;
356  FloatMatrix n;
357 
358  for ( GaussPoint *gp: iRule ) {
359  this->computeEgdeNMatrixAt(n, boundary, gp);
360  dV = this->computeEdgeVolumeAround(gp, boundary);
361 
362  if ( edgeLoad->giveFormulationType() == Load :: FT_Entity ) {
363  edgeLoad->computeValueAt(force, tStep, gp->giveNaturalCoordinates(), mode);
364  } else {
365  //this->computeEdgeIpGlobalCoords(globalIPcoords, gp, boundary);
366  this->giveInterpolation()->boundaryEdgeLocal2Global( globalIPcoords, boundary, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
367  edgeLoad->computeValueAt(force, tStep, globalIPcoords, mode);
368  }
369 
370  // transform force
371  if ( edgeLoad->giveCoordSystMode() == Load :: CST_Global ) {
372  } else {
373  // transform from local boundary to element local c.s
374  if ( this->computeLoadLEToLRotationMatrix(T, boundary, gp) ) {
375  force.rotatedWith(T, 'n');
376  }
377  // then to global c.s
378  if ( this->computeLoadGToLRotationMtrx(T) ) {
379  force.rotatedWith(T, 't');
380  }
381  }
382 
383  ntf.beTProductOf(n, force);
384  answer.add(dV, ntf);
385  }
386 
387 
388  return;
389  } else {
390  OOFEM_ERROR("incompatible load");
391  return;
392  }
393 
394 
395 }
396 
397 
398 
399 /*
400  * double
401  * TrPlanestressRotAllman :: computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
402  * {
403  * // edge with linear geometry -> one can use linear interpolation safely
404  * double detJ = this->interp.edgeGiveTransformationJacobian( iEdge, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
405  * return detJ *gp->giveWeight();
406  * }
407  *
408  *
409  * void
410  * TrPlanestressRotAllman :: computeEdgeIpGlobalCoords(FloatArray &answer, GaussPoint *gp, int iEdge)
411  * {
412  * // edge with linear geometry -> one can use linear interpolation safely
413  * this->interp.edgeLocal2global( answer, iEdge, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
414  * }
415  *
416  *
417  * int
418  * TrPlanestressRotAllman :: computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp)
419  * {
420  * // returns transformation matrix from
421  * // edge local coordinate system
422  * // to element local c.s
423  * // (same as global c.s in this case)
424  * //
425  * // i.e. f(element local) = T * f(edge local)
426  * //
427  * double dx, dy, length;
428  * Node *nodeA, *nodeB;
429  * int aNode = 0, bNode = 0;
430  *
431  * answer.resize(2, 2);
432  * answer.zero();
433  *
434  * if ( iEdge == 1 ) { // edge between nodes 1 2
435  * aNode = 1;
436  * bNode = 2;
437  * } else if ( iEdge == 2 ) { // edge between nodes 2 3
438  * aNode = 2;
439  * bNode = 3;
440  * } else if ( iEdge == 3 ) { // edge between nodes 2 3
441  * aNode = 3;
442  * bNode = 1;
443  * } else {
444  * OOFEM_ERROR("wrong egde number");
445  * }
446  *
447  * nodeA = this->giveNode(aNode);
448  * nodeB = this->giveNode(bNode);
449  *
450  * dx = nodeB->giveCoordinate(1) - nodeA->giveCoordinate(1);
451  * dy = nodeB->giveCoordinate(2) - nodeA->giveCoordinate(2);
452  * length = sqrt(dx * dx + dy * dy);
453  *
454  * answer.at(1, 1) = dx / length;
455  * answer.at(1, 2) = -dy / length;
456  * answer.at(2, 1) = dy / length;
457  * answer.at(2, 2) = dx / length;
458  *
459  * return 1;
460  * }
461  *
462  *
463  *
464  *
465  * double TrPlanestressRotAllman :: computeVolumeAround(GaussPoint *gp)
466  * // Returns the portion of the receiver which is attached to gp.
467  * {
468  * double detJ, weight;
469  *
470  * weight = gp->giveWeight();
471  * // safe to use linear interpolation here (geometry is linear)
472  * detJ = fabs( this->interp.giveTransformationJacobian( gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) ) );
473  *
474  * return detJ *weight *this->giveCrossSection()->give(CS_Thickness), gp;
475  * }
476  */
477 
480 {
482  if ( result != IRRT_OK ) {
483  return result;
484  }
486  return IRRT_OK;
487 }
488 
489 
490 
491 int
493 {
494  return this->numberOfGaussPoints;
495 }
496 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
virtual void edgeEvalN(FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of edge interpolation functions (shape functions) at given point.
Definition: fei2dtrlin.C:184
The element interface required by ZZNodalRecoveryModel.
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Class and object Domain.
Definition: domain.h:115
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 IntegrationRule * giveDefaultIntegrationRulePtr()
Access method for default integration rule.
Definition: element.h:822
Bottom z coordinate.
Definition: crosssection.h:72
virtual void computeLocalEdgeMapping(IntArray &edgeNodes, int iedge)
Definition: fei2dtrquad.C:234
#define Gxy
Definition: matconst.h:72
Wrapper around cell with vertex coordinates stored in FloatArray**.
Definition: feinterpol.h:115
virtual FEInterpolation * giveInterpolation() const
Definition: trplanstrss.C:69
The element interface required by ZZNodalRecoveryModel.
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...
Load is specified in global c.s.
Definition: load.h:70
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
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
int giveInterpolationOrder()
Returns the interpolation order.
Definition: feinterpol.h:159
virtual CoordSystType giveCoordSystMode()
Returns receiver&#39;s coordinate system.
Definition: boundaryload.h:151
Top z coordinate.
Definition: crosssection.h:71
Interface * giveInterface(InterfaceType interface)
Interface requesting service.
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.
void computeStiffnessMatrixZeroEnergyStabilization(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes the stiffness matrix stabilization of zero energy mode (equal rotations) ...
void giveEdgeDofMapping(IntArray &answer, int iEdge) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
void rotatedWith(FloatMatrix &r, char mode)
Returns the receiver a rotated according the change-of-base matrix r.
Definition: floatarray.C:799
virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
Computes volume related to integration point on local edge.
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes the stiffness matrix of receiver.
double giveNaturalCoordinate(int i) const
Returns i-th natural element coordinate of receiver.
Definition: gausspoint.h:136
virtual double giveArea(const FEICellGeometry &cellgeo) const
Computes the exact area.
Definition: fei2dtrlin.C:267
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
Definition: fei2dtrquad.C:43
virtual void computeLocalNodalCoordinates(std::vector< FloatArray > &lxy)
Abstract base class representing a boundary load (force, momentum, ...) that acts directly on a bound...
Definition: boundaryload.h:110
void computeEgdeNMatrixAt(FloatMatrix &answer, int iedge, GaussPoint *gp)
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
void computeBoundaryEdgeLoadVector(FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global)
Computes the contribution of the given load at the given boundary edge.
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
virtual FormulationType giveFormulationType()
Specifies is load should take local or global coordinates.
Definition: load.h:155
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
void times(double f)
Multiplies receiver by factor f.
Definition: floatmatrix.C:1594
virtual void computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep)
Computes full 3D strain vector in element layer.
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
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: fei2dtrquad.C:60
Abstract base class representing an edge load (force, momentum, ...) that acts directly on a edge bou...
Definition: boundaryload.h:200
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
virtual int giveApproxOrder()=0
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix^T and anArray.
Definition: floatarray.C:708
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 void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
Definition: crosssection.C:54
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
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
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes the stiffness matrix of receiver.
CharType
Definition: chartype.h:87
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
void add(const FloatMatrix &a)
Adds matrix to the receiver.
Definition: floatmatrix.C:1023
This class implements an triangular three-node plane-stress elasticity finite element.
Definition: trplanstrss.h:60
Class Interface.
Definition: interface.h:82
void beTProductOf(const FloatMatrix &a, const FloatMatrix &b)
Assigns to the receiver product of .
Definition: floatmatrix.C:367
The spatial localizer element interface associated to spatial localizer.
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
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
The element interface required by LayeredCrossSection.
Second order triangular interpolation in 2D (6 nodes).
Definition: fei2dtrquad.h:44
virtual void boundaryEdgeLocal2Global(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0
Maps the local boundary coordinates to global.
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
virtual void computeValueAt(FloatArray &answer, TimeStep *tStep, const FloatArray &coords, ValueModeType mode)
Computes components values of load at given point - global coordinates (coordinates given)...
Definition: boundaryload.C:58
virtual int SetUpPointsOnLine(int nPoints, MaterialMode mode)
Sets up receiver&#39;s integration points on unit line integration domain.
virtual void edgeEvalN(FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of edge interpolation functions (shape functions) at given point.
Definition: fei2dtrquad.C:146
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
void add(const FloatArray &src)
Adds array src to receiver.
Definition: floatarray.C:156
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
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