OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
tr1_2d_supg2_axi.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_supg2_axi.h"
36 #include "fluidmodel.h"
37 #include "node.h"
38 #include "material.h"
39 #include "gausspoint.h"
40 #include "gaussintegrationrule.h"
41 #include "floatmatrix.h"
42 #include "floatarray.h"
43 #include "intarray.h"
44 #include "domain.h"
45 #include "mathfem.h"
46 #include "engngm.h"
47 #include "fluiddynamicmaterial.h"
48 #include "fluidcrosssection.h"
49 #include "load.h"
50 #include "timestep.h"
51 #include "boundaryload.h"
52 #include "fei2dtrlin.h"
53 #include "fei2dquadlin.h"
54 #include "geotoolbox.h"
55 #include "crosssection.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_SUPG2_AXI);
68 
69 //#define TR1_2D_SUPG2_AXI_DEBUG
70 
72  TR1_2D_SUPG(n, aDomain)
73  // Constructor.
74 {
75  numberOfDofMans = 3;
76 }
77 
79 // Destructor
80 { }
81 
82 
85 {
86  IRResultType result; // Required by IR_GIVE_FIELD macro
87 
88  this->vof = 0.0;
90  if ( vof > 0.0 ) {
92  this->temp_vof = this->vof;
93  } else {
94  this->vof = 0.0;
96  this->temp_vof = this->vof;
97  }
98 
99  this->mat [ 0 ] = this->mat [ 1 ] = this->material;
102  this->material = this->mat [ 0 ];
103 
104  result = SUPGElement :: initializeFrom(ir);
105  if ( result != IRRT_OK ) {
106  return result;
107  }
108  this->initGeometry();
109  return IRRT_OK;
110 }
111 
112 
113 void
115 {
117  if ( this->permanentVofFlag ) {
118  input.setField(this->vof, _IFT_Tr1SUPG_pvof);
119  } else {
120  input.setField(this->vof, _IFT_Tr1SUPG_vof);
121  }
122 
123  input.setField(this->mat [ 0 ], _IFT_Tr1SUPG2_mat0);
124  input.setField(this->mat [ 1 ], _IFT_Tr1SUPG2_mat1);
125 }
126 
127 void
129 // Sets up the array containing the four Gauss points of the receiver.
130 {
131  if ( integrationRulesArray.size() == 0 ) {
132  integrationRulesArray.resize( 2 );
133  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 3, true) );
134  integrationRulesArray [ 1 ].reset( new GaussIntegrationRule(2, this, 1, 3, true) );
135  }
136 }
137 
138 
139 /*
140  * Integration template
141  * // loop over each fluid
142  * for (ifluid = 0; ifluid< 2; ifluid++) {
143  * for (GaussPoint *gp: *integrationRulesArray [ ifluid ] ) {
144  * mapped_gp = sub_mapped_IPRule[ifluid]->getIntegrationPoint(ip) ;
145  * this->computeNMtrx (n, mapped_gp);
146  * dV = this->computeVolumeAroundID(gp,id[ifluid], vcoords[ifluid]) ;
147  * // compute integral here
148  * }
149  * }
150  */
151 
152 
153 
154 void
156 {
157  answer.resize(6, 6);
158  answer.zero();
159  FloatArray n, un;
160 
161  double _val, u, v;
162 
163  this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), un);
164 
165  // consistent mass
166  // loop over each fluid
167  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
168  for ( GaussPoint *gp: *integrationRulesArray [ ifluid ] ) {
169  double rho = this->_giveMaterial(ifluid)->give('d', gp);
170  double dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
171  this->computeNMtrx(n, gp);
172 
173  for ( int i = 1; i <= 3; i++ ) {
174  for ( int j = 1; j <= 3; j++ ) {
175  /* consistent mass */
176  _val = n.at(i) * n.at(j) * rho * dV;
177  answer.at(2 * i - 1, 2 * j - 1) += _val;
178  answer.at(2 * i, 2 * j) += _val;
179 
180  /* SUPG stabilization term */
181  u = n.at(1) * un.at(1) + n.at(2) * un.at(3) + n.at(3) * un.at(5);
182  v = n.at(1) * un.at(2) + n.at(2) * un.at(4) + n.at(3) * un.at(6);
183  answer.at(2 * i - 1, 2 * j - 1) += ( u * b [ i - 1 ] + v * c [ j - 1 ] ) * n.at(j) * rho * t_supg * dV;
184  answer.at(2 * i, 2 * j) += ( u * b [ i - 1 ] + v * c [ j - 1 ] ) * n.at(j) * rho * t_supg * dV;
185  }
186  }
187  }
188  }
189 }
190 
191 
192 void
194 {
195  answer.resize(6);
196  answer.zero();
197 
198  FloatArray n, u, un;
199  double dudx, dudy, dvdx, dvdy, _u, _v;
200  this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), un);
201  this->computeVectorOfVelocities(VM_Total, tStep, u);
202 
203 
204  dudx = b [ 0 ] * u.at(1) + b [ 1 ] * u.at(3) + b [ 2 ] * u.at(5);
205  dudy = c [ 0 ] * u.at(1) + c [ 1 ] * u.at(3) + c [ 2 ] * u.at(5);
206  dvdx = b [ 0 ] * u.at(2) + b [ 1 ] * u.at(4) + b [ 2 ] * u.at(6);
207  dvdy = c [ 0 ] * u.at(2) + c [ 1 ] * u.at(4) + c [ 2 ] * u.at(6);
208 
209  // standard galerkin term
210  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
211  for ( GaussPoint *gp: *integrationRulesArray [ ifluid ] ) {
212  double rho = this->_giveMaterial(ifluid)->give('d', gp);
213  double dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
214  this->computeNMtrx(n, gp);
215 
216  _u = n.at(1) * un.at(1) + n.at(2) * un.at(3) + n.at(3) * un.at(5);
217  _v = n.at(1) * un.at(2) + n.at(2) * un.at(4) + n.at(3) * un.at(6);
218 
219  // standard galerkin term
220  for ( int i = 1; i <= 3; i++ ) {
221  answer.at(2 * i - 1) += rho * n.at(i) * ( _u * dudx + _v * dudy ) * dV;
222  answer.at(2 * i) += rho * n.at(i) * ( _u * dvdx + _v * dvdy ) * dV;
223  }
224 
225  // supg stabilization term
226  for ( int i = 1; i <= 3; i++ ) {
227  answer.at(2 * i - 1) += ( _u * b [ i - 1 ] + _v * c [ i - 1 ] ) * ( _u * dudx + _v * dudy ) * rho * t_supg * dV;
228  answer.at(2 * i) += ( _u * b [ i - 1 ] + _v * c [ i - 1 ] ) * ( _u * dvdx + _v * dvdy ) * rho * t_supg * dV;
229  }
230  }
231  }
232 }
233 
234 
235 void
237 {
238  answer.resize(6, 6);
239  answer.zero();
240 
241  FloatArray u, un, n;
242  this->computeVectorOfVelocities(VM_Total, tStep, u);
243  this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), un);
244  double _u, _v;
245  int w_dof_addr, u_dof_addr, dij;
246 
247  // dN(v)/dv
248  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
249  for ( GaussPoint *gp: *integrationRulesArray [ ifluid ] ) {
250  double rho = this->_giveMaterial(ifluid)->give('d', gp);
251  double dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
252  this->computeNMtrx(n, gp);
253 
254  _u = n.at(1) * un.at(1) + n.at(2) * un.at(3) + n.at(3) * un.at(5);
255  _v = n.at(1) * un.at(2) + n.at(2) * un.at(4) + n.at(3) * un.at(6);
256 
257  // dN(v)/dv
258  for ( int i = 1; i <= 2; i++ ) { // test function index
259  for ( int k = 1; k <= 3; k++ ) { // nodal val of function w
260  for ( int j = 1; j <= 2; j++ ) { // velocity vector component
261  for ( int m = 1; m <= 3; m++ ) { // nodal components
262  w_dof_addr = ( k - 1 ) * 2 + i;
263  u_dof_addr = ( m - 1 ) * 2 + j;
264  //d1j = (j==1); d2j=(j==2);
265  dij = ( i == j );
266  answer.at(w_dof_addr, u_dof_addr) += dV * rho * n.at(k) * ( dij * _u * b [ m - 1 ] + dij * _v * c [ m - 1 ] );
267  }
268  }
269  }
270  }
271 
272  // stabilization term dN_delta/du
273  for ( int i = 1; i <= 2; i++ ) { // test function index
274  for ( int k = 1; k <= 3; k++ ) { // nodal val of function w
275  for ( int j = 1; j <= 2; j++ ) { // velocity vector component
276  for ( int m = 1; m <= 3; m++ ) { // nodal components
277  w_dof_addr = ( k - 1 ) * 2 + i;
278  u_dof_addr = ( m - 1 ) * 2 + j;
279  //d1j = (j==1); d2j=(j==2);
280  dij = ( i == j );
281  answer.at(w_dof_addr, u_dof_addr) += dV * t_supg * rho *
282  ( _u * b [ k - 1 ] + _v * c [ k - 1 ] ) * ( dij * _u * b [ m - 1 ] + dij * _v * c [ m - 1 ] );
283  }
284  }
285  }
286  }
287  }
288  }
289 }
290 
291 
292 void
294 {
295  answer.resize(6);
296  answer.zero();
297  FloatArray u, un, eps, stress;
298  double Re = static_cast< FluidModel * >( domain->giveEngngModel() )->giveReynoldsNumber();
299  //double dudx,dudy,dvdx,dvdy;
300 
301  this->computeVectorOfVelocities(VM_Total, tStep, u);
302  FloatArray n;
303  double _u, _v, _r;
304  this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), un);
305  FloatMatrix _b(4, 6);
306 
307  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
308  FluidDynamicMaterial *mat = static_cast< FluidDynamicMaterial * >( this->_giveMaterial(ifluid) );
309  for ( GaussPoint *gp: *integrationRulesArray [ ifluid ] ) {
310  double dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
311 
312  this->computeBMtrx(_b, gp);
313  eps.beProductOf(_b, u);
314  mat->computeDeviatoricStressVector(stress, gp, eps, tStep);
315  answer.plusProduct(_b, stress, dV / Re);
316 
317 #if 1
318  // stabilization term k_delta
319  _r = this->computeRadiusAt(gp);
320  computeNVector(n, gp);
321  _u = n.at(1) * un.at(1) + n.at(2) * un.at(3) + n.at(3) * un.at(5);
322  _v = n.at(1) * un.at(2) + n.at(2) * un.at(4) + n.at(3) * un.at(6);
323 
324  for ( int i = 1; i <= 3; i++ ) {
325  answer.at(2 * i - 1) -= t_supg * ( _u * b [ i - 1 ] + _v * c [ i - 1 ] ) * ( stress.at(1) / _r ) * dV / Re;
326  answer.at(2 * i) -= t_supg * ( _u * b [ i - 1 ] + _v * c [ i - 1 ] ) * ( stress.at(4) / _r ) * dV / Re;
327  }
328 
329 #endif
330  }
331  }
332 }
333 
334 
335 void
337 {
338  //double dudx, dudy, dvdx, dvdy;
339  answer.resize(6, 6);
340  answer.zero();
341  FloatMatrix _db, _d, _b;
342  //FloatArray un;
343  double Re = static_cast< FluidModel * >( domain->giveEngngModel() )->giveReynoldsNumber();
344  FloatArray un, u, n, eps, stress;
345  double _u, _v, _r;
346 
347  this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), un);
348  this->computeVectorOfVelocities(VM_Total, tStep, u);
349 
350  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
351  FluidDynamicMaterial *mat = static_cast< FluidDynamicMaterial * >( this->_giveMaterial(ifluid) );
352  for ( GaussPoint *gp: *integrationRulesArray [ ifluid ] ) {
353  double dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
354 
355  this->computeBMtrx(_b, gp);
356  mat->giveDeviatoricStiffnessMatrix(_d, mode, gp, tStep);
357  _db.beProductOf(_d, _b);
358  answer.plusProductUnsym(_b, _db, dV);
359  //answer.plusProductSymmUpper (_bs,_db,dV*t_supg);
360  // }
361 
362 #if 1
363  _r = this->computeRadiusAt(gp);
364  computeNVector(n, gp);
365  eps.beProductOf(_b, u);
366  mat->computeDeviatoricStressVector(stress, gp, eps, tStep);
367  //_mu = mat->giveCharacteristicValue(MRM_Viscosity, gp, tStep);
368 
369  _u = n.at(1) * un.at(1) + n.at(2) * un.at(3) + n.at(3) * un.at(5);
370  _v = n.at(1) * un.at(2) + n.at(2) * un.at(4) + n.at(3) * un.at(6);
371 
372  for ( int i = 1; i <= 3; i++ ) {
373  for ( int j = 1; j <= 6; j++ ) {
374  //answer.at(2*i-1,j) -= t_supg*(_u*b[i-1]+_v*c[i-1])*(_d.at(1,1)*_b.at(1,j)/_r - stress.at(1)/_r/_r)*dV;
375  //answer.at(2*i,j) -= t_supg*(_u*b[i-1]+_v*c[i-1])*(_d.at(4,4)*_b.at(4,j)/_r - stress.at(4)/_r/_r)*dV;
376 
377  answer.at(2 * i - 1, j) -= t_supg * ( _u * b [ i - 1 ] + _v * c [ i - 1 ] ) * ( _db.at(1, j) / _r ) * dV;
378  answer.at(2 * i, j) -= t_supg * ( _u * b [ i - 1 ] + _v * c [ i - 1 ] ) * ( _db.at(4, j) / _r ) * dV;
379  }
380  }
381 
382 #endif
383  }
384  }
385 
386  answer.times(1. / Re);
387 }
388 
389 
390 void
392 {
393  answer.resize(6, 3);
394  answer.zero();
395  FloatArray un, n;
396  double _u, _v;
397 
398  this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), un);
399  double dV, _r;
400  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
401  for ( GaussPoint *gp: *integrationRulesArray [ ifluid ] ) {
402  dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
403  _r = this->computeRadiusAt(gp);
404  computeNVector(n, gp);
405  // G matrix
406  for ( int i = 1; i <= 3; i++ ) {
407  for ( int j = 1; j <= 3; j++ ) {
408  answer.at(2 * i - 1, j) -= b [ i - 1 ] * n.at(j) * dV;
409  answer.at(2 * i, j) -= c [ i - 1 ] * n.at(j) * dV;
410 
411  answer.at(2 * i - 1, j) -= n.at(i) * n.at(j) * dV / _r;
412  }
413  }
414 
415  // stabilization term (G_\delta mtrx)
416  _u = n.at(1) * un.at(1) + n.at(2) * un.at(3) + n.at(3) * un.at(5);
417  _v = n.at(1) * un.at(2) + n.at(2) * un.at(4) + n.at(3) * un.at(6);
418  for ( int i = 1; i <= 3; i++ ) {
419  for ( int j = 1; j <= 3; j++ ) {
420  answer.at(2 * i - 1, j) += ( _u * b [ i - 1 ] + _v * c [ i - 1 ] ) * b [ j - 1 ] * dV * t_supg;
421  answer.at(2 * i, j) += ( _u * b [ i - 1 ] + _v * c [ i - 1 ] ) * c [ j - 1 ] * dV * t_supg;
422  }
423  }
424  }
425  }
426 }
427 
428 
429 void
431 {
432  answer.resize(6, 6);
433  answer.zero();
434  //double coeff = area*t_lsic*rho;
435  double n[] = {
436  b [ 0 ], c [ 0 ], b [ 1 ], c [ 1 ], b [ 2 ], c [ 2 ]
437  };
438 
439  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
440  for ( GaussPoint *gp: *integrationRulesArray [ ifluid ] ) {
441  double rho = this->_giveMaterial(ifluid)->give('d', gp);
442  double dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
443 
444  for ( int i = 1; i <= 6; i++ ) {
445  for ( int j = 1; j <= 6; j++ ) {
446  answer.at(i, j) += dV * t_lsic * rho * n [ i - 1 ] * n [ j - 1 ];
447  }
448  }
449  }
450  }
451 }
452 
453 
454 void
456 {
457  answer.resize(3, 6);
458  answer.zero();
459 
460  FloatArray n;
461 
462  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
463  for ( auto &gp : *integrationRulesArray [ ifluid ] ) {
464  double dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
465  double _r = this->computeRadiusAt(gp);
466  computeNVector(n, gp);
467  for ( int i = 1; i <= 3; i++ ) {
468  for ( int j = 1; j <= 3; j++ ) {
469  answer.at(j, 2 * i - 1) += b [ i - 1 ] * n.at(j) * dV;
470  answer.at(j, 2 * i) += c [ i - 1 ] * n.at(j) * dV;
471 
472  answer.at(i, 1 + ( j - 1 ) * 2) += n.at(i) * n.at(j) * dV / _r;
473  }
474  }
475  }
476  }
477 }
478 
479 void
481 {
482  // N_epsilon (due to PSPG stabilization)
483  double dudx, dudy, dvdx, dvdy, _u, _v;
484  FloatArray u, un, n;
485 
486  answer.resize(3);
487  answer.zero();
488 
489  this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), un);
490  this->computeVectorOfVelocities(VM_Total, tStep, u);
491  dudx = b [ 0 ] * u.at(1) + b [ 1 ] * u.at(3) + b [ 2 ] * u.at(5);
492  dudy = c [ 0 ] * u.at(1) + c [ 1 ] * u.at(3) + c [ 2 ] * u.at(5);
493  dvdx = b [ 0 ] * u.at(2) + b [ 1 ] * u.at(4) + b [ 2 ] * u.at(6);
494  dvdy = c [ 0 ] * u.at(2) + c [ 1 ] * u.at(4) + c [ 2 ] * u.at(6);
495 
496  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
497  for ( auto &gp : *integrationRulesArray [ ifluid ] ) {
498  double dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
499  computeNVector(n, gp);
500 
501  _u = n.at(1) * un.at(1) + n.at(2) * un.at(3) + n.at(3) * un.at(5);
502  _v = n.at(1) * un.at(2) + n.at(2) * un.at(4) + n.at(3) * un.at(6);
503 
504  for ( int i = 1; i <= 3; i++ ) {
505  answer.at(i) += t_pspg * dV * ( b [ i - 1 ] * ( _u * dudx + _v * dudy ) + c [ i - 1 ] * ( _u * dvdx + _v * dvdy ) );
506  }
507  }
508  }
509 }
510 
511 
512 void
514 {
515  answer.resize(3, 6);
516  answer.zero();
517  int w_dof_addr, u_dof_addr, d1j, d2j, km1, mm1;
518  FloatArray u, un, n;
519  double _u, _v;
520 
521  this->computeVectorOfVelocities(VM_Total, tStep, u);
522  this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), un);
523 
524  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
525  for ( auto &gp : *integrationRulesArray [ ifluid ] ) {
526  double dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
527  computeNVector(n, gp);
528 
529  _u = n.at(1) * un.at(1) + n.at(2) * un.at(3) + n.at(3) * un.at(5);
530  _v = n.at(1) * un.at(2) + n.at(2) * un.at(4) + n.at(3) * un.at(6);
531 
532  for ( int k = 1; k <= 3; k++ ) { // nodal val of function w
533  km1 = k - 1;
534  for ( int j = 1; j <= 2; j++ ) { // velocity vector component
535  for ( int m = 1; m <= 3; m++ ) { // nodal components
536  w_dof_addr = k;
537  u_dof_addr = ( m - 1 ) * 2 + j;
538  mm1 = m - 1;
539  d1j = ( j == 1 );
540  d2j = ( j == 2 );
541  answer.at(w_dof_addr, u_dof_addr) += t_pspg * dV * ( b [ km1 ] * ( _u * d1j * b [ mm1 ] + _v * d1j * c [ mm1 ] ) + c [ km1 ] * ( _u * d2j * b [ mm1 ] + _v * d2j * c [ mm1 ] ) );
542  }
543  }
544  }
545  }
546  }
547 }
548 
550 {
551  answer.resize(3);
552  answer.zero();
553 
554 #if 1
555  double Re = static_cast< FluidModel * >( domain->giveEngngModel() )->giveReynoldsNumber();
556  FloatArray eps, stress, u;
557  FloatMatrix _b;
558  FluidDynamicMaterial *mat = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial();
559 
560  // stabilization term K_eps
561  this->computeVectorOfVelocities(VM_Total, tStep, u);
562 
563  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
564  for ( auto &gp : *integrationRulesArray [ ifluid ] ) {
565  double rho = this->_giveMaterial(ifluid)->give('d', gp);
566  double dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
567  double _r = this->computeRadiusAt(gp);
568  this->computeBMtrx(_b, gp);
569  eps.beProductOf(_b, u);
570  mat->computeDeviatoricStressVector(stress, gp, eps, tStep);
571  stress.times(1. / Re);
572  for ( int i = 1; i <= 3; i++ ) {
573  answer.at(i) -= t_pspg * ( b [ i - 1 ] * stress.at(1) + c [ i - 1 ] * stress.at(4) ) * dV / rho / _r;
574  }
575  }
576  }
577 
578 #endif
579 }
580 
582 {
583  answer.resize(3, 6);
584  answer.zero();
585 
586 #if 1
587  double Re = static_cast< FluidModel * >( domain->giveEngngModel() )->giveReynoldsNumber();
588  FloatMatrix _d, _b, _db;
589  //FloatArray eps, stress,u;
590 
591  // stabilization term K_eps
592  //this -> computeVectorOfVelocities(VM_Total,tStep, u) ;
593  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
594  FluidDynamicMaterial *mat = static_cast< FluidDynamicMaterial * >( this->_giveMaterial(ifluid) );
595  for ( GaussPoint *gp: *integrationRulesArray [ ifluid ] ) {
596  double dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
597  double rho = mat->give('d', gp);
598  double _r = this->computeRadiusAt(gp);
599  this->computeBMtrx(_b, gp);
600  mat->giveDeviatoricStiffnessMatrix(_d, TangentStiffness, gp, tStep);
601  _db.beProductOf(_d, _b);
602  //eps.beProductOf (_b, u);
603  //mat->computeDeviatoricStressVector (stress,gp,eps,tStep);
604 
605  for ( int i = 1; i <= 3; i++ ) {
606  for ( int j = 1; j <= 6; j++ ) {
607  answer.at(i, j) -= t_pspg * ( b [ i - 1 ] * ( _db.at(1, j) / _r ) + c [ i - 1 ] * ( _db.at(4, j) / _r ) ) * dV / rho;
608  }
609  }
610  }
611  }
612 
613  answer.times(1. / Re);
614 
615 #endif
616 }
617 
618 
619 
620 
621 void
623 {
624  answer.resize(3, 6);
625  answer.zero();
626  FloatArray n;
627  // M_\epsilon
628 
629  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
630  for ( GaussPoint *gp: *integrationRulesArray [ ifluid ] ) {
631  double dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
632  computeNVector(n, gp);
633 
634  for ( int i = 1; i <= 3; i++ ) {
635  for ( int j = 1; j <= 3; j++ ) {
636  answer.at(i, 2 * j - 1) += t_pspg * dV * b [ i - 1 ] * n.at(j);
637  answer.at(i, 2 * j) += t_pspg * dV * c [ i - 1 ] * n.at(j);
638  }
639  }
640  }
641  }
642 }
643 
644 void
646 {
647  answer.resize(3, 3);
648  answer.zero();
649 
650  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
651  FluidDynamicMaterial *mat = static_cast< FluidDynamicMaterial * >( this->_giveMaterial(ifluid) );
652  for ( GaussPoint *gp: *integrationRulesArray [ ifluid ] ) {
653  double rho = mat->give('d', gp);
654  double dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
655 
656 
657  for ( int i = 1; i <= 3; i++ ) {
658  for ( int j = 1; j <= 3; j++ ) {
659  answer.at(i, j) += t_pspg * dV * ( b [ i - 1 ] * b [ j - 1 ] + c [ i - 1 ] * c [ j - 1 ] ) / rho;
660  }
661  }
662  }
663  }
664 
665 #ifdef TR1_2D_SUPG2_DEBUG
666  /* test */
667  FloatMatrix test;
668  std :: swap(integrationRulesArray [ 0 ], integrationRulesArray [ 1 ]);
670  std :: swap(integrationRulesArray [ 0 ], integrationRulesArray [ 1 ]);
671  for ( int i = 1; i <= 3; i++ ) {
672  for ( int j = 1; j <= 3; j++ ) {
673  if ( fabs( ( answer.at(i, j) - test.at(i, j) ) / test.at(i, j) ) >= 1.e-8 ) {
674  OOFEM_ERROR("test failure (err=%e)", ( answer.at(i, j) - test.at(i, j) ) / test.at(i, j));
675  }
676  }
677  }
678 
679 #endif
680 }
681 
683 void
685 {
686  answer.resize(6);
687  answer.zero();
688 
689  int nLoads;
690  FloatArray un, nV, gVector;
691  double u, v;
692 
693  // add body load (gravity) termms
694  this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), un);
695 
696  nLoads = this->giveBodyLoadArray()->giveSize();
697  for ( int i = 1; i <= nLoads; i++ ) {
698  auto load = domain->giveLoad( bodyLoadArray.at(i) );
699  bcGeomType ltype = load->giveBCGeoType();
700  if ( ltype == BodyLoadBGT && load->giveBCValType() == ForceLoadBVT ) {
701  load->computeComponentArrayAt(gVector, tStep, VM_Total);
702  if ( gVector.giveSize() ) {
703  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
704  for ( GaussPoint *gp: *integrationRulesArray [ ifluid ] ) {
705  double rho = this->_giveMaterial(ifluid)->give('d', gp);
706  double dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
707  double coeff = rho * dV;
708  computeNVector(nV, gp);
709  u = nV.at(1) * un.at(1) + nV.at(2) * un.at(3) + nV.at(3) * un.at(5);
710  v = nV.at(1) * un.at(2) + nV.at(2) * un.at(4) + nV.at(3) * un.at(6);
711 
712  answer.at(1) += coeff * ( gVector.at(1) * ( nV.at(1) + t_supg * ( b [ 0 ] * u + c [ 0 ] * v ) ) );
713  answer.at(2) += coeff * ( gVector.at(2) * ( nV.at(1) + t_supg * ( b [ 0 ] * u + c [ 0 ] * v ) ) );
714  answer.at(3) += coeff * ( gVector.at(1) * ( nV.at(2) + t_supg * ( b [ 1 ] * u + c [ 1 ] * v ) ) );
715  answer.at(4) += coeff * ( gVector.at(2) * ( nV.at(2) + t_supg * ( b [ 1 ] * u + c [ 1 ] * v ) ) );
716  answer.at(5) += coeff * ( gVector.at(1) * ( nV.at(3) + t_supg * ( b [ 2 ] * u + c [ 2 ] * v ) ) );
717  answer.at(6) += coeff * ( gVector.at(2) * ( nV.at(3) + t_supg * ( b [ 2 ] * u + c [ 2 ] * v ) ) );
718  }
719  }
720  }
721  }
722  }
723 
724  // loop over sides
725  int n1, n2, code, sid;
726  double tx, ty, l, side_r;
727  //IntArray nodecounter (3);
728  for ( int j = 1; j <= boundarySides.giveSize(); j++ ) {
729  code = boundaryCodes.at(j);
730  sid = boundarySides.at(j);
731  if ( ( code & FMElement_PrescribedTractionBC ) ) {
732  FloatArray t, coords(1);
733  int n, id;
734  // integrate tractions
735  n1 = sid;
736  n2 = ( n1 == 3 ? 1 : n1 + 1 );
737 
738  tx = giveNode(n2)->giveCoordinate(1) - giveNode(n1)->giveCoordinate(1);
739  ty = giveNode(n2)->giveCoordinate(2) - giveNode(n1)->giveCoordinate(2);
740  l = sqrt(tx * tx + ty * ty);
741  // radius at side center
742  side_r = 0.5 * ( giveNode(n2)->giveCoordinate(1) + giveNode(n1)->giveCoordinate(1) );
743 
744  // if no traction bc applied but side marked as with traction load
745  // then zero traction is assumed !!!
746 
747  // loop over boundary load array
748  int numLoads = this->giveBoundaryLoadArray()->giveSize() / 2;
749  for ( int i = 1; i <= numLoads; i++ ) {
750  n = boundaryLoadArray.at(1 + ( i - 1 ) * 2);
751  id = boundaryLoadArray.at(i * 2);
752  if ( id != sid ) {
753  continue;
754  }
755 
756  auto bLoad = dynamic_cast< BoundaryLoad * >( domain->giveLoad(n) );
757  if ( bLoad ) {
758  bLoad->computeValueAt(t, tStep, coords, VM_Total);
759 
760  // here it is assumed constant traction, one point integration only
761  // n1 (u,v)
762  answer.at( ( n1 - 1 ) * 2 + 1 ) += t.at(1) * side_r * l / 2.;
763  answer.at(n1 * 2) += t.at(2) * side_r * l / 2.;
764  // n2 (u,v)
765  answer.at( ( n2 - 1 ) * 2 + 1 ) += t.at(1) * side_r * l / 2.;
766  answer.at(n2 * 2) += t.at(2) * side_r * l / 2.;
767 
768  //answer.at(n1)+= (t.at(1)*nx + t.at(2)*ny) * side_r * l/2.;
769  //answer.at(n2)+= (t.at(1)*nx + t.at(2)*ny) * side_r * l/2.;
770  }
771  }
772  }
773  }
774 }
775 
776 void
778 {
779  int nLoads;
780  FloatArray gVector;
781 
782  answer.resize(3);
783  answer.zero();
784  nLoads = this->giveBodyLoadArray()->giveSize();
785  for ( int i = 1; i <= nLoads; i++ ) {
786  Load *load = domain->giveLoad( bodyLoadArray.at(i) );
787  bcGeomType ltype = load->giveBCGeoType();
788  if ( ( ltype == BodyLoadBGT ) && ( load->giveBCValType() == ForceLoadBVT ) ) {
789  load->computeComponentArrayAt(gVector, tStep, VM_Total);
790  if ( gVector.giveSize() ) {
791  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
792  for ( GaussPoint *gp: *integrationRulesArray [ ifluid ] ) {
793  double dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
794  double coeff = t_pspg * dV;
795 
796  answer.at(1) += coeff * ( b [ 0 ] * gVector.at(1) + c [ 0 ] * gVector.at(2) );
797  answer.at(2) += coeff * ( b [ 1 ] * gVector.at(1) + c [ 1 ] * gVector.at(2) );
798  answer.at(3) += coeff * ( b [ 2 ] * gVector.at(1) + c [ 2 ] * gVector.at(2) );
799  }
800  }
801  }
802  }
803  }
804 }
805 
806 
807 
808 void
810 {
811  //TR1_2D_SUPG :: updateStabilizationCoeffs (tStep);
812 #if 0
813  int i, j, k, l, w_dof_addr, u_dof_addr, ip, ifluid;
814  double __g_norm, __gamma_norm, __gammav_norm, __beta_norm, __betav_norm, __c_norm, __e_norm, __k_norm, __Re;
815  double __t_p1, __t_p2, __t_p3, __t_pv1, __t_pv2, __t_pv3;
816  double nu, nu0, nu1, usum, vsum, rho, dV, u1, u2;
817  FloatArray u, un, a;
818 
819  // compute averaged viscosity based on rule of mixture
820  GaussPoint *gp;
821  if ( integrationRulesArray [ 0 ]->giveNumberOfIntegrationPoints() ) {
822  gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
823  } else {
824  gp = integrationRulesArray [ 1 ]->getIntegrationPoint(0);
825  }
826 
827  nu0 = this->_giveMaterial(0)->giveCharacteristicValue( MRM_Viscosity, gp, tStep->givePreviousStep() );
828  nu1 = this->_giveMaterial(1)->giveCharacteristicValue( MRM_Viscosity, gp, tStep->givePreviousStep() );
829  nu = vof * nu0 + ( 1. - vof ) * nu1;
830 
831  //this -> computeVectorOfVelocities(VM_Total,tStep->givePreviousStep(),un) ;
832  this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), u);
833  this->computeVectorOfVelocities(VM_Acceleration, tStep, a);
834  un = u;
835  usum = un.at(1) + un.at(3) + un.at(5);
836  vsum = un.at(2) + un.at(4) + un.at(6);
837 
838  FloatMatrix __tmp;
839  FloatArray __tmpvec, n;
840  // assemble g matrix
841  __tmp.resize(6, 3);
842  double ar3 = area / 3.0;
843 
844  __tmp.at(1, 1) = __tmp.at(1, 2) = __tmp.at(1, 3) = b [ 0 ] * ar3;
845  __tmp.at(3, 1) = __tmp.at(3, 2) = __tmp.at(3, 3) = b [ 1 ] * ar3;
846  __tmp.at(5, 1) = __tmp.at(5, 2) = __tmp.at(5, 3) = b [ 2 ] * ar3;
847 
848  __tmp.at(2, 1) = __tmp.at(2, 2) = __tmp.at(2, 3) = c [ 0 ] * ar3;
849  __tmp.at(4, 1) = __tmp.at(4, 2) = __tmp.at(4, 3) = c [ 1 ] * ar3;
850  __tmp.at(6, 1) = __tmp.at(6, 2) = __tmp.at(6, 3) = c [ 2 ] * ar3;
851 
852  __g_norm = __tmp.computeFrobeniusNorm();
853 
854  // assemble \gamma matrix (advectionTerm of mass conservation eq)
855  __tmp.resize(3, 6);
856  for ( k = 1; k <= 3; k++ ) {
857  for ( l = 1; l <= 3; l++ ) {
858  __tmp.at(k, l * 2 - 1) = ar3 * b [ k - 1 ] * ( usum * b [ l - 1 ] + vsum * c [ l - 1 ] );
859  __tmp.at(k, l * 2) = ar3 * c [ k - 1 ] * ( usum * b [ l - 1 ] + vsum * c [ l - 1 ] );
860  }
861  }
862 
863  __gamma_norm = __tmp.computeFrobeniusNorm();
864  __tmpvec.beProductOf(__tmp, u);
865  __gammav_norm = __tmpvec.computeNorm();
866  // compute beta mtrx (acceleration term of mass conservation eq)
867  __tmp.resize(3, 6);
868  __tmp.zero();
869  __tmp.at(1, 1) = __tmp.at(1, 3) = __tmp.at(1, 5) = ar3 * b [ 0 ];
870  __tmp.at(1, 2) = __tmp.at(1, 4) = __tmp.at(1, 6) = ar3 * c [ 0 ];
871  __tmp.at(2, 1) = __tmp.at(2, 3) = __tmp.at(2, 5) = ar3 * b [ 1 ];
872  __tmp.at(2, 2) = __tmp.at(2, 4) = __tmp.at(2, 6) = ar3 * c [ 1 ];
873  __tmp.at(3, 1) = __tmp.at(3, 3) = __tmp.at(3, 5) = ar3 * b [ 2 ];
874  __tmp.at(3, 2) = __tmp.at(3, 4) = __tmp.at(3, 6) = ar3 * c [ 2 ];
875  __beta_norm = __tmp.computeFrobeniusNorm();
876  __tmpvec.beProductOf(__tmp, a);
877  __betav_norm = __tmpvec.computeNorm();
878  // compute c mtrx (advection term of momentum balance)
879  // standard galerkin term
880  __tmp.resize(6, 6);
881  __tmp.zero();
882  for ( ifluid = 0; ifluid < 2; ifluid++ ) {
883  for ( GaussPoint *gp: integrationRulesArray [ ifluid ] ) {
884  rho = this->_giveMaterial(ifluid)->give('d', gp);
885  this->computeNMtrx(n, gp);
886  dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
887 
888  u1 = n.at(1) * un.at(1) + n.at(2) * un.at(3) + n.at(3) * un.at(5);
889  u2 = n.at(1) * un.at(2) + n.at(2) * un.at(4) + n.at(3) * un.at(6);
890  for ( i = 1; i <= 2; i++ ) {
891  for ( k = 1; k <= 3; k++ ) {
892  for ( l = 1; l <= 3; l++ ) {
893  w_dof_addr = ( k - 1 ) * 2 + i;
894  u_dof_addr = ( l - 1 ) * 2 + i;
895  __tmp.at(w_dof_addr, u_dof_addr) += rho * dV * n.at(k) * ( u1 * b [ l - 1 ] + u2 * c [ l - 1 ] );
896  }
897  }
898  }
899  }
900  }
901 
902  __c_norm = __tmp.computeFrobeniusNorm();
903  // compute e mtrx (advection term of momentum balance)
904  __tmp.resize(6, 6);
905  __tmp.zero();
906  double __n[] = {
907  b [ 0 ], c [ 0 ], b [ 1 ], c [ 1 ], b [ 2 ], c [ 2 ]
908  };
909  for ( ifluid = 0; ifluid < 2; ifluid++ ) {
910  for ( GaussPoint *gp: *integrationRulesArray [ ifluid ] ) {
911  rho = this->_giveMaterial(ifluid)->give('d', gp);
912  dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
913 
914  for ( i = 1; i <= 6; i++ ) {
915  for ( j = 1; j <= 6; j++ ) {
916  __tmp.at(i, j) += dV * rho * __n [ i - 1 ] * __n [ j - 1 ];
917  }
918  }
919  }
920  }
921 
922  __e_norm = __tmp.computeFrobeniusNorm();
923  // compute element level Reynolds number
924  // compute k norm first
925  __tmp.resize(6, 6);
926  __tmp.zero();
927  for ( ifluid = 0; ifluid < 2; ifluid++ ) {
928  for ( GaussPoint *gp: *integrationRulesArray [ ifluid ] ) {
929  rho = this->_giveMaterial(ifluid)->give('d', gp);
930  this->computeNMtrx(n, gp);
931  dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
932 
933  u1 = n.at(1) * un.at(1) + n.at(2) * un.at(3) + n.at(3) * un.at(5);
934  u2 = n.at(1) * un.at(2) + n.at(2) * un.at(4) + n.at(3) * un.at(6);
935  for ( k = 1; k <= 3; k++ ) {
936  for ( l = 1; l <= 3; l++ ) {
937  __tmp.at(k * 2 - 1, l * 2 - 1) += rho * dV * ( u1 * b [ k - 1 ] + u2 * c [ k - 1 ] ) * ( u1 * b [ l - 1 ] + u2 * c [ l - 1 ] );
938  __tmp.at(k * 2, l * 2) += rho * dV * ( u1 * b [ k - 1 ] + u2 * c [ k - 1 ] ) * ( u1 * b [ l - 1 ] + u2 * c [ l - 1 ] );
939  }
940  }
941  }
942  }
943 
944  __k_norm = __tmp.computeFrobeniusNorm();
945  double u_1, u_2, vnorm = 0.;
946  int im1;
947  for ( i = 1; i <= 3; i++ ) {
948  im1 = i - 1;
949  u_1 = u.at( ( im1 ) * 2 + 1 );
950  u_2 = u.at( ( im1 ) * 2 + 2 );
951  vnorm = max( vnorm, sqrt(u_1 * u_1 + u_2 * u_2) );
952  }
953 
954  if ( vnorm == 0.0 ) {
955  //t_sugn1 = inf;
956  double t_sugn2 = tStep->giveTimeIncrement() / 2.0;
957  //t_sugn3 = inf;
958  this->t_supg = 1. / sqrt( 1. / ( t_sugn2 * t_sugn2 ) );
959  this->t_pspg = this->t_supg;
960  this->t_lsic = 0.0;
961  } else {
962  __Re = vnorm * vnorm * __c_norm / __k_norm / nu;
963 
964  __t_p1 = __g_norm / __gamma_norm;
965  __t_p2 = tStep->giveTimeIncrement() * __g_norm / 2.0 / __beta_norm;
966  __t_p3 = __t_p1 * __Re;
967  this->t_pspg = 1. / sqrt( 1. / ( __t_p1 * __t_p1 ) + 1. / ( __t_p2 * __t_p2 ) + 1. / ( __t_p3 * __t_p3 ) );
968 
969  __t_pv1 = __t_p1;
970  __t_pv2 = __t_pv1 * __gammav_norm / __betav_norm;
971  __t_pv3 = __t_pv1 * __Re;
972  this->t_supg = 1. / sqrt( 1. / ( __t_pv1 * __t_pv1 ) + 1. / ( __t_pv2 * __t_pv2 ) + 1. / ( __t_pv3 * __t_pv3 ) );
973 
974  this->t_lsic = __c_norm / __e_norm;
975  }
976 
977 #else
978  /* UGN-Based Stabilization */
979  double h_ugn, sum = 0.0, vnorm, t_sugn1, t_sugn2, t_sugn3, u_1, u_2, z, Re_ugn;
980  double dscale, uscale, lscale, tscale, dt;
981  //bool zeroFlag = false;
982  int im1;
983  FloatArray u;
984 
989 
990  this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), u);
991  u.times(uscale);
992  double nu, nu0, nu1;
993 
994  // compute averaged viscosity based on rule of mixture
995  GaussPoint *gp;
996  if ( integrationRulesArray [ 0 ]->giveNumberOfIntegrationPoints() ) {
997  gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
998  } else {
999  gp = integrationRulesArray [ 1 ]->getIntegrationPoint(0);
1000  }
1001 
1002  nu0 = static_cast< FluidDynamicMaterial * >( this->_giveMaterial(0) )->giveEffectiveViscosity( gp, tStep->givePreviousStep() );
1003  nu1 = static_cast< FluidDynamicMaterial * >( this->_giveMaterial(1) )->giveEffectiveViscosity( gp, tStep->givePreviousStep() );
1004  nu = vof * nu0 + ( 1. - vof ) * nu1;
1006 
1007  dt = tStep->giveTimeIncrement() * tscale;
1008 
1009  for ( int i = 1; i <= 3; i++ ) {
1010  im1 = i - 1;
1011  sum += fabs(u.at( ( im1 ) * 2 + 1 ) * b [ im1 ] / lscale + u.at(im1 * 2 + 2) * c [ im1 ] / lscale);
1012  }
1013 
1014  /*
1015  * u_1=(u.at(1)+u.at(3)+u.at(5))/3.0;
1016  * u_2=(u.at(2)+u.at(4)+u.at(6))/3.0;
1017  * vnorm=sqrt(u_1*u_1+u_2*u_2);
1018  */
1019  vnorm = 0.;
1020  for ( int i = 1; i <= 3; i++ ) {
1021  im1 = i - 1;
1022  u_1 = u.at( ( im1 ) * 2 + 1 );
1023  u_2 = u.at( ( im1 ) * 2 + 2 );
1024  vnorm = max( vnorm, sqrt(u_1 * u_1 + u_2 * u_2) );
1025  }
1026 
1027  if ( ( vnorm == 0.0 ) || ( sum < vnorm * 1e-10 ) ) {
1028  //t_sugn1 = inf;
1029  t_sugn2 = dt / 2.0;
1030  //t_sugn3 = inf;
1031  this->t_supg = 1. / sqrt( 1. / ( t_sugn2 * t_sugn2 ) );
1032  this->t_pspg = this->t_supg;
1033  this->t_lsic = 0.0;
1034  } else {
1035  h_ugn = 2.0 * vnorm / sum;
1036  t_sugn1 = 1. / sum;
1037  t_sugn2 = dt / 2.0;
1038  t_sugn3 = h_ugn * h_ugn / 4.0 / nu;
1039 
1040  this->t_supg = 1. / sqrt( 1. / ( t_sugn1 * t_sugn1 ) + 1. / ( t_sugn2 * t_sugn2 ) + 1. / ( t_sugn3 * t_sugn3 ) );
1041  this->t_pspg = this->t_supg;
1042 
1043  Re_ugn = vnorm * h_ugn / ( 2. * nu );
1044  z = ( Re_ugn <= 3. ) ? Re_ugn / 3. : 1.0;
1045  this->t_lsic = h_ugn * vnorm * z / 2.0;
1046  }
1047 
1048  // if (this->number == 1) {
1049  // printf ("t_supg %e t_pspg %e t_lsic %e\n", t_supg, t_pspg, t_lsic);
1050  // }
1051 
1052 
1053  this->t_supg *= uscale / lscale;
1054  this->t_pspg *= 1. / ( lscale * dscale );
1055  this->t_lsic *= ( dscale * uscale ) / ( lscale * lscale );
1056 
1057  this->t_lsic = 0.0;
1058 
1059 #endif
1060 
1061  //this->t_lsic=0.0;
1062  //this->t_pspg=0.0;
1063 }
1064 
1065 
1066 
1067 double
1069 {
1070  return 1.e3;
1071 }
1072 
1073 
1074 void
1076 {
1077  /* one computes here average deviatoric stress, based on rule of mixture (this is used only for postprocessing) */
1078  FloatArray u, eps, s0, s1;
1079  FloatMatrix _b;
1080  answer.resize(3);
1081 
1082 
1083  this->computeVectorOfVelocities(VM_Total, tStep, u);
1084  this->computeBMtrx(_b, gp);
1085  eps.beProductOf(_b, u);
1086 
1087  static_cast< FluidDynamicMaterial * >( this->_giveMaterial(0) )->computeDeviatoricStressVector(s0, gp, eps, tStep);
1088  static_cast< FluidDynamicMaterial * >( this->_giveMaterial(1) )->computeDeviatoricStressVector(s1, gp, eps, tStep);
1089 
1090  for ( int i = 1; i <= 3; i++ ) {
1091  answer.at(i) = ( temp_vof ) * s0.at(i) + ( 1. - temp_vof ) * s1.at(i);
1092  }
1093 }
1094 
1095 
1096 /*
1097  * double
1098  * TR1_2D_SUPG2 :: computeCriticalTimeStep (TimeStep* tStep)
1099  * {
1100  * FloatArray u;
1101  * double dt1, dt2, dt;
1102  * double Re = static_cast<FluidModel*>(domain->giveEngngModel())->giveReynoldsNumber();
1103  *
1104  * this -> computeVectorOfVelocities(VM_Total,tStep, u) ;
1105  *
1106  * double vn1 = sqrt(u.at(1)*u.at(1)+u.at(2)*u.at(2));
1107  * double vn2 = sqrt(u.at(3)*u.at(3)+u.at(4)*u.at(4));
1108  * double vn3 = sqrt(u.at(5)*u.at(5)+u.at(6)*u.at(6));
1109  * double veln = max (vn1, max(vn2,vn3));
1110  *
1111  * double l1 = 1.0/(sqrt(b[0]*b[0]+c[0]*c[0]));
1112  * double l2 = 1.0/(sqrt(b[1]*b[1]+c[1]*c[1]));
1113  * double l3 = 1.0/(sqrt(b[2]*b[2]+c[2]*c[2]));
1114  *
1115  * double ln = min (l1, min (l2,l3));
1116  *
1117  * // viscous limit
1118  * dt2 = 0.5*ln*ln*Re;
1119  * if (veln != 0.0) {
1120  * dt1 = ln/veln;
1121  * dt = dt1*dt2/(dt1+dt2);
1122  * } else {
1123  * dt = dt2;
1124  * }
1125  * return dt;
1126  * }
1127  */
1128 
1129 
1130 double
1131 TR1_2D_SUPG2_AXI :: computeLEPLICVolumeFraction(const FloatArray &n, const double p, LEPlic *matInterface, bool updFlag)
1132 {
1133  Polygon pg;
1134  double answer, volume = computeMyVolume(matInterface, updFlag);
1135  this->formVolumeInterfacePoly(pg, matInterface, n, p, updFlag);
1136  answer = fabs(pg.computeVolume() / volume);
1137  if ( answer > 1.000000001 ) {
1138  OOFEM_WARNING("VOF fraction out of bounds, vof = %e\n", answer);
1139  return 1.0;
1140  } else {
1141  return answer;
1142  }
1143 }
1144 
1145 void
1147  const FloatArray &normal, const double p, bool updFlag)
1148 {
1149  double x, y;
1150  Vertex v;
1151 
1152  matvolpoly.clear();
1153 
1154  if ( this->vof <= TRSUPG_ZERO_VOF ) {
1155  return;
1156  } else if ( this->vof >= ( 1 - TRSUPG_ZERO_VOF ) ) {
1157  for ( int i = 1; i <= 3; i++ ) {
1158  if ( updFlag ) {
1159  x = matInterface->giveUpdatedXCoordinate( this->giveNode(i)->giveNumber() );
1160  y = matInterface->giveUpdatedYCoordinate( this->giveNode(i)->giveNumber() );
1161  } else {
1162  x = this->giveNode(i)->giveCoordinate(1);
1163  y = this->giveNode(i)->giveCoordinate(2);
1164  }
1165 
1166  v.setCoords(x, y);
1167  matvolpoly.addVertex(v);
1168  }
1169 
1170  return;
1171  }
1172 
1173  this->formVolumeInterfacePoly(matvolpoly, matInterface, normal, p, updFlag);
1174 }
1175 
1176 
1177 void
1179  const FloatArray &normal, const double p, bool updFlag)
1180 {
1181  int next;
1182  bool nodeIn [ 3 ];
1183  double nx = normal.at(1), ny = normal.at(2), x, y;
1184  double tx, ty;
1185  Vertex v;
1186 
1187  matvolpoly.clear();
1188 
1189  for ( int i = 1; i <= 3; i++ ) {
1190  if ( updFlag ) {
1191  x = matInterface->giveUpdatedXCoordinate( this->giveNode(i)->giveNumber() );
1192  y = matInterface->giveUpdatedYCoordinate( this->giveNode(i)->giveNumber() );
1193  } else {
1194  x = this->giveNode(i)->giveCoordinate(1);
1195  y = this->giveNode(i)->giveCoordinate(2);
1196  }
1197 
1198  if ( ( nx * x + ny * y + p ) >= 0. ) {
1199  nodeIn [ i - 1 ] = true;
1200  } else {
1201  nodeIn [ i - 1 ] = false;
1202  }
1203  }
1204 
1205  if ( nodeIn [ 0 ] && nodeIn [ 1 ] && nodeIn [ 2 ] ) { // all nodes inside
1206  for ( int i = 1; i <= 3; i++ ) {
1207  if ( updFlag ) {
1208  x = matInterface->giveUpdatedXCoordinate( this->giveNode(i)->giveNumber() );
1209  y = matInterface->giveUpdatedYCoordinate( this->giveNode(i)->giveNumber() );
1210  } else {
1211  x = this->giveNode(i)->giveCoordinate(1);
1212  y = this->giveNode(i)->giveCoordinate(2);
1213  }
1214 
1215  v.setCoords(x, y);
1216  matvolpoly.addVertex(v);
1217  }
1218 
1219  return;
1220  } else if ( !( nodeIn [ 0 ] || nodeIn [ 1 ] || nodeIn [ 2 ] ) ) { // all nodes outside
1221  return;
1222  } else {
1223  for ( int i = 1; i <= 3; i++ ) {
1224  next = i < 3 ? i + 1 : 1;
1225  if ( nodeIn [ i - 1 ] ) {
1226  if ( updFlag ) {
1227  v.setCoords( matInterface->giveUpdatedXCoordinate( this->giveNode(i)->giveNumber() ),
1228  matInterface->giveUpdatedYCoordinate( this->giveNode(i)->giveNumber() ) );
1229  } else {
1230  v.setCoords( this->giveNode(i)->giveCoordinate(1),
1231  this->giveNode(i)->giveCoordinate(2) );
1232  }
1233 
1234  matvolpoly.addVertex(v);
1235  }
1236 
1237  if ( nodeIn [ next - 1 ] ^ nodeIn [ i - 1 ] ) {
1238  // compute intersection with (i,next) edge
1239  if ( updFlag ) {
1240  x = matInterface->giveUpdatedXCoordinate( this->giveNode(i)->giveNumber() );
1241  y = matInterface->giveUpdatedYCoordinate( this->giveNode(i)->giveNumber() );
1242  tx = matInterface->giveUpdatedXCoordinate( this->giveNode(next)->giveNumber() ) - x;
1243  ty = matInterface->giveUpdatedYCoordinate( this->giveNode(next)->giveNumber() ) - y;
1244  } else {
1245  x = this->giveNode(i)->giveCoordinate(1);
1246  y = this->giveNode(i)->giveCoordinate(2);
1247  tx = this->giveNode(next)->giveCoordinate(1) - x;
1248  ty = this->giveNode(next)->giveCoordinate(2) - y;
1249  }
1250 
1251  double s, sd = nx * tx + ny * ty;
1252  if ( fabs(sd) > 1.e-10 ) {
1253  s = ( -p - ( nx * x + ny * y ) ) / sd;
1254  v.setCoords(x + tx * s, y + ty * s);
1255  matvolpoly.addVertex(v);
1256  } else {
1257  // pathological case - lines are parallel
1258  if ( nodeIn [ i - 1 ] ) {
1259  if ( updFlag ) {
1260  v.setCoords( matInterface->giveUpdatedXCoordinate( this->giveNode(next)->giveNumber() ),
1261  matInterface->giveUpdatedYCoordinate( this->giveNode(next)->giveNumber() ) );
1262  } else {
1263  v.setCoords( this->giveNode(next)->giveCoordinate(1), this->giveNode(next)->giveCoordinate(2) );
1264  }
1265 
1266  matvolpoly.addVertex(v);
1267  } else {
1268  v.setCoords(x, y);
1269  matvolpoly.addVertex(v);
1270  if ( updFlag ) {
1271  v.setCoords( matInterface->giveUpdatedXCoordinate( this->giveNode(next)->giveNumber() ),
1272  matInterface->giveUpdatedYCoordinate( this->giveNode(next)->giveNumber() ) );
1273  } else {
1274  v.setCoords( this->giveNode(next)->giveCoordinate(1), this->giveNode(next)->giveCoordinate(2) );
1275  }
1276 
1277  matvolpoly.addVertex(v);
1278  }
1279  }
1280  }
1281  } // end loop over elem nodes
1282  }
1283 }
1284 
1285 
1286 void
1287 TR1_2D_SUPG2_AXI :: updateVolumePolygons(Polygon &referenceFluidPoly, Polygon &secondFluidPoly, int &rfPoints, int &sfPoints,
1288  const FloatArray &normal, const double p, bool updFlag)
1289 {
1290  /*
1291  * this method updates two polygons, one filled with reference fluid and second filled with
1292  * other fluid (air). These two polygons are used in integrating element contributions.
1293  */
1294  int next;
1295  bool nodeIn [ 3 ];
1296  double nx = normal.at(1), ny = normal.at(2), x, y;
1297  double tx, ty;
1298  Vertex v;
1299 
1300  rfPoints = sfPoints = 0;
1301  referenceFluidPoly.clear();
1302  secondFluidPoly.clear();
1303 
1304  for ( int i = 1; i <= 3; i++ ) {
1305  x = this->giveNode(i)->giveCoordinate(1);
1306  y = this->giveNode(i)->giveCoordinate(2);
1307 
1308  if ( ( nx * x + ny * y + p ) >= 0. ) {
1309  nodeIn [ i - 1 ] = true;
1310  } else {
1311  nodeIn [ i - 1 ] = false;
1312  }
1313  }
1314 
1315  if ( nodeIn [ 0 ] && nodeIn [ 1 ] && nodeIn [ 2 ] ) { // all nodes inside
1316  for ( int i = 1; i <= 3; i++ ) {
1317  x = this->giveNode(i)->giveCoordinate(1);
1318  y = this->giveNode(i)->giveCoordinate(2);
1319 
1320  v.setCoords(x, y);
1321  referenceFluidPoly.addVertex(v);
1322  rfPoints++;
1323  }
1324 
1325  return;
1326  } else if ( !( nodeIn [ 0 ] || nodeIn [ 1 ] || nodeIn [ 2 ] ) ) { // all nodes outside
1327  for ( int i = 1; i <= 3; i++ ) {
1328  x = this->giveNode(i)->giveCoordinate(1);
1329  y = this->giveNode(i)->giveCoordinate(2);
1330 
1331  v.setCoords(x, y);
1332  secondFluidPoly.addVertex(v);
1333  sfPoints++;
1334  }
1335 
1336  return;
1337  } else {
1338  for ( int i = 1; i <= 3; i++ ) {
1339  next = i < 3 ? i + 1 : 1;
1340  if ( nodeIn [ i - 1 ] ) {
1341  v.setCoords( this->giveNode(i)->giveCoordinate(1),
1342  this->giveNode(i)->giveCoordinate(2) );
1343 
1344  referenceFluidPoly.addVertex(v);
1345  rfPoints++;
1346  } else {
1347  v.setCoords( this->giveNode(i)->giveCoordinate(1),
1348  this->giveNode(i)->giveCoordinate(2) );
1349 
1350  secondFluidPoly.addVertex(v);
1351  sfPoints++;
1352  }
1353 
1354  if ( nodeIn [ next - 1 ] ^ nodeIn [ i - 1 ] ) {
1355  // compute intersection with (i,next) edge
1356  x = this->giveNode(i)->giveCoordinate(1);
1357  y = this->giveNode(i)->giveCoordinate(2);
1358  tx = this->giveNode(next)->giveCoordinate(1) - x;
1359  ty = this->giveNode(next)->giveCoordinate(2) - y;
1360 
1361  double s, sd = nx * tx + ny * ty;
1362  if ( fabs(sd) > 1.e-10 ) {
1363  s = ( -p - ( nx * x + ny * y ) ) / sd;
1364  v.setCoords(x + tx * s, y + ty * s);
1365  referenceFluidPoly.addVertex(v);
1366  secondFluidPoly.addVertex(v);
1367  rfPoints++;
1368  sfPoints++;
1369  } else {
1370  // pathological case - lines are parallel
1371  if ( nodeIn [ i - 1 ] ) {
1372  v.setCoords( this->giveNode(next)->giveCoordinate(1), this->giveNode(next)->giveCoordinate(2) );
1373 
1374  referenceFluidPoly.addVertex(v);
1375  secondFluidPoly.addVertex(v);
1376  rfPoints++;
1377  sfPoints++;
1378  } else {
1379  //v.setCoords(x, y);
1380  //referenceFluidPoly.addVertex(v);
1381  v.setCoords( this->giveNode(next)->giveCoordinate(1), this->giveNode(next)->giveCoordinate(2) );
1382 
1383  referenceFluidPoly.addVertex(v);
1384  secondFluidPoly.addVertex(v);
1385  rfPoints++;
1386  sfPoints++;
1387  }
1388  }
1389  }
1390  } // end loop over elem nodes
1391  }
1392 }
1393 
1394 double
1395 TR1_2D_SUPG2_AXI :: truncateMatVolume(const Polygon &matvolpoly, double &volume)
1396 {
1397  Polygon me, clip;
1398  Graph g;
1399 
1400  this->formMyVolumePoly(me, NULL, false);
1401  g.clip(clip, me, matvolpoly);
1402 #ifdef __OOFEG
1403  EASValsSetColor( gc [ 0 ].getActiveCrackColor() );
1404  //GraphicObj *go = clip.draw(::gc[OOFEG_DEBUG_LAYER],true);
1405  clip.draw(gc [ OOFEG_DEBUG_LAYER ], true);
1406  //EVFastRedraw(myview);
1407 #endif
1408  volume = clip.computeVolume();
1409  return volume / area;
1410 }
1411 
1412 void
1413 TR1_2D_SUPG2_AXI :: formMyVolumePoly(Polygon &me, LEPlic *matInterface, bool updFlag)
1414 {
1415  double x, y;
1416  Vertex v;
1417 
1418  me.clear();
1419 
1420  for ( int i = 1; i <= 3; i++ ) {
1421  if ( updFlag ) {
1422  x = matInterface->giveUpdatedXCoordinate( this->giveNode(i)->giveNumber() );
1423  y = matInterface->giveUpdatedYCoordinate( this->giveNode(i)->giveNumber() );
1424  } else {
1425  x = this->giveNode(i)->giveCoordinate(1);
1426  y = this->giveNode(i)->giveCoordinate(2);
1427  }
1428 
1429  v.setCoords(x, y);
1430  me.addVertex(v);
1431  }
1432 }
1433 
1434 
1435 double
1436 TR1_2D_SUPG2_AXI :: computeMyVolume(LEPlic *matInterface, bool updFlag)
1437 {
1438  double x1, x2, x3, y1, y2, y3;
1439  if ( updFlag ) {
1440  x1 = matInterface->giveUpdatedXCoordinate( this->giveNode(1)->giveNumber() );
1441  x2 = matInterface->giveUpdatedXCoordinate( this->giveNode(2)->giveNumber() );
1442  x3 = matInterface->giveUpdatedXCoordinate( this->giveNode(3)->giveNumber() );
1443 
1444  y1 = matInterface->giveUpdatedYCoordinate( this->giveNode(1)->giveNumber() );
1445  y2 = matInterface->giveUpdatedYCoordinate( this->giveNode(2)->giveNumber() );
1446  y3 = matInterface->giveUpdatedYCoordinate( this->giveNode(3)->giveNumber() );
1447  return 0.5 * ( x2 * y3 + x1 * y2 + y1 * x3 - x2 * y1 - x3 * y2 - x1 * y3 );
1448  } else {
1449  return area;
1450  }
1451 }
1452 
1453 void
1454 TR1_2D_SUPG2_AXI :: giveElementCenter(LEPlic *mat_interface, FloatArray &center, bool upd)
1455 {
1456  FloatArray coords;
1457 
1458  center.resize(2);
1459  center.zero();
1460  if ( upd ) {
1461  for ( int i = 1; i <= 3; i++ ) {
1462  mat_interface->giveUpdatedCoordinate( coords, giveNode(i)->giveNumber() );
1463  center.add(coords);
1464  }
1465  } else {
1466  for ( int i = 1; i <= 3; i++ ) {
1467  center.at(1) += this->giveNode(i)->giveCoordinate(1);
1468  center.at(2) += this->giveNode(i)->giveCoordinate(2);
1469  }
1470  }
1471 
1472  center.times(1. / 3.);
1473 }
1474 
1475 void
1477 {
1478  int c [ 2 ];
1479  int nip;
1480  double x, y;
1481  Vertex v;
1482 
1483  for ( int i = 0; i < 2; i++ ) {
1484  myPoly [ i ].clear();
1485  }
1486 
1487  if ( this->temp_vof <= TRSUPG_ZERO_VOF ) {
1488  for ( int i = 1; i <= 3; i++ ) {
1489  x = this->giveNode(i)->giveCoordinate(1);
1490  y = this->giveNode(i)->giveCoordinate(2);
1491  v.setCoords(x, y);
1492  myPoly [ 1 ].addVertex(v);
1493  c [ 1 ] = 3;
1494  c [ 0 ] = 0;
1495  }
1496  } else if ( this->temp_vof >= ( 1. - TRSUPG_ZERO_VOF ) ) {
1497  for ( int i = 1; i <= 3; i++ ) {
1498  x = this->giveNode(i)->giveCoordinate(1);
1499  y = this->giveNode(i)->giveCoordinate(2);
1500  v.setCoords(x, y);
1501  myPoly [ 0 ].addVertex(v);
1502  c [ 0 ] = 3;
1503  c [ 1 ] = 0;
1504  }
1505  } else {
1506  this->updateVolumePolygons(myPoly [ 0 ], myPoly [ 1 ], c [ 0 ], c [ 1 ], temp_normal, temp_p, false);
1507  }
1508 
1509  integrationRulesArray [ 0 ]->clear();
1510  integrationRulesArray [ 1 ]->clear();
1511 
1512  FloatArray gc, lc;
1513  const Vertex *p;
1514  FEI2dTrLin triaApprox(1, 2);
1515  FEI2dQuadLin quadApprox(1, 2);
1516  FEInterpolation *approx = NULL;
1517  // set up integration points
1518  for ( int i = 0; i < 2; i++ ) {
1519  if ( c [ i ] == 3 ) {
1520  id [ i ] = _Triangle;
1521  approx = & triaApprox;
1522  } else if ( c [ i ] == 4 ) {
1523  id [ i ] = _Square;
1524  approx = & quadApprox;
1525  } else if ( c [ i ] == 0 ) {
1526  continue;
1527  } else {
1528  OOFEM_ERROR("cannot set up integration domain for %d vertex polygon", c [ i ]);
1529  }
1530 
1531  vcoords [ i ].clear();
1532  vcoords [ i ].reserve( c [ i ] );
1533 
1534  // set up vertex coords
1536  while ( it.giveNext(& p) ) {
1537  vcoords [ i ].push_back( *p->getCoords() );
1538  }
1539 
1540  if ( id [ i ] == _Triangle ) {
1541  nip = 4;
1542  } else {
1543  nip = 4;
1544  }
1545 
1546  this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ i ], nip, this);
1547 
1548  // remap ip coords into area coords of receiver
1549  for ( GaussPoint *gp: *integrationRulesArray [ i ] ) {
1550  approx->local2global( gc, gp->giveNaturalCoordinates(), FEIVertexListGeometryWrapper(vcoords [ i ]) );
1551  triaApprox.global2local( lc, gc, FEIElementGeometryWrapper(this) );
1552  // modify original ip coords to target ones
1553  gp->setSubPatchCoordinates( gp->giveNaturalCoordinates() );
1554  gp->setNaturalCoordinates(lc);
1555  //gp->setWeight (gp->giveWeight()*a/area);
1556  }
1557  }
1558 
1559  // internal test -> compute receiver area
1560  double dV, __area = 0.0;
1561  for ( int ifluid = 0; ifluid < 2; ifluid++ ) {
1562  for ( GaussPoint *gp: *integrationRulesArray [ ifluid ] ) {
1563  dV = this->computeVolumeAroundID(gp, id [ ifluid ], vcoords [ ifluid ]);
1564  // compute integral here
1565  __area += dV;
1566  }
1567  }
1568 
1569  /*
1570  * double __err = fabs(__area-area)/area;
1571  * if (__err > 1.e-6) {
1572  * OOFEM_WARNING("volume inconsistency (%5.2f\%)", __err*100);
1573  *
1574  * __area=0.0;
1575  * for (ifluid = 0; ifluid< 2; ifluid++) {
1576  * for (GaussPoint *gp: *integrationRulesArray[ifluid]) {
1577  * dV = this->computeVolumeAroundID(gp,id[ifluid], vcoords[ifluid]) ;
1578  * // compute integral here
1579  * __area += dV;
1580  * }
1581  * }
1582  * }
1583  */
1584 }
1585 
1586 
1587 double
1589 {
1590  double r1, r2, r3;
1591  double n1, n2, n3;
1592 
1593  r1 = this->giveNode(1)->giveCoordinate(1);
1594  r2 = this->giveNode(2)->giveCoordinate(1);
1595  r3 = this->giveNode(3)->giveCoordinate(1);
1596 
1597  n1 = gp->giveNaturalCoordinate(1);
1598  n2 = gp->giveNaturalCoordinate(2);
1599  n3 = 1. - n1 - n2;
1600 
1601  return n1 * r1 + n2 * r2 + n3 * r3;
1602 }
1603 
1604 
1605 
1606 double
1607 TR1_2D_SUPG2_AXI :: computeVolumeAroundID(GaussPoint *gp, integrationDomain id, const std::vector< FloatArray > &idpoly)
1608 {
1609  double weight = gp->giveWeight();
1610  double _r = computeRadiusAt(gp);
1611 
1612  if ( id == _Triangle ) {
1613  FEI2dTrLin __interpolation(1, 2);
1614  return _r *weight *fabs( __interpolation.giveTransformationJacobian ( gp->giveSubPatchCoordinates(), FEIVertexListGeometryWrapper(idpoly) ) );
1615  } else {
1616  FEI2dQuadLin __interpolation(1, 2);
1617  double det = fabs( __interpolation.giveTransformationJacobian( gp->giveSubPatchCoordinates(), FEIVertexListGeometryWrapper(idpoly) ) );
1618  return _r * det * weight;
1619  }
1620 }
1621 
1623 {
1624  _b.resize(4, 6);
1625  double _r = this->computeRadiusAt(gp);
1626 
1627 
1628  _b.at(1, 1) = b [ 0 ];
1629  _b.at(1, 2) = 0.;
1630  _b.at(1, 3) = b [ 1 ];
1631  _b.at(1, 4) = 0.;
1632  _b.at(1, 5) = b [ 2 ];
1633  _b.at(1, 6) = 0.;
1634  _b.at(2, 1) = 1. / _r;
1635  _b.at(2, 2) = 0.;
1636  _b.at(2, 3) = 1. / _r;
1637  _b.at(2, 4) = 0.;
1638  _b.at(2, 5) = 1. / _r;
1639  _b.at(2, 6) = 0.;
1640  _b.at(3, 1) = 0.0;
1641  _b.at(3, 2) = c [ 0 ];
1642  _b.at(3, 3) = 0.0;
1643  _b.at(3, 4) = c [ 1 ];
1644  _b.at(3, 5) = 0.0;
1645  _b.at(3, 6) = c [ 2 ];
1646  _b.at(4, 1) = c [ 0 ];
1647  _b.at(4, 2) = b [ 0 ];
1648  _b.at(4, 3) = c [ 1 ];
1649  _b.at(4, 4) = b [ 1 ];
1650  _b.at(4, 5) = c [ 2 ];
1651  _b.at(4, 6) = b [ 2 ];
1652 
1653  /*
1654  *
1655  * double _c = (1./3.);
1656  * double u1 = _c*(b[0]+1./_r), u2 = _c*(b[1]+1./_r), u3 = _c*(b[2]+1./_r);
1657  * double w1 = _c*c[0], w2 = _c*c[1], w3 = _c*c[2];
1658  *
1659  * //u1=u2=u3=w1=w2=w3=0.0;
1660  *
1661  * _b.at(1,1) = b[0]-u1; _b.at(1,2)=0.-w1; _b.at(1,3)= b[1]-u2; _b.at(1,4)=0.-w2; _b.at(1,5)=b[2]-u3; _b.at(1,6)=0.-w3;
1662  * _b.at(2,1) = 1./_r-u1;_b.at(2,2)=0.-w1; _b.at(2,3)= 1./_r-u2;_b.at(2,4)=0.-w2; _b.at(2,5)=1./_r-u3;_b.at(2,6)=0.-w3;
1663  * _b.at(3,1) = 0.0-u1; _b.at(3,2)=c[0]-w1; _b.at(3,3)= 0.0-u2; _b.at(3,4)=c[1]-w2; _b.at(3,5)=0.0-u3; _b.at(3,6)=c[2]-w3;
1664  * _b.at(4,1) = c[0]; _b.at(4,2)=b[0]; _b.at(4,3)= c[1]; _b.at(4,4)=b[1]; _b.at(4,5)=c[2]; _b.at(4,6)=b[2];
1665  */
1666 }
1667 
1668 /*
1669  * double
1670  * TR1_2D_SUPG2::computeVolumeAroundID(GaussPoint* gp, integrationDomain id, const Polygon& matvolpoly) {
1671  * }
1672  */
1673 void
1675 {
1677 }
1678 
1679 void
1681 // Performs end-of-step operations.
1682 {
1683  SUPGElement :: printOutputAt(file, tStep);
1684 
1685  GaussPoint *gp;
1686  if ( integrationRulesArray [ 0 ]->giveNumberOfIntegrationPoints() ) {
1687  gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
1688  } else {
1689  gp = integrationRulesArray [ 1 ]->getIntegrationPoint(0);
1690  }
1691 
1692  double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
1693  fprintf(file, "VOF %e, density %e\n\n", this->giveVolumeFraction(), rho);
1694 }
1695 
1696 
1697 void
1699  InternalStateType type, TimeStep *tStep)
1700 {
1701  GaussPoint *gp;
1702  if ( integrationRulesArray [ 0 ]->giveNumberOfIntegrationPoints() ) {
1703  gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
1704  } else {
1705  gp = integrationRulesArray [ 1 ]->getIntegrationPoint(0);
1706  }
1707 
1708  this->giveIPValue(answer, gp, type, tStep);
1709 }
1710 
1711 void
1713 {
1714  pap.resize(3);
1715  pap.at(1) = this->giveNode(1)->giveNumber();
1716  pap.at(2) = this->giveNode(2)->giveNumber();
1717  pap.at(3) = this->giveNode(3)->giveNumber();
1718 }
1719 
1720 void
1722 {
1723  answer.resize(1);
1724  if ( ( pap == this->giveNode(1)->giveNumber() ) ||
1725  ( pap == this->giveNode(2)->giveNumber() ) ||
1726  ( pap == this->giveNode(3)->giveNumber() ) ) {
1727  answer.at(1) = pap;
1728  } else {
1729  OOFEM_ERROR("node unknown");
1730  }
1731 }
1732 
1733 int
1735 { return 1; }
1736 
1737 
1740 {
1741  return SPRPatchType_2dxy;
1742 }
1743 
1744 
1745 
1746 #ifdef __OOFEG
1747 int
1749  int node, TimeStep *tStep)
1750 {
1751  /*
1752  * if (type == IST_VOFFraction) {
1753  * answer.resize(1);
1754  * answer.at(1) = this->giveTempVolumeFraction();
1755  * return 1;
1756  * } else if (type == IST_Density) {
1757  * answer.resize(1);
1758  * answer.at(1) = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
1759  * return 1;
1760  *
1761  * } else
1762  */return SUPGElement :: giveInternalStateAtNode(answer, type, mode, node, tStep);
1763 }
1764 
1765 void
1767 {
1768  WCRec p [ 3 ];
1769  GraphicObj *go;
1770 
1771  if ( !gc.testElementGraphicActivity(this) ) {
1772  return;
1773  }
1774 
1775  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
1776  EASValsSetColor( gc.getElementColor() );
1777  EASValsSetEdgeColor( gc.getElementEdgeColor() );
1778  EASValsSetEdgeFlag(true);
1779  EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
1780  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
1781  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
1782  p [ 0 ].z = 0.;
1783  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
1784  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
1785  p [ 1 ].z = 0.;
1786  p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveCoordinate(1);
1787  p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveCoordinate(2);
1788  p [ 2 ].z = 0.;
1789 
1790  go = CreateTriangle3D(p);
1791  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
1792  EGAttachObject(go, ( EObjectP ) this);
1793  EMAddGraphicsToModel(ESIModel(), go);
1794 }
1795 
1797 {
1798  int i, indx, result = 0;
1799  WCRec p [ 3 ];
1800  GraphicObj *tr;
1801  FloatArray v1, v2, v3;
1802  double s [ 3 ];
1803 
1804  if ( !gc.testElementGraphicActivity(this) ) {
1805  return;
1806  }
1807 
1808  EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
1809 
1810  // if ((gc.giveIntVarMode() == ISM_local) && (gc.giveIntVarType() == IST_VOFFraction)) {
1811  if ( ( gc.giveIntVarType() == IST_VOFFraction ) && ( gc.giveIntVarMode() == ISM_local ) ) {
1812  Polygon matvolpoly;
1813  this->formMaterialVolumePoly(matvolpoly, NULL, temp_normal, temp_p, false);
1814  EASValsSetColor( gc.getStandardSparseProfileColor() );
1815  //GraphicObj *go = matvolpoly.draw(gc,true,OOFEG_VARPLOT_PATTERN_LAYER);
1816  matvolpoly.draw(gc, true, OOFEG_VARPLOT_PATTERN_LAYER);
1817  return;
1818  }
1819 
1820  /*
1821  * if ((gc.giveIntVarMode() == ISM_local) && (gc.giveIntVarType() == IST_VOFFraction)) {
1822  * Polygon matvolpoly;
1823  * FloatArray _n;
1824  * double _p;
1825  * this->doCellDLS (_n, this->giveNumber());
1826  * this->findCellLineConstant (_p, _n, this->giveNumber());
1827  * this->formMaterialVolumePoly(matvolpoly, NULL, _n, _p, false);
1828  * GraphicObj *go = matvolpoly.draw(gc,true,OOFEG_VARPLOT_PATTERN_LAYER);
1829  * return;
1830  * }
1831  */
1832 
1833  if ( gc.giveIntVarMode() == ISM_recovered ) {
1834  result += this->giveInternalStateAtNode(v1, gc.giveIntVarType(), gc.giveIntVarMode(), 1, tStep);
1835  result += this->giveInternalStateAtNode(v2, gc.giveIntVarType(), gc.giveIntVarMode(), 2, tStep);
1836  result += this->giveInternalStateAtNode(v3, gc.giveIntVarType(), gc.giveIntVarMode(), 3, tStep);
1837  } else if ( gc.giveIntVarMode() == ISM_local ) {
1838  GaussPoint *gp;
1839  if ( integrationRulesArray [ 0 ]->giveNumberOfIntegrationPoints() ) {
1840  gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
1841  } else {
1842  gp = integrationRulesArray [ 1 ]->getIntegrationPoint(0);
1843  }
1844 
1845  result += giveIPValue(v1, gp, gc.giveIntVarType(), tStep);
1846  v2 = v1;
1847  v3 = v1;
1848  result *= 3;
1849  }
1850 
1851  if ( result != 3 ) {
1852  return;
1853  }
1854 
1855  indx = gc.giveIntVarIndx();
1856 
1857  s [ 0 ] = v1.at(indx);
1858  s [ 1 ] = v2.at(indx);
1859  s [ 2 ] = v3.at(indx);
1860 
1861  if ( gc.getScalarAlgo() == SA_ISO_SURF ) {
1862  for ( i = 0; i < 3; i++ ) {
1863  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
1864  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
1865  p [ i ].z = 0.;
1866  }
1867 
1868  //EASValsSetColor(gc.getYieldPlotColor(ratio));
1869  gc.updateFringeTableMinMax(s, 3);
1870  tr = CreateTriangleWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ]);
1871  EGWithMaskChangeAttributes(LAYER_MASK, tr);
1872  EMAddGraphicsToModel(ESIModel(), tr);
1873  } else if ( ( gc.getScalarAlgo() == SA_ZPROFILE ) || ( gc.getScalarAlgo() == SA_COLORZPROFILE ) ) {
1874  double landScale = gc.getLandScale();
1875 
1876  for ( i = 0; i < 3; i++ ) {
1877  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
1878  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
1879  p [ i ].z = s [ i ] * landScale;
1880  }
1881 
1882  if ( gc.getScalarAlgo() == SA_ZPROFILE ) {
1883  EASValsSetColor( gc.getDeformedElementColor() );
1884  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
1885  EASValsSetFillStyle(FILL_SOLID);
1886  tr = CreateTriangle3D(p);
1887  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | FILL_MASK | LAYER_MASK, tr);
1888  } else {
1889  gc.updateFringeTableMinMax(s, 3);
1890  EASValsSetFillStyle(FILL_SOLID);
1891  tr = CreateTriangleWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ]);
1892  EGWithMaskChangeAttributes(FILL_MASK | LAYER_MASK, tr);
1893  }
1894 
1895  EMAddGraphicsToModel(ESIModel(), tr);
1896  }
1897 }
1898 
1899 
1900 
1901 #endif
1902 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
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 setField(int item, InputFieldType id)
void computeNMtrx(FloatArray &answer, GaussPoint *gp)
Definition: tr1_2d_supg.C:1544
integrationDomain
Used by integrator class to supply integration points for proper domain to be integrated (Area...
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
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
integrationDomain id[2]
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
Fluid cross-section.
ScalarAlgorithmType getScalarAlgo()
virtual void computeDiffusionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes diffusion derivative terms for mass conservation equation.
#define FMElement_PrescribedTractionBC
Definition: fmelement.h:44
virtual double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
Definition: feinterpol.C:43
Abstract base class for all fluid materials.
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition: femcmpnn.h:82
Wrapper around cell with vertex coordinates stored in FloatArray**.
Definition: feinterpol.h:115
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...
double computeVolume() const
Definition: geotoolbox.C:72
const FloatArray & giveSubPatchCoordinates()
Returns local sub-patch coordinates of the receiver.
Definition: gausspoint.h:144
#define _IFT_Tr1SUPG_pvof
Definition: tr1_2d_supg.h:50
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
virtual void computeLSICStabilizationTerm_MB(FloatMatrix &answer, TimeStep *tStep)
Computes SLIC stabilization term for momentum balance equation(s).
virtual void computeAdvectionDerivativeTerm_MB(FloatMatrix &answer, TimeStep *tStep)
Computes the derivative of advection terms for momentum balance equations(s) with respect to nodal ve...
#define OOFEG_RAW_GEOMETRY_LAYER
void setPermanentVolumeFraction(double v)
Definition: leplic.h:106
virtual void formMyVolumePoly(Polygon &myPoly, LEPlic *mat_interface, bool updFlag)
Assembles receiver volume.
Base class for fluid problems.
Definition: fluidmodel.h:46
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
IntArray boundarySides
Array of boundary sides.
Definition: supgelement.h:62
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
double giveUpdatedYCoordinate(int num)
Definition: leplic.h:190
virtual void computeDiffusionTerm_MC(FloatArray &answer, TimeStep *tStep)
Computes diffusion terms for mass conservation equation.
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
EngngModel * giveEngngModel()
Returns engineering model to which receiver is associated.
Definition: domain.C:433
virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
Definition: supgelement.C:80
void computeNVector(FloatArray &answer, GaussPoint *gp)
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
virtual void computeComponentArrayAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Computes boundary condition value - its components values at given time.
Definition: load.C:82
std::vector< FloatArray > vcoords[2]
virtual void computePressureTerm_MB(FloatMatrix &answer, TimeStep *tStep)
Computes pressure terms for momentum balance equations(s).
virtual double giveCoordinate(int i)
Definition: node.C:82
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
MatResponseMode
Describes the character of characteristic material matrix.
virtual void computePressureTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes pressure terms for mass conservation equation.
Definition: tr1_2d_supg.C:557
void updateVolumePolygons(Polygon &referenceFluidPoly, Polygon &secondFluidPoly, int &rfPoints, int &sfPoints, const FloatArray &normal, const double p, bool updFlag)
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Distributed body load.
Definition: bcgeomtype.h:43
const FloatArray * getCoords() const
Definition: geotoolbox.h:81
double giveTimeIncrement()
Returns solution step associated time increment.
Definition: timestep.h:150
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
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).
void plusProductUnsym(const FloatMatrix &a, const FloatMatrix &b, double dV)
Adds to the receiver the product .
Definition: floatmatrix.C:780
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: tr1_2d_supg.C:1917
virtual double give(int aProperty, GaussPoint *gp)
Returns the value of material property &#39;aProperty&#39;.
Definition: material.C:52
InternalStateType giveIntVarType()
virtual void computeAccelerationTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes acceleration terms for mass conservation equation.
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Definition: element.C:759
double giveNaturalCoordinate(int i) const
Returns i-th natural element coordinate of receiver.
Definition: gausspoint.h:136
virtual void updateStabilizationCoeffs(TimeStep *tStep)
Updates the stabilization coefficients used for CBS and SUPG algorithms.
Class representing a 2d triangular linear interpolation based on area coordinates.
Definition: fei2dtrlin.h:44
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
virtual void computeAdvectionTerm_MC(FloatArray &answer, TimeStep *tStep)
Computes advection terms for mass conservation equation.
virtual void computeVectorOfVelocities(ValueModeType mode, TimeStep *tStep, FloatArray &velocities)
Definition: fmelement.C:48
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
virtual double giveVariableScale(VarScaleType varId)
Returns the scale factor for given variable type.
Definition: engngm.h:1087
virtual void initGeometry()
Definition: tr1_2d_supg.C:1504
int material
Number of associated material.
Definition: element.h:153
virtual SPRPatchType SPRNodalRecoveryMI_givePatchType()
#define OOFEG_RAW_GEOMETRY_WIDTH
virtual double giveWeight()
Returns integration weight of receiver.
Definition: gausspoint.h:181
void clip(Polygon &result, const Polygon &a, const Polygon &b)
Definition: geotoolbox.C:575
void times(double f)
Multiplies receiver by factor f.
Definition: floatmatrix.C:1594
virtual void computeAdvectionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes the derivative of advection terms for mass conservation equation with respect to nodal veloc...
double t_supg
Stabilization coefficients, updated for each solution step in updateStabilizationCoeffs() ...
Definition: supgelement.h:70
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
IntArray * giveBodyLoadArray()
Returns array containing load numbers of loads acting on element.
Definition: element.C:372
virtual void computeLinearAdvectionTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes the linear advection term for mass conservation equation.
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
Definition: floatarray.C:676
IntArray bodyLoadArray
Array containing indexes of loads (body loads and boundary loads are kept separately), that apply on receiver.
Definition: element.h:160
virtual void computeDiffusionTerm_MB(FloatArray &answer, TimeStep *tStep)
Computes diffusion terms for momentum balance equations(s).
Class representing the special graph constructed from two polygons that is used to perform boolean op...
Definition: geotoolbox.h:191
virtual int global2local(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Default implementation using Newton&#39;s method to find the local coordinates.
Definition: fei2dtrlin.C:103
TimeStep * givePreviousStep()
Returns pointer to previous solution step.
Definition: timestep.C:114
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 computeBCRhsTerm_MC(FloatArray &answer, TimeStep *tStep)
Computes Rhs terms due to boundary conditions.
virtual double truncateMatVolume(const Polygon &matvolpoly, double &volume)
Truncates given material polygon to receiver.
virtual double computeCriticalTimeStep(TimeStep *tStep)
Computes the critical time increment.
virtual void computeBCRhsTerm_MB(FloatArray &answer, TimeStep *tStep)
Computes Rhs terms due to boundary conditions.
virtual void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
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 computeMyVolume(LEPlic *matInterface, bool updFlag)
Computes the volume of receiver.
virtual void computeDeviatoricStress(FloatArray &answer, GaussPoint *gp, TimeStep *)
#define TRSUPG_ZERO_VOF
double computeRadiusAt(GaussPoint *gp)
Abstract base class representing Lagrangian-Eulerian (moving) material interfaces.
Definition: leplic.h:153
virtual void computeDiffusionDerivativeTerm_MB(FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep)
Computes the derivative of diffusion terms for momentum balance equations(s) with respect to nodal ve...
Class representing vector of real numbers.
Definition: floatarray.h:82
void computeBMtrx(FloatMatrix &answer, GaussPoint *gp)
virtual double computeLEPLICVolumeFraction(const FloatArray &n, const double p, LEPlic *matInterface, bool updFlag)
Computes corresponding volume fraction to given interface position.
Class representing 2D polygon.
Definition: geotoolbox.h:91
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
double computeFrobeniusNorm() const
Computes the Frobenius norm of the receiver.
Definition: floatmatrix.C:1613
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
#define _IFT_Tr1SUPG2_mat0
Definition: tr1_2d_supg.h:52
void setCoords(double x, double y)
Definition: geotoolbox.h:77
static FEI2dTrLin interp
Definition: tr1_2d_supg.h:74
virtual double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
Definition: fei2dtrlin.C:147
#define _IFT_Tr1SUPG2_mat1
Definition: tr1_2d_supg.h:53
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
virtual void computeAdvectionTerm_MB(FloatArray &answer, TimeStep *tStep)
Computes nonlinear advection terms for momentum balance equations(s).
double giveUpdatedXCoordinate(int num)
Definition: leplic.h:189
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.
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
virtual bcGeomType giveBCGeoType() const
Returns geometry character of boundary condition.
Polygon myPoly[2]
myPoly[0] Occupied by reference fluid.
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
IntArray boundaryCodes
Boundary sides codes.
Definition: supgelement.h:64
virtual int SPRNodalRecoveryMI_giveNumberOfIP()
virtual void computeAccelerationTerm_MB(FloatMatrix &answer, TimeStep *tStep)
Computes acceleration terms (generalized mass matrix with stabilization terms) for momentum balance e...
Class representing the a dynamic Input Record.
void times(double s)
Multiplies receiver with scalar.
Definition: floatarray.C:818
int mat[2]
mat[0] reference fluid mat[1] second fluid
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
Material * _giveMaterial(int indx)
Class representing vertex.
Definition: geotoolbox.h:60
double vof
Volume fraction of reference fluid in element.
Definition: leplic.h:63
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
void updateFringeTableMinMax(double *s, int size)
Class representing 2d linear triangular element for solving incompressible fluid with SUPG solver...
Definition: tr1_2d_supg.h:68
#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
Class representing a 2d isoparametric linear interpolation based on natural coordinates for quadrilat...
Definition: fei2dquadlin.h:45
#define IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:78
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
Assigns to the receiver product of .
Definition: floatmatrix.C:337
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: supgelement.C:65
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
TR1_2D_SUPG2_AXI(int n, Domain *d)
double computeVolumeAroundID(GaussPoint *gp, integrationDomain id, const std::vector< FloatArray > &idpoly)
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.
int giveNumber() const
Definition: femcmpnn.h:107
virtual void giveElementCenter(LEPlic *mat_interface, FloatArray &center, bool updFlag)
Computes the receiver center (in updated Lagrangian configuration).
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
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
#define OOFEG_VARPLOT_PATTERN_LAYER
Class representing integration point in finite element program.
Definition: gausspoint.h:93
#define OOFEM_WARNING(...)
Definition: error.h:62
IntArray boundaryLoadArray
Definition: element.h:160
FloatArray normal
Interface segment normal.
Definition: leplic.h:67
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: supgelement.C:366
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
virtual void local2global(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0
Evaluates global coordinates from given local ones.
virtual void computePressureTerm_MC(FloatMatrix &answer, TimeStep *tStep)
Computes pressure terms for mass conservation equation.
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
#define _IFT_Tr1SUPG_vof
Definition: tr1_2d_supg.h:51
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
Adds the product .
Definition: floatarray.C:226
virtual void giveDeviatoricStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
Computes the deviatoric stiffness; .

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