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