/[escript]/trunk/doc/examples/cookbook/example02.py
ViewVC logotype

Annotation of /trunk/doc/examples/cookbook/example02.py

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3259 - (hide annotations)
Mon Oct 11 01:48:14 2010 UTC (9 years ago) by jfenwick
File MIME type: text/x-python
File size: 5227 byte(s)
Merging dudley and scons updates from branches

1 ahallam 2401
2     ########################################################
3     #
4 gross 2904 # Copyright (c) 2009 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 2904 __copyright__="""Copyright (c) 2009 by University of Queensland
15 ahallam 2401 Earth Systems Science Computational Center (ESSCC)
16     http://www.uq.edu.au/esscc
17     Primary Business: Queensland, Australia"""
18     __license__="""Licensed under the Open Software License version 3.0
19     http://www.opensource.org/licenses/osl-3.0.php"""
20     __url__="https://launchpad.net/escript-finley"
21    
22     """
23     Author: Antony Hallam antony.hallam@uqconnect.edu.au
24     """
25    
26 ahallam 2977 ############################################################FILE HEADER
27     # example02.py
28     # Model temperature diffusion along an insulated iron rod with a
29     # heating element at the left hand side.
30    
31     #######################################################EXTERNAL MODULES
32 gross 2904 # To solve the problem it is necessary to import the modules we require.
33     from esys.escript import * # This imports everything from the escript library
34     from esys.escript.unitsSI import *
35     from esys.escript.linearPDEs import LinearPDE # This defines LinearPDE as LinearPDE
36 jfenwick 3259 from esys.finley import Rectangle # This imports the rectangle domain function
37 ahallam 2658 #For interactive use, you can comment out the next two lines
38     import matplotlib
39     matplotlib.use('agg') #It's just here for automated testing
40 ahallam 2606 import pylab as pl #Plotting package.
41     import numpy as np #Array package.
42 gross 2904 import os, sys #This package is necessary to handle saving our data.
43 ahallam 2401
44 ahallam 2977 ########################################################MPI WORLD CHECK
45 ahallam 2658 if getMPISizeWorld() > 1:
46 gross 2904 import sys
47     print "This example will not run in an MPI world."
48     sys.exit(0)
49 jfenwick 2648
50 ahallam 2977 #################################################ESTABLISHING VARIABLES
51 gross 2904 #Domain related.
52     mx = 1*m #meters - model length
53     my = .1*m #meters - model width
54     ndx = 100 # mesh steps in x direction
55     ndy = 1 # mesh steps in y direction - one dimension means one element
56 ahallam 2401 #PDE related
57 gross 2904 rho = 7874. *kg/m**3 #kg/m^{3} density of iron
58     cp = 449.*J/(kg*K) # J/Kg.K thermal capacity
59     rhocp = rho*cp
60     kappa = 80.*W/m/K # watts/m.Kthermal conductivity
61     qH=0 * J/(sec*m**3) # J/(sec.m^{3}) no heat source
62     Tref = 20 * Celsius # base temperature of the rod
63     T0 = 100 * Celsius # temperature at heating element
64 ahallam 2606
65 ahallam 2977 ################################################ESTABLISHING PARAMETERS
66 gross 2904 t=0 * day # our start time, usually zero
67     tend= 0.5 *day # - time to end simulation
68     outputs = 200 # number of time steps required.
69 ahallam 2401 h=(tend-t)/outputs #size of time step
70 ahallam 2606 #user warning statement
71 gross 2904 print "Expected Number of time outputs is: ", (tend-t)/h
72     i=0 #loop counter
73     #the folder to put our outputs in, leave blank "" for script path
74 gross 2949 save_path= os.path.join("data","example02")
75 gross 2904 #ensure the dir exists
76     mkDir(save_path, os.path.join(save_path,"tempT"))
77 ahallam 2401
78 ahallam 2977 ####################################################DOMAIN CONSTRUCTION
79 gross 2904 rod = Rectangle(l0=mx,l1=my,n0=ndx, n1=ndy)
80     x=Solution(rod).getX()
81 ahallam 2977 ###############################################ESCRIPT PDE CONSTRUCTION
82 gross 2904 mypde=LinearPDE(rod)
83     A=zeros((2,2))
84     A[0,0]=kappa
85     q=whereZero(x[0])
86     mypde.setValue(A=A, D=rhocp/h, q=q, r=T0)
87     # ... set initial temperature ....
88     T= T0*whereZero(x[0])+Tref*(1-whereZero(x[0]))
89 ahallam 2401
90 gross 2904 # ... open a collector for the time marks and corresponding total energy
91     t_list=[]
92     E_list=[]
93     # ... convert solution points for plotting
94 ahallam 2606 plx = x.toListOfTuples()
95     plx = np.array(plx) #convert to tuple to numpy array
96     plx = plx[:,0] #extract x locations
97 ahallam 2977 ########################################################START ITERATION
98 gross 2904 while t<tend:
99     i+=1
100     t+=h
101     mypde.setValue(Y=qH+rhocp/h*T)
102     T=mypde.getSolution()
103     totE=integrate(rhocp*T)
104     print "time step %s at t=%e minutes completed. total energy = %e."%(i,t/minute,totE)
105     t_list.append(t)
106     E_list.append(totE)
107 ahallam 2401
108 gross 2904 #establish figure 1 for temperature vs x plots
109     tempT = T.toListOfTuples()
110     pl.figure(1) #current figure
111     pl.plot(plx,tempT) #plot solution
112     # add title
113     pl.axis([0,mx,Tref*.9,T0*1.1])
114     pl.title("Temperature across rod at time %e hours"%(t/hour))
115     #save figure to file
116     pl.savefig(os.path.join(save_path,"tempT", "rodpyplot%03d.png"%i))
117     pl.clf() #clear figure
118 ahallam 2977
119     ###############################################################PLOTTING
120 gross 2904 # plot the total energy over time:
121     pl.figure(2)
122     pl.plot(t_list,E_list)
123     pl.title("Total Energy")
124     # pl.axis([0,max(t_list),0,max(E_list)*1.1])
125     pl.savefig(os.path.join(save_path,"totE.png"))
126     pl.clf()
127    
128 ahallam 2977 #########################################################CREATE A MOVIE
129 gross 2904 # compile the *.png files to create a*.avi video that show T change
130 ahallam 2658 # with time. This opperation uses linux mencoder. For other operating
131 jfenwick 2667 # systems it may be possible to use your favourite video compiler to
132 ahallam 2658 # convert image files to videos. To enable this step uncomment the
133     # following lines.
134    
135 gross 2904 #os.system("mencoder mf://"+save_path+"/tempT"+"/*.png -mf type=png:\
136 ahallam 2658 #w=800:h=600:fps=25 -ovc lavc -lavcopts vcodec=mpeg4 -oac copy -o \
137 gross 2949 #example02tempT.avi")

  ViewVC Help
Powered by ViewVC 1.1.26