OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
tr21_2d_supg.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 "tr21_2d_supg.h"
36 #include "fei2dtrquad.h"
37 #include "fei2dtrlin.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 "mathfem.h"
46 #include "fluiddynamicmaterial.h"
47 #include "fluidcrosssection.h"
48 #include "timestep.h"
49 #include "contextioerr.h"
50 #include "crosssection.h"
51 #include "classfactory.h"
52 #include "engngm.h"
53 
54 #ifdef __OOFEG
55  #include "oofeggraphiccontext.h"
56  #include "connectivitytable.h"
57 #endif
58 
59 namespace oofem {
60 REGISTER_Element(TR21_2D_SUPG);
61 
64 
65 
68  // Constructor.
69 {
70  numberOfDofMans = 6;
71 }
72 
74 // Destructor
75 { }
76 
79 {
80  return & this->velocityInterpolation;
81 }
82 
85 {
86  if ( id == P_f ) {
87  return & this->pressureInterpolation;
88  } else {
89  return & this->velocityInterpolation;
90  }
91 }
92 
93 int
95 {
96  return 15;
97 }
98 
99 void
101 {
102  if ( inode < 4 ) {
103  answer = {V_u, V_v, P_f};
104  } else {
105  answer = {V_u, V_v};
106  }
107 }
108 
109 void
111 // Sets up the array containing the four Gauss points of the receiver.
112 {
113  if ( integrationRulesArray.size() == 0 ) {
114  integrationRulesArray.resize(3);
115 
116  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 3) );
118 
119  //seven point Gauss integration
120  integrationRulesArray [ 1 ].reset( new GaussIntegrationRule(2, this, 1, 3) );
122 
123  integrationRulesArray [ 2 ].reset( new GaussIntegrationRule(3, this, 1, 3) );
125 
126 
127  //integrationRulesArray [ 3 ] = new GaussIntegrationRule(4, this, 1, 3);
128  //this->giveCrossSection()->setupIntegrationPoints( *integrationRulesArray[3], 27, this );
129  }
130 }
131 
132 
133 void
135 {
136  FloatArray n;
137 
139 
140  answer.beNMatrixOf(n, 2);
141 }
142 
143 void
145 {
146  FloatMatrix n, dn;
147  FloatArray u, un;
149  this->computeNuMatrix(n, gp);
150  this->computeVectorOfVelocities(VM_Total, tStep, un);
151 
152  u.beProductOf(n, un);
153 
154  answer.resize(2, 12);
155  answer.zero();
156  for ( int i = 1; i <= 6; i++ ) {
157  answer.at(1, 2 * i - 1) = dn.at(i, 1) * u.at(1) + dn.at(i, 2) * u.at(2);
158  answer.at(2, 2 * i) = dn.at(i, 1) * u.at(1) + dn.at(i, 2) * u.at(2);
159  }
160 }
161 
162 void
164 {
165  FloatMatrix dn;
167 
168  answer.resize(3, 12);
169  answer.zero();
170 
171  for ( int i = 1; i <= 6; i++ ) {
172  answer.at(1, 2 * i - 1) = dn.at(i, 1);
173  answer.at(2, 2 * i) = dn.at(i, 2);
174  answer.at(3, 2 * i - 1) = dn.at(i, 2);
175  answer.at(3, 2 * i) = dn.at(i, 1);
176  }
177 }
178 
179 void
181 {
182  FloatMatrix dn;
184 
185  answer.resize(1, 12);
186  answer.zero();
187 
188  for ( int i = 1; i <= 6; i++ ) {
189  answer.at(1, 2 * i - 1) = dn.at(i, 1);
190  answer.at(1, 2 * i) = dn.at(i, 2);
191  }
192 }
193 
194 void
196 {
197  FloatArray n;
199 
200  answer.resize(1, 3);
201  answer.zero();
202 
203  answer.at(1, 1) = n.at(1);
204  answer.at(1, 2) = n.at(2);
205  answer.at(1, 3) = n.at(3);
206 }
207 
208 
209 void
211 {
212  FloatArray u;
213  FloatMatrix dn, um(2, 6);
214  this->computeVectorOfVelocities(VM_Total, tStep, u);
215 
217  for ( int i = 1; i <= 6; i++ ) {
218  um.at(1, i) = u.at(2 * i - 1);
219  um.at(2, i) = u.at(2 * i);
220  }
221 
222  answer.beProductOf(um, dn);
223 }
224 
225 void
227 {
228  FloatMatrix dn;
230 
231  answer.beTranspositionOf(dn);
232 }
233 
234 
235 
236 void
238 {
239  FloatMatrix D, d2n;
240 
241  answer.resize(2, 12);
242  answer.zero();
243 
244 
246 
247  static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial()->giveDeviatoricStiffnessMatrix(D, TangentStiffness, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
248 
249  for ( int i = 1; i <= 6; i++ ) {
250  answer.at(1, 2 * i - 1) = D.at(1, 1) * d2n.at(i, 1) + D.at(1, 3) * d2n.at(i, 3) + D.at(3, 1) * d2n.at(i, 3) + D.at(3, 3) * d2n.at(i, 2);
251  answer.at(1, 2 * i) = D.at(1, 2) * d2n.at(i, 3) + D.at(1, 3) * d2n.at(i, 1) + D.at(3, 2) * d2n.at(i, 2) + D.at(3, 3) * d2n.at(i, 3);
252  answer.at(2, 2 * i - 1) = D.at(2, 1) * d2n.at(i, 3) + D.at(2, 3) * d2n.at(i, 2) + D.at(3, 1) * d2n.at(i, 1) + D.at(3, 3) * d2n.at(i, 3);
253  answer.at(2, 2 * i) = D.at(2, 2) * d2n.at(i, 2) + D.at(2, 3) * d2n.at(i, 3) + D.at(3, 2) * d2n.at(i, 3) + D.at(3, 3) * d2n.at(i, 1);
254  }
255 }
256 
257 
258 void
260 {
261  double mu_min, norm_N, norm_N_d, norm_M_d, norm_LSIC;
262  FloatMatrix N, N_d, M_d, LSIC;
263  FloatArray u;
264  mu_min = 1;
265  for ( GaussPoint *gp: *integrationRulesArray [ 1 ] ) {
266  double mu = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial()->giveEffectiveViscosity(gp, tStep);
267  if ( mu_min > mu ) {
268  mu_min = mu;
269  }
270  }
271 
272  //this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), un);
273  this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), u);
274 
275  this->computeAdvectionTerm(N, tStep);
276  this->computeAdvectionDeltaTerm(N_d, tStep);
277  this->computeMassDeltaTerm(M_d, tStep);
278  this->computeLSICTerm(LSIC, tStep);
279 
280  norm_N = N.computeFrobeniusNorm();
281  norm_N_d = N_d.computeFrobeniusNorm();
282  norm_M_d = M_d.computeFrobeniusNorm();
283  norm_LSIC = LSIC.computeFrobeniusNorm();
284 
285  if ( ( norm_N == 0 ) || ( norm_N_d == 0 ) || ( norm_M_d == 0 ) ) {
286  t_supg = 0;
287  } else {
288  //t_supg = 1. / sqrt( 1. / ( t_s1 * t_s1 ) + 1. / ( t_s2 * t_s2 ) + 1. / ( t_s3 * t_s3 ) );
289  t_supg = 0;
290  }
291 
292  if ( norm_LSIC == 0 ) {
293  t_lsic = 0;
294  } else {
295  //t_lsic = norm_N / norm_LSIC;
296 
297  t_lsic = 0;
298  }
299 
300  t_pspg = 0;
301 }
302 
303 
304 
305 void
307 {
308  FloatMatrix n, b;
309 
310  answer.clear();
311 
312  /* consistent part + supg stabilization term */
313  for ( GaussPoint *gp: *integrationRulesArray [ 1 ] ) {
314  this->computeNuMatrix(n, gp);
315  this->computeUDotGradUMatrix(b, gp, tStep);
316  double dV = this->computeVolumeAround(gp);
317  double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
318  answer.plusProductUnsym(n, b, rho * dV);
319  }
320 }
321 
322 
323 void
325 {
326  FloatMatrix n, b;
327 
328  answer.clear();
329 
330  /* consistent part + supg stabilization term */
331  for ( GaussPoint *gp: *integrationRulesArray [ 1 ] ) {
332  this->computeNuMatrix(n, gp);
333  this->computeUDotGradUMatrix(b, gp, tStep);
334  double dV = this->computeVolumeAround(gp);
335  double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
336 
337  answer.plusProductUnsym(b, b, rho * dV);
338  }
339 }
340 
341 
342 
343 void
345 {
346  FloatMatrix n, b;
347 
348  answer.clear();
349 
350  /* mtrx for computing t_supg, norm of this mtrx is computed */
351  for ( GaussPoint *gp: *integrationRulesArray [ 1 ] ) {
352  this->computeNuMatrix(n, gp);
353  this->computeUDotGradUMatrix(b, gp, tStep);
354  double dV = this->computeVolumeAround(gp);
355  double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
356 
357  answer.plusProductUnsym(b, n, rho * dV);
358  }
359 }
360 
361 void
363 {
364  FloatMatrix b;
365 
366  answer.clear();
367 
368  for ( GaussPoint *gp: *integrationRulesArray [ 1 ] ) {
369  double dV = this->computeVolumeAround(gp);
370  double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
371  this->computeDivUMatrix(b, gp);
372 
373  answer.plusProductSymmUpper(b, b, dV * rho);
374  }
375 
376  answer.symmetrized();
377 }
378 
379 
380 int
382 {
383  return 2;
384 }
385 
386 
387 double
389 {
390  return 1.e6;
391 }
392 
393 
394 double
396 {
397  double answer = 0.0, vol = 0.0;
398  FloatMatrix n, dn;
399  FloatArray fi(6), u, un, gfi;
400 
401  this->computeVectorOfVelocities(VM_Total, tStep, un);
402 
403  for ( int i = 1; i <= 6; i++ ) {
404  fi.at(i) = ls->giveLevelSetDofManValue( dofManArray.at(i) );
405  }
406 
407  for ( GaussPoint *gp: *integrationRulesArray [ 1 ] ) {
408  double dV = this->computeVolumeAround(gp);
409  velocityInterpolation.evaldNdx( dn, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
410  this->computeNuMatrix(n, gp);
411  u.beProductOf(n, un);
412  gfi.beTProductOf(dn, fi);
413 
414  vol += dV;
415  answer += dV * u.dotProduct(gfi) / gfi.computeNorm();
416  }
417 
418  return answer / vol;
419 }
420 
421 void
423 {
424  FloatMatrix dn;
425 
426  answer.clear();
427 
428  for ( GaussPoint *gp: *integrationRulesArray [ 1 ] ) {
429 
430  velocityInterpolation.evaldNdx( dn, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
431 
432  answer.add(dn);
433  }
434 }
435 
436 void
437 TR21_2D_SUPG :: LS_PCS_computeVolume(double &answer, const FloatArray **coordinates)
438 {
439  //double answer = 0.0;
440  answer = 0.0;
441 
442  for ( GaussPoint *gp: *integrationRulesArray [ 1 ] ) {
443  //answer += this->computeVolumeAround(gp);
444 
445  double determinant, weight, volume;
446 
447  determinant = fabs( this->velocityInterpolation.giveTransformationJacobian( gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) ) );
448 
449  weight = gp->giveWeight();
450  volume = determinant * weight;
451 
452  answer += volume;
453  }
454 }
455 
456 
457 double
459 {
460  double answer = 0.0;
461 
462  for ( GaussPoint *gp: *integrationRulesArray [ 1 ] ) {
463  answer += this->computeVolumeAround(gp);
464  }
465 
466  return answer;
467 }
468 
469 double
471 {
472  FloatArray fi(6), un, n;
473 
474  double vol = 0.0, eps = 0.0, _fi, S = 0.0;
475 
476  for ( int i = 1; i <= 6; i++ ) {
477  fi.at(i) = ls->giveLevelSetDofManValue( dofManArray.at(i) );
478  }
479 
480  for ( GaussPoint *gp: *integrationRulesArray [ 1 ] ) {
481  double dV = this->computeVolumeAround(gp);
482  velocityInterpolation.evalN( n, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
483  vol += dV;
484  _fi = n.dotProduct(fi);
485  S += _fi / ( _fi * _fi + eps * eps ) * dV;
486  }
487 
488  return S / vol;
489 }
490 
491 
492 
493 void
495 {
496  int neg = 0, pos = 0, zero = 0, si = 0, sqi = 0;
497 
498 
499  answer.resize(2);
500 
501  for ( double ifi: fi ) { //comparing values of fi in vertices
502  if ( ifi > 0. ) {
503  pos++;
504  } else if ( ifi < 0.0 ) {
505  neg++;
506  } else {
507  zero++;
508  }
509  }
510 
511  //control !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
512  int first_control = 0;
513  first_control = this->giveNumber();
514  OOFEM_LOG_INFO("TR21_2D_SUPG :: LS_PCS_computeVOFFractions - First control, sign(fi) in vertices, element no. %d", first_control);
515 
516 
517 
518  if ( neg == 0 ) { // all level set values positive
519  answer.at(1) = 1.0;
520  answer.at(2) = 0.0;
521  } else if ( pos == 0 ) { // all level set values negative
522  answer.at(1) = 0.0;
523  answer.at(2) = 1.0;
524  } else if ( zero == 3 ) {
525  // ???????
526  answer.at(1) = 1.0;
527  answer.at(2) = 0.0;
528  } else {
529  // zero level set inside
530  // three main cases of crossing the triangle are possible
531 
532  int inter_case, negq = 0, posq = 0, zeroq = 0; //inter_case is variable that differs type of which fi crosses the triangle
533 
534  for ( double ifi: fi ) { //loop over edge nodes
535  if ( ifi > 0. ) {
536  posq++;
537  } else if ( ifi < 0.0 ) {
538  negq++;
539  } else {
540  zeroq++;
541  }
542  }
543 
544  for ( int i = 1; i <= 3; i++ ) { //loop over vertex nodes, searching for which vertex has different value of fi
545  if ( ( pos > neg ) && ( fi.at(i) < 0.0 ) ) {
546  si = i;
547  break;
548  }
549 
550  if ( ( pos < neg ) && ( fi.at(i) >= 0.0 ) ) {
551  si = i;
552  break;
553  }
554  }
555 
556 
557 
558  for ( int i = 4; i <= 6; i++ ) { //loop over edge nodes, searching for which edge node has different value of fi
559  if ( ( posq > negq ) && ( fi.at(i) < 0.0 ) ) {
560  sqi = i;
561  break;
562  }
563 
564  if ( ( posq < negq ) && ( fi.at(i) >= 0.0 ) ) {
565  sqi = i;
566  break;
567  }
568  }
569 
570  if ( negq == 0 ) {
571  inter_case = 11; //only one node has different level set value from other five
572  } else if ( posq == 0 ) {
573  inter_case = 12; //only one node has different level set value from other five
574  } else { // level set devides element 2:4 or 3:3 nodes
575  if ( fi.at(si) * fi.at(sqi) > 0 ) {
576  inter_case = 2;
577  } else {
578  inter_case = 3;
579  }
580  }
581 
582 
583  int edge1 = 0, edge2 = 0;
584  FloatArray helplcoords(3);
585 
586 
587  if ( inter_case > 3 ) { //only one node (vertex in this case) of triangle has diffenert level set value, than others
588  FloatArray inter1(2), inter2(2);
589 
590  //kontrola!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
591  int second_control1;
592  second_control1 = this->giveNumber();
593  OOFEM_LOG_INFO("TR21_2D_SUPG :: LS_PCS_computeVOFFractions - case 1 - first type of element deviation by LS, element no. %d", second_control1);
594 
595 
596 
597 
598  if ( si == 1 ) { // computing intersection points in order to vertex with different sign of level set funct
599  this->computeIntersection(1, inter1, fi);
600 
601 
602  this->computeIntersection(3, inter2, fi);
603 
604  edge1 = 1;
605  edge2 = 3;
606  } else if ( si == 2 ) {
607  this->computeIntersection(2, inter1, fi);
608 
609 
610  this->computeIntersection(1, inter2, fi);
611 
612  edge1 = 2;
613  edge2 = 1;
614  } else if ( si == 3 ) {
615  this->computeIntersection(3, inter1, fi);
616 
617 
618  this->computeIntersection(2, inter2, fi);
619 
620  edge1 = 3;
621  edge2 = 2;
622  }
623 
624 
625 
626 
627  OOFEM_LOG_INFO("case 1 - after intersections of LS and edges, element no. %d", second_control1);
628 
629 
630 
631  //computing point on zero level set curve: [xM, yM]
632  // this point lies on line from the "si" vertex, the condition is fi([xM, yM]) = 0, M = [xM, yM]
633 
634  FloatArray _l(4), M(2), X_si(2), _Mid(2), line(6), _X1(2), Mid1(2), Mid2(2), Coeff(3), loc_Mid(3), loc_X1(3), N_Mid, N_X1;
635  double x1, xsi, y1, ysi, t, fi_X1, fi_Mid, r1, r11, r12;
636 
637  _l.at(1) = inter1.at(1);
638  _l.at(2) = inter1.at(2);
639  _l.at(3) = inter2.at(1);
640  _l.at(4) = inter2.at(2);
641 
642  this->computeCenterOf(_Mid, _l, 1);
643 
644  xsi = this->giveNode(si)->giveCoordinate(1);
645  ysi = this->giveNode(si)->giveCoordinate(2);
646 
647  X_si.at(1) = xsi;
648  X_si.at(2) = ysi;
649 
650  x1 = xsi + 2 * ( _Mid.at(1) - xsi );
651  y1 = ysi + 2 * ( _Mid.at(2) - ysi );
652 
653  _X1.at(1) = x1;
654 
655  _X1.at(2) = y1;
656 
657  this->velocityInterpolation.global2local( loc_Mid, _Mid, FEIElementGeometryWrapper(this) );
659 
660  this->velocityInterpolation.evalN( N_Mid, loc_Mid, FEIElementGeometryWrapper(this) );
661  this->velocityInterpolation.evalN( N_X1, loc_X1, FEIElementGeometryWrapper(this) );
662 
663  fi_Mid = N_Mid.dotProduct(fi);
664  fi_X1 = N_X1.dotProduct(fi);
665 
666  line.at(1) = 0;
667  line.at(2) = fi.at(si);
668  line.at(3) = sqrt( ( _Mid.at(1) - xsi ) * ( _Mid.at(1) - xsi ) + ( _Mid.at(2) - ysi ) * ( _Mid.at(2) - ysi ) );
669  line.at(4) = fi_Mid;
670  line.at(5) = sqrt( ( _Mid.at(1) - _X1.at(1) ) * ( _Mid.at(1) - _X1.at(1) ) + ( _Mid.at(2) - _X1.at(2) ) * ( _Mid.at(2) - _X1.at(2) ) );
671  line.at(6) = fi_X1;
672 
673  this->computeQuadraticFunct(Coeff, line);
674 
675  this->computeQuadraticRoots(Coeff, r11, r12);
676 
677  if ( r11 > line.at(6) || r11 < line.at(1) ) {
678  r1 = r12;
679  } else {
680  r1 = r11;
681  }
682 
683  t = r1 / sqrt( ( _Mid.at(1) - xsi ) * ( _Mid.at(1) - xsi ) + ( _Mid.at(2) - ysi ) * ( _Mid.at(2) - ysi ) );
684 
685  M.at(1) = xsi + t * ( _Mid.at(1) - xsi );
686  M.at(2) = ysi + t * ( _Mid.at(2) - ysi );
687 
688 
689  OOFEM_LOG_INFO("case 1 - after computing third point on zero level set curve inside element, element no. %d", second_control1);
690 
691 
692 
693  //computing middle points to get six-point triangle
694  FloatArray borderpoints(4);
695 
696  borderpoints.at(1) = xsi;
697  borderpoints.at(2) = ysi;
698  borderpoints.at(3) = inter1.at(1);
699  borderpoints.at(4) = inter1.at(2);
700 
701  this->computeMiddlePointOnParabolicArc(Mid1, edge1, borderpoints);
702 
703  borderpoints.at(1) = xsi;
704  borderpoints.at(2) = ysi;
705  borderpoints.at(3) = inter2.at(1);
706  borderpoints.at(4) = inter2.at(2);
707 
708 
709  this->computeMiddlePointOnParabolicArc(Mid2, edge2, borderpoints);
710 
711 
712  const FloatArray *c1 [ 6 ];
713 
714 
715 
716 
717  double vol_1, vol;
718 
719  c1 [ 0 ] = & X_si;
720  c1 [ 1 ] = & inter1;
721  c1 [ 2 ] = & inter2;
722  c1 [ 3 ] = & Mid1;
723  c1 [ 4 ] = & M;
724  c1 [ 5 ] = & Mid2;
725 
726 
727  //volume of small triangle
728  this->LS_PCS_computeVolume(vol_1, c1);
729  //volume of whole triangle
730  vol = LS_PCS_computeVolume();
731 
732 
733 
734 
735  if ( inter_case == 11 ) {
736  answer.at(2) = vol_1 / vol;
737  answer.at(1) = 1.0 - answer.at(2);
738  } else {
739  answer.at(1) = vol_1 / vol;
740  answer.at(2) = 1.0 - answer.at(1);
741  }
742  } //end case inter_case > 3
743  else if ( inter_case == 2 ) {
744  //kontrola!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
745  int second_control2 = 0;
746  second_control2 = this->giveNumber();
747  OOFEM_LOG_INFO("case 2 - second type of element deviation by LS, element no. %d", second_control2);
748 
749 
750  FloatArray inter1(2), inter2(2), crosssect(4);
751  crosssect.zero();
752 
753  if ( si == 1 ) { // computing intersection points in order to vertex with different sign of level set funct
754  this->computeIntersection(1, inter1, fi);
755 
756  this->computeIntersection(3, inter2, fi);
757 
758  edge1 = 1;
759  edge2 = 3;
760  } else if ( si == 2 ) {
761  this->computeIntersection(2, inter1, fi);
762 
763  this->computeIntersection(1, inter2, fi);
764 
765  edge1 = 2;
766  edge2 = 1;
767  } else if ( si == 3 ) {
768  this->computeIntersection(3, inter1, fi);
769 
770  this->computeIntersection(2, inter2, fi);
771 
772  edge1 = 3;
773  edge2 = 2;
774  }
775 
776  //computing point on zero level set curve: [xM, yM]
777  // this point lies on line from the "si" vertex, the condition is fi([xM, yM]) = 0
778  FloatArray X_qsi(2), X_si(2), _Mid(2), line(6), _X1(2), Mid1(2), Coeff(3), loc_Mid, loc_X1, N_Mid, N_X1, M(2);
779  double x1, xsi, y1, ysi, t, fi_X1, fi_Mid, r1, r11, r12;
780  this->computeCenterOf(_Mid, crosssect, 1);
781 
782  xsi = this->giveNode(si)->giveCoordinate(1);
783  ysi = this->giveNode(si)->giveCoordinate(2);
784 
785  X_si.at(1) = xsi;
786  X_si.at(2) = ysi;
787 
788  X_qsi.at(1) = this->giveNode(sqi)->giveCoordinate(1);
789  X_qsi.at(2) = this->giveNode(sqi)->giveCoordinate(2);
790 
791 
792  x1 = xsi + 1.3 * ( _Mid.at(1) - xsi );
793  y1 = ysi + 1.3 * ( _Mid.at(2) - ysi );
794 
795  _X1.at(1) = x1;
796  _X1.at(2) = y1;
797 
798  this->velocityInterpolation.global2local( loc_Mid, _Mid, FEIElementGeometryWrapper(this) );
800 
801  this->velocityInterpolation.evalN( N_Mid, loc_Mid, FEIElementGeometryWrapper(this) );
802  this->velocityInterpolation.evalN( N_X1, loc_X1, FEIElementGeometryWrapper(this) );
803 
804  fi_Mid = N_Mid.dotProduct(fi);
805  fi_X1 = N_X1.dotProduct(fi);
806 
807  line.at(1) = 0;
808  line.at(2) = fi.at(si);
809  line.at(3) = sqrt( ( _Mid.at(1) - xsi ) * ( _Mid.at(1) - xsi ) + ( _Mid.at(2) - ysi ) * ( _Mid.at(2) - ysi ) );
810  line.at(4) = fi_Mid;
811  line.at(5) = sqrt( ( _Mid.at(1) - _X1.at(1) ) * ( _Mid.at(1) - _X1.at(1) ) + ( _Mid.at(2) - _X1.at(2) ) * ( _Mid.at(2) - _X1.at(2) ) );
812  line.at(6) = fi_X1;
813 
814  this->computeQuadraticFunct(Coeff, line);
815 
816  this->computeQuadraticRoots(Coeff, r11, r12);
817 
818  if ( r11 > line.at(6) || r11 < line.at(1) ) {
819  r1 = r12;
820  } else {
821  r1 = r11;
822  }
823 
824  t = r1 / sqrt( ( _Mid.at(1) - xsi ) * ( _Mid.at(1) - xsi ) + ( _Mid.at(2) - ysi ) * ( _Mid.at(2) - ysi ) );
825 
826  M.at(1) = xsi + t * ( _Mid.at(1) - xsi );
827  M.at(2) = ysi + t * ( _Mid.at(2) - ysi );
828 
829  //computing middle points to get six-point triangle
830  FloatArray borderpoints(4);
831  int sub_case = 0;
832  if ( ( si == 1 ) && ( sqi == 4 ) ) {
833  borderpoints.at(1) = xsi;
834  borderpoints.at(2) = ysi;
835  borderpoints.at(3) = inter2.at(1);
836  borderpoints.at(4) = inter2.at(2);
837  sub_case = 2;
838 
839  this->computeMiddlePointOnParabolicArc(Mid1, edge2, borderpoints);
840  }
841 
842  if ( ( si == 1 ) && ( sqi == 6 ) ) {
843  borderpoints.at(1) = xsi;
844  borderpoints.at(2) = ysi;
845  borderpoints.at(3) = inter1.at(1);
846  borderpoints.at(4) = inter1.at(2);
847  sub_case = 1;
848  this->computeMiddlePointOnParabolicArc(Mid1, edge1, borderpoints);
849  }
850 
851  if ( ( si == 2 ) && ( sqi == 4 ) ) {
852  borderpoints.at(1) = xsi;
853  borderpoints.at(2) = ysi;
854  borderpoints.at(3) = inter1.at(1);
855  borderpoints.at(4) = inter1.at(2);
856  sub_case = 1;
857  this->computeMiddlePointOnParabolicArc(Mid1, edge1, borderpoints);
858  }
859 
860  if ( ( si == 2 ) && ( sqi == 5 ) ) {
861  borderpoints.at(1) = xsi;
862  borderpoints.at(2) = ysi;
863  borderpoints.at(3) = inter2.at(1);
864  borderpoints.at(4) = inter2.at(2);
865  sub_case = 2;
866  this->computeMiddlePointOnParabolicArc(Mid1, edge2, borderpoints);
867  }
868 
869  if ( ( si == 3 ) && ( sqi == 6 ) ) {
870  borderpoints.at(1) = xsi;
871  borderpoints.at(2) = ysi;
872  borderpoints.at(3) = inter2.at(1);
873  borderpoints.at(4) = inter2.at(2);
874  sub_case = 2;
875  this->computeMiddlePointOnParabolicArc(Mid1, edge2, borderpoints);
876  }
877 
878  if ( ( si == 3 ) && ( sqi == 5 ) ) {
879  borderpoints.at(1) = xsi;
880  borderpoints.at(2) = ysi;
881  borderpoints.at(3) = inter1.at(1);
882  borderpoints.at(4) = inter1.at(2);
883  sub_case = 1;
884  this->computeMiddlePointOnParabolicArc(Mid1, edge1, borderpoints);
885  }
886 
887  const FloatArray *c1 [ 6 ];
888 
889 
890 
891  double vol_1, vol;
892 
893 
894 
895  if ( sub_case == 1 ) {
896  c1 [ 0 ] = & X_si;
897  c1 [ 1 ] = & inter1;
898  c1 [ 2 ] = & inter2;
899  c1 [ 3 ] = & Mid1;
900  c1 [ 4 ] = & M;
901  c1 [ 5 ] = & X_qsi;
902  } else {
903  c1 [ 0 ] = & X_si;
904  c1 [ 1 ] = & inter1;
905  c1 [ 2 ] = & inter2;
906  c1 [ 3 ] = & X_qsi;
907  c1 [ 4 ] = & M;
908  c1 [ 5 ] = & Mid1;
909 
910 
911 
912  //volume of small triangle
913  this->LS_PCS_computeVolume(vol_1, c1);
914  //volume of whole triangle
915  vol = LS_PCS_computeVolume();
916 
917 
918 
919  if ( fi(si) < 0 ) {
920  answer.at(2) = vol_1 / vol;
921  answer.at(1) = 1.0 - answer.at(2);
922  } else {
923  answer.at(1) = vol_1 / vol;
924  answer.at(2) = 1.0 - answer.at(1);
925  }
926  } //end case inter_case == 2
927  } else { //inter_case == 3
928  //kontrola!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
929  int second_control3 = 0;
930  second_control3 = this->giveNumber();
931  OOFEM_LOG_INFO("case 3 - third type of element deviation by LS, element no. %d", second_control3);
932 
933 
934  FloatArray inter1(2), inter2(2), crosssect(4);
935  crosssect.zero();
936 
937  if ( si == 1 ) { // computing intersection points in order to vertex with different sign of level set funct
938  this->computeIntersection(1, inter1, fi);
939 
940  this->computeIntersection(3, inter2, fi);
941 
942  //edge1 = 1;
943  //edge2 = 3;
944  } else if ( si == 2 ) {
945  this->computeIntersection(2, inter1, fi);
946 
947 
948  this->computeIntersection(1, inter2, fi);
949 
950  //edge1 = 2;
951  //edge2 = 1;
952  } else if ( si == 3 ) {
953  this->computeIntersection(3, inter1, fi);
954 
955 
956  this->computeIntersection(2, inter2, fi);
957 
958  //edge1 = 3;
959  //edge2 = 2;
960  }
961 
962  //computing point on zero level set curve: [xM, yM]
963  // this point lies on line from the "si" vertex, the condition is fi([xM, yM]) = 0
964  FloatArray X_si(2), M(2), _Mid(2), line(6), _X1(2), Mid1(2), Coeff(3), loc_Mid, loc_X1, N_Mid, N_X1;
965  double x1, xsi, y1, ysi, t, fi_X1, fi_Mid, r1, r11, r12;
966  this->computeCenterOf(_Mid, crosssect, 1);
967 
968  xsi = this->giveNode(si)->giveCoordinate(1);
969  ysi = this->giveNode(si)->giveCoordinate(2);
970 
971  X_si.at(1) = xsi;
972  X_si.at(2) = ysi;
973 
974  x1 = xsi + 0.5 * ( _Mid.at(1) - xsi );
975  y1 = ysi + 0.5 * ( _Mid.at(2) - ysi );
976 
977  _X1.at(1) = x1;
978  _X1.at(2) = y1;
979 
980  this->velocityInterpolation.global2local( loc_Mid, _Mid, FEIElementGeometryWrapper(this) );
982 
983  this->velocityInterpolation.evalN( N_Mid, loc_Mid, FEIElementGeometryWrapper(this) );
984  this->velocityInterpolation.evalN( N_X1, loc_X1, FEIElementGeometryWrapper(this) );
985 
986  fi_Mid = N_Mid.dotProduct(fi);
987  fi_X1 = N_X1.dotProduct(fi);
988 
989  line.at(1) = 0;
990  line.at(2) = fi.at(si);
991  line.at(3) = sqrt( ( _Mid.at(1) - xsi ) * ( _Mid.at(1) - xsi ) + ( _Mid.at(2) - ysi ) * ( _Mid.at(2) - ysi ) );
992  line.at(4) = fi_Mid;
993  line.at(5) = sqrt( ( _Mid.at(1) - _X1.at(1) ) * ( _Mid.at(1) - _X1.at(1) ) + ( _Mid.at(2) - _X1.at(2) ) * ( _Mid.at(2) - _X1.at(2) ) );
994  line.at(6) = fi_X1;
995 
996  this->computeQuadraticFunct(Coeff, line);
997 
998  this->computeQuadraticRoots(Coeff, r11, r12);
999 
1000  if ( r11 > line.at(6) || r11 < line.at(1) ) {
1001  r1 = r12;
1002  } else {
1003  r1 = r11;
1004  }
1005 
1006  t = r1 / sqrt( ( _Mid.at(1) - xsi ) * ( _Mid.at(1) - xsi ) + ( _Mid.at(2) - ysi ) * ( _Mid.at(2) - ysi ) );
1007 
1008  M.at(1) = xsi + t * ( _Mid.at(1) - xsi );
1009  M.at(2) = ysi + t * ( _Mid.at(2) - ysi );
1010 
1011 
1012  FloatArray X_e1(2), X_e2(2); //points on edges close to the vertex si, "quadratic" nodes
1013  double vol_1, vol;
1014 
1015  X_e1.at(1) = this->giveNode(si + 3)->giveCoordinate(1);
1016  X_e1.at(2) = this->giveNode(si + 3)->giveCoordinate(2);
1017 
1018  X_e2.at(1) = this->giveNode( ( ( si + 4 ) % 3 ) + 4 )->giveCoordinate(1);
1019  X_e2.at(2) = this->giveNode( ( ( si + 4 ) % 3 ) + 4 )->giveCoordinate(2);
1020 
1021 
1022 
1023 
1024  const FloatArray *c1 [ 6 ];
1025 
1026 
1027  c1 [ 0 ] = & X_si;
1028  c1 [ 1 ] = & inter1;
1029  c1 [ 2 ] = & inter2;
1030  c1 [ 3 ] = & X_e1;
1031  c1 [ 4 ] = & M;
1032  c1 [ 5 ] = & X_e2;
1033 
1034 
1035  //volume of small triangle
1036  this->LS_PCS_computeVolume(vol_1, c1);
1037  //volume of whole triangle
1038  vol = LS_PCS_computeVolume();
1039 
1040 
1041 
1042 
1043  if ( fi(si) < 0 ) {
1044  answer.at(2) = vol_1 / vol;
1045  answer.at(1) = 1.0 - answer.at(2);
1046  } else {
1047  answer.at(1) = vol_1 / vol;
1048  answer.at(2) = 1.0 - answer.at(1);
1049  }
1050  }
1051  }
1052 }
1053 
1054 void
1056 {
1057  FloatArray Coeff(3), helplcoords(3);
1058  double fi1, fi2, fi3, r1, r11, r12;
1059  IntArray edge(3);
1060  intcoords.resize(2);
1061  intcoords.zero();
1062 
1064  fi1 = fi.at( edge.at(1) );
1065  fi2 = fi.at( edge.at(2) );
1066  fi3 = fi.at( edge.at(3) );
1067 
1068  Coeff.at(1) = fi1 + fi2 - 2 * fi3;
1069  Coeff.at(2) = fi2 - fi1;
1070  Coeff.at(3) = 2 * fi3;
1071 
1072  this->computeQuadraticRoots(Coeff, r11, r12);
1073 
1074  if ( r11 > 1.0 || r11 < -1.0 ) {
1075  r1 = r12;
1076  } else {
1077  r1 = r11;
1078  }
1079 
1080  helplcoords.zero();
1081  helplcoords.at(1) = r1;
1082 
1083  this->velocityInterpolation.edgeLocal2global( intcoords, 3, helplcoords, FEIElementGeometryWrapper(this) );
1084 
1085  //this->velocityInterpolation.evaldNdx(dn, this->giveDomain(), dofManArray, gp->giveNaturalCoordinates(),tStep->giveTime());
1086 }
1087 
1088 
1089 #if 0
1090 void
1091 TR21_2D_SUPG :: computeIntersection(int iedge, FloatArray &intcoords, FloatArray &fi)
1092 {
1093  FloatArray Coeff(3), helplcoords(3);
1094  double fi1, fi2, fi3, r1, r11, r12;
1095 
1096  intcoords.resize(2);
1097  intcoords.zero();
1098 
1100  fi1 = fi.at( edge.at(1) );
1101  fi2 = fi.at( edge.at(2) );
1102  fi3 = fi.at( edge.at(3) );
1103 
1104  Coeff.at(1) = fi1 + fi2 - 2 * fi3;
1105  Coeff.at(2) = fi2 - fi1;
1106  Coeff.at(3) = 2 * fi3;
1107 
1108  this->computeQuadraticRoots(Coeff, r11, r12);
1109 
1110  if ( r11 > 1.0 || r11 < -1.0 ) {
1111  r1 = r12;
1112  } else {
1113  r1 = r11;
1114  }
1115 
1116  helplcoords.zero();
1117  helplcoords.at(1) = r1;
1118 
1119  this->velocityInterpolation.edgeLocal2global( intcoords, iedge, this->giveDomain(), dofManArray, helplcoords, tStep->giveTime() );
1120 
1121  //this->velocityInterpolation.evaldNdx(dn, this->giveDomain(), dofManArray, gp->giveNaturalCoordinates(),tStep->giveTime());
1122 }
1123 #endif
1124 
1125 void
1127 {
1128  double a, b, c;
1129  FloatArray Coeff, C(2);
1130 
1131 
1132 
1133 
1134  this->computeQuadraticFunct(Coeff, iedge);
1135 
1136 
1137  this->computeCenterOf(C, borderpoints, 1);
1138 
1139  a = Coeff.at(1);
1140  b = Coeff.at(2);
1141  c = Coeff.at(3);
1142 
1143  answer.at(1) = C.at(1);
1144  answer.at(2) = a * C.at(1) * C.at(1) + b *C.at(1) + c;
1145 }
1146 
1147 
1148 
1149 
1150 void
1152 {
1153  switch ( dim ) {
1154  case 1:
1155 
1156  C.at(1) = ( c.at(1) + c.at(3) ) / 2;
1157  C.at(2) = ( c.at(2) + c.at(4) ) / 2;
1158 
1159  break;
1160 
1161  case 2:
1162 
1163  C.at(1) = ( c.at(1) + c.at(3) + c.at(5) ) / 3;
1164  C.at(2) = ( c.at(2) + c.at(4) + c.at(6) ) / 3;
1165 
1166  break;
1167  }
1168 }
1169 
1170 
1171 void
1173 {
1174  double a, b, c;
1175 
1176  a = Coeff.at(1);
1177  b = Coeff.at(2);
1178  c = Coeff.at(3);
1179 
1180  r1 = ( -b + sqrt(b * b - 4 * a * c) ) / ( 2 * a );
1181  r2 = ( -b - sqrt(b * b - 4 * a * c) ) / ( 2 * a );
1182 }
1183 
1184 
1185 
1186 void
1188 
1189 {
1190  IntArray edge;
1191 
1193 
1194  answer.at(1) = this->giveNode( edge.at(1) )->giveCoordinate(1);
1195  answer.at(2) = this->giveNode( edge.at(1) )->giveCoordinate(2);
1196 
1197  answer.at(3) = this->giveNode( edge.at(2) )->giveCoordinate(1);
1198  answer.at(4) = this->giveNode( edge.at(2) )->giveCoordinate(2);
1199 
1200  answer.at(5) = this->giveNode( edge.at(3) )->giveCoordinate(1);
1201  answer.at(6) = this->giveNode( edge.at(3) )->giveCoordinate(2);
1202 }
1203 
1204 
1205 
1206 //this function computes coefficients of quadratic function a,b,c in y = a*x^2 + b*x + c
1207 //iedge is number of i-th edge of quadratic triangle
1208 void
1210 {
1211  double x1, y1, x2, y2, x3, y3, detA, detA1;
1212  FloatMatrix A(3, 3), A1(3, 3);
1213  FloatArray edge, crds(6);
1214 
1215 
1216  answer.resize(3);
1217 
1218 
1219 
1220  this->computeCoordsOfEdge(crds, iedge);
1221 
1222 
1223  x1 = crds.at(1);
1224  y1 = crds.at(2);
1225  x2 = crds.at(5);
1226  y2 = crds.at(6);
1227  x3 = crds.at(3);
1228  y3 = crds.at(4);
1229 
1230  A.at(1, 1) = x1 * x1;
1231  A.at(2, 1) = x2 * x2;
1232  A.at(3, 1) = x3 * x3;
1233  A.at(1, 2) = x1;
1234  A.at(2, 2) = x2;
1235  A.at(3, 2) = x3;
1236  A.at(1, 3) = 1.0;
1237  A.at(2, 3) = 1.0;
1238  A.at(3, 3) = 1.0;
1239 
1240  FloatArray b(3);
1241 
1242  detA = A.giveDeterminant();
1243 
1244  b.at(1) = y1;
1245  b.at(2) = y2;
1246  b.at(3) = y3;
1247 
1248  for ( int k = 1; k <= 3; k++ ) {
1249  for ( int i = 1; i <= 3; i++ ) {
1250  for ( int j = 1; j <= 3; j++ ) {
1251  A1.at(i, j) = A.at(i, j);
1252  }
1253 
1254  A1.at(i, k) = b.at(i);
1255  }
1256 
1257  detA1 = A1.giveDeterminant();
1258  answer.at(k) = detA1 / detA;
1259  }
1260 }
1261 
1262 void
1264 {
1265  double x1, y1, x2, y2, x3, y3, detA, detA1;
1266  FloatMatrix A(3, 3), A1(3, 3);
1267  FloatArray edge;
1268 
1269 
1270  answer.resize(3);
1271 
1272 
1273  x1 = line.at(1);
1274  y1 = line.at(2);
1275  x2 = line.at(3);
1276  y2 = line.at(4);
1277  x3 = line.at(5);
1278  y3 = line.at(6);
1279 
1280  A.at(1, 1) = x1 * x1;
1281  A.at(2, 1) = x2 * x2;
1282  A.at(3, 1) = x3 * x3;
1283  A.at(1, 2) = x1;
1284  A.at(2, 2) = x2;
1285  A.at(3, 2) = x3;
1286  A.at(1, 3) = 1.0;
1287  A.at(2, 3) = 1.0;
1288  A.at(3, 3) = 1.0;
1289 
1290  FloatArray b(3);
1291 
1292  detA = A.giveDeterminant();
1293 
1294  b.at(1) = y1;
1295  b.at(2) = y2;
1296  b.at(3) = y3;
1297 
1298  for ( int k = 1; k <= 3; k++ ) {
1299  for ( int i = 1; i <= 3; i++ ) {
1300  for ( int j = 1; j <= 3; j++ ) {
1301  A1.at(i, j) = A.at(i, j);
1302  }
1303 
1304  A1.at(i, k) = b.at(i);
1305  }
1306 
1307  detA1 = A1.giveDeterminant();
1308  answer.at(k) = detA1 / detA;
1309  }
1310 }
1311 
1312 
1313 void
1315  InternalStateType type, TimeStep *tStep)
1316 {
1317  GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
1318  this->giveIPValue(answer, gp, type, tStep);
1319 }
1320 
1321 
1322 void
1324 { }
1325 
1326 
1327 int
1329 {
1331 }
1332 
1333 
1334 void
1336 {
1338 }
1339 
1340 int
1342 {
1343  return SUPGElement2 :: giveIPValue(answer, gp, type, tStep);
1344 }
1345 
1347 //
1348 // saves full element context (saves state variables, that completely describe
1349 // current state)
1350 //
1351 {
1352  contextIOResultType iores;
1353 
1354  if ( ( iores = SUPGElement :: saveContext(stream, mode, obj) ) != CIO_OK ) {
1355  THROW_CIOERR(iores);
1356  }
1357 
1358  return CIO_OK;
1359 }
1360 
1361 
1362 
1364 //
1365 // restores full element context (saves state variables, that completely describe
1366 // current state)
1367 //
1368 {
1369  contextIOResultType iores;
1370 
1371  if ( ( iores = SUPGElement :: restoreContext(stream, mode, obj) ) != CIO_OK ) {
1372  THROW_CIOERR(iores);
1373  }
1374 
1375  return CIO_OK;
1376 }
1377 
1378 
1379 double
1381 // Returns the portion of the receiver which is attached to gp.
1382 {
1383  double determinant, weight, volume;
1384 
1386 
1387 
1388  weight = gp->giveWeight();
1389  volume = determinant * weight;
1390 
1391  return volume;
1392 }
1393 
1394 
1395 Interface *
1397 {
1398  if ( interface == LevelSetPCSElementInterfaceType ) {
1399  return static_cast< LevelSetPCSElementInterface * >(this);
1400  } else if ( interface == ZZNodalRecoveryModelInterfaceType ) {
1401  return static_cast< ZZNodalRecoveryModelInterface * >(this);
1402  } else if ( interface == NodalAveragingRecoveryModelInterfaceType ) {
1403  return static_cast< NodalAveragingRecoveryModelInterface * >(this);
1404  }
1405 
1406  return NULL;
1407 }
1408 
1409 
1410 void
1412 {
1413  map = {1, 2, 4, 5, 7, 8, 10, 11, 12, 13, 14, 15};
1414 }
1415 
1416 void
1418 {
1419  map = {3, 6, 9};
1420 }
1421 
1422 
1423 
1424 
1425 
1426 #ifdef __OOFEG
1427 int
1429  int node, TimeStep *tStep)
1430 {
1431  return SUPGElement :: giveInternalStateAtNode(answer, type, mode, node, tStep);
1432 }
1433 
1434 
1435 
1436 void
1438 {
1439  WCRec p [ 3 ];
1440  GraphicObj *go;
1441 
1442  if ( !gc.testElementGraphicActivity(this) ) {
1443  return;
1444  }
1445 
1446  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
1447  EASValsSetColor( gc.getElementColor() );
1448  EASValsSetEdgeColor( gc.getElementEdgeColor() );
1449  EASValsSetEdgeFlag(true);
1450  EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
1451  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
1452  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
1453  p [ 0 ].z = 0.;
1454  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
1455  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
1456  p [ 1 ].z = 0.;
1457  p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveCoordinate(1);
1458  p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveCoordinate(2);
1459  p [ 2 ].z = 0.;
1460 
1461  go = CreateTriangle3D(p);
1462  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
1463  EGAttachObject(go, ( EObjectP ) this);
1464  EMAddGraphicsToModel(ESIModel(), go);
1465 }
1466 
1468 {
1469  int i, indx, result = 0;
1470  WCRec p [ 3 ];
1471  GraphicObj *tr;
1472  FloatArray v1, v2, v3;
1473  double s [ 3 ];
1474 
1475  if ( !gc.testElementGraphicActivity(this) ) {
1476  return;
1477  }
1478 
1479  EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
1480 
1481  // if ((gc.giveIntVarMode() == ISM_local) && (gc.giveIntVarType() == IST_VOFFraction)) {
1482  if ( ( gc.giveIntVarType() == IST_VOFFraction ) && ( gc.giveIntVarMode() == ISM_local ) ) {
1483  Polygon matvolpoly;
1484  //this->formMaterialVolumePoly(matvolpoly, NULL, temp_normal, temp_p, false);
1485  EASValsSetColor( gc.getStandardSparseProfileColor() );
1486  //GraphicObj *go = matvolpoly.draw(gc,true,OOFEG_VARPLOT_PATTERN_LAYER);
1487  matvolpoly.draw(gc, true, OOFEG_VARPLOT_PATTERN_LAYER);
1488  return;
1489  }
1490 
1491  if ( gc.giveIntVarMode() == ISM_recovered ) {
1492  result += this->giveInternalStateAtNode(v1, gc.giveIntVarType(), gc.giveIntVarMode(), 1, tStep);
1493  result += this->giveInternalStateAtNode(v2, gc.giveIntVarType(), gc.giveIntVarMode(), 2, tStep);
1494  result += this->giveInternalStateAtNode(v3, gc.giveIntVarType(), gc.giveIntVarMode(), 3, tStep);
1495  } else if ( gc.giveIntVarMode() == ISM_local ) {
1496  GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
1497  result += giveIPValue(v1, gp, gc.giveIntVarType(), tStep);
1498  v2 = v1;
1499  v3 = v1;
1500  result *= 3;
1501  }
1502 
1503  if ( result != 3 ) {
1504  return;
1505  }
1506 
1507  indx = gc.giveIntVarIndx();
1508 
1509  s [ 0 ] = v1.at(indx);
1510  s [ 1 ] = v2.at(indx);
1511  s [ 2 ] = v3.at(indx);
1512 
1513  EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
1514 
1515  if ( gc.getScalarAlgo() == SA_ISO_SURF ) {
1516  for ( i = 0; i < 3; i++ ) {
1517  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
1518  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
1519  p [ i ].z = 0.;
1520  }
1521 
1522  //EASValsSetColor(gc.getYieldPlotColor(ratio));
1523  gc.updateFringeTableMinMax(s, 3);
1524  tr = CreateTriangleWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ]);
1525  EGWithMaskChangeAttributes(LAYER_MASK, tr);
1526  EMAddGraphicsToModel(ESIModel(), tr);
1527  } else if ( ( gc.getScalarAlgo() == SA_ZPROFILE ) || ( gc.getScalarAlgo() == SA_COLORZPROFILE ) ) {
1528  double landScale = gc.getLandScale();
1529 
1530  for ( i = 0; i < 3; i++ ) {
1531  p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
1532  p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
1533  p [ i ].z = s [ i ] * landScale;
1534  }
1535 
1536  if ( gc.getScalarAlgo() == SA_ZPROFILE ) {
1537  EASValsSetColor( gc.getDeformedElementColor() );
1538  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
1539  EASValsSetFillStyle(FILL_SOLID);
1540  tr = CreateTriangle3D(p);
1541  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | FILL_MASK | LAYER_MASK, tr);
1542  } else {
1543  gc.updateFringeTableMinMax(s, 3);
1544  EASValsSetFillStyle(FILL_SOLID);
1545  tr = CreateTriangleWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ]);
1546  EGWithMaskChangeAttributes(FILL_MASK | LAYER_MASK, tr);
1547  }
1548 
1549  EMAddGraphicsToModel(ESIModel(), tr);
1550  }
1551 }
1552 
1553 
1554 
1555 #endif
1556 } // end namespace oofem
double giveDeterminant() const
Returns the trace (sum of diagonal components) of the receiver.
Definition: floatmatrix.C:1408
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.
The element interface required by NodalAvergagingRecoveryModel.
IntArray dofManArray
Array containing dofmanager numbers.
Definition: element.h:151
virtual int checkConsistency()
Used to check consistency and initialize some element geometry data (area,b,c).
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
virtual void giveLocalVelocityDofMap(IntArray &map)
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 FEInterpolation * giveInterpolation() const
Definition: tr21_2d_supg.C:78
virtual double evaldNdx(FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the matrix of derivatives of interpolation functions (shape functions) at given point...
Definition: fei2dtrlin.C:53
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 computeDivTauMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)
Definition: tr21_2d_supg.C:237
virtual void computeLocalEdgeMapping(IntArray &edgeNodes, int iedge)
Definition: fei2dtrquad.C:234
virtual double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
Definition: feinterpol.C:43
The element interface required by ZZNodalRecoveryModel.
virtual void computeAdvectionTerm(FloatMatrix &answer, TimeStep *tStep)
Definition: tr21_2d_supg.C:306
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
double giveLevelSetDofManValue(int i)
Returns level set value in specific node.
Definition: levelsetpcs.h:168
#define OOFEG_RAW_GEOMETRY_LAYER
Element interface for LevelSetPCS class representing level-set like material interface.
Definition: levelsetpcs.h:68
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
void computeCenterOf(FloatArray &C, FloatArray c, int dim)
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual void LS_PCS_computeVOFFractions(FloatArray &answer, FloatArray &fi)
Returns VOF fractions for each material on element according to nodal values of level set function (p...
Definition: tr21_2d_supg.C:494
static FEI2dTrLin pressureInterpolation
Definition: tr21_2d_supg.h:61
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: element.C:779
virtual void computeDivUMatrix(FloatMatrix &answer, GaussPoint *gp)
Definition: tr21_2d_supg.C:180
virtual void computeAdvectionDeltaTerm(FloatMatrix &answer, TimeStep *tStep)
Definition: tr21_2d_supg.C:324
#define S(p)
Definition: mdm.C:481
virtual void evald2Ndx2(FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the matrix of second derivatives of interpolation functions (shape functions) at given poin...
Definition: fei2dtrquad.C:81
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
void computeQuadraticFunct(FloatArray &answer, int iedge)
virtual void initGeometry()
#define THROW_CIOERR(e)
Definition: contextioerr.h:61
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
virtual int computeNumberOfDofs()
Computes or simply returns total number of element&#39;s local DOFs.
Definition: tr21_2d_supg.C:94
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
This abstract class represent a general base element class for fluid dynamic problems.
Definition: supgelement2.h:57
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: element.C:970
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
virtual void giveLocalPressureDofMap(IntArray &map)
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
void plusProductUnsym(const FloatMatrix &a, const FloatMatrix &b, double dV)
Adds to the receiver the product .
Definition: floatmatrix.C:780
InternalStateType giveIntVarType()
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
Definition: fei2dtrquad.C:43
#define OOFEM_LOG_INFO(...)
Definition: logger.h:127
virtual void computeVectorOfVelocities(ValueModeType mode, TimeStep *tStep, FloatArray &velocities)
Definition: fmelement.C:48
double dotProduct(const FloatArray &x) const
Computes the dot product (or inner product) of receiver and argument.
Definition: floatarray.C:463
virtual void computeGradUMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)
Definition: tr21_2d_supg.C:210
void computeQuadraticRoots(FloatArray Coeff, double &r1, double &r2)
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
REGISTER_Element(LSpace)
virtual void computeNpMatrix(FloatMatrix &answer, GaussPoint *gp)
Definition: tr21_2d_supg.C:195
DofIDItem
Type representing particular dof type.
Definition: dofiditem.h:86
#define OOFEG_RAW_GEOMETRY_WIDTH
virtual double giveWeight()
Returns integration weight of receiver.
Definition: gausspoint.h:181
static FEI2dTrQuad velocityInterpolation
Definition: tr21_2d_supg.h:60
double t_supg
Stabilization coefficients, updated for each solution step in updateStabilizationCoeffs() ...
Definition: supgelement.h:70
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
virtual double evaldNdx(FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the matrix of derivatives of interpolation functions (shape functions) at given point...
Definition: fei2dtrquad.C:60
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
Definition: floatarray.C:676
#define N(p, q)
Definition: mdm.C:367
void plusProductSymmUpper(const FloatMatrix &a, const FloatMatrix &b, double dV)
Adds to the receiver the product .
Definition: floatmatrix.C:698
TimeStep * givePreviousStep()
Returns pointer to previous solution step.
Definition: timestep.C:114
void computeMiddlePointOnParabolicArc(FloatArray &answer, int iedge, FloatArray borderpoints)
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix^T and anArray.
Definition: floatarray.C:708
void computeCoordsOfEdge(FloatArray &answer, int iedge)
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
virtual void computeNuMatrix(FloatMatrix &answer, GaussPoint *gp)
Definition: tr21_2d_supg.C:134
TR21_2D_SUPG(int n, Domain *aDomain)
Definition: tr21_2d_supg.C:66
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
Definition: crosssection.C:54
InternalStateMode giveIntVarMode()
void beNMatrixOf(const FloatArray &n, int nsd)
Assigns the receiver to be a repeated diagonal matrix.
Definition: floatmatrix.C:505
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual Interface * giveInterface(InterfaceType)
Interface requesting service.
Class representing 2D polygon.
Definition: geotoolbox.h:91
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
virtual void edgeLocal2global(FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates edge global coordinates from given local ones.
Definition: fei2dtrquad.C:215
double computeFrobeniusNorm() const
Computes the Frobenius norm of the receiver.
Definition: floatmatrix.C:1613
Abstract base class representing Level Set representation of material interfaces. ...
Definition: levelsetpcs.h:114
virtual void LS_PCS_computedN(FloatMatrix &answer)
Returns gradient of shape functions.
Definition: tr21_2d_supg.C:422
virtual void computeMassDeltaTerm(FloatMatrix &answer, TimeStep *tStep)
Definition: tr21_2d_supg.C:344
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
void add(const FloatMatrix &a)
Adds matrix to the receiver.
Definition: floatmatrix.C:1023
Class Interface.
Definition: interface.h:82
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: element.C:885
virtual void computeUDotGradUMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)
Definition: tr21_2d_supg.C:144
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
void beTranspositionOf(const FloatMatrix &src)
Assigns to the receiver the transposition of parameter.
Definition: floatmatrix.C:323
virtual double computeCriticalTimeStep(TimeStep *tStep)
Computes the critical time increment.
Definition: tr21_2d_supg.C:388
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
virtual double LS_PCS_computeVolume()
Returns receiver&#39;s volume.
Definition: tr21_2d_supg.C:458
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
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...
Domain * giveDomain() const
Definition: femcmpnn.h:100
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
virtual void updateStabilizationCoeffs(TimeStep *tStep)
Updates the stabilization coefficients used for CBS and SUPG algorithms.
Definition: tr21_2d_supg.C:259
virtual int giveNumberOfSpatialDimensions()
Definition: tr21_2d_supg.C:381
void updateFringeTableMinMax(double *s, int size)
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
Assigns to the receiver product of .
Definition: floatmatrix.C:337
GraphicObj * draw(oofegGraphicContext &, bool filled, int layer=OOFEG_DEBUG_LAYER)
Definition: geotoolbox.C:152
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual ~TR21_2D_SUPG()
Definition: tr21_2d_supg.C:73
void clear()
Sets size of receiver to be an empty matrix. It will have zero rows and zero columns size...
Definition: floatmatrix.h:516
virtual double LS_PCS_computeS(LevelSetPCS *ls, TimeStep *tStep)
Evaluates S in level set equation of the form where .
Definition: tr21_2d_supg.C:470
int giveNumber() const
Definition: femcmpnn.h:107
virtual int checkConsistency()
Performs consistency check.
Definition: supgelement.C:339
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
virtual void computeLSICTerm(FloatMatrix &answer, TimeStep *tStep)
Definition: tr21_2d_supg.C:362
void symmetrized()
Initializes the lower half of the receiver according to the upper half.
Definition: floatmatrix.C:1576
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
Definition: tr21_2d_supg.C:100
#define OOFEG_VARPLOT_PATTERN_LAYER
Class representing integration point in finite element program.
Definition: gausspoint.h:93
virtual int global2local(FloatArray &answer, const FloatArray &gcoords, const FEICellGeometry &cellgeo)
Default implementation using Newton&#39;s method to find the local coordinates.
Definition: feinterpol2d.C:68
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
void computeIntersection(int iedge, FloatArray &intcoords, FloatArray &fi)
Class representing solution step.
Definition: timestep.h:80
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: tr21_2d_supg.C:110
virtual void computeGradPMatrix(FloatMatrix &answer, GaussPoint *gp)
Definition: tr21_2d_supg.C:226
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
InternalStateMode
Determines the mode of internal variable.
virtual void computeBMatrix(FloatMatrix &anwer, GaussPoint *gp)
Definition: tr21_2d_supg.C:163
virtual double LS_PCS_computeF(LevelSetPCS *ls, TimeStep *tStep)
Evaluates F in level set equation of the form where for interface position driven by flow with speed...
Definition: tr21_2d_supg.C:395
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:32 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011