NeoPZ
TPZLinearConvecDiff.h
Go to the documentation of this file.
1 //---------------------------------------------------------------------------
2 
3 #ifndef TPZLinearConvecDiffH
4 #define TPZLinearConvecDiffH
5 #include <iostream>
6 #include "pzfmatrix.h"
7 #include "TPZMaterial.h"
8 
9 
18 
19  protected :
20 
22  STATE fXf;
23 
25  STATE fK;
26 
28  REAL fConvDir[2];
29 
31  STATE fSD;
32 
33 public:
34 
35  TPZLinearConvecDiff(int nummat, REAL k, const TPZVec<REAL> &conv, REAL f, REAL SD);
36 
37  TPZLinearConvecDiff(int matid);
38 
40 
42 
43  virtual ~TPZLinearConvecDiff();
44 
45  virtual TPZMaterial * NewMaterial() override {
46  return new TPZLinearConvecDiff(*this);
47  }
48 
49  virtual int Dimension() const override { return 2;}
50 
51  virtual int NStateVariables() const override { return 1; }
52 
53  virtual void Print(std::ostream & out) override;
54 
55  virtual std::string Name() override { return "TPZLinearConvecDiff"; }
56 
57  virtual void Contribute(TPZMaterialData &data,REAL weight,TPZFMatrix<STATE> &ek,TPZFMatrix<STATE> &ef) override;
58 
59  virtual void ContributeBC(TPZMaterialData &data,REAL weight,
61 
62  virtual int VariableIndex(const std::string &name) override;
63 
64  virtual int NSolutionVariables(int var) override;
65 
66  virtual int NFluxes() override { return 2;}
67 
68  virtual void Solution(TPZMaterialData &data, int var, TPZVec<STATE> &Solout) override;
69 
72  TPZVec<STATE> &u_exact,TPZFMatrix<STATE> &du_exact,TPZVec<REAL> &values) override;
73 
74  virtual int NEvalErrors() override {return 3;}
75 
76  public:
77 virtual int ClassId() const override;
78 
79  virtual void Write(TPZStream &buf, int withclassid) const override {
80  DebugStop();
81  }
82 
83  virtual void Read(TPZStream &buf, void *context) override {
84  DebugStop();
85  }
86 
87 };
88 
89 
90 #endif
virtual void ContributeBC(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc) override
void
virtual int NEvalErrors() override
Returns the number of norm errors. Default is 3: energy, L2 and H1.
virtual TPZMaterial * NewMaterial() override
To create another material of the same type.
clarg::argBool bc("-bc", "binary checkpoints", false)
REAL fConvDir[2]
Convection vector.
virtual std::string Name() override
Returns the name of the material.
virtual int NStateVariables() const override
Returns the number of state variables associated with the material.
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...
virtual int NFluxes() override
Returns the number of components which form the flux function.
STATE fK
Coeficient which multiplies the Laplacian operator.
virtual void Solution(TPZMaterialData &data, int var, TPZVec< STATE > &Solout) override
Returns the solution associated with the var index based on the finite element approximation.
This abstract class defines the behaviour which each derived class needs to implement.
Definition: TPZMaterial.h:39
f
Definition: test.py:287
virtual int ClassId() const override
Define the class id associated with the class.
virtual void Write(TPZStream &buf, int withclassid) const override
Writes this object to the TPZStream buffer. Include the classid if withclassid = true.
Contains TPZMatrixclass which implements full matrix (using column major representation).
#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
virtual int VariableIndex(const std::string &name) override
Returns the variable index associated with the name.
virtual void Print(std::ostream &out) override
Prints out the data associated with the material.
virtual void Read(TPZStream &buf, void *context) override
read objects from the stream
virtual int Dimension() const override
Returns the integrable dimension of the material.
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...
Defines the interface for saving and reading data. Persistency.
Definition: TPZStream.h:50
def values
Definition: rdt.py:119
Convecção-difusão linear 2D.
STATE fXf
Forcing function value.
STATE fSD
Multiplication value for the streamline diffusion term.
virtual int NSolutionVariables(int var) override
Returns the number of variables associated with the variable indexed by var.