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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2957 - (hide annotations)
Tue Mar 2 01:28:19 2010 UTC (9 years, 6 months ago) by artak
File MIME type: application/x-tex
File size: 37823 byte(s)
more corrections
1 ahallam 2401
2     %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
3     %
4 jfenwick 2881 % Copyright (c) 2003-2010 by University of Queensland
5 ahallam 2401 % 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 gross 2905 \begin{figure}[h!]
15     \centerline{\includegraphics[width=4.in]{figures/onedheatdiff001}}
16 gross 2952 \caption{Example 1: Temperature differential along a single interface between two granite blocks.}
17 gross 2905 \label{fig:onedgbmodel}
18     \end{figure}
19 ahallam 2801
20 gross 2949 \section{Example 1: One Dimensional Heat Diffusion in Granite}
21 gross 2905 \label{Sec:1DHDv00}
22 gross 2878
23 gross 2905 The first model consists of two blocks of isotropic material, for instance granite, sitting next to each other.
24 artak 2957 Initial temperature in \textit{Block 1} is \verb|T1| and in \textit{Block 2} is \verb|T2|.
25 gross 2905 We assume that the system is insulated.
26     What would happen to the temperature distribution in each block over time?
27 artak 2957 Intuition tells us that heat will be transported from the hotter block to the cooler until both
28 gross 2905 blocks have the same temperature.
29    
30 ahallam 2495 \subsection{1D Heat Diffusion Equation}
31 artak 2957 We can model the heat distribution of this problem over time using 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}};
32 ahallam 2494 which is defined as:
33 ahallam 2401 \begin{equation}
34     \rho c\hackscore p \frac{\partial T}{\partial t} - \kappa \frac{\partial^{2} T}{\partial x^{2}} = q\hackscore H
35     \label{eqn:hd}
36     \end{equation}
37 gross 2861 where $\rho$ is the material density, $c\hackscore p$ is the specific heat and $\kappa$ is the thermal
38     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
39     parameters are \textbf{constant}.
40     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$.
41 ahallam 2401
42 gross 2861 \subsection{PDEs and the General Form}
43     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
44 gross 2905 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 modeller.
45 gross 2477
46 gross 2861 Firstly, we will discretise the PDE \refEq{eqn:hd} in the time direction which will
47     leave as with a steady linear PDE which is involving spatial derivatives only and needs to be solved in each time
48     step to progress in time - \esc can help us here.
49    
50     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
51     approximation
52     \begin{equation}
53     \frac{\partial T(t)}{\partial t} \approx \frac{T(t)-T(t-h)}{h}
54     \label{eqn:beuler}
55     \end{equation}
56     for $\frac{\partial T}{\partial t}$ at time $t$
57     where $h$ is the time step size. This can also be written as;
58     \begin{equation}
59     \frac{\partial T}{\partial t}(t^{(n)}) \approx \frac{T^{(n)} - T^{(n-1)}}{h}
60     \label{eqn:Tbeuler}
61     \end{equation}
62     where the upper index $n$ denotes the n\textsuperscript{th} time step. So one has
63     \begin{equation}
64     \begin{array}{rcl}
65     t^{(n)} & = & t^{(n-1)}+h \\
66     T^{(n)} & = & T(t^{(n-1)}) \\
67     \end{array}
68     \label{eqn:Neuler}
69     \end{equation}
70     Substituting \refEq{eqn:Tbeuler} into \refEq{eqn:hd} we get;
71     \begin{equation}
72     \frac{\rho c\hackscore p}{h} (T^{(n)} - T^{(n-1)}) - \kappa \frac{\partial^{2} T^{(n)}}{\partial x^{2}} = q\hackscore H
73     \label{eqn:hddisc}
74     \end{equation}
75     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
76     approach is called the \textbf{forward Euler} scheme. This scheme can provide some computational advantages which
77     we are not discussed here but has the major disadvantage that depending on the
78 gross 2905 material parameter as well as the discretization 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
79     that the approximation of the temperature will not grow beyond its initial bounds and becomes non-physical.
80 gross 2861 The backward Euler which we use here is unconditionally stable meaning that under the assumption of
81     physically correct problem set-up the temperature approximation remains physical for all times.
82     The user needs to keep in mind that the discretization error introduced by \refEq{eqn:beuler}
83     is sufficiently small so a good approximation of the true temperature is calculated. It is
84     therefore crucial that the user remains critical about his/her results and for instance compares
85     the results for different time and spatial step sizes.
86    
87     To get the temperature $T^{(n)}$ at time $t^{(n)}$ we need to solve the linear
88     differential equation \refEq{eqn:hddisc} which is only including spatial derivatives. To solve this problem
89     we want to to use \esc.
90    
91     \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
92 gross 2477 $-(A\hackscore{jl} u\hackscore{,l})\hackscore{,j}+D u =Y$}
93 ahallam 2606 is described by;
94 gross 2477 \begin{equation}\label{eqn:commonform nabla}
95 jfenwick 2657 -\nabla\cdot(A\cdot\nabla u) + Du = f
96 ahallam 2411 \end{equation}
97 gross 2861 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
98 ahallam 2495 the spatial derivative of its subject - in this case $u$. Lets assume for a moment that we deal with a one-dimensional problem then ;
99 gross 2477 \begin{equation}
100     \nabla = \frac{\partial}{\partial x}
101     \end{equation}
102 gross 2861 and we can write \refEq{eqn:commonform nabla} as;
103 ahallam 2411 \begin{equation}\label{eqn:commonform}
104 gross 2477 -A\frac{\partial^{2}u}{\partial x^{2}} + Du = f
105 ahallam 2411 \end{equation}
106 gross 2861 if $A$ is constant. To match this simplified general form to our problem \refEq{eqn:hddisc}
107     we rearrange \refEq{eqn:hddisc};
108 ahallam 2411 \begin{equation}
109 ahallam 2645 \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)}
110 ahallam 2494 \label{eqn:hdgenf}
111     \end{equation}
112 ahallam 2775 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
113 ahallam 2495 $t^{(0)}=0$ and $t^{(n)}=t^{(n-1)}+h$ where $h>0$ is the step size and assumed to be constant.
114 gross 2861 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;
115 gross 2862 \begin{equation}\label{ESCRIPT SET}
116 gross 2861 u=T^{(n)};
117 ahallam 2494 A = \kappa; D = \frac{\rho c \hackscore{p}}{h}; f = q \hackscore{H} + \frac{\rho c\hackscore p}{h} T^{(n-1)}
118     \end{equation}
119    
120 ahallam 2495 \subsection{Boundary Conditions}
121 gross 2878 \label{SEC BOUNDARY COND}
122 gross 2862 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 gross 2905 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 the boundary or on the entire boundary of the region of interest.
124     We discuss Dirichlet boundary condition in our second example presented in Section~\ref{Sec:1DHDv0}.
125 ahallam 2495
126 gross 2905 We make the model assumption that the system is insulated so we need
127     to add an appropriate boundary condition to prevent
128     any loss or inflow of energy at boundary of our domain. Mathematically this is expressed by prescribing
129     the heat flux $\kappa \frac{\partial T}{\partial x}$ to zero. In our simplified one dimensional model this is expressed
130 gross 2862 in the form;
131 ahallam 2494 \begin{equation}
132 gross 2862 \kappa \frac{\partial T}{\partial x} = 0
133 ahallam 2494 \end{equation}
134 gross 2862 or in a more general case as
135     \begin{equation}\label{NEUMAN 1}
136     \kappa \nabla T \cdot n = 0
137     \end{equation}
138     where $n$ is the outer normal field \index{outer normal field} at the surface of the domain.
139 gross 2905 The $\cdot$ (dot) refers to the dot product of the vectors $\nabla T$ and $n$. In fact, the term $\nabla T \cdot n$ is the normal derivative of
140 gross 2862 the temperature $T$. Other notations which are used are\footnote{The \esc notation for the normal
141     derivative is $T\hackscore{,i} n\hackscore i$.};
142 ahallam 2645 \begin{equation}
143 gross 2862 \nabla T \cdot n = \frac{\partial T}{\partial n} \; .
144 ahallam 2645 \end{equation}
145 gross 2862 A condition of the type \refEq{NEUMAN 1} defines a \textbf{Neuman boundary condition} for the PDE.
146 ahallam 2494
147 gross 2905 The PDE \refEq{eqn:hdgenf}
148     and the Neuman boundary condition~\ref{eqn:hdgenf} (potentially together with the Dirichlet boundary condition set) define a \textbf{boundary value problem}.
149 gross 2862 It is a nature of a boundary value problem that it allows to make statements on the solution in the
150     interior of the domain from information known on the boundary only. In most cases
151     we use the term partial differential equation but in fact mean a boundary value problem.
152     It is important to keep in mind that boundary conditions need to be complete and consistent in the sense that
153     at any point on the boundary either a Dirichlet or a Neuman boundary condition must be set.
154    
155 gross 2878 Conveniently, \esc makes default assumption on the boundary conditions which the user may modify where appropriate.
156 gross 2862 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
157     $n\hackscore{j}A\hackscore{jl} u\hackscore{,l}=0$.} is;
158     \begin{equation}\label{NEUMAN 2}
159 gross 2948 -n\cdot A \cdot\nabla u = 0
160 gross 2862 \end{equation}
161     which is used everywhere on the boundary. Again $n$ denotes the outer normal field.
162     Notice that the coefficient $A$ is the same as in the \esc PDE~\ref{eqn:commonform nabla}.
163     With the settings for the coefficients we have already identified in \refEq{ESCRIPT SET} this
164     condition translates into
165 gross 2867 \begin{equation}\label{NEUMAN 2b}
166 gross 2862 \kappa \frac{\partial T}{\partial x} = 0
167     \end{equation}
168 gross 2905 for the boundary of the domain. 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.
169 gross 2862
170 gross 2870 \subsection{Outline of the Implementation}
171     \label{sec:outline}
172     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
173 gross 2905 calculate the temperature. For our problem this is the joint blocks of granite which has a rectangular shape. Secondly we need to define the PDE
174 gross 2878 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. As a work flow this takes the form;
175     \begin{enumerate}
176     \item create domain
177     \item create PDE
178     \item while end time not reached:
179     \begin{enumerate}
180     \item set PDE coefficients
181     \item solve PDE
182     \item update time marker
183     \end{enumerate}
184     \item end of calculation
185     \end{enumerate}
186 gross 2870 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
187 gross 2905 rather than its actual representation. So we will create a domain object to describe the geometry of the two
188     granite blocks. The main feature
189 gross 2870 of the object we will use is the fact that we can define PDEs and spatially distributed values such as the temperature
190     on a domain. In fact the domain object has many more features - most of them you will
191     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
192 gross 2878 later stage you may use more advanced features of the PDE class but you need to worry about them only at the point when you use them.
193 gross 2870
194    
195     \begin{figure}[t]
196     \centering
197     \includegraphics[width=6in]{figures/functionspace.pdf}
198     \label{fig:fs}
199     \caption{\esc domain construction overview}
200     \end{figure}
201    
202     \subsection{The Domain Constructor in \esc}
203     \label{ss:domcon}
204     It is helpful to have a better understanding how spatially distributed value such as the temperature or PDE coefficients are interpreted in \esc. Again
205     from the user's point of view the representation of these spatially distributed values is not relevant.
206    
207     There are various ways to construct domain objects. The simplest form is as rectangular shaped region with a length and height. There is
208     a ready to use function call for this. Besides the spatial dimensions the function call will require you to specify the number
209 gross 2905 elements or cells to be used along the length and height, see \reffig{fig:fs}. Any spatially distributed value
210     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} for details.}. Therefore we will have access to an approximation of the true PDE solution only.
211 gross 2870 The quality of the approximation depends - besides other factors- mainly on the number of elements being used. In fact, the
212     approximation becomes better the more elements are used. However, computational costs and compute time grow with the number of
213 gross 2878 elements being used. It therefore important that you find the right balance between the demand in accuracy and acceptable resource usage.
214 gross 2870
215     In general, one can thinks about a domain object as a composition of nodes and elements.
216 gross 2905 As shown in \reffig{fig:fs}, an element is defined by the nodes used to describe its vertices.
217 gross 2870 To represent spatial distributed values the user can use
218     the values at the nodes, at the elements in the interior of the domain or at elements located at the surface of the domain.
219     The different approach used to represent values is called \textbf{function space} and is attached to all objects
220     in \esc representing a spatial distributed value such as the solution of a PDE. The three
221     function spaces we will use at the moment are;
222     \begin{enumerate}
223     \item the nodes, called by \verb|ContinuousFunction(domain)| ;
224     \item the elements/cells, called by \verb|Function(domain)| ; and
225     \item the boundary, called by \verb|FunctionOnBoundary(domain)| .
226     \end{enumerate}
227     A function space object such as \verb|ContinuousFunction(domain)| has the method \verb|getX| attached to it. This method returns the
228     location of the so-called \textbf{sample points} used to represent values with the particular function space attached to it. So the
229     call \verb|ContinuousFunction(domain).getX()| will return the coordinates of the nodes used to describe the domain while
230     the \verb|Function(domain).getX()| returns the coordinates of numerical integration points within elements, see
231 gross 2905 \reffig{fig:fs}.
232 gross 2870
233 gross 2878 This distinction between different representations of spatial distributed values
234 gross 2870 is important in order to be able to vary the degrees of smoothness in a PDE problem.
235     The coefficients of a PDE need not be continuous thus this qualifies as a \verb|Function()| type.
236     On the other hand a temperature distribution must be continuous and needs to be represented with a \verb|ContinuousFunction()| function space.
237     An influx may only be defined at the boundary and is therefore a \verb FunctionOnBoundary() object.
238     \esc allows certain transformations of the function spaces. A \verb ContinuousFunction() can be transformed into a \verb|FunctionOnBoundary()|
239     or \verb|Function()|. On the other hand there is not enough information in a \verb FunctionOnBoundary() to transform it to a \verb ContinuousFunction() .
240 gross 2878 These transformations, which are called \textbf{interpolation} are invoked automatically by \esc if needed.
241 gross 2870
242     Later in this introduction we will discuss how
243     to define specific areas of geometry with different materials which are represented by different material coefficients such the
244 gross 2878 thermal conductivities $kappa$. A very powerful technique to define these types of PDE
245 gross 2870 coefficients is tagging. Blocks of materials and boundaries can be named and values can be defined on subregions based on their names.
246 gross 2878 This is simplifying PDE coefficient and flux definitions. It makes for much easier scripting. We will discuss this technique in Section~\ref{STEADY-STATE HEAT REFRACTION}.
247 gross 2870
248    
249     \subsection{A Clarification for the 1D Case}
250 gross 2931 \label{SEC: 1D CLARIFICATION}
251 gross 2861 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;
252 gross 2477 \begin{equation}\label{eqn:commonform2D}
253     -A\hackscore{00}\frac{\partial^{2}u}{\partial x^{2}}
254     -A\hackscore{01}\frac{\partial^{2}u}{\partial x\partial y}
255     -A\hackscore{10}\frac{\partial^{2}u}{\partial y\partial x}
256     -A\hackscore{11}\frac{\partial^{2}u}{\partial y^{2}}
257     + Du = f
258     \end{equation}
259 ahallam 2606 Notice that for the higher dimensional case $A$ becomes a matrix. It is also
260 ahallam 2495 important to notice that the usage of the Nabla operator creates
261     a compact formulation which is also independent from the spatial dimension.
262 gross 2861 So to make the general PDE \refEq{eqn:commonform2D} one dimensional as
263     shown in \refEq{eqn:commonform} we need to set
264 ahallam 2606 \begin{equation}
265 ahallam 2494 A\hackscore{00}=A; A\hackscore{01}=A\hackscore{10}=A\hackscore{11}=0
266 gross 2477 \end{equation}
267    
268 gross 2867
269 ahallam 2495 \subsection{Developing a PDE Solution Script}
270 ahallam 2801 \label{sec:key}
271 gross 2949 \sslist{example01a.py}
272 gross 2870 We will write a simple \pyt script which uses the \modescript, \modfinley and \modmpl modules.
273 gross 2905 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.}
274 ahallam 2495 that we will require.
275 ahallam 2775 \begin{python}
276 ahallam 2495 from esys.escript import *
277 ahallam 2606 # This defines the LinearPDE module as LinearPDE
278     from esys.escript.linearPDEs import LinearPDE
279     # This imports the rectangle domain function from finley.
280     from esys.finley import Rectangle
281     # A useful unit handling package which will make sure all our units
282     # match up in the equations under SI.
283     from esys.escript.unitsSI import *
284 ahallam 2775 \end{python}
285 gross 2878 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.
286 gross 2477
287 gross 2905 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 granite blocks, we will simply use a rectangular model.
288 ahallam 2401
289 gross 2905 Using a rectangular model simplifies our granite blocks which would in reality be a \textit{3D} object, into a single dimension. The granite blocks will have a lengthways cross section that looks like a rectangle. As a result we do not need to model the volume of the block. 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.
290 ahallam 2775 \begin{python}
291 gross 2905 mx = 500.*m #meters - model length
292     my = 100.*m #meters - model width
293     ndx = 50 # mesh steps in x direction
294     ndy = 1 # mesh steps in y direction
295     boundloc = mx/2 # location of boundary between the two blocks
296 ahallam 2775 \end{python}
297 gross 2905 The material constants and the temperature variables must also be defined. For the granite in the model they are defined as:
298 ahallam 2775 \begin{python}
299 ahallam 2495 #PDE related
300 gross 2905 rho = 2750. *kg/m**3 #kg/m^{3} density of iron
301     cp = 790.*J/(kg*K) # J/Kg.K thermal capacity
302 gross 2878 rhocp = rho*cp
303 gross 2905 kappa = 2.2*W/m/K # watts/m.Kthermal conductivity
304 gross 2878 qH=0 * J/(sec*m**3) # J/(sec.m^{3}) no heat source
305 gross 2905 T1=20 * Celsius # initial temperature at Block 1
306     T2=2273. * Celsius # base temperature at Block 2
307 ahallam 2775 \end{python}
308 ahallam 2495 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:
309 ahallam 2775 \begin{python}
310 gross 2878 t=0 * day #our start time, usually zero
311     tend=1. * day # - time to end simulation
312 ahallam 2495 outputs = 200 # number of time steps required.
313     h=(tend-t)/outputs #size of time step
314 ahallam 2606 #user warning statement
315     print "Expected Number of time outputs is: ", (tend-t)/h
316     i=0 #loop counter
317 ahallam 2775 \end{python}
318 gross 2905 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 model as:
319 ahallam 2775 \begin{python}
320 ahallam 2606 #generate domain using rectangle
321 gross 2905 blocks = Rectangle(l0=mx,l1=my,n0=ndx, n1=ndy)
322 ahallam 2775 \end{python}
323 gross 2905 \verb blocks now describes a domain in the manner of Section \ref{ss:domcon}. T
324 ahallam 2658
325 ahallam 2775 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.
326     \begin{python}
327 gross 2905 mypde=LinearPDE(blocks)
328 gross 2878 A=zeros((2,2)))
329     A[0,0]=kappa
330 gross 2905 mypde.setValue(A=A, D=rhocp/h)
331 ahallam 2775 \end{python}
332 gross 2878 In a many cases it may be possible to decrease the computational time of the solver if the PDE is symmetric.
333     Symmetry of a PDE is defined by;
334 ahallam 2495 \begin{equation}\label{eqn:symm}
335     A\hackscore{jl}=A\hackscore{lj}
336     \end{equation}
337 gross 2878 Symmetry is only dependent on the $A$ coefficient in the general form and the other coefficients $D$ as well as the right hand side $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;
338 ahallam 2775 \begin{python}
339 ahallam 2495 myPDE.setSymmetryOn()
340 ahallam 2775 \end{python}
341 gross 2905 Next we need to establish the initial temperature distribution \verb|T|. We need to
342     assign the value \verb|T1| to all sample points left to the contact interface at $x\hackscore{0}=\frac{mx}{2}$
343     and the value \verb|T2| right to the contact interface. \esc
344     provides the \verb|whereNegative| function to construct this. In fact,
345     \verb|whereNegative| returns the value $1$ at those sample points where the argument
346     has a negative value. Otherwise zero is returned. If \verb|x| are the $x\hackscore{0}$
347     coordinates of the sample points used to represent the temperature distribution
348     then \verb|x[0]-boundloc| gives us a negative value for
349     all sample points left to the interface and non-negative value to
350     the right of the interface. So with;
351 gross 2878 \begin{python}
352     # ... set initial temperature ....
353 gross 2905 T= T1*whereNegative(x[0]-boundloc)+T2*(1-whereNegative(x[0]-boundloc))
354 gross 2878 \end{python}
355 gross 2905 we get the desired temperature distribution. To get the actual sample points \verb|x| we use
356     the \verb|getX()| method of the function space \verb|Solution(blocks)|
357     which is used to represent the solution of a PDE;
358 gross 2878 \begin{python}
359 gross 2905 x=Solution(blocks).getX()
360     \end{python}
361     As \verb|x| are the sample points for the function space \verb|Solution(blocks)|
362     the initial temperature \verb|T| is using these sample points for representation.
363     Although \esc is trying to be forgiving with the choice of sample points and to convert
364     where necessary the adjustment of the function space is not always possible. So it is
365     advisable to make a careful choice on the function space used.
366    
367     Finally we will initialise an iteration loop to solve our PDE for all the time steps we specified in the variable section. As the right hand side 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.
368     \begin{python}
369 gross 2878 while t < tend:
370     i+=1 #increment the counter
371     t+=h #increment the current time
372     mypde.setValue(Y=qH+rhocp/h*T) # set variable PDE coefficients
373     T=mypde.getSolution() #get the PDE solution
374     totE = integrate(rhocp*T) #get the total heat (energy) in the system
375     \end{python}
376     The last statement in this script calculates the total energy in the system as volume integral
377 gross 2905 of $\rho c\hackscore{p} T$ over the block. As the blocks are insulated no energy should be get lost or added.
378     The total energy should stay constant for the example discussed here.
379 ahallam 2401
380 gross 2905 \subsection{Running the Script}
381     The script presented so for is available under
382 gross 2949 \verb|example01a.py|. You can edit this file with your favourite text editor.
383 jfenwick 2911 On most operating systems\footnote{The you can use \texttt{run-escript} launcher is not supported under {\it MS Windows} yet.} you can use the \program{run-escript} command
384 gross 2905 to launch {\it escript} scripts. For the example script use;
385     \begin{verbatim}
386 gross 2949 run-escript example01a.py
387 gross 2905 \end{verbatim}
388     The program will print a progress report. Alternatively, you can use
389     the python interpreter directly;
390     \begin{verbatim}
391 gross 2949 python example01a.py
392 gross 2905 \end{verbatim}
393     if the system is configured correctly (Please talk to your system administrator).
394    
395     \begin{figure}
396     \begin{center}
397     \includegraphics[width=4in]{figures/ttblockspyplot150}
398 gross 2952 \caption{Example 1b: Total Energy in the Blocks over Time (in seconds).}
399 gross 2905 \label{fig:onedheatout1}
400     \end{center}
401     \end{figure}
402    
403 gross 2878 \subsection{Plotting the Total Energy}
404 gross 2949 \sslist{example01b.py}
405 gross 2878
406 gross 2905 \esc does not include its own plotting capabilities. However, it is possible to use a variety of free \pyt packages for visualisation.
407 gross 2878 Two types will be demonstrated in this cookbook; \mpl\footnote{\url{http://matplotlib.sourceforge.net/}} and \verb VTK \footnote{\url{http://www.vtk.org/}} visualisation.
408     The \mpl package is a component of SciPy\footnote{\url{http://www.scipy.org}} and is good for basic graphs and plots.
409     For more complex visualisation tasks in particular when it comes to two and three dimensional problems it is recommended to us more advanced tools for instance \mayavi \footnote{\url{http://code.enthought.com/projects/mayavi/}}
410 gross 2905 which bases on the \verb|VTK| toolkit. We will discuss the usage of \verb|VTK| based
411 gross 2878 visualization in Chapter~\ref{Sec:2DHD} where will discuss a two dimensional PDE.
412    
413     For our simple problem we have two plotting tasks: Firstly we are interested in showing the
414 gross 2905 behaviour of the total energy over time and secondly in how the temperature distribution within the block is
415 gross 2878 developing over time. Lets start with the first task.
416    
417     The trick is to create a record of the time marks and the corresponding total energies observed.
418     \pyt provides the concept of lists for this. Before
419     the time loop is opened we create empty lists for the time marks \verb|t_list| and the total energies \verb|E_list|.
420     After the new temperature as been calculated by solving the PDE we append the new time marker and total energy
421     to the corresponding list using the \verb|append| method. With these modifications the script looks as follows:
422 ahallam 2775 \begin{python}
423 gross 2878 t_list=[]
424     E_list=[]
425     # ... start iteration:
426     while t<tend:
427     t+=h
428     mypde.setValue(Y=qH+rhocp/h*T) # set variable PDE coefficients
429     T=mypde.getSolution() #get the PDE solution
430     totE=integrate(rhocp*T)
431     t_list.append(t) # add current time mark to record
432     E_list.append(totE) # add current total energy to record
433 ahallam 2775 \end{python}
434 gross 2878 To plot $t$ over $totE$ we use the \mpl a module contained within \pylab which needs to be loaded before used;
435 ahallam 2775 \begin{python}
436 gross 2878 import pylab as pl # plotting package.
437 ahallam 2775 \end{python}
438 gross 2905 Here we are not using the \verb|from pylab import *| in order to avoid name clashes for function names
439     within \esc.
440 ahallam 2401
441 gross 2878 The following statements are added to the script after the time loop has been completed;
442 ahallam 2775 \begin{python}
443 gross 2878 pl.plot(t_list,E_list)
444     pl.title("Total Energy")
445 gross 2905 pl.axis([0,max(t_list),0,max(E_list)*1.1])
446 gross 2878 pl.savefig("totE.png")
447 ahallam 2775 \end{python}
448 gross 2878 The first statement hands over the time marks and corresponding total energies to the plotter.
449 gross 2905 The second statment is setting the title for the plot. The third statement
450     sets the axis ranges. In most cases these are set appropriately by the plotter.
451     The last statement renders the plot and writes the
452     result into the file \verb|totE.png| which can be displayed by (almost) any image viewer.
453     As expected the total energy is constant over time, see \reffig{fig:onedheatout1}.
454 ahallam 2401
455 gross 2878 \subsection{Plotting the Temperature Distribution}
456 gross 2931 \label{sec: plot T}
457 gross 2949 \sslist{example01c.py}
458 gross 2878 For plotting the spatial distribution of the temperature we need to modify the strategy we have used
459     for the total energy. Instead of producing a final plot at the end we will generate a
460     picture at each time step which can be browsed as slide show or composed to a movie.
461     The first problem we encounter is that if we produce an image in each time step we need
462     to make sure that the images previously generated are not overwritten.
463    
464     To develop an incrementing file name we can use the following convention. It is convenient to
465     put all image file showing the same variable - in our case the temperature distribution -
466     into a separate directory. As part of the \verb|os| module\footnote{The \texttt{os} module provides
467     a powerful interface to interact with the operating system, see \url{http://docs.python.org/library/os.html}.} \pyt
468     provides the \verb|os.path.join| command to build file and
469     directory names in a platform independent way. Assuming that
470     \verb|save_path| is name of directory we want to put the results the command is;
471 ahallam 2775 \begin{python}
472 gross 2878 import os
473     os.path.join(save_path, "tempT%03d.png"%i )
474 ahallam 2775 \end{python}
475 gross 2878 where \verb|i| is the time step counter.
476     There are two arguments to the \verb join command. The \verb save_path variable is a predefined string pointing to the directory we want to save our data in, for example a single sub-folder called \verb data would be defined by;
477     \begin{verbatim}
478     save_path = "data"
479     \end{verbatim}
480 gross 2949 while a sub-folder of \verb data called \verb example01 would be defined by;
481 gross 2878 \begin{verbatim}
482 gross 2949 save_path = os.path.join("data","example01")
483 gross 2878 \end{verbatim}
484 gross 2905 The second argument of \verb join \xspace contains a string which is the file name or subdirectory name. We can use the operator \verb|%| to increment our file names with the value \verb|i| denoting a incrementing counter. The sub-string \verb %03d does this by defining the following parameters;
485 gross 2878 \begin{itemize}
486     \item \verb 0 becomes the padding number;
487     \item \verb 3 tells us the amount of padding numbers that are required; and
488     \item \verb d indicates the end of the \verb % operator.
489     \end{itemize}
490     To increment the file name a \verb %i is required directly after the operation the string is involved in. When correctly implemented the output files from this command would be place in the directory defined by \verb save_path as;
491     \begin{verbatim}
492 gross 2905 blockspyplot.png
493     blockspyplot.png
494     blockspyplot.png
495 gross 2878 ...
496     \end{verbatim}
497     and so on.
498    
499     A sub-folder check/constructor is available in \esc. The command;
500     \begin{verbatim}
501     mkDir(save_path)
502     \end{verbatim}
503     will check for the existence of \verb save_path and if missing, make the required directories.
504    
505     We start by modifying our solution script from before.
506 gross 2905 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 we create the node coordinates of the sample points used to represent
507     the temperature as a \pyt list of tuples or a \numpy array as requested by the plotting function.
508     We need to convert thearray \verb|x| previously set as \verb|Solution(blocks).getX()| into a \pyt list
509     and then to a \numpy array. The $x\hackscore{0}$ component is then extracted via an array slice to the variable \verb|plx|;
510 ahallam 2775 \begin{python}
511 gross 2878 import numpy as np # array package.
512     #convert solution points for plotting
513 gross 2905 plx = x.toListOfTuples()
514 gross 2878 plx = np.array(plx) # convert to tuple to numpy array
515     plx = plx[:,0] # extract x locations
516 ahallam 2775 \end{python}
517 gross 2878
518 ahallam 2645 \begin{figure}
519     \begin{center}
520 gross 2905 \includegraphics[width=4in]{figures/blockspyplot001}
521     \includegraphics[width=4in]{figures/blockspyplot050}
522     \includegraphics[width=4in]{figures/blockspyplot200}
523 gross 2952 \caption{Example 1c: Temperature ($T$) distribution in the blocks at time steps $1$, $50$ and $200$.}
524 ahallam 2645 \label{fig:onedheatout}
525     \end{center}
526     \end{figure}
527    
528 gross 2878 For each time step we will generate a plot of the temperature distribution and save each to a file. We use the same
529     techniques provided by \mpl as we have used to plot the total energy over time.
530     The following is appended to the end of the \verb while loop and creates one figure of the temperature distribution. We start by converting the solution to a tuple and then plotting this against our \textit{x coordinates} \verb plx we have generated before. We add a title to the diagram before it is rendered into a file.
531     Finally, the figure is saved to a \verb|*.png| file and cleared for the following iteration.
532 ahallam 2775 \begin{python}
533 gross 2878 # ... start iteration:
534     while t<tend:
535     ....
536     T=mypde.getSolution() #get the PDE solution
537     tempT = T.toListOfTuples() # convert to a tuple
538     pl.plot(plx,tempT) # plot solution
539     # set scale (Temperature should be between Tref and T0)
540     pl.axis([0,mx,Tref*.9,T0*1.1])
541     # add title
542 gross 2905 pl.title("Temperature across the blocks at time %e minutes"%(t/day))
543 gross 2878 #save figure to file
544 gross 2905 pl.savefig(os.path.join(save_path,"tempT","blockspyplot%03d.png") %i)
545 ahallam 2775 \end{python}
546 gross 2905 Some results are shown in \reffig{fig:onedheatout}.
547 jfenwick 2657
548 ahallam 2606 \subsection{Make a video}
549 gross 2878 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.
550 ahallam 2775 \begin{python}
551 gross 2878 # compile the *.png files to create a *.avi videos that show T change
552     # with time. This operation uses Linux mencoder. For other operating
553 gross 2905 # systems it is possible to use your favourite video compiler to
554 ahallam 2606 # convert image files to videos.
555 gross 2477
556 ahallam 2606 os.system("mencoder mf://"+save_path+"/tempT"+"/*.png -mf type=png:\
557 gross 2878 w=800:h=600:fps=25 -ovc lavc -lavcopts vcodec=mpeg4 -oac copy -o \
558 gross 2949 example01tempT.avi")
559 ahallam 2775 \end{python}
560 gross 2477

  ViewVC Help
Powered by ViewVC 1.1.26