1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829
|
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
% Copyright (c) 2003-2018 by The University of Queensland
% http://www.uq.edu.au
%
% Primary Business: Queensland, Australia
% Licensed under the Apache License, version 2.0
% http://www.apache.org/licenses/LICENSE-2.0
%
% Development until 2012 by Earth Systems Science Computational Center (ESSCC)
% Development 2012-2013 by School of Earth Sciences
% Development from 2014 by Centre for Geoscience Computing (GeoComp)
%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\section{Example 1: One Dimensional Heat Diffusion in Granite}
\label{Sec:1DHDv00}
The first model consists of two blocks of isotropic material, for instance
granite, sitting next to each other (\autoref{fig:onedgbmodel}).
Initial temperature in \textit{Block 1} is \verb|T1| and in \textit{Block 2} is
\verb|T2|.
We assume that the system is insulated.
What would happen to the temperature distribution in each block over time?
Intuition tells us that heat will be transported from the hotter block to the
cooler one until both
blocks have the same temperature.
\begin{figure}[ht]
\centerline{\includegraphics[width=4.in]{figures/onedheatdiff001}}
\caption{Example 1: Temperature differential along a single interface between
two granite blocks.}
\label{fig:onedgbmodel}
\end{figure}
\subsection{1D Heat Diffusion Equation}
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}};
which is defined as:
\begin{equation}
\rho c_p \frac{\partial T}{\partial t} - \kappa \frac{\partial^{2}
T}{\partial x^{2}} = q_H
\label{eqn:hd}
\end{equation}
where $\rho$ is the material density, $c_p$ is the specific heat and
$\kappa$ is the thermal
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
parameters are \textbf{constant}.
The heat source is defined by the right hand side of \refEq{eqn:hd} as
$q_{H}$; this can take the form of a constant or a function of time and
space. For example $q_{H} = q_{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 signed distance from the block-block interface $x$.
\subsection{PDEs and the General Form}
It is possible to solve PDE \refEq{eqn:hd} analytically and obtain an exact
solution to our problem. However, it is not always practical to solve the
problem this way. Alternatively, computers can be used to find the solution. To
do this, a numerical approach is required to discretise
the PDE \refEq{eqn:hd} across time and space, this reduces the problem to a
finite number of equations for a finite number of spatial points and time steps.
These parameters together define the model. While discretisation introduces
approximations and a degree of error, a sufficiently sampled model is generally
accurate enough to satisfy the accuracy requirements for the final solution.
Firstly, we discretise the PDE \refEq{eqn:hd} in time. This leaves us with a
steady linear PDE which involves spatial derivatives only and needs to be solved
in each time step to progress in time. \esc can help us here.
For time discretisation we use the Backward Euler approximation
scheme\footnote{see \url{http://en.wikipedia.org/wiki/Euler_method}}. It is
based on the approximation
\begin{equation}
\frac{\partial T(t)}{\partial t} \approx \frac{T(t)-T(t-h)}{h}
\label{eqn:beuler}
\end{equation}
for $\frac{\partial T}{\partial t}$ at time $t$
where $h$ is the time step size. This can also be written as;
\begin{equation}
\frac{\partial T}{\partial t}(t^{(n)}) \approx \frac{T^{(n)} - T^{(n-1)}}{h}
\label{eqn:Tbeuler}
\end{equation}
where the upper index $n$ denotes the n\textsuperscript{th} time step. So one
has
\begin{equation}
\begin{array}{rcl}
t^{(n)} & = & t^{(n-1)}+h \\
T^{(n)} & = & T(t^{(n-1)}) \\
\end{array}
\label{eqn:Neuler}
\end{equation}
Substituting \refEq{eqn:Tbeuler} into \refEq{eqn:hd} we get;
\begin{equation}
\frac{\rho c_p}{h} (T^{(n)} - T^{(n-1)}) - \kappa \frac{\partial^{2}
T^{(n)}}{\partial x^{2}} = q_H
\label{eqn:hddisc}
\end{equation}
Notice that we evaluate the spatial derivative term at the current time
$t^{(n)}$ - therefore the name \textbf{backward Euler} scheme. Alternatively,
one can evaluate the spatial derivative term at the previous time $t^{(n-1)}$.
This approach is called the \textbf{forward Euler} scheme. This scheme can
provide some computational advantages, which
are not discussed here. However, the \textbf{forward Euler} scheme has a major
disadvantage. Namely, depending on the
material parameters as well as the domain 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. Stability is achieved if
the temperature does not grow beyond its initial bounds and becomes
non-physical.
The backward Euler scheme, which we use here, is unconditionally stable meaning
that under the assumption of a
physically correct problem set-up the temperature approximation remains physical
for all time steps.
The user needs to keep in mind that the discretisation error introduced by
\refEq{eqn:beuler}
is sufficiently small, thus a good approximation of the true temperature is
computed. It is therefore very important that any results are viewed with
caution. For example, one may compare the results for different time and
spatial step sizes.
To get the temperature $T^{(n)}$ at time $t^{(n)}$ we need to solve the linear
differential equation \refEq{eqn:hddisc} which only includes spatial
derivatives. To solve this problem we want to use \esc.
In \esc any given PDE can be described by the general form. For the purpose of
this introduction we illustrate a simpler version of the general form for full
linear PDEs which is available in the \esc user's guide. A simplified form that
suits our heat diffusion problem\footnote{The form in the \esc users guide which
uses the Einstein convention is written as
$-(A_{jl} u_{,l})_{,j}+D u =Y$}
is described by;
\begin{equation}\label{eqn:commonform nabla}
-\nabla\cdot(A\cdot\nabla u) + Du = f
\end{equation}
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
the spatial derivative of its subject - in this case $u$. Lets assume for a
moment that we deal with a one-dimensional problem then ;
\begin{equation}
\nabla = \frac{\partial}{\partial x}
\end{equation}
and we can write \refEq{eqn:commonform nabla} as;
\begin{equation}\label{eqn:commonform}
-A\frac{\partial^{2}u}{\partial x^{2}} + Du = f
\end{equation}
if $A$ is constant. To match this simplified general form to our problem
\refEq{eqn:hddisc}
we rearrange \refEq{eqn:hddisc};
\begin{equation}
\frac{\rho c_p}{h} T^{(n)} - \kappa \frac{\partial^2 T^{(n)}}{\partial
x^2} = q_H + \frac{\rho c_p}{h} T^{(n-1)}
\label{eqn:hdgenf}
\end{equation}
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
$t^{(0)}=0$ and $t^{(n)}=t^{(n-1)}+h$ where $h>0$ is the step size and assumed
to be constant.
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} one can see
that;
\begin{equation}\label{ESCRIPT SET}
u=T^{(n)};
A = \kappa; D = \frac{\rho c _{p}}{h}; f = q _{H} + \frac{\rho
c_p}{h} T^{(n-1)}
\end{equation}
\subsection{Boundary Conditions}
\label{SEC BOUNDARY COND}
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.
A \textbf{Dirichlet boundary condition} is conceptually simpler and is used to
prescribe a known value to the unknown solution (in our example the temperature)
on parts of the boundary or on the entire boundary of the region of interest.
We discuss the Dirichlet boundary condition in our second example presented in
Section~\ref{Sec:1DHDv0}.
However, for this example we have made the model assumption that the system is
insulated, so we need to add an appropriate boundary condition to prevent
any loss or inflow of energy at the boundary of our domain. Mathematically this
is expressed by prescribing
the heat flux $\kappa \frac{\partial T}{\partial x}$ to zero. In our simplified
one dimensional model this is expressed
in the form;
\begin{equation}
\kappa \frac{\partial T}{\partial x} = 0
\end{equation}
or in a more general case as
\begin{equation}\label{NEUMAN 1}
\kappa \nabla T \cdot n = 0
\end{equation}
where $n$ is the outer normal field \index{outer normal field} at the surface
of the domain.
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
the temperature $T$. Other notations used here are\footnote{The \esc notation
for the normal
derivative is $T_{,i} n_i$.};
\begin{equation}
\nabla T \cdot n = \frac{\partial T}{\partial n} \; .
\end{equation}
A condition of the type \refEq{NEUMAN 1} defines a \textbf{Neumann boundary
condition} for the PDE.
The PDE \refEq{eqn:hdgenf}
and the Neumann boundary condition~\ref{eqn:hdgenf} (potentially together with
the Dirichlet boundary conditions) define a \textbf{boundary value problem}.
It is the nature of a boundary value problem to allow making statements about
the solution in the
interior of the domain from information known on the boundary only. In most
cases we use the term partial differential equation but in fact it is a
boundary value problem.
It is important to keep in mind that boundary conditions need to be complete and
consistent in the sense that
at any point on the boundary either a Dirichlet or a Neumann boundary condition
must be set.
Conveniently, \esc makes a default assumption on the boundary conditions which
the user may modify where appropriate.
For a problem of the form in~\refEq{eqn:commonform nabla} the default
condition\footnote{In the \esc user guide which uses the Einstein convention
this is written as
$n_{j}A_{jl} u_{,l}=0$.} is;
\begin{equation}\label{NEUMAN 2}
-n\cdot A \cdot\nabla u = 0
\end{equation}
which is used everywhere on the boundary. Again $n$ denotes the outer normal
field.
Notice that the coefficient $A$ is the same as in the \esc
PDE~\ref{eqn:commonform nabla}.
With the settings for the coefficients we have already identified in
\refEq{ESCRIPT SET} this
condition translates into
\begin{equation}\label{NEUMAN 2b}
\kappa \frac{\partial T}{\partial x} = 0
\end{equation}
for the boundary of the domain. This is identical to the Neumann boundary
condition we want to set. \esc will take care of this condition for us. We
discuss the Dirichlet boundary condition later.
\subsection{Outline of the Implementation}
\label{sec:outline}
To solve the heat diffusion equation (\refEq{eqn:hd}) we 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 (discussed in \refSec{sec:key}) has
four major steps. Firstly, we need to define the domain where we want to
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 to solve in
each time step to get the updated temperature. Thirdly, we need to define 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. The work flow
is described in \reffig{fig:wf}.
% \begin{enumerate}
% \item create domain
% \item create PDE
% \item while end time not reached:
% \begin{enumerate}
% \item set PDE coefficients
% \item solve PDE
% \item update time marker
% \end{enumerate}
% \item end of calculation
% \end{enumerate}
\begin{figure}[h!]
\centering
\includegraphics[width=1in]{figures/workflow.png}
\caption{Workflow for developing an \esc model and solution}
\label{fig:wf}
\end{figure}
In the terminology of \pyt, the domain and PDE are represented by
\textbf{objects}. The nice feature of an object is that it is defined by its
usage and features
rather than its actual representation. So we will create a domain object to
describe the geometry of the two
granite blocks. Then we define PDEs and spatially distributed values such as the
temperature
on this domain. Similarly, to define a PDE object we use the fact that one needs
only to define the coefficients of the PDE and solve the PDE. The PDE object has
advanced features, but these are not required in simple cases.
\begin{figure}[htp]
\centering
\includegraphics[width=6in]{figures/functionspace.pdf}
\caption{\esc domain construction overview}
\label{fig:fs}
\end{figure}
\subsection{The Domain Constructor in \esc}
\label{ss:domcon}
Whilst it is not strictly relevant or necessary, a better understanding of
how values are spatially distributed (\textit{e.g.} Temperature) and how PDE
coefficients are interpreted in \esc can be helpful.
There are various ways to construct domain objects. The simplest form is a
rectangular shaped region with a length and height. There is
a ready to use function for this named \verb|rectangle()|. Besides the spatial
dimensions this function requires to specify the number of
elements or cells to be used along the length and height, see \reffig{fig:fs}.
Any spatially distributed value
and the PDE is represented in discrete form using this element
representation\footnote{We 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.
The quality of the approximation depends - besides other factors - mainly on the
number of elements being used. In fact, the
approximation becomes better when more elements are used. However, computational
cost grows with the number of
elements being used. It is therefore important that you find the right balance
between the demand in accuracy and acceptable resource usage.
In general, one can think about a domain object as a composition of nodes and
elements.
As shown in \reffig{fig:fs}, an element is defined by the nodes that are used to
describe its vertices.
To represent spatially distributed values the user can use
the values at the nodes, at the elements in the interior of the domain or at the
elements located on the surface of the domain.
The different approach used to represent values is called \textbf{function
space} and is attached to all objects
in \esc representing a spatially distributed value such as the solution of
a PDE. The three function spaces we use at the moment are;
\begin{enumerate}
\item the nodes, called by \verb|ContinuousFunction(domain)| ;
\item the elements/cells, called by \verb|Function(domain)| ; and
\item the boundary, called by \verb|FunctionOnBoundary(domain)|.
\end{enumerate}
A function space object such as \verb|ContinuousFunction(domain)| has the method
\verb|getX| attached to it. This method returns the
location of the so-called \textbf{sample points} used to represent values of the
particular function space. So the
call \verb|ContinuousFunction(domain).getX()| will return the coordinates of the
nodes used to describe the domain while
\verb|Function(domain).getX()| returns the coordinates of numerical
integration points within elements, see \reffig{fig:fs}.
This distinction between different representations of spatially distributed
values
is important in order to be able to vary the degrees of smoothness in a PDE
problem.
The coefficients of a PDE do not need to be continuous, thus this qualifies as a
\verb|Function()| type.
On the other hand a temperature distribution must be continuous and needs to be
represented with a \verb|ContinuousFunction()| function space.
An influx may only be defined at the boundary and is therefore a
\verb|FunctionOnBoundary()| object.
\esc allows certain transformations of the function spaces. A
\verb|ContinuousFunction()| can be transformed into a
\verb|FunctionOnBoundary()| or \verb|Function()|. On the other hand there is
not enough information in a \verb|FunctionOnBoundary()| to transform it to a
\verb|ContinuousFunction()|.
These transformations, which are called \textbf{interpolation} are invoked
automatically by \esc if needed.
Later in this introduction we discuss how
to define specific areas of geometry with different materials which are
represented by different material coefficients such as the
thermal conductivities $\kappa$. A very powerful technique to define these types
of PDE
coefficients is tagging. Blocks of materials and boundaries can be named and
values can be defined on subregions based on their names.
This is a method for simplifying PDE coefficient and flux definitions. It makes
scripting much easier and we will discuss this technique in
Section~\ref{STEADY-STATE HEAT REFRACTION}.
\subsection{A Clarification for the 1D Case}
\label{SEC: 1D CLARIFICATION}
It is necessary for clarification that we revisit our general PDE from
\refeq{eqn:commonform nabla} for a two dimensional domain. \esc is inherently
designed to solve problems that are multi-dimensional 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})$. Assuming the coefficient $A$ is constant, the \refEq{eqn:commonform nabla}
takes the following form;
\begin{equation}\label{eqn:commonform2D}
-A_{00}\frac{\partial^{2}u}{\partial x^{2}}
-A_{01}\frac{\partial^{2}u}{\partial x\partial y}
-A_{10}\frac{\partial^{2}u}{\partial y\partial x}
-A_{11}\frac{\partial^{2}u}{\partial y^{2}}
+ Du = f
\end{equation}
Notice that for the higher dimensional case $A$ becomes a matrix. It is also
important to notice that the usage of the Nabla operator creates
a compact formulation which is also independent from the spatial dimension.
To make the general PDE \refEq{eqn:commonform2D} one dimensional as
shown in \refEq{eqn:commonform} we need to set
\begin{equation}
A_{00}=A; A_{01}=A_{10}=A_{11}=0
\end{equation}
\subsection{Developing a PDE Solution Script}
\label{sec:key}
\sslist{example01a.py}
We write a simple \pyt script which uses the \modescript, \modfinley and \modmpl
modules.
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.}
that we will require.
\begin{python}
from esys.escript import *
# This defines the LinearPDE module as LinearPDE
from esys.escript.linearPDEs import LinearPDE
# This imports the rectangle domain function from finley.
from esys.finley import Rectangle
# A useful unit handling package which will make sure all our units
# match up in the equations under SI.
from esys.escript.unitsSI import *
\end{python}
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
domain. The module \verb|unitsSI| provides support for SI unit definitions with
our variables.
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
need values. Firstly, the domain upon which we wish to solve our problem needs
to be defined. There are different types of domains in \modescript which we
demonstrate in later tutorials but for our granite blocks, we simply use a
rectangular domain.
Using a rectangular domain 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 due to symmetry. There are four arguments
we must consider when we decide to create a rectangular domain, the domain
\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 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| needs 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.
\begin{python}
mx = 500.*m #meters - model length
my = 100.*m #meters - model width
ndx = 50 # mesh steps in x direction
ndy = 1 # mesh steps in y direction
boundloc = mx/2 # location of boundary between the two blocks
\end{python}
The material constants and the temperature variables must also be defined. For
the granite in the model they are defined as:
\begin{python}
#PDE related
rho = 2750. *kg/m**3 #kg/m^{3} density of iron
cp = 790.*J/(kg*K) # J/Kg.K thermal capacity
rhocp = rho*cp
kappa = 2.2*W/m/K # watts/m.Kthermal conductivity
qH=0 * J/(sec*m**3) # J/(sec.m^{3}) no heat source
T1=20 * Celsius # initial temperature at Block 1
T2=2273. * Celsius # base temperature at Block 2
\end{python}
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:
\begin{python}
t=0 * day # our start time, usually zero
tend=50 * yr # - time to end simulation
outputs = 200 # number of time steps required.
h=(tend-t)/outputs #size of time step
#user warning statement
print("Expected Number of time outputs is: ", (tend-t)/h)
i=0 #loop counter
\end{python}
Now that we know our inputs we will build a domain using the
\verb|Rectangle()| function from \FINLEY. The four arguments allow us to
define our domain \verb|model| as:
\begin{python}
#generate domain using rectangle
blocks = Rectangle(l0=mx,l1=my,n0=ndx, n1=ndy)
\end{python}
\verb|blocks| now describes a domain in the manner of Section \ref{ss:domcon}.
With a domain and all the 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 contrast to a system of PDEs which we discuss
later.}. We also need to state the values of our general form variables.
\begin{python}
mypde=LinearPDE(blocks)
A=zeros((2,2)))
A[0,0]=kappa
mypde.setValue(A=A, D=rhocp/h)
\end{python}
In many cases it may be possible to decrease the computational time of the
solver if the PDE is symmetric.
Symmetry of a PDE is defined by;
\begin{equation}\label{eqn:symm}
A_{jl}=A_{lj}
\end{equation}
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$. 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 enable symmetry via;
\begin{python}
myPDE.setSymmetryOn()
\end{python}
Next we need to establish the initial temperature distribution \verb|T|. We need
to
assign the value \verb|T1| to all sample points left to the contact interface at
$x_{0}=\frac{mx}{2}$
and the value \verb|T2| right to the contact interface. \esc
provides the \verb|whereNegative| function to construct this. More
specifically, \verb|whereNegative| returns the value $1$ at those sample points
where the argument has a negative value. Otherwise zero is returned.
If \verb|x| are the $x_{0}$
coordinates of the sample points used to represent the temperature distribution
then \verb|x[0]-boundloc| gives us a negative value for
all sample points left to the interface and non-negative value to
the right of the interface. So with;
\begin{python}
# ... set initial temperature ....
T= T1*whereNegative(x[0]-boundloc)+T2*(1-whereNegative(x[0]-boundloc))
\end{python}
we get the desired temperature distribution. To get the actual sample points
\verb|x| we use the \verb|getX()| method of the function space
\verb|Solution(blocks)| which is used to represent the solution of a PDE;
\begin{python}
x=Solution(blocks).getX()
\end{python}
As \verb|x| are the sample points for the function space
\verb|Solution(blocks)|
the initial temperature \verb|T| is using these sample points for
representation.
Although \esc is trying to be forgiving with the choice of sample points and to
convert
where necessary the adjustment of the function space is not always possible. So
it is advisable to make a careful choice on the function space used.
Finally we 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.
\begin{python}
while t < tend:
i+=1 #increment the counter
t+=h #increment the current time
mypde.setValue(Y=qH+rhocp/h*T) # set variable PDE coefficients
T=mypde.getSolution() #get the PDE solution
totE = integrate(rhocp*T) #get the total heat (energy) in the system
\end{python}
The last statement in this script calculates the total energy in the system as
the volume integral of $\rho c_{p} T$ over the block.
As the blocks are insulated no energy should be lost or added.
The total energy should stay constant for the example discussed here.
\subsection{Running the Script}
The script presented so far is available under
\verb|example01a.py|. You can edit this file with your favourite text editor.
On most operating systems\footnote{The \texttt{run-escript} launcher is not
supported under {\it MS Windows}.} you can use the
\program{run-escript} command
to launch {\it escript} scripts. For the example script use;
\begin{verbatim}
run-escript example01a.py
\end{verbatim}
The program will print a progress report. Alternatively, you can use
the python interpreter directly;
\begin{verbatim}
python example01a.py
\end{verbatim}
if the system is configured correctly (please talk to your system
administrator).
\subsection{Plotting the Total Energy}
\sslist{example01b.py}
\esc does not include its own plotting capabilities. However, it is possible to
use a variety of free \pyt packages for visualisation.
Two types will be demonstrated in this cookbook;
\mpl\footnote{\url{http://matplotlib.sourceforge.net/}} and
\verb|VTK|\footnote{\url{http://www.vtk.org/}}.
The \mpl package is a component of SciPy\footnote{\url{http://www.scipy.org}}
and is good for basic graphs and plots.
For more complex visualisation tasks, in particular two and three dimensional
problems we recommend the use of more advanced tools. For instance, \mayavi
\footnote{\url{http://code.enthought.com/projects/mayavi/}}
which is based upon the \verb|VTK| toolkit. The usage of \verb|VTK| based
visualisation is discussed in Chapter~\ref{Sec:2DHD} which focuses on a two
dimensional PDE.
For our simple granite block problem, we have two plotting tasks. Firstly, we
are interested in showing the
behaviour of the total energy over time and secondly, how the temperature
distribution within the block is developing over time.
Let us start with the first task.
The idea is to create a record of the time marks and the corresponding total
energies observed.
\pyt provides the concept of lists for this. Before
the time loop is opened we create empty lists for the time marks \verb|t_list|
and the total energies \verb|E_list|.
After the new temperature has been calculated by solving the PDE we append the
new time marker and the total energy value for that time
to the corresponding list using the \verb|append| method. With these
modifications our script looks as follows:
\begin{python}
t_list=[]
E_list=[]
# ... start iteration:
while t<tend:
t+=h
mypde.setValue(Y=qH+rhocp/h*T) # set variable PDE coefficients
T=mypde.getSolution() #get the PDE solution
totE=integrate(rhocp*T)
t_list.append(t) # add current time mark to record
E_list.append(totE) # add current total energy to record
\end{python}
To plot $t$ over $totE$ we use \mpl a module contained within \pylab which needs
to be loaded before use;
\begin{python}
import pylab as pl # plotting package.
\end{python}
Here we are not using \verb|from pylab import *| in order to avoid name
clashes for function names within \esc.
The following statements are added to the script after the time loop has been
completed;
\begin{python}
pl.plot(t_list,E_list)
pl.title("Total Energy")
pl.axis([0,max(t_list),0,max(E_list)*1.1])
pl.savefig("totE.png")
\end{python}
The first statement hands over the time marks and corresponding total energies
to the plotter.
The second statement sets the title for the plot. The third statement
sets the axis ranges. In most cases these are set appropriately by the plotter.
The last statement generates the plot and writes the result into the file
\verb|totE.png| which can be displayed by (almost) any image viewer.
As expected the total energy is constant over time, see
\reffig{fig:onedheatout1}.
\begin{figure}[ht]
\begin{center}
\includegraphics[width=4in]{figures/ttblockspyplot150}
\caption{Example 1b: Total Energy in the Blocks over Time (in seconds)}
\label{fig:onedheatout1}
\end{center}
\end{figure}
\clearpage
\subsection{Plotting the Temperature Distribution}
\label{sec: plot T}
\sslist{example01c.py}
For plotting the spatial distribution of the temperature we need to modify the
strategy we have used for the total energy.
Instead of producing a final plot at the end we will generate a
picture at each time step which can be browsed as a slide show or composed into
a movie.
The first problem we encounter is that if we produce an image at each time step
we need to make sure that the images previously generated are not overwritten.
To develop an incrementing file name we can use the following convention. It is
convenient to put all image files showing the same variable - in our case the
temperature distribution - into a separate directory.
As part of the \verb|os| module\footnote{The \texttt{os} module provides
a powerful interface to interact with the operating system, see
\url{http://docs.python.org/library/os.html}.} \pyt
provides the \verb|os.path.join| command to build file and directory names in a
platform independent way. Assuming that
\verb|save_path| is the name of the directory we want to put the results in the
command is;
\begin{python}
import os
os.path.join(save_path, "tempT%03d.png"%i )
\end{python}
where \verb|i| is the time step counter.
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, for example a single sub-folder called \verb|data| would be defined by;
\begin{verbatim}
save_path = "data"
\end{verbatim}
while a sub-folder of \verb|data| called \verb|example01| would be defined by;
\begin{verbatim}
save_path = os.path.join("data","example01")
\end{verbatim}
The second argument of \verb|join| contains a string which is the file
name or subdirectory name. We can use the operator \verb|%| to use the value of
\verb|i| as part of our filename. The sub-string \verb|%03d| indicates that we
want to substitute a value into the name;
\begin{itemize}
\item \verb|0| means that small numbers should have leading zeroes;
\item \verb|3| means that numbers should be written using at least 3 digits;
and
\item \verb|d| means that the value to substitute will be a decimal integer.
\end{itemize}
To actually substitute the value of \verb|i| into the name write \verb|%i| after
the string.
When done correctly, the output files from this command will be placed in the
directory defined by \verb|save_path| as;
\begin{verbatim}
blockspyplot001.png
blockspyplot002.png
blockspyplot003.png
...
\end{verbatim}
and so on.
A sub-folder check/constructor is available in \esc. The command;
\begin{verbatim}
mkDir(save_path)
\end{verbatim}
will check for the existence of \verb|save_path| and if missing, create the
required directories.
We start by modifying our solution script.
Prior to the \verb|while| loop we 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
the temperature as a \pyt list of tuples or a \numpy array as requested by the
plotting function.
We need to convert the array \verb|x| previously set as
\verb|Solution(blocks).getX()| into a \pyt list
and then to a \numpy array. The $x_{0}$ component is then extracted via
an array slice to the variable \verb|plx|;
\begin{python}
import numpy as np # array package.
#convert solution points for plotting
plx = x.toListOfTuples()
plx = np.array(plx) # convert to tuple to numpy array
plx = plx[:,0] # extract x locations
\end{python}
\begin{figure}
\begin{center}
\includegraphics[width=4in]{figures/blockspyplot001}
\includegraphics[width=4in]{figures/blockspyplot050}
\includegraphics[width=4in]{figures/blockspyplot200}
\caption{Example 1c: Temperature ($T$) distribution in the blocks at time steps
$1$, $50$ and $200$}
\label{fig:onedheatout}
\end{center}
\end{figure}
\clearpage
We use the same techniques provided by \mpl as we have used to plot the total
energy over time.
For each time step we generate a plot of the temperature distribution and save
each to a file.
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.
Finally, the figure is saved to a \verb|*.png| file and cleared for the
following iteration.
\begin{python}
# ... start iteration:
while t<tend:
i+=1
t+=h
mypde.setValue(Y=qH+rhocp/h*T)
T=mypde.getSolution()
totE=integrate(rhocp*T)
print("time step %s at t=%e days completed. total energy = %e."%(i,t/day,totE))
t_list.append(t)
E_list.append(totE)
#establish figure 1 for temperature vs x plots
tempT = T.toListOfTuples()
pl.figure(1) #current figure
pl.plot(plx,tempT) #plot solution
# add title
pl.axis([0,mx,T1*.9,T2*1.1])
pl.title("Temperature across blocks at time %d days"%(t/day))
#save figure to file
pl.savefig(os.path.join(save_path,"tempT", "blockspyplot%03d.png"%i))
pl.clf() #clear figure
\end{python}
Some results are shown in \reffig{fig:onedheatout}.
\subsection{Making a Video}
Our saved plots from the previous section can be cast into a video using the
following command appended to the end of the script. The \verb mencoder command
is not available on every platform, so some users need to use an alternative
video encoder.
\begin{python}
# compile the *.png files to create a *.avi video that shows T change
# with time. This operation uses Linux mencoder. For other operating
# systems it is possible to use your favourite video compiler to
# convert image files to videos.
os.system("mencoder mf://"+save_path+"/tempT"+"/*.png -mf type=png:\
w=800:h=600:fps=25 -ovc lavc -lavcopts vcodec=mpeg4 -oac copy -o \
example01tempT.avi")
\end{python}
|