OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
tr_shell01.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/Shells/tr_shell01.h"
36 #include "fei2dtrlin.h"
37 #include "contextioerr.h"
38 #include "gaussintegrationrule.h"
39 #include "gausspoint.h"
40 #include "classfactory.h"
41 #include "node.h"
42 
43 #ifdef __OOFEG
44  #include "node.h"
45  #include "oofeggraphiccontext.h"
46  #include "connectivitytable.h"
47 #endif
48 
49 
50 namespace oofem {
51 REGISTER_Element(TR_SHELL01);
52 
53 IntArray TR_SHELL01 :: loc_plate = {3, 4, 5, 9, 10, 11, 15, 16, 17};
54 IntArray TR_SHELL01 :: loc_membrane = {1, 2, 6, 7, 8, 12, 13, 14, 18};
55 
57  plate(new CCTPlate3d(-1, aDomain)), membrane(new TrPlaneStrRot3d(-1, aDomain))
58 {
59  numberOfDofMans = 3;
60 }
61 
62 
65 {
66  // proc tady neni return = this... ??? termitovo
68  if ( result != IRRT_OK ) {
69  return result;
70  }
71 
72 #if 0
74  if ( val != -1 ) {
75  OOFEM_WARNING("key word NIP is not allowed for element TR_SHELL01");
76  return result;
77  }
79  if ( val != -1 ) {
80  OOFEM_WARNING("key word NIProt is not allowed for element TR_SHELL01");
81  return result;
82  }
83 #endif
84 
85  result = plate->initializeFrom(ir);
86  if ( result != IRRT_OK ) {
87  return result;
88  }
89  result = membrane->initializeFrom(ir);
90  if ( result != IRRT_OK ) {
91  return result;
92  }
93 
94  return IRRT_OK;
95 }
96 
97 void
99 {
101 
102  if ( plate->giveDefaultIntegrationRulePtr()->giveNumberOfIntegrationPoints() != membrane->giveDefaultIntegrationRulePtr()->giveNumberOfIntegrationPoints() ) {
103  OOFEM_ERROR("incompatible integration rules detected %d vs %d IPs",
104  plate->giveDefaultIntegrationRulePtr()->giveNumberOfIntegrationPoints(),
105  membrane->giveDefaultIntegrationRulePtr()->giveNumberOfIntegrationPoints());
106  }
107 }
108 
109 void
111 {
113  plate->updateLocalNumbering(f);
114  membrane->updateLocalNumbering(f);
115 }
116 
118 {
120  plate->setCrossSection(csIndx);
121  membrane->setCrossSection(csIndx);
122 }
123 
124 
125 void
127 //
128 // returns characteristics vector of receiver accordind to mtrx
129 //
130 {
131  FloatArray aux;
132 
133  answer.resize(18);
134  answer.zero();
135 
136  plate->giveCharacteristicVector(aux, mtrx, mode, tStep);
137  if ( !aux.isEmpty() ) answer.assemble(aux, loc_plate);
138 
139  membrane->giveCharacteristicVector(aux, mtrx, mode, tStep);
140  if ( !aux.isEmpty() ) answer.assemble(aux, loc_membrane);
141 }
142 
143 void
145 //
146 // returns characteristics matrix of receiver accordind to mtrx
147 //
148 {
149  FloatMatrix aux;
150 
151  answer.resize(18, 18);
152  answer.zero();
153 
154  plate->giveCharacteristicMatrix(aux, mtrx, tStep);
155  if ( aux.isNotEmpty() ) answer.assemble(aux, loc_plate);
156 
157  membrane->giveCharacteristicMatrix(aux, mtrx, tStep);
158  if ( aux.isNotEmpty() ) answer.assemble(aux, loc_membrane);
159 }
160 
162 {
163  FloatArray aux;
164 
165  answer.resize(18);
166  answer.zero();
167 
168  plate->computeBodyLoadVectorAt(aux, forLoad, tStep, mode);
169  if ( !aux.isEmpty() ) answer.assemble(aux, loc_plate);
170 
171  membrane->computeBodyLoadVectorAt(aux, forLoad, tStep, mode);
172  if ( !aux.isEmpty() ) answer.assemble(aux, loc_membrane);
173 
174 
175 
176 }
177 
178 bool
180 {
181  FloatMatrix aux1, aux2;
182  int ncol;
183 
184  bool t1 = plate->giveRotationMatrix(aux1);
185  bool t2 = membrane->giveRotationMatrix(aux2);
186 
187  if ( t1 != t2 ) {
188  OOFEM_ERROR("Transformation demand mismatch");
189  }
190 
191  if ( t1 ) {
192  ncol = aux1.giveNumberOfColumns();
193  answer.resize(18, ncol);
194 
195  for ( int i = 1; i <= 9; i++ ) { // row index
196  for ( int j = 1; j <= ncol; j++ ) {
197  answer.at(loc_plate.at(i), j) = aux1.at(i, j);
198  }
199  }
200 
201  for ( int i = 1; i <= 9; i++ ) { // row index
202  for ( int j = 1; j <= ncol; j++ ) {
203  answer.at(loc_membrane.at(i), j) = aux2.at(i, j);
204  }
205  }
206  }
207 
208  return t1;
209 }
210 
211 void
213 // Updates the receiver at end of step.
214 {
215  plate->updateInternalState(tStep);
216  membrane->updateInternalState(tStep);
217 }
218 
219 void
221 {
223 
224  plate->updateYourself(tStep);
225  membrane->updateYourself(tStep);
226 }
227 
228 
229 Interface *
231 {
232  if ( interface == ZZNodalRecoveryModelInterfaceType ) {
233  return static_cast< ZZNodalRecoveryModelInterface * >(this);
234  } else if ( interface == NodalAveragingRecoveryModelInterfaceType ) {
235  return static_cast< NodalAveragingRecoveryModelInterface * >(this);
236  } else if ( interface == ZZErrorEstimatorInterfaceType ) {
237  return static_cast< ZZErrorEstimatorInterface * >(this);
238  } else if ( interface == SpatialLocalizerInterfaceType ) {
239  return static_cast< SpatialLocalizerInterface * >(this);
240  }
241 
242  return NULL;
243 }
244 
245 double
247 {
248  return plate->computeVolumeAround(gp);
249 }
250 
251 int
253 {
254  if ( type == IST_ShellForceTensor || type == IST_ShellStrainTensor ||
255  type == IST_ShellMomentTensor || type == IST_CurvatureTensor ) {
256  FloatArray aux;
257  GaussPoint *membraneGP = membrane->giveDefaultIntegrationRulePtr()->getIntegrationPoint(gp->giveNumber() - 1);
258  GaussPoint *plateGP = plate->giveDefaultIntegrationRulePtr()->getIntegrationPoint(gp->giveNumber() - 1);
259 
260  plate->giveIPValue(answer, plateGP, type, tStep);
261  membrane->giveIPValue(aux, membraneGP, type, tStep);
262  answer.add(aux);
263  return 1;
264  } else {
265  return StructuralElement :: giveIPValue(answer, gp, type, tStep);
266  }
267 }
268 
269 
270 void
272  InternalStateType type, TimeStep *tStep)
273 {
274  this->giveIPValue(answer, this->giveDefaultIntegrationRulePtr()->getIntegrationPoint(0), type, tStep);
275 }
276 
277 
278 void
280 {
281  FloatArray v;
283  fprintf( file, "element %d (%8d) :\n", this->giveLabel(), this->giveNumber() );
284 
285  for ( auto &gp: *iRule ) {
286  fprintf( file, " GP %2d.%-2d :", iRule->giveNumber(), gp->giveNumber() );
287  // Strain - Curvature
288  this->giveIPValue(v, gp, IST_ShellStrainTensor, tStep);
289 
290  fprintf(file, " strains ");
291  // eps_x, eps_y, eps_z, eps_yz, eps_xz, eps_xy
292  fprintf(file, " %.4e %.4e %.4e %.4e %.4e %.4e ",
293  v.at(1), v.at(2), v.at(3), v.at(4), v.at(5), v.at(6) );
294 
295  this->giveIPValue(v, gp, IST_CurvatureTensor, tStep);
296 
297  fprintf(file, "\n curvatures ");
298  // k_x, k_y, k_z, k_yz, k_xz, k_xy
299  fprintf(file, " %.4e %.4e %.4e %.4e %.4e %.4e ",
300  v.at(1), v.at(2), v.at(3), v.at(4), v.at(5), v.at(6) );
301 
302  // Forces - Moments
303  this->giveIPValue(v, gp, IST_ShellForceTensor, tStep);
304 
305  fprintf(file, "\n stresses ");
306  // n_x, n_y, n_z, v_yz, v_xz, v_xy
307  fprintf(file, " %.4e %.4e %.4e %.4e %.4e %.4e ",
308  v.at(1), v.at(2), v.at(3), v.at(4), v.at(5), v.at(6) );
309 
310  this->giveIPValue(v, gp, IST_ShellMomentTensor, tStep);
311 
312  fprintf(file, "\n moments ");
313  // m_x, m_y, m_z, m_yz, m_xz, m_xy
314  fprintf(file, " %.4e %.4e %.4e %.4e %.4e %.4e ",
315  v.at(1), v.at(2), v.at(3), v.at(4), v.at(5), v.at(6) );
316 
317  fprintf(file, "\n");
318  }
319 }
320 
321 
324 {
325  contextIOResultType iores;
326  if ( ( iores = StructuralElement :: saveContext(stream, mode, obj) ) != CIO_OK ) {
327  THROW_CIOERR(iores);
328  }
329  if ( ( iores = this->plate->saveContext(stream, mode, obj) ) != CIO_OK ) {
330  THROW_CIOERR(iores);
331  }
332  if ( ( iores = this->membrane->saveContext(stream, mode, obj) ) != CIO_OK ) {
333  THROW_CIOERR(iores);
334  }
335  return iores;
336 }
337 
340 {
341  contextIOResultType iores;
342  if ( ( iores = StructuralElement :: restoreContext(stream, mode, obj) ) != CIO_OK ) {
343  THROW_CIOERR(iores);
344  }
345  if ( ( iores = this->plate->restoreContext(stream, mode, obj) ) != CIO_OK ) {
346  THROW_CIOERR(iores);
347  }
348  if ( ( iores = this->membrane->restoreContext(stream, mode, obj) ) != CIO_OK ) {
349  THROW_CIOERR(iores);
350  }
351  return iores;
352 }
353 
356 {
357  if ( !this->compositeIR ) {
358  this->compositeIR.reset( new GaussIntegrationRule(1, this, 1, 12) );
359  this->compositeIR->SetUpPointsOnTriangle(plate->giveDefaultIntegrationRulePtr()->giveNumberOfIntegrationPoints(), _3dShell);
360  }
361  return this->compositeIR.get();
362 }
363 
364 void
366 {
367  // sig is global ShellForceMomentTensor
368  FloatMatrix globTensor(3, 3);
369  const FloatMatrix *GtoLRotationMatrix = plate->computeGtoLRotationMatrix();
370  FloatMatrix LtoGRotationMatrix;
371 
372  answer.resize(8); // reduced, local form
373  LtoGRotationMatrix.beTranspositionOf(* GtoLRotationMatrix);
374 
375  // Forces
376  globTensor.at(1, 1) = sig.at(1); //sxForce
377  globTensor.at(1, 2) = sig.at(6); //qxyForce
378  globTensor.at(1, 3) = sig.at(5); //qxzForce
379 
380  globTensor.at(2, 1) = sig.at(6); //qxyForce
381  globTensor.at(2, 2) = sig.at(2); //syForce
382  globTensor.at(2, 3) = sig.at(4); //syzForce
383 
384  globTensor.at(3, 1) = sig.at(5); //qxzForce
385  globTensor.at(3, 2) = sig.at(4); //syzForce
386  globTensor.at(3, 3) = sig.at(3); //szForce
387 
388  globTensor.rotatedWith(LtoGRotationMatrix);
389  // Forces: now globTensoris transformed into local c.s
390 
391  // answer should be in reduced, local form
392  answer.at(1) = globTensor.at(1, 1); //sxForce
393  answer.at(2) = globTensor.at(2, 2); //syForce
394  answer.at(3) = globTensor.at(1, 2); //qxyForce
395  answer.at(7) = globTensor.at(2, 3); //syzForce
396  answer.at(8) = globTensor.at(1, 3); //qxzForce
397 
398 
399  // Moments:
400  globTensor.at(1, 1) = sig.at(7); //mxForce
401  globTensor.at(1, 2) = sig.at(12); //mxyForce
402  globTensor.at(1, 3) = sig.at(11); //mxzForce
403 
404  globTensor.at(2, 1) = sig.at(12); //mxyForce
405  globTensor.at(2, 2) = sig.at(8); //myForce
406  globTensor.at(2, 3) = sig.at(10); //myzForce
407 
408  globTensor.at(3, 1) = sig.at(11); //mxzForce
409  globTensor.at(3, 2) = sig.at(10); //myzForce
410  globTensor.at(3, 3) = sig.at(9); //mzForce
411 
412  globTensor.rotatedWith(LtoGRotationMatrix);
413  // now globTensoris transformed into local c.s
414 
415  answer.at(4) = globTensor.at(1, 1); //mxForce
416  answer.at(5) = globTensor.at(2, 2); //myForce
417  answer.at(6) = globTensor.at(1, 2); //mxyForce
418 }
419 
420 
421 void
423 {
424  FloatArray lt3, gt3; // global vector in the element thickness direction of lenght thickeness/2
425  const FloatMatrix *GtoLRotationMatrix = plate->computeGtoLRotationMatrix();
426 
427  // setup vector in the element local cs. perpendicular to element plane of thickness/2 length
428  lt3 = {0., 0., 1.}; //this->giveCrossSection()->give(CS_Thickness)/2.0; // HUHU
429  // transform it to globa cs
430  gt3.beTProductOf(* GtoLRotationMatrix, lt3);
431 
432  // use gt3 to construct element bounding box respecting true element volume
433 
434  FloatArray _c;
435 
436  for ( int i = 1; i <= this->giveNumberOfNodes(); ++i ) {
437  FloatArray *coordinates = this->giveNode(i)->giveCoordinates();
438 
439  _c = * coordinates;
440  _c.add(gt3);
441  if ( i == 1 ) {
442  bb0 = bb1 = _c;
443  } else {
444  bb0.beMinOf(bb0, _c);
445  bb1.beMaxOf(bb1, _c);
446  }
447 
448  _c = * coordinates;
449  _c.subtract(gt3);
450  bb0.beMinOf(bb0, _c);
451  bb1.beMaxOf(bb1, _c);
452  }
453 }
454 
455 //
456 // io routines
457 //
458 #ifdef __OOFEG
459 
460 void
462 {
463  WCRec p [ 3 ];
464  GraphicObj *go;
465 
466  if ( !gc.testElementGraphicActivity(this) ) {
467  return;
468  }
469 
470  if ( this->giveMaterial()->isActivated(tStep) ) {
471  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
472  EASValsSetColor( gc.getElementColor() );
473  EASValsSetEdgeColor( gc.getElementEdgeColor() );
474  EASValsSetEdgeFlag(true);
475  EASValsSetFillStyle(FILL_SOLID);
476  EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
477  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
478  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
479  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
480  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
481  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
482  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
483  p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveCoordinate(1);
484  p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveCoordinate(2);
485  p [ 2 ].z = ( FPNum ) this->giveNode(3)->giveCoordinate(3);
486 
487  go = CreateTriangle3D(p);
488  EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
489  EGAttachObject(go, ( EObjectP ) this);
490  EMAddGraphicsToModel(ESIModel(), go);
491  }
492 }
493 
494 void
496 {
497  WCRec p [ 3 ];
498  GraphicObj *go;
499  double defScale = gc.getDefScale();
500 
501  if ( !gc.testElementGraphicActivity(this) ) {
502  return;
503  }
504 
505  if ( this->giveMaterial()->isActivated(tStep) ) {
506  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
507  EASValsSetColor( gc.getDeformedElementColor() );
508  EASValsSetEdgeColor( gc.getElementEdgeColor() );
509  EASValsSetEdgeFlag(true);
510  EASValsSetFillStyle(FILL_SOLID);
511  EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
512  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
513  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
514  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
515  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
516  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
517  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
518  p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveUpdatedCoordinate(1, tStep, defScale);
519  p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveUpdatedCoordinate(2, tStep, defScale);
520  p [ 2 ].z = ( FPNum ) this->giveNode(3)->giveUpdatedCoordinate(3, tStep, defScale);
521 
522  go = CreateTriangle3D(p);
523  EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
524  EMAddGraphicsToModel(ESIModel(), go);
525  }
526 }
527 
528 void
530 {
531  int indx, result = 0;
532  WCRec p [ 3 ];
533  GraphicObj *tr;
534  FloatArray v1, v2, v3;
535  double s [ 3 ], defScale;
536  if ( !gc.testElementGraphicActivity(this) ) {
537  return;
538  }
539 
540  if ( !this->giveMaterial()->isActivated(tStep) ) {
541  return;
542  }
543 
544  if ( gc.giveIntVarMode() == ISM_recovered ) {
545  result += this->giveInternalStateAtNode(v1, gc.giveIntVarType(), gc.giveIntVarMode(), 1, tStep);
546  result += this->giveInternalStateAtNode(v2, gc.giveIntVarType(), gc.giveIntVarMode(), 2, tStep);
547  result += this->giveInternalStateAtNode(v3, gc.giveIntVarType(), gc.giveIntVarMode(), 3, tStep);
548  } else if ( gc.giveIntVarMode() == ISM_local ) {
549  double tot_w = 0.;
550  FloatArray a, v;
551  for ( GaussPoint *gp: *plate->giveDefaultIntegrationRulePtr() ) {
552  this->giveIPValue(a, gp, IST_ShellMomentTensor, tStep);
553  v.add(gp->giveWeight(), a);
554  tot_w += gp->giveWeight();
555  }
556  v.times(1. / tot_w);
557  v1 = v;
558  v2 = v;
559  v3 = v;
560  }
561 
562  indx = gc.giveIntVarIndx();
563 
564  s [ 0 ] = v1.at(indx);
565  s [ 1 ] = v2.at(indx);
566  s [ 2 ] = v3.at(indx);
567  EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
568  if ( gc.getScalarAlgo() == SA_ISO_SURF ) {
569  for ( int i = 0; i < 3; i++ ) {
570  if ( gc.getInternalVarsDefGeoFlag() ) {
571  // use deformed geometry
572  defScale = gc.getDefScale();
573  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(1, tStep, defScale);
574  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(2, tStep, defScale);
575  p [ i ].z = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(3, tStep, defScale);
576  } else {
577  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
578  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
579  p [ i ].z = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(3);
580  }
581  }
582  // //EASValsSetColor(gc.getYieldPlotColor(ratio));
583  gc.updateFringeTableMinMax(s, 3);
584  tr = CreateTriangleWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ]);
585  EGWithMaskChangeAttributes(LAYER_MASK, tr);
586  EMAddGraphicsToModel(ESIModel(), tr);
587  }
588 }
589 
590 #endif
591 } // end namespace oofem
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.
int giveNumberOfColumns() const
Returns number of columns of receiver.
Definition: floatmatrix.h:158
void subtract(const FloatArray &src)
Subtracts array src to receiver.
Definition: floatarray.C:258
virtual void postInitialize()
Performs post initialization steps.
Definition: element.C:752
virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType mtrx, TimeStep *tStep)
Computes characteristic matrix of receiver of requested type in given time step.
Definition: tr_shell01.C:144
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
ScalarAlgorithmType getScalarAlgo()
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: tr_shell01.C:339
The element interface required by ZZNodalRecoveryModel.
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
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...
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: tr_shell01.C:220
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
#define OOFEG_RAW_GEOMETRY_LAYER
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual bool giveRotationMatrix(FloatMatrix &answer)
Transformation matrices updates rotation matrix between element-local and primary DOFs...
Definition: tr_shell01.C:179
This class represent triangular plane stress element with rotational degree of freedom around normal ...
Definition: trplanrot3d.h:68
virtual double giveCoordinate(int i)
Definition: node.C:82
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: tr_shell01.C:252
#define THROW_CIOERR(e)
Definition: contextioerr.h:61
virtual void SpatialLocalizerI_giveBBox(FloatArray &bb0, FloatArray &bb1)
Creates a bounding box of the nodes and checks if it includes the given coordinate.
Definition: tr_shell01.C:422
#define OOFEG_DEFORMED_GEOMETRY_LAYER
void beMaxOf(const FloatArray &a, const FloatArray &b)
Sets receiver to maximum of a or b&#39;s respective elements.
Definition: floatarray.C:288
static IntArray loc_plate
Definition: tr_shell01.h:71
Abstract base class representing integration rule.
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: element.C:970
virtual int giveNumberOfNodes() const
Returns number of nodes of receiver.
Definition: element.h:662
virtual void updateInternalState(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: tr_shell01.C:212
bool isNotEmpty() const
Tests for empty matrix.
Definition: floatmatrix.h:162
InternalStateType giveIntVarType()
Abstract base class for all "structural" finite elements.
std::unique_ptr< IntegrationRule > compositeIR
Element integraton rule (plate and membrane parts have their own integration rules) this one used to ...
Definition: tr_shell01.h:69
void beMinOf(const FloatArray &a, const FloatArray &b)
Sets receiver to be minimum of a or b&#39;s respective elements.
Definition: floatarray.C:315
This class represent CCT plate element that can be arbitrary oriented in space, in contrast to base C...
Definition: cct3d.h:72
The element interface corresponding to ZZErrorEstimator.
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
static IntArray loc_membrane
Definition: tr_shell01.h:72
int giveNumber()
Returns number of receiver.
Definition: gausspoint.h:184
#define OOFEG_RAW_GEOMETRY_WIDTH
UnknownType
Type representing particular unknown (its physical meaning).
Definition: unknowntype.h:55
virtual void ZZErrorEstimatorI_computeLocalStress(FloatArray &answer, FloatArray &sig)
Returns stress vector in global c.s.
Definition: tr_shell01.C:365
void rotatedWith(const FloatMatrix &r, char mode= 'n')
Returns the receiver &#39;a&#39; transformed using give transformation matrix r.
Definition: floatmatrix.C:1557
void updateLocalNumbering(EntityRenumberingFunctor &f)
Local renumbering support.
Definition: tr_shell01.C:110
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: tr_shell01.C:64
virtual void setCrossSection(int csIndx)
Sets the cross section model of receiver.
Definition: element.h:653
virtual double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Returns updated ic-th coordinate of receiver.
Definition: node.C:245
TR_SHELL01(int n, Domain *d)
Constructor.
Definition: tr_shell01.C:56
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Definition: tr_shell01.C:529
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix^T and anArray.
Definition: floatarray.C:708
bool isEmpty() const
Returns true if receiver is empty.
Definition: floatarray.h:222
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
virtual IntegrationRule * ZZErrorEstimatorI_giveIntegrationRule()
Returns element integration rule used to evaluate error.
Definition: tr_shell01.C:355
InternalStateMode giveIntVarMode()
virtual void giveCharacteristicVector(FloatArray &answer, CharType mtrx, ValueModeType mode, TimeStep *tStep)
Computes characteristic vector of receiver of requested type in given time step.
Definition: tr_shell01.C:126
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: tr_shell01.C:461
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 void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
void assemble(const FloatArray &fe, const IntArray &loc)
Assembles the array fe (typically, the load vector of a finite element) into the receiver, using loc as location array.
Definition: floatarray.C:551
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
Class Interface.
Definition: interface.h:82
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type)
Definition: tr_shell01.C:495
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: element.C:885
std::unique_ptr< CCTPlate3d > plate
Pointer to plate element.
Definition: tr_shell01.h:61
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: tr_shell01.C:230
virtual void postInitialize()
Performs post initialization steps.
Definition: tr_shell01.C:98
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
void times(double s)
Multiplies receiver with scalar.
Definition: floatarray.C:818
virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode)
Computes the load vector due to body load acting on receiver, at given time step. ...
Definition: tr_shell01.C:161
void beTranspositionOf(const FloatMatrix &src)
Assigns to the receiver the transposition of parameter.
Definition: floatmatrix.C:323
The spatial localizer element interface associated to spatial localizer.
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
Definition: tr_shell01.C:271
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
#define new
virtual bool isActivated(TimeStep *tStep)
Definition: material.h:161
virtual FloatArray * giveCoordinates()
Definition: node.h:114
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: tr_shell01.C:323
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
void updateFringeTableMinMax(double *s, int size)
#define _IFT_TrPlaneStrRot_niprot
Definition: trplanrot.h:43
Load is base abstract class for all loads.
Definition: load.h:61
#define IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:78
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Access method for default integration rule.
Definition: tr_shell01.h:110
std::unique_ptr< TrPlaneStrRot3d > membrane
Pointer to membrane (plane stress) element.
Definition: tr_shell01.h:63
int giveLabel() const
Definition: element.h:1055
the oofem namespace is to define a context or scope in which all oofem names are defined.
void assemble(const FloatMatrix &src, const IntArray &loc)
Assembles the contribution using localization array into receiver.
Definition: floatmatrix.C:251
virtual void updateLocalNumbering(EntityRenumberingFunctor &f)
Local renumbering support.
Definition: element.C:1511
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
#define OOFEG_VARPLOT_PATTERN_LAYER
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
#define OOFEM_WARNING(...)
Definition: error.h:62
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
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Definition: tr_shell01.C:279
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: tr_shell01.C:246
void setCrossSection(int csIndx)
Sets the cross section model of receiver.
Definition: tr_shell01.C:117
#define _IFT_Element_nip
Definition: element.h:72
virtual Material * giveMaterial()
Definition: element.C:484
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:32 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011