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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 4853 - (show annotations)
Wed Apr 9 05:41:57 2014 UTC (5 years, 6 months ago) by jfenwick
File MIME type: text/x-python
File size: 5631 byte(s)
Added import division and changed a few / to //

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

  ViewVC Help
Powered by ViewVC 1.1.26