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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3148 - (show annotations)
Fri Sep 3 02:09:47 2010 UTC (10 years, 10 months ago) by jfenwick
File MIME type: text/x-python
File size: 6159 byte(s)
Another attempt to patch the X issue

1
2 ########################################################
3 #
4 # Copyright (c) 2009-2010 by University of Queensland
5 # 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 __copyright__="""Copyright (c) 2009-2010 by University of Queensland
15 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 ############################################################FILE HEADER
23 # example09.py
24 # Antony Hallam
25 # Seismic Wave Equation Simulation using acceleration solution.
26 # 3D model with multiple layers.
27
28 #######################################################EXTERNAL MODULES
29 from esys.escript import *
30 from esys.finley import Rectangle
31 import os
32 # smoothing operator
33 from esys.escript.pdetools import Projector, Locator
34 from esys.escript.unitsSI import *
35 import numpy as np
36 import matplotlib
37 matplotlib.use('agg') #It's just here for automated testing
38
39 import pylab as pl
40 import matplotlib.cm as cm
41 from esys.escript.linearPDEs import LinearPDE
42 from esys.finley import ReadMesh
43
44 ########################################################MPI WORLD CHECK
45 if getMPISizeWorld() > 1:
46 import sys
47 print "This example will not run in an MPI world."
48 sys.exit(0)
49
50 #################################################ESTABLISHING VARIABLES
51 # where to save output data
52 savepath = "data/example09b2"
53 meshpath = "data/example09m"
54 mkDir(savepath)
55 #Geometric and material property related variables.
56 mx = 200. # model lenght
57 my = 200. # model width
58 mz=100.0
59 step=4.0 # the element size
60 ndx = int(mx/step) # steps in x direction
61 ndy = int(my/step) # steps in y direction
62 ndz = int(mz/step)
63
64 vel2=1800.; vel1=3000.
65 rho2=2300.; rho1=3100. #density
66 mu2=vel2**2.*rho2/4.; mu1=vel1**2.*rho1/4. #bulk modulus
67 lam2=vel2**2.*rho2/2.; lam1=vel1**2.*rho1/2. #lames constant
68
69 # Time related variables.
70 tend=0.1 # end time
71 h=0.00005 # time step
72 # data recording times
73 rtime=0.0 # first time to record
74 rtime_inc=tend/200.0 # time increment to record
75 #Check to make sure number of time steps is not too large.
76 print "Time step size= ",h, "Expected number of outputs= ",tend/h
77
78 U0=0.1 # amplitude of point source
79 ls=500 # length of the source
80 source=np.zeros(ls,'float') # source array
81 decay1=np.zeros(ls,'float') # decay curve one
82 decay2=np.zeros(ls,'float') # decay curve two
83 time=np.zeros(ls,'float') # time values
84 g=np.log(0.01)/ls
85
86 dfeq=50 #Dominant Frequency
87 a = 2.0 * (np.pi * dfeq)**2.0
88 t0 = 5.0 / (2.0 * np.pi * dfeq)
89 srclength = 5. * t0
90 ls = int(srclength/h)
91 print 'source length',ls
92 source=np.zeros(ls,'float') # source array
93 ampmax=0
94 for it in range(0,ls):
95 t = it*h
96 tt = t-t0
97 dum1 = np.exp(-a * tt * tt)
98 source[it] = -2. * a * tt * dum1
99 if (abs(source[it]) > ampmax):
100 ampmax = abs(source[it])
101 time[t]=t*h
102 #tdecay=decay1*decay2*U0
103 #decay1=decay1*U0; decay2=decay2*U0
104 #pl.clf();
105 #pl.plot(source)
106 #pl.plot(time,decay1);pl.plot(time,decay2);
107 #pl.plot(time,tdecay)
108 #pl.savefig(os.path.join(savepath,'source.png'))
109
110 # will introduce a spherical source at middle left of bottom face
111 xc=[mx/2,my/2,0]
112
113 ####################################################DOMAIN CONSTRUCTION
114 domain=ReadMesh(os.path.join(meshpath,'example09m.fly')) # create the domain
115 x=domain.getX() # get the locations of the nodes in the domain
116
117 lam=Scalar(0,Function(domain))
118 lam.setTaggedValue("vintfa",lam1)
119 lam.setTaggedValue("vintfb",lam2)
120 mu=Scalar(0,Function(domain))
121 mu.setTaggedValue("vintfa",mu1)
122 mu.setTaggedValue("vintfb",mu2)
123 rho=Scalar(0,Function(domain))
124 rho.setTaggedValue("vintfa",rho1)
125 rho.setTaggedValue("vintfb",rho2)
126
127 ##########################################################ESTABLISH PDE
128 mypde=LinearPDE(domain) # create pde
129 mypde.setSymmetryOn() # turn symmetry on
130 # turn lumping on for more efficient solving
131 #mypde.getSolverOptions().setSolverMethod(mypde.getSolverOptions().LUMPING)
132 kmat = kronecker(domain) # create the kronecker delta function of the domain
133 mypde.setValue(D=rho*kmat) #set the general form value D
134
135 ############################################FIRST TIME STEPS AND SOURCE
136 # define small radius around point xc
137 src_length = 20; print "src_length = ",src_length
138 # set initial values for first two time steps with source terms
139 xb=FunctionOnBoundary(domain).getX()
140 yx=(cos(length(xb-xc)*3.1415/src_length)+1)*whereNegative(length(xb-xc)-src_length)
141 stop=Scalar(0.0,FunctionOnBoundary(domain))
142 stop.setTaggedValue("stop",1.0)
143 src_dir=numpy.array([0.,1.,0.0]) # defines direction of point source as down
144
145 mypde.setValue(y=source[0]*yx*src_dir*stop) #set the source as a function on the boundary
146 # initial value of displacement at point source is constant (U0=0.01)
147 # for first two time steps
148 u=[0.0,0.0,0.0]*wherePositive(x)
149 u_m1=u
150
151 ####################################################ITERATION VARIABLES
152 n=0 # iteration counter
153 t=0 # time counter
154 ##############################################################ITERATION
155 while t<tend:
156 # get current stress
157 g=grad(u); stress=lam*trace(g)*kmat+mu*(g+transpose(g))#*abc
158 mypde.setValue(X=-stress) # set PDE values
159 accel = mypde.getSolution() #get PDE solution for accelleration
160 u_p1=(2.*u-u_m1)+h*h*accel #calculate displacement
161 u_p1=u_p1#*abc # apply boundary conditions
162 u_m1=u; u=u_p1 # shift values by 1
163 # save current displacement, acceleration and pressure
164 if (t >= rtime):
165 saveVTK(os.path.join(savepath,"ex09b.%05d.vtu"%n),displacement=length(u),\
166 acceleration=length(accel),tensor=stress)
167 rtime=rtime+rtime_inc #increment data save time
168 # increment loop values
169 t=t+h; n=n+1
170 if (n < ls):
171 mypde.setValue(y=source[n]*yx*src_dir*stop) #set the source as a function on the boundary
172 print n,"-th time step t ",t

  ViewVC Help
Powered by ViewVC 1.1.26