OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
ltrspace.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/3D/ltrspace.h"
36 #include "../sm/CrossSections/structuralcrosssection.h"
37 #include "node.h"
38 #include "material.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 "fei3dtetlin.h"
46 #include "classfactory.h"
47 
48 #ifdef __OOFEG
49  #include "oofeggraphiccontext.h"
50  #include "oofegutils.h"
51  #include "Materials/rcm2.h"
52 
53  #include <Etetrawd.h>
54 #endif
55 
56 namespace oofem {
57 REGISTER_Element(LTRSpace);
58 
59 FEI3dTetLin LTRSpace :: interpolation;
60 
61 LTRSpace :: LTRSpace(int n, Domain *aDomain) :
66 
67 {
68  numberOfDofMans = 4;
70 }
71 
72 
73 Interface *
75 {
76  if ( interface == ZZNodalRecoveryModelInterfaceType ) {
77  return static_cast< ZZNodalRecoveryModelInterface * >(this);
78  } else if ( interface == NodalAveragingRecoveryModelInterfaceType ) {
79  return static_cast< NodalAveragingRecoveryModelInterface * >(this);
80  } else if ( interface == SPRNodalRecoveryModelInterfaceType ) {
81  return static_cast< SPRNodalRecoveryModelInterface * >(this);
82  } else if ( interface == SpatialLocalizerInterfaceType ) {
83  return static_cast< SpatialLocalizerInterface * >(this);
84  } else if ( interface == ZZErrorEstimatorInterfaceType ) {
85  return static_cast< ZZErrorEstimatorInterface * >(this);
86  } else if ( interface == HuertaErrorEstimatorInterfaceType ) {
87  return static_cast< HuertaErrorEstimatorInterface * >(this);
88  }
89 
90  return NULL;
91 }
92 
93 
96 {
97  return & interpolation;
98 }
99 
100 
101 
102 void
104 // Returns the lumped mass matrix of the receiver.
105 {
106  GaussPoint *gp;
107  double dV, mss1;
108 
109  answer.resize(12, 12);
110  answer.zero();
111  gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
112 
113  dV = this->computeVolumeAround(gp);
114  double density = this->giveStructuralCrossSection()->give('d', gp);
115  mss1 = dV * density / 4.;
116 
117  for ( int i = 1; i <= 12; i++ ) {
118  answer.at(i, i) = mss1;
119  }
120 }
121 
122 
123 void
125  InternalStateType type, TimeStep *tStep)
126 {
127  GaussPoint *gp;
128 
129  if ( numberOfGaussPoints != 1 ) {
130  answer.clear(); // for more gp's need to be refined
131  return;
132  }
133 
134  gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
135  giveIPValue(answer, gp, type, tStep);
136 }
137 
138 
139 void
141 {
142  pap.resize(4);
143  pap.at(1) = this->giveNode(1)->giveNumber();
144  pap.at(2) = this->giveNode(2)->giveNumber();
145  pap.at(3) = this->giveNode(3)->giveNumber();
146  pap.at(4) = this->giveNode(4)->giveNumber();
147 }
148 
149 
150 void
152 {
153  int found = 0;
154  answer.resize(1);
155 
156  for ( int i = 1; i <= 4; i++ ) {
157  if ( this->giveNode(i)->giveNumber() == pap ) {
158  found = 1;
159  }
160  }
161 
162  if ( found ) {
163  answer.at(1) = pap;
164  } else {
165  OOFEM_ERROR("node unknown");
166  }
167 }
168 
169 
170 int
173 
174 
177 {
178  return SPRPatchType_3dBiLin;
179 }
180 
181 
182 void
184  IntArray &localNodeIdArray, IntArray &globalNodeIdArray,
186  int &localNodeId, int &localElemId, int &localBcId,
187  IntArray &controlNode, IntArray &controlDof,
189 {
190  FloatArray *corner [ 4 ], midSide [ 6 ], midFace [ 4 ], midNode;
191  double x = 0.0, y = 0.0, z = 0.0;
192  int inode, nodes = 4, iside, sides = 6, iface, faces = 4, nd, nd1, nd2;
193 
194  static int sideNode [ 6 ] [ 2 ] = { { 1, 2 }, { 2, 3 }, { 3, 1 }, { 1, 4 }, { 2, 4 }, { 3, 4 } };
195  static int faceNode [ 4 ] [ 3 ] = { { 1, 2, 3 }, { 1, 2, 4 }, { 2, 3, 4 }, { 3, 1, 4 } };
196 
197  /* ordering of hexa nodes must be compatible with refinedElement connectivity ordering;
198  * generally the ordering is: corner side side face side face face center;
199  * however the concrete ordering is element dependent (see refineMeshGlobally source if in doubts) */
200 
201  int hexaSideNode [ 4 ] [ 3 ] = { { 1, 3, 4 }, { 2, 1, 5 }, { 3, 2, 6 }, { 4, 6, 5 } };
202  int hexaFaceNode [ 4 ] [ 3 ] = { { 1, 2, 4 }, { 1, 3, 2 }, { 1, 4, 3 }, { 4, 2, 3 } };
203 
206  for ( inode = 0; inode < nodes; inode++ ) {
207  corner [ inode ] = this->giveNode(inode + 1)->giveCoordinates();
208 
209  x += corner [ inode ]->at(1);
210  y += corner [ inode ]->at(2);
211  z += corner [ inode ]->at(3);
212  }
213 
214  for ( iside = 0; iside < sides; iside++ ) {
215  midSide [ iside ].resize(3);
216 
217  nd1 = sideNode [ iside ] [ 0 ] - 1;
218  nd2 = sideNode [ iside ] [ 1 ] - 1;
219 
220  midSide [ iside ].at(1) = ( corner [ nd1 ]->at(1) + corner [ nd2 ]->at(1) ) / 2.0;
221  midSide [ iside ].at(2) = ( corner [ nd1 ]->at(2) + corner [ nd2 ]->at(2) ) / 2.0;
222  midSide [ iside ].at(3) = ( corner [ nd1 ]->at(3) + corner [ nd2 ]->at(3) ) / 2.0;
223  }
224 
225  midNode.resize(3);
226 
227  midNode.at(1) = x / nodes;
228  midNode.at(2) = y / nodes;
229  midNode.at(3) = z / nodes;
230 
231  for ( iface = 0; iface < faces; iface++ ) {
232  x = y = z = 0.0;
233  for ( inode = 0; inode < 3; inode++ ) {
234  nd = faceNode [ iface ] [ inode ] - 1;
235  x += corner [ nd ]->at(1);
236  y += corner [ nd ]->at(2);
237  z += corner [ nd ]->at(3);
238  }
239 
240  midFace [ iface ].resize(3);
241 
242  midFace [ iface ].at(1) = x / 3;
243  midFace [ iface ].at(2) = y / 3;
244  midFace [ iface ].at(3) = z / 3;
245  }
246  }
247 
248  this->setupRefinedElementProblem3D(this, refinedElement, level, nodeId, localNodeIdArray, globalNodeIdArray,
249  sMode, tStep, nodes, corner, midSide, midFace, midNode,
250  localNodeId, localElemId, localBcId, hexaSideNode, hexaFaceNode,
251  controlNode, controlDof, aMode, "LSpace");
252 }
253 
255 {
257 }
258 
259 #ifdef __OOFEG
260  #define TR_LENGHT_REDUCT 0.3333
261 
263 {
264  WCRec p [ 4 ];
265  GraphicObj *go;
266 
267  if ( !gc.testElementGraphicActivity(this) ) {
268  return;
269  }
270 
271  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
272  EASValsSetColor( gc.getElementColor() );
273  EASValsSetEdgeColor( gc.getElementEdgeColor() );
274  EASValsSetEdgeFlag(true);
275  EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
276  EASValsSetFillStyle(FILL_SOLID);
277 #ifdef __PARALLEL_MODE
278  if (this->giveParallelMode() == Element_remote) {
279  EASValsSetColor( gc.getRemoteElementColor() );
280  EASValsSetEdgeColor( gc.getRemoteElementEdgeColor() );
281  }
282 #endif
283  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
284  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
285  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
286  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
287  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
288  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
289  p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveCoordinate(1);
290  p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveCoordinate(2);
291  p [ 2 ].z = ( FPNum ) this->giveNode(3)->giveCoordinate(3);
292  p [ 3 ].x = ( FPNum ) this->giveNode(4)->giveCoordinate(1);
293  p [ 3 ].y = ( FPNum ) this->giveNode(4)->giveCoordinate(2);
294  p [ 3 ].z = ( FPNum ) this->giveNode(4)->giveCoordinate(3);
295 
296  go = CreateTetra(p);
297  EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
298  EGAttachObject(go, ( EObjectP ) this);
299  EMAddGraphicsToModel(ESIModel(), go);
300 }
301 
302 
304 {
305  WCRec p [ 4 ];
306  GraphicObj *go;
307  double defScale = gc.getDefScale();
308 
309  if ( !gc.testElementGraphicActivity(this) ) {
310  return;
311  }
312 
313  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
314  EASValsSetColor( gc.getDeformedElementColor() );
315  EASValsSetEdgeColor( gc.getElementEdgeColor() );
316  EASValsSetEdgeFlag(true);
317  EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
318  EASValsSetFillStyle(FILL_SOLID);
319  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
320  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
321  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
322  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
323  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
324  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
325  p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveUpdatedCoordinate(1, tStep, defScale);
326  p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveUpdatedCoordinate(2, tStep, defScale);
327  p [ 2 ].z = ( FPNum ) this->giveNode(3)->giveUpdatedCoordinate(3, tStep, defScale);
328  p [ 3 ].x = ( FPNum ) this->giveNode(4)->giveUpdatedCoordinate(1, tStep, defScale);
329  p [ 3 ].y = ( FPNum ) this->giveNode(4)->giveUpdatedCoordinate(2, tStep, defScale);
330  p [ 3 ].z = ( FPNum ) this->giveNode(4)->giveUpdatedCoordinate(3, tStep, defScale);
331 
332  go = CreateTetra(p);
333  EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
334  EMAddGraphicsToModel(ESIModel(), go);
335 }
336 
337 
339 {
340  int i, indx, result = 0;
341  WCRec p [ 4 ];
342  GraphicObj *tr;
343  FloatArray v [ 4 ];
344  double s [ 4 ], defScale = 0.0;
345 
346  if ( !gc.testElementGraphicActivity(this) ) {
347  return;
348  }
349 
350  if ( gc.giveIntVarMode() == ISM_recovered ) {
351  for ( i = 1; i <= 4; i++ ) {
352  result += this->giveInternalStateAtNode(v [ i - 1 ], gc.giveIntVarType(), gc.giveIntVarMode(), i, tStep);
353  }
354 
355  if ( result != 4 ) {
356  return;
357  }
358  } else if ( gc.giveIntVarMode() == ISM_local ) {
359  if ( integrationRulesArray [ 0 ]->giveNumberOfIntegrationPoints() != 1 ) return;
361  if ( giveIPValue(v [ 0 ], gp, gc.giveIntVarType(), tStep) == 0 ) {
362  return;
363  }
364  v[1]=v[0];
365  v[2]=v[0];
366  v[3]=v[0];
367  }
368 
369  indx = gc.giveIntVarIndx();
370 
371  for ( i = 1; i <= 4; i++ ) {
372  s [ i - 1 ] = v [ i - 1 ].at(indx);
373  }
374 
375  EASValsSetEdgeColor( gc.getElementEdgeColor() );
376  EASValsSetEdgeFlag(true);
377  EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
378  if ( gc.getScalarAlgo() == SA_ISO_SURF ) {
379  for ( i = 0; i < 4; i++ ) {
380  if ( gc.getInternalVarsDefGeoFlag() ) {
381  // use deformed geometry
382  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(1, tStep, defScale);
383  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(2, tStep, defScale);
384  p [ i ].z = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(3, tStep, defScale);
385  } else {
386  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
387  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
388  p [ i ].z = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(3);
389  }
390  }
391 
392  gc.updateFringeTableMinMax(s, 4);
393  tr = CreateTetraWD(p, s);
394  EGWithMaskChangeAttributes(LAYER_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK, tr);
395  EMAddGraphicsToModel(ESIModel(), tr);
396  }
397 }
398 
399 void
401 {
402  int i, j, k;
403  WCRec q [ 4 ];
404  GraphicObj *tr;
405  double defScale = gc.getDefScale();
406  FloatArray crackStatuses, cf;
407 
408  if ( !gc.testElementGraphicActivity(this) ) {
409  return;
410  }
411 
412  if ( gc.giveIntVarType() == IST_CrackState ) {
413  int crackStatus;
414  double xc, yc, zc, length;
415  FloatArray crackDir;
416 
417  if ( numberOfGaussPoints != 1 ) {
418  return;
419  }
420 
421  // for (GaussPoint *gp: *integrationRulesArray [ 0 ] ) {
422  {
423  GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
424  if ( this->giveIPValue(cf, gp, IST_CrackedFlag, tStep) == 0 ) {
425  return;
426  }
427 
428  if ( ( int ) cf.at(1) == 0 ) {
429  return;
430  }
431 
432  //
433  // obtain gp global coordinates - here only one exists
434  // it is in centre of gravity.
435  xc = yc = zc = 0.;
436  for ( i = 0; i < 4; i++ ) {
437  if ( gc.getInternalVarsDefGeoFlag() ) {
438  // use deformed geometry
439  xc += ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(1, tStep, defScale);
440  yc += ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(2, tStep, defScale);
441  zc += ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(3, tStep, defScale);
442  } else {
443  xc += ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
444  yc += ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
445  zc += ( FPNum ) this->giveNode(i + 1)->giveCoordinate(3);
446  }
447  }
448 
449  xc = xc / 4.;
450  yc = yc / 4.;
451  zc = zc / 4.;
452  length = TR_LENGHT_REDUCT * pow(this->computeVolumeAround(gp), 1. / 3.) / 2.0;
453  if ( this->giveIPValue(crackDir, gp, IST_CrackDirs, tStep) ) {
454  this->giveIPValue(crackStatuses, gp, IST_CrackStatuses, tStep);
455 
456 
457  for ( i = 1; i <= 3; i++ ) {
458  crackStatus = ( int ) crackStatuses.at(i);
459  if ( ( crackStatus != pscm_NONE ) && ( crackStatus != pscm_CLOSED ) ) {
460  // draw a crack
461  // this element is 3d element
462 
463  if ( i == 1 ) {
464  j = 2;
465  k = 3;
466  } else if ( i == 2 ) {
467  j = 3;
468  k = 1;
469  } else {
470  j = 1;
471  k = 2;
472  }
473 
474  q [ 0 ].x = ( FPNum ) xc + 0.5 * crackDir.at(0 + j) * length + 0.5 * crackDir.at(0 + k) * length;
475  q [ 0 ].y = ( FPNum ) yc + 0.5 * crackDir.at(3 + j) * length + 0.5 * crackDir.at(3 + k) * length;
476  q [ 0 ].z = ( FPNum ) zc + 0.5 * crackDir.at(6 + j) * length + 0.5 * crackDir.at(6 + k) * length;
477  q [ 1 ].x = ( FPNum ) xc + 0.5 * crackDir.at(0 + j) * length - 0.5 * crackDir.at(0 + k) * length;
478  q [ 1 ].y = ( FPNum ) yc + 0.5 * crackDir.at(3 + j) * length - 0.5 * crackDir.at(3 + k) * length;
479  q [ 1 ].z = ( FPNum ) zc + 0.5 * crackDir.at(6 + j) * length - 0.5 * crackDir.at(6 + k) * length;
480  q [ 2 ].x = ( FPNum ) xc - 0.5 * crackDir.at(0 + j) * length - 0.5 * crackDir.at(0 + k) * length;
481  q [ 2 ].y = ( FPNum ) yc - 0.5 * crackDir.at(3 + j) * length - 0.5 * crackDir.at(3 + k) * length;
482  q [ 2 ].z = ( FPNum ) zc - 0.5 * crackDir.at(6 + j) * length - 0.5 * crackDir.at(6 + k) * length;
483  q [ 3 ].x = ( FPNum ) xc - 0.5 * crackDir.at(0 + j) * length + 0.5 * crackDir.at(0 + k) * length;
484  q [ 3 ].y = ( FPNum ) yc - 0.5 * crackDir.at(3 + j) * length + 0.5 * crackDir.at(3 + k) * length;
485  q [ 3 ].z = ( FPNum ) zc - 0.5 * crackDir.at(6 + j) * length + 0.5 * crackDir.at(6 + k) * length;
486 
487  EASValsSetLayer(OOFEG_CRACK_PATTERN_LAYER);
488  EASValsSetLineWidth(OOFEG_CRACK_PATTERN_WIDTH);
489  if ( ( crackStatus == pscm_SOFTENING ) || ( crackStatus == pscm_OPEN ) ) {
490  EASValsSetColor( gc.getActiveCrackColor() );
491  } else {
492  EASValsSetColor( gc.getCrackPatternColor() );
493  }
494 
495  // EASValsSetFillStyle (FILL_HOLLOW);
496  tr = CreateQuad3D(q);
497  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, tr);
498  EMAddGraphicsToModel(ESIModel(), tr);
499  }
500  }
501  }
502  }
503  }
504 }
505 
506 #endif
507 
508 
511 {
512  IntegrationRule *iRule = new GaussIntegrationRule(1, this, 1, 1);
513  int npoints = iRule->getRequiredNumberOfIntegrationPoints(_Triangle, approxOrder);
514  iRule->SetUpPointsOnTriangle(npoints, _Unknown);
515  return iRule;
516 }
517 
518 
519 
520 } // end namespace oofem
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Definition: ltrspace.C:103
virtual IntegrationRule * GetSurfaceIntegrationRule(int)
Definition: ltrspace.C:510
int testElementGraphicActivity(Element *)
Test if particular element passed fulfills various filtering criteria for its graphics output...
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
The element interface required by NodalAvergagingRecoveryModel.
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
virtual void drawSpecial(oofegGraphicContext &gc, TimeStep *tStep)
Definition: ltrspace.C:400
ScalarAlgorithmType getScalarAlgo()
The element interface required by ZZNodalRecoveryModel.
const FloatArray & giveSubPatchCoordinates()
Returns local sub-patch coordinates of the receiver.
Definition: gausspoint.h:144
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
Definition: ltrspace.C:124
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
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
void setupRefinedElementProblem3D(Element *element, RefinedElement *refinedElement, int level, int nodeId, IntArray &localNodeIdArray, IntArray &globalNodeIdArray, HuertaErrorEstimatorInterface::SetupMode mode, TimeStep *tStep, int nodes, FloatArray **corner, FloatArray *midSide, FloatArray *midFace, FloatArray &midNode, int &localNodeId, int &localElemId, int &localBcId, int hexaSideNode[1][3], int hexaFaceNode[1][3], IntArray &controlNode, IntArray &controlDof, HuertaErrorEstimator::AnalysisMode aMode, const char *hexatype)
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: ltrspace.C:262
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 void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
Definition: ltrspace.C:151
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: ltrspace.C:183
#define OOFEG_DEFORMED_GEOMETRY_LAYER
LTRSpace(int n, Domain *d)
Definition: ltrspace.C:61
Abstract base class representing integration rule.
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
InternalStateType giveIntVarType()
#define OOFEG_CRACK_PATTERN_LAYER
The element interface corresponding to ZZErrorEstimator.
#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
virtual int SPRNodalRecoveryMI_giveNumberOfIP()
Definition: ltrspace.C:171
UnknownType
Type representing particular unknown (its physical meaning).
Definition: unknowntype.h:55
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
#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
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 Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: ltrspace.C:74
virtual SPRPatchType SPRNodalRecoveryMI_givePatchType()
Definition: ltrspace.C:176
int numberOfGaussPoints
Number of integration points as specified by nip.
Definition: element.h:188
InternalStateMode giveIntVarMode()
#define TR_LENGHT_REDUCT
Definition: ltrspace.C:260
Class representing vector of real numbers.
Definition: floatarray.h:82
elementParallelMode giveParallelMode() const
Return elementParallelMode of receiver.
Definition: element.h:1069
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.
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1358
virtual int getRequiredNumberOfIntegrationPoints(integrationDomain dType, int approxOrder)
Abstract service.
Class Interface.
Definition: interface.h:82
int giveNumberOfIntegrationPoints() const
Returns number of integration points of receiver.
virtual FEInterpolation * giveInterpolation() const
Definition: ltrspace.C:95
Base class 3D elements.
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
Definition: ltrspace.C:303
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 int SetUpPointsOnTriangle(int, MaterialMode mode)
Sets up receiver&#39;s integration points on triangular (area coords) integration domain.
Element in active domain is only mirror of some remote element.
Definition: element.h:102
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
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
void updateFringeTableMinMax(double *s, int size)
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
Definition: ltrspace.C:140
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Definition: ltrspace.C:338
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
static FEI3dTetLin interpolation
Definition: ltrspace.h:64
#define pscm_SOFTENING
Definition: fcm.h:56
#define OOFEG_VARPLOT_PATTERN_LAYER
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
virtual void HuertaErrorEstimatorI_computeNmatrixAt(GaussPoint *gp, FloatMatrix &answer)
Definition: ltrspace.C:254
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:29 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011