1 |
ksteube |
1811 |
|
2 |
|
|
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
3 |
ksteube |
1316 |
% |
4 |
ksteube |
1811 |
% Copyright (c) 2003-2008 by University of Queensland |
5 |
|
|
% Earth Systems Science Computational Center (ESSCC) |
6 |
|
|
% http://www.uq.edu.au/esscc |
7 |
gross |
625 |
% |
8 |
ksteube |
1811 |
% 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 |
gross |
625 |
% |
12 |
ksteube |
1811 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
13 |
jgs |
102 |
|
14 |
ksteube |
1811 |
|
15 |
ksteube |
1318 |
\chapter{The Module \linearPDEs} |
16 |
jgs |
102 |
|
17 |
|
|
|
18 |
gross |
999 |
|
19 |
|
|
\section{Linear Partial Differential Equations} |
20 |
jgs |
102 |
\label{SEC LinearPDE} |
21 |
|
|
|
22 |
|
|
The \LinearPDE class is used to define a general linear, steady, second order PDE |
23 |
|
|
for an unknown function $u$ on a given $\Omega$ defined through a \Domain object. |
24 |
|
|
In the following $\Gamma$ denotes the boundary of the domain $\Omega$. $n$ denotes |
25 |
gross |
660 |
the outer normal field on $\Gamma$. |
26 |
jgs |
102 |
|
27 |
gross |
660 |
For a single PDE with a solution with a single component the linear PDE is defined in the |
28 |
jgs |
102 |
following form: |
29 |
|
|
\begin{equation}\label{LINEARPDE.SINGLE.1} |
30 |
gross |
660 |
-(A\hackscore{jl} u\hackscore{,l})\hackscore{,j}-(B\hackscore{j} u)\hackscore{,j}+C\hackscore{l} u\hackscore{,l}+D u =-X\hackscore{j,j}+Y \; . |
31 |
jgs |
102 |
\end{equation} |
32 |
gross |
660 |
$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. |
33 |
|
|
The coefficients $A$, $B$, $C$, $D$, $X$ and $Y$ have to be specified through \Data objects in the |
34 |
|
|
\Function on the PDE or objects that can be converted into such \Data objects. |
35 |
|
|
$A$ is a \RankTwo, $B$, $C$ and $X$ are \RankOne and $D$ and $Y$ are scalar. |
36 |
jgs |
102 |
The following natural |
37 |
|
|
boundary conditions are considered \index{boundary condition!natural} on $\Gamma$: |
38 |
|
|
\begin{equation}\label{LINEARPDE.SINGLE.2} |
39 |
|
|
n\hackscore{j}(A\hackscore{jl} u\hackscore{,l}+B\hackscore{j} u)+d u=n\hackscore{j}X\hackscore{j} + y \;. |
40 |
|
|
\end{equation} |
41 |
gross |
660 |
Notice that the coefficients $A$, $B$ and $X$ are defined in the PDE. The coefficients $d$ and $y$ are |
42 |
|
|
each a \Scalar in the \FunctionOnBoundary. Constraints \index{constraint} for the solution prescribing the value of the |
43 |
jgs |
102 |
solution at certain locations in the domain. They have the form |
44 |
|
|
\begin{equation}\label{LINEARPDE.SINGLE.3} |
45 |
|
|
u=r \mbox{ where } q>0 |
46 |
|
|
\end{equation} |
47 |
|
|
$r$ and $q$ are each \Scalar where $q$ is the characteristic function |
48 |
|
|
\index{characteristic function} defining where the constraint is applied. |
49 |
|
|
The constraints defined by \eqn{LINEARPDE.SINGLE.3} override any other condition set by \eqn{LINEARPDE.SINGLE.1} |
50 |
gross |
660 |
or \eqn{LINEARPDE.SINGLE.2}. |
51 |
gross |
625 |
|
52 |
jgs |
102 |
For a system of PDEs and a solution with several components the PDE has the form |
53 |
|
|
\begin{equation}\label{LINEARPDE.SYSTEM.1} |
54 |
gross |
660 |
-(A\hackscore{ijkl} u\hackscore{k,l})\hackscore{,j}-(B\hackscore{ijk} u\hackscore{k})\hackscore{,j}+C\hackscore{ikl} u\hackscore{k,l}+D\hackscore{ik} u\hackscore{k} =-X\hackscore{ij,j}+Y\hackscore{i} \; . |
55 |
jgs |
102 |
\end{equation} |
56 |
gross |
660 |
$A$ is a \RankFour, $B$ and $C$ are each a \RankThree, $D$ and $X$ are each a \RankTwo and $Y$ is a \RankOne. |
57 |
jgs |
102 |
The natural boundary conditions \index{boundary condition!natural} take the form: |
58 |
|
|
\begin{equation}\label{LINEARPDE.SYSTEM.2} |
59 |
gross |
625 |
n\hackscore{j}(A\hackscore{ijkl} u\hackscore{k,l}+B\hackscore{ijk} u\hackscore{k})+d\hackscore{ik} u\hackscore{k}=n\hackscore{j}X\hackscore{ij}+y\hackscore{i} \;. |
60 |
jgs |
102 |
\end{equation} |
61 |
gross |
660 |
The coefficient $d$ is a \RankTwo and $y$ is a |
62 |
jgs |
102 |
\RankOne both in the \FunctionOnBoundary. Constraints \index{constraint} take the form |
63 |
|
|
\begin{equation}\label{LINEARPDE.SYSTEM.3} |
64 |
|
|
u\hackscore{i}=r\hackscore{i} \mbox{ where } q\hackscore{i}>0 |
65 |
|
|
\end{equation} |
66 |
gross |
660 |
$r$ and $q$ are each \RankOne. Notice that not necessarily all components must |
67 |
gross |
625 |
have a constraint at all locations. |
68 |
|
|
|
69 |
jgs |
102 |
\LinearPDE also supports solution discontinuities \index{discontinuity} over contact region $\Gamma^{contact}$ |
70 |
|
|
in the domain $\Omega$. To specify the conditions across the discontinuity we are using the |
71 |
|
|
generalised flux $J$ which is in the case of a systems of PDEs and several components of the solution |
72 |
gross |
660 |
defined as |
73 |
jgs |
102 |
\begin{equation}\label{LINEARPDE.SYSTEM.5} |
74 |
|
|
J\hackscore{ij}=A\hackscore{ijkl}u\hackscore{k,l}+B\hackscore{ijk}u\hackscore{k}-X\hackscore{ij} |
75 |
|
|
\end{equation} |
76 |
|
|
For the case of single solution component and single PDE $J$ is defined |
77 |
|
|
\begin{equation}\label{LINEARPDE.SINGLE.5} |
78 |
|
|
J\hackscore{j}=A\hackscore{jl}u\hackscore{,l}+B\hackscore{j}u\hackscore{k}-X\hackscore{j} |
79 |
|
|
\end{equation} |
80 |
gross |
660 |
In the context of discontinuities \index{discontinuity} $n$ denotes the normal on the |
81 |
jgs |
102 |
discontinuity pointing from side 0 towards side 1. For a system of PDEs |
82 |
|
|
the contact condition takes the form |
83 |
|
|
\begin{equation}\label{LINEARPDE.SYSTEM.6} |
84 |
|
|
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} \; . |
85 |
|
|
\end{equation} |
86 |
|
|
where $J^{0}$ and $J^{1}$ are the fluxes on side $0$ and side $1$ of the |
87 |
|
|
discontinuity $\Gamma^{contact}$, respectively. $[u]$, which is the difference |
88 |
|
|
of the solution at side 1 and at side 0, denotes the jump of $u$ across $\Gamma^{contact}$. |
89 |
gross |
660 |
The coefficient $d^{contact}$ is a \RankTwo and $y^{contact}$ is a |
90 |
jgs |
102 |
\RankOne both in the \FunctionOnContactZero or \FunctionOnContactOne. |
91 |
|
|
In case of a single PDE and a single component solution the contact condition takes the form |
92 |
|
|
\begin{equation}\label{LINEARPDE.SINGLE.6} |
93 |
|
|
n\hackscore{j} J^{0}\hackscore{j}=n\hackscore{j} J^{1}\hackscore{j}=y^{contact} - d^{contact}[u] |
94 |
|
|
\end{equation} |
95 |
ksteube |
1316 |
In this case the the coefficient $d^{contact}$ and $y^{contact}$ are each \Scalar |
96 |
jgs |
102 |
both in the \FunctionOnContactZero or \FunctionOnContactOne. |
97 |
gross |
625 |
|
98 |
|
|
The PDE is symmetrical \index{symmetrical} if |
99 |
|
|
\begin{equation}\label{LINEARPDE.SINGLE.4} |
100 |
|
|
A\hackscore{jl}=A\hackscore{lj} \mbox{ and } B\hackscore{j}=C\hackscore{j} |
101 |
|
|
\end{equation} |
102 |
|
|
The system of PDEs is symmetrical \index{symmetrical} if |
103 |
|
|
\begin{eqnarray} |
104 |
|
|
\label{LINEARPDE.SYSTEM.4} |
105 |
jfenwick |
1959 |
A\hackscore{ijkl}&=&A\hackscore{klij} \\ |
106 |
|
|
B\hackscore{ijk}&=&C\hackscore{kij} \\ |
107 |
|
|
D\hackscore{ik}&=&D\hackscore{ki} \\ |
108 |
|
|
d\hackscore{ik}&=&d\hackscore{ki} \\ |
109 |
|
|
d^{contact}\hackscore{ik}&=&d^{contact}\hackscore{ki} |
110 |
gross |
625 |
\end{eqnarray} |
111 |
jfenwick |
1959 |
Note that in contrast with the scalar case~\eqn{LINEARPDE.SINGLE.4} now the coefficients $D$, $d$ abd $d^{contact}$ |
112 |
gross |
625 |
have to be inspected. |
113 |
|
|
|
114 |
gross |
999 |
|
115 |
|
|
\subsection{Classes} |
116 |
|
|
\declaremodule{extension}{esys.escript.linearPDEs} |
117 |
ksteube |
1316 |
\modulesynopsis{Linear partial differential equation handler} |
118 |
gross |
999 |
The module \linearPDEs provides an interface to define and solve linear partial |
119 |
jfenwick |
1959 |
differential equations within \escript. The module \linearPDEs does not provide any |
120 |
gross |
999 |
solver capabilities in itself but hands the PDE over to |
121 |
|
|
the PDE solver library defined through the \Domain of the PDE. |
122 |
|
|
The general interface is provided through the \LinearPDE class. The |
123 |
|
|
\AdvectivePDE which is derived from the \LinearPDE class |
124 |
jfenwick |
1959 |
provides an interface to a PDE dominated by its advective terms. The \Poisson |
125 |
gross |
999 |
class which is also derived form the \LinearPDE class should be used |
126 |
|
|
to define the Poisson equation \index{Poisson}. |
127 |
|
|
|
128 |
|
|
\subsection{\LinearPDE class} |
129 |
ksteube |
1316 |
This is the general class to define a linear PDE in \escript. We list a selection of the most |
130 |
jfenwick |
1959 |
important methods of the class. For a complete list, see the reference at \ReferenceGuide. |
131 |
gross |
625 |
|
132 |
jgs |
102 |
\begin{classdesc}{LinearPDE}{domain,numEquations=0,numSolutions=0} |
133 |
|
|
opens a linear, steady, second order PDE on the \Domain \var{domain}. \var{numEquations} |
134 |
ksteube |
1316 |
and \var{numSolutions} gives the number of equations and the number of solution components. |
135 |
gross |
660 |
If \var{numEquations} and \var{numSolutions} is non-positive, the number of equations |
136 |
ksteube |
1316 |
and the number solutions, respectively, stay undefined until a coefficient is |
137 |
gross |
660 |
defined. |
138 |
jgs |
102 |
\end{classdesc} |
139 |
|
|
|
140 |
jfenwick |
1959 |
\subsubsection{\LinearPDE methods} |
141 |
|
|
|
142 |
gross |
625 |
\begin{methoddesc}[LinearPDE]{setValue}{ |
143 |
gross |
660 |
\optional{A}\optional{, B}, |
144 |
|
|
\optional{, C}\optional{, D} |
145 |
|
|
\optional{, X}\optional{, Y} |
146 |
|
|
\optional{, d}\optional{, y} |
147 |
|
|
\optional{, d_contact}\optional{, y_contact} |
148 |
|
|
\optional{, q}\optional{, r}} |
149 |
ksteube |
1316 |
assigns new values to coefficients. By default all values are assumed to be zero\footnote{ |
150 |
gross |
660 |
In fact it is assumed they are not present by assigning the value \code{escript.Data()}. The |
151 |
|
|
can by used by the solver library to reduce computational costs. |
152 |
|
|
} |
153 |
gross |
625 |
If the new coefficient value is not a \Data object, it is converted into a \Data object in the |
154 |
jgs |
102 |
appropriate \FunctionSpace. |
155 |
|
|
\end{methoddesc} |
156 |
|
|
|
157 |
|
|
\begin{methoddesc}[LinearPDE]{getCoefficient}{name} |
158 |
gross |
660 |
return the value assigned to coefficient \var{name}. If \var{name} is not a valid name |
159 |
|
|
an exception is raised. |
160 |
jgs |
102 |
\end{methoddesc} |
161 |
|
|
|
162 |
|
|
\begin{methoddesc}[LinearPDE]{getShapeOfCoefficient}{name} |
163 |
|
|
returns the shape of coefficient \var{name} even if no value has been assigned to it. |
164 |
|
|
\end{methoddesc} |
165 |
|
|
|
166 |
gross |
625 |
\begin{methoddesc}[LinearPDE]{getFunctionSpaceForCoefficient}{name} |
167 |
jgs |
102 |
returns the \FunctionSpace of coefficient \var{name} even if no value has been assigned to it. |
168 |
|
|
\end{methoddesc} |
169 |
|
|
|
170 |
|
|
\begin{methoddesc}[LinearPDE]{setDebugOn}{} |
171 |
jfenwick |
1959 |
switches on debug mode. |
172 |
jgs |
102 |
\end{methoddesc} |
173 |
|
|
|
174 |
|
|
\begin{methoddesc}[LinearPDE]{setDebugOff}{} |
175 |
jfenwick |
1959 |
switches off debug mode. |
176 |
jgs |
102 |
\end{methoddesc} |
177 |
|
|
|
178 |
gross |
625 |
\begin{methoddesc}[LinearPDE]{isUsingLumping}{} |
179 |
ksteube |
1316 |
returns \True if \LUMPING is set as the solver for the system of linear equations. |
180 |
gross |
653 |
Otherwise \False is returned. |
181 |
jgs |
102 |
\end{methoddesc} |
182 |
|
|
|
183 |
gross |
653 |
\begin{methoddesc}[LinearPDE]{setSolverMethod}{\optional{solver=LinearPDE.DEFAULT}\optional{, preconditioner=LinearPDE.DEFAULT}} |
184 |
jfenwick |
1959 |
sets the solver method and preconditioner to be used. It should be noted that a PDE solver library |
185 |
gross |
660 |
may not know the specified solver method but may choose a similar method and preconditioner. |
186 |
jgs |
102 |
\end{methoddesc} |
187 |
|
|
|
188 |
gross |
653 |
\begin{methoddesc}[LinearPDE]{getSolverMethodName}{} |
189 |
jfenwick |
1959 |
returns the name of the solver method and preconditioner which is in use. |
190 |
gross |
653 |
\end{methoddesc} |
191 |
|
|
|
192 |
|
|
\begin{methoddesc}[LinearPDE]{getSolverMethod}{} |
193 |
jfenwick |
1959 |
returns the solver method and preconditioner which is in use. |
194 |
gross |
653 |
\end{methoddesc} |
195 |
|
|
|
196 |
|
|
\begin{methoddesc}[LinearPDE]{setSolverPackage}{\optional{package=LinearPDE.DEFAULT}} |
197 |
jfenwick |
1959 |
sets the solver package to be used by PDE library to solve the linear systems of equations. The |
198 |
ksteube |
1316 |
specified package may not be supported by the PDE solver library. In this case, depending on |
199 |
|
|
the PDE solver, the default solver is used or an exception is thrown. |
200 |
gross |
660 |
If \var{package} is not specified, the default package of the PDE solver library is used. |
201 |
gross |
653 |
\end{methoddesc} |
202 |
|
|
|
203 |
|
|
\begin{methoddesc}[LinearPDE]{getSolverPackage}{} |
204 |
|
|
returns the linear solver package currently by the PDE solver library |
205 |
|
|
\end{methoddesc} |
206 |
|
|
|
207 |
|
|
|
208 |
jfenwick |
1959 |
\begin{methoddesc}[LinearPDE]{setTolerance}{\optional{tol=1.e-8}} |
209 |
|
|
resets the tolerance for solution. The actually meaning of tolerance depends |
210 |
|
|
on the underlying PDE library. In most cases, the tolerance |
211 |
ksteube |
1316 |
will only consider the error from solving the discrete problem but will |
212 |
gross |
625 |
not consider any discretization error. |
213 |
|
|
\end{methoddesc} |
214 |
jgs |
102 |
|
215 |
lgraham |
1700 |
\begin{methoddesc}[LinearPDE]{setToleranceReductionFactor}{TOL} |
216 |
jfenwick |
1959 |
lowers the tolerance by a factor of TOL. |
217 |
lgraham |
1700 |
\end{methoddesc} |
218 |
|
|
|
219 |
gross |
625 |
\begin{methoddesc}[LinearPDE]{getTolerance}{} |
220 |
|
|
returns the current tolerance of the solution |
221 |
jgs |
102 |
\end{methoddesc} |
222 |
|
|
|
223 |
gross |
625 |
\begin{methoddesc}[LinearPDE]{getDomain}{} |
224 |
|
|
returns the \Domain of the PDE. |
225 |
jgs |
102 |
\end{methoddesc} |
226 |
|
|
|
227 |
gross |
625 |
\begin{methoddesc}[LinearPDE]{getDim}{} |
228 |
|
|
returns the spatial dimension of the PDE. |
229 |
jgs |
102 |
\end{methoddesc} |
230 |
|
|
|
231 |
gross |
625 |
\begin{methoddesc}[LinearPDE]{getNumEquations}{} |
232 |
|
|
returns the number of equations. |
233 |
|
|
\end{methoddesc} |
234 |
jgs |
102 |
|
235 |
gross |
625 |
\begin{methoddesc}[LinearPDE]{getNumSolutions}{} |
236 |
|
|
returns the number of components of the solution. |
237 |
jgs |
102 |
\end{methoddesc} |
238 |
|
|
|
239 |
gross |
625 |
\begin{methoddesc}[LinearPDE]{checkSymmetry}{verbose=\False} |
240 |
gross |
660 |
returns \True if the PDE is symmetric and \False otherwise. |
241 |
jfenwick |
1959 |
The method is very computationally expensive and should only be |
242 |
gross |
625 |
called for testing purposes. The symmetry flag is not altered. |
243 |
|
|
If \var{verbose}=\True information about where symmetry is violated |
244 |
|
|
are printed. |
245 |
jgs |
102 |
\end{methoddesc} |
246 |
|
|
|
247 |
gross |
625 |
\begin{methoddesc}[LinearPDE]{getFlux}{u} |
248 |
|
|
returns the flux $J\hackscore{ij}$ \index{flux} for given solution \var{u} |
249 |
|
|
defined by \eqn{LINEARPDE.SYSTEM.5} and \eqn{LINEARPDE.SINGLE.5}, respectively. |
250 |
jgs |
102 |
\end{methoddesc} |
251 |
|
|
|
252 |
gross |
625 |
|
253 |
jgs |
102 |
\begin{methoddesc}[LinearPDE]{isSymmetric}{} |
254 |
|
|
returns \True if the PDE has been indicated to be symmetric. |
255 |
|
|
Otherwise \False is returned. |
256 |
|
|
\end{methoddesc} |
257 |
|
|
|
258 |
|
|
\begin{methoddesc}[LinearPDE]{setSymmetryOn}{} |
259 |
|
|
indicates that the PDE is symmetric. |
260 |
|
|
\end{methoddesc} |
261 |
|
|
|
262 |
|
|
\begin{methoddesc}[LinearPDE]{setSymmetryOff}{} |
263 |
|
|
indicates that the PDE is not symmetric. |
264 |
|
|
\end{methoddesc} |
265 |
|
|
|
266 |
|
|
\begin{methoddesc}[LinearPDE]{setReducedOrderOn}{} |
267 |
gross |
660 |
switches on the reduction of polynomial order for the solution and equation evaluation even if |
268 |
|
|
a quadratic or higher interpolation order is defined in the \Domain. This feature may not |
269 |
gross |
625 |
be supported by all PDE libraries. |
270 |
jgs |
102 |
\end{methoddesc} |
271 |
|
|
|
272 |
|
|
\begin{methoddesc}[LinearPDE]{setReducedOrderOff}{} |
273 |
gross |
660 |
switches off the reduction of polynomial order for the solution and |
274 |
jgs |
102 |
equation evaluation. |
275 |
|
|
\end{methoddesc} |
276 |
|
|
|
277 |
|
|
\begin{methoddesc}[LinearPDE]{getOperator}{} |
278 |
|
|
returns the \Operator of the PDE. |
279 |
|
|
\end{methoddesc} |
280 |
|
|
|
281 |
gross |
625 |
\begin{methoddesc}[LinearPDE]{getRightHandSide}{} |
282 |
jgs |
102 |
returns the right hand side of the PDE as a \Data object. If |
283 |
jfenwick |
1959 |
\var{ignoreConstraint}=\True, then the constraints are not considered |
284 |
jgs |
102 |
when building up the right hand side. |
285 |
|
|
\end{methoddesc} |
286 |
|
|
|
287 |
|
|
\begin{methoddesc}[LinearPDE]{getSystem}{} |
288 |
|
|
returns the \Operator and right hand side of the PDE. |
289 |
|
|
\end{methoddesc} |
290 |
|
|
|
291 |
gross |
625 |
\begin{methoddesc}[LinearPDE]{getSolution}{ |
292 |
|
|
\optional{verbose=False} |
293 |
|
|
\optional{, reordering=LinearPDE.NO_REORDERING} |
294 |
|
|
\optional{, iter_max=1000} |
295 |
|
|
\optional{, drop_tolerance=0.01} |
296 |
|
|
\optional{, drop_storage=1.20} |
297 |
|
|
\optional{, truncation=-1} |
298 |
|
|
\optional{, restart=-1} |
299 |
|
|
} |
300 |
jfenwick |
1959 |
returns (an approximation of) the solution of the PDE. If \code{verbose=\True}, then some information is printed during the solution process. |
301 |
gross |
660 |
\var{reordering} selects a reordering methods that is applied before or during the solution process |
302 |
jfenwick |
1959 |
(=\NOREORDERING, \MINIMUMFILLIN, \NESTEDDESCTION). |
303 |
gross |
660 |
\var{iter_max} specifies the maximum number of iteration steps that are allowed to reach the specified tolerance. |
304 |
gross |
625 |
\var{drop_tolerance} specifies a relative tolerance for small elements to be dropped when building a preconditioner |
305 |
gross |
653 |
(eg. in \ILUT). \var{drop_storage} limits the extra storage allowed when building a preconditioner |
306 |
gross |
660 |
(eg. in \ILUT). The extra storage is given relative to the size of the stiffness matrix, eg. |
307 |
|
|
\var{drop_storage=1.2} will allow the preconditioner to use the $1.2$ fold storage space than used |
308 |
|
|
for the stiffness matrix. \var{truncation} defines the truncation. |
309 |
jgs |
102 |
\end{methoddesc} |
310 |
|
|
|
311 |
jfenwick |
1959 |
\subsubsection{\LinearPDE symbols/members} |
312 |
|
|
|
313 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{DEFAULT} |
314 |
gross |
660 |
default method, preconditioner or package to be used to solve the PDE. An appropriate method should be |
315 |
gross |
625 |
chosen by the used PDE solver library. |
316 |
|
|
\end{memberdesc} |
317 |
jgs |
102 |
|
318 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{SCSL} |
319 |
gross |
660 |
the SCSL library by SGI,~\Ref{SCSL}\footnote{The SCSL library will only be available on SGI systems} |
320 |
gross |
625 |
\end{memberdesc} |
321 |
jgs |
102 |
|
322 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{MKL} |
323 |
lgraham |
1700 |
the MKL library by Intel,~\Ref{MKL}\footnote{The MKL library will only be available when the Intel compilation environment is used.}. |
324 |
gross |
625 |
\end{memberdesc} |
325 |
jgs |
102 |
|
326 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{UMFPACK} |
327 |
gross |
653 |
the UMFPACK,~\Ref{UMFPACK}. Remark: UMFPACK is not parallelized. |
328 |
gross |
625 |
\end{memberdesc} |
329 |
jgs |
102 |
|
330 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{PASO} |
331 |
gross |
653 |
the solver library of \finley, see \Sec{CHAPTER ON FINLEY}. |
332 |
gross |
625 |
\end{memberdesc} |
333 |
jgs |
102 |
|
334 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{ITERATIVE} |
335 |
gross |
653 |
the default iterative method and preconditioner. The actually used method depends on the |
336 |
ksteube |
1316 |
PDE solver library and the solver package been chosen. Typically, \PCG is used for symmetric PDEs |
337 |
gross |
660 |
and \BiCGStab otherwise, both with \JACOBI preconditioner. |
338 |
gross |
625 |
\end{memberdesc} |
339 |
jgs |
102 |
|
340 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{DIRECT} |
341 |
gross |
660 |
the default direct linear solver. |
342 |
gross |
625 |
\end{memberdesc} |
343 |
jgs |
102 |
|
344 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{CHOLEVSKY} |
345 |
gross |
660 |
direct solver based on Cholevsky factorization (or similar), see~\Ref{Saad}. The solver will require a symmetric PDE. |
346 |
gross |
625 |
\end{memberdesc} |
347 |
jgs |
110 |
|
348 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{PCG} |
349 |
gross |
653 |
preconditioned conjugate gradient method, see~\Ref{WEISS}\index{linear solver!PCG}\index{PCG}. The solver will require a symmetric PDE. |
350 |
gross |
625 |
\end{memberdesc} |
351 |
jgs |
110 |
|
352 |
artak |
1978 |
\begin{memberdesc}[LinearPDE]{TFQMR} |
353 |
|
|
transpose-free quasi-minimal residual method, see~\Ref{WEISS}\index{linear solver!TFQMR}\index{TFQMR}. \end{memberdesc} |
354 |
|
|
|
355 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{GMRES} |
356 |
gross |
653 |
the GMRES method, see~\Ref{WEISS}\index{linear solver!GMRES}\index{GMRES}. Truncation and restart are controlled by the parameters |
357 |
gross |
625 |
\var{truncation} and \var{restart} of \method{getSolution}. |
358 |
|
|
\end{memberdesc} |
359 |
jgs |
102 |
|
360 |
artak |
1978 |
\begin{memberdesc}[LinearPDE]{MINRES} |
361 |
|
|
minimal residual method method, \index{linear solver!MINRES}\index{MINRES} \end{memberdesc} |
362 |
|
|
|
363 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{LUMPING} |
364 |
gross |
660 |
uses lumping to solve the system of linear equations~\index{linear solver!lumping}\index{lumping}. This solver technique |
365 |
|
|
condenses the stiffness matrix to a diagonal matrix so the solution of the linear systems becomes very cheap. It can be used when |
366 |
gross |
653 |
only \var{D} is present but in any case has to applied with care. The difference in the solutions with and without lumping can be significant |
367 |
jfenwick |
1959 |
but is expected to converge to zero when the mesh gets finer. |
368 |
gross |
660 |
Lumping does not use the linear system solver library. |
369 |
gross |
625 |
\end{memberdesc} |
370 |
jgs |
107 |
|
371 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{PRES20} |
372 |
gross |
653 |
the GMRES method with truncation after five residuals and |
373 |
gross |
625 |
restart after 20 steps, see~\Ref{WEISS}. |
374 |
gross |
999 |
\end{memberdesc} |
375 |
gross |
625 |
|
376 |
|
|
\begin{memberdesc}[LinearPDE]{CGS} |
377 |
|
|
conjugate gradient squared method, see~\Ref{WEISS}. |
378 |
jgs |
107 |
\end{memberdesc} |
379 |
|
|
|
380 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{BICGSTAB} |
381 |
gross |
660 |
stabilized bi-conjugate gradients methods, see~\Ref{WEISS}. |
382 |
jgs |
107 |
\end{memberdesc} |
383 |
|
|
|
384 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{SSOR} |
385 |
gross |
653 |
symmetric successive over-relaxation method, see~\Ref{WEISS}. Typically used as preconditioner but some linear solver libraries support |
386 |
gross |
660 |
this as a solver. |
387 |
gross |
625 |
\end{memberdesc} |
388 |
|
|
\begin{memberdesc}[LinearPDE]{ILU0} |
389 |
gross |
660 |
the incomplete LU factorization preconditioner with no fill-in, see~\Ref{Saad}. |
390 |
gross |
653 |
\end{memberdesc} |
391 |
|
|
|
392 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{ILUT} |
393 |
gross |
653 |
the incomplete LU factorization preconditioner with fill-in, see~\Ref{Saad}. During the LU-factorization element with |
394 |
|
|
relative size less then \var{drop_tolerance} are dropped. Moreover, the size of the LU-factorization is restricted to the |
395 |
gross |
660 |
\var{drop_storage}-fold of the stiffness matrix. \var{drop_tolerance} and \var{drop_storage} are both set in the |
396 |
|
|
\method{getSolution} call. |
397 |
gross |
653 |
\end{memberdesc} |
398 |
|
|
|
399 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{JACOBI} |
400 |
gross |
653 |
the Jacobi preconditioner, see~\Ref{Saad}. |
401 |
|
|
\end{memberdesc} |
402 |
|
|
|
403 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{AMG} |
404 |
gross |
660 |
the algebraic--multi grid method, see~\Ref{AMG}. This method can be used as linear solver method but is more robust when used |
405 |
gross |
653 |
in a preconditioner. |
406 |
|
|
\end{memberdesc} |
407 |
|
|
|
408 |
artak |
1978 |
\begin{memberdesc}[LinearPDE]{GS} |
409 |
|
|
the symmetric Gauss-Seidel preconditioner, see~\Ref{Saad}. |
410 |
|
|
\end{memberdesc} |
411 |
|
|
|
412 |
gross |
625 |
\begin{memberdesc}[LinearPDE]{RILU} |
413 |
gross |
653 |
recursive incomplete LU factorization preconditioner, see~\Ref{RILU}. This method is similar to \ILUT but uses smoothing |
414 |
|
|
between levels. During the LU-factorization element with |
415 |
|
|
relative size less then \var{drop_tolerance} are dropped. Moreover, the size of the LU-factorization is restricted to the |
416 |
gross |
660 |
\var{drop_storage}-fold of the stiffness matrix. \var{drop_tolerance} and \var{drop_storage} are both set in the |
417 |
|
|
\method{getSolution} call. |
418 |
gross |
653 |
\end{memberdesc} |
419 |
jgs |
107 |
|
420 |
gross |
653 |
\begin{memberdesc}[LinearPDE]{NO_REORDERING} |
421 |
|
|
no ordering is used during factorization. |
422 |
|
|
\end{memberdesc} |
423 |
gross |
625 |
|
424 |
gross |
653 |
\begin{memberdesc}[LinearPDE]{MINIMUM_FILL_IN} |
425 |
|
|
applies reordering before factorization using a fill-in minimization strategy. You have to check with the particular solver library or |
426 |
|
|
linear solver package if this is supported. In any case, it is advisable to apply reordering on the mesh to minimize fill-in. |
427 |
|
|
\end{memberdesc} |
428 |
gross |
625 |
|
429 |
|
|
\begin{memberdesc}[LinearPDE]{NESTED_DISSECTION} |
430 |
gross |
653 |
applies reordering before factorization using a nested dissection strategy. You have to check with the particular solver library or |
431 |
|
|
linear solver package if this is supported. In any case, it is advisable to apply reordering on the mesh to minimize fill-in. |
432 |
|
|
\end{memberdesc} |
433 |
gross |
625 |
|
434 |
gross |
999 |
\subsection{The \Poisson Class} |
435 |
jgs |
102 |
The \Poisson class provides an easy way to define and solve the Poisson |
436 |
|
|
equation |
437 |
|
|
\begin{equation}\label{POISSON.1} |
438 |
|
|
-u\hackscore{,ii}=f\; . |
439 |
|
|
\end{equation} |
440 |
|
|
with homogeneous boundary conditions |
441 |
|
|
\begin{equation}\label{POISSON.2} |
442 |
|
|
n\hackscore{i}u\hackscore{,i}=0 |
443 |
|
|
\end{equation} |
444 |
|
|
and homogeneous constraints |
445 |
|
|
\begin{equation}\label{POISSON.3} |
446 |
|
|
u=0 \mbox{ where } q>0 |
447 |
|
|
\end{equation} |
448 |
|
|
$f$ has to be a \Scalar in the \Function and $q$ must be |
449 |
gross |
660 |
a \Scalar in the \SolutionFS. |
450 |
jgs |
102 |
|
451 |
|
|
\begin{classdesc}{Poisson}{domain} |
452 |
|
|
opens a Poisson equation on the \Domain domain. \Poisson is derived from \LinearPDE. |
453 |
|
|
\end{classdesc} |
454 |
|
|
\begin{methoddesc}[Poisson]{setValue}{f=escript.Data(),q=escript.Data()} |
455 |
|
|
assigns new values to \var{f} and \var{q}. |
456 |
|
|
\end{methoddesc} |
457 |
gross |
625 |
|
458 |
gross |
999 |
\subsection{The \Helmholtz Class} |
459 |
gross |
660 |
The \Helmholtz class defines the Helmholtz problem |
460 |
|
|
\begin{equation}\label{HZ.1} |
461 |
|
|
\omega \; u - (k\; u\hackscore{,j})\hackscore{,j} = f |
462 |
|
|
\end{equation} |
463 |
ksteube |
1316 |
with natural boundary conditions |
464 |
gross |
660 |
\begin{equation}\label{HZ.2} |
465 |
|
|
k\; u\hackscore{,j} n\hackscore{,j} = g- \alpha \; u |
466 |
|
|
\end{equation} |
467 |
|
|
and constraints: |
468 |
|
|
\begin{equation}\label{HZ.3} |
469 |
|
|
u=r \mbox{ where } q>0 |
470 |
|
|
\end{equation} |
471 |
|
|
$\omega$, $k$, $f$ have to be a \Scalar in the \Function, |
472 |
|
|
$g$ and $\alpha$ must be a \Scalar in the \FunctionOnBoundary, |
473 |
|
|
and $q$ and $r$ must be a \Scalar in the \SolutionFS or must be mapped or interpolated into the particular \FunctionSpace. |
474 |
gross |
625 |
|
475 |
gross |
660 |
\begin{classdesc}{Helmholtz}{domain} |
476 |
|
|
opens a Helmholtz equation on the \Domain domain. \Helmholtz is derived from \LinearPDE. |
477 |
|
|
\end{classdesc} |
478 |
|
|
\begin{methoddesc}[Helmholtz]{setValue}{ \optional{omega} \optional{, k} \optional{, f} \optional{, alpha} \optional{, g} \optional{, r} \optional{, q}} |
479 |
|
|
assigns new values to \var{omega}, \var{k}, \var{f}, \var{alpha}, \var{g}, \var{r}, \var{q}. By default all values are set to be zero. |
480 |
|
|
\end{methoddesc} |
481 |
|
|
|
482 |
gross |
999 |
\subsection{The \Lame Class} |
483 |
gross |
660 |
The \Lame class defines a Lame equation problem: |
484 |
|
|
\begin{equation}\label{LE.1} |
485 |
|
|
-\mu (u\hackscore{i,j}+u\hackscore{j,i})+\lambda u\hackscore{k,k})\hackscore{j} = F\hackscore{i}-\sigma\hackscore{ij,j} |
486 |
|
|
\end{equation} |
487 |
ksteube |
1316 |
with natural boundary conditions: |
488 |
gross |
660 |
\begin{equation}\label{LE.2} |
489 |
|
|
n\hackscore{j}(\mu \; (u\hackscore{i,j}+u\hackscore{j,i})+\lambda*u\hackscore{k,k}) = f\hackscore{i}+n\hackscore{j}\sigma\hackscore{ij} |
490 |
|
|
\end{equation} |
491 |
|
|
and constraint |
492 |
|
|
\begin{equation}\label{LE.3} |
493 |
|
|
u\hackscore{i}=r\hackscore{i} \mbox{ where } q\hackscore{i}>0 |
494 |
|
|
\end{equation} |
495 |
|
|
$\mu$, $\lambda$ have to be a \Scalar in the \Function, |
496 |
|
|
$F$ has to be a \Vector in the \Function, |
497 |
|
|
$\sigma$ has to be a \Tensor in the \Function, |
498 |
|
|
$f$ must be a \Vector in the \FunctionOnBoundary, |
499 |
|
|
and $q$ and $r$ must be a \Vector in the \SolutionFS or must be mapped or interpolated into the particular \FunctionSpace. |
500 |
gross |
625 |
|
501 |
gross |
660 |
\begin{classdesc}{Lame}{domain} |
502 |
|
|
opens a Lame equation on the \Domain domain. \Lame is derived from \LinearPDE. |
503 |
|
|
\end{classdesc} |
504 |
|
|
\begin{methoddesc}[Lame]{setValue}{ \optional{lame_lambda} \optional{, lame_mu} \optional{, F} \optional{, sigma} \optional{, f} \optional{, r} \optional{, q}} |
505 |
|
|
assigns new values to |
506 |
|
|
\var{lame_lambda}, |
507 |
|
|
\var{lame_mu}, |
508 |
|
|
\var{F}, |
509 |
|
|
\var{sigma}, |
510 |
|
|
\var{f}, |
511 |
|
|
\var{r} and |
512 |
|
|
\var{q} |
513 |
|
|
By default all values are set to be zero. |
514 |
|
|
\end{methoddesc} |
515 |
|
|
|
516 |
gross |
2208 |
% \section{Transport Problems} |
517 |
|
|
% \label{SEC Transport} |