NeoPZ
pzplaca.cpp
Go to the documentation of this file.
1 
6 #include "pzplaca.h"
7 #include "TPZMaterial.h"
8 
9 #include "pzbndcond.h"
10 #include <math.h>
11 #include "pzvec.h"
12 #include "pzerror.h"
13 
14 using namespace std;
15 
16 TPZPlaca::TPZPlaca(int num, STATE h, STATE f, STATE E1 , STATE E2 ,
17  STATE ni1 , STATE ni2 , STATE G12 , STATE G13 ,
18  STATE G23 , TPZFMatrix<STATE> &naxes, TPZVec<STATE> &xf) :
19 TPZRegisterClassId(&TPZPlaca::ClassId),
20 TPZMaterial(num), fnaxes(naxes),
21 fE1(E1), fE2(E2), fG12(G12), fG13(G13), fG23(G23),
22 fh(h),ff(f),fmi(1./(-1.+ni1*ni2)),fni1(ni1),fni2(ni2),
23 fRmat(6,6,0.),fRmatT(6,6,0.),
24 fKxxR(6,6,0.),fKyyR(6,6,0.),fKxyR(6,6,0.),fKyxR(6,6,0.),
25 fBx0R(6,6,0.),fXF(xf) {
26 
27 
28  TPZFMatrix<STATE> Kxx(6,6,0.),Kxy(6,6,0.),Kyx(6,6,0.),Kyy(6,6,0.),
29  Bx0(6,6,0.),B0x(6,6,0.),By0(6,6,0.),B0y(6,6,0.),B00(6,6,0.),
30  B0xR(6,6,0.),By0R(6,6,0.),B0yR(6,6,0.),B00R(6,6,0.);
31  TPZFMatrix<STATE> Kn1n1(6,6,0.),Kn1n2(6,6,0.),Kn2n1(6,6,0.),Kn2n2(6,6,0.),
32  Bn10(6,6,0.),B0n1(6,6,0.),Bn20(6,6,0.),B0n2(6,6,0.),B000(6,6,0.);
33  // TPZFMatrix<STATE> fRmat(6,6),fRmatT(6,6);
34  STATE Small , k, mi;
35  Small = 1.E-5;
36  k = 5./6.; // coeficiente de cisalhamento
37  mi = 1./(-1.0 + ni1 * ni2);
38 
39  fRmat(0,0) = fnaxes(0,0); fRmat(0,1) = fnaxes(0,1); fRmat(0,2) = fnaxes(0,2);
40  fRmat(1,0) = fnaxes(1,0); fRmat(1,1) = fnaxes(1,1); fRmat(1,2) = fnaxes(1,2);
41  fRmat(2,0) = fnaxes(2,0); fRmat(2,1) = fnaxes(2,1); fRmat(2,2) = fnaxes(2,2);
42 
43  fRmat(3,3) = fnaxes(0,0); fRmat(3,4) = fnaxes(0,1); fRmat(3,5) = fnaxes(0,2);
44  fRmat(4,3) = fnaxes(1,0); fRmat(4,4) = fnaxes(1,1); fRmat(4,5) = fnaxes(1,2);
45  fRmat(5,3) = fnaxes(2,0); fRmat(5,4) = fnaxes(2,1); fRmat(5,5) = fnaxes(2,2);
46 
48 
49 
50 
51  Kxx(0,0) = -E1*h*mi;
52  Kxx(0,4) = -E1*f*h*mi;
53  Kxx(1,1) = G12*h/2. + h*Small/4.;
54  Kxx(1,3) = -f*G12*h/2.;
55  Kxx(2,2) = G13*h*k/2.;
56  Kxx(3,1) = -f*G12*h/2.;
57  Kxx(3,3) = f*f*G12*h/2. + G12*h*h*h/24.;
58  Kxx(4,0) = -E1*f*h*mi;
59  Kxx(4,4) = -E1*f*f*h*mi - E1*h*h*h*mi/12.;
60 
61  Kyx(0,1) = G12*h/2. - h*Small/4.;
62  Kyx(0,3) = -f*G12*h/2.;
63  Kyx(1,0) = -E1*h*mi*ni2;
64  Kyx(1,4) = -E1*f*h*mi*ni2;
65  Kyx(3,0) = E1*f*h*mi*ni2;
66  Kyx(3,4) = E1*f*f*h*mi*ni2 + E1*h*h*h*mi*ni2/12.;
67  Kyx(4,1) = f*G12*h/2.;
68  Kyx(4,3) = -f*f*G12*h/2. - G12*h*h*h/24.;
69 
70  Kyx.Transpose(&Kxy);
71 
72  Kyy(0,0) = G12*h/2. + h*Small/4.;
73  Kyy(0,4) = f*G12*h/2.;
74  Kyy(1,1) = -E2*h*mi;
75  Kyy(1,3) = E2*f*h*mi;
76  Kyy(2,2) = G23*h*k/2.;
77  Kyy(3,1) = E2*f*h*mi;
78  Kyy(3,3) = -E2*f*f*h*mi - E2*h*h*h*mi/12.;
79  Kyy(4,0) = f*G12*h/2.;
80  Kyy(4,4) = f*f*G12*h/2. + G12*h*h*h/24.;
81 
82  B0x(4,2) = G13*h*k/2.;
83  B0x(5,1) = -h*Small/2.;
84 
85  B0y(3,2) = -G23*h*k/2.;
86  B0y(5,0) = h*Small/2.;
87 
88 
89  B0x.Transpose(&Bx0);
90 
91  B0y.Transpose(&By0);
92 
93  B00(3,3) = G23*h*k/2.;
94  B00(4,4) = G13*h*k/2.;
95  B00(5,5) = h*Small;
96 
97  fKxxR = fRmatT * (Kxx * fRmat);
98  fKyxR = fRmatT * (Kyx * fRmat);
99  fKxyR = fRmatT * (Kxy * fRmat);
100  fKyyR = fRmatT * (Kyy * fRmat);
101  fB0xR = fRmatT * (B0x * fRmat);
102  fB0yR = fRmatT * (B0y * fRmat);
103  fBx0R = fRmatT * (Bx0 * fRmat);
104  fBy0R = fRmatT * (By0 * fRmat);
105  fB00R = fRmatT * (B00 * fRmat);
106 
107 }
108 
110  REAL weight,
111  TPZFMatrix<STATE> &ek,
112  TPZFMatrix<STATE> &ef) {
113 
114  TPZFMatrix<REAL> &dphi = data.dphix;
115  TPZFMatrix<REAL> &phi = data.phi;
116  TPZManVector<REAL,3> &x = data.x;
117  TPZFMatrix<REAL> &axes=data.axes;
118  // this method adds the contribution of the material to the stiffness
119  // matrix and right hand side
120 
121  // check on the validity of the arguments
122  //rows x cols
123  if(phi.Cols() != 1 || dphi.Rows() != 2 || phi.Rows() != dphi.Cols()){
124  PZError << "TPZPlaca.contr, inconsistent input data : phi.Cols() = "
125  << phi.Cols() << " dphi.Cols + " << dphi.Cols() <<
126  " phi.Rows = " << phi.Rows() << " dphi.Rows = " <<
127  dphi.Rows() << "\n";
128  }
129  if(fForcingFunction) {
130  fForcingFunction->Execute(x,fXF);//fXf = xfloat
131  }
132 
133  STATE Dax1n1, Dax1n2, Dax2n1, Dax2n2;
134 
135  Dax1n1 = axes(0,0)* fnaxes(0,0) + axes(0,1)* fnaxes(0,1) + axes(0,2)* fnaxes(0,2);
136  Dax1n2 = axes(0,0)* fnaxes(1,0) + axes(0,1)* fnaxes(1,1) + axes(0,2)* fnaxes(1,2);
137  Dax2n1 = axes(1,0)* fnaxes(0,0) + axes(1,1)* fnaxes(0,1) + axes(1,2)* fnaxes(0,2);
138  Dax2n2 = axes(1,0)* fnaxes(1,0) + axes(1,1)* fnaxes(1,1) + axes(1,2)* fnaxes(1,2);
139 
140  TPZFMatrix<STATE> Kn1n1(6,6),Kn1n2(6,6),Kn2n1(6,6),Kn2n2(6,6),
141  Bn10(6,6),B0n1(6,6),Bn20(6,6),B0n2(6,6),B000(6,6);
142 
143 
144  Kn1n1 = fKxxR * Dax1n1 * Dax1n1 + fKyyR * Dax1n2 * Dax1n2 +
145  fKxyR * Dax1n1 * Dax1n2 + fKyxR * Dax1n1 * Dax1n2 ;
146 
147  Kn2n1 = fKxxR * Dax2n1 * Dax1n1 + fKyyR * Dax2n2 * Dax1n2 +
148  fKyxR * Dax2n2 * Dax1n1 + fKxyR * Dax2n1 * Dax1n2 ;
149 
150  Kn1n2 = fKxyR * Dax1n1 * Dax2n2 + fKxxR * Dax1n1 * Dax2n1 +
151  fKyyR * Dax1n2 * Dax2n2 + fKyxR * Dax1n2 * Dax2n1;
152 
153  Kn2n2 = fKyyR * Dax2n2 * Dax2n2 + fKxxR * Dax2n1 * Dax2n1 +
154  fKxyR * Dax2n1 * Dax2n2 + fKyxR * Dax2n2 * Dax2n1;
155 
156  B0n1 = fB0xR * Dax1n1 + fB0yR * Dax1n2 ;
157 
158  B0n2 = fB0yR * Dax2n2 + fB0xR * Dax2n1 ;
159 
160  Bn10 = fBx0R * Dax1n1 + fBy0R * Dax1n2 ;
161 
162  Bn20 = fBy0R * Dax2n2 + fBx0R * Dax2n1 ;
163 
164  B000 = fB00R ;
165 
166  int idf,jdf,i,j;
167  int nshape = phi.Rows();
168  TPZFMatrix<STATE> KIJ(6,6);
169 
170  for(i=0; i<nshape; i++) {
171  STATE dphi_0i(dphi(0,i)),dphi_1i(dphi(1,i)),phi_i(phi(i,0));
172  if(fExactFunction) {
173  TPZFMatrix<STATE> u(6,1),du(2,6),Fun(6,1),dux(6,1),duy(6,1);
174  fExactFunction(axes,x,u,du);
175  for(int k=0;k<6;k++){
176  dux(k,0) = du(0,k);
177  duy(k,0) = du(1,k);
178  }
179  Fun = (dphi_0i*(Kn1n1*dux)+
180  dphi_0i*(Kn1n2*duy)+
181  dphi_1i*(Kn2n1*dux)+
182  dphi_1i*(Kn2n2*duy)+
183  dphi_0i*( Bn10* u )+
184  dphi_1i*( Bn20* u )+
185  phi_i *(B0n1 *dux)+
186  phi_i *(B0n2 *duy)+
187  phi_i *( B000* u ));
188  for(idf=0; idf<6; idf++) ef(6*i+idf,0) += weight*Fun(idf,0)*phi(i,0);
189  }
190 
191  for(idf=0; idf<6; idf++) ef(6*i+idf,0) += weight*fXF[idf]*phi(i,0);
192 
193  for(j=0; j<nshape; j++) {
194  STATE dphi_0j(dphi(0,j)),dphi_1j(dphi(1,j)),phi_j(phi(j,0));
195  KIJ = (STATE(weight))*(dphi_0i*dphi_0j*Kn1n1+
196  dphi_0i*dphi_1j*Kn1n2+
197  dphi_1i*dphi_0j*Kn2n1+
198  dphi_1i*dphi_1j*Kn2n2+
199  dphi_0i*phi_j *Bn10 +
200  dphi_1i*phi_j *Bn20 +
201  phi_i *dphi_0j*B0n1 +
202  phi_i *dphi_1j*B0n2 +
203  phi_i *phi_j *B000 );
204  for(idf=0; idf<6; idf++) for(jdf=0; jdf<6; jdf++)
205  ek(i*6+idf,j*6+jdf) += KIJ(idf,jdf);
206  }
207  }
208 }
209 
211  REAL weight,
212  TPZFMatrix<STATE> &ek,
213  TPZFMatrix<STATE> &ef,
214  TPZBndCond &bc) {
215 
216  TPZFMatrix<REAL> &phi = data.phi;
217 
218  if(bc.Material() != this){
219  PZError << "TPZMat1dLin.apply_bc warning : this material didn't create the boundary condition!\n";
220  }
221 
222  if(bc.Type() < 0 && bc.Type() > 2){
223  PZError << "TPZMat1dLin.aplybc, unknown boundary condition type :" <<
224  bc.Type() << " boundary condition ignored\n";
225  }
226 
227  int numdof = NStateVariables();
228  int numnod = ek.Rows()/numdof;
229  int r = numdof;
230 
231  int idf,jdf,in,jn;
232  switch(bc.Type()){
233 
234  case 0:
235  for(in=0 ; in<numnod ; ++in){
236  for(idf = 0;idf<r;idf++) {
237  (ef)(in*r+idf,0) += gBigNumber*phi(in,0)*bc.Val2()(idf,0)*weight;
238  }
239  for(jn=0 ; jn<numnod ; ++jn) {
240  for(idf = 0;idf<r;idf++) {
241  ek(in*r+idf,jn*r+idf) += gBigNumber*phi(in,0)*phi(jn,0)*weight;
242  }
243  }
244  }
245  break;
246 
247  case 1:
248  for(in=0 ; in<numnod ; ++in){
249  for(idf = 0;idf<r;idf++) {
250  (ef)(in*r+idf,0) += phi(in,0)*bc.Val2()(idf,0)*weight;
251  }
252  }
253  break;
254 
255  case 2:
256  for(in=0 ; in<numnod ; ++in){
257  for(idf = 0;idf<r;idf++) {
258  (ef)(in*r+idf,0) += phi(in,0)*bc.Val2()(idf,0)*weight;
259  }
260  for(jn=0 ; jn<numnod ; ++jn) {
261  for(idf = 0;idf<r;idf++) {
262  for(jdf = 0;jdf<r;jdf++) {
263  ek(in*r+idf,jn*r+jdf) += bc.Val1()(idf,jdf)*phi(in,0)*phi(jn,0)*weight;
264  }
265  }
266  }
267  }//fim switch
268  }
269 }
270 
271 int TPZPlaca::NFluxes() {return 1;}
272 
274  PZError << "TPZPlaca::Flux is called\n";
275 }
276 
277 void TPZPlaca::Print(std::ostream & out) {
278  //out << "Material type TPZPlaca -- number = " << Id() << "\n";
279  //out << "Matrix xk -> "; fXk.Print("fXk",out);
280  //out << "Matrix xC -> "; fCk.Print("fCf",out);
281  //out << "Matrix xf -> "; fXf.Print("fXf",out);
282 }
283 
285 int TPZPlaca::VariableIndex(const std::string &name){
286  if(!strcmp(name.c_str(),"Deslocx")) return 2;// Desloc. eixo x global
287  if(!strcmp(name.c_str(),"Deslocy")) return 3;// Desloc. eixo y global
288  if(!strcmp(name.c_str(),"Deslocz")) return 4;// Desloc. eixo z global
289  if(!strcmp(name.c_str(),"Mn1")) return 5;// Mom. fletor eixo n1 da fibra
290  if(!strcmp(name.c_str(),"Mn2")) return 6;// Mom. fletor eixo n2 da fibra
291  if(!strcmp(name.c_str(),"Mn1n2")) return 7;// Mom. volvente eixos n1 e n2 da fibra
292  if(!strcmp(name.c_str(),"Sign1")) return 8;// tens� normal na dire�o n1
293  if(!strcmp(name.c_str(),"Sign2")) return 9;// tens� normal na dire�o n2
294  if(!strcmp(name.c_str(),"Taun1n2")) return 10;// tens� cisalhamento eixos n1 e n2
295  if(!strcmp(name.c_str(),"Taun1n3")) return 11;// tens� cisalhamento eixos n1 e n3
296  if(!strcmp(name.c_str(),"Taun2n3")) return 12;// tens� cisalhamento eixos n2 e n2
297  if(!strcmp(name.c_str(),"Displacement")) return 13;// deslocamento x,y,z
298  int var;
299  cout << "TPZPlaca name not found " << name << endl;
300  cout.flush();
301  cin >> var;
302 
303  return TPZMaterial::VariableIndex(name);
304 }
305 
309  if(var == 2) return 1;
310  if(var == 3) return 1;
311  if(var == 4) return 1;
312  if(var == 5) return 1;
313  if(var == 6) return 1;
314  if(var == 7) return 1;
315  if(var == 8) return 1;
316  if(var == 9) return 1;
317  if(var == 10) return 1;
318  if(var == 11) return 1;
319  if(var == 12) return 1;
320  if(var == 13) return 3;
321 
323 }
324 
327  TPZFMatrix<REAL> &axes,int var,TPZVec<STATE> &Solout){
328 
329 
330  if(var == 2) {
331  Solout.Resize(1);
332  Solout[0] = Sol[0];
333  return;
334  }
335  if(var == 3) {
336  Solout.Resize(1);
337  Solout[0] = Sol[1];
338  return;
339  }
340  if(var == 4) {
341  Solout.Resize(1);
342  Solout[0] = Sol[2];
343  return;
344  }
345  if(var == 13) {
346  Solout.Resize(3);
347  Solout[0]=Sol[0];
348  Solout[1]=Sol[1];
349  Solout[2]=Sol[2];
350  return;
351  }
352  if(var > 4) {
353  Solout.Resize(1);
354  TPZVec<REAL> Soln(6);
355  TPZFMatrix<REAL> DSolnax(2,6),DSolnn(2,6);
356  int idf,jdf;
357  for(idf=0; idf<6; idf++) {
358  Soln[idf] = 0;
359  DSolnax(0,idf) = 0;
360  DSolnax(1,idf) = 0.;
361  for(jdf=0; jdf<6; jdf++) {
362  Soln[idf] += fRmat(idf,jdf)*Sol[jdf];
363  DSolnax(0,idf) += fRmat(idf,jdf)*DSol(0,jdf);
364  DSolnax(1,idf) += fRmat(idf,jdf)*DSol(1,jdf);
365  }
366  }
367  REAL Dax1n1, Dax1n2, Dax2n1, Dax2n2;
368 
369  Dax1n1 = axes(0,0)* fnaxes(0,0) + axes(0,1)* fnaxes(0,1) + axes(0,2)* fnaxes(0,2);
370  Dax1n2 = axes(0,0)* fnaxes(1,0) + axes(0,1)* fnaxes(1,1) + axes(0,2)* fnaxes(1,2);
371  Dax2n1 = axes(1,0)* fnaxes(0,0) + axes(1,1)* fnaxes(0,1) + axes(1,2)* fnaxes(0,2);
372  Dax2n2 = axes(1,0)* fnaxes(1,0) + axes(1,1)* fnaxes(1,1) + axes(1,2)* fnaxes(1,2);
373 
374  for(idf=0;idf<6;idf++) {
375  DSolnn(0,idf) = Dax1n1*DSolnax(0,idf)+Dax2n1*DSolnax(1,idf);
376  DSolnn(1,idf) = Dax1n2*DSolnax(0,idf)+Dax2n2*DSolnax(1,idf);
377  }
378 
379  if (var == 5) {
380  REAL Mn1;
381  Mn1 = -fE1*fh*fh*fh*fmi*(-fni2*DSolnn(1,3)+DSolnn(0,4))/12.;
382  Mn1 += fE1*ff*fh*fmi*(ff*fni2*DSolnn(1,3)-fni2*DSolnn(1,1)
383  -ff*DSolnn(0,4)-DSolnn(0,0));
384  Solout[0] = -Mn1;
385  return;
386  }
387 
388  if (var == 6) {
389  REAL Mn2;
390  Mn2 =-(fE2*fh*fh*fh*fmi*(-DSolnn(1,3) + fni1*DSolnn(0,4)))/12.;
391  Mn2 += fE2*ff*fh*fmi*(ff*DSolnn(1,3) - DSolnn(1,1) -
392  ff*fni1*DSolnn(0,4) - fni1*DSolnn(0,0));
393  Solout[0] = -Mn2;
394  return;
395  }
396 
397  if (var == 7 ) {
398  REAL Mn1n2;
399  Mn1n2 = fG12*fh*fh*fh*(DSolnn(1,4) - DSolnn(0,3))/24.;
400  Mn1n2 += (ff*fG12*fh*(ff*DSolnn(1,4) + DSolnn(1,0) -
401  ff*DSolnn(0,3) + DSolnn(0,1)))/2.;
402  Solout[0] = -Mn1n2;
403  return;
404  }
405  if (var == 8 ) {
406  REAL Sign1;
407  REAL z;
408  z=0.;
409  Sign1 = -fE1*fmi*(fni2*(-(ff + z)*DSolnn(1,3) + DSolnn(1,1))+
410  (ff + z)*DSolnn(0,4) + DSolnn(0,0));
411  Solout[0] = Sign1;
412  return;
413  }
414  if (var == 9 ) {
415  REAL Sign2;
416  REAL z;
417  z=0.;
418  Sign2=-fE2*fmi*(-(ff + z)*DSolnn(1,3) + DSolnn(1,1) +
419  fni1*((ff + z)*DSolnn(0,4) + DSolnn(0,0)));
420  Solout[0] = Sign2;
421  return;
422  }
423  if (var == 10 ) {
424  REAL Taun1n2;
425  REAL z;
426  z=0.;
427  Taun1n2=(fG12*((ff + z)*DSolnn(1,4) + DSolnn(1,0) -
428  (ff + z)*DSolnn(0,3) + DSolnn(0,1)))/2.;
429  Solout[0] = Taun1n2;
430  return;
431  }
432  if (var == 11 ) {
433  REAL Taun1n3;
434  Taun1n3=fG13*(Soln[4] + DSolnn(0,2))/2.;
435  Solout[0] = Taun1n3;
436  return;
437  }
438  if (var == 12 ) {
439  REAL Taun2n3;
440  Taun2n3=fG23*(-Soln[3] + DSolnn(1,2))/2.;
441  Solout[0] = Taun2n3;
442  return;
443  }
444 
445  }
446 
447  TPZMaterial::Solution(Sol,DSol,axes,var,Solout);
448 }
449 
451  TPZVec<STATE> &u_exact,TPZFMatrix<STATE> &du_exact,TPZVec<REAL> &values) {
452 
453  //ENERGY NORM
454  //solution error:
455  int ndof = NStateVariables();
456  TPZFMatrix<STATE> err(ndof,1),dxerr(ndof,1),dyerr(ndof,1);
457  int i;
458  // u - uh
459  for(i=0;i<6;i++) err(i,0) = u_exact[i] - u[i];
460  // du/dx - duh/dx
461  for(i=0;i<6;i++){
462  dxerr(i,0) = du_exact(0,i) - dudx(0,i);
463  dyerr(i,0) = du_exact(1,i) - dudx(1,i);
464  }
465  // error transpose
466  TPZFMatrix<STATE> errt(1,ndof),dxerrt(1,ndof),dyerrt(1,ndof);
467  err.Transpose(&errt);
468  dxerr.Transpose(&dxerrt);
469  dyerr.Transpose(&dyerrt);
470  TPZFMatrix<STATE> ENERGY(1,1);
471 
472  //ENERGY norm calculation
473  ENERGY = dxerrt * (fKxxR * dxerr) + dyerrt * (fKyyR * dyerr);
474  ENERGY += dxerrt * (fKxyR * dyerr) + dyerrt * (fKyxR * dxerr);
475  ENERGY += errt * (fB0xR * dxerr) + dxerrt * (fBx0R * err);
476  ENERGY += errt * (fB0yR * dyerr) + dyerrt * (fBy0R * err);
477  ENERGY += errt * (fB00R * err);
478  values[0] = ENERGY(0,0);
479 
480 }
481 
482 int TPZPlaca::ClassId() const{
483  return Hash("TPZPlaca") ^ TPZMaterial::ClassId() << 1;
484 }
virtual void Solution(TPZMaterialData &data, int var, TPZVec< STATE > &Solout)
Returns the solution associated with the var index based on the finite element approximation.
virtual void Execute(const TPZVec< REAL > &x, TPZVec< TVar > &f, TPZFMatrix< TVar > &df)
Performs function computation.
Definition: pzfunction.h:38
STATE fmi
Definition: pzplaca.h:24
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: pzplaca.cpp:210
TPZFMatrix< STATE > fBy0R
Definition: pzplaca.h:26
STATE fE1
Definition: pzplaca.h:24
virtual int NFluxes() override
Returns the number of components which form the flux function.
Definition: pzplaca.cpp:271
TPZManVector< REAL, 3 > x
value of the coordinate at the integration point
clarg::argBool bc("-bc", "binary checkpoints", false)
STATE fG12
Definition: pzplaca.h:24
Contains the TPZPlaca class.
TPZFMatrix< STATE > fKyyR
Definition: pzplaca.h:26
int Type()
Definition: pzbndcond.h:249
Templated vector implementation.
TPZFMatrix< STATE > fKxxR
Definition: pzplaca.h:26
virtual int VariableIndex(const std::string &name)
Returns the variable index associated with the name.
Defines PZError.
TPZFMatrix< STATE > fB0xR
Definition: pzplaca.h:26
DESCRIBE PLEASE.
Definition: pzplaca.h:21
STATE fG23
Definition: pzplaca.h:24
clarg::argBool h("-h", "help message", false)
TPZFMatrix< STATE > fKxyR
Definition: pzplaca.h:26
STATE fni1
Definition: pzplaca.h:24
virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef) override
It computes a contribution to the stiffness matrix and load vector at one integration point...
Definition: pzplaca.cpp:109
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)
TPZFNMatrix< 660, REAL > dphix
values of the derivative of the shape functions
TPZPlaca(int num, STATE h, STATE f, STATE E1, STATE E2, STATE ni1, STATE ni2, STATE G12, STATE G13, STATE G23, TPZFMatrix< STATE > &naxes, TPZVec< STATE > &xf)
Definition: pzplaca.cpp:16
TPZFMatrix< STATE > fBx0R
Definition: pzplaca.h:26
This abstract class defines the behaviour which each derived class needs to implement.
Definition: TPZMaterial.h:39
Contains the TPZBndCond class which implements a boundary condition for TPZMaterial objects...
STATE ff
Definition: pzplaca.h:24
virtual void Resize(const int64_t newsize, const T &object)
Resizes the vector object reallocating the necessary storage, copying the existing objects to the new...
Definition: pzvec.h:373
f
Definition: test.py:287
STATE fni2
Definition: pzplaca.h:24
virtual int VariableIndex(const std::string &name) override
Definition: pzplaca.cpp:285
This class defines the boundary condition for TPZMaterial objects.
Definition: pzbndcond.h:29
virtual void Errors(TPZVec< REAL > &x, TPZVec< STATE > &u, TPZFMatrix< STATE > &dudx, TPZFMatrix< REAL > &axes, TPZVec< STATE > &flux, TPZVec< STATE > &u_exact, TPZFMatrix< STATE > &du_exact, TPZVec< REAL > &values) override
Computes the error due to the difference between the interpolated flux and the flux computed based o...
Definition: pzplaca.cpp:450
TPZFMatrix< STATE > fB0yR
Definition: pzplaca.h:26
virtual void Solution(TPZVec< STATE > &Sol, TPZFMatrix< STATE > &DSol, TPZFMatrix< REAL > &axes, int var, TPZVec< STATE > &Solout) override
Definition: pzplaca.cpp:326
int64_t Rows() const
Returns number of rows.
Definition: pzmatrix.h:803
virtual int NSolutionVariables(int var) override
Definition: pzplaca.cpp:308
TPZFNMatrix< 9, REAL > axes
axes indicating the directions of the derivatives of the shapefunctions
TPZVec< STATE > fXF
Definition: pzplaca.h:28
STATE fG13
Definition: pzplaca.h:24
TPZFMatrix< STATE > & Val1()
Definition: pzbndcond.h:253
virtual int NStateVariables() const override
Returns the number of state variables associated with the material.
Definition: pzplaca.h:35
TPZFMatrix< STATE > fRmatT
Definition: pzplaca.h:25
virtual void Flux(TPZVec< REAL > &x, TPZVec< STATE > &u, TPZFMatrix< STATE > &dudx, TPZFMatrix< REAL > &axes, TPZVec< STATE > &fl) override
Computes the value of the flux function to be used by ZZ error estimator.
Definition: pzplaca.cpp:273
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
int ClassId() const override
Unique identifier for serialization purposes.
TPZFMatrix< STATE > fnaxes
Definition: pzplaca.h:23
static int idf[6]
virtual int NSolutionVariables(int var)
Returns the number of variables associated with the variable indexed by var.
void Print(std::ostream &out) override
Prints out the data associated with the material.
Definition: pzplaca.cpp:277
virtual int ClassId() const override
Define the class id associated with the class.
Definition: pzplaca.cpp:482
TPZFMatrix< STATE > fKyxR
Definition: pzplaca.h:26
TPZAutoPointer< TPZFunction< STATE > > fForcingFunction
Pointer to forcing function, it is the right member at differential equation.
Definition: TPZMaterial.h:47
int64_t Cols() const
Returns number of cols.
Definition: pzmatrix.h:809
void(* fExactFunction)(TPZFMatrix< REAL > &axes, TPZVec< REAL > &x, TPZFMatrix< STATE > &uexact, TPZFMatrix< STATE > &duexact)
Definition: pzplaca.h:99
TPZFMatrix< STATE > fB00R
Definition: pzplaca.h:26
STATE fh
Definition: pzplaca.h:24
def values
Definition: rdt.py:119
void Transpose(TPZMatrix< TVar > *const T) const override
It makes *T the transpose of current matrix.
Definition: pzfmatrix.cpp:1073
#define PZError
Defines the output device to error messages and the DebugStop() function.
Definition: pzerror.h:15
STATE fE2
Definition: pzplaca.h:24
TPZMaterial * Material() const
Definition: pzbndcond.h:263
TPZFMatrix< STATE > fRmat
Definition: pzplaca.h:25