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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2645 - (show annotations)
Thu Sep 3 02:20:33 2009 UTC (9 years, 7 months ago) by ahallam
File MIME type: application/x-tex
File size: 23618 byte(s)
More revisions to cookbook, fixed cblib.py for Scons, new figures.
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 \section{One Dimensional Heat Diffusion in an Iron Rod}
15 %\label{Sec:1DHDv0}
16 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 \ESCRIPT and 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 \ESCRIPT and implementation should become a clearer as we develop our understanding further into the cookbook.}
17
18 \begin{figure}[h!]
19 \centerline{\includegraphics[width=4.in]{figures/onedheatdiff}}
20 \caption{One dimensional model of an Iron bar.}
21 \label{fig:onedhdmodel}
22 \end{figure}
23 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 the heat source.
24
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 conductivity constant for a given material\footnote{A list of some common thermal conductivities is available from Wikipedia \url{http://en.wikipedia.org/wiki/List_of_thermal_conductivities}}.
33 The heat source is defined by the right hand side of \ref{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} = Te^{-\gamma t}$ where we have the output of our heat source decaying with time. There are also two partial derivatives in \ref{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$.
34
35 \subsection{Escript, PDEs and The General Form}
36 Potentially, it is now possible to solve \ref{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 when a large number of sums or a more complex visualisation is required. To do this, a numerical approach is required - \ESCRIPT can help us here - and it becomes necessary to discretize the equation so that 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 modeller.
37
38 \ESCRIPT interfaces with any given PDE via a general form. In this example we will illustrate a simpler version of the full linear PDE general form which is available in the \ESCRIPT user's guide. A simplified form that suits our heat diffusion problem\footnote{In the form of the \ESCRIPT users guide which using the Einstein convention is written as
39 $-(A\hackscore{jl} u\hackscore{,l})\hackscore{,j}+D u =Y$}
40 is described by;
41 \begin{equation}\label{eqn:commonform nabla}
42 -\nabla.(A.\nabla u) + Du = f
43 \end{equation}
44 where $A$, $D$ and $f$ are known values. The symbol $\nabla$ which is called the \textit{Nabla operator} or \textit{del operator} represents
45 the spatial derivative of its subject - in this case $u$. Lets assume for a moment that we deal with a one-dimensional problem then ;
46 \begin{equation}
47 \nabla = \frac{\partial}{\partial x}
48 \end{equation}
49 and we can write \ref{eqn:commonform nabla} as;
50 \begin{equation}\label{eqn:commonform}
51 -A\frac{\partial^{2}u}{\partial x^{2}} + Du = f
52 \end{equation}
53 if $A$ is constant then \ref{eqn:commonform} is consistent with our heat diffusion problem in \ref{eqn:hd} with the exception of $u$. Thus when comparing equations \ref{eqn:hd} and \ref{eqn:commonform} we see that;
54 \begin{equation}
55 A = \kappa; D = \rho c \hackscore{p}; f = q \hackscore{H}
56 \end{equation}
57
58 We can write the partial $\frac{\partial T}{\partial t}$ in terms of $u$ by discretising the time of our solution. The method we will use is the Backwards Euler approximation, which states;
59 \begin{equation}
60 \frac{\partial f(x)}{\partial x} \approx \frac{f(x+h)-f(x)}{h}
61 \label{eqn:beuler}
62 \end{equation}
63 where h is the the discrete step size $\Delta x$.
64 Now if the temperature $T(t)$ is substituted in by letting $f(x) = T(t)$ then from \ref{eqn:beuler} we see that;
65 \begin{equation}
66 T'(t) \approx \frac{T(t+h) - T(t)}{h}
67 \end{equation}
68 which can also be written as;
69 \begin{equation}
70 \frac{\partial T^{(n)}}{\partial t} \approx \frac{T^{(n)} - T^{(n-1)}}{h}
71 \label{eqn:Tbeuler}
72 \end{equation}
73 where $n$ denotes the n\textsuperscript{th} time step. Substituting \ref{eqn:Tbeuler} into \ref{eqn:hd} we get;
74 \begin{equation}
75 \frac{\rho c\hackscore p}{h} (T^{(n)} - T^{(n-1)}) - \kappa \frac{\partial^{2} T}{\partial x^{2}} = q\hackscore H
76 \label{eqn:hddisc}
77 \end{equation}
78 To fit our simplified general form we can rearrange \ref{eqn:hddisc};
79 \begin{equation}
80 \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)}
81 \label{eqn:hdgenf}
82 \end{equation}
83 The PDE is now in a form that statisfies \refEq{eqn:commonform nabla} which is required for \ESCRIPT to solve our PDE. This can be done by generating a solution for sucessive increments in the time nodes $t^{(n)}$ where
84 $t^{(0)}=0$ and $t^{(n)}=t^{(n-1)}+h$ where $h>0$ is the step size and assumed to be constant.
85 In the following the upper index ${(n)}$ refers to a value at time $t^{(n)}$. Finally, by comparing \ref{eqn:hdgenf} with \ref{eqn:commonform} it can be seen that;
86 \begin{equation}
87 A = \kappa; D = \frac{\rho c \hackscore{p}}{h}; f = q \hackscore{H} + \frac{\rho c\hackscore p}{h} T^{(n-1)}
88 \end{equation}
89
90 Now that the general form has been established, it can be submitted to \ESCRIPT. 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 \ref{eqn:Tbeuler}. Our model stipulates a starting temperature in the iron bar of 0\textcelsius. Thus the temperature distribution is simply;
91 \begin{equation}
92 T(x,0) = T\hackscore{ref} = 0
93 \end{equation}
94 for all $x$ in the domain.
95
96 \subsection{Boundary Conditions}
97 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 Neumann and Dirichlet\footnote{More information on Boundary Conditions is available at Wikipedia \url{http://en.wikipedia.org/wiki/Boundary_conditions}}. In this example, we have utilised both conditions. Dirichlet is conceptually simpler and is used to prescribe a known value to the model on its boundary. This is like holding a snake by the tail; we know where the tail will be as we hold it however, we have no control over the rest of the snake. Dirichlet boundary conditions exist where we have applied our heat source. As the heat source is a constant, we can simulate its presence on that boundary. This is done by continuously resetting the temperature of the boundary, so that is is the same as the heat source.
98
99 Neumann boundary conditions describe the radiation or flux normal to the boundayor surface. This aptly describes our insulation conditions as we do not want to exert a constant temperature as with the heat source. However, we do want to prevent any loss of energy from the system. These natural boundary conditions can be described by specifying a radiation condition which prescribes the normal component of the flux $\kappa T\hackscore{,i}$ to be proportional
100 to the difference of the current temperature to the surrounding temperature $T\hackscore{ref}$; in general terms this is;
101 \begin{equation}
102 \kappa T\hackscore{,i} n\hackscore i = \eta (T\hackscore{ref}-T)
103 \label{eqn:hdbc}
104 \end{equation}
105 and simplified to our one dimensional model we have;
106 \begin{equation}
107 \kappa \frac{\partial T}{\partial dx} n\hackscore x = \eta (T\hackscore{ref}-T)
108 \end{equation}
109 where $\eta$ is a given material coefficient depending on the material of the block and the surrounding medium and $n\hackscore i$ is the $i$-th component of the outer normal field \index{outer normal field} at the surface of the domain. These two conditions form a boundary value problem that has to be solved for each time step. Due to the perfect insulation in our model we can set $\eta = 0$ which results in zero flux - no energy in or out - we do not need to worry about the Neumann terms of the general form for this example.
110
111 \subsection{A \textit{1D} Clarification}
112 It is necessary for clarification that we revisit the general PDE from \refeq{eqn:commonform nabla} under the light of a two dimensional domain. \ESCRIPT is inherently designed to solve problems that are greater than one dimension and so \ref{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, \ref{eqn:commonform nabla} assuming a constant coefficient $A$, takes the form;
113 \begin{equation}\label{eqn:commonform2D}
114 -A\hackscore{00}\frac{\partial^{2}u}{\partial x^{2}}
115 -A\hackscore{01}\frac{\partial^{2}u}{\partial x\partial y}
116 -A\hackscore{10}\frac{\partial^{2}u}{\partial y\partial x}
117 -A\hackscore{11}\frac{\partial^{2}u}{\partial y^{2}}
118 + Du = f
119 \end{equation}
120 Notice that for the higher dimensional case $A$ becomes a matrix. It is also
121 important to notice that the usage of the Nabla operator creates
122 a compact formulation which is also independent from the spatial dimension.
123 So to make the general PDE~\ref{eqn:commonform2D} one dimensional as
124 shown in~\ref{eqn:commonform} we need to set
125 \begin{equation}
126 A\hackscore{00}=A; A\hackscore{01}=A\hackscore{10}=A\hackscore{11}=0
127 \end{equation}
128
129 \subsection{Developing a PDE Solution Script}
130 To solve \ref{eqn:hd} we will write a simple python script which uses the \modescript, \modfinley and \modmpl modules. At this point we assume that you have some basic understanding of the python programming language. If not there are some pointers and links available in Section \ref{sec:escpybas} .
131
132 Our goal here is to develop a script for \ESCRIPT that will solve the heat equation at successive time steps for a predefined period using our general form \ref{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 sin and cos functions or more complicated like those from our \ESCRIPT library.}
133 that we will require.
134 \begin{verbatim}
135 from esys.escript import *
136 # This defines the LinearPDE module as LinearPDE
137 from esys.escript.linearPDEs import LinearPDE
138 # This imports the rectangle domain function from finley.
139 from esys.finley import Rectangle
140 # A useful unit handling package which will make sure all our units
141 # match up in the equations under SI.
142 from esys.escript.unitsSI import *
143 import pylab as pl #Plotting package.
144 import numpy as np #Array package.
145 import os #This package is necessary to handle saving our data.
146 \end{verbatim}
147 It is generally a good idea to import all of the \modescript library, although if you know the packages you need you can specify them individually. The function \verb|LinearPDE| has been imported for ease of use later in the script. \verb|Rectangle| is going to be our type of domain. The package \verb unitsSI is a module of \esc that provides support for units definitions with our variables; and the \verb|os| package 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.
148
149 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 \ESCRIPT 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 domain upon which we wish to solve our problem needs to be defined. There are many different types of domains in \modescript which we will demonstrate in later tutorials but for our iron rod, we will simply use a rectangular domain.
150
151 Using a rectangular domain 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 domain, 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 domain 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 our \textit{1D} problem we will define our bar as being 1 metre long. An appropriate \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 so the temperature does not vary with $y$. Thus the domain 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.
152 \begin{verbatim}
153 #Domain related.
154 mx = 1*m #meters - model length
155 my = .1*m #meters - model width
156 ndx = 100 # mesh steps in x direction
157 ndy = 1 # mesh steps in y direction - one dimension means one element
158 \end{verbatim}
159 The material constants and the temperature variables must also be defined. For the iron rod in the model they are defined as:
160 \begin{verbatim}
161 #PDE related
162 q=200. * Celsius #Kelvin - our heat source temperature
163 Tref = 0. * Celsius #Kelvin - starting temp of iron bar
164 rho = 7874. *kg/m**3 #kg/m^{3} density of iron
165 cp = 449.*J/(kg*K) #j/Kg.K thermal capacity
166 rhocp = rho*cp
167 kappa = 80.*W/m/K #watts/m.Kthermal conductivity
168 \end{verbatim}
169 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:
170 \begin{verbatim}
171 t=0 #our start time, usually zero
172 tend=5.*minute #seconds - time to end simulation
173 outputs = 200 # number of time steps required.
174 h=(tend-t)/outputs #size of time step
175 #user warning statement
176 print "Expected Number of time outputs is: ", (tend-t)/h
177 i=0 #loop counter
178 #the folder to put our outputs in, leave blank "" for script path
179 save_path="data/onedheatdiff001"
180 \end{verbatim}
181 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:
182 \begin{verbatim}
183 #generate domain using rectangle
184 rod = Rectangle(l0=mx,l1=my,n0=ndx, n1=ndy)
185 \end{verbatim}
186 In this form \verb rod does not represent any discrete points, but rather an area of \verb ndx*ndy cells that fit into a rectangular space with opposing vertices at the origin and the point \verb [mx,my] . Our domain is constructed this way to allow the user to determine where the discrete points of each model will be located. Discretisation may be at the corners of each cell, the middle point of a cell, halfway along each side of the cell and so on. Depending on the PDE or the model there may be advantages and disadvantages for each case. Fortunately \modescript offers 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;
187 \begin{verbatim}
188 #extract finite points - the solution points
189 x=rod.getX()
190 \end{verbatim}
191 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 \ESCRIPT. 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.
192 \begin{verbatim}
193 mypde=LinearSinglePDE(rod)
194 mypde.setValue(A=kappa*kronecker(rod),D=rhocp/h)
195 \end{verbatim}
196
197 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;
198 \begin{equation}\label{eqn:symm}
199 A\hackscore{jl}=A\hackscore{lj}
200 \end{equation}
201 Symmetry is only dependent on the $A$ coefficient in the general form and the others $D$ and $d$ as well as the RHS coefficients $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;
202 \begin{verbatim}
203 myPDE.setSymmetryOn()
204 \end{verbatim}
205
206 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 domain. We will set our bar to:
207 \begin{verbatim}
208 T = Tref
209 \end{verbatim}
210 Boundary conditions are a little more difficult. Fortunately the escript 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$ . \ESCRIPT makes this easy by letting us define areas in our domain. The finite points in the domain 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 \ESCRIPT. 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.
211 \begin{verbatim}
212 # ... set heat source: ....
213 qH=q*whereZero(x[0])
214 \end{verbatim}
215
216 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 timestep is \verb T the heat distribution and \verb totT the total heat in the system.
217 \begin{verbatim}
218 while t<=tend:
219 i+=1 #increment the counter
220 t+=h #increment the current time
221 mypde.setValue(Y=qH+rhocp/h*T) #set variable PDE coefficients
222 T=mypde.getSolution() #get the PDE solution
223 totT = rhocp*T #get the total heat solution in the system
224 \end{verbatim}
225
226 \subsection{Plotting the heat solutions}
227 Visualisation of the solution can be achieved using \mpl a module contained with \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 .
228 \begin{verbatim}
229 #convert solution points for plotting
230 plx = x.toListOfTuples()
231 plx = np.array(plx) #convert to tuple to numpy array
232 plx = plx[:,0] #extract x locations
233 \end{verbatim}
234 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 labeling. 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. The figure is then saved to a *.png file and cleared for the following iteration.
235 \begin{verbatim}
236 #establish figure 1 for temperature vs x plots
237 tempT = T.toListOfTuples(scalarastuple=False)
238 pl.figure(1) #current figure
239 pl.plot(plx,tempT) #plot solution
240 #define axis extents and title
241 pl.axis([0,1.0,273.14990+0.00008,0.004+273.1499])
242 pl.title("Temperature accross Rod")
243 #save figure to file
244 pl.savefig(os.path.join(save_path+"/tempT","rodpyplot%03d.png") %i)
245 pl.clf() #clear figure
246
247 #establish figure 2 for total temperature vs x plots and repeat
248 tottempT = totT.toListOfTuples(scalarastuple=False)
249 pl.figure(2)
250 pl.plot(plx,tottempT)
251 pl.axis([0,1.0,9.657E08,12000+9.657E08])
252 pl.title("Total temperature accross Rod")
253 pl.savefig(os.path.join(save_path+"/totT","ttrodpyplot%03d.png")%i)
254 pl.clf()
255 \end{verbatim}
256 \begin{figure}
257 \begin{center}
258 \includegraphics[width=4in]{figures/ttrodpyplot150}
259 \caption{Total temperature ($T$) distribution in rod at $t=150$}
260 \label{fig:onedheatout}
261 \end{center}
262 \end{figure}
263
264 \subsection{Make a video}
265 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.
266 \begin{verbatim}
267 # compile the *.png files to create two *.avi videos that show T change
268 # with time. This opperation uses linux mencoder. For other operating
269 # systems it is possible to use your favourite video compiler to
270 # convert image files to videos.
271
272 os.system("mencoder mf://"+save_path+"/tempT"+"/*.png -mf type=png:\
273 w=800:h=600:fps=25 -ovc lavc -lavcopts vcodec=mpeg4 -oac copy -o \
274 onedheatdiff001tempT.avi")
275
276 os.system("mencoder mf://"+save_path+"/totT"+"/*.png -mf type=png:\
277 w=800:h=600:fps=25 -ovc lavc -lavcopts vcodec=mpeg4 -oac copy -o \
278 onedheatdiff001totT.avi")
279 \end{verbatim}
280

  ViewVC Help
Powered by ViewVC 1.1.26