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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 625 - (show annotations)
Thu Mar 23 00:41:25 2006 UTC (13 years, 1 month ago) by gross
File MIME type: application/x-tex
File size: 14593 byte(s)
some updates and linearPDE.tex switched off
1 % $Id$
2 %
3 % Copyright © 2006 by ACcESS MNRF
4 % \url{http://www.access.edu.au
5 % Primary Business: Queensland, Australia.
6 % Licensed under the Open Software License version 3.0
7 % http://www.opensource.org/licenses/osl-3.0.php
8 %
9
10
11 \section{The First Steps}
12 \label{FirstSteps}
13
14 \begin{figure}
15 \centerline{\includegraphics[width=\figwidth]{figures/FirstStepDomain}}
16 \caption{Domain $\Omega=[0,1]^2$ with outer normal field $n$.}
17 \label{fig:FirstSteps.1}
18 \end{figure}
19
20 In this chapter we will give an introduction how to use \escript to solve
21 a partial differential equation \index{partial differential equation} (PDE \index{partial differential equation!PDE}). The reader should be familiar with Python. The knowledge presented at the Python tutorial at \url{http://docs.python.org/tut/tut.html}
22 is sufficient. It is helpful if the reader has some basic knowledge of PDEs \index{partial differential equation}.
23
24 The PDE \index{partial differential equation} we wish to solve is the Poisson equation \index{Poisson equation}
25 \begin{equation}
26 -\Delta u =f
27 \label{eq:FirstSteps.1}
28 \end{equation}
29 for the solution $u$. The function $f$ is the given right hand side. The domain of interest, denoted by $\Omega$
30 is the unit square
31 \begin{equation}
32 \Omega=[0,1]^2=\{ (x\hackscore 0;x\hackscore 1) | 0\le x\hackscore{0} \le 1 \mbox{ and } 0\le x\hackscore{1} \le 1 \}
33 \label{eq:FirstSteps.1b}
34 \end{equation}
35 The domain is shown in \fig{fig:FirstSteps.1}.
36
37 $\Delta$ denotes the Laplace operator\index{Laplace operator} which is defined by
38 \begin{equation}
39 \Delta u = (u\hackscore {,0})\hackscore{,0}+(u\hackscore{,1})\hackscore{,1}
40 \label{eq:FirstSteps.1.1}
41 \end{equation}
42 where, for any function $w$ and any direction $i$, $u\hackscore{,i}$
43 denotes the partial derivative \index{partial derivative} of $u$ with respect to $i$.
44 \footnote{Some readers
45 may be more familiar with the Laplace operator\index{Laplace operator} being written
46 as $\nabla^2$, and written in the form
47 \begin{equation*}
48 \nabla^2 u = \nabla^t \cdot \nabla u = \frac{\partial^2 u}{\partial x\hackscore 0^2}
49 + \frac{\partial^2 u}{\partial x\hackscore 1^2}
50 \end{equation*}
51 and \eqn{eq:FirstSteps.1} as
52 \begin{equation*}
53 -\nabla^2 u = f
54 \end{equation*}
55 }
56 Basically, in the subindex of a function, any index to the left of the comma denotes a spatial derivative with respect
57 to the index. To get a more compact form we will write $w\hackscore{,ij}=(w\hackscore {,i})\hackscore{,j}$
58 which leads to
59 \begin{equation}
60 \Delta u = u\hackscore{,00}+u\hackscore{,11}=\sum\hackscore{i=0}^2 u\hackscore{,ii}
61 \label{eq:FirstSteps.1.1b}
62 \end{equation}
63 In some cases, and we will see examples for this in the next chapter,
64 the usage of the nested $\sum$ symbols blows up the formulas and therefore
65 it is convenient to use the Einstein summation convention \index{summation convention}. This
66 drops the $\sum$ sign and assumes that a summation over a repeated index is performed
67 ("repeated index means summation"). For instance we write
68 \begin{eqnarray}
69 x\hackscore{i}y\hackscore{i}=\sum\hackscore{i=0}^2 x\hackscore{i}y\hackscore{i} \\
70 x\hackscore{i}u\hackscore{,i}=\sum\hackscore{i=0}^2 x\hackscore{i}u\hackscore{,i} \\
71 u\hackscore{,ii}=\sum\hackscore{i=0}^2 u\hackscore{,ii} \\
72 x\hackscore{ij}u\hackscore{i,j}=\sum\hackscore{j=0}^2\sum\hackscore{i=0}^2 x\hackscore{ij}u\hackscore{i,j} \\
73 \label{eq:FirstSteps.1.1c}
74 \end{eqnarray}
75 With the summation convention we can write the Poisson equation \index{Poisson equation} as
76 \begin{equation}
77 - u\hackscore{,ii} =1
78 \label{eq:FirstSteps.1.sum}
79 \end{equation}
80 where $f=1$ in this example.
81
82 On the boundary of the domain $\Omega$ the normal derivative $n\hackscore{i} u\hackscore{,i}$
83 of the solution $u$ shall be zero, ie. $u$ shall fulfill
84 the homogeneous Neumann boundary condition\index{Neumann
85 boundary condition!homogeneous}
86 \begin{equation}
87 n\hackscore{i} u\hackscore{,i}= 0 \;.
88 \label{eq:FirstSteps.2}
89 \end{equation}
90 $n=(n\hackscore{i})$ denotes the outer normal field
91 of the domain, see \fig{fig:FirstSteps.1}. Remember that we
92 are applying the Einstein summation convention \index{summation convention}, i.e
93 $n\hackscore{i} u\hackscore{,i}= n\hackscore{0} u\hackscore{,0} +
94 n\hackscore{1} u\hackscore{,1}$.
95 \footnote{Some readers may familiar with the notation
96 \begin{equation*}
97 \frac{\partial u}{\partial n} = n\hackscore{i} u\hackscore{,i}
98 \end{equation*}
99 for the normal derivative.}
100 The Neumann boundary condition of \eqn{eq:FirstSteps.2} should be fulfilled on the
101 set $\Gamma^N$ which is the top and right edge of the domain:
102 \begin{equation}
103 \Gamma^N=\{(x\hackscore 0;x\hackscore 1) \in \Omega | x\hackscore{0}=1 \mbox{ or } x\hackscore{1}=1 \}
104 \label{eq:FirstSteps.2b}
105 \end{equation}
106 On the bottom and the left edge of the domain which is defined
107 as
108 \begin{equation}
109 \Gamma^D=\{(x\hackscore 0;x\hackscore 1) \in \Omega | x\hackscore{0}=0 \mbox{ or } x\hackscore{1}=0 \}
110 \label{eq:FirstSteps.2c}
111 \end{equation}
112 the solution shall be identically zero:
113 \begin{equation}
114 u=0 \; .
115 \label{eq:FirstSteps.2d}
116 \end{equation}
117 This kind of boundary condition is called a homogeneous Dirichlet boundary condition
118 \index{Dirichlet boundary condition!homogeneous}. The partial differential equation in \eqn{eq:FirstSteps.1.sum} together
119 with the Neumann boundary condition \eqn{eq:FirstSteps.2} and
120 Dirichlet boundary condition in \eqn{eq:FirstSteps.2d} form a so
121 called boundary value
122 problem\index{boundary value problem} (BVP\index{boundary value problem!BVP}) for
123 the unknown
124 function $u$.
125
126
127 \begin{figure}
128 \centerline{\includegraphics[width=\figwidth]{figures/FirstStepMesh.eps}}
129 \caption{Mesh of $4 \time 4$ elements on a rectangular domain. Here
130 each element is a quadrilateral and described by four nodes, namely
131 the corner points. The solution is interpolated by a bi-linear
132 polynomial.}
133 \label{fig:FirstSteps.2}
134 \end{figure}
135
136 In general the BVP\index{boundary value problem!BVP} cannot be solved analytically and numerical
137 methods have to be used construct an approximation of the solution
138 $u$. Here we will use the finite element method\index{finite element
139 method} (FEM\index{finite element
140 method!FEM}). The basic idea is to fill the domain with a
141 set of points called nodes. The solution is approximated by its
142 values on the nodes\index{finite element
143 method!nodes}. Moreover, the domain is subdivided into smaller
144 sub-domains called elements \index{finite element
145 method!element}. On each element the solution is
146 represented by a polynomial of a certain degree through its values at
147 the nodes located in the element. The nodes and its connection through
148 elements is called a mesh\index{finite element
149 method!mesh}. \fig{fig:FirstSteps.2} shows an
150 example of a FEM mesh with four elements in the $x_0$ and four elements
151 in the $x_1$ direction over the unit square.
152 For more details we refer the reader to the literature, for instance
153 \Ref{Zienc,NumHand}.
154
155 \escript provides the class \Poisson to define a Poisson equation \index{Poisson equation}.
156 (We will discuss a more general form of a PDE \index{partial differential equation!PDE}
157 that can be defined through the \LinearPDE class later). The instantiation of
158 a \Poisson class object requires the specification of the domain $\Omega$. In \escript
159 the \Domain class objects are used to describe the geometry of a domain but it also
160 contains information about the discretization methods and the actual solver which is used
161 to solve the PDE. Here we are using the FEM library \finley \index{finite element
162 method}. The following statements create the \Domain object \var{mydomain} from the
163 \finley method \method{Rectangle}
164 \begin{python}
165 from esys.finley import Rectangle
166 mydomain = Rectangle(l0=1.,l1=1.,n0=40, n1=20)
167 \end{python}
168 In this case the domain is a rectangle with the lower, left corner at point $(0,0)$ and
169 the right, upper corner at $(\var{l0},\var{l1})=(1,1)$.
170 The arguments \var{n0} and \var{n1} define the number of elements in $x\hackscore{0}$ and
171 $x\hackscore{1}$-direction respectively. For more details on \method{Rectangle} and
172 other \Domain generators within the \finley module,
173 see \Chap{CHAPTER ON FINLEY}.
174
175 The following statements define the \Poisson class object \var{mypde} with domain \var{mydomain} and
176 the right hand side $f$ of the PDE to constant $1$:
177 \begin{python}
178 from esys.escript.linearPDEs import Poisson
179 mypde = Poisson(mydomain)
180 mypde.setValue(f=1)
181 \end{python}
182 We have not specified any boundary condition but the
183 \Poisson class implicitly assumes homogeneous Neuman boundary conditions \index{Neumann
184 boundary condition!homogeneous} defined by \eqn{eq:FirstSteps.2}. With this boundary
185 condition the BVP\index{boundary value problem!BVP} we have defined has no unique solution. In fact, with any solution $u$
186 and any constant $C$ the function $u+C$ becomes a solution as well. We have to add
187 a Dirichlet boundary condition \index{Dirichlet boundary condition}. This is done
188 by defining a characteristic function \index{characteristic function}
189 which has positive values at locations $x=(x\hackscore{0},x\hackscore{1})$ where Dirichlet boundary condition is set
190 and $0$ elsewhere. In our case of $\Gamma^D$ defined by \eqn{eq:FirstSteps.2c},
191 we need to construct a function \var{gammaD} which is positive for the cases $x\hackscore{0}=0$ or $x\hackscore{1}=0$. To get
192 an object \var{x} which represents locations in the domain one uses
193 \begin{python}
194 x=mydomain.getX()
195 \end{python}
196 In fact \var{x} is a \Data object which we will learn more about in Chapter~\ref{X}. At this stage we only have to know
197 that \var{x} has a
198
199 In the first statement, the method \method{getX} of the \Domain \var{mydomain}
200 gives access to locations
201 in the domain defined by \var{mydomain}. The object \var{x} is actually a \Data object which is
202 discussed in Chpater\ref{X} in more details. What we need to know here is that
203 \var{x} has \Rank (=number of dimensions) and a \Shape (=tuple of dimensions) which can be checked by
204 calling the \method{getRank} and \method{getShape} methods:
205 \begin{python}
206 print "rank ",x.getRank(),", shape ",x.getShape()
207 \end{python}
208 will print something like
209 \begin{python}
210 rank 1, shape (2,)
211 \end{python}
212 The \Data object also maintains type information which is represented by the
213 \FunctionSpace of the object. For instance
214 \begin{python}
215 print x.getFunctionSpace()
216 \end{python}
217 will print
218 \begin{python}
219 Function space type: Finley_Nodes on FinleyMesh
220 \end{python}
221 which tells us that the coordinates are stored on the nodes of a \finley mesh.
222 To get the $x\hackscore{0}$ coordinates of the locations we use the
223 statement
224 \begin{python}
225 x0=x[0]
226 \end{python}
227 Object \var{x0}
228 is again a \Data object now with \Rank $0$ and
229 \Shape $()$. It inherits the \FunctionSpace from \var{x}:
230 \begin{python}
231 print x0.getRank(),x0.getShape(),x0.getFunctionSpace()
232 \end{python}
233 will print
234 \begin{python}
235 0 () Function space type: Finley_Nodes on FinleyMesh
236 \end{python}
237 We can now construct the function \var{gammaD} by
238 \begin{python}
239 from esys.escript import whereZero
240 gammaD=whereZero(x[0])+whereZero(x[1])
241 \end{python}
242 where
243 \code{whereZero(x[0])} creates function which equals $1$ where \code{x[0]} is (allmost) equal to zero
244 and $0$ elsewhere.
245 Similarly, \code{whereZero(x[1])} creates function which equals $1$ where \code{x[1]} is
246 equal to zero and $0$ elsewhere.
247 The sum of the results of \code{whereZero(x[0])} and \code{whereZero(x[1])}
248 gives a function on the domain \var{mydomain} which is exactly positive where $x\hackscore{0}$ or $x\hackscore{1}$ is equal to zero.
249 Note that \var{gammaD} has the same \Rank, \Shape and \FunctionSpace like \var{x0} used to define it. So from
250 \begin{python}
251 print gammaD.getRank(),gammaD.getShape(),gammaD.getFunctionSpace()
252 \end{python}
253 one gets
254 \begin{python}
255 0 () Function space type: Finley_Nodes on FinleyMesh
256 \end{python}
257 The additional parameter \var{q} of the \code{setValue} method of the \Poisson class defines the
258 characteristic function \index{characteristic function} of the locations
259 of the domain where homogeneous Dirichlet boundary condition \index{Dirichlet boundary condition!homogeneous}
260 are set. The complete definition of our example is now:
261 \begin{python}
262 from esys.linearPDEs import Poisson
263 x = mydomain.getX()
264 gammaD = whereZero(x[0])+whereZero(x[1])
265 mypde = Poisson(domain=mydomain)
266 mypde.setValue(f=1,q=gammaD)
267 \end{python}
268 The first statement imports the \Poisson class definition from the \linearPDEs module \escript package.
269 To get the solution of the Poisson equation defined by \var{mypde} we just have to call its
270 \method{getSolution}.
271
272 Now we can write the script to solve our test problem (Remember that
273 lines starting with '\#' are comment lines in Python) (available as \file{poisson.py}
274 in the \ExampleDirectory):
275 \begin{python}
276 from esys.escript import *
277 from esys.escript.linearPDEs import Poisson
278 from esys.finley import Rectangle
279 # generate domain:
280 mydomain = Rectangle(l0=1.,l1=1.,n0=40, n1=20)
281 # define characteristic function of Gamma^D
282 x = mydomain.getX()
283 gammaD = whereZero(x[0])+whereZero(x[1])
284 # define PDE and get its solution u
285 mypde = Poisson(domain=mydomain)
286 mypde.setValue(f=1,q=gammaD)
287 u = mypde.getSolution()
288 # write u to an external file
289 saveVTK("u.xml",sol=u)
290 \end{python}
291 The last statement writes the solution tagged with the name "sol" to the external file \file{u.xml} in
292 \VTK file format. \VTK is a software library
293 for the visualization of scientific, engineering and analytical data and is freely available
294 from \url{http://www.vtk.org}. There are a variety of graphical user interfaces
295 for \VTK available, for instance \mayavi which can be downloaded from \url{http://mayavi.sourceforge.net/} but is also available on most
296 \LINUX distributions.
297
298 \begin{figure}
299 \centerline{\includegraphics[width=\figwidth]{figures/FirstStepResult.eps}}
300 \caption{Visualization of the Poisson Equation Solution for $f=1$}
301 \label{fig:FirstSteps.3}
302 \end{figure}
303
304 You can edit the script file using your favourite text editor (or the Integrated DeveLopment Environment IDLE
305 for Python, see \url{http://idlefork.sourceforge.net}). If the script file has the name \file{poisson.py} \index{scripts!\file{poisson.py}} you can run the
306 script from any shell using the command:
307 \begin{python}
308 python poisson.py
309 \end{python}
310 After the script has (hopefully successfully) been completed you will find the file \file{u.xml} in the current
311 directory. An easy way to visualize the results is the command
312 \begin{python}
313 mayavi -d u.xml -m SurfaceMap &
314 \end{python}
315 to show the results, see \fig{fig:FirstSteps.3}.

Properties

Name Value
svn:eol-style native
svn:keywords Author Date Id Revision

  ViewVC Help
Powered by ViewVC 1.1.26