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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 5288 - (show annotations)
Tue Dec 2 23:18:40 2014 UTC (4 years, 9 months ago) by sshaw
File MIME type: text/x-python
File size: 5660 byte(s)
fixing tests for cases where required domains not built
1 from __future__ import division, print_function
2 ##############################################################################
3 #
4 # Copyright (c) 2009-2014 by University of Queensland
5 # http://www.uq.edu.au
6 #
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 # Development until 2012 by Earth Systems Science Computational Center (ESSCC)
12 # Development 2012-2013 by School of Earth Sciences
13 # Development from 2014 by Centre for Geoscience Computing (GeoComp)
14 #
15 ##############################################################################
16
17 __copyright__="""Copyright (c) 2009-2014 by University of Queensland
18 http://www.uq.edu.au
19 Primary Business: Queensland, Australia"""
20 __license__="""Licensed under the Open Software License version 3.0
21 http://www.opensource.org/licenses/osl-3.0.php"""
22 __url__="https://launchpad.net/escript-finley"
23
24 """
25 Author: Antony Hallam antony.hallam@uqconnect.edu.au
26 """
27
28 ############################################################FILE HEADER
29 # example05a.py
30 # Create either a 2D syncline or anticline model using pycad meshing
31 # tools. Then model steady state heat solution.
32
33 #######################################################EXTERNAL MODULES
34 import matplotlib
35 matplotlib.use('agg') #It's just here for automated testing
36 from esys.pycad import * #domain constructor
37 from esys.pycad.gmsh import Design #Finite Element meshing package
38 import os #file path tool
39 from math import * # math package
40 from esys.escript import *
41 from esys.escript.unitsSI import *
42 from esys.escript.linearPDEs import LinearPDE
43 from cblib import toRegGrid
44 import pylab as pl #Plotting package
45
46 try:
47 # This imports the rectangle domain function
48 from esys.finley import MakeDomain#Converter for escript
49 HAVE_FINLEY = True
50 except ImportError:
51 print("Finley module not available")
52 HAVE_FINLEY = False
53 ########################################################MPI WORLD CHECK
54 if getMPISizeWorld() > 1:
55 import sys
56 print("This example will not run in an MPI world.")
57 sys.exit(0)
58
59 if HAVE_FINLEY:
60 #################################################ESTABLISHING VARIABLES
61 #set modal to 1 for a syncline or -1 for an anticline structural
62 #configuration
63 modal=-1
64
65 # the folder to put our outputs in, leave blank "" for script path -
66 # note this folder path must exist to work
67 save_path= os.path.join("data","example05")
68 mkDir(save_path)
69
70 ################################################ESTABLISHING PARAMETERS
71 #Model Parameters
72 width=5000.0*m #width of model
73 depth=-6000.0*m #depth of model
74 Ttop=20*K # top temperature
75 qin=70*Milli*W/(m*m) # bottom heat influx
76
77 sspl=51 #number of discrete points in spline
78 dsp=width/(sspl-1) #dx of spline steps for width
79 dep_sp=2500.0*m #avg depth of spline
80 h_sp=1500.0*m #heigh of spline
81 orit=-1.0 #orientation of spline 1.0=>up -1.0=>down
82
83 ####################################################DOMAIN CONSTRUCTION
84 # Domain Corners
85 p0=Point(0.0, 0.0, 0.0)
86 p1=Point(0.0, depth, 0.0)
87 p2=Point(width, depth, 0.0)
88 p3=Point(width, 0.0, 0.0)
89
90 # Generate Material Boundary
91 x=[ Point(i*dsp\
92 ,-dep_sp+modal*orit*h_sp*cos(pi*i*dsp/dep_sp+pi))\
93 for i in range(0,sspl)\
94 ]
95 mysp = Spline(*tuple(x))
96 # Start and end of material boundary.
97 x1=mysp.getStartPoint()
98 x2=mysp.getEndPoint()
99
100 # Create TOP BLOCK
101 # lines
102 tbl1=Line(p0,x1)
103 tbl2=mysp
104 tbl3=Line(x2,p3)
105 l30=Line(p3, p0)
106 # curve
107 tblockloop = CurveLoop(tbl1,tbl2,tbl3,l30)
108 # surface
109 tblock = PlaneSurface(tblockloop)
110
111
112 # Create BOTTOM BLOCK
113 # lines
114 bbl1=Line(x1,p1)
115 bbl3=Line(p2,x2)
116 bbl4=-mysp
117 l12=Line(p1, p2)
118 # curve
119 bblockloop = CurveLoop(bbl1,l12,bbl3,bbl4)
120 # surface
121 bblock = PlaneSurface(bblockloop)
122
123 ################################################CREATE MESH FOR ESCRIPT
124 # Create a Design which can make the mesh
125 d=Design(dim=2, element_size=200)
126 # Add the subdomains and flux boundaries.
127 d.addItems(PropertySet("top",tblock),PropertySet("bottom",bblock),\
128 PropertySet("linebottom",l12))
129 # Create the geometry, mesh and Escript domain
130 d.setScriptFileName(os.path.join(save_path,"example05.geo"))
131 d.setMeshFileName(os.path.join(save_path,"example05.msh"))
132 domain=MakeDomain(d, optimizeLabeling=True)
133 print("Domain has been generated ...")
134 ##############################################################SOLVE PDE
135 mypde=LinearPDE(domain)
136 mypde.getSolverOptions().setVerbosityOn()
137 mypde.setSymmetryOn()
138 kappa=Scalar(0,Function(domain))
139 kappa.setTaggedValue("top",2.0*W/m/K)
140 kappa.setTaggedValue("bottom",4.0*W/m/K)
141 mypde.setValue(A=kappa*kronecker(domain))
142 x=Solution(domain).getX()
143 mypde.setValue(q=whereZero(x[1]-sup(x[1])),r=Ttop)
144 qS=Scalar(0,FunctionOnBoundary(domain))
145 qS.setTaggedValue("linebottom",qin)
146 mypde.setValue(y=qS)
147 print("PDE has been generated ...")
148 ###########################################################GET SOLUTION
149 T=mypde.getSolution()
150 print("PDE has been solved ...")
151
152 #######################################################################
153 xi, yi, zi = toRegGrid(T, nx=50, ny=50)
154 pl.matplotlib.pyplot.autumn()
155 pl.contourf(xi,yi,zi,10)
156 pl.xlabel("Horizontal Displacement (m)")
157 pl.ylabel("Depth (m)")
158 pl.savefig(os.path.join(save_path,"Tcontour.png"))
159 print("Solution has been plotted ...")

  ViewVC Help
Powered by ViewVC 1.1.26