/[escript]/trunk/doc/cookbook/example01.tex
ViewVC logotype

Contents of /trunk/doc/cookbook/example01.tex

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2870 - (show annotations)
Mon Jan 25 08:34:15 2010 UTC (9 years, 7 months ago) by gross
Original Path: trunk/doc/cookbook/onedheatdiff001.tex
File MIME type: application/x-tex
File size: 32859 byte(s)
more review on the cookbook.
1
2 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
3 %
4 % Copyright (c) 2003-2009 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 We will start by examining a simple one dimensional heat diffusion example. This problem will provide a good launch pad to build our knowledge of \esc and demonstrate how to solve simple partial differential equations (PDEs)\footnote{Wikipedia provides an excellent and comprehensive introduction to \textit{Partial Differential Equations} \url{http://en.wikipedia.org/wiki/Partial_differential_equation}, however their relevance to \esc and implementation should become a clearer as we develop our understanding further into the cookbook.}
15
16 \section{One Dimensional Heat Diffusion in an Iron Rod}
17 \sslist{onedheatdiff001.py and cblib.py}
18 %\label{Sec:1DHDv0}
19 The first model consists of a simple cold iron bar at a constant temperature of zero \reffig{fig:onedhdmodel}. The bar is perfectly insulated on all sides with a heating element at one end. Intuition tells us that as heat is applied; energy will disperse along the bar via conduction. With time the bar will reach a constant temperature equivalent to that of the heat source.
20 \begin{figure}[h!]
21 \centerline{\includegraphics[width=4.in]{figures/onedheatdiff}}
22 \caption{One dimensional model of an Iron bar.}
23 \label{fig:onedhdmodel}
24 \end{figure}
25 \subsection{1D Heat Diffusion Equation}
26 We can model the heat distribution of this problem over time using the one dimensional heat diffusion equation\footnote{A detailed discussion on how the heat diffusion equation is derived can be found at \url{http://online.redwoods.edu/instruct/darnold/DEProj/sp02/AbeRichards/paper.pdf}};
27 which is defined as:
28 \begin{equation}
29 \rho c\hackscore p \frac{\partial T}{\partial t} - \kappa \frac{\partial^{2} T}{\partial x^{2}} = q\hackscore H
30 \label{eqn:hd}
31 \end{equation}
32 where $\rho$ is the material density, $c\hackscore p$ is the specific heat and $\kappa$ is the thermal
33 conductivity\footnote{A list of some common thermal conductivities is available from Wikipedia \url{http://en.wikipedia.org/wiki/List_of_thermal_conductivities}}. Here we assume that these material
34 parameters are \textbf{constant}.
35 The heat source is defined by the right hand side of \refEq{eqn:hd} as $q\hackscore{H}$; this can take the form of a constant or a function of time and space. For example $q\hackscore{H} = q\hackscore{0}e^{-\gamma t}$ where we have the output of our heat source decaying with time. There are also two partial derivatives in \refEq{eqn:hd}; $\frac{\partial T}{\partial t}$ describes the change in temperature with time while $\frac{\partial ^2 T}{\partial x^2}$ is the spatial change of temperature. As there is only a single spatial dimension to our problem, our temperature solution $T$ is only dependent on the time $t$ and our position along the iron bar $x$.
36
37 \subsection{PDEs and the General Form}
38 Potentially, it is now possible to solve PDE \refEq{eqn:hd} analytically and this would produce an exact solution to our problem. However, it is not always possible or practical to solve a problem this way. Alternatively, computers can be used to solve these kinds of problems. To do this, a numerical approach is required to discretised
39 the PDE \refEq{eqn:hd} in time and space so finally we are left with a finite number of equations for a finite number of spatial and time steps in the model. While discretization introduces approximations and a degree of error, we find that a sufficiently sampled model is generally accurate enough for the requirements of the modeler.
40
41 Firstly, we will discretise the PDE \refEq{eqn:hd} in the time direction which will
42 leave as with a steady linear PDE which is involving spatial derivatives only and needs to be solved in each time
43 step to progress in time - \esc can help us here.
44
45 For the discretization in time we will use is the Backwards Euler approximation scheme\footnote{see \url{http://en.wikipedia.org/wiki/Euler_method}}. It bases on the
46 approximation
47 \begin{equation}
48 \frac{\partial T(t)}{\partial t} \approx \frac{T(t)-T(t-h)}{h}
49 \label{eqn:beuler}
50 \end{equation}
51 for $\frac{\partial T}{\partial t}$ at time $t$
52 where $h$ is the time step size. This can also be written as;
53 \begin{equation}
54 \frac{\partial T}{\partial t}(t^{(n)}) \approx \frac{T^{(n)} - T^{(n-1)}}{h}
55 \label{eqn:Tbeuler}
56 \end{equation}
57 where the upper index $n$ denotes the n\textsuperscript{th} time step. So one has
58 \begin{equation}
59 \begin{array}{rcl}
60 t^{(n)} & = & t^{(n-1)}+h \\
61 T^{(n)} & = & T(t^{(n-1)}) \\
62 \end{array}
63 \label{eqn:Neuler}
64 \end{equation}
65 Substituting \refEq{eqn:Tbeuler} into \refEq{eqn:hd} we get;
66 \begin{equation}
67 \frac{\rho c\hackscore p}{h} (T^{(n)} - T^{(n-1)}) - \kappa \frac{\partial^{2} T^{(n)}}{\partial x^{2}} = q\hackscore H
68 \label{eqn:hddisc}
69 \end{equation}
70 Notice that we evaluate the spatial derivative term at current time $t^{(n)}$ - therefore the name \textbf{backward Euler} scheme. Alternatively, one can use evaluate the spatial derivative term at the previous time $t^{(n-1)}$. This
71 approach is called the \textbf{forward Euler} scheme. This scheme can provide some computational advantages which
72 we are not discussed here but has the major disadvantage that depending on the
73 material parameter as well as the discretiztion of the spatial derivative term the time step size $h$ needs to be chosen sufficiently small to achieve a stable temperature when progressing in time. The term \textit{stable} means
74 that the approximation of the temperature will not grow beyond its initial bounds and becomes unphysical.
75 The backward Euler which we use here is unconditionally stable meaning that under the assumption of
76 physically correct problem set-up the temperature approximation remains physical for all times.
77 The user needs to keep in mind that the discretization error introduced by \refEq{eqn:beuler}
78 is sufficiently small so a good approximation of the true temperature is calculated. It is
79 therefore crucial that the user remains critical about his/her results and for instance compares
80 the results for different time and spatial step sizes.
81
82 To get the temperature $T^{(n)}$ at time $t^{(n)}$ we need to solve the linear
83 differential equation \refEq{eqn:hddisc} which is only including spatial derivatives. To solve this problem
84 we want to to use \esc.
85
86 \esc interfaces with any given PDE via a general form. For the purpose of this introduction we will illustrate a simpler version of the full linear PDE general form which is available in the \esc user's guide. A simplified form that suits our heat diffusion problem\footnote{In the form of the \esc users guide which using the Einstein convention is written as
87 $-(A\hackscore{jl} u\hackscore{,l})\hackscore{,j}+D u =Y$}
88 is described by;
89 \begin{equation}\label{eqn:commonform nabla}
90 -\nabla\cdot(A\cdot\nabla u) + Du = f
91 \end{equation}
92 where $A$, $D$ and $f$ are known values and $u$ is the unknown solution. The symbol $\nabla$ which is called the \textit{Nabla operator} or \textit{del operator} represents
93 the spatial derivative of its subject - in this case $u$. Lets assume for a moment that we deal with a one-dimensional problem then ;
94 \begin{equation}
95 \nabla = \frac{\partial}{\partial x}
96 \end{equation}
97 and we can write \refEq{eqn:commonform nabla} as;
98 \begin{equation}\label{eqn:commonform}
99 -A\frac{\partial^{2}u}{\partial x^{2}} + Du = f
100 \end{equation}
101 if $A$ is constant. To match this simplified general form to our problem \refEq{eqn:hddisc}
102 we rearrange \refEq{eqn:hddisc};
103 \begin{equation}
104 \frac{\rho c\hackscore p}{h} T^{(n)} - \kappa \frac{\partial^2 T^{(n)}}{\partial x^2} = q\hackscore H + \frac{\rho c\hackscore p}{h} T^{(n-1)}
105 \label{eqn:hdgenf}
106 \end{equation}
107 The PDE is now in a form that satisfies \refEq{eqn:commonform nabla} which is required for \esc to solve our PDE. This can be done by generating a solution for successive increments in the time nodes $t^{(n)}$ where
108 $t^{(0)}=0$ and $t^{(n)}=t^{(n-1)}+h$ where $h>0$ is the step size and assumed to be constant.
109 In the following the upper index ${(n)}$ refers to a value at time $t^{(n)}$. Finally, by comparing \refEq{eqn:hdgenf} with \refEq{eqn:commonform} it can be seen that;
110 \begin{equation}\label{ESCRIPT SET}
111 u=T^{(n)};
112 A = \kappa; D = \frac{\rho c \hackscore{p}}{h}; f = q \hackscore{H} + \frac{\rho c\hackscore p}{h} T^{(n-1)}
113 \end{equation}
114
115 % Now that the general form has been established, it can be submitted to \esc. Note that it is necessary to establish the state of our system at time zero or $T^{(n=0)}$. This is due to the time derivative approximation we have used from \refEq{eqn:Tbeuler}. Our model stipulates a starting temperature in the iron bar of 0\textcelsius. Thus the temperature distribution is simply;
116 % \begin{equation}
117 % T(x,0) = \left
118 % \end{equation}
119 % for all $x$ in the domain.
120
121 \subsection{Boundary Conditions}
122 With the PDE sufficiently modified, consideration must now be given to the boundary conditions of our model. Typically there are two main types of boundary conditions known as \textbf{Neumann} and \textbf{Dirichlet} boundary conditions\footnote{More information on Boundary Conditions is available at Wikipedia \url{http://en.wikipedia.org/wiki/Boundary_conditions}}, respectively.
123 A \textbf{Dirichlet boundary condition} is conceptually simpler and is used to prescribe a known value to the unknown - in our example the temperature - on parts of or the entire boundary of the region of interest.
124 For our model problem we want to keep the initial temperature setting on the left side of the
125 iron bar over time. This defines a Dirichlet boundary condition for the PDE \refEq{eqn:hddisc} to be solved at each time step.
126
127 On the other end of the iron rod we want to add an appropriate boundary condition to define insolation to prevent
128 any loss or inflow of energy at the right end of the rod. Mathematically this is expressed by prescribing
129 the heat flux $\kappa \frac{\partial T}{\partial x}$ to zero on the right end of the rod
130 In our simplified one dimensional model this is expressed
131 in the form;
132 \begin{equation}
133 \kappa \frac{\partial T}{\partial x} = 0
134 \end{equation}
135 or in a more general case as
136 \begin{equation}\label{NEUMAN 1}
137 \kappa \nabla T \cdot n = 0
138 \end{equation}
139 where $n$ is the outer normal field \index{outer normal field} at the surface of the domain.
140 For the iron rod the outer normal field on the right hand side is the vector $(1,0)$. The $\cdot$ (dot) refers to the
141 dot product of the vectors $\nabla T$ and $n$. In fact, the term $\nabla T \cdot n$ is the normal derivative of
142 the temperature $T$. Other notations which are used are\footnote{The \esc notation for the normal
143 derivative is $T\hackscore{,i} n\hackscore i$.};
144 \begin{equation}
145 \nabla T \cdot n = \frac{\partial T}{\partial n} \; .
146 \end{equation}
147 A condition of the type \refEq{NEUMAN 1} defines a \textbf{Neuman boundary condition} for the PDE.
148
149 The PDE \refEq{eqn:hdgenf} together with the Dirichlet boundary condition set on the left face of the rod
150 and the Neuman boundary condition~\ref{eqn:hdgenf} define a \textbf{boundary value problem}.
151 It is a nature of a boundary value problem that it allows to make statements on the solution in the
152 interior of the domain from information known on the boundary only. In most cases
153 we use the term partial differential equation but in fact mean a boundary value problem.
154 It is important to keep in mind that boundary conditions need to be complete and consistent in the sense that
155 at any point on the boundary either a Dirichlet or a Neuman boundary condition must be set.
156
157 Conviniently, \esc makes default assumption on the boundary conditions which the user may modify where appropriate.
158 For a problem of the form in~\refEq{eqn:commonform nabla} the default condition\footnote{In the form of the \esc users guide which is using the Einstein convention is written as
159 $n\hackscore{j}A\hackscore{jl} u\hackscore{,l}=0$.} is;
160 \begin{equation}\label{NEUMAN 2}
161 n\cdot A \cdot\nabla u = 0
162 \end{equation}
163 which is used everywhere on the boundary. Again $n$ denotes the outer normal field.
164 Notice that the coefficient $A$ is the same as in the \esc PDE~\ref{eqn:commonform nabla}.
165 With the settings for the coefficients we have already identified in \refEq{ESCRIPT SET} this
166 condition translates into
167 \begin{equation}\label{NEUMAN 2b}
168 \kappa \frac{\partial T}{\partial x} = 0
169 \end{equation}
170 for the right hand side of the rod. This is identical to the Neuman boundary condition we want to set. \esc will take care of this condition for us. We will discuss the Dirichlet boundary condition later.
171
172 \subsection{Outline of the Implementation}
173 \label{sec:outline}
174 To solve the heat diffusion equation (equation \refEq{eqn:hd}) we will write a simple \pyt script. At this point we assume that you have some basic understanding of the \pyt programming language. If not there are some pointers and links available in Section \ref{sec:escpybas}. The script we will discuss later in details will have four major steps. Firstly we need to define the domain where we want to
175 calculate the temperature. For our problem this is the iron rod which has a rectangular shape. Secondly we need to define the PDE
176 we need to solve in each time step to get the updated temperature. Thirdly we need to define the the coefficients of the PDE and finally we need to solve the PDE. The last two steps need to be repeated until the final time marker has been reached.
177
178 In the terminology of \pyt the domain and PDE are represented by \textbf{objects}. The nice feature of an object is that it defined by it usage and features
179 rather than its actual representation. So we will create a domain object to decribe the geometry of our iron rod. The main feature
180 of the object we will use is the fact that we can define PDEs and spatially distributed values such as the temperature
181 on a domain. In fact the domain object has many more features - most of them you will
182 never use and do not need to understand. Similar a PDE object is defined by the fact that we can define the coefficients of the PDE and solve the PDE. At a
183 later stage you may use more advanved features of the PDE class but you need to worry about them only at the point when you use them.
184
185
186 \begin{figure}[t]
187 \centering
188 \includegraphics[width=6in]{figures/functionspace.pdf}
189 \label{fig:fs}
190 \caption{\esc domain construction overview}
191 \end{figure}
192
193 \subsection{The Domain Constructor in \esc}
194 \label{ss:domcon}
195 It is helpful to have a better understanding how spatially distributed value such as the temperature or PDE coefficients are interpreted in \esc. Again
196 from the user's point of view the representation of these spatially distributed values is not relevant.
197
198 There are various ways to construct domain objects. The simplest form is as rectangular shaped region with a length and height. There is
199 a ready to use function call for this. Besides the spatial dimensions the function call will require you to specify the number
200 elements or cells to be used along the length and height, see Figure~\ref{fig:fs}. Any spacially distributed value
201 and the PDE is represented in discrete form using this element representation\footnote{We will use the finite element method (FEM), see \url{http://en.wikipedia.org/wiki/Finite_element_method}i for details.}. Therefore we will have access to an approximation of the true PDE solution only.
202 The quality of the approximation depends - besides other factors- mainly on the number of elements being used. In fact, the
203 approximation becomes better the more elements are used. However, computational costs and compute time grow with the number of
204 elements being used. It therefore important that you find the right balance between the demand in accuarcy and acceptable resource usage.
205
206 In general, one can thinks about a domain object as a composition of nodes and elements.
207 As shown in Figure~\ref{fig:fs}, an element is defined by the nodes used to describe its vertices.
208 To represent spatial distributed values the user can use
209 the values at the nodes, at the elements in the interior of the domain or at elements located at the surface of the domain.
210 The different approach used to represent values is called \textbf{function space} and is attached to all objects
211 in \esc representing a spatial distributed value such as the solution of a PDE. The three
212 function spaces we will use at the moment are;
213 \begin{enumerate}
214 \item the nodes, called by \verb|ContinuousFunction(domain)| ;
215 \item the elements/cells, called by \verb|Function(domain)| ; and
216 \item the boundary, called by \verb|FunctionOnBoundary(domain)| .
217 \end{enumerate}
218 A function space object such as \verb|ContinuousFunction(domain)| has the method \verb|getX| attached to it. This method returns the
219 location of the so-called \textbf{sample points} used to represent values with the particular function space attached to it. So the
220 call \verb|ContinuousFunction(domain).getX()| will return the coordinates of the nodes used to describe the domain while
221 the \verb|Function(domain).getX()| returns the coordinates of numerical integration points within elements, see
222 Figure~\ref{fig:fs}.
223
224 This distiction between different representations of spatial distributed values
225 is important in order to be able to vary the degrees of smoothness in a PDE problem.
226 The coefficients of a PDE need not be continuous thus this qualifies as a \verb|Function()| type.
227 On the other hand a temperature distribution must be continuous and needs to be represented with a \verb|ContinuousFunction()| function space.
228 An influx may only be defined at the boundary and is therefore a \verb FunctionOnBoundary() object.
229 \esc allows certain transformations of the function spaces. A \verb ContinuousFunction() can be transformed into a \verb|FunctionOnBoundary()|
230 or \verb|Function()|. On the other hand there is not enough information in a \verb FunctionOnBoundary() to transform it to a \verb ContinuousFunction() .
231 These transformations, which ar ecalled \textbf{interpolation} are invoked autmatically by \esc if needed.
232
233 Later in this introduction we will discuss how
234 to define specific areas of geometry with different materials which are represented by different material coefficients such the
235 thermal conductivities $kappa$. A very powerful techique to define these types of PDE
236 coefficients is tagging. Blocks of materials and boundaries can be named and values can be defined on subregions based on their names.
237 This is simplifing PDE coefficient and flux definitions. It makes for much easier scripting. We will discuss this technique in Section~\ref{STEADY-STATE HEAT REFRACTION}.
238
239
240 \subsection{A Clarification for the 1D Case}
241 It is necessary for clarification that we revisit the general PDE from \refeq{eqn:commonform nabla} under the light of a two dimensional domain. \esc is inherently designed to solve problems that are greater than one dimension and so \refEq{eqn:commonform nabla} needs to be read as a higher dimensional problem. In the case of two spatial dimensions the \textit{Nabla operator} has in fact two components $\nabla = (\frac{\partial}{\partial x}, \frac{\partial}{\partial y})$. In full, \refEq{eqn:commonform nabla} assuming a constant coefficient $A$, takes the form;
242 \begin{equation}\label{eqn:commonform2D}
243 -A\hackscore{00}\frac{\partial^{2}u}{\partial x^{2}}
244 -A\hackscore{01}\frac{\partial^{2}u}{\partial x\partial y}
245 -A\hackscore{10}\frac{\partial^{2}u}{\partial y\partial x}
246 -A\hackscore{11}\frac{\partial^{2}u}{\partial y^{2}}
247 + Du = f
248 \end{equation}
249 Notice that for the higher dimensional case $A$ becomes a matrix. It is also
250 important to notice that the usage of the Nabla operator creates
251 a compact formulation which is also independent from the spatial dimension.
252 So to make the general PDE \refEq{eqn:commonform2D} one dimensional as
253 shown in \refEq{eqn:commonform} we need to set
254 \begin{equation}
255 A\hackscore{00}=A; A\hackscore{01}=A\hackscore{10}=A\hackscore{11}=0
256 \end{equation}
257
258
259 \subsection{Developing a PDE Solution Script}
260 \label{sec:key}
261 We will write a simple \pyt script which uses the \modescript, \modfinley and \modmpl modules.
262 By developing a script for \esc, the heat diffusion equation can be solved at successive time steps for a predefined period using our general form \refEq{eqn:hdgenf}. Firstly it is necessary to import all the libraries\footnote{The libraries contain predefined scripts that are required to solve certain problems, these can be simple like $sine$ and $cosine$ functions or more complicated like those from our \esc library.}
263 that we will require.
264 \begin{python}
265 from esys.escript import *
266 # This defines the LinearPDE module as LinearPDE
267 from esys.escript.linearPDEs import LinearPDE
268 # This imports the rectangle domain function from finley.
269 from esys.finley import Rectangle
270 # A useful unit handling package which will make sure all our units
271 # match up in the equations under SI.
272 from esys.escript.unitsSI import *
273 import pylab as pl #Plotting package.
274 import numpy as np #Array package.
275 import os #This package is necessary to handle saving our data.
276 \end{python}
277 It is generally a good idea to import all of the \modescript library, although if the functions and classes required are known they can be specified individually. The function \verb|LinearPDE| has been imported explicitly for ease of use later in the script. \verb|Rectangle| is going to be our type of model. The module \verb unitsSI provides support for SI unit definitions with our variables; and the \verb|os| module is needed to handle file outputs once our PDE has been solved. \verb pylab and \verb numpy are modules developed independently of \esc. They are used because they have efficient plotting and array handling capabilities.
278
279 Once our library dependencies have been established, defining the problem specific variables is the next step. In general the number of variables needed will vary between problems. These variables belong to two categories. They are either directly related to the PDE and can be used as inputs into the \esc solver, or they are script variables used to control internal functions and iterations in our problem. For this PDE there are a number of constants which will need values. Firstly, the model upon which we wish to solve our problem needs to be defined. There are many different types of models in \modescript which we will demonstrate in later tutorials but for our iron rod, we will simply use a rectangular model.
280
281 Using a rectangular model simplifies our rod which would be a \textit{3D} object, into a single dimension. The iron rod will have a lengthways cross section that looks like a rectangle. As a result we do not need to model the volume of the rod because a cylinder is symmetrical about its centre. There are four arguments we must consider when we decide to create a rectangular model, the model \textit{length}, \textit{width} and \textit{step size} in each direction. When defining the size of our problem it will help us determine appropriate values for our model arguments. If we make our dimensions large but our step sizes very small we will to a point, increase the accuracy of our solution. Unfortunately we also increase the number of calculations that must be solved per time step. This means more computational time is required to produce a solution. In this \textit{1D} problem, the bar is defined as being 1 metre long. An appropriate step size \verb|ndx| would be 1 to 10\% of the length. Our \verb|ndy| need only be 1, this is because our problem stipulates no partial derivatives in the $y$ direction. Thus the temperature does not vary with $y$. Hence, the model parameters can be defined as follows; note we have used the \verb unitsSI convention to make sure all our input units are converted to SI.
282 \begin{python}
283 #Domain related.
284 mx = 1*m #meters - model length
285 my = .1*m #meters - model width
286 ndx = 100 # mesh steps in x direction
287 ndy = 1 # mesh steps in y direction - one dimension means one element
288 \end{python}
289 The material constants and the temperature variables must also be defined. For the iron rod in the model they are defined as:
290 \begin{python}
291 #PDE related
292 q=200. * Celsius #Kelvin - our heat source temperature
293 Tref = 0. * Celsius #Kelvin - starting temp of iron bar
294 rho = 7874. *kg/m**3 #kg/m^{3} density of iron
295 cp = 449.*J/(kg*K) #j/Kg.K thermal capacity
296 rhocp = rho*cp
297 kappa = 80.*W/m/K #watts/m.Kthermal conductivity
298 \end{python}
299 Finally, to control our script we will have to specify our timing controls and where we would like to save the output from the solver. This is simple enough:
300 \begin{python}
301 t=0 #our start time, usually zero
302 tend=5.*minute #seconds - time to end simulation
303 outputs = 200 # number of time steps required.
304 h=(tend-t)/outputs #size of time step
305 #user warning statement
306 print "Expected Number of time outputs is: ", (tend-t)/h
307 i=0 #loop counter
308 #the folder to put our outputs in, leave blank "" for script path
309 save_path="data/onedheatdiff001"
310 \end{python}
311 Now that we know our inputs we will build a domain using the \verb Rectangle() function from \verb finley . The four arguments allow us to define our domain \verb rod as:
312 \begin{python}
313 #generate domain using rectangle
314 rod = Rectangle(l0=mx,l1=my,n0=ndx, n1=ndy)
315 \end{python}
316 \verb rod now describes a domain in the manner of Section \ref{ss:domcon}. As we define our variables, various function spaces will be created to accommodate them. There is an easy way to extract finite points from the domain \verb|rod| using the domain property function \verb|getX()| . This function sets the vertices of each cell as finite points to solve in the solution. If we let \verb|x| be these finite points, then;
317 \begin{python}
318 #extract finite points - the solution points
319 x=rod.getX()
320 \end{python}
321 The data locations of specific function spaces can be returned in a similar manner by extracting the relevant function space from the domain followed by the \verb .getX() operator.
322
323 With a domain and all our required variables established, it is now possible to set up our PDE so that it can be solved by \esc. The first step is to define the type of PDE that we are trying to solve in each time step. In this example it is a single linear PDE\footnote{in comparison to a system of PDEs which will be discussed later.}. We also need to state the values of our general form variables.
324 \begin{python}
325 mypde=LinearSinglePDE(rod)
326 mypde.setValue(A=kappa*kronecker(rod),D=rhocp/h)
327 \end{python}
328
329 In a few special cases it may be possible to decrease the computational time of the solver if our PDE is symmetric. Symmetry of a PDE is defined by;
330 \begin{equation}\label{eqn:symm}
331 A\hackscore{jl}=A\hackscore{lj}
332 \end{equation}
333 Symmetry is only dependent on the $A$ coefficient in the general form and the other coefficients $D$ and $d$ as well as the RHS $Y$ and $y$ may take any value. From the above definition we can see that our PDE is symmetric. The \verb LinearPDE class provides the method \method{checkSymmetry} to check if the given PDE is symmetric. As our PDE is symmetrical we will enable symmetry via;
334 \begin{python}
335 myPDE.setSymmetryOn()
336 \end{python}
337
338 We now need to specify our boundary conditions and initial values. The initial values required to solve this PDE are temperatures for each discrete point in our model. We will set our bar to:
339 \begin{python}
340 T = Tref
341 \end{python}
342 Boundary conditions are a little more difficult. Fortunately the \esc solver will handle our insulated boundary conditions by default with a zero flux operator. However, we will need to apply our heat source $q_{H}$ to the end of the bar at $x=0$ . \esc makes this easy by letting us define areas in our model. The finite points in the model were previously defined as \verb x and it is possible to set all of points that satisfy $x=0$ to \verb q via the \verb whereZero() function. There are a few \verb where functions available in \esc. They will return a value \verb 1 where they are satisfied and \verb 0 where they are not. In this case our \verb qH is only applied to the far LHS of our model as required.
343 \begin{python}
344 # ... set heat source: ....
345 qH=q*whereZero(x[0])
346 \end{python}
347
348 Finally we will initialise an iteration loop to solve our PDE for all the time steps we specified in the variable section. As the RHS of the general form is dependent on the previous values for temperature \verb T across the bar this must be updated in the loop. Our output at each time step is \verb T the heat distribution and \verb totT the total heat in the system.
349 \begin{python}
350 while t<=tend:
351 i+=1 #increment the counter
352 t+=h #increment the current time
353 mypde.setValue(Y=qH+rhocp/h*T) #set variable PDE coefficients
354 T=mypde.getSolution() #get the PDE solution
355 totT = rhocp*T #get the total heat solution in the system
356 \end{python}
357
358 \subsection{Plotting the heat solutions}
359 Visualisation of the solution can be achieved using \mpl a module contained within \pylab. We start by modifying our solution script from before. Prior to the \verb while loop we will need to extract our finite solution points to a data object that is compatible with \mpl. First it is necessary to convert \verb x to a list of tuples. These are then converted to a \numpy array and the $x$ locations extracted via an array slice to the variable \verb plx .
360 \begin{python}
361 #convert solution points for plotting
362 plx = x.toListOfTuples()
363 plx = np.array(plx) #convert to tuple to numpy array
364 plx = plx[:,0] #extract x locations
365 \end{python}
366 As there are two solution outputs, we will generate two plots and save each to a file for every time step in the solution. The following is appended to the end of the \verb while loop and creates two figures. The first figure is for the temperature distribution, and the second the total temperature in the bar. Both cases are similar with a few minor changes for scale and labelling. We start by converting the solution to a tuple and then plotting this against our \textit{x coordinates} \verb plx from before. The axis is then standardised and a title applied. Finally, the figure is saved to a *.png file and cleared for the following iteration.
367 \begin{python}
368 #establish figure 1 for temperature vs x plots
369 tempT = T.toListOfTuples(scalarastuple=False)
370 pl.figure(1) #current figure
371 pl.plot(plx,tempT) #plot solution
372 #define axis extents and title
373 pl.axis([0,1.0,273.14990+0.00008,0.004+273.1499])
374 pl.title("Temperature accross Rod")
375 #save figure to file
376 pl.savefig(os.path.join(save_path+"/tempT","rodpyplot%03d.png") %i)
377 pl.clf() #clear figure
378
379 #establish figure 2 for total temperature vs x plots and repeat
380 tottempT = totT.toListOfTuples(scalarastuple=False)
381 pl.figure(2)
382 pl.plot(plx,tottempT)
383 pl.axis([0,1.0,9.657E08,12000+9.657E08])
384 pl.title("Total temperature across Rod")
385 pl.savefig(os.path.join(save_path+"/totT","ttrodpyplot%03d.png")%i)
386 pl.clf()
387 \end{python}
388 \begin{figure}
389 \begin{center}
390 \includegraphics[width=4in]{figures/ttrodpyplot150}
391 \caption{Total temperature ($T$) distribution in rod at $t=150$}
392 \label{fig:onedheatout}
393 \end{center}
394 \end{figure}
395
396 \subsubsection{Parallel scripts (MPI)}
397 In some of the example files for this cookbook the plotting commands are a little different.
398 For example,
399 \begin{python}
400 if getMPIRankWorld() == 0:
401 pl.savefig(os.path.join(save_path+"/totT","ttrodpyplot%03d.png")%i)
402 pl.clf()
403 \end{python}
404
405 The additional \verb if statement is not necessary for normal desktop use.
406 It becomes important for scripts run on parallel computers.
407 Its purpose is to ensure that only one copy of the file is written.
408 For more details on writing scripts for parallel computing please consult the \emph{user's guide}.
409
410 \subsection{Make a video}
411 Our saved plots from the previous section can be cast into a video using the following command appended to the end of the script. \verb mencoder is linux only however, and other platform users will need to use an alternative video encoder.
412 \begin{python}
413 # compile the *.png files to create two *.avi videos that show T change
414 # with time. This operation uses linux mencoder. For other operating
415 # systems it is possible to use your favourite video compiler to
416 # convert image files to videos.
417
418 os.system("mencoder mf://"+save_path+"/tempT"+"/*.png -mf type=png:\
419 w=800:h=600:fps=25 -ovc lavc -lavcopts vcodec=mpeg4 -oac copy -o \
420 onedheatdiff001tempT.avi")
421
422 os.system("mencoder mf://"+save_path+"/totT"+"/*.png -mf type=png:\
423 w=800:h=600:fps=25 -ovc lavc -lavcopts vcodec=mpeg4 -oac copy -o \
424 onedheatdiff001totT.avi")
425 \end{python}
426

  ViewVC Help
Powered by ViewVC 1.1.26