OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
tr1_2d_cbs.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_cbs.h"
36 #include "fei2dtrlin.h"
37 #include "cbs.h"
38 #include "node.h"
39 #include "material.h"
40 #include "gausspoint.h"
41 #include "gaussintegrationrule.h"
42 #include "floatmatrix.h"
43 #include "floatarray.h"
44 #include "intarray.h"
45 #include "domain.h"
46 #include "mathfem.h"
47 #include "engngm.h"
48 #include "fluiddynamicmaterial.h"
49 #include "fluidcrosssection.h"
50 #include "load.h"
51 #include "timestep.h"
52 #include "boundaryload.h"
53 #include "geotoolbox.h"
54 #include "crosssection.h"
55 #include "contextioerr.h"
56 #include "dynamicinputrecord.h"
57 #include "classfactory.h"
58 
59 #ifdef __OOFEG
60  #include "oofeggraphiccontext.h"
61  #include "connectivitytable.h"
62 #endif
63 
64 namespace oofem {
65 #define TRSUPG_ZERO_VOF 1.e-8
66 
67 REGISTER_Element(TR1_2D_CBS);
68 
69 FEI2dTrLin TR1_2D_CBS :: interp(1, 2);
70 
71 TR1_2D_CBS :: TR1_2D_CBS(int n, Domain *aDomain) :
72  CBSElement(n, aDomain),
75  //<RESTRICTED_SECTION>
77  //</RESTRICTED_SECTION>
78 {
79  // Constructor.
80  numberOfDofMans = 3;
81 }
82 
84 { }
85 
88 
89 int
91 {
92  return 9;
93 }
94 
95 void
97 {
98  answer = {V_u, V_v, P_f};
99 }
100 
101 
104 {
105  //<RESTRICTED_SECTION>
106  IRResultType result; // Required by IR_GIVE_FIELD macro
107 
108  this->vof = 0.0;
110  if ( vof > 0.0 ) {
112  this->temp_vof = vof;
113  } else {
114  this->vof = 0.0;
116  this->temp_vof = this->vof;
117  }
118 
119  //</RESTRICTED_SECTION>
120 
121  return CBSElement :: initializeFrom(ir);
122 }
123 
124 
125 void
127 {
129  if ( this->permanentVofFlag ) {
130  input.setField(this->vof, _IFT_Tr1CBS_pvof);
131  } else {
132  input.setField(this->vof, _IFT_Tr1CBS_vof);
133  }
134 }
135 
136 
137 void
139 // Sets up the array containing the four Gauss points of the receiver.
140 {
141  if ( integrationRulesArray.size() == 0 ) {
142  integrationRulesArray.resize(1);
143  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 3) );
145  }
146 }
147 
148 void
150 {
151  answer.resize(9, 9);
152  answer.zero();
153  double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity( integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
154 
155  double ar6 = rho * area / 6.0;
156  double ar12 = rho * area / 12.0;
157 
158  answer.at(1, 1) = answer.at(2, 2) = ar6;
159  answer.at(4, 4) = answer.at(5, 5) = ar6;
160  answer.at(7, 7) = answer.at(8, 8) = ar6;
161 
162  answer.at(1, 4) = answer.at(1, 7) = ar12;
163  answer.at(4, 1) = answer.at(4, 7) = ar12;
164  answer.at(7, 1) = answer.at(7, 4) = ar12;
165 
166  answer.at(2, 5) = answer.at(2, 8) = ar12;
167  answer.at(5, 2) = answer.at(5, 8) = ar12;
168  answer.at(8, 2) = answer.at(8, 5) = ar12;
169 }
170 
171 
172 void
174 {
175  answer.resize(9);
176  answer.zero();
177 
178  double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity( integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
179  double mm = rho * this->area / 3.0;
180  for ( int i = 0; i < 3; i++ ) {
181  answer.at(i * 3 + 1) = mm;
182  answer.at(i * 3 + 2) = mm;
183  }
184 }
185 
186 void
188 {
189  // calculates advection component for (*) velocities
190 
191  double ar12, ar3, dudx, dudy, dvdx, dvdy;
192  double usum, vsum;
193  double adu11, adu21, adu31, adv11, adv21, adv31;
194  double adu12, adu22, adu32, adv12, adv22, adv32;
195  double dt = tStep->giveTimeIncrement();
196  //double rho =static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
197  double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity( integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
198  int nLoads;
199  bcGeomType ltype;
200  Load *load;
201  FloatArray gVector;
202 
203  FloatArray u;
204 
205  ar12 = area / 12.;
206  ar3 = area / 3.0;
207 
208  this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), u);
209 
210  dudx = b [ 0 ] * u.at(1) + b [ 1 ] * u.at(3) + b [ 2 ] * u.at(5);
211  dudy = c [ 0 ] * u.at(1) + c [ 1 ] * u.at(3) + c [ 2 ] * u.at(5);
212  dvdx = b [ 0 ] * u.at(2) + b [ 1 ] * u.at(4) + b [ 2 ] * u.at(6);
213  dvdy = c [ 0 ] * u.at(2) + c [ 1 ] * u.at(4) + c [ 2 ] * u.at(6);
214 
215  usum = u.at(1) + u.at(3) + u.at(5);
216  vsum = u.at(2) + u.at(4) + u.at(6);
217 
218  // compute Cu*U term
219 
220  adu11 = ar12 * ( dudx * ( usum + u.at(1) ) + dudy * ( vsum + u.at(2) ) );
221  adu21 = ar12 * ( dudx * ( usum + u.at(3) ) + dudy * ( vsum + u.at(4) ) );
222  adu31 = ar12 * ( dudx * ( usum + u.at(5) ) + dudy * ( vsum + u.at(6) ) );
223  adv11 = ar12 * ( dvdx * ( usum + u.at(1) ) + dvdy * ( vsum + u.at(2) ) );
224  adv21 = ar12 * ( dvdx * ( usum + u.at(3) ) + dvdy * ( vsum + u.at(4) ) );
225  adv31 = ar12 * ( dvdx * ( usum + u.at(5) ) + dvdy * ( vsum + u.at(6) ) );
226 
227 
228  // compute Ku*U term
229  double uu = ( usum * usum + u.at(1) * u.at(1) + u.at(3) * u.at(3) + u.at(5) * u.at(5) );
230  double uv = ( usum * vsum + u.at(1) * u.at(2) + u.at(3) * u.at(4) + u.at(5) * u.at(6) );
231  double vv = ( vsum * vsum + u.at(2) * u.at(2) + u.at(4) * u.at(4) + u.at(6) * u.at(6) );
232 
233  adu12 = dt * 0.5 * ar12 * ( b [ 0 ] * dudx * uu + b [ 0 ] * dudy * uv + c [ 0 ] * dudx * uv + c [ 0 ] * dudy * vv );
234  adu22 = dt * 0.5 * ar12 * ( b [ 1 ] * dudx * uu + b [ 1 ] * dudy * uv + c [ 1 ] * dudx * uv + c [ 1 ] * dudy * vv );
235  adu32 = dt * 0.5 * ar12 * ( b [ 2 ] * dudx * uu + b [ 2 ] * dudy * uv + c [ 2 ] * dudx * uv + c [ 2 ] * dudy * vv );
236  adv12 = dt * 0.5 * ar12 * ( b [ 0 ] * dvdx * uu + b [ 0 ] * dvdy * uv + c [ 0 ] * dvdx * uv + c [ 0 ] * dvdy * vv );
237  adv22 = dt * 0.5 * ar12 * ( b [ 1 ] * dvdx * uu + b [ 1 ] * dvdy * uv + c [ 1 ] * dvdx * uv + c [ 1 ] * dvdy * vv );
238  adv32 = dt * 0.5 * ar12 * ( b [ 2 ] * dvdx * uu + b [ 2 ] * dvdy * uv + c [ 2 ] * dvdx * uv + c [ 2 ] * dvdy * vv );
239 
240  answer.resize(9);
241  answer.zero();
242 
243  answer.at(1) = -adu11 - adu12;
244  answer.at(2) = -adv11 - adv12;
245  answer.at(4) = -adu21 - adu22;
246  answer.at(5) = -adv21 - adv22;
247  answer.at(7) = -adu31 - adu32;
248  answer.at(8) = -adv31 - adv32;
249 
250 
251  // body load (gravity effects)
252  nLoads = this->giveBodyLoadArray()->giveSize();
253  for ( int i = 1; i <= nLoads; i++ ) {
254  load = domain->giveLoad( bodyLoadArray.at(i) );
255  ltype = load->giveBCGeoType();
256  if ( ( ltype == BodyLoadBGT ) && ( load->giveBCValType() == ForceLoadBVT ) ) {
257  load->computeComponentArrayAt(gVector, tStep, VM_Total);
258  if ( gVector.giveSize() ) {
259  answer.at(1) -= dt * 0.5 * ar3 * ( b [ 0 ] * usum + c [ 0 ] * vsum ) * gVector.at(1);
260  answer.at(2) -= dt * 0.5 * ar3 * ( b [ 0 ] * usum + c [ 0 ] * vsum ) * gVector.at(2);
261  answer.at(4) -= dt * 0.5 * ar3 * ( b [ 1 ] * usum + c [ 1 ] * vsum ) * gVector.at(1);
262  answer.at(5) -= dt * 0.5 * ar3 * ( b [ 1 ] * usum + c [ 1 ] * vsum ) * gVector.at(2);
263  answer.at(7) -= dt * 0.5 * ar3 * ( b [ 2 ] * usum + c [ 2 ] * vsum ) * gVector.at(1);
264  answer.at(8) -= dt * 0.5 * ar3 * ( b [ 2 ] * usum + c [ 2 ] * vsum ) * gVector.at(2);
265  }
266  }
267  }
268 
269 
270  answer.times(rho);
271 }
272 
273 void
275 {
276  FloatArray stress;
277  double Re = static_cast< FluidModel * >( domain->giveEngngModel() )->giveReynoldsNumber();
278  GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
279  double coeff, rho = static_cast< FluidCrossSection* >(this->giveCrossSection())->giveDensity( gp );
280  FloatArray gVector;
281  double ar3 = area / 3.0;
282 
283  stress = static_cast< FluidDynamicMaterialStatus * >( gp->giveMaterialStatus() )->giveDeviatoricStressVector();
284  stress.times(1. / Re);
285 
286  // \int dNu/dxj \Tau_ij
287  answer.resize(9);
288  answer.zero();
289  for ( int i = 0; i < 3; i++ ) {
290  answer.at(i * 3 + 1) = -area * ( stress.at(1) * b [ i ] + stress.at(3) * c [ i ] );
291  answer.at(i * 3 + 2) = -area * ( stress.at(3) * b [ i ] + stress.at(2) * c [ i ] );
292  }
293 
294  // add boundary termms
295  coeff = ar3 * rho;
296  // body load (gravity effects)
297  int nLoads = this->giveBodyLoadArray()->giveSize();
298  for ( int i = 1; i <= nLoads; i++ ) {
299  auto load = domain->giveLoad( bodyLoadArray.at(i) );
300  bcGeomType ltype = load->giveBCGeoType();
301  if ( ltype == BodyLoadBGT && load->giveBCValType() == ForceLoadBVT ) {
302  load->computeComponentArrayAt(gVector, tStep, VM_Total);
303  if ( gVector.giveSize() ) {
304  answer.at(1) += coeff * gVector.at(1);
305  answer.at(2) += coeff * gVector.at(2);
306  answer.at(4) += coeff * gVector.at(1);
307  answer.at(5) += coeff * gVector.at(2);
308  answer.at(7) += coeff * gVector.at(1);
309  answer.at(8) += coeff * gVector.at(2);
310  }
311  }
312  }
313 
314  // terms now computed even for boundary nodes that have prescribed velocities! (but they will not be localized)
315  // loop over sides
316  int n1, n2;
317  double tx, ty, l, nx, ny;
318  nLoads = boundarySides.giveSize();
319  for ( int j = 1; j <= nLoads; j++ ) {
320  int code = boundaryCodes.at(j);
321  if ( ( code & FMElement_PrescribedTractionBC ) ) {
322  //printf ("TR1_2D_CBS :: computeDiffusionTermsI tractions detected\n");
323 
324  //_error("computeDiffusionTermsI: traction bc not supported");
325  FloatArray t, coords(1);
326 
327  BoundaryLoad *load;
328  // integrate tractions
329  n1 = boundarySides.at(j);
330  n2 = ( n1 == 3 ? 1 : n1 + 1 );
331 
332  tx = giveNode(n2)->giveCoordinate(1) - giveNode(n1)->giveCoordinate(1);
333  ty = giveNode(n2)->giveCoordinate(2) - giveNode(n1)->giveCoordinate(2);
334  l = sqrt(tx * tx + ty * ty);
335  nx = ty / l;
336  ny = -tx / l;
337 
338  // loop over boundary load array
339  int numLoads = this->giveBoundaryLoadArray()->giveSize() / 2;
340  for ( int i = 1; i <= numLoads; i++ ) {
341  int n = boundaryLoadArray.at(1 + ( i - 1 ) * 2);
342  load = dynamic_cast< BoundaryLoad * >( domain->giveLoad(n) );
343  if ( load ) {
344  load->computeValueAt(t, tStep, coords, VM_Total);
345 
346  //printf ("TR1_2D_CBS :: computeDiffusionTermsI traction (%e,%e) detected\n", t.at(1), t.at(2));
347 
348  answer.at( ( n1 - 1 ) * 3 + 1 ) += 0.5 * l * ( t.at(1) * nx );
349  answer.at( ( n1 - 1 ) * 3 + 2 ) += 0.5 * l * ( t.at(2) * ny );
350 
351  answer.at( ( n2 - 1 ) * 3 + 1 ) += 0.5 * l * ( t.at(1) * nx );
352  answer.at( ( n2 - 1 ) * 3 + 2 ) += 0.5 * l * ( t.at(2) * ny );
353  }
354  }
355  } else if ( !( ( code & FMElement_PrescribedUnBC ) && ( code & FMElement_PrescribedUsBC ) ) ) {
356  /*
357  * n1 = i;
358  * n2 = (n1==3?n2=1:n2=n1+1);
359  *
360  * //if (giveNode(n1)->isBoundary() && giveNode(n2)->isBoundary()) {
361  *
362  * tx = giveNode(n2)->giveCoordinate(1)-giveNode(n1)->giveCoordinate(1);
363  * ty = giveNode(n2)->giveCoordinate(2)-giveNode(n1)->giveCoordinate(2);
364  * l = sqrt(tx*tx+ty*ty);
365  * nx = ty/l; ny = -tx/l;
366  *
367  * // normal displacement precribed
368  * answer.at((n1-1)*3+1)+=l*0.5*(stress.at(1)*nx + stress.at(3)*ny);
369  * answer.at((n1-1)*3+2)+=l*0.5*(stress.at(3)*nx + stress.at(2)*ny);
370  *
371  * answer.at((n2-1)*3+1)+=l*0.5*(stress.at(1)*nx + stress.at(3)*ny);
372  * answer.at((n2-1)*3+2)+=l*0.5*(stress.at(3)*nx + stress.at(2)*ny);
373  * //}
374  */
375  }
376  }
377 }
378 
379 void
381 {
382  // computes velocity terms on RHS for density equation
383  double velu = 0.0, velv = 0.0; // dudx=0.0, dvdy=0.0;
384  //double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
385  double theta1 = static_cast< CBS * >( domain->giveEngngModel() )->giveTheta1();
386  double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity( integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
387  FloatArray u, ustar;
388 
389  answer.resize(9);
390  answer.zero();
391 
392  this->computeVectorOfVelocities(VM_Total, tStep, u);
394  // Should we really want to add this term?
395  // This will produce velocities that differs from their actual Dirichlet b.c.s;
396  // The patch05.in test has a Dirichlet V_v = 1.0 on node 2, but this will produce the value V_v 1.5 for that dof.
397  // It also requires knowing theta1, which we could otherwise do without.
398  // It is added for now to mirror the old code below.
399  this->computeVectorOfPrescribed({V_u, V_v}, VM_Incremental, tStep, ustar);
400  u.add(theta1, ustar);
401 
402  //this->computeVectorOfVelocities(VM_Incremental, tStep, ustar);
403  //u.add((theta1-1, ustar);
404  /*
405  * printf("u_new = "); u.printYourself();
406  * this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), u);
407  * this->computeVectorOfAuxVelocities(VM_Incremental, tStep, ustar);
408  * u.add(theta1, ustar);
409  * printf("u_corr = "); u.printYourself();
410  */
411 
412  for ( int i = 0; i < 3; i++ ) {
413  velu += u.at(i * 2 + 1);
414  velv += u.at(i * 2 + 2);
415  }
416 
417  for ( int i = 0; i < 3; i++ ) {
418  answer.at( ( i + 1 ) * 3 ) = area * ( ( b [ i ] * velu + c [ i ] * velv ) ) / 3.0;
419  }
420 
421  /* account for normal prescribed velocity on boundary*/
422  /* on the rest of the boundary the tractions or pressure is prescribed -> to be implemented later */
423  int n1, n2;
424  double tx, ty, l, nx, ny, un1, un2;
425 
426 
427  // loop over sides
428  int code;
429  for ( int j = 1; j <= boundarySides.giveSize(); j++ ) {
430  code = boundaryCodes.at(j);
431  if ( ( code & FMElement_PrescribedPressureBC ) ) {
432  continue;
433  } else if ( ( code & FMElement_PrescribedUnBC ) ) {
434  this->computeVectorOfPrescribed({V_u, V_v}, VM_Total, tStep, u);
435 
436  n1 = boundarySides.at(j);
437  n2 = ( n1 == 3 ? 1 : n1 + 1 );
438 
439  //if (giveNode(n1)->isBoundary() && giveNode(n2)->isBoundary()) {
440 
441  tx = giveNode(n2)->giveCoordinate(1) - giveNode(n1)->giveCoordinate(1);
442  ty = giveNode(n2)->giveCoordinate(2) - giveNode(n1)->giveCoordinate(2);
443  l = sqrt(tx * tx + ty * ty);
444  nx = ty / l;
445  ny = -tx / l;
446  un1 = nx * u.at( ( n1 - 1 ) * 2 + 1 ) + ny *u.at(n1 * 2);
447  un2 = nx * u.at( ( n2 - 1 ) * 2 + 1 ) + ny *u.at(n2 * 2);
448  //if ((un1 != 0.) && (un2 != 0.)) {
449  // normal displacement precribed
450  answer.at(n1 * 3) -= ( un1 * l / 3. + un2 * l / 6. );
451  answer.at(n2 * 3) -= ( un2 * l / 3. + un1 * l / 6. );
452  //}
453  //}
454  } else if ( ( code & FMElement_PrescribedTractionBC ) ) {
455  continue;
456  }
457  }
458 
459  answer.times(rho);
460 }
461 
462 
463 void
465 {
466  /*
467  * this method computes the prescribed pressure due to applied traction
468  * p = tau(i,j)*n(i)*n(j) - traction(i)*n(i)
469  * this pressure is enforced as dirichlet bc in density/pressure equation
470  */
471  FloatArray stress;
472  double Re = static_cast< FluidModel * >( domain->giveEngngModel() )->giveReynoldsNumber();
473  GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
474  stress = static_cast< FluidDynamicMaterialStatus * >( gp->giveMaterialStatus() )->giveDeviatoricStressVector();
475  stress.times(1. / Re);
476 
477  answer.resize(9);
478  answer.zero();
479 
480  // loop over sides
481  int n1, n2, code, sid;
482  double tx, ty, l, nx, ny, pcoeff;
483  //IntArray nodecounter (3);
484  for ( int j = 1; j <= boundarySides.giveSize(); j++ ) {
485  code = boundaryCodes.at(j);
486  sid = boundarySides.at(j);
487  if ( ( code & FMElement_PrescribedTractionBC ) ) {
488  FloatArray t, coords(1);
489  int nLoads, n, id;
490  BoundaryLoad *load;
491  // integrate tractions
492  n1 = sid;
493  n2 = ( n1 == 3 ? 1 : n1 + 1 );
494 
495  tx = giveNode(n2)->giveCoordinate(1) - giveNode(n1)->giveCoordinate(1);
496  ty = giveNode(n2)->giveCoordinate(2) - giveNode(n1)->giveCoordinate(2);
497  l = sqrt(tx * tx + ty * ty);
498  nx = ty / l;
499  ny = -tx / l;
500 
501 
502  //nodecounter.at(n1)++;
503  //nodecounter.at(n2)++;
504  pcoeff = stress.at(1) * nx * nx + stress.at(2) * ny * ny + 2.0 * stress.at(3) * nx * ny;
505  answer.at(n1 * 3) += pcoeff;
506  answer.at(n2 * 3) += pcoeff;
507 
508  // if no traction bc applied but side marked as with traction load
509  // then zero traction is assumed !!!
510 
511  // loop over boundary load array
512  nLoads = this->giveBoundaryLoadArray()->giveSize() / 2;
513  for ( int i = 1; i <= nLoads; i++ ) {
514  n = boundaryLoadArray.at(1 + ( i - 1 ) * 2);
515  id = boundaryLoadArray.at(i * 2);
516  if ( id != sid ) {
517  continue;
518  }
519 
520  load = dynamic_cast< BoundaryLoad * >( domain->giveLoad(n) );
521  if ( load ) {
522  load->computeValueAt(t, tStep, coords, VM_Total);
523 
524  answer.at(n1 * 3) -= t.at(1) * nx + t.at(2) * ny;
525  answer.at(n2 * 3) -= t.at(1) * nx + t.at(2) * ny;
526  }
527  }
528  }
529  }
530 
531  //for (i=1; i<=3; i++) {
532  // if (nodecounter.at(i)) answer.at(i) /= nodecounter.at(i);
533  //}
534 }
535 
536 
537 void
539 {
540  /*
541  * this method computes the prescribed pressure due to applied traction
542  * p = tau(i,j)*n(i)*n(j) - traction(i)*n(i)
543  * this pressure is enforced as dirichlet bc in density/pressure equation
544  */
545  answer.resize(9);
546  answer.zero();
547 
548  // loop over sides
549  int n1, n2, code;
550  IntArray nodecounter(3);
551  for ( int j = 1; j <= boundarySides.giveSize(); j++ ) {
552  code = boundaryCodes.at(j);
553  if ( ( code & FMElement_PrescribedTractionBC ) ) {
554  n1 = boundarySides.at(j);
555  n2 = ( n1 == 3 ? 1 : n1 + 1 );
556  answer.at(n1 * 3)++;
557  answer.at(n2 * 3)++;
558  }
559  }
560 }
561 
562 
563 
564 
565 void
567 {
568  // computes pressure terms on RHS for density equation
569  FloatArray p;
570  double theta1 = static_cast< CBS * >( domain->giveEngngModel() )->giveTheta1();
571 
572  this->computeVectorOfPressures(VM_Total, tStep->givePreviousStep(), p);
573  answer.resize(9);
574  answer.zero();
575 
576  double dpdx = 0.0, dpdy = 0.0;
577  for ( int i = 0; i < 3; i++ ) {
578  dpdx += b [ i ] * p.at(i + 1);
579  dpdy += c [ i ] * p.at(i + 1);
580  }
581 
582  for ( int i = 0; i < 3; i++ ) {
583  answer.at( ( i + 1 ) * 3 ) = -theta1 *tStep->giveTimeIncrement() * area * ( b [ i ] * dpdx + c [ i ] * dpdy );
584  }
585 }
586 
587 
588 void
590 {
591  // calculates the pressure LHS
592 
593  answer.resize(3, 3);
594 
595  answer.at(1, 1) = area * ( b [ 0 ] * b [ 0 ] + c [ 0 ] * c [ 0 ] );
596  answer.at(2, 2) = area * ( b [ 1 ] * b [ 1 ] + c [ 1 ] * c [ 1 ] );
597  answer.at(3, 3) = area * ( b [ 2 ] * b [ 2 ] + c [ 2 ] * c [ 2 ] );
598 
599  answer.at(1, 2) = answer.at(2, 1) = area * ( b [ 0 ] * b [ 1 ] + c [ 0 ] * c [ 1 ] );
600  answer.at(1, 3) = answer.at(3, 1) = area * ( b [ 0 ] * b [ 2 ] + c [ 0 ] * c [ 2 ] );
601  answer.at(2, 3) = answer.at(3, 2) = area * ( b [ 1 ] * b [ 2 ] + c [ 1 ] * c [ 2 ] );
602 }
603 
604 
605 void
607 {
608  //Evaluates the RHS of velocity correction step
609  FloatArray p, u;
610  double pn1, ar3;
611  double usum, vsum, coeff;
612 
613 
614  this->computeVectorOfPressures(VM_Total, tStep, p);
615 
616  double dpdx = 0.0, dpdy = 0.0;
617  for ( int i = 0; i < 3; i++ ) {
618  pn1 = p.at(i + 1);
619  dpdx += b [ i ] * pn1;
620  dpdy += c [ i ] * pn1;
621  }
622 
623  answer.resize(9);
624  answer.zero();
625 
626  ar3 = area / 3.0;
627  answer.at(1) = answer.at(4) = answer.at(7) = -ar3 * dpdx;
628  answer.at(2) = answer.at(5) = answer.at(8) = -ar3 * dpdy;
629 
630 
631  this->computeVectorOfPressures(VM_Total, tStep->givePreviousStep(), p);
632  this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), u);
633  dpdx = 0.0, dpdy = 0.0;
634  for ( int i = 0; i < 3; i++ ) {
635  pn1 = p.at(i + 1);
636  dpdx += b [ i ] * pn1;
637  dpdy += c [ i ] * pn1;
638  }
639 
640  usum = u.at(1) + u.at(3) + u.at(5);
641  vsum = u.at(2) + u.at(4) + u.at(6);
642  coeff = ar3 * tStep->giveTimeIncrement() / 2.0;
643  answer.at(1) -= coeff * dpdx * ( b [ 0 ] * usum + c [ 0 ] * vsum );
644  answer.at(4) -= coeff * dpdx * ( b [ 1 ] * usum + c [ 1 ] * vsum );
645  answer.at(7) -= coeff * dpdx * ( b [ 2 ] * usum + c [ 2 ] * vsum );
646 
647  answer.at(2) -= coeff * dpdy * ( b [ 0 ] * usum + c [ 0 ] * vsum );
648  answer.at(5) -= coeff * dpdy * ( b [ 1 ] * usum + c [ 1 ] * vsum );
649  answer.at(8) -= coeff * dpdy * ( b [ 2 ] * usum + c [ 2 ] * vsum );
650 }
651 
652 
653 Interface *
655 {
656  if ( interface == ZZNodalRecoveryModelInterfaceType ) {
657  return static_cast< ZZNodalRecoveryModelInterface * >(this);
658  } else if ( interface == NodalAveragingRecoveryModelInterfaceType ) {
659  return static_cast< NodalAveragingRecoveryModelInterface * >(this);
660  } else if ( interface == SPRNodalRecoveryModelInterfaceType ) {
661  return static_cast< SPRNodalRecoveryModelInterface * >(this);
662  } else if ( interface == SpatialLocalizerInterfaceType ) {
663  return static_cast< SpatialLocalizerInterface * >(this);
664  } else if ( interface == EIPrimaryFieldInterfaceType ) {
665  return static_cast< EIPrimaryFieldInterface * >(this);
666  }
667  //<RESTRICTED_SECTION>
668  else if ( interface == LEPlicElementInterfaceType ) {
669  return static_cast< LEPlicElementInterface * >(this);
670  }
671 
672  //</RESTRICTED_SECTION>
673  return NULL;
674 }
675 
676 
677 void
679 {
680  /* one should call material driver instead */
681  FloatArray u, eps(3);
682 
683  this->computeVectorOfVelocities(VM_Total, tStep, u);
684 
685  eps.at(1) = ( b [ 0 ] * u.at(1) + b [ 1 ] * u.at(3) + b [ 2 ] * u.at(5) );
686  eps.at(2) = ( c [ 0 ] * u.at(2) + c [ 1 ] * u.at(4) + c [ 2 ] * u.at(6) );
687  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) );
688  static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial()->computeDeviatoricStressVector(answer, gp, eps, tStep);
689 }
690 
691 int
693 {
694  Node *node1, *node2, *node3;
695  double x1, x2, x3, y1, y2, y3;
696 
697  node1 = giveNode(1);
698  node2 = giveNode(2);
699  node3 = giveNode(3);
700 
701  // init geometry data
702  x1 = node1->giveCoordinate(1);
703  x2 = node2->giveCoordinate(1);
704  x3 = node3->giveCoordinate(1);
705 
706  y1 = node1->giveCoordinate(2);
707  y2 = node2->giveCoordinate(2);
708  y3 = node3->giveCoordinate(2);
709 
710  this->area = 0.5 * ( x2 * y3 + x1 * y2 + y1 * x3 - x2 * y1 - x3 * y2 - x1 * y3 );
711 
712  b [ 0 ] = ( y2 - y3 ) / ( 2. * area );
713  c [ 0 ] = ( x3 - x2 ) / ( 2. * area );
714  b [ 1 ] = ( y3 - y1 ) / ( 2. * area );
715  c [ 1 ] = ( x1 - x3 ) / ( 2. * area );
716  b [ 2 ] = ( y1 - y2 ) / ( 2. * area );
717  c [ 2 ] = ( x2 - x1 ) / ( 2. * area );
718 
720 }
721 
722 double
724 {
725  FloatArray u;
726  double dt1, dt2, dt;
727  double Re = static_cast< FluidModel * >( domain->giveEngngModel() )->giveReynoldsNumber();
728 
729  this->computeVectorOfVelocities(VM_Total, tStep, u);
730 
731  double vn1 = sqrt( u.at(1) * u.at(1) + u.at(2) * u.at(2) );
732  double vn2 = sqrt( u.at(3) * u.at(3) + u.at(4) * u.at(4) );
733  double vn3 = sqrt( u.at(5) * u.at(5) + u.at(6) * u.at(6) );
734  double veln = max( vn1, max(vn2, vn3) );
735 
736  double l1 = 1.0 / ( sqrt(b [ 0 ] * b [ 0 ] + c [ 0 ] * c [ 0 ]) );
737  double l2 = 1.0 / ( sqrt(b [ 1 ] * b [ 1 ] + c [ 1 ] * c [ 1 ]) );
738  double l3 = 1.0 / ( sqrt(b [ 2 ] * b [ 2 ] + c [ 2 ] * c [ 2 ]) );
739 
740  double ln = min( l1, min(l2, l3) );
741 
742  // viscous limit
743  dt2 = 0.5 * ln * ln * Re;
744  if ( veln != 0.0 ) {
745  dt1 = ln / veln;
746  dt = dt1 * dt2 / ( dt1 + dt2 );
747  } else {
748  dt = dt2;
749  }
750 
751  return dt;
752 }
753 
754 //<RESTRICTED_SECTION>
755 double
756 TR1_2D_CBS :: computeLEPLICVolumeFraction(const FloatArray &n, const double p, LEPlic *matInterface, bool updFlag)
757 {
758  Polygon pg;
759  double answer, volume = computeMyVolume(matInterface, updFlag);
760  this->formVolumeInterfacePoly(pg, matInterface, n, p, updFlag);
761  answer = fabs(pg.computeVolume() / volume);
762  if ( answer > 1.000000001 ) {
763  return 1.0;
764  } else {
765  return answer;
766  }
767 }
768 
769 void
771  const FloatArray &normal, const double p, bool updFlag)
772 {
773  double x, y;
774  Vertex v;
775 
776  matvolpoly.clear();
777 
778  if ( this->vof <= TRSUPG_ZERO_VOF ) {
779  return;
780  } else if ( this->vof >= ( 1 - TRSUPG_ZERO_VOF ) ) {
781  for ( int i = 1; i <= 3; i++ ) {
782  if ( updFlag ) {
783  x = matInterface->giveUpdatedXCoordinate( this->giveNode(i)->giveNumber() );
784  y = matInterface->giveUpdatedYCoordinate( this->giveNode(i)->giveNumber() );
785  } else {
786  x = this->giveNode(i)->giveCoordinate(1);
787  y = this->giveNode(i)->giveCoordinate(2);
788  }
789 
790  v.setCoords(x, y);
791  matvolpoly.addVertex(v);
792  }
793 
794  return;
795  }
796 
797  this->formVolumeInterfacePoly(matvolpoly, matInterface, normal, p, updFlag);
798 }
799 
800 
801 void
803  const FloatArray &normal, const double p, bool updFlag)
804 {
805  int next;
806  bool nodeIn [ 3 ];
807  double nx = normal.at(1), ny = normal.at(2), x, y;
808  double tx, ty;
809  Vertex v;
810 
811  matvolpoly.clear();
812 
813  for ( int i = 1; i <= 3; i++ ) {
814  if ( updFlag ) {
815  x = matInterface->giveUpdatedXCoordinate( this->giveNode(i)->giveNumber() );
816  y = matInterface->giveUpdatedYCoordinate( this->giveNode(i)->giveNumber() );
817  } else {
818  x = this->giveNode(i)->giveCoordinate(1);
819  y = this->giveNode(i)->giveCoordinate(2);
820  }
821 
822  if ( ( nx * x + ny * y + p ) >= 0. ) {
823  nodeIn [ i - 1 ] = true;
824  } else {
825  nodeIn [ i - 1 ] = false;
826  }
827  }
828 
829  if ( nodeIn [ 0 ] && nodeIn [ 1 ] && nodeIn [ 2 ] ) { // all nodes inside
830  for ( int i = 1; i <= 3; i++ ) {
831  if ( updFlag ) {
832  x = matInterface->giveUpdatedXCoordinate( this->giveNode(i)->giveNumber() );
833  y = matInterface->giveUpdatedYCoordinate( this->giveNode(i)->giveNumber() );
834  } else {
835  x = this->giveNode(i)->giveCoordinate(1);
836  y = this->giveNode(i)->giveCoordinate(2);
837  }
838 
839  v.setCoords(x, y);
840  matvolpoly.addVertex(v);
841  }
842 
843  return;
844  } else if ( !( nodeIn [ 0 ] || nodeIn [ 1 ] || nodeIn [ 2 ] ) ) { // all nodes outside
845  return;
846  } else {
847  for ( int i = 1; i <= 3; i++ ) {
848  next = i < 3 ? i + 1 : 1;
849  if ( nodeIn [ i - 1 ] ) {
850  if ( updFlag ) {
851  v.setCoords( matInterface->giveUpdatedXCoordinate( this->giveNode(i)->giveNumber() ),
852  matInterface->giveUpdatedYCoordinate( this->giveNode(i)->giveNumber() ) );
853  } else {
854  v.setCoords( this->giveNode(i)->giveCoordinate(1),
855  this->giveNode(i)->giveCoordinate(2) );
856  }
857 
858  matvolpoly.addVertex(v);
859  }
860 
861  if ( nodeIn [ next - 1 ] ^ nodeIn [ i - 1 ] ) {
862  // compute intersection with (i,next) edge
863  if ( updFlag ) {
864  x = matInterface->giveUpdatedXCoordinate( this->giveNode(i)->giveNumber() );
865  y = matInterface->giveUpdatedYCoordinate( this->giveNode(i)->giveNumber() );
866  tx = matInterface->giveUpdatedXCoordinate( this->giveNode(next)->giveNumber() ) - x;
867  ty = matInterface->giveUpdatedYCoordinate( this->giveNode(next)->giveNumber() ) - y;
868  } else {
869  x = this->giveNode(i)->giveCoordinate(1);
870  y = this->giveNode(i)->giveCoordinate(2);
871  tx = this->giveNode(next)->giveCoordinate(1) - x;
872  ty = this->giveNode(next)->giveCoordinate(2) - y;
873  }
874 
875  double s, sd = nx * tx + ny * ty;
876  if ( fabs(sd) > 1.e-10 ) {
877  s = ( -p - ( nx * x + ny * y ) ) / sd;
878  v.setCoords(x + tx * s, y + ty * s);
879  matvolpoly.addVertex(v);
880  } else {
881  // pathological case - lines are parallel
882  if ( nodeIn [ i - 1 ] ) {
883  if ( updFlag ) {
884  v.setCoords( matInterface->giveUpdatedXCoordinate( this->giveNode(next)->giveNumber() ),
885  matInterface->giveUpdatedYCoordinate( this->giveNode(next)->giveNumber() ) );
886  } else {
887  v.setCoords( this->giveNode(next)->giveCoordinate(1), this->giveNode(next)->giveCoordinate(2) );
888  }
889 
890  matvolpoly.addVertex(v);
891  } else {
892  v.setCoords(x, y);
893  matvolpoly.addVertex(v);
894  if ( updFlag ) {
895  v.setCoords( matInterface->giveUpdatedXCoordinate( this->giveNode(next)->giveNumber() ),
896  matInterface->giveUpdatedYCoordinate( this->giveNode(next)->giveNumber() ) );
897  } else {
898  v.setCoords( this->giveNode(next)->giveCoordinate(1), this->giveNode(next)->giveCoordinate(2) );
899  }
900 
901  matvolpoly.addVertex(v);
902  }
903  }
904  }
905  } // end loop over elem nodes
906  }
907 }
908 
909 
910 double
911 TR1_2D_CBS :: truncateMatVolume(const Polygon &matvolpoly, double &volume)
912 {
913  Polygon me, clip;
914  Graph g;
915 
916  this->formMyVolumePoly(me, NULL, false);
917  g.clip(clip, me, matvolpoly);
918 #ifdef __OOFEG
919  EASValsSetColor( gc [ 0 ].getActiveCrackColor() );
920  clip.draw(gc [ OOFEG_DEBUG_LAYER ], true);
921  //EVFastRedraw(myview);
922 #endif
923  volume = clip.computeVolume();
924  return volume / area;
925 }
926 
927 void
928 TR1_2D_CBS :: formMyVolumePoly(Polygon &me, LEPlic *matInterface, bool updFlag)
929 {
930  double x, y;
931  Vertex v;
932 
933  me.clear();
934 
935  for ( int i = 1; i <= 3; i++ ) {
936  if ( updFlag ) {
937  x = matInterface->giveUpdatedXCoordinate( this->giveNode(i)->giveNumber() );
938  y = matInterface->giveUpdatedYCoordinate( this->giveNode(i)->giveNumber() );
939  } else {
940  x = this->giveNode(i)->giveCoordinate(1);
941  y = this->giveNode(i)->giveCoordinate(2);
942  }
943 
944  v.setCoords(x, y);
945  me.addVertex(v);
946  }
947 }
948 
949 
950 double
951 TR1_2D_CBS :: computeMyVolume(LEPlic *matInterface, bool updFlag)
952 {
953  double x1, x2, x3, y1, y2, y3;
954  if ( updFlag ) {
955  x1 = matInterface->giveUpdatedXCoordinate( this->giveNode(1)->giveNumber() );
956  x2 = matInterface->giveUpdatedXCoordinate( this->giveNode(2)->giveNumber() );
957  x3 = matInterface->giveUpdatedXCoordinate( this->giveNode(3)->giveNumber() );
958 
959  y1 = matInterface->giveUpdatedYCoordinate( this->giveNode(1)->giveNumber() );
960  y2 = matInterface->giveUpdatedYCoordinate( this->giveNode(2)->giveNumber() );
961  y3 = matInterface->giveUpdatedYCoordinate( this->giveNode(3)->giveNumber() );
962  return 0.5 * ( x2 * y3 + x1 * y2 + y1 * x3 - x2 * y1 - x3 * y2 - x1 * y3 );
963  } else {
964  return area;
965  }
966 }
967 
968 void
969 TR1_2D_CBS :: giveElementCenter(LEPlic *mat_interface, FloatArray &center, bool upd)
970 {
971  FloatArray coords;
972 
973  center.resize(2);
974  center.zero();
975  if ( upd ) {
976  for ( int i = 1; i <= 3; i++ ) {
977  mat_interface->giveUpdatedCoordinate( coords, giveNode(i)->giveNumber() );
978  center.add(coords);
979  }
980  } else {
981  for ( int i = 1; i <= 3; i++ ) {
982  center.at(1) += this->giveNode(i)->giveCoordinate(1);
983  center.at(2) += this->giveNode(i)->giveCoordinate(2);
984  }
985  }
986 
987  center.times(1. / 3.);
988 }
989 
990 //</RESTRICTED_SECTION>
991 
992 
993 int
995  const FloatArray &coords, IntArray &dofId, ValueModeType mode,
996  TimeStep *tStep)
997 {
998 #if 0
999  int es = dofId.giveSize();
1001  this->computeVectorOf(pf, dofId, mode, tStep, elemvector, true);
1002  answer.resize( es );
1003  answer.zero();
1004  for ( int i = 1; i <= es; i++ ) {
1005  for ( int j = 1; j <= lc.giveSize(); j++ ) {
1006  answer.at(i) += lc.at(j) * elemvector.at( es * ( j - 1 ) + i );
1007  }
1008  }
1009 #endif
1010  int indx, es;
1011  double sum;
1012  FloatArray elemvector, f, lc;
1013  //FloatMatrix n;
1014  IntArray elemdofs;
1015  // determine element dof ids
1016  this->giveElementDofIDMask(elemdofs);
1017  es = elemdofs.giveSize();
1018  // first evaluate element unknown vector
1019  this->computeVectorOf(pf, elemdofs, mode, tStep, elemvector);
1020  // determine corresponding local coordinates
1021  if ( this->computeLocalCoordinates(lc, coords) ) {
1022  // compute interpolation matrix
1023  // this->computeNmatrixAt(n, &lc);
1024  // compute answer
1025  answer.resize( dofId.giveSize() );
1026  for ( int i = 1; i <= dofId.giveSize(); i++ ) {
1027  if ( ( indx = elemdofs.findFirstIndexOf( dofId.at(i) ) ) ) {
1028  sum = 0.0;
1029  for ( int j = 1; j <= 3; j++ ) {
1030  sum += lc.at(j) * elemvector.at(es * ( j - 1 ) + indx);
1031  }
1032 
1033  answer.at(i) = sum;
1034  } else {
1035  //_error("EIPrimaryFieldI_evaluateFieldVectorAt: unknown dof id encountered");
1036  answer.at(i) = 0.0;
1037  }
1038  }
1039 
1040  return 0; // ok
1041  } else {
1042  OOFEM_ERROR("target point not in receiver volume");
1043  return 1; // fail
1044  }
1045 }
1046 
1047 
1048 void
1050 {
1052  //<RESTRICTED_SECTION>
1054  //</RESTRICTED_SECTION>
1055 }
1056 
1057 int
1059 {
1060  if ( type == IST_VOFFraction ) {
1061  answer.resize(1);
1062  answer.at(1) = this->giveTempVolumeFraction();
1063  return 1;
1064  } else {
1065  return CBSElement :: giveIPValue(answer, gp, type, tStep);
1066  }
1067 }
1068 
1069 void
1071  InternalStateType type, TimeStep *tStep)
1072 {
1073  GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
1074  this->giveIPValue(answer, gp, type, tStep);
1075 }
1076 
1077 void
1079 {
1080  pap.resize(3);
1081  pap.at(1) = this->giveNode(1)->giveNumber();
1082  pap.at(2) = this->giveNode(2)->giveNumber();
1083  pap.at(3) = this->giveNode(3)->giveNumber();
1084 }
1085 
1086 void
1088 {
1089  answer.resize(1);
1090  if ( ( pap == this->giveNode(1)->giveNumber() ) ||
1091  ( pap == this->giveNode(2)->giveNumber() ) ||
1092  ( pap == this->giveNode(3)->giveNumber() ) ) {
1093  answer.at(1) = pap;
1094  } else {
1095  OOFEM_ERROR("node unknown");
1096  }
1097 }
1098 
1099 int
1101 { return 1; }
1102 
1103 
1106 {
1107  return SPRPatchType_2dxy;
1108 }
1109 
1110 
1111 
1112 void
1114 // Performs end-of-step operations.
1115 {
1116  CBSElement :: printOutputAt(file, tStep);
1117  //<RESTRICTED_SECTION>
1118  double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity( integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
1119  fprintf(file, "VOF %e, density %e\n\n", this->giveVolumeFraction(), rho);
1120  //</RESTRICTED_SECTION>
1121 }
1122 
1123 
1124 
1126 //
1127 // saves full element context (saves state variables, that completely describe
1128 // current state)
1129 //
1130 {
1131  contextIOResultType iores;
1132 
1133  if ( ( iores = CBSElement :: saveContext(stream, mode, obj) ) != CIO_OK ) {
1134  THROW_CIOERR(iores);
1135  }
1136 
1137  //<RESTRICTED_SECTION>
1138  if ( ( iores = LEPlicElementInterface :: saveContext(stream, mode, obj) ) != CIO_OK ) {
1139  THROW_CIOERR(iores);
1140  }
1141 
1142  //</RESTRICTED_SECTION>
1143 
1144  return CIO_OK;
1145 }
1146 
1147 
1148 
1150 //
1151 // restores full element context (saves state variables, that completely describe
1152 // current state)
1153 //
1154 {
1155  contextIOResultType iores;
1156 
1157  if ( ( iores = CBSElement :: restoreContext(stream, mode, obj) ) != CIO_OK ) {
1158  THROW_CIOERR(iores);
1159  }
1160 
1161  //<RESTRICTED_SECTION>
1162  if ( ( iores = LEPlicElementInterface :: restoreContext(stream, mode, obj) ) != CIO_OK ) {
1163  THROW_CIOERR(iores);
1164  }
1165 
1166  //</RESTRICTED_SECTION>
1167 
1168 
1169  return CIO_OK;
1170 }
1171 
1172 
1173 #ifdef __OOFEG
1174 int
1176  int node, TimeStep *tStep)
1177 {
1178  //<RESTRICTED_SECTION>
1179  if ( type == IST_VOFFraction ) {
1180  answer.resize(1);
1181  answer.at(1) = this->giveTempVolumeFraction();
1182  return 1;
1183  } else
1184  //</RESTRICTED_SECTION>
1185  if ( type == IST_Density ) {
1186  answer.resize(1);
1187  answer.at(1) = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity( integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
1188  return 1;
1189  } else {
1190  return CBSElement :: giveInternalStateAtNode(answer, type, mode, node, tStep);
1191  }
1192 }
1193 
1194 void
1196 {
1197  WCRec p [ 3 ];
1198  GraphicObj *go;
1199 
1200  if ( !gc.testElementGraphicActivity(this) ) {
1201  return;
1202  }
1203 
1204  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
1205  EASValsSetColor( gc.getElementColor() );
1206  EASValsSetEdgeColor( gc.getElementEdgeColor() );
1207  EASValsSetEdgeFlag(true);
1208  EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
1209  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
1210  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
1211  p [ 0 ].z = 0.;
1212  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
1213  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
1214  p [ 1 ].z = 0.;
1215  p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveCoordinate(1);
1216  p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveCoordinate(2);
1217  p [ 2 ].z = 0.;
1218 
1219  go = CreateTriangle3D(p);
1220  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
1221  EGAttachObject(go, ( EObjectP ) this);
1222  EMAddGraphicsToModel(ESIModel(), go);
1223 }
1224 
1226 {
1227  int i, indx, result = 0;
1228  WCRec p [ 3 ];
1229  GraphicObj *tr;
1230  FloatArray v1, v2, v3;
1231  double s [ 3 ];
1232 
1233  if ( !gc.testElementGraphicActivity(this) ) {
1234  return;
1235  }
1236 
1237  if ( gc.giveIntVarMode() == ISM_recovered ) {
1238  result += this->giveInternalStateAtNode(v1, gc.giveIntVarType(), gc.giveIntVarMode(), 1, tStep);
1239  result += this->giveInternalStateAtNode(v2, gc.giveIntVarType(), gc.giveIntVarMode(), 2, tStep);
1240  result += this->giveInternalStateAtNode(v3, gc.giveIntVarType(), gc.giveIntVarMode(), 3, tStep);
1241  } else if ( gc.giveIntVarMode() == ISM_local ) {
1242  GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
1243  result += giveIPValue(v1, gp, gc.giveIntVarType(), tStep);
1244  v2 = v1;
1245  v3 = v1;
1246  result *= 3;
1247  }
1248 
1249  if ( result != 3 ) {
1250  return;
1251  }
1252 
1253  indx = gc.giveIntVarIndx();
1254 
1255  s [ 0 ] = v1.at(indx);
1256  s [ 1 ] = v2.at(indx);
1257  s [ 2 ] = v3.at(indx);
1258 
1259  EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
1260 
1261  if ( gc.getScalarAlgo() == SA_ISO_SURF ) {
1262  for ( i = 0; i < 3; i++ ) {
1263  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
1264  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
1265  p [ i ].z = 0.;
1266  }
1267 
1268  //EASValsSetColor(gc.getYieldPlotColor(ratio));
1269  gc.updateFringeTableMinMax(s, 3);
1270  tr = CreateTriangleWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ]);
1271  EGWithMaskChangeAttributes(LAYER_MASK, tr);
1272  EMAddGraphicsToModel(ESIModel(), tr);
1273  } else if ( ( gc.getScalarAlgo() == SA_ZPROFILE ) || ( gc.getScalarAlgo() == SA_COLORZPROFILE ) ) {
1274  double landScale = gc.getLandScale();
1275 
1276  for ( i = 0; i < 3; i++ ) {
1277  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
1278  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
1279  p [ i ].z = s [ i ] * landScale;
1280  }
1281 
1282  if ( gc.getScalarAlgo() == SA_ZPROFILE ) {
1283  EASValsSetColor( gc.getDeformedElementColor() );
1284  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
1285  EASValsSetFillStyle(FILL_SOLID);
1286  tr = CreateTriangle3D(p);
1287  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | FILL_MASK | LAYER_MASK, tr);
1288  } else {
1289  gc.updateFringeTableMinMax(s, 3);
1290  EASValsSetFillStyle(FILL_SOLID);
1291  tr = CreateTriangleWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ]);
1292  EGWithMaskChangeAttributes(FILL_MASK | LAYER_MASK, tr);
1293  }
1294 
1295  EMAddGraphicsToModel(ESIModel(), tr);
1296  }
1297 }
1298 
1299 
1300 
1301 #endif
1302 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
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...
Definition: cbselement.C:168
virtual void computeDiffusionTermsI(FloatArray &answer, TimeStep *)
Calculates contribution from diffusion terms for (*) velocities.
Definition: tr1_2d_cbs.C:274
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: tr1_2d_cbs.C:103
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.
void setField(int item, InputFieldType id)
virtual int checkConsistency()
Used to check consistency and initialize some element geometry data (area,b,c)
Definition: tr1_2d_cbs.C:692
virtual ~TR1_2D_CBS()
Definition: tr1_2d_cbs.C:83
The element interface required by ZZNodalRecoveryModel.
contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores context of receiver into given stream.
Definition: leplic.C:86
virtual int computeNumberOfDofs()
Computes or simply returns total number of element&#39;s local DOFs.
Definition: tr1_2d_cbs.C:90
void giveUpdatedCoordinate(FloatArray &answer, int num)
Returns updated nodal positions.
Definition: leplic.h:183
IntArray * giveBoundaryLoadArray()
Returns array containing load numbers of boundary loads acting on element.
Definition: element.C:381
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 void formMyVolumePoly(Polygon &myPoly, LEPlic *mat_interface, bool updFlag)
Assembles receiver volume.
Definition: tr1_2d_cbs.C:928
Element interface for LEPlic class representing Lagrangian-Eulerian (moving) material interface...
Definition: leplic.h:58
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
Definition: cbselement.C:72
Class and object Domain.
Definition: domain.h:115
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Returns local vector of unknowns.
Definition: element.C:86
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: cbselement.C:58
Fluid cross-section.
ScalarAlgorithmType getScalarAlgo()
#define FMElement_PrescribedTractionBC
Definition: fmelement.h:44
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition: femcmpnn.h:82
virtual int SPRNodalRecoveryMI_giveNumberOfIP()
Definition: tr1_2d_cbs.C:1100
double computeVolume() const
Definition: geotoolbox.C:72
The element interface required by ZZNodalRecoveryModel.
IntArray boundarySides
Array of boundary sides.
Definition: cbselement.h:60
Abstract class representing field of primary variables (those, which are unknown and are typically as...
Definition: primaryfield.h:104
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
virtual double truncateMatVolume(const Polygon &matvolpoly, double &volume)
Truncates given material polygon to receiver.
Definition: tr1_2d_cbs.C:911
bcGeomType
Type representing the geometric character of loading.
Definition: bcgeomtype.h:40
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
int max(int i, int j)
Returns bigger value form two given decimals.
Definition: mathfem.h:71
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
IntArray boundaryCodes
Boundary sides codes.
Definition: cbselement.h:62
#define OOFEG_RAW_GEOMETRY_LAYER
virtual void computePressureLhs(FloatMatrix &answer, TimeStep *tStep)
Calculates the pressure LHS.
Definition: tr1_2d_cbs.C:589
void setPermanentVolumeFraction(double v)
Definition: leplic.h:106
Base class for fluid problems.
Definition: fluidmodel.h:46
virtual void giveElementCenter(LEPlic *mat_interface, FloatArray &center, bool upd)
Computes the receiver center (in updated Lagrangian configuration).
Definition: tr1_2d_cbs.C:969
virtual SPRPatchType SPRNodalRecoveryMI_givePatchType()
Definition: tr1_2d_cbs.C:1105
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
double giveUpdatedYCoordinate(int num)
Definition: leplic.h:190
EngngModel * giveEngngModel()
Returns engineering model to which receiver is associated.
Definition: domain.C:433
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: element.C:779
virtual void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
Definition: tr1_2d_cbs.C:1078
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: tr1_2d_cbs.C:1049
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: tr1_2d_cbs.C:138
virtual void computeComponentArrayAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Computes boundary condition value - its components values at given time.
Definition: load.C:82
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_cbs.C:678
virtual double giveCoordinate(int i)
Definition: node.C:82
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: tr1_2d_cbs.C:1125
virtual void computeNumberOfNodalPrescribedTractionPressureContributions(FloatArray &answer, TimeStep *tStep)
Computes number of edges/sides with prescribed traction contributing to node with prescribed pressure...
Definition: tr1_2d_cbs.C:538
#define FMElement_PrescribedPressureBC
Definition: fmelement.h:47
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
Definition: tr1_2d_cbs.C:126
virtual int EIPrimaryFieldI_evaluateFieldVectorAt(FloatArray &answer, PrimaryField &pf, const FloatArray &coords, IntArray &dofId, ValueModeType mode, TimeStep *tStep)
Evaluates the value of field at given point of interest (should be located inside receiver&#39;s volume) ...
Definition: tr1_2d_cbs.C:994
virtual FEInterpolation * giveInterpolation() const
Definition: tr1_2d_cbs.C:87
#define THROW_CIOERR(e)
Definition: contextioerr.h:61
Distributed body load.
Definition: bcgeomtype.h:43
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: element.C:970
void computeVectorOfPrescribed(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Returns local vector of prescribed unknowns.
Definition: element.C:181
double giveTimeIncrement()
Returns solution step associated time increment.
Definition: timestep.h:150
virtual void computeDensityRhsVelocityTerms(FloatArray &answer, TimeStep *tStep)
Computes velocity terms on RHS for density equation.
Definition: tr1_2d_cbs.C:380
void updateYourself(TimeStep *tStep)
Definition: leplic.h:122
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
#define FMElement_PrescribedUnBC
Definition: fmelement.h:45
Element interface class.
Definition: primaryfield.h:58
virtual void computeVectorOfPressures(ValueModeType mode, TimeStep *tStep, FloatArray &pressures)
Definition: fmelement.C:55
InternalStateType giveIntVarType()
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Definition: element.C:759
#define _IFT_Tr1CBS_vof
Definition: tr1_2d_cbs.h:51
This class implements a transport material status information.
Abstract base class representing a boundary load (force, momentum, ...) that acts directly on a bound...
Definition: boundaryload.h:110
TR1_2D_CBS(int n, Domain *aDomain)
Definition: tr1_2d_cbs.C:71
virtual void computeVectorOfVelocities(ValueModeType mode, TimeStep *tStep, FloatArray &velocities)
Definition: fmelement.C:48
virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
Definition: tr1_2d_cbs.C:1087
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
#define OOFEG_RAW_GEOMETRY_WIDTH
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: tr1_2d_cbs.C:1058
void clip(Polygon &result, const Polygon &a, const Polygon &b)
Definition: geotoolbox.C:575
virtual void computeConsistentMassMtrx(FloatMatrix &answer, TimeStep *)
Calculates consistent mass matrix.
Definition: tr1_2d_cbs.C:149
IntArray * giveBodyLoadArray()
Returns array containing load numbers of loads acting on element.
Definition: element.C:372
virtual double computeMyVolume(LEPlic *matInterface, bool updFlag)
Computes the volume of receiver.
Definition: tr1_2d_cbs.C:951
IntArray bodyLoadArray
Array containing indexes of loads (body loads and boundary loads are kept separately), that apply on receiver.
Definition: element.h:160
#define _IFT_Tr1CBS_pvof
Definition: tr1_2d_cbs.h:52
Class representing the special graph constructed from two polygons that is used to perform boolean op...
Definition: geotoolbox.h:191
virtual void giveElementDofIDMask(IntArray &answer) const
Returns element dof mask for node.
Definition: element.h:498
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Definition: tr1_2d_cbs.C:1225
TimeStep * givePreviousStep()
Returns pointer to previous solution step.
Definition: timestep.C:114
virtual Interface * giveInterface(InterfaceType)
Interface requesting service.
Definition: tr1_2d_cbs.C:654
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
void resize(int n)
Checks size of receiver towards requested bounds.
Definition: intarray.C:124
virtual void formMaterialVolumePoly(Polygon &matvolpoly, LEPlic *matInterface, const FloatArray &normal, const double p, bool updFlag)
Assembles the true element material polygon (takes receiver vof into accout).
Definition: tr1_2d_cbs.C:770
This class represents CBS algorithm for solving incompressible Navier-Stokes equations.
Definition: cbs.h:171
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
Definition: crosssection.C:54
InternalStateMode giveIntVarMode()
virtual double computeLEPLICVolumeFraction(const FloatArray &n, const double p, LEPlic *matInterface, bool updFlag)
Computes corresponding volume fraction to given interface position.
Definition: tr1_2d_cbs.C:756
Abstract base class representing Lagrangian-Eulerian (moving) material interfaces.
Definition: leplic.h:153
virtual void computePrescribedTractionPressure(FloatArray &answer, TimeStep *tStep)
Computes prescribed pressure due to applied tractions.
Definition: tr1_2d_cbs.C:464
Class representing vector of real numbers.
Definition: floatarray.h:82
Class representing 2D polygon.
Definition: geotoolbox.h:91
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Definition: tr1_2d_cbs.C:1113
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
Definition: tr1_2d_cbs.C:1070
virtual void computeDensityRhsPressureTerms(FloatArray &answer, TimeStep *tStep)
Computes pressure terms on RHS for density equation.
Definition: tr1_2d_cbs.C:566
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
void setCoords(double x, double y)
Definition: geotoolbox.h:77
IntegrationPointStatus * giveMaterialStatus()
Returns reference to associated material status (NULL if not defined).
Definition: gausspoint.h:205
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: tr1_2d_cbs.C:1149
virtual void computeDiagonalMassMtrx(FloatArray &answer, TimeStep *)
Calculates diagonal mass matrix.
Definition: tr1_2d_cbs.C:173
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
double giveUpdatedXCoordinate(int num)
Definition: leplic.h:189
virtual double computeCriticalTimeStep(TimeStep *tStep)
Calculates critical time step.
Definition: tr1_2d_cbs.C:723
Class Interface.
Definition: interface.h:82
#define FMElement_PrescribedUsBC
Definition: fmelement.h:46
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
virtual bcGeomType giveBCGeoType() const
Returns geometry character of boundary condition.
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: element.C:885
This abstract class represent a general base element class for fluid dynamic problems solved using CB...
Definition: cbselement.h:56
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
Class representing the a dynamic Input Record.
void times(double s)
Multiplies receiver with scalar.
Definition: floatarray.C:818
The spatial localizer element interface associated to spatial localizer.
virtual bool computeLocalCoordinates(FloatArray &answer, const FloatArray &gcoords)
Computes the element local coordinates from given global coordinates.
Definition: element.C:1222
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
contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores context of receiver from given stream.
Definition: leplic.C:107
Class representing vertex.
Definition: geotoolbox.h:60
double vof
Volume fraction of reference fluid in element.
Definition: leplic.h:63
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: tr1_2d_cbs.C:1195
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
int min(int i, int j)
Returns smaller value from two given decimals.
Definition: mathfem.h:59
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
Definition: tr1_2d_cbs.C:96
virtual int checkConsistency()
Performs consistency check.
Definition: cbselement.C:140
void updateFringeTableMinMax(double *s, int size)
virtual void computeConvectionTermsI(FloatArray &answer, TimeStep *)
Calculates convection component for (*) velocities.
Definition: tr1_2d_cbs.C:187
#define OOFEG_DEBUG_LAYER
Load is base abstract class for all loads.
Definition: load.h:61
double p
Line constant of line segment representing interface.
Definition: leplic.h:65
#define IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:78
GraphicObj * draw(oofegGraphicContext &, bool filled, int layer=OOFEG_DEBUG_LAYER)
Definition: geotoolbox.C:152
int giveSize() const
Definition: intarray.h:203
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
the oofem namespace is to define a context or scope in which all oofem names are defined.
void addVertex(Vertex v)
Definition: geotoolbox.h:96
virtual bcValType giveBCValType() const
Returns receiver load type.
Class implementing node in finite element mesh.
Definition: node.h:87
#define TRSUPG_ZERO_VOF
Definition: tr1_2d_cbs.C:65
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 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
Load * giveLoad(int n)
Service for accessing particular domain load.
Definition: domain.C:222
#define OOFEG_VARPLOT_PATTERN_LAYER
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...
Definition: tr1_2d_cbs.C:1175
Class representing integration point in finite element program.
Definition: gausspoint.h:93
IntArray boundaryLoadArray
Definition: element.h:160
FloatArray normal
Interface segment normal.
Definition: leplic.h:67
virtual void formVolumeInterfacePoly(Polygon &matvolpoly, LEPlic *matInterface, const FloatArray &normal, const double p, bool updFlag)
Assembles receiver material polygon based solely on given interface line.
Definition: tr1_2d_cbs.C:802
Class representing solution step.
Definition: timestep.h:80
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
InternalStateMode
Determines the mode of internal variable.
void add(const FloatArray &src)
Adds array src to receiver.
Definition: floatarray.C:156
int findFirstIndexOf(int value) const
Finds index of first occurrence of given value in array.
Definition: intarray.C:331
static FEI2dTrLin interp
Definition: tr1_2d_cbs.h:77
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631
virtual void computeCorrectionRhs(FloatArray &answer, TimeStep *tStep)
Calculates the RHS of velocity correction step.
Definition: tr1_2d_cbs.C:606

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