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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 4087 - (hide annotations)
Thu Nov 22 22:28:01 2012 UTC (6 years, 10 months ago) by caltinay
File MIME type: text/x-python
File size: 5494 byte(s)
Moved matplotlib imports in test scripts before escript since there is an
import chain which pulls it so the use() function has no effect.

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

  ViewVC Help
Powered by ViewVC 1.1.26