NeoPZ
pzburger.cpp
Go to the documentation of this file.
1 
6 #include "pzburger.h"
7 #include "pzbndcond.h"
8 
9 using namespace std;
10 
12 
13 TPZBurger::TPZBurger(int nummat, int dim):TPZRegisterClassId(&TPZBurger::ClassId),
14 TPZMatPoisson3dReferred(nummat, dim){
15  this->fIsReferred = true;
16  this->fSolRef = 1.;
18 }
19 
22  this->fIsReferred = cp.fIsReferred;
23  this->fSolRef = cp.fSolRef;
24 }
25 
27 
28 }
29 
33  if (this->IsReferred()){
34  this->SetConvectionTerm(dsol, axes);
35  }
36 
37  int phr = phi.Rows();
38 
39 
40 
41  if(fForcingFunction) { // phi(in, 0) = phi_in
43  fForcingFunction->Execute(x,res); // dphi(i,j) = dphi_j/dxi
44  fXf = res[0];
45  }
46  REAL delx = 0.;
47  REAL ConvDirAx[3] = {0.};
48  if(fC != 0.0) {
49  int di,dj;
50  delx = 0.;
51  for(di=0; di<fDim; di++) {
52  for(dj=0; dj<fDim; dj++) {
53  delx = (delx<fabs(jacinv(di,dj))) ? fabs(jacinv(di,dj)) : delx;
54  }
55  }
56  delx = 2./delx;
57 
58 
59  switch(fDim) {
60  case 1:
61  ConvDirAx[0] = axes(0,0)*fConvDir[0]+axes(0,1)*fConvDir[1]+axes(0,2)*fConvDir[2];
62  break;
63  case 2:
64  ConvDirAx[0] = axes(0,0)*fConvDir[0]+axes(0,1)*fConvDir[1]+axes(0,2)*fConvDir[2];
65  ConvDirAx[1] = axes(1,0)*fConvDir[0]+axes(1,1)*fConvDir[1]+axes(1,2)*fConvDir[2];
66  break;
67  case 3:
68  ConvDirAx[0] = axes(0,0)*fConvDir[0]+axes(0,1)*fConvDir[1]+axes(0,2)*fConvDir[2];
69  ConvDirAx[1] = axes(1,0)*fConvDir[0]+axes(1,1)*fConvDir[1]+axes(1,2)*fConvDir[2];
70  ConvDirAx[2] = axes(2,0)*fConvDir[0]+axes(2,1)*fConvDir[1]+axes(2,2)*fConvDir[2];
71  break;
72  default:
73  PZError << "TPZBurger::Contribute dimension error " << fDim << endl;
74  }
75  }
76 
77  REAL signsol = (sol[0] < 0) ? -1. : +1.;
78  for( int in = 0; in < phr; in++ ) {
79  int kd;
80  REAL dphiic = 0;
81 
82  ef(in, 0) += - weight * ( fXf*phi(in,0) + 0.5*fSD*delx*fC*dphiic*fXf * fabs(sol[0])/fSolRef );
83  for(kd = 0; kd < fDim; kd++){
84  ef(in, 0) += -1. * weight * ( +fK * ( dphi(kd,in) * dsol(kd,0) )
85  -fC * ( ConvDirAx[kd]* dphi(kd,in) * sol[0]*sol[0]/fSolRef )
86  +0.5 * fSD * delx * fC * dphiic * dsol(kd,0) * ConvDirAx[kd] * fabs(sol[0])/fSolRef );
87  }//kd
88 
89  for( int jn = 0; jn < phr; jn++ ) {
90  ek(in, jn) += weight * (0.5*fSD*delx*fC*dphiic*fXf*phi(jn,0)*signsol/fSolRef + 0.5*fSD*delx*fC*dphiic*fXf * fabs(sol[0])/fSolRef);
91  for(kd=0; kd<fDim; kd++) {
92  ek(in,jn) += weight * (
93  +fK * ( dphi(kd,in) * dphi(kd,jn) )
94  -fC * ( ConvDirAx[kd]* dphi(kd,in) * phi(jn) * 2.*sol[0]/fSolRef )
95  +0.5 * fSD * delx * fC * dphiic * ConvDirAx[kd] * signsol * (dphi(kd,jn)*sol[0]+dsol(kd,0)*phi(jn,0))/fSolRef
96  );
97  }
98  }
99  }//in
100 
101  if (this->fC == 0.){
102  if ( !ek.VerifySymmetry() ) cout << __PRETTY_FUNCTION__ << "\nMATRIZ NAO SIMETRICA" << endl;
103  }
104 }
105 
108  if (this->IsReferred()){
109  this->SetConvectionTerm(dsol, axes);
110  }
111 
112  int phr = phi.Rows();
113 
114  if(fForcingFunction) { // phi(in, 0) = phi_in
116  fForcingFunction->Execute(x,res); // dphi(i,j) = dphi_j/dxi
117  fXf = res[0];
118  }
119  REAL delx = 0.;
120  REAL ConvDirAx[3] = {0.};
121  if(fC != 0.0) {
122  int di,dj;
123  delx = 0.;
124  for(di=0; di<fDim; di++) {
125  for(dj=0; dj<fDim; dj++) {
126  delx = (delx<fabs(jacinv(di,dj))) ? fabs(jacinv(di,dj)) : delx;
127  }
128  }
129  delx = 2./delx;
130 
131 
132  switch(fDim) {
133  case 1:
134  ConvDirAx[0] = axes(0,0)*fConvDir[0]+axes(0,1)*fConvDir[1]+axes(0,2)*fConvDir[2];
135  break;
136  case 2:
137  ConvDirAx[0] = axes(0,0)*fConvDir[0]+axes(0,1)*fConvDir[1]+axes(0,2)*fConvDir[2];
138  ConvDirAx[1] = axes(1,0)*fConvDir[0]+axes(1,1)*fConvDir[1]+axes(1,2)*fConvDir[2];
139  break;
140  case 3:
141  ConvDirAx[0] = axes(0,0)*fConvDir[0]+axes(0,1)*fConvDir[1]+axes(0,2)*fConvDir[2];
142  ConvDirAx[1] = axes(1,0)*fConvDir[0]+axes(1,1)*fConvDir[1]+axes(1,2)*fConvDir[2];
143  ConvDirAx[2] = axes(2,0)*fConvDir[0]+axes(2,1)*fConvDir[1]+axes(2,2)*fConvDir[2];
144  break;
145  default:
146  PZError << "TPZBurger::Contribute dimension error " << fDim << endl;
147  }
148  }
149 
150  REAL signsol = (sol[0] < 0) ? -1. : +1.;
151  for( int in = 0; in < phr; in++ ) {
152  int kd;
153  REAL dphiic = 0;
154  for(kd = 0; kd<fDim; kd++) dphiic += ConvDirAx[kd]*dphi(kd,in);
155  REAL norm = 0.;
156  for(kd = 0; kd<fDim; kd++) norm += (ConvDirAx[kd]*fC) * (ConvDirAx[kd]*fC);
157  norm = sqrt(norm);
158 
159  ef(in, 0) += - weight * ( fXf*phi(in,0) + 0.5*fSD*delx*fC*dphiic*fXf * fabs(sol[0])/fSolRef );
160  for(kd = 0; kd < fDim; kd++){
161  ef(in, 0) += -1. * weight * ( +fK * ( dphi(kd,in) * dsol(kd,0) )
162  -fC * ( ConvDirAx[kd]* dphi(kd,in) * sol[0]*sol[0]/fSolRef )
163  +0.5 * fSD * delx * fC * dphiic * dsol(kd,0) * ConvDirAx[kd] * fabs(sol[0])/fSolRef );
164  }//kd
165 
166  for( int jn = 0; jn < phr; jn++ ) {
167  ek(in, jn) += weight * (0.5*fSD*delx*fC*dphiic*fXf*phi(jn,0)*signsol/fSolRef);
168  for(kd=0; kd<fDim; kd++) {
169  ek(in,jn) += weight * (
170  +fK * ( dphi(kd,in) * dphi(kd,jn) )
171  -fC * ( ConvDirAx[kd]* dphi(kd,in) * phi(jn) * 2.*sol[0]/fSolRef )
172  +0.5 * fSD * delx * fC * dphiic * ConvDirAx[kd] * signsol * (dphi(kd,jn)*sol[0]+dsol(kd,0)*phi(jn,0))/fSolRef
173  );
174  }
175  }
176  }//in
177 
178  if (fC == 0.){
179  if ( !ek.VerifySymmetry() ) cout << __PRETTY_FUNCTION__ << "\nMATRIZ NAO SIMETRICA" << endl;
180  }
181 }
182 
184  REAL weight,
185  TPZFMatrix<STATE> &ek,
186  TPZFMatrix<STATE> &ef,
187  TPZBndCond &bc){
188  int numbersol = data.sol.size();
189  if(numbersol != 1)
190  {
191  DebugStop();
192  }
193  TPZFMatrix<REAL> &phi = data.phi;
194  TPZVec<STATE> &sol=data.sol[0];
195  TPZFMatrix<REAL> &axes=data.axes;
196 
197  int phr = phi.Rows();
198  short in,jn;
199  REAL v2[1];
200  v2[0] = bc.Val2()(0,0);
201 
202  switch (bc.Type()) {
203  case 0 : { // Dirichlet condition
204  for(in = 0 ; in < phr; in++) {
205  ef(in,0) += weight * ( gBigNumber * v2[0] * phi(in,0) - gBigNumber * phi(in,0) * sol[0] );
206  for (jn = 0 ; jn < phr; jn++) {
207  ek(in,jn) += gBigNumber * phi(in,0) * phi(jn,0) * weight;
208  }
209  }
210  }
211  break;
212 
213  case 1 : { // Neumann condition
214  for(in = 0 ; in < phi.Rows(); in++) {
215  ef(in,0) += v2[0] * phi(in,0) * weight;
216  }
217  }
218  break;
219 
220  case 2 :{ // condicao mista
221  cout << __PRETTY_FUNCTION__ << " at line " << __LINE__ << " not implemented\n";
222  }
223  break;
224 
225  case 3 : { // outflow condition
226 
227  int id, il, jl;
228  REAL normal[3];
229  if (fDim == 1) PZError << __PRETTY_FUNCTION__ << " - ERROR! The normal vector is not available for 1D TPZInterpolatedElement\n";
230  if (fDim == 2){
231  normal[0] = axes(0,1);
232  normal[1] = axes(1,1);
233  }
234  if (fDim == 3){
235  normal[0] = axes(0,2);
236  normal[1] = axes(1,2);
237  normal[2] = axes(2,2);
238  }
239  REAL ConvNormal = 0.;
240  for(id=0; id<fDim; id++) ConvNormal += fC*fConvDir[id]*normal[id];
241  if(ConvNormal > 0.) {
242  for(il=0; il<phr; il++) {
243  for(jl=0; jl<phr; jl++) {
244  ek(il,jl) += weight * ConvNormal * phi(il)*phi(jl)*2.*sol[0]/fSolRef;
245  }
246  ef(il,0) += -1. * weight * ConvNormal * phi(il) * sol[0]*sol[0]/fSolRef;
247  }
248  }
249  else{
250  if (ConvNormal < 0.) std::cout << "Boundary condition error: inflow detected in outflow boundary condition: ConvNormal = " << ConvNormal << "\n";
251  }
252  }
253  break;
254 
255  default :{
256  PZError << __PRETTY_FUNCTION__ << " at line " << __LINE__ << " - Error! Wrong boundary condition type\n";
257  }
258  break;
259  }
260 
261  if (this->IsSymetric()) {//only 1.e-3 because of bignumbers.
262  if ( !ek.VerifySymmetry( 1.e-3 ) ) cout << __PRETTY_FUNCTION__ << "\nMATRIZ NAO SIMETRICA" << endl;
263  }
264 }
265 
267  REAL weight,
268  TPZFMatrix<STATE> &ek,
269  TPZFMatrix<STATE> &ef){
270  int numbersol = dataleft.sol.size();
271  if (numbersol != 1) {
272  DebugStop();
273  }
274  TPZFMatrix<REAL> &dphiL = dataleft.dphix;
275  TPZFMatrix<REAL> &dphiR = dataright.dphix;
276  TPZFMatrix<REAL> &phiL = dataleft.phi;
277  TPZFMatrix<REAL> &phiR = dataright.phi;
278  TPZManVector<REAL,3> &normal = data.normal;
279  TPZVec<STATE> &solL=dataleft.sol[0];
280  TPZVec<STATE> &solR=dataright.sol[0];
281  TPZFMatrix<STATE> &dsolL=dataleft.dsol[0];
282  TPZFMatrix<STATE> &dsolR=dataright.dsol[0];
283 
284  if (this->IsReferred()){
285  this->SetConvectionTermInterface(dsolL, dsolR);
286  }
287 
288  int nrowl = phiL.Rows();
289  int nrowr = phiR.Rows();
290  int il,jl,ir,jr,id;
291 
292  //Convection term
293  REAL ConvNormal = 0.;
294  for(id=0; id<fDim; id++) ConvNormal += fC * fConvDir[id]*normal[id];
295  if(ConvNormal > 0.) {
296  for(il=0; il<nrowl; il++) {
297  ef(il, 0) += -1. * weight * ConvNormal * phiL(il) * solL[0]*solL[0]/fSolRef;
298  for(jl=0; jl<nrowl; jl++) {
299  ek(il,jl) += weight * ConvNormal * phiL(il)*phiL(jl) *2.*solL[0]/fSolRef;
300  }
301  }
302  for(ir=0; ir<nrowr; ir++) {
303  ef(ir+nrowl,0) += -1. * (-1. * weight * ConvNormal * phiR(ir) * solL[0]*solL[0]/fSolRef);
304  for(jl=0; jl<nrowl; jl++) {
305  ek(ir+nrowl,jl) -= weight * ConvNormal * phiR(ir) * phiL(jl) * 2.*solL[0]/fSolRef;
306  }
307  }
308  } else {
309  for(ir=0; ir<nrowr; ir++) {
310  ef(ir+nrowl,0) += -1. * (-1. * weight * ConvNormal * phiR(ir) * solR[0]*solR[0]/fSolRef );
311  for(jr=0; jr<nrowr; jr++) {
312  ek(ir+nrowl,jr+nrowl) -= weight * ConvNormal * phiR(ir) * phiR(jr) *2.*solR[0]/fSolRef;
313  }
314  }
315  for(il=0; il<nrowl; il++) {
316  ef(il,0) += -1. * weight * ConvNormal * phiL(il) * solR[0]*solR[0]/fSolRef;
317  for(jr=0; jr<nrowr; jr++) {
318  ek(il,jr+nrowl) += weight * ConvNormal * phiL(il) * phiR(jr) * 2.*solR[0]/fSolRef;
319  }
320  }
321  }
322 
323 
324  //diffusion term
325  REAL leftK, rightK;
326  leftK = this->fK;
327  rightK = this->fK;
328 
329  //Compute GradSol . normal
330  REAL DSolLNormal = 0.;
331  REAL DSolRNormal = 0.;
332  for(id=0; id<fDim; id++) {
333  DSolLNormal += dsolL(id,0)*normal[id];
334  DSolRNormal += dsolR(id,0)*normal[id];
335  }//for
336 
337  // 1) phi_I_left, phi_J_left
338  for(il=0; il<nrowl; il++) {
339  REAL dphiLinormal = 0.;
340  for(id=0; id<fDim; id++) {
341  dphiLinormal += dphiL(id,il)*normal[id];
342  }
343 
344  ef(il,0) += -1. * (weight * leftK * (this->fSymmetry * 0.5 * dphiLinormal*solL[0]-0.5*DSolLNormal*phiL(il,0)));
345 
346  for(jl=0; jl<nrowl; jl++) {
347  REAL dphiLjnormal = 0.;
348  for(id=0; id<fDim; id++) {
349  dphiLjnormal += dphiL(id,jl)*normal[id];
350  }
351  ek(il,jl) += weight * leftK * (
352  this->fSymmetry * 0.5*dphiLinormal*phiL(jl,0)-0.5*dphiLjnormal*phiL(il,0)
353  );
354  }
355  }
356 
357  // 2) phi_I_right, phi_J_right
358  for(ir=0; ir<nrowr; ir++) {
359  REAL dphiRinormal = 0.;
360  for(id=0; id<fDim; id++) {
361  dphiRinormal += dphiR(id,ir)*normal[id];
362  }
363 
364  //ef = F - K u
365  ef(ir+nrowl,0) += -1. * weight * rightK * ( this->fSymmetry * (-0.5 * dphiRinormal * solR[0] ) + 0.5 * DSolRNormal * phiR(ir) );
366 
367  for(jr=0; jr<nrowr; jr++) {
368  REAL dphiRjnormal = 0.;
369  for(id=0; id<fDim; id++) {
370  dphiRjnormal += dphiR(id,jr)*normal[id];
371  }
372  ek(ir+nrowl,jr+nrowl) += weight * rightK * (
373  this->fSymmetry * (-0.5 * dphiRinormal * phiR(jr) ) + 0.5 * dphiRjnormal * phiR(ir)
374  );
375  }
376  }
377 
378  // 3) phi_I_left, phi_J_right
379  for(il=0; il<nrowl; il++) {
380  REAL dphiLinormal = 0.;
381  for(id=0; id<fDim; id++) {
382  dphiLinormal += dphiL(id,il)*normal[id];
383  }
384 
385  //ef = F - K u
386  ef(il,0) += -1. * weight * ( this->fSymmetry * (-0.5 * dphiLinormal * leftK * solR[0] ) - 0.5 * DSolRNormal * rightK * phiL(il) );
387 
388  for(jr=0; jr<nrowr; jr++) {
389  REAL dphiRjnormal = 0.;
390  for(id=0; id<fDim; id++) {
391  dphiRjnormal += dphiR(id,jr)*normal[id];
392  }
393  ek(il,jr+nrowl) += weight * (
394  this->fSymmetry * (-0.5 * dphiLinormal * leftK * phiR(jr) ) - 0.5 * dphiRjnormal * rightK * phiL(il)
395  );
396  }
397  }
398 
399  // 4) phi_I_right, phi_J_left
400  for(ir=0; ir<nrowr; ir++) {
401  REAL dphiRinormal = 0.;
402  for(id=0; id<fDim; id++) {
403  dphiRinormal += dphiR(id,ir)*normal[id];
404  }
405 
406  //ef = F - K u
407  ef(ir+nrowl,0) += -1. * weight * (this->fSymmetry * 0.5 * dphiRinormal * rightK * solL[0] + 0.5 * DSolLNormal * leftK * phiR(ir));
408 
409  for(jl=0; jl<nrowl; jl++) {
410  REAL dphiLjnormal = 0.;
411  for(id=0; id<fDim; id++) {
412  dphiLjnormal += dphiL(id,jl)*normal[id];
413  }
414  ek(ir+nrowl,jl) += weight * (
415  this->fSymmetry * 0.5 * dphiRinormal * rightK * phiL(jl) + 0.5 * dphiLjnormal * leftK * phiR(ir)
416  );
417  }
418  }
419 
420  if (this->IsSymetric()){
421  if ( !ek.VerifySymmetry() ) cout << __PRETTY_FUNCTION__ << "\nMATRIZ NAO SIMETRICA" << endl;
422  }
423 }
424 
426  REAL weight,
427  TPZFMatrix<STATE> &ek,
428  TPZFMatrix<STATE> &ef,
429  TPZBndCond &bc) {
430  int numbersol = dataleft.sol.size();
431  if (numbersol != 1) {
432  DebugStop();
433  }
434  TPZFMatrix<REAL> &dphiL = dataleft.dphix;
435  TPZFMatrix<REAL> &phiL = dataleft.phi;
436  TPZManVector<REAL,3> &normal = data.normal;
437  TPZVec<STATE> &solL=dataleft.sol[0];
438  TPZFMatrix<STATE> &dsolL=dataleft.dsol[0];
439 
440  if (this->IsReferred()){
441  this->SetConvectionTermInterface(dsolL, dsolL);
442  }
443 
444  int il,jl,nrowl,id;
445  nrowl = phiL.Rows();
446  REAL ConvNormal = 0.;
447  for(id=0; id<fDim; id++) ConvNormal += fC * fConvDir[id]*normal[id];
448 
449  //Compute GradSol . normal
450  REAL DSolLNormal = 0.;
451  for(id=0; id<fDim; id++) {
452  DSolLNormal += dsolL(id,0)*normal[id];
453  }//for
454 
455  switch(bc.Type()) {
456  case 0: // DIRICHLET
457 
458  //Diffusion
459  for(il=0; il<nrowl; il++) {
460  REAL dphiLinormal = 0.;
461  for(id=0; id<fDim; id++) {
462  dphiLinormal += dphiL(id,il)*normal[id];
463  }
464  ef(il,0) += weight*fK*dphiLinormal*bc.Val2()(0,0) * this->fSymmetry;
465 
466  //ef = F - K u
467  ef(il,0) += -1. * weight*fK*(this->fSymmetry * dphiLinormal * solL[0] - DSolLNormal * phiL(il,0));
468 
469  for(jl=0; jl<nrowl; jl++) {
470  REAL dphiLjnormal = 0.;
471  for(id=0; id<fDim; id++) {
472  dphiLjnormal += dphiL(id,jl)*normal[id];
473  }
474  ek(il,jl) += weight*fK*(this->fSymmetry * dphiLinormal * phiL(jl,0) - dphiLjnormal * phiL(il,0));
475  }
476  }
477 
478  //Convection
479  if(ConvNormal > 0.) {
480  for(il=0; il<nrowl; il++) {
481  ef(il,0) += -1. * weight * ConvNormal * phiL(il) * solL[0]*solL[0]/fSolRef;
482  for(jl=0; jl<nrowl; jl++) {
483  ek(il,jl) += weight * ConvNormal * phiL(il)*phiL(jl) * 2.*solL[0]/fSolRef;
484  }
485  }
486  } else {
487  for(il=0; il<nrowl; il++) {
488  ef(il,0) -= weight * ConvNormal * bc.Val2()(0,0) * phiL(il);
489  }
490  }
491 
492  break;
493  case 1: // Neumann
494  for(il=0; il<nrowl; il++) {
495  ef(il,0) += weight*phiL(il,0)*bc.Val2()(0,0);
496  }
497  break;
498 
499  case 3: // outflow condition
500  if(ConvNormal > 0.) {
501  for(il=0; il<nrowl; il++) {
502  for(jl=0; jl<nrowl; jl++) {
503  ek(il,jl) += weight * ConvNormal * phiL(il)*phiL(jl) *2.*solL[0]/fSolRef;
504  }
505  ef(il,0) += -1. * weight * ConvNormal * phiL(il) * solL[0]*solL[0]/fSolRef;
506  }
507  }
508  else {
509  if (ConvNormal < 0.){
510  std::cout << "Boundary condition error: inflow detected in outflow boundary condition: ConvNormal = " << ConvNormal << "\n";
511  }
512  }
513  break;
514 
515  default:
516  PZError << __PRETTY_FUNCTION__ << " - Wrong boundary condition type\n";
517  break;
518  }
519  if (this->IsSymetric()){
520  if ( !ek.VerifySymmetry() ) cout << __PRETTY_FUNCTION__ << "\nMATRIZ NAO SIMETRICA" << endl;
521  }
522 }
523 
524 int TPZBurger::ClassId() const{
525  return Hash("TPZBurger") ^ TPZMatPoisson3dReferred::ClassId() << 1;
526 }
virtual void Execute(const TPZVec< REAL > &x, TPZVec< TVar > &f, TPZFMatrix< TVar > &df)
Performs function computation.
Definition: pzfunction.h:38
expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ fabs
Definition: tfadfunc.h:140
virtual void ContributeInterface(TPZMaterialData &data, TPZMaterialData &dataleft, TPZMaterialData &dataright, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef) override
It computes a contribution to stiffness matrix and load vector at one integration point...
Definition: pzburger.cpp:266
TPZManVector< REAL, 3 > normal
normal to the element at the integration point
Implements a vector class which allows to use external storage provided by the user. Utility.
Definition: pzquad.h:16
bool IsReferred()
Definition: pzburger.h:36
virtual void ContributeBCInterface(TPZMaterialData &data, TPZMaterialData &dataleft, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc) override
Definition: pzburger.cpp:425
clarg::argBool bc("-bc", "binary checkpoints", false)
STATE fSD
Multiplication value for the streamline diffusion term.
Definition: pzpoisson3d.h:53
int Type()
Definition: pzbndcond.h:249
int ClassId() const override
Unique identifier for serialization purposes.
Definition: pzburger.cpp:524
bool fIsReferred
Definition: pzburger.h:119
TPZGradSolVec dsol
vector of the derivatives of the solution at the integration point
TPZFMatrix< STATE > & Val2(int loadcase=0)
Definition: pzbndcond.h:255
TPZFNMatrix< 220, REAL > phi
vector of shapefunctions (format is dependent on the value of shapetype)
void SetConvectionTermInterface(TPZFMatrix< STATE > &dsolL, TPZFMatrix< STATE > &dsolR)
Sets convection term for ContributeInterface methods.
REAL fC
Variable which multiplies the convection term of the equation.
Definition: pzpoisson3d.h:40
TPZFNMatrix< 660, REAL > dphix
values of the derivative of the shape functions
static int gStabilizationScheme
Definition: pzburger.h:28
Contains the TPZBndCond class which implements a boundary condition for TPZMaterial objects...
int64_t size() const
Returns the number of elements of the vector.
Definition: pzvec.h:196
STATE fXf
Forcing function value.
Definition: pzpoisson3d.h:31
REAL fSolRef
Definition: pzburger.h:40
STATE fK
Coeficient which multiplies the Laplacian operator.
Definition: pzpoisson3d.h:37
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
This class defines the boundary condition for TPZMaterial objects.
Definition: pzbndcond.h:29
int64_t Rows() const
Returns number of rows.
Definition: pzmatrix.h:803
virtual void ContributeBC(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc) override
It computes a contribution to the stiffness matrix and load vector at one BC integration point...
Definition: pzburger.cpp:183
expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ sqrt
Definition: tfadfunc.h:120
TPZFNMatrix< 9, REAL > axes
axes indicating the directions of the derivatives of the shapefunctions
string res
Definition: test.py:151
Contains the TPZBurger class which implements a linear convection equation using a burger flux...
static REAL gBigNumber
Big number to penalization method, used for Dirichlet conditions.
Definition: TPZMaterial.h:83
int32_t Hash(std::string str)
Definition: TPZHash.cpp:10
virtual ~TPZBurger()
Destructor.
Definition: pzburger.cpp:26
virtual int VerifySymmetry(REAL tol=1.e-13) const
Checks if current matrix value is symmetric.
Definition: pzmatrix.cpp:1501
TPZBurger(int nummat, int dim)
Constructor with id of material and dimension of the space.
Definition: pzburger.cpp:13
TPZAutoPointer< TPZFunction< STATE > > fForcingFunction
Pointer to forcing function, it is the right member at differential equation.
Definition: TPZMaterial.h:47
void ContributeGradStab(TPZVec< REAL > &x, TPZFMatrix< REAL > &jacinv, TPZVec< STATE > &sol, TPZFMatrix< STATE > &dsol, REAL weight, TPZFMatrix< REAL > &axes, TPZFMatrix< REAL > &phi, TPZFMatrix< REAL > &dphi, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef)
Definition: pzburger.cpp:30
REAL fConvDir[3]
Direction of the convection operator.
Definition: pzpoisson3d.h:43
int fDim
Problem dimension.
Definition: pzpoisson3d.h:34
void SetConvectionTerm(TPZFMatrix< STATE > &dsol, TPZFMatrix< REAL > &axes)
Sets convection term.
void ContributeSUPG(TPZVec< REAL > &x, TPZFMatrix< REAL > &jacinv, TPZVec< STATE > &sol, TPZFMatrix< STATE > &dsol, REAL weight, TPZFMatrix< REAL > &axes, TPZFMatrix< REAL > &phi, TPZFMatrix< REAL > &dphi, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef)
Definition: pzburger.cpp:106
REAL fSymmetry
Symmetry coefficient of elliptic term.
Definition: pzpoisson3d.h:50
TPZSolVec sol
vector of the solutions at the integration point
This class implements a version of TPZMatPoisson3d where the convection term is given at each integr...
#define PZError
Defines the output device to error messages and the DebugStop() function.
Definition: pzerror.h:15
This class implements a linear convection equation using a burger flux instead of the linear flux...
Definition: pzburger.h:22
virtual int ClassId() const override
Unique identifier for serialization purposes.