OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
quad1planestrain.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 
36 #include "fei2dquadlin.h"
37 #include "node.h"
38 #include "crosssection.h"
39 #include "gausspoint.h"
40 #include "gaussintegrationrule.h"
41 #include "floatmatrix.h"
42 #include "floatarray.h"
43 #include "intarray.h"
44 #include "mathfem.h"
45 #include "classfactory.h"
46 
47 #ifdef __OOFEG
48  #include "oofeggraphiccontext.h"
49  #include "oofegutils.h"
50  #include "Materials/rcm2.h"
51 #endif
52 
53 namespace oofem {
54 REGISTER_Element(Quad1PlaneStrain);
55 
56 FEI2dQuadLin Quad1PlaneStrain :: interp(1, 2);
57 
62 {
63  numberOfDofMans = 4;
65 }
66 
67 
69 { }
70 
71 
73 
74 
75 void
77 //
78 // Returns the [4x8] strain-displacement matrix {B} of the receiver,
79 // evaluated at gp.
80 // (epsilon_x,epsilon_y,epsilon_z,gamma_xy) = B . r
81 // r = ( u1,v1,u2,v2,u3,v3,u4,v4)
82 {
83  FloatMatrix dN;
84 
86 
87  // Reshape
88  answer.resize(4, 8);
89  answer.zero();
90 
91 #ifdef Quad1PlaneStrain_reducedVolumetricIntegration
92  FloatMatrix dN_red;
93  FloatArray lcoords_red(2);
94  lcoords_red.zero();
95  this->interp.evaldNdx( dN_red, lcoords_red, FEIElementGeometryWrapper(this) );
96 
97  for ( int i = 1; i <= 4; i++ ) {
98  answer.at(1, 2 * i - 1) = ( dN.at(i, 1) + dN_red.at(i, 1) ) / 2.;
99  answer.at(1, 2 * i - 0) = (-dN.at(i, 2) + dN_red.at(i, 2) ) / 2.;
100  answer.at(2, 2 * i - 1) = (-dN.at(i, 1) + dN_red.at(i, 1) ) / 2.;
101  answer.at(2, 2 * i - 0) = ( dN.at(i, 2) + dN_red.at(i, 2) ) / 2.;
102  answer.at(4, 2 * i - 1) = dN.at(i, 2);
103  answer.at(4, 2 * i - 0) = dN.at(i, 1);
104  }
105 #else
106 #ifdef Quad1PlaneStrain_reducedShearIntegration
107  FloatMatrix dN_red;
108  FloatArray lcoords_red(2);
109  lcoords_red.zero();
110  this->interp.evaldNdx( dN_red, lcoords_red, FEIElementGeometryWrapper(this) );
111 
112  for ( int i = 1; i <= 4; i++ ) {
113  answer.at(1, 2 * i - 1) = dN.at(i, 1);
114  answer.at(2, 2 * i - 0) = dN.at(i, 2);
115  answer.at(4, 2 * i - 1) = dN_red.at(i, 2);
116  answer.at(4, 2 * i - 0) = dN_red.at(i, 1);
117  }
118 
119 #else
120  for ( int i = 1; i <= 4; i++ ) {
121  answer.at(1, 2 * i - 1) = dN.at(i, 1);
122  answer.at(2, 2 * i - 0) = dN.at(i, 2);
123  answer.at(4, 2 * i - 1) = dN.at(i, 2);
124  answer.at(4, 2 * i - 0) = dN.at(i, 1);
125  }
126 
127 #endif
128 #endif
129 }
130 
131 void
133 // Returns the [5x8] displacement gradient matrix {BH} of the receiver,
134 // evaluated at gp.
135 {
136  FloatMatrix dnx;
138 
139  answer.resize(5, 8);
140  answer.zero();
141  // 3rd row is zero -> dw/dz = 0
142  for ( int i = 1; i <= 4; i++ ) {
143  answer.at(1, 2 * i - 1) = dnx.at(i, 1); // du/dx -1
144  answer.at(2, 2 * i - 0) = dnx.at(i, 2); // dv/dy -2
145  answer.at(4, 2 * i - 1) = dnx.at(i, 2); // du/dy -6
146  answer.at(5, 2 * i - 0) = dnx.at(i, 1); // dv/dx -9
147  }
148 }
149 
150 
153 {
156  if ( result != IRRT_OK ) {
157  return result;
158  }
159 
160  if ( !( ( numberOfGaussPoints == 4 ) ||
161  ( numberOfGaussPoints == 1 ) ||
162  ( numberOfGaussPoints == 9 ) ||
163  ( numberOfGaussPoints == 16 ) ) ) {
165  }
166 
167  return IRRT_OK;
168 }
169 
170 
171 
172 Interface *
174 {
175  if ( interface == ZZNodalRecoveryModelInterfaceType ) {
176  return static_cast< ZZNodalRecoveryModelInterface * >(this);
177  } else if ( interface == SPRNodalRecoveryModelInterfaceType ) {
178  return static_cast< SPRNodalRecoveryModelInterface * >(this);
179  } else if ( interface == SpatialLocalizerInterfaceType ) {
180  return static_cast< SpatialLocalizerInterface * >(this);
181  } else if ( interface == HuertaErrorEstimatorInterfaceType ) {
182  return static_cast< HuertaErrorEstimatorInterface * >(this);
183  }
184 
185  return NULL;
186 }
187 
188 
189 void
191  IntArray &localNodeIdArray, IntArray &globalNodeIdArray,
193  int &localNodeId, int &localElemId, int &localBcId,
194  IntArray &controlNode, IntArray &controlDof,
196 {
197  int inode, nodes = 4, iside, sides = 4, nd1, nd2;
198  FloatArray *corner [ 4 ], midSide [ 4 ], midNode, cor [ 4 ];
199  double x = 0.0, y = 0.0;
200 
201  static int sideNode [ 4 ] [ 2 ] = { { 1, 2 }, { 2, 3 }, { 3, 4 }, { 4, 1 } };
202 
205  for ( inode = 0; inode < nodes; inode++ ) {
206  corner [ inode ] = this->giveNode(inode + 1)->giveCoordinates();
207  if ( corner [ inode ]->giveSize() != 3 ) {
208  cor [ inode ].resize(3);
209  cor [ inode ].at(1) = corner [ inode ]->at(1);
210  cor [ inode ].at(2) = corner [ inode ]->at(2);
211  cor [ inode ].at(3) = 0.0;
212 
213  corner [ inode ] = & ( cor [ inode ] );
214  }
215 
216  x += corner [ inode ]->at(1);
217  y += corner [ inode ]->at(2);
218  }
219 
220  for ( iside = 0; iside < sides; iside++ ) {
221  midSide [ iside ].resize(3);
222 
223  nd1 = sideNode [ iside ] [ 0 ] - 1;
224  nd2 = sideNode [ iside ] [ 1 ] - 1;
225 
226  midSide [ iside ].at(1) = ( corner [ nd1 ]->at(1) + corner [ nd2 ]->at(1) ) / 2.0;
227  midSide [ iside ].at(2) = ( corner [ nd1 ]->at(2) + corner [ nd2 ]->at(2) ) / 2.0;
228  midSide [ iside ].at(3) = 0.0;
229  }
230 
231  midNode.resize(3);
232 
233  midNode.at(1) = x / nodes;
234  midNode.at(2) = y / nodes;
235  midNode.at(3) = 0.0;
236  }
237 
238  this->setupRefinedElementProblem2D(this, refinedElement, level, nodeId, localNodeIdArray, globalNodeIdArray,
239  sMode, tStep, nodes, corner, midSide, midNode,
240  localNodeId, localElemId, localBcId,
241  controlNode, controlDof, aMode, "Quad1PlaneStrain");
242 }
243 
244 
246 {
248 }
249 
250 
251 
252 #ifdef __OOFEG
253  #define TR_LENGHT_REDUCT 0.3333
254 
256 {
257  WCRec p [ 4 ];
258  GraphicObj *go;
259 
260  if ( !gc.testElementGraphicActivity(this) ) {
261  return;
262  }
263 
264  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
265  EASValsSetColor( gc.getElementColor() );
266  EASValsSetEdgeColor( gc.getElementEdgeColor() );
267  EASValsSetEdgeFlag(true);
268  EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
269  EASValsSetFillStyle(FILL_HOLLOW);
270  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
271  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
272  p [ 0 ].z = 0.;
273  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
274  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
275  p [ 1 ].z = 0.;
276  p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveCoordinate(1);
277  p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveCoordinate(2);
278  p [ 2 ].z = 0.;
279  p [ 3 ].x = ( FPNum ) this->giveNode(4)->giveCoordinate(1);
280  p [ 3 ].y = ( FPNum ) this->giveNode(4)->giveCoordinate(2);
281  p [ 3 ].z = 0.;
282 
283  go = CreateQuad3D(p);
284  EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
285  EGAttachObject(go, ( EObjectP ) this);
286  EMAddGraphicsToModel(ESIModel(), go);
287 }
288 
289 
291 {
292  WCRec p [ 4 ];
293  GraphicObj *go;
294  double defScale = gc.getDefScale();
295 
296  if ( !gc.testElementGraphicActivity(this) ) {
297  return;
298  }
299 
300  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
301  EASValsSetColor( gc.getDeformedElementColor() );
302  EASValsSetEdgeColor( gc.getElementEdgeColor() );
303  EASValsSetEdgeFlag(true);
304  EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
305  EASValsSetFillStyle(FILL_HOLLOW);
306  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
307  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
308  p [ 0 ].z = 0.;
309  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
310  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
311  p [ 1 ].z = 0.;
312  p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveUpdatedCoordinate(1, tStep, defScale);
313  p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveUpdatedCoordinate(2, tStep, defScale);
314  p [ 2 ].z = 0.;
315  p [ 3 ].x = ( FPNum ) this->giveNode(4)->giveUpdatedCoordinate(1, tStep, defScale);
316  p [ 3 ].y = ( FPNum ) this->giveNode(4)->giveUpdatedCoordinate(2, tStep, defScale);
317  p [ 3 ].z = 0.;
318 
319  go = CreateQuad3D(p);
320  EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
321  EMAddGraphicsToModel(ESIModel(), go);
322 }
323 
324 
325 
327 {
328  int i, indx, result = 0;
329  WCRec p [ 4 ];
330  GraphicObj *tr;
331  FloatArray v [ 4 ];
332  double s [ 4 ], defScale;
333 
334  if ( !gc.testElementGraphicActivity(this) ) {
335  return;
336  }
337 
338  EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
339  if ( gc.giveIntVarMode() == ISM_recovered ) {
340  for ( i = 1; i <= 4; i++ ) {
341  result += this->giveInternalStateAtNode(v [ i - 1 ], gc.giveIntVarType(), gc.giveIntVarMode(), i, tStep);
342  }
343 
344  if ( result != 4 ) {
345  return;
346  }
347 
348  indx = gc.giveIntVarIndx();
349 
350  for ( i = 1; i <= 4; i++ ) {
351  s [ i - 1 ] = v [ i - 1 ].at(indx);
352  }
353 
354  if ( gc.getScalarAlgo() == SA_ISO_SURF ) {
355  for ( i = 0; i < 4; i++ ) {
356  if ( gc.getInternalVarsDefGeoFlag() ) {
357  // use deformed geometry
358  defScale = gc.getDefScale();
359  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(1, tStep, defScale);
360  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(2, tStep, defScale);
361  p [ i ].z = 0.;
362  } else {
363  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
364  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
365  p [ i ].z = 0.;
366  }
367  }
368 
369  //EASValsSetColor(gc.getYieldPlotColor(ratio));
370  gc.updateFringeTableMinMax(s, 4);
371  tr = CreateQuadWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ], s [ 3 ]);
372  EGWithMaskChangeAttributes(LAYER_MASK, tr);
373  EMAddGraphicsToModel(ESIModel(), tr);
374  } else if ( ( gc.getScalarAlgo() == SA_ZPROFILE ) || ( gc.getScalarAlgo() == SA_COLORZPROFILE ) ) {
375  double landScale = gc.getLandScale();
376 
377  for ( i = 0; i < 4; i++ ) {
378  if ( gc.getInternalVarsDefGeoFlag() ) {
379  // use deformed geometry
380  defScale = gc.getDefScale();
381  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(1, tStep, defScale);
382  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(2, tStep, defScale);
383  p [ i ].z = s [ i ] * landScale;
384  } else {
385  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
386  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
387  p [ i ].z = s [ i ] * landScale;
388  }
389 
390  // this fixes a bug in ELIXIR
391  if ( fabs(s [ i ]) < 1.0e-6 ) {
392  s [ i ] = 1.0e-6;
393  }
394  }
395 
396  if ( gc.getScalarAlgo() == SA_ZPROFILE ) {
397  EASValsSetColor( gc.getDeformedElementColor() );
398  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
399  tr = CreateQuad3D(p);
400  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, tr);
401  } else {
402  gc.updateFringeTableMinMax(s, 4);
403  tr = CreateQuadWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ], s [ 3 ]);
404  EGWithMaskChangeAttributes(LAYER_MASK, tr);
405  }
406 
407  EMAddGraphicsToModel(ESIModel(), tr);
408  }
409  } else if ( gc.giveIntVarMode() == ISM_local ) {
410  if ( numberOfGaussPoints != 4 ) {
411  return;
412  }
413 
414  IntArray ind(4);
415  WCRec pp [ 9 ];
416 
417  for ( i = 0; i < 4; i++ ) {
418  if ( gc.getInternalVarsDefGeoFlag() ) {
419  // use deformed geometry
420  defScale = gc.getDefScale();
421  pp [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(1, tStep, defScale);
422  pp [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(2, tStep, defScale);
423  pp [ i ].z = 0.;
424  } else {
425  pp [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
426  pp [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
427  pp [ i ].z = 0.;
428  }
429  }
430 
431  for ( i = 0; i < 3; i++ ) {
432  pp [ i + 4 ].x = 0.5 * ( pp [ i ].x + pp [ i + 1 ].x );
433  pp [ i + 4 ].y = 0.5 * ( pp [ i ].y + pp [ i + 1 ].y );
434  pp [ i + 4 ].z = 0.5 * ( pp [ i ].z + pp [ i + 1 ].z );
435  }
436 
437  pp [ 7 ].x = 0.5 * ( pp [ 3 ].x + pp [ 0 ].x );
438  pp [ 7 ].y = 0.5 * ( pp [ 3 ].y + pp [ 0 ].y );
439  pp [ 7 ].z = 0.5 * ( pp [ 3 ].z + pp [ 0 ].z );
440 
441  pp [ 8 ].x = 0.25 * ( pp [ 0 ].x + pp [ 1 ].x + pp [ 2 ].x + pp [ 3 ].x );
442  pp [ 8 ].y = 0.25 * ( pp [ 0 ].y + pp [ 1 ].y + pp [ 2 ].y + pp [ 3 ].y );
443  pp [ 8 ].z = 0.25 * ( pp [ 0 ].z + pp [ 1 ].z + pp [ 2 ].z + pp [ 3 ].z );
444 
445  for ( GaussPoint *gp: *this->giveDefaultIntegrationRulePtr() ) {
446  const FloatArray& gpCoords = gp->giveNaturalCoordinates();
447  if ( ( gpCoords.at(1) > 0. ) && ( gpCoords.at(2) > 0. ) ) {
448  ind.at(1) = 0;
449  ind.at(2) = 4;
450  ind.at(3) = 8;
451  ind.at(4) = 7;
452  } else if ( ( gpCoords.at(1) < 0. ) && ( gpCoords.at(2) > 0. ) ) {
453  ind.at(1) = 4;
454  ind.at(2) = 1;
455  ind.at(3) = 5;
456  ind.at(4) = 8;
457  } else if ( ( gpCoords.at(1) < 0. ) && ( gpCoords.at(2) < 0. ) ) {
458  ind.at(1) = 5;
459  ind.at(2) = 2;
460  ind.at(3) = 6;
461  ind.at(4) = 8;
462  } else {
463  ind.at(1) = 6;
464  ind.at(2) = 3;
465  ind.at(3) = 7;
466  ind.at(4) = 8;
467  }
468 
469  if ( giveIPValue(v [ 0 ], gp, gc.giveIntVarType(), tStep) == 0 ) {
470  return;
471  }
472 
473  indx = gc.giveIntVarIndx();
474 
475  for ( i = 1; i <= 4; i++ ) {
476  s [ i - 1 ] = v [ 0 ].at(indx);
477  }
478 
479  for ( i = 0; i < 4; i++ ) {
480  p [ i ].x = pp [ ind.at(i + 1) ].x;
481  p [ i ].y = pp [ ind.at(i + 1) ].y;
482  p [ i ].z = pp [ ind.at(i + 1) ].z;
483  }
484 
485  gc.updateFringeTableMinMax(s, 4);
486  tr = CreateQuadWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ], s [ 3 ]);
487  EGWithMaskChangeAttributes(LAYER_MASK, tr);
488  EMAddGraphicsToModel(ESIModel(), tr);
489  }
490  }
491 }
492 
493 
494 void
496 {
497  int i;
498  WCRec l [ 2 ];
499  GraphicObj *tr;
500  double defScale = gc.getDefScale();
501  FloatArray crackStatuses, cf;
502 
503  if ( !gc.testElementGraphicActivity(this) ) {
504  return;
505  }
506 
507  if ( gc.giveIntVarType() == IST_CrackState ) {
508  // ask if any active crack exist
509  int crackStatus;
510  double ax, ay, bx, by, norm, xc, yc, length;
511  FloatArray crackDir;
512  FloatArray gpglobalcoords;
513 
514  for ( GaussPoint *gp: *this->giveDefaultIntegrationRulePtr() ) {
515 
516  if ( this->giveIPValue(cf, gp, IST_CrackedFlag, tStep) == 0 ) {
517  return;
518  }
519 
520  if ( ( int ) cf.at(1) == 0 ) {
521  return;
522  }
523 
524  if ( this->giveIPValue(crackDir, gp, IST_CrackDirs, tStep) ) {
525  this->giveIPValue(crackStatuses, gp, IST_CrackStatuses, tStep);
526  for ( i = 1; i <= 3; i++ ) {
527  crackStatus = ( int ) crackStatuses.at(i);
528  if ( ( crackStatus != pscm_NONE ) && ( crackStatus != pscm_CLOSED ) ) {
529  // draw a crack
530  // this element is 2d element in x-y plane
531  //
532  // compute perpendicular line to normal in xy plane
533  ax = crackDir.at(i);
534  ay = crackDir.at(3 + i);
535  if ( fabs(ax) > 1.e-6 ) {
536  by = 1.;
537  bx = -ay * by / ax;
538  norm = sqrt(bx * bx + by * by);
539  bx = bx / norm; // normalize to obtain unit vector
540  by = by / norm;
541  } else {
542  by = 0.0;
543  bx = 1.0;
544  }
545 
546  // obtain gp global coordinates
547  if ( gc.getInternalVarsDefGeoFlag() ) {
548  double ksi, eta, n1, n2, n3, n4;
549  ksi = gp->giveNaturalCoordinate(1);
550  eta = gp->giveNaturalCoordinate(2);
551 
552  n1 = ( 1. + ksi ) * ( 1. + eta ) * 0.25;
553  n2 = ( 1. - ksi ) * ( 1. + eta ) * 0.25;
554  n3 = ( 1. - ksi ) * ( 1. - eta ) * 0.25;
555  n4 = ( 1. + ksi ) * ( 1. - eta ) * 0.25;
556 
557  gpglobalcoords.resize(2);
558 
559  gpglobalcoords.at(1) = ( n1 * this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale) +
560  n2 * this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale) +
561  n3 * this->giveNode(3)->giveUpdatedCoordinate(1, tStep, defScale) +
562  n4 * this->giveNode(4)->giveUpdatedCoordinate(1, tStep, defScale) );
563  gpglobalcoords.at(2) = ( n1 * this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale) +
564  n2 * this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale) +
565  n3 * this->giveNode(3)->giveUpdatedCoordinate(2, tStep, defScale) +
566  n4 * this->giveNode(4)->giveUpdatedCoordinate(2, tStep, defScale) );
567  } else {
568  computeGlobalCoordinates( gpglobalcoords, gp->giveNaturalCoordinates() );
569  }
570 
571  xc = gpglobalcoords.at(1);
572  yc = gpglobalcoords.at(2);
573  length = sqrt( computeVolumeAround(gp) / this->giveCrossSection()->give(CS_Thickness, gp) ) / 2.;
574  l [ 0 ].x = ( FPNum ) xc + bx * length;
575  l [ 0 ].y = ( FPNum ) yc + by * length;
576  l [ 0 ].z = 0.;
577  l [ 1 ].x = ( FPNum ) xc - bx * length;
578  l [ 1 ].y = ( FPNum ) yc - by * length;
579  l [ 1 ].z = 0.;
580  EASValsSetLayer(OOFEG_CRACK_PATTERN_LAYER);
581  EASValsSetLineWidth(OOFEG_CRACK_PATTERN_WIDTH);
582  if ( ( crackStatus == pscm_SOFTENING ) || ( crackStatus == pscm_OPEN ) ) {
583  EASValsSetColor( gc.getActiveCrackColor() );
584  } else {
585  EASValsSetColor( gc.getCrackPatternColor() );
586  }
587 
588  tr = CreateLine3D(l);
589  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, tr);
590  EMAddGraphicsToModel(ESIModel(), tr);
591  }
592  }
593  }
594  }
595  }
596 }
597 
598 #endif
599 
600 
601 void
603 {
604  pap.resize(4);
605  for ( int i = 1; i < 5; i++ ) {
606  pap.at(i) = this->giveNode(i)->giveNumber();
607  }
608 }
609 
610 
611 void
613 {
614  int found = 0;
615  answer.resize(1);
616 
617  for ( int i = 1; i < 5; i++ ) {
618  if ( pap == this->giveNode(i)->giveNumber() ) {
619  found = 1;
620  }
621  }
622 
623  if ( found ) {
624  answer.at(1) = pap;
625  } else {
626  OOFEM_ERROR("node unknown");
627  }
628 }
629 
630 
631 int
633 {
635 }
636 
637 
640 {
641  return SPRPatchType_2dxy;
642 }
643 
644 } // end namespace oofem
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.
Class and object Domain.
Definition: domain.h:115
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Access method for default integration rule.
Definition: element.h:822
ScalarAlgorithmType getScalarAlgo()
virtual int SPRNodalRecoveryMI_giveNumberOfIP()
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
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
#define pscm_NONE
Definition: fcm.h:54
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
virtual double giveCoordinate(int i)
Definition: node.C:82
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
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)
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
virtual SPRPatchType SPRNodalRecoveryMI_givePatchType()
#define OOFEG_DEFORMED_GEOMETRY_LAYER
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
InternalStateType giveIntVarType()
virtual void computeBHmatrixAt(GaussPoint *gp, FloatMatrix &answer)
Computes a matrix which, multiplied by the column matrix of nodal displacements, gives the displaceme...
#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)
#define OOFEG_RAW_GEOMETRY_WIDTH
UnknownType
Type representing particular unknown (its physical meaning).
Definition: unknowntype.h:55
#define OOFEG_CRACK_PATTERN_WIDTH
virtual void HuertaErrorEstimatorI_computeNmatrixAt(GaussPoint *gp, FloatMatrix &answer)
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
virtual double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Returns updated ic-th coordinate of receiver.
Definition: node.C:245
Quad1PlaneStrain(int n, Domain *d)
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 SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
int numberOfGaussPoints
Number of integration points as specified by nip.
Definition: element.h:188
InternalStateMode giveIntVarMode()
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
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)
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
static FEI2dQuadLin interp
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 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)
virtual FEInterpolation * giveInterpolation() const
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
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
virtual void drawSpecial(oofegGraphicContext &gc, TimeStep *tStep)
#define OOFEG_VARPLOT_PATTERN_LAYER
Class representing integration point in finite element program.
Definition: gausspoint.h:93
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
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
virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
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:31 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011