1 |
% $Id$ |
2 |
|
3 |
\chapter{The module \linearPDEsPack} |
4 |
|
5 |
\declaremodule{extension}{linearPDEs} \modulesynopsis{Linear partial pifferential equation handler} |
6 |
The module \linearPDEsPack provides an interface to define and solve linear partial |
7 |
differential equations within \escript. \linearPDEsPack does not provide any |
8 |
solver capabilities in itself but hands the PDE over to |
9 |
the PDE solver library defined through the \Domain of the PDE. |
10 |
The general interface is provided through the \LinearPDE class. The |
11 |
\AdvectivePDE which is derived from the \LinearPDE class |
12 |
provides an interface to PDE dominated by its advective terms. The \Poisson |
13 |
class which is also derived form the \LinearPDE class should be used |
14 |
to define the Poisson equation \index{Poisson}. |
15 |
|
16 |
\section{\LinearPDE Class} |
17 |
\label{SEC LinearPDE} |
18 |
|
19 |
The \LinearPDE class is used to define a general linear, steady, second order PDE |
20 |
for an unknown function $u$ on a given $\Omega$ defined through a \Domain object. |
21 |
In the following $\Gamma$ denotes the boundary of the domain $\Omega$. $n$ denotes |
22 |
the outer normal field on $\Gamma$. |
23 |
|
24 |
For a single PDE with a solution with a single component the linear PDE is defined in the |
25 |
following form: |
26 |
\begin{equation}\label{LINEARPDE.SINGLE.1} |
27 |
-(A\hackscore{jl} u\hackscore{,l}){,j}+(B\hackscore{j} u)\hackscore{,j}+C\hackscore{l} u\hackscore{,l}+D u =-X\hackscore{j,j}+Y \; . |
28 |
\end{equation} |
29 |
$u_{,j}$ denotes the derivative of $u$ with respect to the $j$-th spatial direction. Einstein's summation convention, ie. summation over indexes appearing twice in a term of a sum is performed, is used. |
30 |
The coefficients $A$, $B$, $C$, $D$, $X$ and $Y$ have to be specified through \Data objects in the |
31 |
\Function on the PDE or objects that can be converted into such \Data objects. |
32 |
$A$ is a \RankTwo, $B$, $C$ and $X$ are \RankOne and $D$ and $Y$ are scalar. |
33 |
The following natural |
34 |
boundary conditions are considered \index{boundary condition!natural} on $\Gamma$: |
35 |
\begin{equation}\label{LINEARPDE.SINGLE.2} |
36 |
n\hackscore{j}(A\hackscore{jl} u\hackscore{,l}+B\hackscore{j} u)+d u=n\hackscore{j}X\hackscore{j} + y \;. |
37 |
\end{equation} |
38 |
Notice that the coefficients $A$, $B$ and $X$ are defined in the PDE. The coefficients $d$ and $y$ are |
39 |
each a \Scalar in the \FunctionOnBoundary. Constraints \index{constraint} for the solution prescribing the value of the |
40 |
solution at certain locations in the domain. They have the form |
41 |
\begin{equation}\label{LINEARPDE.SINGLE.3} |
42 |
u=r \mbox{ where } q>0 |
43 |
\end{equation} |
44 |
$r$ and $q$ are each \Scalar where $q$ is the characteristic function |
45 |
\index{characteristic function} defining where the constraint is applied. |
46 |
The constraints defined by \eqn{LINEARPDE.SINGLE.3} override any other condition set by \eqn{LINEARPDE.SINGLE.1} |
47 |
or \eqn{LINEARPDE.SINGLE.2}. The PDE is symmetrical \index{symmetrical} if |
48 |
\begin{equation}\label{LINEARPDE.SINGLE.4} |
49 |
A\hackscore{jl}=A\hackscore{lj} \mbox{ and } B\hackscore{j}=C\hackscore{j} |
50 |
\end{equation} |
51 |
For a system of PDEs and a solution with several components the PDE has the form |
52 |
\begin{equation}\label{LINEARPDE.SYSTEM.1} |
53 |
-(A\hackscore{ijkl} u\hackscore{k,l}){,j}+(B\hackscore{ijk} u_k)\hackscore{,j}+C\hackscore{ikl} u\hackscore{k,l}+D\hackscore{ik} u_k =-X\hackscore{ij,j}+Y\hackscore{i} \; . |
54 |
\end{equation} |
55 |
$A$ is a \RankFour, $B$ and $C$ are each a \RankThree, $D$ and $X$ are each a \RankTwo and $Y$ is a \RankOne. |
56 |
The natural boundary conditions \index{boundary condition!natural} take the form: |
57 |
\begin{equation}\label{LINEARPDE.SYSTEM.2} |
58 |
n\hackscore{j}(A\hackscore{ijkl} u\hackscore{k,l}){,j}+(B\hackscore{ijk} u_k)+d\hackscore{ik} u_k=n\hackscore{j}-X\hackscore{ij}+y\hackscore{i} \;. |
59 |
\end{equation} |
60 |
The coefficient $d$ is a \RankTwo and $y$ is a |
61 |
\RankOne both in the \FunctionOnBoundary. Constraints \index{constraint} take the form |
62 |
\begin{equation}\label{LINEARPDE.SYSTEM.3} |
63 |
u\hackscore{i}=r\hackscore{i} \mbox{ where } q\hackscore{i}>0 |
64 |
\end{equation} |
65 |
$r$ and $q$ are each \RankOne. Notice that at some locations not necessarily all components must |
66 |
have a constraint. The system of PDEs is symmetrical \index{symmetrical} if |
67 |
\begin{eqnarray}\label{LINEARPDE.SYSTEM.4} |
68 |
A\hackscore{ijkl}=A\hackscore{klij} \\ |
69 |
B\hackscore{ijk}=C\hackscore{kij} \\ |
70 |
D\hackscore{ik}=D\hackscore{ki} \\ |
71 |
d\hackscore{ik}=d\hackscore{ki} \ |
72 |
\end{eqnarray} |
73 |
\LinearPDE also supports solution discontinuities \index{discontinuity} over contact region $\Gamma^{contact}$ |
74 |
in the domain $\Omega$. To specify the conditions across the discontinuity we are using the |
75 |
generalised flux $J$ which is in the case of a systems of PDEs and several components of the solution |
76 |
defined as |
77 |
\begin{equation}\label{LINEARPDE.SYSTEM.5} |
78 |
J\hackscore{ij}=A\hackscore{ijkl}u\hackscore{k,l}+B\hackscore{ijk}u\hackscore{k}-X\hackscore{ij} |
79 |
\end{equation} |
80 |
For the case of single solution component and single PDE $J$ is defined |
81 |
\begin{equation}\label{LINEARPDE.SINGLE.5} |
82 |
J\hackscore{j}=A\hackscore{jl}u\hackscore{,l}+B\hackscore{j}u\hackscore{k}-X\hackscore{j} |
83 |
\end{equation} |
84 |
In the context of discontinuities \index{discontinuity} $n$ denotes the normal on the |
85 |
discontinuity pointing from side 0 towards side 1. For a system of PDEs |
86 |
the contact condition takes the form |
87 |
\begin{equation}\label{LINEARPDE.SYSTEM.6} |
88 |
n\hackscore{j} J^{0}\hackscore{ij}=n\hackscore{j} J^{1}\hackscore{ij}=y^{contact}\hackscore{i} - d^{contact}\hackscore{ik} [u]\hackscore{k} \; . |
89 |
\end{equation} |
90 |
where $J^{0}$ and $J^{1}$ are the fluxes on side $0$ and side $1$ of the |
91 |
discontinuity $\Gamma^{contact}$, respectively. $[u]$, which is the difference |
92 |
of the solution at side 1 and at side 0, denotes the jump of $u$ across $\Gamma^{contact}$. |
93 |
The coefficient $d^{contact}$ is a \RankTwo and $y^{contact}$ is a |
94 |
\RankOne both in the \FunctionOnContactZero or \FunctionOnContactOne. |
95 |
In case of a single PDE and a single component solution the contact condition takes the form |
96 |
\begin{equation}\label{LINEARPDE.SINGLE.6} |
97 |
n\hackscore{j} J^{0}\hackscore{j}=n\hackscore{j} J^{1}\hackscore{j}=y^{contact} - d^{contact}[u] |
98 |
\end{equation} |
99 |
In this case the the coefficient $d^{contact}$ and $y^{contact}$ are eaach \Scalar |
100 |
both in the \FunctionOnContactZero or \FunctionOnContactOne. |
101 |
|
102 |
\begin{classdesc}{LinearPDE}{domain,numEquations=0,numSolutions=0} |
103 |
opens a linear, steady, second order PDE on the \Domain \var{domain}. \var{numEquations} |
104 |
and \var{numSolutions} gives the number of equations and the number of solutiopn components. |
105 |
If \var{numEquations} and \var{numSolutions} is non-positive, the number of equations |
106 |
and the number solutions, respctively, stay undefined until a coefficient is |
107 |
defined. |
108 |
\end{classdesc} |
109 |
|
110 |
\begin{methoddesc}[LinearPDE]{setValues}{arg1,arg2,...,argN} |
111 |
assigns new values to coefficients \var{arg1}, \var{arg2}, $\ldots$ \var{argN}. |
112 |
The coefficients must be one of the valid coefficient names |
113 |
\var{A}, |
114 |
\var{B}, |
115 |
\var{C}, |
116 |
\var{D}, |
117 |
\var{X}, |
118 |
\var{Y}, |
119 |
\var{d}, |
120 |
\var{y}, |
121 |
\var{dcontact}, |
122 |
\var{ycontact}, |
123 |
\var{r}, |
124 |
or \var{q}. |
125 |
If the coefficient is not a \Data object, it is converted into a \Data object in the |
126 |
appropriate \FunctionSpace. |
127 |
\end{methoddesc} |
128 |
|
129 |
\begin{methoddesc}[LinearPDE]{createCoefficient}{name} |
130 |
facilitates the creation of a \Data object corresponding to coefficient of name |
131 |
\var{name}. If \var{name} is not a valid name an exception is raised. The |
132 |
returned \Data object is initialised to zero and is not set as a value of the |
133 |
LinearPDE. |
134 |
\end{methoddesc} |
135 |
|
136 |
\begin{methoddesc}[LinearPDE]{getCoefficient}{name} |
137 |
return the value assigned to coefficient \var{name}. If \var{name} is not a valid name |
138 |
an exception is raised. |
139 |
\end{methoddesc} |
140 |
|
141 |
\begin{methoddesc}[LinearPDE]{cleanCoefficients}{} |
142 |
resets all coefficients to their initialization values. This method is useful to call when trying to save memory |
143 |
as is releaces eqnerences to coefficients but keeping the differential operator. |
144 |
\end{methoddesc} |
145 |
\begin{methoddesc}[LinearPDE]{createNewCoefficient}{name} |
146 |
returns a \Data object which has the \FunctionSpace and \Shape of coefficient \var{name}. The initial value is $0$. |
147 |
\end{methoddesc} |
148 |
|
149 |
\begin{methoddesc}[LinearPDE]{getShapeOfCoefficient}{name} |
150 |
returns the shape of coefficient \var{name} even if no value has been assigned to it. |
151 |
\end{methoddesc} |
152 |
|
153 |
|
154 |
\begin{methoddesc}[LinearPDE]{getFunctionSpaceOfCoefficient}{name} |
155 |
returns the \FunctionSpace of coefficient \var{name} even if no value has been assigned to it. |
156 |
\end{methoddesc} |
157 |
|
158 |
\begin{methoddesc}[LinearPDE]{hasCoefficient}{name} |
159 |
returns \True if \var{name} is valid name of a coefficient |
160 |
\end{methoddesc} |
161 |
|
162 |
|
163 |
\begin{methoddesc}[LinearPDE]{getFunctionSpaceForEquation}{} |
164 |
returns \FunctionSpace of the right hand side which is identical to \FunctionSpace of |
165 |
the result of the PDE operator. |
166 |
\end{methoddesc} |
167 |
|
168 |
\begin{methoddesc}[LinearPDE]{getFunctionSpaceForSolution}{} |
169 |
returns \FunctionSpace of the solution of the PDE hich is identical to \FunctionSpace of |
170 |
the result of argument of the PDE operator. |
171 |
\end{methoddesc} |
172 |
|
173 |
\begin{methoddesc}[LinearPDE]{setDebugOn}{} |
174 |
switches the debug mode to on. |
175 |
\end{methoddesc} |
176 |
|
177 |
\begin{methoddesc}[LinearPDE]{setDebugOff}{} |
178 |
switches the debug mode to on. |
179 |
\end{methoddesc} |
180 |
|
181 |
\begin{methoddesc}[LinearPDE]{debug}{} |
182 |
returns \True if the debug mode is switched on. Otherwise it return \False. |
183 |
\end{methoddesc} |
184 |
|
185 |
\begin{methoddesc}[LinearPDE]{setLumpingOn}{} |
186 |
switches on lumping \index{lumping}. If lumping is switched on |
187 |
the operator is |
188 |
\end{methoddesc} |
189 |
|
190 |
|
191 |
\begin{methoddesc}[LinearPDE]{setLumpingOff} |
192 |
switches lumping off \index{lumping}. |
193 |
\end{methoddesc} |
194 |
|
195 |
\begin{methoddesc}[LinearPDE]{setLumping}{flag=\False} |
196 |
switches on lumping if \var{flag} is \True. Otherwise lumping is swiched off. |
197 |
\end{methoddesc} |
198 |
|
199 |
\begin{methoddesc}[LinearPDE]{isUsingLumping}{} |
200 |
returns \True if lumping is switched on. Otherwise \False is returned. |
201 |
\end{methoddesc} |
202 |
|
203 |
\begin{memberdesc}[LinearPDE]{DEFAULT_METHOD} |
204 |
default method to be used to solve the PDE. An appropriate method should be |
205 |
chosen by the used PDE solver library. |
206 |
\end{memberdesc} |
207 |
|
208 |
\begin{memberdesc}[LinearPDE]{DIRECT} |
209 |
direct linear solver~\Ref{SAAD} |
210 |
\end{memberdesc} |
211 |
|
212 |
\begin{memberdesc}[LinearPDE]{CHOLEVSKY} |
213 |
direct solver based on Cholevsky factorization (or similar), see~\Ref{SAAD}. The solver will require a symmetric PDE. |
214 |
\end{memberdesc} |
215 |
|
216 |
\begin{memberdesc}[LinearPDE]{PCG} |
217 |
preconditioned conjugate gradient method, see~\Ref{WEISS}. The solver will require a symmetric PDE. |
218 |
\end{memberdesc} |
219 |
|
220 |
\begin{memberdesc}[LinearPDE]{GMRES} |
221 |
the GMRES method, see~\Ref{WEISS}. Truncation and restart ar econtrolled by the parameters |
222 |
\var{truncation} and \var{restart} of \method{getSolution}. |
223 |
\end{memberdesc} |
224 |
|
225 |
\begin{memberdesc}[LinearPDE]{PRES20} |
226 |
the GMRES method with trunction after five residuals and |
227 |
restart after 20 steps, see~\Ref{WEISS}. |
228 |
\end{memberdesc} |
229 |
|
230 |
\begin{memberdesc}[LinearPDE]{CR} |
231 |
conjugate residual method, see~\Ref{WEISS}. |
232 |
\end{memberdesc} |
233 |
|
234 |
\begin{memberdesc}[LinearPDE]{CGS} |
235 |
conjugate gradient squared method, see~\Ref{WEISS}. |
236 |
\end{memberdesc} |
237 |
|
238 |
\begin{memberdesc}[LinearPDE]{BICGSTAB} |
239 |
stabilzed bi-conjugate gradients methods, see~\Ref{WEISS}. |
240 |
\end{memberdesc} |
241 |
|
242 |
\begin{memberdesc}[LinearPDE]{SSOR} |
243 |
symmetric successive overrelaxtion method, see~\Ref{WEISS}. |
244 |
\end{memberdesc} |
245 |
|
246 |
\begin{methoddesc}[LinearPDE]{setSolverMethod}{solver=linearPDE.DEFAULT_METHOD} |
247 |
sets the solver method to be used. It is pointed out that the PDE solver library |
248 |
does not know the specified solver method but may choose a similar method. |
249 |
\end{methoddesc} |
250 |
|
251 |
\begin{methoddesc}[LinearPDE]{setTolerance}{tol=1.e-8} |
252 |
resets the tolerance for solution. The actually meaning of tolerance is |
253 |
depending on the underlying PDE library. In most cases, the tolerance |
254 |
will only consider the error from solving the discerete problem but will |
255 |
not consider any discretization error. |
256 |
\end{methoddesc} |
257 |
|
258 |
\begin{methoddesc}[LinearPDE]{getTolerance}{} |
259 |
returns the current tolerance of the solution |
260 |
\end{methoddesc} |
261 |
|
262 |
\begin{methoddesc}[LinearPDE]{isSymmetric}{} |
263 |
returns \True if the PDE has been indicated to be symmetric. |
264 |
Otherwise \False is returned. |
265 |
\end{methoddesc} |
266 |
|
267 |
\begin{methoddesc}[LinearPDE]{setSymmetryOn}{} |
268 |
indicates that the PDE is symmetric. |
269 |
\end{methoddesc} |
270 |
|
271 |
\begin{methoddesc}[LinearPDE]{setSymmetryOff}{} |
272 |
indicates that the PDE is not symmetric. |
273 |
\end{methoddesc} |
274 |
|
275 |
\begin{methoddesc}[LinearPDE]{setSymmetryTo}{flag=\False} |
276 |
indicates that the PDE is symmetric if \var{flag}=\True |
277 |
and indicates a non-symmetric PDE is \var{flag}=\False. |
278 |
\end{methoddesc} |
279 |
|
280 |
\begin{methoddesc}[LinearPDE]{setReducedOrderOn}{} |
281 |
switches on the reduction of polynomial order for the solution and |
282 |
equation evaluation. |
283 |
\end{methoddesc} |
284 |
|
285 |
\begin{methoddesc}[LinearPDE]{setReducedOrderOff}{} |
286 |
switches off the reduction of polynomial order for the solution and |
287 |
equation evaluation. |
288 |
\end{methoddesc} |
289 |
|
290 |
\begin{methoddesc}[LinearPDE]{setReducedOrderTo}{flag=\False} |
291 |
switches on the reduction of polynomial order for the solution and |
292 |
equation evaluation if \var{flag}=\True. Otherwise |
293 |
the order reduction is switched off. |
294 |
\end{methoddesc} |
295 |
|
296 |
\begin{methoddesc}[LinearPDE]{setReducedOrderForSolutionOn}{} |
297 |
switches on reduction of polynomial order for the solution. |
298 |
\end{methoddesc} |
299 |
|
300 |
\begin{methoddesc}[LinearPDE]{setReducedOrderForSolutionOff}{} |
301 |
switches off reduction of polynomial order for the solution. |
302 |
\end{methoddesc} |
303 |
|
304 |
\begin{methoddesc}[LinearPDE]{setReducedOrderForSolutionTo}{flag=\False} |
305 |
switches on the reduction of polynomial order for the solution |
306 |
if \var{flag}=\True. Otherwise |
307 |
the order reduction is switched off. |
308 |
\end{methoddesc} |
309 |
|
310 |
\begin{methoddesc}[LinearPDE]{setReducedOrderForEquationOn}{} |
311 |
switches on reduction of polynomial order for the equation evaluation. |
312 |
\end{methoddesc} |
313 |
|
314 |
\begin{methoddesc}[LinearPDE]{setReducedOrderForEquationOff}{} |
315 |
switches off reduction of polynomial order for the equation evaluation. |
316 |
\end{methoddesc} |
317 |
|
318 |
\begin{methoddesc}[LinearPDE]{setReducedOrderForEquationTo}{flag=\False} |
319 |
switches on the reduction of polynomial order for the equation |
320 |
evaluation if \var{flag}=\True. Otherwise |
321 |
the order reduction is switched off. |
322 |
\end{methoddesc} |
323 |
|
324 |
\begin{methoddesc}[LinearPDE]{getOperator}{} |
325 |
returns the \Operator of the PDE. |
326 |
\end{methoddesc} |
327 |
|
328 |
\begin{methoddesc}[LinearPDE]{getRightHandSide}{ignoreConstraint=\False} |
329 |
returns the right hand side of the PDE as a \Data object. If |
330 |
\var{ignoreConstraint}=\True the constraints are not considered |
331 |
when building up the right hand side. |
332 |
\end{methoddesc} |
333 |
|
334 |
\begin{methoddesc}[LinearPDE]{getSystem}{} |
335 |
returns the \Operator and right hand side of the PDE. |
336 |
\end{methoddesc} |
337 |
|
338 |
\begin{methoddesc}[LinearPDE]{getSolution}{option1,option2,...,optionN} |
339 |
returns (an approximation of) the solution of the PDE. \var{options1}, \var{options2} |
340 |
$\ldots$ \var{optionsN} are options handed over to the underlying PDE solver library. |
341 |
\end{methoddesc} |
342 |
|
343 |
\begin{methoddesc}[LinearPDE]{getDomain}{} |
344 |
returns the \Domain of the PDE. |
345 |
\end{methoddesc} |
346 |
|
347 |
\begin{methoddesc}[LinearPDE]{getDim}{} |
348 |
returns the spatial dimension of the PDE. |
349 |
\end{methoddesc} |
350 |
|
351 |
\begin{methoddesc}[LinearPDE]{getNumEquations}{} |
352 |
returns the number of equations. |
353 |
\end{methoddesc} |
354 |
|
355 |
\begin{methoddesc}[LinearPDE]{getNumSolutions}{} |
356 |
returns the number of components of the solution. |
357 |
\end{methoddesc} |
358 |
|
359 |
\begin{methoddesc}[LinearPDE]{checkSymmetry}{verbose=\False} |
360 |
returns \True if the PDE is symmetric and \False otherwise. |
361 |
The method is very computational expensive and should only be |
362 |
called for testing purposes. The symmetry flag is not altered. |
363 |
If \var{verbose}=\True information about where symmetry is violated |
364 |
are printed. |
365 |
\end{methoddesc} |
366 |
|
367 |
\begin{methoddesc}[LinearPDE]{getFlux}{u} |
368 |
returns the flux $J\hackscore{ij}$ \index{flux} for given solution \var{u} |
369 |
defined by \eqn{LINEARPDE.SYSTEM.5} and \eqn{LINEARPDE.SINGLE.5}, respectively. |
370 |
\end{methoddesc} |
371 |
|
372 |
\begin{methoddesc}[LinearPDE]{applyOperator}{u} |
373 |
applies the PDE operator to \var{u} |
374 |
\end{methoddesc} |
375 |
|
376 |
\begin{methoddesc}[LinearPDE]{getResidual}{u} |
377 |
returns the residual when insering \var{u} into the PDE |
378 |
\end{methoddesc} |
379 |
|
380 |
\section{\AdvectivePDE Class} |
381 |
In cases of PDEs dominated by the advection terms $B$ and $C$ against the diffusion term $A$ |
382 |
up-winding has been used. |
383 |
The \AdvectivePDE class applies upwinding to the advective terms, see \Ref{SUPG}. |
384 |
To measure the dominance of the advective terms over the diffusive term $A$ the |
385 |
Pelclet number is used \index{Pelclet number}. The are defined as |
386 |
\begin{eqnarray}\label{LINEARPDE.Peclet.single} |
387 |
P^{B}=\frac{h\|B\hackscore{:}\|}{2\|A\hackscore{::}\|} |
388 |
\mbox{ and } |
389 |
P^{C}=\frac{h\|C\hackscore{:}\|}{2\|A\hackscore{::}\|} |
390 |
\end{eqnarray} |
391 |
\begin{eqnarray}\label{LINEARPDE.Peclet.system} |
392 |
P^{B}_{ik}=\frac{h\|B\hackscore{i:k}\|}{2\|A\hackscore{i:k:}\|} |
393 |
\mbox{ and } |
394 |
P^{C}_{ik}=\frac{h\|C\hackscore{ik:}\|}{2\|A\hackscore{i:k:}\|} |
395 |
\end{eqnarray} |
396 |
where $h$ is the local cell size and |
397 |
\begin{eqnarray}\label{LINEARPDE.ADVECTIVE.1b} |
398 |
\|C\hackscore{:}\|^2=C\hackscore{j}C\hackscore{j} \\ |
399 |
\|A\hackscore{::}\|^2=A\hackscore{jl}A\hackscore{jl} \\ |
400 |
\|C\hackscore{i:k}\|^2=C\hackscore{ijk}C\hackscore{ijk} \\ |
401 |
\|A\hackscore{i:k:}\|^2=A\hackscore{ijkl}A\hackscore{ijkl} \; . |
402 |
\end{eqnarray} |
403 |
From the Pelclet number the stabilization parameters $\Xi^{B}$ and $\Xi^{C}$ are calculated: |
404 |
\begin{eqnarray}\label{LINEARPDE.Peclet.2} |
405 |
\Xi^{B}=\frac{\xi(P^{B}) h}{\|B\hackscore{:}\|} |
406 |
\mbox{ and } |
407 |
\Xi^{C}=\frac{\xi(P^{C}) h}{\|C\hackscore{:}\|} \\ |
408 |
\mbox{ or } |
409 |
\Xi^{B}\hackscore{ik}=\frac{\xi(P^{B}\hackscore{ik}) h}{\|B\hackscore{i:k}\|} |
410 |
\mbox{ and } |
411 |
\Xi^{C}\hackscore{ik}=\frac{\xi(P^{C}\hackscore{ik}) h}{\|C\hackscore{ik:}\|} |
412 |
\end{eqnarray} |
413 |
where $\xi$ is a suitable function of the Peclet number. |
414 |
In the case of a single PDE the coefficient are up-dated in the following way: |
415 |
\begin{eqnarray}\label{LINEARPDE.ADVECTIVE.1} |
416 |
A\hackscore{jl} \leftarrow A\hackscore{jl} + \Xi^{B} B\hackscore{j} B\hackscore{l} + \Xi^{C} C\hackscore{j} C\hackscore{l} \\ |
417 |
B\hackscore{j} \leftarrow B\hackscore{j} + \Xi^{C} C\hackscore{j} D \\ |
418 |
C\hackscore{j} \leftarrow C\hackscore{j} + \Xi^{B} B\hackscore{j} D \\ |
419 |
X\hackscore{j} \leftarrow X\hackscore{j} + (Xi^{B} B\hackscore{j} + \Xi^{C} C\hackscore{j}) Y \\ |
420 |
\end{eqnarray} |
421 |
Similar for the case of a systems of PDEs: |
422 |
\begin{eqnarray}\label{LINEARPDE.ADVECTIVE.SYSTEM} |
423 |
A\hackscore{ijkl} \leftarrow A\hackscore{ijl} + \Xi^{B}\hackscore{ik} B\hackscore{ijk} B\hackscore{ilk} + |
424 |
\Xi^{C}\hackscore{ik} C\hackscore{ikj} C\hackscore{ikl} \\ |
425 |
B\hackscore{ijk} \leftarrow B\hackscore{ijk} + \Xi^{C}\hackscore{ij} C\hackscore{ikj} D\hackscore{ik} \\ |
426 |
C\hackscore{ikl} \leftarrow C\hackscore{ikl} + \Xi^{B}\hackscore{ik} B\hackscore{ilk} D\hackscore{ik} \\ |
427 |
X\hackscore{ij} \leftarrow X\hackscore{ij} + (Xi^{B}\hackscore{ik} B\hackscore{ij} + \Xi^{C}\hackscore{ik} C\hackscore{ij}) Y\hackscore{i} \\ |
428 |
\end{eqnarray} |
429 |
Using upwinding in this form, introduces an additonal error which is proprtional to the cell size $h$ |
430 |
but with the intension to stabilize the solution. |
431 |
|
432 |
\begin{classdesc}{AdvectivePDE}{domain,numEquations=0,numSolutions=0,xi=AdvectivePDE.ELMAN_RAMAGE} |
433 |
opens a linear, steady, second order PDE on the \Domain \var{domain}. \var{numEquations} |
434 |
and \var{numSolutions} gives the number of equations and the number of solutiopn components. |
435 |
If \var{numEquations} and \var{numSolutions} is non-positive, the number of equations |
436 |
and the number solutions, respectively, stay undefined until a coefficient is |
437 |
defined. \var{xi} defines a function which returns for any given Preclet number $P\ge 0$ the |
438 |
$\xi$-value used to define the stabilization parameters $\Xi^{B}$ and $\Xi^{C}$. |
439 |
\AdvectivePDE is derived from \LinearPDE. |
440 |
\end{classdesc} |
441 |
|
442 |
\begin{memberdesc}[AdvectivePDE]{SIMPLIFIED_BROOKS_HUGHES}{} |
443 |
Predefined function to set a values for $\xi$ from a Preclet number $P$. |
444 |
This function uses the method suggested in \Ref{SUPG1} |
445 |
where $\xi(P)$ is given by |
446 |
\begin{equation}\label{LINEARPDE.ADVECTIVE.1d} |
447 |
\xi(P)=coth(P)-\frac{1}{P} \;. |
448 |
\end{equation} |
449 |
As the evaluation of $coth$ is expensive we are using the approximation: |
450 |
The function $\xi$ is approximated by |
451 |
\begin{equation}\label{LINEARPDE.ADVECTIVE.23} |
452 |
\xi(P)= |
453 |
\left\{ |
454 |
\begin{array}{lc} |
455 |
\frac{P}{6} & P<3 \\ |
456 |
\frac{1}{2} & \mbox{otherwise} |
457 |
\end{array} |
458 |
\right. |
459 |
\end{equation} |
460 |
\end{memberdesc} |
461 |
|
462 |
\begin{memberdesc}[AdvectivePDE]{ELMAN_RAMAGE}{} |
463 |
Predefined function to set a values for $\xi$ from a Preclet number $P$. |
464 |
This function uses the method suggested in \Ref{SUPG2} |
465 |
where $\xi(P)$ is given by |
466 |
\begin{equation}\label{LINEARPDE.ADVECTIVE.23b} |
467 |
\xi(P)= |
468 |
\left\{ |
469 |
\begin{array}{lc} |
470 |
0 & P<1 \\ |
471 |
\frac{1}{2}(1-\frac{1}{P}) & \mbox{otherwise} |
472 |
\end{array} |
473 |
\right. |
474 |
\end{equation} |
475 |
\end{memberdesc} |
476 |
|
477 |
|
478 |
\section{The \Poisson Class} |
479 |
|
480 |
The \Poisson class provides an easy way to define and solve the Poisson |
481 |
equation |
482 |
\begin{equation}\label{POISSON.1} |
483 |
-u\hackscore{,ii}=f\; . |
484 |
\end{equation} |
485 |
with homogeneous boundary conditions |
486 |
\begin{equation}\label{POISSON.2} |
487 |
n\hackscore{i}u\hackscore{,i}=0 |
488 |
\end{equation} |
489 |
and homogeneous constraints |
490 |
\begin{equation}\label{POISSON.3} |
491 |
u=0 \mbox{ where } q>0 |
492 |
\end{equation} |
493 |
$f$ has to be a \Scalar in the \Function and $q$ must be |
494 |
a \Scalar in the \SolutionFS. |
495 |
|
496 |
\begin{classdesc}{Poisson}{domain} |
497 |
opens a Poisson equation on the \Domain domain. \Poisson is derived from \LinearPDE. |
498 |
\end{classdesc} |
499 |
\begin{methoddesc}[Poisson]{setValue}{f=escript.Data(),q=escript.Data()} |
500 |
assigns new values to \var{f} and \var{q}. |
501 |
\end{methoddesc} |