OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
planstrss.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/planstrss.h"
36 #include "../sm/Materials/structuralms.h"
37 #include "fei2dquadlin.h"
38 #include "node.h"
39 #include "crosssection.h"
40 #include "gausspoint.h"
41 #include "gaussintegrationrule.h"
42 #include "floatmatrix.h"
43 #include "floatarray.h"
44 #include "intarray.h"
45 #include "domain.h"
46 #include "mathfem.h"
47 #include "strainvector.h"
48 #include "classfactory.h"
49 
50 #ifdef __OOFEG
51  #include "oofeggraphiccontext.h"
52  #include "oofegutils.h"
53  #include "connectivitytable.h"
54  #include "Materials/rcm2.h"
55 #endif
56 
57 namespace oofem {
58 REGISTER_Element(PlaneStress2d);
59 
60 FEI2dQuadLin PlaneStress2d :: interpolation(1, 2);
61 
66  // Constructor.
67 {
68  numberOfDofMans = 4;
70 }
71 
73 // Destructor
74 { }
75 
77 
78 void
80 //
81 // Returns the [3x8] strain-displacement matrix {B} of the receiver,
82 // evaluated at gp.
83 // (epsilon_x,epsilon_y,gamma_xy) = B . r
84 // r = ( u1,v1,u2,v2,u3,v3,u4,v4)
85 {
86  FloatMatrix dnx;
87 
89 
90  answer.resize(3, 8);
91  answer.zero();
92 
93  for ( int i = 1; i <= 4; i++ ) {
94  answer.at(1, 2 * i - 1) = dnx.at(i, 1);
95  answer.at(2, 2 * i - 0) = dnx.at(i, 2);
96  }
97 
98 #ifdef PlaneStress2d_reducedShearIntegration
99  this->interpolation.evaldNdx( dnx, {0., 0.}, *this->giveCellGeometryWrapper() );
100 #endif
101 
102  for ( int i = 1; i <= 4; i++ ) {
103  answer.at(3, 2 * i - 1) = dnx.at(i, 2);
104  answer.at(3, 2 * i - 0) = dnx.at(i, 1);
105  }
106 }
107 
108 
109 void
111 //
112 // Returns the [4x8] displacement gradient matrix {BH} of the receiver,
113 // evaluated at gp.
114 // @todo not checked if correct
115 {
116  FloatMatrix dnx;
117 
119 
120  answer.resize(4, 8);
121 
122  for ( int i = 1; i <= 4; i++ ) {
123  answer.at(1, 2 * i - 1) = dnx.at(i, 1); // du/dx -1
124  answer.at(2, 2 * i - 0) = dnx.at(i, 2); // dv/dy -2
125  }
126 
127 #ifdef PlaneStress2d_reducedShearIntegration
128  this->interpolation.evaldNdx( dnx, {0., 0.}, *this->giveCellGeometryWrapper() );
129 #endif
130 
131  for ( int i = 1; i <= 4; i++ ) {
132  answer.at(3, 2 * i - 1) = dnx.at(i, 2); // du/dy -6
133  answer.at(4, 2 * i - 0) = dnx.at(i, 1); // dv/dx -9
134  }
135 }
136 
137 
140 {
143  if ( result != IRRT_OK ) {
144  return result;
145  }
146 
149  OOFEM_WARNING("Number of Gauss points enforced to 4");
150  }
151 
152  return IRRT_OK;
153 }
154 
155 
156 
157 
158 double
160 //
161 // returns receiver's characteristic size at gp (for some material models)
162 // for crack formed in plane with normal normalToCrackPlane
163 // using the selected method
164 //
165 {
166  if ( method == ECSM_SquareRootOfArea ) {
167  // square root of element area
168  return this->computeMeanSize();
169  }
170 
171  if ( method == ECSM_Projection ) {
172  // standard projection method
173  return this->giveCharacteristicLength(normalToCrackPlane);
174  }
175 
176  // evaluate average strain and its maximum principal direction
177  FloatArray sumstrain, averageNormal;
178  for ( GaussPoint *gpi: *this->giveDefaultIntegrationRulePtr() ) {
179  StructuralMaterialStatus *matstatus = dynamic_cast< StructuralMaterialStatus * >( gpi->giveMaterialStatus() );
180  if ( matstatus ) {
181  sumstrain.add( matstatus->giveTempStrainVector() );
182  }
183  }
184 
185  StrainVector sumstrainvec(sumstrain, _PlaneStress);
186  sumstrainvec.computeMaxPrincipalDir(averageNormal);
187 
188  if ( method == ECSM_ProjectionCentered ) {
189  // projection method based on principal direction of average strain
190  normalToCrackPlane = averageNormal;
191  return this->giveLengthInDir(normalToCrackPlane);
192  }
193 
194  if ( method == ECSM_Oliver1 || method == ECSM_Oliver1modified ) {
195  // method based on derivative of auxiliary function phi at each Gauss point
196  // in the maximum principal strain direction determined at
197  // ECSM_Oliver1 ... at each Gauss point
198  // ECSM_Oliver1modified ... at element center (from average strain)
199 
200  // coordinates of the element center
201  FloatArray center(2);
202  double cx = 0., cy = 0.;
203  for ( int i = 1; i <= 4; i++ ) {
204  cx += giveNode(i)->giveCoordinate(1);
205  cy += giveNode(i)->giveCoordinate(2);
206  }
207 
208  cx /= 4.;
209  cy /= 4.;
210 
211  // nodal values of function phi (0 or 1)
212  FloatArray phi(4);
213  for ( int i = 1; i <= 4; i++ ) {
214  if ( ( ( giveNode(i)->giveCoordinate(1) - cx ) * averageNormal.at(1) + ( giveNode(i)->giveCoordinate(2) - cy ) * averageNormal.at(2) ) > 0. ) {
215  phi.at(i) = 1.;
216  } else {
217  phi.at(i) = 0.;
218  }
219  }
220 
221  // gradient of function phi at the current GP
222  FloatMatrix dnx;
224  FloatArray gradPhi(2);
225  gradPhi.zero();
226  for ( int i = 1; i <= 4; i++ ) {
227  gradPhi.at(1) += phi.at(i) * dnx.at(i, 1);
228  gradPhi.at(2) += phi.at(i) * dnx.at(i, 2);
229  }
230 
231  // scalar product of the gradient with crack normal (at GP)
232  double dPhidN = 0.;
233  if ( method == ECSM_Oliver1modified ) {
234  normalToCrackPlane = averageNormal;
235  }
236 
237  for ( int i = 1; i <= 2; i++ ) {
238  dPhidN += gradPhi.at(i) * normalToCrackPlane.at(i);
239  }
240 
241  if ( dPhidN == 0. ) {
242  OOFEM_ERROR("Zero value of dPhidN");
243  }
244 
245  return 1. / fabs(dPhidN);
246  }
247 
248  OOFEM_ERROR("invalid method");
249  return 0.;
250 }
251 
252 
253 
254 Interface *
256 {
257  if ( interface == ZZNodalRecoveryModelInterfaceType ) {
258  return static_cast< ZZNodalRecoveryModelInterface * >(this);
259  } else if ( interface == SPRNodalRecoveryModelInterfaceType ) {
260  return static_cast< SPRNodalRecoveryModelInterface * >(this);
261  } else if ( interface == SpatialLocalizerInterfaceType ) {
262  return static_cast< SpatialLocalizerInterface * >(this);
263  } else if ( interface == HuertaErrorEstimatorInterfaceType ) {
264  return static_cast< HuertaErrorEstimatorInterface * >(this);
265  }
266 
267  return NULL;
268 }
269 
270 
271 void
273  IntArray &localNodeIdArray, IntArray &globalNodeIdArray,
275  int &localNodeId, int &localElemId, int &localBcId,
276  IntArray &controlNode, IntArray &controlDof,
278 {
279  int inode, nodes = 4, iside, sides = 4, nd1, nd2;
280  FloatArray *corner [ 4 ], midSide [ 4 ], midNode, cor [ 4 ];
281  double x = 0.0, y = 0.0;
282 
283  static int sideNode [ 4 ] [ 2 ] = { { 1, 2 }, { 2, 3 }, { 3, 4 }, { 4, 1 } };
284 
287  for ( inode = 0; inode < nodes; inode++ ) {
288  corner [ inode ] = this->giveNode(inode + 1)->giveCoordinates();
289  if ( corner [ inode ]->giveSize() != 3 ) {
290  cor [ inode ].resize(3);
291  cor [ inode ].at(1) = corner [ inode ]->at(1);
292  cor [ inode ].at(2) = corner [ inode ]->at(2);
293  cor [ inode ].at(3) = 0.0;
294 
295  corner [ inode ] = & ( cor [ inode ] );
296  }
297 
298  x += corner [ inode ]->at(1);
299  y += corner [ inode ]->at(2);
300  }
301 
302  for ( iside = 0; iside < sides; iside++ ) {
303  midSide [ iside ].resize(3);
304 
305  nd1 = sideNode [ iside ] [ 0 ] - 1;
306  nd2 = sideNode [ iside ] [ 1 ] - 1;
307 
308  midSide [ iside ].at(1) = ( corner [ nd1 ]->at(1) + corner [ nd2 ]->at(1) ) / 2.0;
309  midSide [ iside ].at(2) = ( corner [ nd1 ]->at(2) + corner [ nd2 ]->at(2) ) / 2.0;
310  midSide [ iside ].at(3) = 0.0;
311  }
312 
313  midNode.resize(3);
314 
315  midNode.at(1) = x / nodes;
316  midNode.at(2) = y / nodes;
317  midNode.at(3) = 0.0;
318  }
319 
320  this->setupRefinedElementProblem2D(this, refinedElement, level, nodeId, localNodeIdArray, globalNodeIdArray,
321  sMode, tStep, nodes, corner, midSide, midNode,
322  localNodeId, localElemId, localBcId,
323  controlNode, controlDof, aMode, "PlaneStress2d");
324 }
325 
327 {
329 }
330 
331 
332 #ifdef __OOFEG
333  #define TR_LENGHT_REDUCT 0.3333
334 
336 {
337  WCRec p [ 4 ];
338  GraphicObj *go;
339 
340  if ( !gc.testElementGraphicActivity(this) ) {
341  return;
342  }
343 
344  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
345  EASValsSetColor( gc.getElementColor() );
346  EASValsSetEdgeColor( gc.getElementEdgeColor() );
347  EASValsSetEdgeFlag(true);
348  EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
349  EASValsSetFillStyle(FILL_HOLLOW);
350  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
351  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
352  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
353  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
354  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
355  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
356  p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveCoordinate(1);
357  p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveCoordinate(2);
358  p [ 2 ].z = ( FPNum ) this->giveNode(3)->giveCoordinate(3);
359  p [ 3 ].x = ( FPNum ) this->giveNode(4)->giveCoordinate(1);
360  p [ 3 ].y = ( FPNum ) this->giveNode(4)->giveCoordinate(2);
361  p [ 3 ].z = ( FPNum ) this->giveNode(4)->giveCoordinate(3);
362 
363  go = CreateQuad3D(p);
364  EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
365  EGAttachObject(go, ( EObjectP ) this);
366  EMAddGraphicsToModel(ESIModel(), go);
367 }
368 
369 
371 {
372  WCRec p [ 4 ];
373  GraphicObj *go;
374  double defScale = gc.getDefScale();
375 
376  if ( !gc.testElementGraphicActivity(this) ) {
377  return;
378  }
379 
380  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
381  EASValsSetColor( gc.getDeformedElementColor() );
382  EASValsSetEdgeColor( gc.getElementEdgeColor() );
383  EASValsSetEdgeFlag(true);
384  EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
385  EASValsSetFillStyle(FILL_HOLLOW);
386  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
387  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
388  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
389  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
390  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
391  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
392  p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveUpdatedCoordinate(1, tStep, defScale);
393  p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveUpdatedCoordinate(2, tStep, defScale);
394  p [ 2 ].z = ( FPNum ) this->giveNode(3)->giveUpdatedCoordinate(3, tStep, defScale);
395  p [ 3 ].x = ( FPNum ) this->giveNode(4)->giveUpdatedCoordinate(1, tStep, defScale);
396  p [ 3 ].y = ( FPNum ) this->giveNode(4)->giveUpdatedCoordinate(2, tStep, defScale);
397  p [ 3 ].z = ( FPNum ) this->giveNode(4)->giveUpdatedCoordinate(3, tStep, defScale);
398 
399  go = CreateQuad3D(p);
400  EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
401  EMAddGraphicsToModel(ESIModel(), go);
402 }
403 
404 
405 
407 {
408  int i, indx, result = 0;
409  WCRec p [ 4 ];
410  GraphicObj *tr;
411  FloatArray v [ 4 ];
412  double s [ 4 ], defScale;
413 
414  if ( !gc.testElementGraphicActivity(this) ) {
415  return;
416  }
417 
418  EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
419  if ( gc.giveIntVarMode() == ISM_recovered ) {
420  for ( i = 1; i <= 4; i++ ) {
421  result += this->giveInternalStateAtNode(v [ i - 1 ], gc.giveIntVarType(), gc.giveIntVarMode(), i, tStep);
422  }
423 
424  if ( result != 4 ) {
425  return;
426  }
427 
428  indx = gc.giveIntVarIndx();
429 
430  for ( i = 1; i <= 4; i++ ) {
431  s [ i - 1 ] = v [ i - 1 ].at(indx);
432  }
433 
434  if ( gc.getScalarAlgo() == SA_ISO_SURF ) {
435  for ( i = 0; i < 4; i++ ) {
436  if ( gc.getInternalVarsDefGeoFlag() ) {
437  // use deformed geometry
438  defScale = gc.getDefScale();
439  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(1, tStep, defScale);
440  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(2, tStep, defScale);
441  p [ i ].z = 0.;
442  } else {
443  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
444  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
445  p [ i ].z = 0.;
446  }
447  }
448 
449  //EASValsSetColor(gc.getYieldPlotColor(ratio));
450  gc.updateFringeTableMinMax(s, 4);
451  tr = CreateQuadWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ], s [ 3 ]);
452  EGWithMaskChangeAttributes(LAYER_MASK, tr);
453  EMAddGraphicsToModel(ESIModel(), tr);
454 
455  /* } else if (gc.getScalarAlgo() == SA_ISO_LINE) {
456  *
457  * EASValsSetColor(context.getActiveCrackColor());
458  * EASValsSetLineWidth(OOFEG_ISO_LINE_WIDTH);
459  *
460  * for (i=0; i< 4; i++) {
461  * if (gc.getInternalVarsDefGeoFlag()) {
462  * // use deformed geometry
463  * defScale = gc.getDefScale();
464  * p[i].x = (FPNum) this->giveNode(i+1)->giveUpdatedCoordinate(1,tStep,defScale);
465  * p[i].y = (FPNum) this->giveNode(i+1)->giveUpdatedCoordinate(2,tStep,defScale);
466  * p[i].z = 0.;
467  *
468  * } else {
469  * p[i].x = (FPNum) this->giveNode(i+1)->giveCoordinate(1);
470  * p[i].y = (FPNum) this->giveNode(i+1)->giveCoordinate(2);
471  * p[i].z = 0.;
472  * }
473  * }
474  *
475  * // isoline implementation
476  * oofeg_drawIsoLinesOnQuad (p, s);
477  */
478  } else if ( ( gc.getScalarAlgo() == SA_ZPROFILE ) || ( gc.getScalarAlgo() == SA_COLORZPROFILE ) ) {
479  double landScale = gc.getLandScale();
480 
481  for ( i = 0; i < 4; i++ ) {
482  if ( gc.getInternalVarsDefGeoFlag() ) {
483  // use deformed geometry
484  defScale = gc.getDefScale();
485  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(1, tStep, defScale);
486  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(2, tStep, defScale);
487  p [ i ].z = s [ i ] * landScale;
488  } else {
489  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
490  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
491  p [ i ].z = s [ i ] * landScale;
492  }
493 
494  // this fixes a bug in ELIXIR
495  if ( fabs(s [ i ]) < 1.0e-6 ) {
496  s [ i ] = 1.0e-6;
497  }
498  }
499 
500  if ( gc.getScalarAlgo() == SA_ZPROFILE ) {
501  EASValsSetColor( gc.getDeformedElementColor() );
502  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
503  tr = CreateQuad3D(p);
504  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, tr);
505  } else {
506  gc.updateFringeTableMinMax(s, 4);
507  tr = CreateQuadWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ], s [ 3 ]);
508  EGWithMaskChangeAttributes(LAYER_MASK, tr);
509  }
510 
511  EMAddGraphicsToModel(ESIModel(), tr);
512  }
513  } else if ( gc.giveIntVarMode() == ISM_local ) {
514  if ( integrationRulesArray [ 0 ]->giveNumberOfIntegrationPoints() != 4 ) {
515  return;
516  }
517 
518  IntArray ind(4);
519  WCRec pp [ 9 ];
520 
521  for ( i = 0; i < 4; i++ ) {
522  if ( gc.getInternalVarsDefGeoFlag() ) {
523  // use deformed geometry
524  defScale = gc.getDefScale();
525  pp [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(1, tStep, defScale);
526  pp [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(2, tStep, defScale);
527  pp [ i ].z = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(3, tStep, defScale);
528  } else {
529  pp [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
530  pp [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
531  pp [ i ].z = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(3);
532  }
533  }
534 
535  for ( i = 0; i < 3; i++ ) {
536  pp [ i + 4 ].x = 0.5 * ( pp [ i ].x + pp [ i + 1 ].x );
537  pp [ i + 4 ].y = 0.5 * ( pp [ i ].y + pp [ i + 1 ].y );
538  pp [ i + 4 ].z = 0.5 * ( pp [ i ].z + pp [ i + 1 ].z );
539  }
540 
541  pp [ 7 ].x = 0.5 * ( pp [ 3 ].x + pp [ 0 ].x );
542  pp [ 7 ].y = 0.5 * ( pp [ 3 ].y + pp [ 0 ].y );
543  pp [ 7 ].z = 0.5 * ( pp [ 3 ].z + pp [ 0 ].z );
544 
545  pp [ 8 ].x = 0.25 * ( pp [ 0 ].x + pp [ 1 ].x + pp [ 2 ].x + pp [ 3 ].x );
546  pp [ 8 ].y = 0.25 * ( pp [ 0 ].y + pp [ 1 ].y + pp [ 2 ].y + pp [ 3 ].y );
547  pp [ 8 ].z = 0.25 * ( pp [ 0 ].z + pp [ 1 ].z + pp [ 2 ].z + pp [ 3 ].z );
548 
549  for ( GaussPoint *gp: *this->giveDefaultIntegrationRulePtr() ) {
550  const FloatArray& gpCoords = gp->giveNaturalCoordinates();
551  if ( ( gpCoords.at(1) > 0. ) && ( gpCoords.at(2) > 0. ) ) {
552  ind.at(1) = 0;
553  ind.at(2) = 4;
554  ind.at(3) = 8;
555  ind.at(4) = 7;
556  } else if ( ( gpCoords.at(1) < 0. ) && ( gpCoords.at(2) > 0. ) ) {
557  ind.at(1) = 4;
558  ind.at(2) = 1;
559  ind.at(3) = 5;
560  ind.at(4) = 8;
561  } else if ( ( gpCoords.at(1) < 0. ) && ( gpCoords.at(2) < 0. ) ) {
562  ind.at(1) = 5;
563  ind.at(2) = 2;
564  ind.at(3) = 6;
565  ind.at(4) = 8;
566  } else {
567  ind.at(1) = 6;
568  ind.at(2) = 3;
569  ind.at(3) = 7;
570  ind.at(4) = 8;
571  }
572 
573  if ( giveIPValue(v [ 0 ], gp, gc.giveIntVarType(), tStep) == 0 ) {
574  return;
575  }
576 
577  indx = gc.giveIntVarIndx();
578 
579  for ( i = 1; i <= 4; i++ ) {
580  s [ i - 1 ] = v [ 0 ].at(indx);
581  }
582 
583  for ( i = 0; i < 4; i++ ) {
584  p [ i ].x = pp [ ind.at(i + 1) ].x;
585  p [ i ].y = pp [ ind.at(i + 1) ].y;
586  p [ i ].z = pp [ ind.at(i + 1) ].z;
587  }
588 
589  gc.updateFringeTableMinMax(s, 4);
590  tr = CreateQuadWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ], s [ 3 ]);
591  EGWithMaskChangeAttributes(LAYER_MASK, tr);
592  EMAddGraphicsToModel(ESIModel(), tr);
593  }
594  }
595 }
596 
597 
598 
599 
600 void
602 {
603  int i;
604  WCRec l [ 2 ];
605  GraphicObj *tr;
606  double defScale = gc.getDefScale();
607  FloatArray crackStatuses, cf;
608 
609  if ( !gc.testElementGraphicActivity(this) ) {
610  return;
611  }
612 
613  if ( gc.giveIntVarType() == IST_CrackState ) {
614  // ask if any active crack exist
615  int crackStatus;
616  double ax, ay, bx, by, norm, xc, yc, length;
617  FloatArray crackDir;
618  FloatArray gpglobalcoords;
619 
620  for ( GaussPoint *gp: *this->giveDefaultIntegrationRulePtr() ) {
621 
622  if ( this->giveIPValue(cf, gp, IST_CrackedFlag, tStep) == 0 ) {
623  return;
624  }
625 
626  if ( ( int ) cf.at(1) == 0 ) {
627  return;
628  }
629 
630  if ( this->giveIPValue(crackDir, gp, IST_CrackDirs, tStep) ) {
631  this->giveIPValue(crackStatuses, gp, IST_CrackStatuses, tStep);
632  for ( i = 1; i <= 3; i++ ) {
633  crackStatus = ( int ) crackStatuses.at(i);
634  if ( ( crackStatus != pscm_NONE ) && ( crackStatus != pscm_CLOSED ) ) {
635  // draw a crack
636  // this element is 2d element in x-y plane
637  //
638  // compute perpendicular line to normal in xy plane
639  ax = crackDir.at(i);
640  ay = crackDir.at(3 + i);
641  if ( fabs(ax) > 1.e-6 ) {
642  by = 1.;
643  bx = -ay * by / ax;
644  norm = sqrt(bx * bx + by * by);
645  bx = bx / norm; // normalize to obtain unit vector
646  by = by / norm;
647  } else {
648  by = 0.0;
649  bx = 1.0;
650  }
651 
652  // obtain gp global coordinates
653  if ( gc.getInternalVarsDefGeoFlag() ) {
654  double ksi, eta, n1, n2, n3, n4;
655  ksi = gp->giveNaturalCoordinate(1);
656  eta = gp->giveNaturalCoordinate(2);
657 
658  n1 = ( 1. + ksi ) * ( 1. + eta ) * 0.25;
659  n2 = ( 1. - ksi ) * ( 1. + eta ) * 0.25;
660  n3 = ( 1. - ksi ) * ( 1. - eta ) * 0.25;
661  n4 = ( 1. + ksi ) * ( 1. - eta ) * 0.25;
662 
663  gpglobalcoords.resize(2);
664 
665  gpglobalcoords.at(1) = ( n1 * this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale) +
666  n2 * this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale) +
667  n3 * this->giveNode(3)->giveUpdatedCoordinate(1, tStep, defScale) +
668  n4 * this->giveNode(4)->giveUpdatedCoordinate(1, tStep, defScale) );
669  gpglobalcoords.at(2) = ( n1 * this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale) +
670  n2 * this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale) +
671  n3 * this->giveNode(3)->giveUpdatedCoordinate(2, tStep, defScale) +
672  n4 * this->giveNode(4)->giveUpdatedCoordinate(2, tStep, defScale) );
673  } else {
674  computeGlobalCoordinates( gpglobalcoords, gp->giveNaturalCoordinates() );
675  }
676 
677  xc = gpglobalcoords.at(1);
678  yc = gpglobalcoords.at(2);
679  length = sqrt( computeVolumeAround(gp) / this->giveCrossSection()->give(CS_Thickness, gp) ) / 2.;
680  l [ 0 ].x = ( FPNum ) xc + bx * length;
681  l [ 0 ].y = ( FPNum ) yc + by * length;
682  l [ 0 ].z = 0.;
683  l [ 1 ].x = ( FPNum ) xc - bx * length;
684  l [ 1 ].y = ( FPNum ) yc - by * length;
685  l [ 1 ].z = 0.;
686  EASValsSetLayer(OOFEG_CRACK_PATTERN_LAYER);
687  EASValsSetLineWidth(OOFEG_CRACK_PATTERN_WIDTH);
688  if ( ( crackStatus == pscm_SOFTENING ) || ( crackStatus == pscm_OPEN ) ) {
689  EASValsSetColor( gc.getActiveCrackColor() );
690  } else {
691  EASValsSetColor( gc.getCrackPatternColor() );
692  }
693 
694  tr = CreateLine3D(l);
695  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, tr);
696  EMAddGraphicsToModel(ESIModel(), tr);
697  }
698  }
699  }
700  }
701  }
702 }
703 
704 #endif
705 
706 void
708 {
709  pap.resize(4);
710  for ( int i = 1; i < 5; i++ ) {
711  pap.at(i) = this->giveNode(i)->giveNumber();
712  }
713 }
714 
715 void
717 {
718  int found = 0;
719  answer.resize(1);
720 
721  for ( int i = 1; i < 5; i++ ) {
722  if ( pap == this->giveNode(i)->giveNumber() ) {
723  found = 1;
724  }
725  }
726 
727  if ( found ) {
728  answer.at(1) = pap;
729  } else {
730  OOFEM_ERROR("node unknown");
731  }
732 }
733 
734 int
736 {
738 }
739 
740 
743 {
744  return SPRPatchType_2dxy;
745 }
746 
747 } // end namespace oofem
virtual FEICellGeometry * giveCellGeometryWrapper()
Returns the Cell Geometry Wrapper.
CrossSection * giveCrossSection()
Definition: element.C:495
int testElementGraphicActivity(Element *)
Test if particular element passed fulfills various filtering criteria for its graphics output...
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 void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: planstrss.C:335
Class and object Domain.
Definition: domain.h:115
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Access method for default integration rule.
Definition: element.h:822
virtual void drawSpecial(oofegGraphicContext &gc, TimeStep *tStep)
Definition: planstrss.C:601
ScalarAlgorithmType getScalarAlgo()
double computeMeanSize()
Computes the size of the element defined as its length.
Definition: element.C:1078
void computeMaxPrincipalDir(FloatArray &answer) const
Computes the principal direction of the receiver associated with the maximum principal value...
Definition: strainvector.C:182
The element interface required by ZZNodalRecoveryModel.
const FloatArray & giveSubPatchCoordinates()
Returns local sub-patch coordinates of the receiver.
Definition: gausspoint.h:144
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
virtual int giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *tStep)
Returns internal state variable (like stress,strain) at node of element in Reduced form...
#define pscm_CLOSED
Definition: fcm.h:58
#define OOFEG_RAW_GEOMETRY_LAYER
This class implements a structural material status information.
Definition: structuralms.h:65
#define pscm_NONE
Definition: fcm.h:54
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual ~PlaneStress2d()
Definition: planstrss.C:72
Specialization of a floating point array for representing a strain state.
Definition: strainvector.h:45
virtual double giveCoordinate(int i)
Definition: node.C:82
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
virtual double giveCharacteristicLength(const FloatArray &normalToCrackPlane)
Returns the size of element in the given direction, in some cases adjusted (e.g.
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Definition: planstrss.C:406
#define OOFEG_DEFORMED_GEOMETRY_LAYER
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
InternalStateType giveIntVarType()
#define OOFEG_CRACK_PATTERN_LAYER
#define pscm_OPEN
Definition: rcm2.h:61
The element interface corresponding to HuertaErrorEstimator.
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
Definition: planstrss.C:716
#define OOFEG_RAW_GEOMETRY_WIDTH
UnknownType
Type representing particular unknown (its physical meaning).
Definition: unknowntype.h:55
#define OOFEG_CRACK_PATTERN_WIDTH
virtual double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Returns updated ic-th coordinate of receiver.
Definition: node.C:245
PlaneStress2d(int n, Domain *d)
Definition: planstrss.C:62
const FloatArray & giveTempStrainVector() const
Returns the const pointer to receiver&#39;s temporary strain vector.
Definition: structuralms.h:115
virtual double giveLengthInDir(const FloatArray &normalToCrackPlane)
Default implementation returns length of element projection into specified direction.
Definition: element.C:1143
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
static FEI2dQuadLin interpolation
Definition: planstrss.h:61
virtual void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
Definition: planstrss.C:707
int numberOfGaussPoints
Number of integration points as specified by nip.
Definition: element.h:188
InternalStateMode giveIntVarMode()
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: planstrss.C:139
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: planstrss.C:255
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
void setupRefinedElementProblem2D(Element *element, RefinedElement *refinedElement, int level, int nodeId, IntArray &localNodeIdArray, IntArray &globalNodeIdArray, HuertaErrorEstimatorInterface::SetupMode mode, TimeStep *tStep, int nodes, FloatArray **corner, FloatArray *midSide, FloatArray &midNode, int &localNodeId, int &localElemId, int &localBcId, IntArray &controlNode, IntArray &controlDof, HuertaErrorEstimator::AnalysisMode aMode, const char *quadtype)
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
Definition: planstrss.C:370
virtual int SPRNodalRecoveryMI_giveNumberOfIP()
Definition: planstrss.C:735
double norm(const FloatArray &x)
Definition: floatarray.C:985
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
Class Interface.
Definition: interface.h:82
int giveNumberOfIntegrationPoints() const
Returns number of integration points of receiver.
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
The spatial localizer element interface associated to spatial localizer.
virtual FEInterpolation * giveInterpolation() const
Definition: planstrss.C:76
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 void computeBHmatrixAt(GaussPoint *gp, FloatMatrix &answer)
Computes a matrix which, multiplied by the column matrix of nodal displacements, gives the displaceme...
Definition: planstrss.C:110
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
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
void updateFringeTableMinMax(double *s, int size)
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.
int giveNumber() const
Definition: femcmpnn.h:107
virtual SPRPatchType SPRNodalRecoveryMI_givePatchType()
Definition: planstrss.C:742
virtual void HuertaErrorEstimatorI_setupRefinedElementProblem(RefinedElement *refinedElement, int level, int nodeId, IntArray &localNodeIdArray, IntArray &globalNodeIdArray, HuertaErrorEstimatorInterface::SetupMode sMode, TimeStep *tStep, int &localNodeId, int &localElemId, int &localBcId, IntArray &controlNode, IntArray &controlDof, HuertaErrorEstimator::AnalysisMode aMode)
Definition: planstrss.C:272
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
#define pscm_SOFTENING
Definition: fcm.h:56
#define OOFEG_VARPLOT_PATTERN_LAYER
Class representing integration point in finite element program.
Definition: gausspoint.h:93
#define OOFEM_WARNING(...)
Definition: error.h:62
Class representing solution step.
Definition: timestep.h:80
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element&#39;s local coordinates.
Definition: element.C:1207
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
virtual void HuertaErrorEstimatorI_computeNmatrixAt(GaussPoint *gp, FloatMatrix &answer)
Definition: planstrss.C:326
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Definition: planstrss.C:79
virtual double giveCharacteristicSize(GaussPoint *gp, FloatArray &normalToCrackPlane, ElementCharSizeMethod method)
Returns characteristic element size for a given integration point and given direction.
Definition: planstrss.C:159

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