NeoPZ
swelling.h
Go to the documentation of this file.
1 
6 #ifndef SWELLINGHPP
7 #define SWELLINGHPP
8 
9 #include "TPZMaterial.h"
10 #include "pzfmatrix.h"
11 
12 #ifdef _AUTODIFF
13 #include "fadType.h"
14 #endif
15 
21 class TPZSwelling : public TPZMaterial {
22 
32  STATE fLambda;
34  STATE fGamma;
36  STATE fShear;
38  STATE fAlfa;
40  STATE fM;
42  STATE fDPlus;
44  STATE fDMinus;
46  STATE frHinder;
48  STATE fInitDeform;
50  STATE fCfc;
52  STATE fNf0;
54  STATE fNPlus0;
56  STATE fNMinus0;
57 
59  STATE fTheta;
61  STATE fDelt;
63  static STATE gExtConc;
64 
65 #ifdef _AUTODIFF
66 
67  static REAL gFaraday;
69  static REAL gVPlus;
71  static REAL gVMinus;
73  static REAL gRGas;
75  static REAL gTemp;
77  static REAL gMuRef[3];
78 
79 #else
80 
81  static STATE gFaraday;
83  static STATE gVPlus;
85  static STATE gVMinus;
87  static STATE gRGas;
89  static STATE gTemp;
91  static STATE gMuRef[3];
92 
93 #endif
94 
95  public :
96 
97 
115  TPZSwelling(int matindex, STATE lambda, STATE shear, STATE alfa, STATE M, STATE Gamma, STATE Kperm, STATE DPlus, STATE DMinus,
116  STATE rHinder, STATE Cfc, STATE Nf0, STATE NPlus0, STATE NMinus0);
117 
118  virtual ~TPZSwelling();
119 
121  virtual int Dimension() const override { return 3;}
122 
127  virtual int NStateVariables() const override{ return 8;}
128 
129  virtual void Print(std::ostream & out) override;
130 
131  virtual std::string Name() override { return "TPZSwelling"; }
132 
133  void SetComputationMode(int mode) {
134  switch(mode) {
135  case 0:
136  fComputationMode = mode;
137  break;
138  case 1:
139  fComputationMode = mode;
140  break;
141  default:
142  std::cout << "ComputationMode illegal mode = " << mode << std::endl;
143  }
144  }
145 
151  virtual void Contribute(TPZMaterialData &data,
152  REAL weight,
153  TPZFMatrix<STATE> &ek,
154  TPZFMatrix<STATE> &ef) override {
155  std::cout << "TPZSwelling::Contribute not implemented\n";
156  }
157  virtual void Contribute(TPZMaterialData &data,
158  REAL weight,
159  TPZFMatrix<STATE> &ef) override {
160  std::cout << "TPZSwelling::Contribute not implemented\n";
161  }
162 
163  virtual void ContributeBC(TPZMaterialData &data,
164  REAL weight,
165  TPZFMatrix<STATE> &ek,
166  TPZFMatrix<STATE> &ef,
167  TPZBndCond &bc) override;
168  virtual void ContributeBC(TPZMaterialData &data,
169  REAL weight,
170  TPZFMatrix<STATE> &ef,
171  TPZBndCond &bc) override
172  {
173  TPZMaterial::ContributeBC(data,weight,ef,bc);
174  }
175 
176 
177 #ifdef _AUTODIFF
178 
180  virtual void ContributeElastEnergy(
181  TPZVec<FADFADREAL> &dsol,
182  FADFADREAL &U,
183  REAL weight);
184 
186  virtual void ContributeResidual(TPZVec<REAL> & x,
187  TPZVec<FADREAL> & sol,
188  TPZVec<FADREAL> &dsol,
189  TPZFMatrix<REAL> &phi,
190  TPZFMatrix<REAL> &dphi,
191  TPZVec<FADREAL> &RES,
192  REAL weight);
193 
195  virtual void ContributePrevResidual(TPZVec<REAL> & x,
196  TPZVec<FADREAL> & sol,
197  TPZVec<FADREAL> &dsol,
198  TPZFMatrix<REAL> &phi,
199  TPZFMatrix<REAL> &dphi,
200  TPZVec<FADREAL> &RES,
201  REAL weight);
202 
204  virtual void ContributeBCEnergy(TPZVec<REAL> & x,
205  TPZVec<FADFADREAL> & sol, FADFADREAL &U,
206  REAL weight, TPZBndCond &bc) {
207  std::cout << "TPZSwelling::ContributeBCEnergy is not implemented\n";
208  }
209 
210 #endif
211 
214 private:
215 
216  void ExactSolution(TPZVec<STATE> &mu, STATE ksi, STATE pres, TPZVec<STATE> &N);
217 
220  void ComputeN(TPZVec<STATE> &N, TPZVec<STATE> &mu, STATE ksi);
221 
224  void ComputeInitialGuess(TPZVec<STATE> &mu, STATE J, STATE &pres, STATE &ksi, TPZVec<STATE> &N);
225 
226 #ifdef _AUTODIFF
227 
229  void ComputeW(FADFADREAL &W, TPZVec<STATE> &N);
230 
231 #endif
232 
235  void ComputeN(TPZVec<STATE> &mu, STATE ksi, STATE pressure, TPZVec<STATE> &N);
236 
237 #ifdef _AUTODIFF
238 
241 
244  void ComputeN(TPZVec<FADREAL> &sol, TPZVec<REAL> &N);
245 
246 #endif
247 
249  void NResidual(TPZVec<STATE> &mu, STATE ksi, STATE pressure, TPZVec<STATE> &N, TPZFMatrix<STATE> &res, TPZFMatrix<STATE> &tangent);
250 
251 #ifdef _AUTODIFF
252 
257 
263  static void Solve(TPZVec<TPZVec<FADREAL> > &tangent, TPZVec<FADREAL> &res);
264 public:
265 
266  static int main();
267 
268 #endif
269 
271  int NumCases();
273  void LoadState(TPZFMatrix<STATE> &state);
275  void ComputeTangent(TPZFMatrix<STATE> &tangent,TPZVec<REAL> &coefs, int cases);
277  void Residual(TPZFMatrix<STATE> &res, int cases);
282 
283 public:
284 
285  virtual int VariableIndex(const std::string &name) override;
288  virtual int NSolutionVariables(int var) override;
289 
290 protected:
292  virtual void Solution(TPZVec<STATE> &Sol,TPZFMatrix<STATE> &DSol,TPZFMatrix<REAL> &axes,int var,TPZVec<STATE> &Solout) override;
293 public:
294  virtual void Solution(TPZMaterialData &data,int var,TPZVec<STATE> &Solout) override
295  {
296  int numbersol = data.sol.size();
297  if (numbersol != 1) {
298  DebugStop();
299  }
300  Solution(data.sol[0],data.dsol[0],data.axes,var,Solout);
301  }
302  public:
303 virtual int ClassId() const override;
304 
305 };
306 
307 #endif
STATE fDPlus
Diffusion coeficient for cations [mm^2/s].
Definition: swelling.h:42
virtual ~TPZSwelling()
Definition: swelling.cpp:72
TPZFMatrix< STATE > fKperm
Hydraulic permeability .
Definition: swelling.h:30
virtual void ContributeBC(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ef, TPZBndCond &bc) override
It computes a contribution to the stiffness matrix and load vector at one BC integration point...
Definition: swelling.h:168
STATE fNMinus0
Initial anion volume fraction (no dimension)
Definition: swelling.h:56
virtual void Solution(TPZVec< STATE > &Sol, TPZFMatrix< STATE > &DSol, TPZFMatrix< REAL > &axes, int var, TPZVec< STATE > &Solout) override
Computes a post-processed solution variable corresponding to the variable index.
Definition: swelling.cpp:120
STATE fM
Storage modulus [N/mm^2].
Definition: swelling.h:40
STATE fDMinus
Diffusion coeficient for anions [mm^2/s].
Definition: swelling.h:44
clarg::argBool bc("-bc", "binary checkpoints", false)
int NumCases()
Methods needed to perform convergence checks.
virtual void Print(std::ostream &out) override
Prints out the data associated with the material.
Definition: swelling.cpp:76
static STATE gExtConc
External concentration (used as reference value for pressure) [mmol/mm^3].
Definition: swelling.h:63
STATE fShear
Shear modulus [N/mm^2].
Definition: swelling.h:36
TPZGradSolVec dsol
vector of the derivatives of the solution at the integration point
The TPZSwelling class implements a numerical model of swelling material coupling flow through porous ...
Definition: swelling.h:21
STATE fDelt
Timestep [s].
Definition: swelling.h:61
virtual int NStateVariables() const override
Number of state variables, in this case: 3 displacements, 1 pressure, 3 eletrochemical potencials...
Definition: swelling.h:127
This class implements a simple vector storage scheme for a templated class T. Utility.
Definition: pzgeopoint.h:19
static TPZFMatrix< REAL > gdphi
Definition: swelling.h:281
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.
Definition: swelling.h:294
void ExactSolution(TPZVec< STATE > &mu, STATE ksi, STATE pres, TPZVec< STATE > &N)
Definition: swelling.cpp:779
void ComputeTangent(TPZFMatrix< STATE > &tangent, TPZVec< REAL > &coefs, int cases)
Computes the tangent matrix for a given loadcase.
STATE fAlfa
Biot coupling coeficient (no dimension)
Definition: swelling.h:38
This abstract class defines the behaviour which each derived class needs to implement.
Definition: TPZMaterial.h:39
static STATE gTemp
Absolute temperature [K].
Definition: swelling.h:89
static STATE gRGas
gas constant [Nmm/(mmol K)]
Definition: swelling.h:87
STATE fCfc
Fixed charge density [mmoleq/mm^3].
Definition: swelling.h:50
int64_t size() const
Returns the number of elements of the vector.
Definition: pzvec.h:196
STATE fInitDeform
Initial deformation (assuming everything occurs isotropically, a constant is suficient (no dimension)...
Definition: swelling.h:48
STATE fNf0
Initial fluid volume fraction (do dimension)
Definition: swelling.h:52
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 std::string Name() override
Returns the name of the material.
Definition: swelling.h:131
STATE fTheta
Timestepping parameter theta (no dimension)
Definition: swelling.h:59
static STATE gFaraday
Faraday constant [C/mmol].
Definition: swelling.h:81
TPZFNMatrix< 9, REAL > axes
axes indicating the directions of the derivatives of the shapefunctions
void ComputeN(TPZVec< STATE > &N, TPZVec< STATE > &mu, STATE ksi)
Computes the value of the N coeficients in function of ksi and mus, iterative method, inverting the Hessian of W.
string res
Definition: test.py:151
int fComputationMode
Computation mode: 0->residual wrt to time 1->residual and tangent wrt to time .
Definition: swelling.h:28
static TPZFMatrix< REAL > gphi
Variables which holds the state variables used in the check convergence procedure.
Definition: swelling.h:281
STATE fGamma
Osmotic coeficient (no dimension)
Definition: swelling.h:34
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...
TPZSwelling(int matindex, STATE lambda, STATE shear, STATE alfa, STATE M, STATE Gamma, STATE Kperm, STATE DPlus, STATE DMinus, STATE rHinder, STATE Cfc, STATE Nf0, STATE NPlus0, STATE NMinus0)
Constructor of the class, where the user needs to specify the most important parameters.
Definition: swelling.cpp:45
void LoadState(TPZFMatrix< STATE > &state)
Loads the state within the current object, to be used when computing the tangent matrix.
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: swelling.h:151
STATE fNPlus0
Initial cation volume fraction (no dimension)
Definition: swelling.h:54
void ComputeInitialGuess(TPZVec< STATE > &mu, STATE J, STATE &pres, STATE &ksi, TPZVec< STATE > &N)
Computes the aproximate values of the pressure, ksi and N based on mu and J by direct inversion of th...
Definition: swelling.cpp:750
virtual int VariableIndex(const std::string &name) override
Definition: swelling.cpp:111
static STATE gMuRef[3]
Reference chemical potentials (order f,plus,minus) [mV].
Definition: swelling.h:91
virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ef) override
It computes a contribution to the residual vector at one integration point.
Definition: swelling.h:157
void SetComputationMode(int mode)
Definition: swelling.h:133
static STATE gVPlus
Molar volume cation [mm^3/mmol].
Definition: swelling.h:83
virtual int ClassId() const override
Define the class id associated with the class.
Definition: swelling.cpp:842
STATE frHinder
Hindrance factor (no dimension)
Definition: swelling.h:46
STATE fLambda
Compression modulus [N/mm^2].
Definition: swelling.h:32
virtual int Dimension() const override
Dimension of the problem.
Definition: swelling.h:121
static TPZFMatrix< STATE > gState
Variables which holds the state variables used in the check convergence procedure.
Definition: swelling.h:279
int main()
Definition: benchadolc.cc:42
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: swelling.cpp:128
TPZSolVec sol
vector of the solutions at the integration point
static STATE gVMinus
Molar volume anions [mm^3/mmol].
Definition: swelling.h:85
virtual int NSolutionVariables(int var) override
Returns the number of solution variables associated with a variable index.
Definition: swelling.cpp:115
void Residual(TPZFMatrix< STATE > &res, int cases)
Computes the residual for the given state variable.
void NResidual(TPZVec< STATE > &mu, STATE ksi, STATE pressure, TPZVec< STATE > &N, TPZFMatrix< STATE > &res, TPZFMatrix< STATE > &tangent)
Computes the residual and tangent vector of the system of equations which determines N...