NeoPZ
pzmat2dlin.cpp
Go to the documentation of this file.
1 
6 #include "pzmat2dlin.h"
7 #include "TPZMaterial.h"
8 #include "pzconnect.h"
9 #include "pzbndcond.h"
10 #include <math.h>
11 #include "pzvec.h"
12 #include "pzerror.h"
13 #include "TPZStream.h"
14 
15 #include "pzpoisson3d.h"
16 
17 using namespace std;
18 
19 void TPZMat2dLin::Contribute( TPZMaterialData &data,REAL weight,
21 {
22  // this method adds the contribution of the material to the stiffness
23  // matrix and right hand side
24 
25  // check on the validity of the arguments
26  //rows x cols
27  TPZFMatrix<REAL> &dphi = data.dphix;
28  TPZFMatrix<REAL> &phi = data.phi;
29  TPZManVector<REAL,3> &x = data.x;
30  TPZFMatrix<REAL> &axes=data.axes;
31 
32  if(phi.Cols() != 1 || dphi.Rows() != 2 || phi.Rows() != dphi.Cols())
33  {
34  PZError << "TPZMat2dLin.contr, inconsistent input data : phi.Cols() = "
35  << phi.Cols() << " dphi.Cols + " << dphi.Cols()
36  << " phi.Rows = " << phi.Rows() << " dphi.Rows = " << dphi.Rows() << "\n";
37  DebugStop();
38  }
39  if(fForcingFunction)
40  {
41  TPZManVector<STATE> xfloat(fXf.Rows());
42  fForcingFunction->Execute(x,xfloat);//fXf = xfloat
43  int i;
44  for(i=0; i<fXf.Rows(); i++) fXf(i,0) = xfloat[i];
45  }
46  int r = fKxx.Rows();
47  int ic,jc;
48  for(int in=0 ; in < phi.Rows() ; ++in)
49  {
50  for(ic=0; ic< r; ic++)
51  {
52  ef(in*r+ic, 0) += (STATE)(weight*phi(in,0))*fXf(ic,0);
53  }
54  for(int jn=0 ; jn<phi.Rows() ; ++jn)
55  {
56  for(ic=0; ic<r; ic++) for(jc=0; jc<r; jc++)
57  {
58  REAL dphix = dphi(0,in)*axes(0,0)+axes(1,0)*dphi(1,in);
59  REAL dphiy = dphi(0,in)*axes(0,1)+axes(1,1)*dphi(1,in);
60  REAL dphjx = dphi(0,jn)*axes(0,0)+axes(1,0)*dphi(1,jn);
61  REAL dphjy = dphi(0,jn)*axes(0,1)+axes(1,1)*dphi(1,jn);
62  ek(in*r+ic,jn*r+jc) += (STATE)(weight*phi(in,0)*phi(jn,0))*fK00(ic,jc) +
63  ((STATE)(dphix*dphjx)*fKxx(ic,jc) +
64  (STATE)(dphiy*dphjy)*fKyy(ic,jc) +
65  (STATE)(dphix*phi(jn,0))*fKx0(ic,jc) +
66  (STATE)(dphiy*phi(jn,0))*fKy0(ic,jc) +
67  (STATE)(phi(in,0)*dphjx)*fK0x(ic,jc) +
68  (STATE)(phi(in,0)*dphjy)*fK0y(ic,jc) ) * (STATE)weight;
69  }
70  }
71  }
72 }
73 
76  TPZFMatrix<REAL> &phi = data.phi;
77  TPZFMatrix<REAL> &axes=data.axes;
78 
79  if(bc.Material() != this){
80  PZError << "TPZMat1dLin.apply_bc warning : this material didn't create the boundary condition!\n";
81  }
82 
83  if(bc.Type() < 0 && bc.Type() > 2){
84  PZError << "TPZMat1dLin.aplybc, unknown boundary condition type :" <<
85  bc.Type() << " boundary condition ignored\n";
86  }
87 
88  int r = fKxx.Rows();
89  int numnod = ek.Rows()/r;
90 
91  int idf,jdf,in,jn;
92  switch(bc.Type()){
93 
94  case 0:
95  for(in=0 ; in<numnod ; ++in){
96  for(idf = 0;idf<r;idf++) {
97  (ef)(in*r+idf,0) += (STATE)(weight*gBigNumber*phi(in,0))*bc.Val2()(idf,0);
98  }
99  for(jn=0 ; jn<numnod ; ++jn) {
100  for(idf = 0;idf<r;idf++) {
101  ek(in*r+idf,jn*r+idf) += gBigNumber*phi(in,0)*phi(jn,0)*weight;
102  }
103  }
104  }
105  break;
106 
107  case 1:
108  for(in=0 ; in<numnod ; ++in){
109  for(idf = 0;idf<r;idf++) {
110  (ef)(in*r+idf,0) += (STATE)(weight*phi(in,0))*bc.Val2()(idf,0);
111  }
112  }
113  break;
114 
115  case 2:
116  for(in=0 ; in<numnod ; ++in){
117  for(idf = 0;idf<r;idf++) {
118  (ef)(in*r+idf,0) += (STATE)(weight*phi(in,0))*bc.Val2()(idf,0);
119  }
120  for(jn=0 ; jn<numnod ; ++jn) {
121  for(idf = 0;idf<r;idf++) {
122  for(jdf = 0;jdf<r;jdf++) {
123  ek(in*r+idf,jn*r+jdf) += bc.Val1()(idf,jdf)*(STATE)(phi(in,0)*phi(jn,0)*weight);
124  }
125  }
126  }
127  }
128  break;
129  case 3:
130  for(in=0; in<numnod; in++) {
131  for(idf=0; idf<r; idf++) {
132  for(jdf=0; jdf<r; jdf++) {
133  ef(in*r+idf,0) += (fKx0(idf,jdf)*(STATE)axes(0,1)-fKy0(idf,jdf)*(STATE)axes(0,0))*(STATE)(phi(in,0)*weight);
134  }
135  }
136  }
137  break;
138  }//fim switch
139 }
140 
141 int TPZMat2dLin::NFluxes() {return 1;}
142 
144  PZError << "TPZMat2dLin::Flux is called\n";
145 }
146 
147 void TPZMat2dLin::Print(std::ostream & out) {
148  out << "Material type TPZMat2dLin -- number = " << Id() << "\n";
149  out << "Matrix Kxx -> "; fKxx.Print("fKxx",out);
150  out << "Matrix Kyy -> "; fKyy.Print("fKyy",out);
151  out << "Matrix Kxy -> "; fKxy.Print("fKxy",out);
152  out << "Matrix Kyx -> "; fKyx.Print("fKyx",out);
153  out << "Matrix Kx0 -> "; fKx0.Print("fKx0",out);
154  out << "Matrix K0x -> "; fK0x.Print("fK0x",out);
155  out << "Matrix Ky0 -> "; fKy0.Print("fKy0",out);
156  out << "Matrix K0y -> "; fK0y.Print("fK0y",out);
157  out << "Matrix K00 -> "; fK00.Print("fK00",out);
158  out << "Matrix xf -> "; fXf.Print("fXf",out);
159 }
160 
161 int TPZMat2dLin::VariableIndex(const std::string &name) {
162  if(!strcmp(name.c_str(),"displacement")) return 1;
163  if(!strcmp(name.c_str(),"Solution")) return 1;
164  if(!strcmp("Derivate",name.c_str())) return 2;
165  if(!strcmp("Pressure",name.c_str())) return 3;
166  return TPZMaterial::VariableIndex(name);
167 }
168 
170  if(index == 1) return 3;
171  if(index == 2) return 2;
172  if(index == 3) return 1;
173  return TPZMaterial::NSolutionVariables(index);
174 }
175 
176 
178 
179  if(var == 0) {
180 #ifndef STATE_COMPLEX
181  Solout.Resize(Sol.size());
182  for (int i = 0; i<Sol.size(); i++) {
183  Solout[i] = Sol[i];
184  }
185 #else
186  for (int i=0; i<Sol.size(); i++) {
187  Solout[i] = Sol[i].real();
188  }
189 #endif
190  } else if(var == 1) {
191  Solout.Resize(3,0.);
192 #ifndef STATE_COMPLEX
193  Solout[0] = Sol[0];
194 #else
195  Solout[0] = Sol[0].real();
196 #endif
197  } else if(var == 2) {
198 #ifndef STATE_COMPLEX
199  Solout[0] = DSol(0,0);//derivate
200  Solout[1] = DSol(1,0);//derivate
201 #else
202  Solout[0] = DSol(0,0).real();//derivate
203  Solout[1] = DSol(1,0).real();//derivate
204 #endif
205  return;
206  } else if(var == 3) {
207  Solout.Resize(1,0.);
208 #ifndef STATE_COMPLEX
209  Solout[0] = Sol[0];
210 #else
211  Solout[0] = Sol[0].real();
212 #endif
213  }
214 
215  else TPZMaterial::Solution(Sol,DSol,axes,var,Solout);
216 }
217 
219  TPZVec<STATE> &u_exact,TPZFMatrix<STATE> &du_exact,TPZVec<REAL> &values) {
220 
221  STATE dx = dudx(0,0)*(STATE)axes(0,0)+dudx(1,0)*(STATE)axes(1,0);
222  STATE dy = dudx(0,0)*(STATE)axes(0,1)+dudx(1,0)*(STATE)axes(1,1);
223  //Norma Energia
224  STATE parc1 = dx-du_exact(0,0) ;
225  STATE parc2 = dy-du_exact(1,0) ;
226  values[0] = abs(parc1*parc1 + parc2*parc2);/*pow(parc1,2.)+pow(parc2,2.);*/
227  //Norma L2
228  values[1] = pow((REAL)fabs(u[0] - u_exact[0]),(REAL)2.0);
229  values[2] = 0.;
230 }
231 
233  return new TPZMat2dLin(*this);
234 }
235 
237  fKxx.Redim(1,1);
238  fKxy.Redim(1,1);
239  fKyx.Redim(1,1);
240  fKyy.Redim(1,1);
241  fKx0.Redim(1,1);
242  fK0x.Redim(1,1);
243  fKy0.Redim(1,1);
244  fK0y.Redim(1,1);
245  fK00.Redim(1,1);
246  fXf.Redim(1,1);
247  REAL cosa = cos(angle);
248  REAL sina = sin(angle);
249  fKxx(0,0) = cosa*cosa;
250  fKxy(0,0) = cosa*sina;
251  fKyy(0,0) = sina*sina;
252  fKyx(0,0) = cosa*sina;
253  fKx0(0,0) = -cosa;
254  fKy0(0,0) = -sina;
255 }
256 
258 
259  int nstate = fKxx.Rows();
260  TPZFMatrix<STATE> val1(nstate,nstate),val2(nstate,1);
261  return TPZMaterial::CreateBC(reference,bc,3,val1,val2);
262 }
263 
268  return Hash("TPZMat2dLin") ^ TPZMaterial::ClassId() << 1;
269 }
270 
272 void TPZMat2dLin::Write(TPZStream &buf, int withclassid) const
273 {
274 #ifdef PZDEBUG2
275  if (logger->isDebugEnabled())
276  {
277  std::stringstream sout;
278  sout << __PRETTY_FUNCTION__ << " before write material ";
279  LOGPZ_DEBUG( logger,sout.str().c_str() );
280  }
281 #endif
282  TPZMaterial::Write(buf,withclassid);
283 #ifdef PZDEBUG2
284  if (logger->isDebugEnabled())
285  {
286  std::stringstream sout;
287  sout << __PRETTY_FUNCTION__ << " after write material ";
288  LOGPZ_DEBUG( logger,sout.str().c_str() );
289  }
290 #endif
291 
292  fKxx.Write(buf,0);
293  fKxy.Write(buf,0);
294  fKyx.Write(buf,0);
295  fKyy.Write(buf,0);
296  fKx0.Write(buf,0);
297  fK0x.Write(buf,0);
298  fKy0.Write(buf,0);
299  fK0y.Write(buf,0);
300  fK00.Write(buf,0);
301  fXf.Write(buf,0);
302 }
303 
305 void TPZMat2dLin::Read(TPZStream &buf, void *context)
306 {
307  TPZMaterial::Read(buf,context);
308  fKxx.Read(buf,0);
309  fKxy.Read(buf,0);
310  fKyx.Read(buf,0);
311  fKyy.Read(buf,0);
312  fKx0.Read(buf,0);
313  fK0x.Read(buf,0);
314  fKy0.Read(buf,0);
315  fK0y.Read(buf,0);
316  fK00.Read(buf,0);
317  fXf.Read(buf,0);
318 }
319 
320 #ifndef BORLAND
321 template class TPZRestoreClass<TPZMat2dLin>;
322 #endif
void Solution(TPZVec< STATE > &Sol, TPZFMatrix< STATE > &DSol, TPZFMatrix< REAL > &axes, int var, TPZVec< STATE > &Solout) override
Definition: pzmat2dlin.cpp:177
virtual void Solution(TPZMaterialData &data, int var, TPZVec< STATE > &Solout)
Returns the solution associated with the var index based on the finite element approximation.
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
TPZManVector< REAL, 3 > x
value of the coordinate at the integration point
virtual int VariableIndex(const std::string &name) override
Returns the variable index associated with the name.
Definition: pzmat2dlin.cpp:161
clarg::argBool bc("-bc", "binary checkpoints", false)
int Type()
Definition: pzbndcond.h:249
Templated vector implementation.
virtual int VariableIndex(const std::string &name)
Returns the variable index associated with the name.
Defines PZError.
void Write(TPZStream &buf, int withclassid) const override
Saves the element data to a stream.
TinyFad< 8, T > abs(const TinyFad< 8, T > &in)
Definition: tinyfadeight.h:846
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: pzmat2dlin.cpp:19
Contains the TPZMatPoisson3d class.
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 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: pzmat2dlin.cpp:218
TPZFNMatrix< 660, REAL > dphix
values of the derivative of the shape functions
virtual int ClassId() const override
returns the unique identifier for reading/writing objects to streams
Definition: pzmat2dlin.cpp:267
TPZBndCond * OutflowFlux(TPZMaterial *&reference, int bc)
Definition: pzmat2dlin.cpp:257
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...
virtual void Write(TPZStream &buf, int withclassid) const override
Saves the element data to a stream.
Definition: pzmat2dlin.cpp:272
int64_t size() const
Returns the number of elements of the vector.
Definition: pzvec.h:196
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
sin
Definition: fadfunc.h:63
void Read(TPZStream &buf, void *context) override
Reads the element data from a stream.
expr_ dx(i) *cos(expr_.val())
virtual int NSolutionVariables(int index) override
Returns the number of variables associated with the variable indexed by var.
Definition: pzmat2dlin.cpp:169
virtual TPZMaterial * NewMaterial() override
Creates a copy of the material object.
Definition: pzmat2dlin.cpp:232
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
virtual void Read(TPZStream &buf, void *context) override
Reads the element data from a stream.
Definition: pzmat2dlin.cpp:305
#define LOGPZ_DEBUG(A, B)
Define log for debug info.
Definition: pzlog.h:87
This class defines the boundary condition for TPZMaterial objects.
Definition: pzbndcond.h:29
Contains declaration of TPZConnect class which represents a set of shape functions associated with a ...
int64_t Rows() const
Returns number of rows.
Definition: pzmatrix.h:803
TPZFNMatrix< 9, REAL > axes
axes indicating the directions of the derivatives of the shapefunctions
TPZFMatrix< STATE > & Val1()
Definition: pzbndcond.h:253
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: pzmat2dlin.cpp:74
void ConvectionDiffusion(REAL angle, REAL diff)
Definition: pzmat2dlin.cpp:236
Contains the TPZMat2dLin class which implements a bi-dimensional linear problem.
void Print(std::ostream &out=std::cout) override
Prints out the data associated with the material.
Definition: pzmat2dlin.cpp:147
int32_t Hash(std::string str)
Definition: TPZHash.cpp:10
static REAL angle
Angle in radians to test.
Definition: pzsubcmesh.cpp:53
int ClassId() const override
Unique identifier for serialization purposes.
virtual TPZBndCond * CreateBC(TPZMaterial *reference, int id, int typ, TPZFMatrix< STATE > &val1, TPZFMatrix< STATE > &val2)
Creates an object TPZBndCond derived of TPZMaterial.
TPZFlopCounter pow(const TPZFlopCounter &orig, const TPZFlopCounter &xp)
Returns the power and increments the counter of the power.
Definition: pzreal.h:487
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: pzmat2dlin.cpp:143
static int idf[6]
virtual int NSolutionVariables(int var)
Returns the number of variables associated with the variable indexed by var.
Implements a bi-dimensional linear problem.
Definition: pzmat2dlin.h:22
int64_t Cols() const
Returns number of cols.
Definition: pzmatrix.h:809
Contains declaration of the abstract TPZStream class. TPZStream defines the interface for saving and ...
Defines the interface for saving and reading data. Persistency.
Definition: TPZStream.h:50
def values
Definition: rdt.py:119
virtual int NFluxes() override
Returns the number of components which form the flux function.
Definition: pzmat2dlin.cpp:141
TPZFlopCounter cos(const TPZFlopCounter &orig)
Returns the cosine in radians and increments the counter of the Cosine.
Definition: pzreal.h:514
Implements an interface to register a class id and a restore function. Persistence.
Definition: TPZSavable.h:150
#define PZError
Defines the output device to error messages and the DebugStop() function.
Definition: pzerror.h:15
TPZMaterial * Material() const
Definition: pzbndcond.h:263