NeoPZ
TPZMaterial.h
Go to the documentation of this file.
1 
7 #ifndef PZMATERIALHPP
8 #define PZMATERIALHPP
9 
10 #include "pzreal.h"
11 #include "pzvec.h"
12 #include "pzmatrix.h"
13 #include "pzfmatrix.h"
14 #include "pzadmchunk.h"
15 #include "tpzautopointer.h"
16 #include "TPZSavable.h"
17 #include "pzmaterialdata.h"
18 #include "pzfunction.h"
19 #include "pzcompel.h"
20 
21 #include <iostream>
22 #include <string>
23 
24 class TPZBndCond;
25 class TPZMaterial;
26 class TPZMaterialData;
27 class TPZIntPoints;
28 
39 class TPZMaterial : public virtual TPZSavable
40 {
41 private:
42  int fId;
43 
44 protected:
45 
48 
51 
54 
57 
60 
63 
64 
77 
80 
81 public:
83  static REAL gBigNumber;
84 
87  TPZMaterial(int id);
88 
90  TPZMaterial();
91 
94  TPZMaterial(const TPZMaterial &mat);
95 
97  TPZMaterial &operator=(const TPZMaterial &copy);
98 
100  virtual ~TPZMaterial();
101 
110  virtual void FillDataRequirements(TPZMaterialData &data);
111 
117  virtual void FillDataRequirements(TPZVec<TPZMaterialData > &datavec);
118 
121  {
122  // default is no specific data requirements
123  if(type == 50)
124  {
125  data.fNeedsSol = true;
126  }
127  }
128 
131  {
132  // default is no specific data requirements
133  int nref = datavec.size();
134  if(type == 50)
135  {
136  for(int iref = 0; iref<nref; iref++){
137  datavec[iref].fNeedsSol = true;
138  }
139  }
140  }
141 
144  {
145  data.fNeedsNormal = false;
146  }
147 
150  {
151  data.SetAllRequirements(false);
152  int nref_left = datavec_left.size();
153  for(int iref = 0; iref<nref_left; iref++){
154  datavec_left[iref].SetAllRequirements(false);
155  }
156  int nref_right = datavec_right.size();
157  for(int iref = 0; iref<nref_right; iref++){
158  datavec_right[iref].SetAllRequirements(false);
159  }
160 
161  }
162 
163 
165  virtual std::string Name() { return "no_name"; }
166 
168  virtual int Dimension() const = 0;
169 
170  int Id() const { return fId; }
171  void SetId(int id) {
172 /* if(id == 0) {
173  std::cout << "\n*** Material Id can't be ZERO! ***\n";
174  std::cout << "*** This Will Be a Disaster!!! ***\n";
175  DebugStop();
176  }*/
177  fId = id; }
178 
180  virtual int NStateVariables() const = 0;
181 
183  virtual int NFluxes() {return 0;}
184 
187  {
188  return fNumLoadCases;
189  }
190 
193  {
194  return 1;
195  }
196 
198  void SetNumLoadCases(int numloadcases)
199  {
200  if(numloadcases <= 0)
201  {
202  std::cout << __PRETTY_FUNCTION__ << " numloadcases " << numloadcases << " cannot be less or equal to zero\n";
203  DebugStop();
204  }
205  fNumLoadCases = numloadcases;
206  }
207 
209  void SetPostProcessIndex(int index)
210  {
211 #ifdef PZDEBUG
212  if (index < 0 || index >= fNumLoadCases)
213  {
214  DebugStop();
215  }
216 #endif
217  fPostProcIndex = index;
218  }
219 
221  virtual void Print(std::ostream &out = std::cout);
222 
228  virtual int VariableIndex(const std::string &name);
229 
234  virtual int NSolutionVariables(int var);
235 
237  virtual void Solution(TPZMaterialData &data, int var, TPZVec<STATE> &Solout);
238 
240  virtual void Solution(TPZVec<TPZMaterialData> &datavec, int var, TPZVec<STATE> &Solout);
241 
243  virtual void Solution(TPZMaterialData &data, TPZVec<TPZMaterialData> &dataleftvec, TPZVec<TPZMaterialData> &datarightvec, int var, TPZVec<STATE> &Solout);
244 
246  virtual void Solution(TPZMaterialData &data, TPZVec<TPZMaterialData> &dataleftvec, TPZVec<TPZMaterialData> &datarightvec, int var, TPZVec<STATE> &Solout, TPZCompEl * left, TPZCompEl * ritgh);
247 
248 protected:
250  virtual void Solution(TPZVec<STATE> &Sol,TPZFMatrix<STATE> &DSol,TPZFMatrix<REAL> &axes,int var,TPZVec<STATE> &Solout);
251 
252 public:
253 
255  virtual void Flux(TPZVec<REAL> &x, TPZVec<STATE> &Sol,
256  TPZFMatrix<STATE> &DSol, TPZFMatrix<REAL> &axes,
257  TPZVec<STATE> &flux) {}
258 
262  virtual TPZBndCond *CreateBC(TPZMaterial *reference, int id, int typ, TPZFMatrix<STATE> &val1,
263  TPZFMatrix<STATE> &val2);
264 
277  virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef) = 0;
278 
286  virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix<STATE> &ef);
287 
288 
296  virtual void Contribute(TPZVec<TPZMaterialData> &datavec, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef);
297 
304  virtual void Contribute(TPZVec<TPZMaterialData> &datavec, REAL weight, TPZFMatrix<STATE> &ef)
305  {
306  TPZFMatrix<STATE> ek(ef.Rows(),ef.Rows(),0.);
307  Contribute(datavec, weight, ek, ef);
308  }
309 
319  virtual void ContributeBC(TPZMaterialData &data, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef, TPZBndCond &bc) = 0;
320 
331  virtual void ContributeBC(TPZVec<TPZMaterialData> &datavec, REAL weight, TPZFMatrix<STATE> &ef, TPZBndCond &bc)
332  {
333  TPZFMatrix<STATE> ek(ef.Rows(),ef.Rows(),0.);
334  ContributeBC(datavec, weight, ek, ef, bc);
335 
336  }
337 
348  virtual void ContributeBC(TPZVec<TPZMaterialData> &datavec, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef, TPZBndCond &bc);
349 
358  virtual void ContributeBC(TPZMaterialData &data, REAL weight, TPZFMatrix<STATE> &ef, TPZBndCond &bc);
359 
369  {
370  fForcingFunction = fp;
371  }
372  void SetForcingFunction(void (*fp)(const TPZVec<REAL> &loc, TPZVec<STATE> &result), int porder )
373  {
374  if(fp)
375  {
377  loc->SetPolynomialOrder(porder);
378  fForcingFunction = loc;
379  }
380 
381  else fForcingFunction = NULL;
382  }
383 
384  void SetForcingFunction(void (*fp)(const TPZVec<REAL> &loc, TPZVec<STATE> &result, TPZFMatrix<STATE> &gradu), int porder )
385  {
386  if(fp)
387  {
389  loc->SetPolynomialOrder(porder);
390  fForcingFunction = loc;
391  }
392 
393  else fForcingFunction = NULL;
394  }
395 
398  return fForcingFunction;
399  }
400 
406  {
407  fForcingFunctionExact = fp;
408  }
409 
412  return fForcingFunctionExact;
413  }
414 
420  {
421  fTimeDependentForcingFunction = fp;
422  }
423 
427  }
428 
434  {
435  fTimedependentFunctionExact = fp;
436  }
437 
441  }
442 
448  {
449  fBCForcingFunction = fp;
450  }
451 
454  return fBCForcingFunction;
455  }
456 
462  {
463  fTimedependentBCForcingFunction = fp;
464  }
465 
469  }
470 
472  virtual int HasForcingFunction() {return (fForcingFunction != 0);}
473 
475  virtual int HasForcingFunctionExact() {return (fForcingFunctionExact != 0);}
476 
478  virtual int HasBCForcingFunction() {return (fBCForcingFunction != 0);}
479 
481  virtual int HasTimedependentFunctionExact() {return (fTimedependentFunctionExact != 0);}
482 
484  virtual int HasTimedependentForcingFunction() {return (fTimeDependentForcingFunction != 0);}
485 
487  virtual int HasTimedependentBCForcingFunction() {return (fTimedependentBCForcingFunction != 0);}
488 
489 
491  virtual int IntegrationRuleOrder(int elPMaxOrder) const;
492 
494  virtual int IntegrationRuleOrder(TPZVec<int> &elPMaxOrder) const;
495 
496  virtual void Errors(TPZMaterialData &data, TPZVec<STATE> &u_exact, TPZFMatrix<STATE> &du_exact, TPZVec<REAL> &errors)
497  {
499  Flux(data.x, data.sol[0], data.dsol[0], data.axes, flux);
500  Errors(data.x, data.sol[0], data.dsol[0], data.axes, flux, u_exact, du_exact, errors );
501  }
502  virtual void Errors(TPZVec<TPZMaterialData> &data, TPZVec<STATE> &u_exact, TPZFMatrix<STATE> &du_exact, TPZVec<REAL> &errors)
503  {
504  DebugStop();
505  }
506 
511  virtual void Errors(TPZVec<REAL> &x, TPZVec<STATE> &sol, TPZFMatrix<STATE> &dsol,
512  TPZFMatrix<REAL> &axes, TPZVec<STATE> &flux,
513  TPZVec<STATE> &uexact, TPZFMatrix<STATE> &duexact,
514  TPZVec<REAL> &val) {
515  PZError << __PRETTY_FUNCTION__ << std::endl;
516  PZError << "Method not implemented! Error comparison not available. Please, implement it." << std::endl;
517  }
519  PZError << __PRETTY_FUNCTION__ << std::endl;
520  PZError << "Nao sei o q fazer." << std::endl;
521 
522  }
524  virtual int NEvalErrors() {return 3;}
525 
527  virtual TPZMaterial * NewMaterial();
528 
530  virtual void SetData(std::istream &data);
531 
533  virtual void Clone(std::map<int, TPZMaterial * > &matvec);
534 
536  virtual int FluxType() { return 2; }
537 
538  virtual void ContributeErrors(TPZMaterialData &data,
539  REAL weight,
540  TPZVec<REAL> &nk,
541  int &errorid){
542  PZError << "Error at " << __PRETTY_FUNCTION__ << " - Method not implemented\n";
543  }
544 
552  PZError << "Error at " << __PRETTY_FUNCTION__ << " - Method not implemented\n";
553  return -1.;
554  }
555 
561  virtual int PushMemItem(int sourceIndex = -1){ return -1; }
562 
564  virtual void FreeMemItem(int index){ return; }
565 
567  void SetLinearContext(bool IsLinear);
568 
570  bool GetLinearContext() const {
571  return fLinearContext;
572  }
573 
579  public:
580 int ClassId() const override;
581 
582 
584  void Write(TPZStream &buf, int withclassid) const override;
585 
587  void Read(TPZStream &buf, void *context) override;
588 
591 };
592 
595 
596 #endif
597 
bool fLinearContext
Defines whether the equation context is linear solver or non linear.
Definition: TPZMaterial.h:70
virtual void Solution(TPZMaterialData &data, int var, TPZVec< STATE > &Solout)
Returns the solution associated with the var index based on the finite element approximation.
Contains declaration of the TPZSavable class which defines the interface to save and restore objects ...
int fNumLoadCases
Defines the number of load cases generated by this material.
Definition: TPZMaterial.h:76
virtual void Errors(TPZMaterialData &data, TPZVec< STATE > &u_exact, TPZFMatrix< STATE > &du_exact, TPZVec< REAL > &errors)
Definition: TPZMaterial.h:496
TPZManVector< REAL, 3 > x
value of the coordinate at the integration point
TPZAutoPointer< TPZFunction< STATE > > & TimedependentBCForcingFunction()
Returns a procedure as time variable boundary condition.
Definition: TPZMaterial.h:467
clarg::argBool bc("-bc", "binary checkpoints", false)
virtual int HasBCForcingFunction()
Directive that gives true if the material has a bc forcing function exact.
Definition: TPZMaterial.h:478
TPZAutoPointer< TPZFunction< STATE > > fBCForcingFunction
Pointer to bc forcing function, it is a variable boundary condition at differential equation...
Definition: TPZMaterial.h:59
virtual REAL ComputeSquareResidual(TPZVec< REAL > &X, TPZVec< STATE > &sol, TPZFMatrix< STATE > &dsol)
Computes square of residual of the differential equation at one integration point.
Definition: TPZMaterial.h:551
virtual void ErrorsHdiv(TPZMaterialData &data, TPZVec< STATE > &u_exact, TPZFMatrix< STATE > &du_exact, TPZVec< REAL > &values)
Definition: TPZMaterial.h:518
virtual TPZMaterial * NewMaterial()
To create another material of the same type.
Contains declaration of TPZCompEl class which defines the interface of a computational element...
void SetForcingFunction(void(*fp)(const TPZVec< REAL > &loc, TPZVec< STATE > &result, TPZFMatrix< STATE > &gradu), int porder)
Definition: TPZMaterial.h:384
Templated vector implementation.
TPZMaterial & operator=(const TPZMaterial &copy)
operator =
Definition: TPZMaterial.cpp:62
void SetBCForcingFunction(TPZAutoPointer< TPZFunction< STATE > > fp)
Sets a procedure as variable boundary condition.
Definition: TPZMaterial.h:447
virtual int VariableIndex(const std::string &name)
Returns the variable index associated with the name.
void SetAllRequirements(bool set)
Set all flags at once.
bool GetLinearContext() const
Returns fLinearContext attribute.
Definition: TPZMaterial.h:570
virtual int Dimension() const =0
Returns the integrable dimension of the material.
virtual int NStateVariables() const =0
Returns the number of state variables associated with the material.
void Write(TPZStream &buf, int withclassid) const override
Saves the element data to a stream.
TPZGradSolVec dsol
vector of the derivatives of the solution at the integration point
Declarates the TPZBlock<REAL>class which implements block matrices.
virtual void Errors(TPZVec< REAL > &x, TPZVec< STATE > &sol, TPZFMatrix< STATE > &dsol, TPZFMatrix< REAL > &axes, TPZVec< STATE > &flux, TPZVec< STATE > &uexact, TPZFMatrix< STATE > &duexact, TPZVec< REAL > &val)
Computes the error due to the difference between the interpolated flux and the flux computed based o...
Definition: TPZMaterial.h:511
virtual int MinimumNumberofLoadCases()
returns the minimum number of load cases for this material
Definition: TPZMaterial.h:192
This class implements a simple vector storage scheme for a templated class T. Utility.
Definition: pzgeopoint.h:19
REAL val(STATE &number)
Returns value of the variable.
Definition: pzartdiff.h:23
TPZMaterial()
Default constructor.
Definition: TPZMaterial.cpp:31
virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef)=0
It computes a contribution to the stiffness matrix and load vector at one integration point...
Abstract class defining integration rules. Numerical Integration.
Definition: tpzintpoints.h:19
virtual void Print(std::ostream &out=std::cout)
Prints out the data associated with the material.
int fPostProcIndex
indicates which solution should be used for post processing
Definition: TPZMaterial.h:79
This abstract class defines the behaviour which each derived class needs to implement.
Definition: TPZMaterial.h:39
virtual void ContributeBC(TPZVec< TPZMaterialData > &datavec, REAL weight, TPZFMatrix< STATE > &ef, TPZBndCond &bc)
It computes a contribution to the stiffness matrix and load vector at one BC integration point to mul...
Definition: TPZMaterial.h:331
virtual int FluxType()
To return a numerical flux type to apply over the interfaces of the elements.
Definition: TPZMaterial.h:536
TPZAutoPointer< TPZFunction< STATE > > fTimedependentBCForcingFunction
Pointer to time dependent bc forcing function, it is a variable boundary condition at differential eq...
Definition: TPZMaterial.h:62
TPZAutoPointer< TPZFunction< STATE > > & ForcingFunction()
Returns a procedure as source function for the material.
Definition: TPZMaterial.h:397
int64_t size() const
Returns the number of elements of the vector.
Definition: pzvec.h:196
void SetForcingFunctionExact(TPZAutoPointer< TPZFunction< STATE > > fp)
Sets a procedure as exact solution for the problem.
Definition: TPZMaterial.h:405
virtual void Flux(TPZVec< REAL > &x, TPZVec< STATE > &Sol, TPZFMatrix< STATE > &DSol, TPZFMatrix< REAL > &axes, TPZVec< STATE > &flux)
Computes the value of the flux function to be used by ZZ error estimator.
Definition: TPZMaterial.h:255
virtual int IntegrationRuleOrder(int elPMaxOrder) const
Gets the order of the integration rule necessary to integrate an element with polinomial order p...
virtual int NEvalErrors()
Returns the number of norm errors. Default is 3: energy, L2 and H1.
Definition: TPZMaterial.h:524
void Read(TPZStream &buf, void *context) override
Reads the element data from a stream.
virtual int NFluxes()
Returns the number of components which form the flux function.
Definition: TPZMaterial.h:183
Contains TPZMatrixclass which implements full matrix (using column major representation).
virtual void Errors(TPZVec< TPZMaterialData > &data, TPZVec< STATE > &u_exact, TPZFMatrix< STATE > &du_exact, TPZVec< REAL > &errors)
Definition: TPZMaterial.h:502
virtual void FillDataRequirements(TPZMaterialData &data)
Fill material data parameter with necessary requirements for the.
Definition: TPZMaterial.cpp:81
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
TPZAutoPointer< TPZFunction< STATE > > fForcingFunctionExact
Pointer to exact solution function, needed to calculate exact error.
Definition: TPZMaterial.h:50
virtual int PushMemItem(int sourceIndex=-1)
Pushes a new entry in the context of materials with memory, returning its index at the internal stora...
Definition: TPZMaterial.h:561
This class defines the boundary condition for TPZMaterial objects.
Definition: pzbndcond.h:29
void SetTimeDependentForcingFunction(TPZAutoPointer< TPZFunction< STATE > > fp)
Sets a procedure as time dependent source function for the material.
Definition: TPZMaterial.h:419
virtual int HasForcingFunctionExact()
Directive that gives true if the material has a function exact.
Definition: TPZMaterial.h:475
virtual void FreeMemItem(int index)
Frees an entry in the material with memory internal history storage.
Definition: TPZMaterial.h:564
virtual int HasTimedependentFunctionExact()
Directive that gives true if the material has a time dependent function exact.
Definition: TPZMaterial.h:481
virtual void ContributeErrors(TPZMaterialData &data, REAL weight, TPZVec< REAL > &nk, int &errorid)
Definition: TPZMaterial.h:538
int64_t Rows() const
Returns number of rows.
Definition: pzmatrix.h:803
virtual void FillBoundaryConditionDataRequirement(int type, TPZMaterialData &data)
This method defines which parameters need to be initialized in order to compute the contribution of t...
Definition: TPZMaterial.h:120
virtual int HasTimedependentForcingFunction()
Directive that gives true if the material has a time dependent forcing function.
Definition: TPZMaterial.h:484
TPZFNMatrix< 9, REAL > axes
axes indicating the directions of the derivatives of the shapefunctions
Contains declaration of the TPZAutoPointer class which has Increment and Decrement actions are mutexe...
void SetNumLoadCases(int numloadcases)
changes the number of load cases for this material
Definition: TPZMaterial.h:198
Contains the TPZMaterialData class which implements an interface between TPZCompEl::CalcStiff and TPZ...
virtual void FillDataRequirementsInterface(TPZMaterialData &data, TPZVec< TPZMaterialData > &datavec_left, TPZVec< TPZMaterialData > &datavec_right)
This method defines which parameters need to be initialized in order to compute the contribution of i...
Definition: TPZMaterial.h:149
virtual int HasForcingFunction()
Directive that gives true if the material has a forcing function.
Definition: TPZMaterial.h:472
virtual ~TPZMaterial()
Default destructor.
Definition: TPZMaterial.cpp:43
virtual void ContributeBC(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc)=0
It computes a contribution to the stiffness matrix and load vector at one BC integration point...
void SetTimedependentBCForcingFunction(TPZAutoPointer< TPZFunction< STATE > > fp)
Sets a procedure as time variable boundary condition.
Definition: TPZMaterial.h:461
TPZAutoPointer< TPZFunction< STATE > > & TimeDependentForcingFunction()
Returns a procedure as time dependent source function for the material.
Definition: TPZMaterial.h:425
static REAL gBigNumber
Big number to penalization method, used for Dirichlet conditions.
Definition: TPZMaterial.h:83
TPZAutoPointer< TPZFunction< STATE > > fTimedependentFunctionExact
Pointer to time dependent exact solution function, needed to calculate exact error.
Definition: TPZMaterial.h:56
Contains TPZMatrix<TVar>class, root matrix class.
TPZAutoPointer< TPZFunction< STATE > > & BCForcingFunction()
Returns a procedure as variable boundary condition.
Definition: TPZMaterial.h:453
virtual int HasTimedependentBCForcingFunction()
Directive that gives true if the material has a time dependent bc forcing function.
Definition: TPZMaterial.h:487
int ClassId() const override
Unique identifier for serialization purposes.
virtual void Contribute(TPZVec< TPZMaterialData > &datavec, REAL weight, TPZFMatrix< STATE > &ef)
It computes a contribution to the stiffness matrix and load vector at one integration point to multip...
Definition: TPZMaterial.h:304
virtual TPZBndCond * CreateBC(TPZMaterial *reference, int id, int typ, TPZFMatrix< STATE > &val1, TPZFMatrix< STATE > &val2)
Creates an object TPZBndCond derived of TPZMaterial.
virtual int NSolutionVariables(int var)
Returns the number of variables associated with the variable indexed by var.
void SetPostProcessIndex(int index)
indicates which variable should be post processed
Definition: TPZMaterial.h:209
void SetForcingFunction(TPZAutoPointer< TPZFunction< STATE > > fp)
Sets a procedure as source function for the material.
Definition: TPZMaterial.h:368
void SetPolynomialOrder(int porder)
Definition: pzfunction.h:226
void SetForcingFunction(void(*fp)(const TPZVec< REAL > &loc, TPZVec< STATE > &result), int porder)
Definition: TPZMaterial.h:372
virtual void FillBoundaryConditionDataRequirement(int type, TPZVec< TPZMaterialData > &datavec)
This method defines which parameters need to be initialized in order to compute the contribution of t...
Definition: TPZMaterial.h:130
virtual std::string Name()
Returns the name of the material.
Definition: TPZMaterial.h:165
virtual void FillDataRequirementsInterface(TPZMaterialData &data)
This method defines which parameters need to be initialized in order to compute the contribution of i...
Definition: TPZMaterial.h:143
virtual void Clone(std::map< int, TPZMaterial * > &matvec)
Creates a copy of the material object and put it in the vector which is passed on.
TPZAutoPointer< TPZFunction< STATE > > fForcingFunction
Pointer to forcing function, it is the right member at differential equation.
Definition: TPZMaterial.h:47
int NumLoadCases()
returns the number of load cases for this material object
Definition: TPZMaterial.h:186
int Id() const
Definition: TPZMaterial.h:170
void SetId(int id)
Definition: TPZMaterial.h:171
void SetTimeDependentFunctionExact(TPZAutoPointer< TPZFunction< STATE > > fp)
Sets a procedure as time dependent exact solution for the problem.
Definition: TPZMaterial.h:433
TPZAutoPointer< TPZFunction< STATE > > & TimedependentFunctionExact()
Returns a procedure as time dependent exact solution for the problem.
Definition: TPZMaterial.h:439
Defines the interface for saving and reading data. Persistency.
Definition: TPZStream.h:50
Contains the declaration of TPZFlopCounter class and TPZCounter struct.
void SetLinearContext(bool IsLinear)
Sets fLinearContext attribute.
Definition: TPZMaterial.cpp:77
def values
Definition: rdt.py:119
clarg::argInt porder("-porder", "polinomial order", 1)
This class defines the interface to save and restore objects from TPZStream objects. Persistency.
Definition: TPZSavable.h:67
Defines the interface of a computational element. Computational Element.
Definition: pzcompel.h:59
TPZSolVec sol
vector of the solutions at the integration point
TPZAutoPointer< TPZFunction< STATE > > & ForcingFunctionExact()
Returns a procedure as exact solution for the problem.
Definition: TPZMaterial.h:411
#define PZError
Defines the output device to error messages and the DebugStop() function.
Definition: pzerror.h:15
virtual void SetData(std::istream &data)
Reads data of the material from a istream (file data)
TPZVec< void(*)(const TPZVec< REAL > &, TPZVec< STATE > &) > GFORCINGVEC
Extern variable - Vector of force values.
Definition: TPZMaterial.cpp:25
TPZAutoPointer< TPZFunction< STATE > > fTimeDependentForcingFunction
Pointer to time dependent forcing function, it is the right member at differential equation...
Definition: TPZMaterial.h:53
This class implements a reference counter mechanism to administer a dynamically allocated object...