NeoPZ
TPZMixedElasticityMaterial.h
Go to the documentation of this file.
1 
6 #ifndef MIXEDELASMATHPP
7 #define MIXEDELASMATHPP
8 
9 #include <iostream>
10 
11 #include "TPZMaterial.h"
12 #include "pzdiscgal.h"
13 
19 
21  {
23  REAL fE;
24 
26  REAL fnu;
27 
29  REAL flambda;
30 
32  REAL fmu;
33 
34  TElasticityAtPoint(REAL E, REAL nu) : fE(E), fnu(nu)
35  {
36  flambda = (E * nu) / ((1 + nu)*(1 - 2 * nu));
37  fmu = E / (2 * (1 + nu));
38  }
39  };
40 
41 
42 
43 public:
44 
45  enum MVoigt {
47  };
48 
60  TPZMixedElasticityMaterial(int id, REAL E, REAL nu, REAL fx, REAL fy, int planestress = 1, int fDimension = 1);
61 
63 
65 
68 
69 
70  virtual int VariableIndex(const std::string &name) override;
71 
72 
73  virtual int NSolutionVariables(int var) override;
74 
76  int SIndex() {
77  return 0;
78  }
79 
81  int UIndex() {
82  return 1;
83  }
84 
86  int PIndex() {
87  return 2;
88  }
89 
90  virtual int NEvalErrors() override;
91 
93 
94  void ComputeDeformationVector(TPZVec<STATE> &PhiStress, TPZVec<STATE> &APhiStress, TElasticityAtPoint &elastb);
95 
96  void ComputeStressVector(TPZVec<STATE> &Deformation, TPZVec<STATE> &Stress, TElasticityAtPoint &elast);
97 
99 
101  virtual TPZMaterial * NewMaterial() override {
102  return new TPZMixedElasticityMaterial(*this);
103  }
104 
106  virtual ~TPZMixedElasticityMaterial();
107 
115  void SetElasticParameters(REAL Eyoung, REAL nu) {
116  this->SetElasticity(Eyoung, nu);
117  }
118 
120  void SetElasticity(REAL E, REAL nu) {
121  fE_const = E; // Young modulus
122  fnu_const = nu; // poisson coefficient
123  flambda_const = (E * nu) / ((1 + nu)*(1 - 2 * nu));
124  fmu_const = E / (2 * (1 + nu));
125 
126  }
127 
130  {
131  fElasticity = func;
132  }
133 
135  void SetLameParameters(REAL Lambda, REAL mu) {
136  fE_const = (mu * (3.0 * Lambda + 2.0 * mu)) / (Lambda + mu);
137  fnu_const = (Lambda) / (2 * (Lambda + mu));
138 
139  flambda_const = Lambda;
140  fmu_const = mu;
141  }
142 
144 
146  fAxisSymmetric = 1;
147  }
149 
150  void SetPlaneStrain() {
151  fPlaneStress = 0;
152  }
153 
155 
156  void SetPlaneStress() {
157  fPlaneStress = 1;
158  }
159 
161  void SetBodyForce(REAL fx, REAL fy) {
162  fForce[0] = fx;
163  fForce[1] = fy;
164  fForce[2] = 0.;
165  }
166 
168  int Dimension() const override {
169  return 2;
170  }
171 
173  virtual int NStateVariables() const override;
174 
176  virtual void Print(std::ostream & out = std::cout) override;
177 
179  virtual std::string Name() override {
180  return "TPZMixedElasticityMaterial";
181  }
182 
184  virtual short NumberOfFluxes() {
185  return 3;
186  }
187 
189  virtual int NFluxes() override {
190  return 3;
191  }
192 
197  virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef) override;
198 
200  virtual void Contribute(TPZVec<TPZMaterialData> &data, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef) override;
201 
203  // virtual void Contribute(TPZVec<TPZMaterialData> &data, REAL weight, TPZFMatrix<STATE> &ef)
204  // {
205  // DebugStop();
206  // }
207 
209  virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix<STATE> &ef) override {
210  TPZDiscontinuousGalerkin::Contribute(data, weight, ef);
211  }
212 
213 
215  virtual void ContributeBC(TPZVec<TPZMaterialData> &datavec, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef, TPZBndCond &bc) override;
216 
217 
219  virtual void ContributeBC(TPZMaterialData &data, REAL weight,
220  TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef, TPZBndCond &bc) override;
221 
223  virtual void ContributeBC(TPZMaterialData &data, REAL weight,
224  TPZFMatrix<STATE> &ef, TPZBndCond &bc) override {
225  TPZDiscontinuousGalerkin::ContributeBC(data, weight, ef, bc);
226  }
227 
228  //virtual void FillDataRequirements(TPZMaterialData &data);
229  virtual void FillDataRequirements(TPZMaterialData &data) override;
230  virtual void FillDataRequirements(TPZVec<TPZMaterialData > &datavec) override;
231 
232  virtual void FillBoundaryConditionDataRequirement(int type, TPZMaterialData &data) override;
233 
234  virtual void ContributeInterface(TPZMaterialData &data, TPZMaterialData &dataleft, TPZMaterialData &dataright, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef) override {
235  PZError << "\nFATAL ERROR - Method not implemented: " << __PRETTY_FUNCTION__ << "\n";
236  }
237 
238  virtual void ContributeBCInterface(TPZMaterialData &data, TPZMaterialData &dataleft, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef, TPZBndCond &bc) override {
239  PZError << "\nFATAL ERROR - Method not implemented: " << __PRETTY_FUNCTION__ << "\n";
240  }
241 
242  virtual void ContributeInterface(TPZMaterialData &data, TPZMaterialData &dataleft, TPZMaterialData &dataright, REAL weight, TPZFMatrix<STATE> &ef) override {
243  PZError << "\nFATAL ERROR - Method not implemented: " << __PRETTY_FUNCTION__ << "\n";
244  }
245 
246  virtual void ContributeBCInterface(TPZMaterialData &data, TPZMaterialData &left, REAL weight, TPZFMatrix<STATE> &ef, TPZBndCond &bc) override {
247  PZError << "\nFATAL ERROR - Method not implemented: " << __PRETTY_FUNCTION__ << "\n";
248  }
249 
252 
254  static void ToVoigt(TPZFMatrix<STATE> &S, TPZVec<STATE> &Svoigt);
255 
257  static void FromVoigt(TPZVec<STATE> &Svoigt, TPZFMatrix<STATE> &S);
258 
260  template<class TVar>
261  TVar InnerVec(const TPZVec<TVar> &S, const TPZVec<TVar> &T);
262 
264  STATE Tr(TPZFMatrix<REAL> &GradU);
265 
267  STATE Transpose(TPZFMatrix<REAL> &GradU);
268 
270  void FillGradPhi(TPZMaterialData &dataV, TPZVec< TPZFMatrix<STATE> > &GradPhi);
271 
274 
275 public:
276 
278  virtual void Solution(TPZMaterialData &data, int var, TPZVec<STATE> &Solout) override;
279 
281  virtual void Solution(TPZVec<TPZMaterialData> &data, int var, TPZVec<STATE> &Solout) override;
282 
284  virtual void SolutionDisc(TPZMaterialData &data, TPZMaterialData &dataleft, TPZMaterialData &dataright, int var, TPZVec<STATE> &Solout) {
285  TPZDiscontinuousGalerkin::SolutionDisc(data, dataleft, dataright, var, Solout) ;
286  }
287 
289  virtual void Flux(TPZVec<REAL> &x, TPZVec<STATE> &Sol, TPZFMatrix<STATE> &DSol, TPZFMatrix<REAL> &axes, TPZVec<STATE> &flux) override;
290 
295  void Errors(TPZVec<REAL> &x, TPZVec<STATE> &u,
297  TPZVec<STATE> &u_exact, TPZFMatrix<STATE> &du_exact, TPZVec<REAL> &values) override; //Cedric
298 
299  virtual void Errors(TPZVec<TPZMaterialData> &data, TPZVec<STATE> &u_exact, TPZFMatrix<STATE> &du_exact, TPZVec<REAL> &errors) override;
300 
301  virtual int ClassId() const override;
302 
303  virtual void Read(TPZStream &buf, void *context) override;
304 
305  virtual void Write(TPZStream &buf, int withclassid) const override;
306 
307 
308 
309 protected:
312 
314  int fPlaneStress = 1;
315 
317  int fDimension = 2;
318 
321 
323  REAL fE_const;
324 
326  REAL fnu_const;
327 
330 
332  REAL fmu_const;
333 
334 
335  // Matrix A
337 
339  bool fAxisSymmetric = false;
340 
341 };
342 
343 #endif
Defines the interface which material objects need to implement for discontinuous Galerkin formulation...
Definition: pzdiscgal.h:20
virtual void FillDataRequirements(TPZMaterialData &data) override
Fill material data parameter with necessary requirements for the.
void SolutionDisc(TPZMaterialData &data, TPZMaterialData &dataleft, TPZMaterialData &dataright, int var, TPZVec< STATE > &Solout)
Returns the solution associated with the var index based on the finite element approximation.
Definition: pzdiscgal.h:176
STATE Inner(TPZFMatrix< STATE > &S, TPZFMatrix< STATE > &T)
void SetLameParameters(REAL Lambda, REAL mu)
Set elasticity parameters.
void FillVecShapeIndex(TPZMaterialData &data)
transform a H1 data structure to a vector data structure
clarg::argBool bc("-bc", "binary checkpoints", false)
bool fAxisSymmetric
flag indicates axix-AxisSymmetric
void SetElasticityFunction(TPZAutoPointer< TPZFunction< STATE > > func)
Set a variable elasticity and poisson coefficient.
TPZMixedElasticityMaterial()
Default constructor.
virtual void SolutionDisc(TPZMaterialData &data, TPZMaterialData &dataleft, TPZMaterialData &dataright, int var, TPZVec< STATE > &Solout)
Returns the solution associated with the var index based on the finite element approximation.
TPZAutoPointer< TPZFunction< STATE > > fElasticity
This class implements a simple vector storage scheme for a templated class T. Utility.
Definition: pzgeopoint.h:19
virtual void FillBoundaryConditionDataRequirement(int type, TPZMaterialData &data) override
This method defines which parameters need to be initialized in order to compute the contribution of t...
void ElasticityModulusTensor(TPZFMatrix< STATE > &MatrixElast, TElasticityAtPoint &elast)
virtual void Flux(TPZVec< REAL > &x, TPZVec< STATE > &Sol, TPZFMatrix< STATE > &DSol, TPZFMatrix< REAL > &axes, TPZVec< STATE > &flux) override
Computes the value of the flux function to be used by ZZ error estimator.
This class implements a two dimensional elastic material in plane stress or strain.
void SetElasticity(REAL E, REAL nu)
Set elasticity parameters.
virtual short NumberOfFluxes()
Returns the number of components which form the flux function.
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...
This abstract class defines the behaviour which each derived class needs to implement.
Definition: TPZMaterial.h:39
virtual void ContributeBC(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ef, TPZBndCond &bc) override
Applies the element boundary conditions.
static void FromVoigt(TPZVec< STATE > &Svoigt, TPZFMatrix< STATE > &S)
Transform a Voigt notation to a tensor.
STATE Tr(TPZFMatrix< REAL > &GradU)
virtual int NSolutionVariables(int var) override
int Dimension() const override
Returns the model dimension.
REAL flambda_const
first Lame Parameter
virtual void Read(TPZStream &buf, void *context) override
Reads the element data from a stream.
virtual void Write(TPZStream &buf, int withclassid) const override
Saves the element data to a stream.
virtual int ClassId() const override
Unique identifier for serialization purposes.
void ComputeDivergenceOnDeformed(TPZVec< TPZMaterialData > &datavec, TPZFMatrix< STATE > &DivergenceofPhi)
virtual void ContributeBCInterface(TPZMaterialData &data, TPZMaterialData &dataleft, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc) override
It computes a contribution to stiffness matrix and load vector at one BC integration point...
virtual void ContributeBC(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc) override=0
It computes a contribution to the stiffness matrix and load vector at one BC integration point...
virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ef) override
Calculates the element stiffness matrix - simulate compaction as aditional variable.
void SetElasticParameters(REAL Eyoung, REAL nu)
Set parameters of elastic material:
This class defines the boundary condition for TPZMaterial objects.
Definition: pzbndcond.h:29
virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef) override
Calculates the element stiffness matrix.
virtual void ContributeInterface(TPZMaterialData &data, TPZMaterialData &dataleft, TPZMaterialData &dataright, REAL weight, TPZFMatrix< STATE > &ef) override
It computes a contribution to residual vector at one integration point.
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...
int fDimension
dimension of the material
void ComputeDeformationVector(TPZVec< STATE > &PhiStress, TPZVec< STATE > &APhiStress, TElasticityAtPoint &elastb)
virtual ~TPZMixedElasticityMaterial()
Default destructor.
virtual TPZMaterial * NewMaterial() override
Creates a new material from the current object ??
void SetAxisSymmetric()
set the material configuration to AxisSymmetric
virtual void ContributeBCInterface(TPZMaterialData &data, TPZMaterialData &left, REAL weight, TPZFMatrix< STATE > &ef, TPZBndCond &bc) override
It computes a contribution to residual vector at one BC integration point.
void SetPlaneStrain()
Set the material configuration to plane strain.
STATE Transpose(TPZFMatrix< REAL > &GradU)
void SetPlaneStress()
Set the material configuration to plane stress.
virtual void ContributeBC(TPZVec< TPZMaterialData > &datavec, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc) override
Applies the element boundary conditions Mixed.
void SetBodyForce(REAL fx, REAL fy)
Set forcing function.
virtual int NFluxes() override
Returns the number of components which form the flux function.
virtual std::string Name() override
Returns the material name.
virtual int NEvalErrors() override
Returns the number of norm errors. Default is 3: energy, L2 and H1.
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.
static void ToVoigt(TPZFMatrix< STATE > &S, TPZVec< STATE > &Svoigt)
Transform a tensor to a Voigt notation.
REAL fmu_const
Second Lame Parameter.
virtual int NStateVariables() const override
Returns the number of state variables associated with the material.
void ComputeStressVector(TPZVec< STATE > &Deformation, TPZVec< STATE > &Stress, TElasticityAtPoint &elast)
TVar InnerVec(const TPZVec< TVar > &S, const TPZVec< TVar > &T)
Defines the interface for saving and reading data. Persistency.
Definition: TPZStream.h:50
def values
Definition: rdt.py:119
REAL fnu_const
Poison coefficient.
virtual void Print(std::ostream &out=std::cout) override
Print the material data.
Contains the TPZDiscontinuousGalerkin class which implements the interface for discontinuous Galerkin...
virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef) override=0
It computes a contribution to the stiffness matrix and load vector at one integration point...
virtual int VariableIndex(const std::string &name) override
TPZManVector< REAL, 3 > fForce
Forcing vector.
#define PZError
Defines the output device to error messages and the DebugStop() function.
Definition: pzerror.h:15
void FillGradPhi(TPZMaterialData &dataV, TPZVec< TPZFMatrix< STATE > > &GradPhi)
This class implements a reference counter mechanism to administer a dynamically allocated object...