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

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011