/[escript]/trunk/doc/user/Models.tex
ViewVC logotype

Contents of /trunk/doc/user/Models.tex

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2251 - (show annotations)
Fri Feb 6 06:50:39 2009 UTC (11 years ago) by gross
File MIME type: application/x-tex
File size: 26922 byte(s)
new compressible flow solver
1
2 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
3 %
4 % Copyright (c) 2003-2008 by University of Queensland
5 % Earth Systems Science Computational Center (ESSCC)
6 % http://www.uq.edu.au/esscc
7 %
8 % Primary Business: Queensland, Australia
9 % Licensed under the Open Software License version 3.0
10 % http://www.opensource.org/licenses/osl-3.0.php
11 %
12 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
13
14
15 \chapter{Models}
16 \label{MODELS CHAPTER}
17
18 The following sections give a breif overview of the model classes and their corresponding methods.
19
20 \section{Stokes Problem}
21 The velocity \index{velocity} field $v$ and pressure $p$ of an incompressible fluid \index{incompressible fluid} is given as the solution of the Stokes problem\index{Stokes problem}
22 \begin{equation}\label{Stokes 1}
23 -\left(\eta(v\hackscore{i,j}+ v\hackscore{i,j})\right)\hackscore{,j}+p\hackscore{,i}=f\hackscore{i}-\sigma\hackscore{ij,j}
24 \end{equation}
25 where $\eta$ is the viscosity, $F\hackscore{i}$ defines an internal force \index{force, internal} and $\sigma\hackscore{ij}$ is an intial stress \index{stress, initial}. We assume an incompressible media:
26 \begin{equation}\label{Stokes 2}
27 -v\hackscore{i,i}=0
28 \end{equation}
29 Natural boundary conditions are taken in the form
30 \begin{equation}\label{Stokes Boundary}
31 \left(\eta(v\hackscore{i,j}+ v\hackscore{i,j})\right)n\hackscore{j}-n\hackscore{i}p=s\hackscore{i}+\sigma\hackscore{ij} n\hackscore{i}
32 \end{equation}
33 which can be overwritten by constraints of the form
34 \begin{equation}\label{Stokes Boundary0}
35 v\hackscore{i}(x)=v^D\hackscore{i}(x)
36 \end{equation}
37 at some locations $x$ at the boundary of the domain. The index $i$ may depend on the location $x$ on the boundary.
38 $v^D$ is a given function on the domain.
39
40 \subsection{Solution Method \label{STOKES SOLVE}}
41 In block form equation equations~\ref{Stokes 1} and~\ref{Stokes 2} takes the form of a saddle point problem
42 \index{saddle point problem}
43 \begin{equation}
44 \left[ \begin{array}{cc}
45 A & B^{*} \\
46 B & 0 \\
47 \end{array} \right]
48 \left[ \begin{array}{c}
49 v \\
50 p \\
51 \end{array} \right]
52 =\left[ \begin{array}{c}
53 G \\
54 0 \\
55 \end{array} \right]
56 \label{SADDLEPOINT}
57 \end{equation}
58 where $A$ is coercive, self-adjoint linear operator in a suitable Hilbert space, $B$ is the $(-1) \cdot$ divergence operator and $B^{*}$ is it adjoint operator (=gradient operator). For more details on the mathematics see references \cite{AAMIRBERKYAN2008,MBENZI2005}.
59 We use iterative techniques to solve this problem. To make sure that the incomressibilty condition holds
60 with sufficient accuracy we check for
61 \begin{equation}
62 \|v\hackscore{k,k}\| \hackscore \le \epsilon
63 \|\sqrt{v\hackscore{j,k}v\hackscore{j,k}}\|
64 \label{STOKES STOP}
65 \end{equation}
66 where $\epsilon$ is the desired relative accuracy and
67 \begin{equation}
68 \|p\|^2= \int\hackscore{\Omega} p^2 \; dx
69 \label{PRESSURE NORM}
70 \end{equation}
71 defines the $L^2$-norm. We use the Uzawa scheme \index{Uzawa scheme} to solve the problem.
72
73 In fact the first equation in~\ref{SADDLEPOINT} gives for a known pressure
74 \begin{equation}
75 v=A^{-1}(G-B^{*}p)
76 \label{V CALC}
77 \end{equation}
78 which is inserted into the second equation leading to
79 \begin{equation}
80 S p = B A^{-1} G
81 \end{equation}
82 with the Schur complement \index{Schur complement} $S=BA^{-1}B^{*}$. This problem can be solved iteratively
83 with the preconditioner $\hat{S}$ defined as $q=\hat{S}^{-1}p$ by solving
84 \begin{equation}
85 \frac{1}{\eta}q = p
86 \end{equation}
87 see~\cite{ELMAN} for more details. Note that the residual for the current approximation $p$ is given as
88 \begin{equation}
89 r=B A^{-1} (G - B^* p) = Bv
90 \end{equation}
91 where $v$ is given by~\ref{V CALC}.
92
93 If one uses the generalized minimal residual method (GMRES) \index{generalized minimal residual method!GMRES}
94 the method is directly applied to the preconditioned system
95 \begin{equation}
96 \hat{S}^{-1} S p = \hat{S}^{-1} B A^{-1} G
97 \end{equation}
98 We use the norm
99 \begin{equation}
100 \|p\|\hackscore{GMRES} = \|\hat{S} p \|
101 \end{equation}
102 Notice that for the residual $\hat{r}=\hat{S}^{-1} r$ one has
103 \begin{equation}
104 \
105 \end{equation}
106 If $p^{0}$ provides an initial guess for the pressure we use~\ref{V CALC} to get a first initial guess for the
107 velocity $v^{0}$ which we use to set an absolute tolerance $ATOL =\epsilon \|\sqrt{v^{0}\hackscore{j,k}v^{0}\hackscore{j,k}}\|$.
108 The GMRES is terminated when
109 \begin{equation}
110 \|\hat{r}\|\hackscore{GMRES} \le ATOL
111 \end{equation}
112 Notice that $\|\hat{r}\|\hackscore{GMRES}= \|r \| = \|Bv\| = \|v\hackscore{k,k}\|$ so we we can expect that
113 the target stopping criterion~\ref{STOKES STOP} is fullfilled. However, if $v$ is very different from the
114 initial choice of $v^{0}$ the value of $ATOL$ is corrected and GMRES is restarted with a new tolerance. For time dependend problems this apprach works well as value for $p$ form a previous time step provides a good initial guess.
115
116 Alternatively, as $S$ is symmetric and positive definite one can apply the preconditioned conjugate gradient method (PCG) \index{preconditioned conjugate gradient method!PCG}. PCG use the norm
117 \begin{equation}
118 \|r\|\hackscore{PCG}^2 = \int\hackscore{\Omega} r \hat{S}^{-1}r \; dx = \int\hackscore{\Omega} \eta r^2 \; dx
119 \end{equation}
120 To take the extra factor $\eta$ into consideration when checking the stopping criterion we use the following
121 definition for $ATOL$:
122 \begin{equation}
123 ATOL = \epsilon \frac{\|\sqrt{v^{0}\hackscore{j,k}v^{0}\hackscore{j,k}}\| }{\|v^{0}\hackscore{k,k}\|}
124 \|v^{0}\hackscore{k,k}\|\hackscore{PCG}
125 \end{equation}
126
127
128
129 \subsection{Functions}
130
131 \begin{classdesc}{StokesProblemCartesian}{domain}
132 opens the Stokes problem\index{Stokes problem} on the \Domain domain. The approximation
133 order needs to be two.
134 \end{classdesc}
135
136 \begin{methoddesc}[StokesProblemCartesian]{initialize}{\optional{f=Data(), \optional{fixed_u_mask=Data(), \optional{eta=1, \optional{surface_stress=Data(), \optional{stress=Data()}}}}}}
137 assigns values to the model parameters. In any call all values must be set.
138 \var{f} defines the external force $f$, \var{eta} the viscosity $\eta$,
139 \var{surface_stress} the surface stress $s$ and \var{stress} the initial stress $\sigma$.
140 The locations and compontents where the velocity is fixed are set by
141 the values of \var{fixed_u_mask}. The method will try to cast the given values to appropriate
142 \Data class objects.
143 \end{methoddesc}
144
145 \begin{methoddesc}[StokesProblemCartesian]{solve}{v,p,
146 \optional{max_iter=20, \optional{verbose=False, \optional{usePCG=True}}}}
147 solves the problem and return approximations for velocity and pressure.
148 The arguments \var{v} and \var{p} define initial guess. The values of \var{v} marked
149 by \var{fixed_u_mask} remain unchanged.
150 If \var{usePCG} is set to \True
151 reconditioned conjugate gradient method (PCG) \index{preconditioned conjugate gradient method!PCG} scheme is used. Otherwise the problem is solved generalized minimal residual method (GMRES) \index{generalized minimal residual method!GMRES}. In most cases
152 the PCG scheme is more efficient.
153 \var{max_iter} defines the maximum number of iteration steps.
154 If \var{verbose} is set to \True informations on the progress of of the solver are printed.
155 \end{methoddesc}
156
157
158 \begin{methoddesc}[StokesProblemCartesian]{setTolerance}{\optional{tolerance=1.e-4}}
159 sets the tolerance in an appropriate norm relative to the right hand side. The tolerance must be non-negative and less than 1.
160 \end{methoddesc}
161 \begin{methoddesc}[StokesProblemCartesian]{getTolerance}{}
162 returns the current relative tolerance.
163 \end{methoddesc}
164 \begin{methoddesc}[StokesProblemCartesian]{setAbsoluteTolerance}{\optional{tolerance=0.}}
165 sets the absolute tolerance for the error in the relevant norm. The tolerance must be non-negative. Typically the
166 absolute talerance is set to 0.
167 \end{methoddesc}
168 \begin{methoddesc}[StokesProblemCartesian]{getAbsoluteTolerance}{}
169 sreturns the current absolute tolerance.
170 \end{methoddesc}
171 \begin{methoddesc}[StokesProblemCartesian]{setSubProblemTolerance}{\optional{rtol=None}}
172 sets the tolerance to solve the involved PDEs. The subtolerance \var{rtol} should not be choosen to large
173 in order to avoid feed back of errors in the subproblem solution into the outer iteration.
174 On the otherhand is choosen to small compute time is wasted.
175 If \var{rtol} is set to \var{None} the sub-tolerance is set automatically depending on the
176 tolerance choosen for the oter iteration.
177 \end{methoddesc}
178 \begin{methoddesc}[StokesProblemCartesian]{getSubProblemTolerance}{}
179 return the tolerance for the involved PDEs.
180 \end{methoddesc}
181
182 \subsection{Example: Lit Driven Cavity}
183 The following script \file{lit\hackscore driven\hackscore cavity.py}
184 \index{scripts!\file{helmholtz.py}} which is available in the \ExampleDirectory
185 illustrates the usage of the \class{StokesProblemCartesian} class to solve
186 the lit driven cavity problem~\cite{LITDRIVENCAVITY}:
187 \begin{python}
188 from esys.escript import *
189 from esys.finley import Rectangle
190 from esys.escript.models import StokesProblemCartesian
191 NE=25
192 dom = Rectangle(NE,NE,order=2)
193 x = dom.getX()
194 sc=StokesProblemCartesian(dom)
195 mask= (whereZero(x[0])*[1.,0]+whereZero(x[0]-1))*[1.,0] + \
196 (whereZero(x[1])*[0.,1.]+whereZero(x[1]-1))*[1.,1]
197 sc.initialize(eta=.1, fixed_u_mask= mask)
198 v=Vector(0.,Solution(dom))
199 v[0]+=whereZero(x[1]-1.)
200 p=Scalar(0.,ReducedSolution(dom))
201 v,p=sc.solve(v,p, verbose=True)
202 saveVTK("u.xml",velocity=v,pressure=p)
203 \end{python}
204
205 \section{Darcy Flux}
206 We want to calculate the velocity $u$ and pressure $p$ on a domain $\Omega$ solving
207 the Darcy flux problem \index{Darcy flux}\index{Darcy flow}
208 \begin{equation}\label{DARCY PROBLEM}
209 \begin{array}{rcl}
210 u\hackscore{i} + \kappa\hackscore{ij} p\hackscore{,j} & = & g\hackscore{i} \\
211 u\hackscore{k,k} & = & f
212 \end{array}
213 \end{equation}
214 with the boundary conditions
215 \begin{equation}\label{DARCY BOUNDARY}
216 \begin{array}{rcl}
217 u\hackscore{i} \; n\hackscore{i} = u^{N}\hackscore{i} \; n\hackscore{i} & \mbox{ on } & \Gamma\hackscore{N} \\
218 p = p^{D} & \mbox{ on } & \Gamma\hackscore{D} \\
219 \end{array}
220 \end{equation}
221 where $\Gamma\hackscore{N}$ and $\Gamma\hackscore{D}$ are a partition of the boundary of $\Omega$ with $\Gamma\hackscore{D}$ non empty, $n\hackscore{i}$ is the outer normal field of the boundary of $\Omega$, $u^{N}\hackscore{i}$ and $p^{D}$ are given functions on $\Omega$, $g\hackscore{i}$ and $f$ are given source terms and $\kappa\hackscore{ij}$ is the given permability. We assume that $\kappa\hackscore{ij}$ is symmetric (which is not really required) and positive definite, i.e there are positive constants $\alpha\hackscore{0}$ and $\alpha\hackscore{1}$ wich are independent from the location in $\Omega$ such that
222 \begin{equation}
223 \alpha\hackscore{0} \; x\hackscore{i} x\hackscore{i} \le \kappa\hackscore{ij} x\hackscore{i} x\hackscore{j} \le \alpha\hackscore{1} \; x\hackscore{i} x\hackscore{i}
224 \end{equation}
225 for all $x\hackscore{i}$.
226
227
228 \subsection{Solution Method \label{DARCY SOLVE}}
229 In practical applications it is an advantage to calculate the pressure $p$ as a correction of a 'static' pressure $p^{ref}$ which is the solution of
230 \begin{equation}
231 -(\kappa\hackscore{ki}\kappa\hackscore{kj} p\hackscore{,j}^{ref})\hackscore{,i} = - (\kappa\hackscore{ki} (g\hackscore{k}- u^{N}\hackscore{k}))\hackscore{,i}
232 \mbox{ with }
233 p^{ref} = p^{D} \mbox{ on } \Gamma\hackscore{D}
234 \end{equation}
235 With setting $u \leftarrow u-u^{N}$ and $p \leftarrow p-p^{ref}$ and
236 \begin{equation}
237 \begin{array}{rcl}
238 g\hackscore{i} & \leftarrow & g\hackscore{i} - u^{N}\hackscore{i} - \kappa\hackscore{ij} p^{ref}\hackscore{,j }\\
239 f & \leftarrow & f - u^{N}\hackscore{k,k}
240 \end{array}
241 \end{equation}
242 we can assume that $u^{N}\hackscore{i} \; n\hackscore{i}=0$ and
243 $p^{D}=0$.
244 We set
245 \begin{equation}
246 V = \{ q \in H^{1}(\Omega) : q=0 \mbox{ on } \Gamma\hackscore{D} \}
247 \end{equation}
248 and
249 \begin{equation}
250 W = \{ v \in (L^2(\Omega))^{d} : v\hackscore{k,k} \in L^2(\Omega) \mbox{ and } u\hackscore{i} \; n\hackscore{i} =0 \mbox{ on } \Gamma\hackscore{N} \}
251 \end{equation}
252 and define the operator $Q: V \rightarrow (L^2(\Omega))^{d}$ defined by
253 \begin{equation}
254 (Qp)\hackscore{i} = \kappa\hackscore{ij} p\hackscore{,j}
255 \end{equation}
256 and the operator $D: W \rightarrow L^2(\Omega)$ defined by
257 \begin{equation}
258 Dv = v\hackscore{k,k}
259 \end{equation}
260 In operator notation the Darcy problem~\ref{DARCY PROBLEM} is written in the form
261 \begin{equation}
262 \begin{array}{rcl}
263 u + Qp & = & g \\
264 Du & = & f
265 \end{array}
266 \end{equation}
267 We solve this equation by minimising the functional
268 \begin{equation}
269 J(u,p):=\|u + Qp - g\|^2\hackscore{0} + \|Du-f\|\hackscore{0}^2
270 \end{equation}
271 over $W \times V$ where $\|.\|\hackscore{0}$ denotes the norm in $L^2(\Omega)$. A simple calculation shows that
272 one has to solve
273 \begin{equation}
274 ( v + Qq , u + Qp - g) + (Dv,Du-f) =0
275 \end{equation}
276 for all $v\in W$ and $q \in V$.which translates back into operator notation
277 \begin{equation}
278 \begin{array}{rcl}
279 (I+D^*D)u + Qp & = & D^*f + g \\
280 Q^*u + Q^*Q p & = & Q^*g \\
281 \end{array}
282 \end{equation}
283 where $D^*$ and $Q^*$ denote the adjoint operators.
284 In~\cite{XXX} it has been shown that this problem is continuous and coercive in $W \times V$ and therefore has a unique solution. Also standart FEM methods can be used for discretization. It is also possible
285 to solve the problem is coupled form, however this approach leads in some cases to a very ill-conditioned stiffness matrix in particular in the case of a very small or large permability ($\alpha\hackscore{1} \ll 1$ or $\alpha\hackscore{0} \gg 1$)
286
287 The approach we are taking is to eliminate the velocity $u$ from the problem. Assuming that $p$ is known we have
288 \begin{equation}
289 v= (I+D^*D)^{-1} (D^*f + g - Qp)
290 \end{equation}
291 (notice that $(I+D^*D)$ is coercive in $W$) which is inserted into the second equation
292 \begin{equation}
293 Q^* (I+D^*D)^{-1} (D^*f + g - Qp) + Q^* Q p = Q^*g
294 \end{equation}
295 which is
296 \begin{equation}
297 Q^* ( I - (I+D^*D)^{-1} ) Q p = Q^* (g-(I+D^*D)^{-1} (D^*f + g) )
298 \end{equation}
299 We use the PCG \index{linear solver!PCG}\index{PCG} method to solve this. The residual $r$ ($\in V^*$) is given as
300 \begin{equation}
301 \begin{array}{rcl}
302 r & = & Q^* \left( g -(I+D^*D)^{-1} (D^*f + g) - Qp + (I+D^*D)^{-1}Q p \right)\\
303 & =& Q^* \left( - Qp - (I+D^*D)^{-1} (D^*f + g - Qp) \right) \\
304 & =& Q^* \left( g - Qp - v \right)
305 \end{array}
306 \end{equation}
307 So in a partical implementation we use the pair $(Qp,v)$ to represent the residual. This will save the
308 reconstruction of the velocity $v$. In this notation the right hand side is given as
309 $(0,(I+D^*D)^{-1} (D^*f + g))$. The evaluation of the iteration operator for a given $p$ is then
310 returning $(Qp,w)$ where $w$ is the solution of
311 \begin{equation}\label{UPDATE W}
312 (I+D^*D)w = Qp
313 \end{equation}
314 We use $Q^*Q$ as a a preconditioner for the iteration operator $Q^* ( I - (I+D^*D)^{-1} ) Q$.
315
316 The iteration PCG \index{linear solver!PCG}\index{PCG} is terminated if
317 \begin{equation}\label{DARCY STOP}
318 \int r \cdot (Q^*Q)^{-1} r \; dx \le \mbox{ATOL}^2
319 \end{equation}
320 where ATOL is a given absolute tolerance.
321 The initial residual $r\hackscore{0}$ is
322 \begin{equation}\label{DARCY STOP 2}
323 r\hackscore{0}=Q^* \left( g - v\hackscore{ref} \right) \mbox{ with } v\hackscore{ref} = (I+D^*D)^{-1} (D^*f + g)
324 \end{equation}
325 so the
326 \begin{equation}\label{DARCY NORM 0}
327 \int r\hackscore0 \cdot (Q^*Q)^{-1} r\hackscore0 \; dx = \int \left( g - v\hackscore{ref} \right) \cdot Q p\hackscore{ref} \; dx \mbox{ with }p\hackscore{ref} = (Q^*Q)^{-1} Q^* \left( g - v\hackscore{ref} \right)
328 \end{equation}
329 So we set
330 \begin{equation}\label{DARCY NORM 1}
331 ATOL = atol + rtol \cdot \max(|g - v\hackscore{ref}|\hackscore{0}, |Q p\hackscore{ref} |\hackscore{0} )
332 \end{equation}
333 where atol and rtol a given absolute and relative tolerances, respectively. The reference flux $v\hackscore{ref}$
334 and reference pressure $p\hackscore{ref}$ may be calcualated from their definition which would require to solve to
335 PDEs but in a practical application estimates can be used for instance solutions from previous time steps or for simplified scenarious (e.g. constant permability).
336
337 \subsection{Functions}
338 \begin{classdesc}{DarcyFlow}{domain}
339 opens the Darcy flux problem\index{Darcy flux} on the \Domain domain.
340 \end{classdesc}
341 \begin{methoddesc}[DarcyFlow]{setValue}{\optional{f=None, \optional{g=None, \optional{location_of_fixed_pressure=None, \optional{location_of_fixed_flux=None, \optional{permeability=None}}}}}}
342 assigns values to the model parameters. Values can be assigned using various calls - in particular
343 in a time dependend problem only values that change over time needs to be reset. The permability can be defined as scalar (isotropic), a vector (orthotropic) or a matrix (anisotropic).
344 \var{f} and \var{g} are the corresponding parameters in~\ref{DARCY PROBLEM}.
345 The locations and compontents where the flux is prescribed are set by positive values in
346 \var{location_of_fixed_flux}.
347 The locations where the pressure is prescribed are set by
348 by positive values of \var{location_of_fixed_pressure}.
349 The values of the pressure and flux are defined by the initial guess.
350 Notice that at any point on the boundary of the domain the pressure or the normal component of
351 the flux must be defined. There must be at least one point where the pressure is prescribed.
352 The method will try to cast the given values to appropriate
353 \Data class objects.
354 \end{methoddesc}
355
356 \begin{methoddesc}[DarcyFlow]{setTolerance}{\optional{atol=0,\optional{rtol=1e-8,\optional{p_ref=None,\optional{v_ref=None}}}}}
357 sets the absolute tolerance ATOL according to~\ref{DARCY NORM 1}. If \var{p_ref} is not present $0$ is used.
358 If \var{v_ref} is not present $0$ is used. If the final result ATOL is not positive an exception is thrown.
359 \end{methoddesc}
360
361
362
363 \begin{methoddesc}[DarcyFlow]{solve}{u0,p0, \optional{max_iter=100, \optional{verbose=False \optional{sub_rtol=1.e-8}}}}
364 solves the problem. and returns approximations for the flux $v$ and the pressure $p$.
365 \var{u0} and \var{p0} define initial guess for flux and pressure. Values marked
366 by positive values \var{location_of_fixed_flux} and \var{location_of_fixed_pressure}, respectively, are kept unchanged.
367 \var{sub_rtol} defines the tolerance used to solve the involved PDEs. \var{sub_rtol} needs to be choosen sufficiently small to ensure convergence but users need to keep in mind that a very small value for \var{sub_rtol} will result in a long compute time. Typically $\var{sub_rtol}=\var{rtol}^2$ is a good choice if $\var{rtol}$ is not choosen too small.
368 \end{methoddesc}
369
370
371 \subsection{Example: Gravity Flow}
372 later
373
374 %================================================
375 % \section{Temperature Advection Diffusion\label{TEMP ADV DIFF}}
376
377 %\begin{equation}
378 % \rho c\hackscore{p} \left (\frac{\partial T}{\partial t} + \vec{v} \cdot \nabla T \right ) = k \nabla^{2}T
379 % \label{HEAT EQUATION}
380 % \end{equation}
381
382 % where $\vec{v}$ is the velocity vector, $T$ is the temperature, $\rho$ is the density, $\eta$ is the viscosity, % % $c\hackscore{p}$ is the specific heat at constant pressure and $k$ is the thermal conductivity.
383
384 % \subsection{Description}
385
386 % \subsection{Method}
387 %
388 % \begin{classdesc}{TemperatureCartesian}{dom,theta=THETA,useSUPG=SUPG}
389 % \end{classdesc}
390
391 % \subsection{Benchmark Problem}
392 %===============================================================================================================
393
394 %=========================================================
395 \input{levelsetmodel}
396
397 % \section{Drucker Prager Model}
398
399 \section{Isotropic Kelvin Material \label{IKM}}
400 As proposed by Kelvin~\ref{KELVN} material strain $D\hackscore{ij}=v\hackscore{i,j}+v\hackscore{j,i}$ can be decomposed into
401 an elastic part $D\hackscore{ij}^{el}$ and visco-plastic part $D\hackscore{ij}^{vp}$:
402 \begin{equation}\label{IKM-EQU-2}
403 D\hackscore{ij}=D\hackscore{ij}^{el}+D\hackscore{ij}^{vp}
404 \end{equation}
405 with the elastic strain given as
406 \begin{equation}\label{IKM-EQU-3}
407 D\hackscore{ij}'^{el}=\frac{1}{2 \mu} \dot{\sigma}\hackscore{ij}'
408 \end{equation}
409 where $\sigma'\hackscore{ij}$ is the deviatoric stress (Notice that $\sigma'\hackscore{ii}=0$).
410 If the material is composed by materials $q$ the visco-plastic strain can be decomposed as
411 \begin{equation}\label{IKM-EQU-4}
412 D\hackscore{ij}'^{vp}=\sum\hackscore{q} D\hackscore{ij}'^{q}
413 \end{equation}
414 where $D\hackscore{ij}^{q}$ is the strain in material $q$ given as
415 \begin{equation}\label{IKM-EQU-5}
416 D\hackscore{ij}'^{q}=\frac{1}{2 \eta^{q}} \sigma'\hackscore{ij}
417 \end{equation}
418 where $\eta^{q}$ is the viscosity of material $q$. We assume the following
419 betwee the the strain in material $q$
420 \begin{equation}\label{IKM-EQU-5b}
421 \eta^{q}=\eta^{q}\hackscore{N} \left(\frac{\tau}{\tau\hackscore{t}^q}\right)^{\frac{1}{n^{q}}-1}
422 \mbox{ with } \tau=\sqrt{\frac{1}{2}\sigma'\hackscore{ij} \sigma'\hackscore{ij}}
423 \end{equation}
424 for a given power law coefficients $n^{q}$ and transition stresses $\tau\hackscore{t}^q$, see~\ref{POERLAW}.
425 Notice that $n^{q}=1$ gives a constant viscosity.
426 After inserting equation~\ref{IKM-EQU-5} into equation \ref{IKM-EQU-4} one gets:
427 \begin{equation}\label{IKM-EQU-6}
428 D\hackscore{ij}'^{vp}=\frac{1}{2 \eta^{vp}} \sigma'\hackscore{ij} \mbox{ with } \frac{1}{\eta^{vp}} = \sum\hackscore{q} \frac{1}{\eta^{q}} \;.
429 \end{equation}
430 With
431 \begin{equation}\label{IKM-EQU-8}
432 \dot{\gamma}=\sqrt{2 D\hackscore{ij} D\hackscore{ij}}
433 \end{equation}
434 one gets
435 \begin{equation}\label{IKM-EQU-8b}
436 \tau = \eta^{vp} \dot{\gamma}^{vp} \;.
437 \end{equation}
438 With the Drucker-Prager cohesion factor $\tau\hackscore{Y}$, Drucker-Prager friction $\beta$ and total pressure $p$ we want to achieve
439 \begin{equation}\label{IKM-EQU-8c}
440 \tau \le \tau\hackscore{Y} + \beta \; p
441 \end{equation}
442 which leads to the condition
443 \begin{equation}\label{IKM-EQU-8d}
444 \eta^{vp} \le \frac{\tau\hackscore{Y} + \beta \; p}{ \dot{\gamma}^{vp}} \; .
445 \end{equation}
446 Therefore we modify the definition of $\eta^{vp}$ to the form
447 \begin{equation}\label{IKM-EQU-6b}
448 \frac{1}{\eta^{vp}}=\max(\sum\hackscore{q} \frac{1}{\eta^{q}}, \frac{\dot{\gamma}^{vp}} {\tau\hackscore{Y} + \beta \; p})
449 \end{equation}
450 The deviatoric stress needs to fullfill the equilibrion equation
451 \begin{equation}\label{IKM-EQU-1}
452 -\sigma'\hackscore{ij,j}+p\hackscore{,i}=F\hackscore{i}
453 \end{equation}
454 where $F\hackscore{j}$ is a given external fource. We assume an incompressible media:
455 \begin{equation}\label{IKM-EQU-2}
456 -v\hackscore{i,i}=0
457 \end{equation}
458 Natural boundary conditions are taken in the form
459 \begin{equation}\label{IKM-EQU-Boundary}
460 \sigma'\hackscore{ij}n\hackscore{j}-n\hackscore{i}p=f
461 \end{equation}
462 which can be overwritten by a constraint
463 \begin{equation}\label{IKM-EQU-Boundary0}
464 v\hackscore{i}(x)=0
465 \end{equation}
466 where the index $i$ may depend on the location $x$ on the bondary.
467
468 \subsection{Solution Method \label{IKM-SOLVE}}
469 By using a first order finite difference approximation wit step size $dt>0$~\ref{IKM-EQU-3} get the form
470 \begin{equation}\label{IKM-EQU-3b}
471 D\hackscore{ij}'^{el}=\frac{1}{2 \mu dt } \left( \sigma\hackscore{ij}' - \sigma\hackscore{ij}^{'-} \right)
472 \end{equation}
473 where $\sigma\hackscore{ij}^{'-}$ is the deviatoric stress at the precious time step.
474 Now we can combine equations~\ref{IKM-EQU-2}, \ref{IKM-EQU-3b} and~\ref{IKM-EQU-6b} to get
475 \begin{equation}\label{IKM-EQU-10}
476 \sigma\hackscore{ij}' = 2 \eta\hackscore{eff} \left( D\hackscore{ij}' +
477 \frac{1}{ 2 \mu \; dt} \sigma\hackscore{ij}^{'-}\right) \mbox{ with }
478 \frac{1}{\eta\hackscore{eff}}=\frac{1}{\mu \; dt}+\frac{1}{\eta^{vp}}
479 \end{equation}
480 Notice that $\eta\hackscore{eff}$ is a function of diatoric stress $\sigma\hackscore{ij}'$.
481 After inserting~\ref{IKM-EQU-10} into~\ref{IKM-EQU-1} we get
482 \begin{equation}\label{IKM-EQU-1ib}
483 -\left(\eta\hackscore{eff} (v\hackscore{i,j}+ v\hackscore{i,j})
484 \right)\hackscore{,j}+p\hackscore{,i}=F\hackscore{i}+
485 \frac{\eta\hackscore{eff}}{\mu dt } \sigma\hackscore{ij,j}^{'-}
486 \end{equation}
487 Together with the incomressibilty condition~\ref{IKM-EQU-2} we need to solve a problem with a form almost identical
488 to the Stokes problem discussed in section~\ref{STOKES SOLVE} but with the difference that $\eta\hackscore{eff}$ is depending on the solution. Analog to the iteration scheme~\ref{SADDLEPOINT iteration 2} we can run
489 \begin{equation}
490 \begin{array}{rcl}
491 -\left(\eta\hackscore{eff}(dv\hackscore{i,j}+ dv\hackscore{i,j}
492 )\right)\hackscore{,j} & = & F\hackscore{i}+ \sigma\hackscore{ij,j}'-p\hackscore{,i} \\
493 \frac{1}{\eta\hackscore{eff}} dp & = & - v\hackscore{i,i}^+
494 \end{array}
495 \label{IKM iteration 2}
496 \end{equation}
497 where $v^+=v+dv$. As this problem is non-linear the Jacobi-free Newton-GMRES method is used with the norm
498 \begin{equation}
499 \|(v, p)\|^2= \int\hackscore{\Omega} v\hackscore{i,j}^2 + \frac{1}{\bar{\eta}^2\hackscore{eff}} p^2 \; dx
500 \label{IKM iteration 3}
501 \end{equation}
502 where $\bar{\eta}\hackscore{eff}$ is the caracteristic viscosity, for instance:
503 \begin{equation}
504 \frac{1}{\bar{\eta}\hackscore{eff}} = \frac{1}{\tau^{-}}+\sum\hackscore{q} \frac{1}{\eta^{q}\hackscore{N}}
505 \label{IKM iteration 4}
506 \end{equation}
507 In oder to perform step~\ref{IKM iteration 2} we need to calculate the $\eta\hackscore{eff}$ as well as $\sigma\hackscore{ij}'$ while via $\tau$ the first is a function of the latter. The priority is the
508 calculation of $\eta\hackscore{eff}$ with the Newton-Raphson scheme. This value can then be used to calculate
509 $\sigma\hackscore{ij}'$ via~\ref{IKM-EQU-10}. We need to solve
510 \begin{equation}
511 \tau = \eta\hackscore{eff} \cdot \epsilon \mbox{ with }
512 \epsilon = \sqrt{ 2 \left( D\hackscore{ij}' +
513 \frac{1}{ 2 \mu \; dt} \sigma\hackscore{ij}^{'-}\right)^2}
514 \label{IKM iteration 5}
515 \end{equation}
516 The Newton scheme takes the form
517 \begin{equation}
518 \tau\hackscore{n+1} = \min(\tau\hackscore{n} - \frac{\tau\hackscore{n} - \eta\hackscore{eff} \cdot \epsilon}{1 - \eta\hackscore{eff}' \cdot \epsilon}, \tau\hackscore{Y} + \beta \; p)
519 = \min(\frac{\eta\hackscore{eff} - \tau\hackscore{n} \eta\hackscore{eff}'}
520 {1 - \eta\hackscore{eff}' \cdot \epsilon}, \frac{\tau\hackscore{Y} + \beta \; p}{\epsilon}) \epsilon
521 \label{IKM iteration 6}
522 \end{equation}
523 where $\eta\hackscore{eff}'$ denotes the derivative of $\eta\hackscore{eff}$ with respect of $\tau$. The second term in $\min$ is droped of $\tau\hackscore{Y} + \beta \; p<0$ or $\epsilon=0$. In fact we have
524 \begin{equation}
525 \eta\hackscore{eff}' = - \eta\hackscore{eff}^2 \left(\frac{1}{\eta\hackscore{eff}}\right)'
526 \mbox{ with }
527 \left(\frac{1}{\eta\hackscore{eff}}\right)' = \sum\hackscore{q} \left(\frac{1}{\eta^{q}} \right)'
528 \label{IKM iteration 7}
529 \end{equation}
530 \begin{equation}\label{IKM-EQU-5XX}
531 \left(\frac{1}{\eta^{q}} \right)'
532 = \frac{1-\frac{1}{n^{q}}}{\eta^{q}\hackscore{N}} \frac{\tau^{-\frac{1}{n^{q}}}}{(\tau\hackscore{t}^q)^{1-\frac{1}{n^{q}}}}
533 = \frac{1-\frac{1}{n^{q}}}{ \tau \eta^{q}}
534 \end{equation}
535 Notice that allways $\eta\hackscore{eff}'\le 0$ which makes the denomionator in~\ref{IKM iteration 6}
536 positive.
537
538
539

  ViewVC Help
Powered by ViewVC 1.1.26