OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
tr1_2d_pfem.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 "tr1_2d_pfem.h"
36 #include "node.h"
37 #include "material.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 "domain.h"
45 #include "mathfem.h"
46 #include "engngm.h"
47 #include "pfem.h"
48 #include "fluiddynamicmaterial.h"
49 #include "fluidcrosssection.h"
50 #include "load.h"
51 #include "bodyload.h"
52 #include "timestep.h"
53 #include "boundaryload.h"
54 #include "mathfem.h"
55 #include "contextioerr.h"
56 
57 #ifdef __OOFEG
58  #include "oofeggraphiccontext.h"
59  #include "connectivitytable.h"
60 #endif
61 
62 namespace oofem {
65 
66 IntArray TR1_2D_PFEM :: edge_ordering [ 3 ] = {
67  IntArray(6), IntArray(6), IntArray(6)
68 };
69 
71  1, 2, 4, 5, 7, 8
72 };
74  3, 6, 9
75 };
76 
77 
78 TR1_2D_PFEM :: TR1_2D_PFEM(int n, Domain *aDomain, int particle1, int particle2, int particle3, int mat, int cs) :
79  PFEMElement2d(n, aDomain)
80 {
81  // Constructor.
82  numberOfDofMans = 3;
83  IntArray aBodyLoadArry(1);
84  aBodyLoadArry.at(1) = 3;
85  this->setDofManagers({ particle1, particle2, particle3 });
86  // CHECK THIS - NOT NICE
87  this->setBodyLoads(aBodyLoadArry);
88  this->material = mat;
89  this->setCrossSection(cs);
90  this->postInitialize();
91 }
92 
94 // Destructor
95 { }
96 
97 int
99 {
100  return 9;
101 }
102 
103 void
105 {
106  answer = {
107  V_u, V_v, P_f
108  };
109 }
110 
111 // NOT IN USE
112 void
114 {
115  this->giveDofManDofIDMask(1, answer);
116 }
117 
118 
121 {
123 
124  this->computeGaussPoints();
125  return ret;
126 }
127 
128 void
130 // Sets up the array containing the four Gauss points of the receiver.
131 {
132  if ( integrationRulesArray.size() == 0 ) {
133  integrationRulesArray.resize(1);
134  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 3) );
135  integrationRulesArray [ 0 ]->setUpIntegrationPoints(_Triangle, 3, _2dFlow);
136  }
137 }
138 
139 //copied from tr21stokes
140 void
142 {
143  IntegrationRule *iRule = this->integrationRulesArray [ 0 ].get();
144  FloatArray N, gVector;
145 
146  answer.resize(6);
147  answer.zero();
148 
149  load->computeComponentArrayAt(gVector, atTime, mode);
150  if ( gVector.giveSize() ) {
151  for ( auto &gp : *iRule ) {
152  double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity( integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
153  double detJ = this->pressureInterpolation.giveTransformationJacobian( gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
154  double dA = detJ * gp->giveWeight();
155  this->pressureInterpolation.evalN( N, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
156  for ( int j = 0; j < 3; j++ ) {
157  answer(2 * j) += N(j) * rho * gVector(0) * dA;
158  answer(2 * j + 1) += N(j) * rho * gVector(1) * dA;
159  }
160  }
161  }
162 }
163 
164 
165 // NOT IN USE
166 //copied from tr21stokes
167 void
169 {
170  answer.resize(15);
171  answer.zero();
172 
173  if ( load->giveType() == TransmissionBC ) { // Neumann boundary conditions (traction)
174  BoundaryLoad *boundaryLoad = ( BoundaryLoad * ) load;
175 
176  int numberOfEdgeIPs = ( int ) ceil( ( boundaryLoad->giveApproxOrder() + 1. ) / 2. ) * 2;
177 
178  GaussIntegrationRule iRule(1, this, 1, 1);
179  FloatArray N, t, f(6), f2(6);
180  IntArray edge_mapping;
181 
182  f.zero();
183  iRule.setUpIntegrationPoints(_Line, numberOfEdgeIPs, _Unknown);
184 
185  for ( auto &gp : iRule ) {
186  const FloatArray &lcoords = gp->giveNaturalCoordinates();
187  this->pressureInterpolation.edgeEvalN( N, iEdge, lcoords, FEIElementGeometryWrapper(this) );
188  double dS = gp->giveWeight() * this->pressureInterpolation.edgeGiveTransformationJacobian( iEdge, lcoords, FEIElementGeometryWrapper(this) );
189 
190  if ( boundaryLoad->giveFormulationType() == Load :: FT_Entity ) { // Edge load in xi-eta system
191  boundaryLoad->computeValueAt(t, atTime, lcoords, VM_Total);
192  } else { // Edge load in x-y system
193  FloatArray gcoords;
194  this->pressureInterpolation.edgeLocal2global( gcoords, iEdge, lcoords, FEIElementGeometryWrapper(this) );
195  boundaryLoad->computeValueAt(t, atTime, gcoords, VM_Total);
196  }
197 
198  // Reshape the vector
199  for ( int j = 0; j < 3; j++ ) {
200  f(2 * j) += N(j) * t(0) * dS;
201  f(2 * j + 1) += N(j) * t(1) * dS;
202  }
203  }
204 
205  answer.assemble(f, this->edge_ordering [ iEdge - 1 ]);
206  } else {
207  OOFEM_ERROR("Strange boundary condition type");
208  }
209 }
210 
211 void
213 {
214  answer.resize(6);
215  answer.zero();
216 
217  double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity( integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
218  double mm = rho * this->area / 3.0;
219  for ( int i = 1; i <= 6; i++ ) {
220  answer.at(i) = mm;
221  }
222 }
223 
224 void
226 {
227  answer.resize(6, 6);
228  answer.zero();
229  double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity( integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
230  double mm = rho * this->area / 3.0;
231  for ( int i = 1; i <= 6; i++ ) {
232  answer.at(i, i) = mm;
233  }
234 }
235 
236 Interface *
238 {
239  return NULL;
240 }
241 
242 void
244 {
245  /* one should call material driver instead */
246  FloatArray u(6), eps(3);
247  answer.resize(3);
248 
249 
250  this->computeVectorOf(VM_Total, tStep, u);
251 
252  eps.at(1) = ( b [ 0 ] * u.at(1) + b [ 1 ] * u.at(3) + b [ 2 ] * u.at(5) );
253  eps.at(2) = ( c [ 0 ] * u.at(2) + c [ 1 ] * u.at(4) + c [ 2 ] * u.at(6) );
254  eps.at(3) = ( b [ 0 ] * u.at(2) + b [ 1 ] * u.at(4) + b [ 2 ] * u.at(6) + c [ 0 ] * u.at(1) + c [ 1 ] * u.at(3) + c [ 2 ] * u.at(5) );
255  FluidDynamicMaterial *mat = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial();
256  mat->computeDeviatoricStressVector(answer, gp, eps, tStep);
257 }
258 
259 void
261 {
262  answer.resize(6);
263  answer.zero();
264  FloatArray stress;
265 
266  this->computeDeviatoricStress(stress, integrationRulesArray [ 0 ]->getIntegrationPoint(0), atTime);
267 
268  // \int dNu/dxj \Tau_ij
269  for ( int i = 0; i < 3; i++ ) {
270  answer.at( ( i ) * 2 + 1 ) = area * ( stress.at(1) * b [ i ] + stress.at(3) * c [ i ] );
271  answer.at( ( i + 1 ) * 2 ) = area * ( stress.at(3) * b [ i ] + stress.at(2) * c [ i ] );
272  }
273 }
274 
275 int
277 {
278  Node *node1, *node2, *node3;
279  double x1, x2, x3, y1, y2, y3;
280 
281  node1 = giveNode(1);
282  node2 = giveNode(2);
283  node3 = giveNode(3);
284 
285  // init geometry data
286  x1 = node1->giveCoordinate(1);
287  x2 = node2->giveCoordinate(1);
288  x3 = node3->giveCoordinate(1);
289 
290  y1 = node1->giveCoordinate(2);
291  y2 = node2->giveCoordinate(2);
292  y3 = node3->giveCoordinate(2);
293 
294  this->area = 0.5 * ( x2 * y3 + x1 * y2 + y1 * x3 - x2 * y1 - x3 * y2 - x1 * y3 );
295 
296  b [ 0 ] = ( y2 - y3 ) / ( 2. * area );
297  c [ 0 ] = ( x3 - x2 ) / ( 2. * area );
298  b [ 1 ] = ( y3 - y1 ) / ( 2. * area );
299  c [ 1 ] = ( x1 - x3 ) / ( 2. * area );
300  b [ 2 ] = ( y1 - y2 ) / ( 2. * area );
301  c [ 2 ] = ( x2 - x1 ) / ( 2. * area );
302 
304 }
305 
306 double
308 {
310 
311  // assuming plane x-y !!!!!!!!!!!!!!!!!!!!!!
312 
313  FloatArray u;
314  this->computeVectorOf(VM_Total, tStep, u);
315  FloatArray u1(2), u2(2), u3(2);
316  u1.at(1) = u.at(1);
317  u1.at(2) = u.at(2);
318 
319  u2.at(1) = u.at(4);
320  u2.at(2) = u.at(5);
321 
322  u3.at(1) = u.at(7);
323  u3.at(2) = u.at(8);
324 
325  Node *node1 = giveNode(1);
326  Node *node2 = giveNode(2);
327  Node *node3 = giveNode(3);
328 
329  double x1, x2, x3, y1, y2, y3;
330 
331  x1 = node1->giveCoordinate(1);
332  x2 = node2->giveCoordinate(1);
333  x3 = node3->giveCoordinate(1);
334 
335  y1 = node1->giveCoordinate(2);
336  y2 = node2->giveCoordinate(2);
337  y3 = node3->giveCoordinate(2);
338 
339  FloatArray p1(2), p2(2), p3(2);
340  p1.at(1) = x1;
341  p1.at(2) = y1;
342 
343  p2.at(1) = x2;
344  p2.at(2) = y2;
345 
346  p3.at(1) = x3;
347  p3.at(2) = y3;
348 
349  FloatArray s12(p2), s23(p3), s31(p1);
350 
351  s12.subtract(p1);
352  FloatArray s21(s12);
353  s21.negated();
354 
355  s23.subtract(p2);
356  FloatArray s32(s23);
357  s32.negated();
358 
359  s31.subtract(p3);
360  FloatArray s13(s31);
361  s13.negated();
362 
363  FloatArray foot1(p2);
364  FloatArray foot2(p3);
365  FloatArray foot3(p1);
366 
367  double factor = s21.dotProduct(s23) / s23.dotProduct(s23);
368  foot1.add(factor, s23);
369 
370  FloatArray altitude1(foot1);
371  altitude1.subtract(p1);
372 
373  factor = s32.dotProduct(s31) / s31.dotProduct(s31);
374  foot2.add(factor, s31);
375 
376  FloatArray altitude2(foot2);
377  altitude2.subtract(p2);
378 
379  factor = s13.dotProduct(s12) / s12.dotProduct(s12);
380  foot3.add(factor, s12);
381 
382  FloatArray altitude3(foot3);
383  altitude3.subtract(p3);
384 
385  double dt1(deltaT), dt2(deltaT), dt3(deltaT);
386 
387  double u1_proj = u1.dotProduct(altitude1) / altitude1.computeNorm();
388  if ( u1_proj > 1.e-6 ) {
389  dt1 = altitude1.computeNorm() / u1_proj;
390  }
391 
392  double u2_proj = u2.dotProduct(altitude2) / altitude2.computeNorm();
393  if ( u2_proj > 1.e-6 ) {
394  dt2 = altitude1.computeNorm() / u2_proj;
395  }
396 
397  double u3_proj = u3.dotProduct(altitude3) / altitude3.computeNorm();
398  if ( u3_proj > 1.e-6 ) {
399  dt3 = altitude3.computeNorm() / u3_proj;
400  }
401 
402  double dt_min = min( dt1, min(dt2, dt3) );
403 
404  return dt_min;
405 }
406 
407 
409 {
410  contextIOResultType iores;
411 
412  if ( ( iores = PFEMElement2d :: saveContext(stream, mode, obj) ) != CIO_OK ) {
413  THROW_CIOERR(iores);
414  }
415 
416  return CIO_OK;
417 }
418 
419 
421 {
422  contextIOResultType iores;
423 
424  if ( ( iores = PFEMElement2d :: restoreContext(stream, mode, obj) ) != CIO_OK ) {
425  THROW_CIOERR(iores);
426  }
427 
428 
429  return CIO_OK;
430 }
431 
432 
433 double
435 // Returns the portion of the receiver which is attached to gp.
436 {
438  return determinant * gp->giveWeight();
439 }
440 
441 #ifdef __OOFEG
442 int
444  int node, TimeStep *atTime)
445 {
446  //<RESTRICTED_SECTION>
447  if ( type == IST_VOFFraction ) {
448  answer.resize(1);
449  // answer.at(1) = this->giveTempVolumeFraction();
450  return 1;
451  } else
452  //</RESTRICTED_SECTION>
453  if ( type == IST_Density ) {
454  answer.resize(1);
455  answer.at(1) = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity( integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
456  return 1;
457  } else {
458  return PFEMElement :: giveInternalStateAtNode(answer, type, mode, node, atTime);
459  }
460 }
461 
462 void
464 {
465  WCRec p [ 3 ];
466  GraphicObj *go;
467 
468  if ( !gc.testElementGraphicActivity(this) ) {
469  return;
470  }
471 
472  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
473  EASValsSetColor( gc.getElementColor() );
474  EASValsSetEdgeColor( gc.getElementEdgeColor() );
475  EASValsSetEdgeFlag(TRUE);
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 = 0.;
480  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
481  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
482  p [ 1 ].z = 0.;
483  p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveCoordinate(1);
484  p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveCoordinate(2);
485  p [ 2 ].z = 0.;
486 
487  go = CreateTriangle3D(p);
488  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
489  EGAttachObject(go, ( EObjectP ) this);
490  EMAddGraphicsToModel(ESIModel(), go);
491 }
492 
494 {
495  int i, indx, result = 0;
496  WCRec p [ 3 ];
497  GraphicObj *tr;
498  TimeStep *tStep = this->giveDomain()->giveEngngModel()->giveCurrentStep();
499  FloatArray v1, v2, v3;
500  double s [ 3 ];
501  IntArray map;
502 
503  if ( !context.testElementGraphicActivity(this) ) {
504  return;
505  }
506 
507  if ( context.giveIntVarMode() == ISM_recovered ) {
508  result += this->giveInternalStateAtNode(v1, context.giveIntVarType(), context.giveIntVarMode(), 1, tStep);
509  result += this->giveInternalStateAtNode(v2, context.giveIntVarType(), context.giveIntVarMode(), 2, tStep);
510  result += this->giveInternalStateAtNode(v3, context.giveIntVarType(), context.giveIntVarMode(), 3, tStep);
511  } else if ( context.giveIntVarMode() == ISM_local ) {
512  GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
513  result += giveIPValue(v1, gp, context.giveIntVarType(), tStep);
514  v2 = v1;
515  v3 = v1;
516  result *= 3;
517  }
518 
519  if ( result != 3 ) {
520  return;
521  }
522 
523  indx = context.giveIntVarIndx();
524 
525  s [ 0 ] = v1.at(indx);
526  s [ 1 ] = v2.at(indx);
527  s [ 2 ] = v3.at(indx);
528 
529  EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
530 
531  if ( context.getScalarAlgo() == SA_ISO_SURF ) {
532  for ( i = 0; i < 3; i++ ) {
533  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
534  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
535  p [ i ].z = 0.;
536  }
537 
538  //EASValsSetColor(gc.getYieldPlotColor(ratio));
539  context.updateFringeTableMinMax(s, 3);
540  tr = CreateTriangleWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ]);
541  EGWithMaskChangeAttributes(LAYER_MASK, tr);
542  EMAddGraphicsToModel(ESIModel(), tr);
543  } else if ( ( context.getScalarAlgo() == SA_ZPROFILE ) || ( context.getScalarAlgo() == SA_COLORZPROFILE ) ) {
544  double landScale = context.getLandScale();
545 
546  for ( i = 0; i < 3; i++ ) {
547  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
548  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
549  p [ i ].z = s [ i ] * landScale;
550  }
551 
552  if ( context.getScalarAlgo() == SA_ZPROFILE ) {
553  EASValsSetColor( context.getDeformedElementColor() );
554  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
555  EASValsSetFillStyle(FILL_SOLID);
556  tr = CreateTriangle3D(p);
557  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | FILL_MASK | LAYER_MASK, tr);
558  } else {
559  context.updateFringeTableMinMax(s, 3);
560  EASValsSetFillStyle(FILL_SOLID);
561  tr = CreateTriangleWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ]);
562  EGWithMaskChangeAttributes(FILL_MASK | LAYER_MASK, tr);
563  }
564 
565  EMAddGraphicsToModel(ESIModel(), tr);
566  }
567 }
568 
569 
570 
571 #endif
572 } // end namespace oofem
virtual void edgeLocal2global(FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates edge global coordinates from given local ones.
Definition: fei2dtrlin.C:213
CrossSection * giveCrossSection()
Definition: element.C:495
virtual void drawRawGeometry(oofegGraphicContext &)
Definition: tr1_2d_pfem.C:463
virtual void edgeEvalN(FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of edge interpolation functions (shape functions) at given point.
Definition: fei2dtrlin.C:184
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.
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 int checkConsistency()
Performs consistency check.
Definition: pfemelement2d.C:77
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: element.C:1257
virtual int computeNumberOfDofs()
Computes or simply returns total number of element&#39;s local DOFs.
Definition: tr1_2d_pfem.C:98
Class and object Domain.
Definition: domain.h:115
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
Definition: fei2dtrlin.C:43
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Returns local vector of unknowns.
Definition: element.C:86
Fluid cross-section.
ScalarAlgorithmType getScalarAlgo()
Abstract base class for all fluid materials.
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition: femcmpnn.h:82
virtual void computeEdgeBCSubVectorAt(FloatArray &answer, Load *load, int iEdge, TimeStep *tStep)
Definition: tr1_2d_pfem.C:168
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
IRResultType initializeFrom(InputRecord *ir)
Initializes receiver acording to object description stored in input record.
Definition: pfemelement.C:69
void zero()
Sets all component to zero.
Definition: intarray.C:52
TR1_2D_PFEM(int n, Domain *aDomain, int particle1, int particle2, int particle3, int mat, int cs)
Constructor of TR1_2D_PFEM Element.
Definition: tr1_2d_pfem.C:78
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
Class implementing element body load, acting over whole element volume (e.g., the dead weight)...
Definition: bodyload.h:49
#define OOFEG_RAW_GEOMETRY_LAYER
static IntArray velocityDofMask
Mask of velocity Dofs.
Definition: tr1_2d_pfem.h:83
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: tr1_2d_pfem.C:129
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
EngngModel * giveEngngModel()
Returns engineering model to which receiver is associated.
Definition: domain.C:433
virtual double edgeGiveTransformationJacobian(int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the edge Jacobian of transformation between local and global coordinates.
Definition: feinterpol2d.C:175
virtual void computeDiagonalMassMtrx(FloatArray &answer, TimeStep *tStep)
Calculates diagonal mass matrix as vector.
Definition: tr1_2d_pfem.C:212
virtual void computeDeviatoricStress(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Computes deviatoric stress vector in given integration point and solution step from given total strai...
Definition: tr1_2d_pfem.C:243
virtual void computeComponentArrayAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Computes boundary condition value - its components values at given time.
Definition: load.C:82
virtual double giveCoordinate(int i)
Definition: node.C:82
virtual int giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *atTime)
Returns internal state variable (like stress,strain) at node of element in Reduced form...
Definition: pfemelement.C:235
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
#define THROW_CIOERR(e)
Definition: contextioerr.h:61
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
double giveTimeIncrement()
Returns solution step associated time increment.
Definition: timestep.h:150
virtual void giveElementDofIDMask(IntArray &answer) const
Returns element dof mask for node.
Definition: tr1_2d_pfem.C:113
InternalStateType giveIntVarType()
static FEI2dTrLin velocityInterpolation
Interpolation for velocity unknowns.
Definition: tr1_2d_pfem.h:76
virtual void computeDeviatoricStressVector(FloatArray &stress_dev, double &epsp_vol, GaussPoint *gp, const FloatArray &eps, double pressure, TimeStep *tStep)
Computes the deviatoric stress vector and volumetric strain rate from given deviatoric strain and pre...
Abstract base class representing a boundary load (force, momentum, ...) that acts directly on a bound...
Definition: boundaryload.h:110
double dotProduct(const FloatArray &x) const
Computes the dot product (or inner product) of receiver and argument.
Definition: floatarray.C:463
#define OOFEM_ERROR(...)
Definition: error.h:61
int material
Number of associated material.
Definition: element.h:153
#define OOFEG_RAW_GEOMETRY_WIDTH
virtual double giveWeight()
Returns integration weight of receiver.
Definition: gausspoint.h:181
virtual FormulationType giveFormulationType()
Specifies is load should take local or global coordinates.
Definition: load.h:155
Neumann type (prescribed flux).
Definition: bctype.h:43
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
virtual void setCrossSection(int csIndx)
Sets the cross section model of receiver.
Definition: element.h:653
#define N(p, q)
Definition: mdm.C:367
virtual int giveApproxOrder()=0
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
static FEI2dTrLin pressureInterpolation
Interpolation for pressure unknowns.
Definition: tr1_2d_pfem.h:78
void setDofManagers(const IntArray &dmans)
Sets receiver dofManagers.
Definition: element.C:550
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: tr1_2d_pfem.C:434
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj)
Stores receiver state to output stream.
Definition: tr1_2d_pfem.C:408
InternalStateMode giveIntVarMode()
~TR1_2D_PFEM()
Destructor.
Definition: tr1_2d_pfem.C:93
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual double computeCriticalTimeStep(TimeStep *tStep)
Calculates critical time step.
Definition: tr1_2d_pfem.C:307
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 double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
Definition: fei2dtrlin.C:147
virtual void computeDeviatoricStressDivergence(FloatArray &answer, TimeStep *tStep)
Calculates the divergence of the deviatoric stress.
Definition: tr1_2d_pfem.C:260
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
double computeNorm() const
Computes the norm (or length) of the vector.
Definition: floatarray.C:840
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
void setBodyLoads(const IntArray &bodyLoads)
Sets receiver bodyLoadArray.
Definition: element.C:556
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
Definition: tr1_2d_pfem.C:104
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
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
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
int setUpIntegrationPoints(integrationDomain intdomain, int nPoints, MaterialMode matMode)
Initializes the receiver.
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
Domain * giveDomain() const
Definition: femcmpnn.h:100
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver acording to object description stored in input record.
Definition: tr1_2d_pfem.C:120
virtual TimeStep * giveSolutionStepWhenIcApply(bool force=false)
Returns the solution step when Initial Conditions (IC) apply.
Definition: engngm.h:720
int min(int i, int j)
Returns smaller value from two given decimals.
Definition: mathfem.h:59
void updateFringeTableMinMax(double *s, int size)
Load is base abstract class for all loads.
Definition: load.h:61
static IntArray pressureDofMask
Mask of pressure Dofs.
Definition: tr1_2d_pfem.h:85
virtual TimeStep * giveCurrentStep(bool force=false)
Returns current time step.
Definition: engngm.h:683
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
virtual void computeBodyLoadVectorAt(FloatArray &answer, BodyLoad *load, TimeStep *tStep, ValueModeType mode)
Computes the load vector due to body load acting on receiver, at given time step. ...
Definition: tr1_2d_pfem.C:141
the oofem namespace is to define a context or scope in which all oofem names are defined.
Class implementing node in finite element mesh.
Definition: node.h:87
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj)
Restores the receiver state previously written in stream.
Definition: tr1_2d_pfem.C:420
void negated()
Switches the sign of every coefficient of receiver.
Definition: floatarray.C:739
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
virtual void computeValueAt(FloatArray &answer, TimeStep *tStep, const FloatArray &coords, ValueModeType mode)
Computes components values of load at given point - global coordinates (coordinates given)...
Definition: boundaryload.C:58
virtual void drawScalar(oofegGraphicContext &context)
Definition: tr1_2d_pfem.C:493
#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
This class is the implementation of general 2d element with arbitrary interpolation of velocity and p...
Definition: pfemelement2d.h:68
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
virtual int checkConsistency()
Performs consistency check.
Definition: tr1_2d_pfem.C:276
virtual int giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *atTime)
Returns internal state variable (like stress,strain) at node of element in Reduced form...
Definition: tr1_2d_pfem.C:443
InternalStateMode
Determines the mode of internal variable.
void add(const FloatArray &src)
Adds array src to receiver.
Definition: floatarray.C:156
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
virtual Interface * giveInterface(InterfaceType)
Interface requesting service.
Definition: tr1_2d_pfem.C:237
static IntArray edge_ordering[3]
Definition: tr1_2d_pfem.h:80
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:31 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011