/[escript]/trunk/doc/user/finley.tex
ViewVC logotype

Contents of /trunk/doc/user/finley.tex

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1459 - (show annotations)
Thu Mar 27 01:49:10 2008 UTC (11 years, 1 month ago) by ksteube
File MIME type: application/x-tex
File size: 19578 byte(s)
Implemented ParMETIS for mesh partition optimization under MPI.
Documented optimize=True in User Guide for Brick() and Rectangle().

1 %
2 % $Id$
3 %
4 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5 %
6 % Copyright 2003-2007 by ACceSS MNRF
7 % Copyright 2007 by University of Queensland
8 %
9 % http://esscc.uq.edu.au
10 % Primary Business: Queensland, Australia
11 % Licensed under the Open Software License version 3.0
12 % http://www.opensource.org/licenses/osl-3.0.php
13 %
14 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
15 %
16
17 \chapter{ The Module \finley}
18 \label{CHAPTER ON FINLEY}
19
20 \begin{figure}
21 \centerline{\includegraphics[width=\figwidth]{figures/FinleyMesh.eps}}
22 \caption{Subdivision of an Ellipse into triangles order 1 (\finleyelement{Tri3})}
23 \label{FINLEY FIG 0}
24 \end{figure}
25
26 \begin{figure}
27 \centerline{\includegraphics[width=\figwidth]{figures/FinleyContact.eps}}
28 \caption{Mesh around a contact region (\finleyelement{Rec4})}
29 \label{FINLEY FIG 01}
30 \end{figure}
31
32 \declaremodule{extension}{finley} \modulesynopsis{Solving linear, steady partial differential equations using
33 finite elements}
34
35 {\it finley} is a library of C functions solving linear, steady partial differential equations
36 \index{partial differential equations} (PDEs) or systems of PDEs using isoparametrical finite
37 elements \index{FEM!isoparametrical}.
38 It supports unstructured, 1D, 2D and 3D meshes. The module \finley provides an access to the
39 library through the \LinearPDE class of \escript supporting its full functionality. {\it finley}
40 is parallelized using the OpenMP \index{OpenMP} paradigm.
41
42 \section{Formulation}
43
44 For a single PDE with a solution with a single component the linear PDE is defined in the
45 following form:
46 \begin{equation}\label{FINLEY.SINGLE.1}
47 \begin{array}{cl} &
48 \displaystyle{
49 \int\hackscore{\Omega}
50 A\hackscore{jl} \cdot v\hackscore{,j}u\hackscore{,l}+ B\hackscore{j} \cdot v\hackscore{,j} u+ C\hackscore{l} \cdot v u\hackscore{,l}+D \cdot vu \; d\Omega } \\
51 + & \displaystyle{\int\hackscore{\Gamma} d \cdot vu \; d{\Gamma} }
52 + \displaystyle{\int\hackscore{\Gamma^{contact}} d^{contact} \cdot [v][u] \; d{\Gamma} } \\
53 = & \displaystyle{\int\hackscore{\Omega} X\hackscore{j} \cdot v\hackscore{,j}+ Y \cdot v \; d\Omega }\\
54 + & \displaystyle{\int\hackscore{\Gamma} y \cdot v \; d{\Gamma}} +
55 \displaystyle{\int\hackscore{\Gamma^{contact}} y^{contact}\cdot [v] \; d{\Gamma}} \\
56 \end{array}
57 \end{equation}
58
59 \section{Meshes}
60 To understand the usage of \finley one needs to have an understanding of how the finite element meshes
61 \index{FEM!mesh} are defined. \fig{FINLEY FIG 0} shows an example of the
62 subdivision of an ellipse into so called elements \index{FEM!elements} \index{element}.
63 In this case, triangles have been used but other forms of subdivisions
64 can be constructed, e.g. into quadrilaterals or, in the three dimensional case, into tetrahedrons
65 and hexahedrons. The idea of the finite element method is to approximate the solution by a function
66 which is a polynomial of a certain order and is continuous across it boundary to neighbour elements.
67 In the example of \fig{FINLEY FIG 0} a linear polynomial is used on each triangle. As one can see, the triangulation
68 is quite a poor approximation of the ellipse. It can be improved by introducing a midpoint on each element edge then
69 positioning those nodes located on an edge expected to describe the boundary, onto the boundary.
70 In this case the triangle gets a curved edge which requires a parametrization of the triangle using a
71 quadratic polynomial. For this case, the solution is also approximated by a piecewise quadratic polynomial
72 (which explains the name isoparametrical elements), see \Ref{Zienc,NumHand} for more details.
73
74 The union of all elements defines the domain of the PDE.
75 Each element is defined by the nodes used to describe its shape. In \fig{FINLEY FIG 0} the element,
76 which has type \finleyelement{Tri3},
77 with element reference number $19$ \index{element!reference number} is defined by the nodes
78 with reference numbers $9$, $11$ and $0$ \index{node!reference number}. Notice that the order is counterclockwise.
79 The coefficients of the PDE are evaluated at integration nodes with each individual element.
80 For quadrilateral elements a Gauss quadrature scheme is used. In the case of triangular elements a
81 modified form is applied. The boundary of the domain is also subdivided into elements. \index{element!face} In \fig{FINLEY FIG 0}
82 line elements with two nodes are used. The elements are also defined by their describing nodes, e.g.
83 the face element reference number $20$ which has type \finleyelement{Line2} is defined by the nodes
84 with the reference numbers $11$ and $0$. Again the order is crucial, if moving from the first
85 to second node the domain has to lie on the left hand side (in the case of a two dimension surface element
86 the domain has to lie on the left hand side when moving counterclockwise). If the gradient on the
87 surface of the domain is to be calculated rich face elements face to be used. Rich elements on a face
88 are identical to interior elements but with a modified order of nodes such that the 'first' face of the element aligns
89 with the surface of the domain. In \fig{FINLEY FIG 0}
90 elements of the type \finleyelement{Tri3Face} are used.
91 The face element reference number $20$ as a rich face element is defined by the nodes
92 with reference numbers $11$, $0$ and $9$. Notice that the face element $20$ is identical to the
93 interior element $19$ except that, in this case, the order of the node is different to align the first
94 edge of the triangle (which is the edge starting with the first node) with the boundary of the domain.
95
96 Be aware that face elements and elements in the interior of the domain must match, i.e. a face element must be the face
97 of an interior element or, in case of a rich face element, it must be identical to an interior element.
98 If no face elements are specified
99 \finley implicitly assumes homogeneous natural boundary conditions \index{natural boundary conditions!homogeneous},
100 i.e. \var{d}=$0$ and \var{y}=$0$, on the entire boundary of the domain. For
101 inhomogeneous natural boundary conditions \index{natural boundary conditions!inhomogeneous},
102 the boundary must be described by face elements.
103
104 If discontinuities of the PDE solution are considered contact elements
105 \index{element!contact}\index{contact conditions} are introduced to describe the contact region $\Gamma^{contact}$
106 even if $d^{contact}$ and $y^{contact}$ are zero. \fig{FINLEY FIG 01} shows a simple example of a mesh
107 of rectangular elements around a contact region $\Gamma^{contact}$ \index{element!contact}.
108 The contact region is described by the
109 elements $4$, $3$ and $6$. Their element type is \finleyelement{Line2_Contact}.
110 The nodes $9$, $12$, $6$, $5$ define contact element $4$, where the coordinates of nodes $12$ and $5$ and
111 nodes $4$ and $6$ are identical with the idea that nodes $12$ and $9$ are located above and
112 nodes $5$ and $6$ below the contact region.
113 Again, the order of the nodes within an element is crucial. There is also the option of using rich elements
114 if the gradient is to be calculated on the contact region. Similarly to the rich face elements
115 these are constructed from two interior elements by reordering the nodes such that
116 the 'first' face of the element above and the 'first' face of the element below the
117 contact regions line up. The rich version of element
118 $4$ is of type \finleyelement{Rec4Face_Contact} and is defined by the nodes $9$, $12$, $16$, $18$, $6$, $5$, $0$ and
119 $2$.
120
121 \tab{FINLEY TAB 1} shows the interior element types and the corresponding element types to be used
122 on the face and contacts. \fig{FINLEY.FIG:1}, \fig{FINLEY.FIG:2} and \fig{FINLEY.FIG:4} show the ordering of
123 the nodes within an element.
124
125 \begin{table}
126 \begin{tablev}{l|llll}{textrm}{interior}{face}{rich face}{contact}{rich contact}
127 \linev{\finleyelement{Line2}}{\finleyelement{Point1}}{\finleyelement{Line2Face}}{\finleyelement{Point1_Contact}}{\finleyelement{Line2Face_Contact}}
128 \linev{\finleyelement{Line3}}{\finleyelement{Point1}}{\finleyelement{Line3Face}}{\finleyelement{Point1_Contact}}{\finleyelement{Line3Face_Contact}}
129 \linev{\finleyelement{Tri3}}{\finleyelement{Line2}}{\finleyelement{Tri3Face}}{\finleyelement{Line2_Contact}}{\finleyelement{Tri3Face_Contact}}
130 \linev{\finleyelement{Tri6}}{\finleyelement{Line3}}{\finleyelement{Tri6Face}}{\finleyelement{Line3_Contact}}{\finleyelement{Tri6Face_Contact}}
131 \linev{\finleyelement{Rec4}}{\finleyelement{Line2}}{\finleyelement{Rec4Face}}{\finleyelement{Line2_Contact}}{\finleyelement{Rec4Face_Contact}}
132 \linev{\finleyelement{Rec8}}{\finleyelement{Line3}}{\finleyelement{Rec8Face}}{\finleyelement{Line3_Contact}}{\finleyelement{Rec8Face_Contact}}
133 \linev{\finleyelement{Rec9}}{\finleyelement{Line3}}{\finleyelement{Rec9Face}}{\finleyelement{Line3_Contact}}{\finleyelement{Rec9Face_Contact}}
134 \linev{\finleyelement{Tet4}}{\finleyelement{Tri6}}{\finleyelement{Tet4Face}}{\finleyelement{Tri6_Contact}}{\finleyelement{Tet4Face_Contact}}
135 \linev{\finleyelement{Tet10}}{\finleyelement{Tri9}}{\finleyelement{Tet10Face}}{\finleyelement{Tri9_Contact}}{\finleyelement{Tet10Face_Contact}}
136 \linev{\finleyelement{Hex8}}{\finleyelement{Rec4}}{\finleyelement{Hex8Face}}{\finleyelement{Rec4_Contact}}{\finleyelement{Hex8Face_Contact}}
137 \linev{\finleyelement{Hex20}}{\finleyelement{Rec8}}{\finleyelement{Hex20Face}}{\finleyelement{Rec8_Contact}}{\finleyelement{Hex20Face_Contact}}
138 \end{tablev}
139 \caption{Finley elements and corresponding elements to be used on domain faces and contacts.
140 The rich types have to be used if the gradient of function is to be calculated on faces and contacts, respectively.}
141 \label{FINLEY TAB 1}
142 \end{table}
143
144 The native \finley file format is defined as follows.
145 Each node \var{i} has \var{dim} spatial coordinates \var{Node[i]}, a reference number
146 \var{Node_ref[i]}, a degree of freedom \var{Node_DOF[i]} and tag \var{Node_tag[i]}.
147 In most cases \var{Node_DOF[i]}=\var{Node_ref[i]} however, for periodic boundary conditions,
148 \var{Node_DOF[i]} is chosen differently, see example below. The tag can be used to mark nodes sharing
149 the same properties. Element \var{i} is defined by the \var{Element_numNodes} nodes \var{Element_Nodes[i]}
150 which is a list of node reference numbers. The order is crucial.
151 It has a reference number \var{Element_ref[i]} and a tag \var{Element_tag[i]}. The tag
152 can be used to mark elements sharing the same properties. For instance elements above
153 a contact region are marked with $2$ and elements below a contact region are marked with $1$.
154 \var{Element_Type} and \var{Element_Num} give the element type and the number of elements in the mesh.
155 Analogue notations are used for face and contact elements. The following Python script
156 prints the mesh definition in the \finley file format:
157 \begin{python}
158 print "%s\n"%mesh_name
159 # node coordinates:
160 print "%dD-nodes %d\n"%(dim,numNodes)
161 for i in range(numNodes):
162 print "%d %d %d"%(Node_ref[i],Node_DOF[i],Node_tag[i])
163 for j in range(dim): print " %e"%Node[i][j]
164 print "\n"
165 # interior elements
166 print "%s %d\n"%(Element_Type,Element_Num)
167 for i in range(Element_Num):
168 print "%d %d"%(Element_ref[i],Element_tag[i])
169 for j in range(Element_numNodes): print " %d"%Element_Nodes[i][j]
170 print "\n"
171 # face elements
172 print "%s %d\n"%(FaceElement_Type,FaceElement_Num)
173 for i in range(FaceElement_Num):
174 print "%d %d"%(FaceElement_ref[i],FaceElement_tag[i])
175 for j in range(FaceElement_numNodes): print " %d"%FaceElement_Nodes[i][j]
176 print "\n"
177 # contact elements
178 print "%s %d\n"%(ContactElement_Type,ContactElement_Num)
179 for i in range(ContactElement_Num):
180 print "%d %d"%(ContactElement_ref[i],ContactElement_tag[i])
181 for j in range(ContactElement_numNodes): print " %d"%ContactElement_Nodes[i][j]
182 print "\n"
183 # point sources (not supported yet)
184 write("Point1 0",face_element_type,numFaceElements)
185 \end{python}
186
187 The following example of a mesh file defines the mesh shown in \fig{FINLEY FIG 01}:
188 \begin{verbatim}
189 Example 1
190 2D Nodes 16
191 0 0 0 0. 0.
192 2 2 0 0.33 0.
193 3 3 0 0.66 0.
194 7 4 0 1. 0.
195 5 5 0 0. 0.5
196 6 6 0 0.33 0.5
197 8 8 0 0.66 0.5
198 10 10 0 1.0 0.5
199 12 12 0 0. 0.5
200 9 9 0 0.33 0.5
201 13 13 0 0.66 0.5
202 15 15 0 1.0 0.5
203 16 16 0 0. 1.0
204 18 18 0 0.33 1.0
205 19 19 0 0.66 1.0
206 20 20 0 1.0 1.0
207 Rec4 6
208 0 1 0 2 6 5
209 1 1 2 3 8 6
210 2 1 3 7 10 8
211 5 2 12 9 18 16
212 7 2 13 19 18 9
213 10 2 20 19 13 15
214 Line2 0
215 Line2_Contact 3
216 4 0 9 12 6 5
217 3 0 13 9 8 6
218 6 0 15 13 10 8
219 Point1 0
220 \end{verbatim}
221 Notice that the order in which the nodes and elements are given is arbitrary.
222 In the case that rich contact elements are used the contact element section gets
223 the form
224 \begin{verbatim}
225 Rec4Face_Contact 3
226 4 0 9 12 16 18 6 5 0 2
227 3 0 13 9 18 19 8 6 2 3
228 6 0 15 13 19 20 10 8 3 7
229 \end{verbatim}
230 Periodic boundary condition \index{boundary conditions!periodic} can be introduced by altering \var{Node_DOF}.
231 It allows identification of nodes even if they have different physical locations. For instance, to
232 enforce periodic boundary conditions at the face $x_0=0$ and $x_0=1$ one identifies
233 the degrees of freedom for nodes $0$, $5$, $12$ and $16$ with the degrees of freedom for
234 $7$, $10$, $15$ and $20$, respectively. The node section of the \finley mesh gets now the form:
235 \begin{verbatim}
236 2D Nodes 16
237 0 0 0 0. 0.
238 2 2 0 0.33 0.
239 3 3 0 0.66 0.
240 7 0 0 1. 0.
241 5 5 0 0. 0.5
242 6 6 0 0.33 0.5
243 8 8 0 0.66 0.5
244 10 5 0 1.0 0.5
245 12 12 0 0. 0.5
246 9 9 0 0.33 0.5
247 13 13 0 0.66 0.5
248 15 12 0 1.0 0.5
249 16 16 0 0. 1.0
250 18 18 0 0.33 1.0
251 19 19 0 0.66 1.0
252 20 16 0 1.0 1.0
253 \end{verbatim}
254
255
256 \include{finleyelements}
257
258 \subsection{Linear Solvers in \LinearPDE}
259 Currently \finley supports the linear solvers \PCG, \GMRES, \PRESTWENTY and \BiCGStab.
260 For \GMRES the options \var{truncation} and \var{restart} of the \method{getSolution} can be
261 used to control the truncation and restart during iteration. Default values are
262 \var{truncation}=5 and \var{restart}=20.
263 The default solver is \BiCGStab but if the symmetry flag is set \PCG is the default solver.
264 \finley supports the solver options \var{iter_max} which specifies the maximum number of iterations steps,
265 \var{verbose}=\True or \False and \var{preconditioner}=\constant{JACOBI} or \constant {ILU0}.
266 In some installations \finley supports the \Direct solver and the
267 solver options \var{reordering}=\constant{util.NO_REORDERING},
268 \constant{util.MINIMUM_FILL_IN} or \constant{util.NESTED_DISSECTION} (default is \constant{util.NO_REORDERING}),
269 \var{drop_tolerance} specifying the threshold for values to be dropped in the
270 incomplete elimination process (default is 0.01) and \var{drop_storage} specifying the maximum increase
271 in storage allowed in the
272 incomplete elimination process (default is 1.20).
273
274 \subsection{Functions}
275 \begin{funcdesc}{Mesh}{fileName,integrationOrder=-1}
276 creates a \Domain object form the FEM mesh defined in
277 file \var{fileName}. The file must be given the \finley file format.
278 If \var{integrationOrder} is positive, a numerical integration scheme
279 chosen which is accurate on each element up to a polynomial of
280 degree \var{integrationOrder} \index{integration order}. Otherwise
281 an appropriate integration order is chosen independently.
282 \end{funcdesc}
283
284 \begin{funcdesc}{Rectangle}{n0,n1,order=1,l0=1.,l1=1., integrationOrder=-1, \\
285 periodic0=\False,periodic1=\False,useElementsOnFace=\False,optimize=\False}
286 Generates a \Domain object representing a two dimensional rectangle between
287 $(0,0)$ and $(l0,l1)$ with orthogonal edges. The rectangle is filled with
288 \var{n0} elements along the $x_0$-axis and
289 \var{n1} elements along the $x_1$-axis.
290 For \var{order}=1 and \var{order}=2
291 \finleyelement{Rec4} and
292 \finleyelement{Rec8} are used, respectively.
293 In the case of \var{useElementsOnFace}=\False,
294 \finleyelement{Line2} and
295 \finleyelement{Line3} are used to subdivide the edges of the rectangle, respectively.
296 In the case of \var{useElementsOnFace}=\True (this option should be used if gradients
297 are calculated on domain faces),
298 \finleyelement{Rec4Face} and
299 \finleyelement{Rec8Face} are used on the edges, respectively.
300 If \var{integrationOrder} is positive, a numerical integration scheme
301 chosen which is accurate on each element up to a polynomial of
302 degree \var{integrationOrder} \index{integration order}. Otherwise
303 an appropriate integration order is chosen independently. If
304 \var{periodic0}=\True, periodic boundary conditions \index{periodic boundary conditions}
305 along the $x_0$-directions are enforced. That means when for any solution of a PDE solved by \finley
306 the value on the line $x_0=0$ will be identical to the values on $x_0=\var{l0}$.
307 Correspondingly,
308 \var{periodic1}=\False sets periodic boundary conditions
309 in $x_1$-direction.
310 If \var{optimize}=\True mesh node relabeling will be attempted to reduce the computation and also ParMETIS will be used to improve the mesh partition if running on multiple CPUs with MPI.
311 \end{funcdesc}
312
313 \begin{funcdesc}{Brick}{n0,n1,n2,order=1,l0=1.,l1=1.,l2=1., integrationOrder=-1, \\
314 periodic0=\False,periodic1=\False,periodic2=\False,useElementsOnFace=\False,optimize=\False}
315 Generates a \Domain object representing a three dimensional brick between
316 $(0,0,0)$ and $(l0,l1,l2)$ with orthogonal faces. The brick is filled with
317 \var{n0} elements along the $x_0$-axis,
318 \var{n1} elements along the $x_1$-axis and
319 \var{n2} elements along the $x_2$-axis.
320 For \var{order}=1 and \var{order}=2
321 \finleyelement{Hex8} and
322 \finleyelement{Hex20} are used, respectively.
323 In the case of \var{useElementsOnFace}=\False,
324 \finleyelement{Rec4} and
325 \finleyelement{Rec8} are used to subdivide the faces of the brick, respectively.
326 In the case of \var{useElementsOnFace}=\True (this option should be used if gradients
327 are calculated on domain faces),
328 \finleyelement{Hex8Face} and
329 \finleyelement{Hex20Face} are used on the brick faces, respectively.
330 If \var{integrationOrder} is positive, a numerical integration scheme
331 chosen which is accurate on each element up to a polynomial of
332 degree \var{integrationOrder} \index{integration order}. Otherwise
333 an appropriate integration order is chosen independently. If
334 \var{periodic0}=\True, periodic boundary conditions \index{periodic boundary conditions}
335 along the $x_0$-directions are enforced. That means when for any solution of a PDE solved by \finley
336 the value on the plane $x_0=0$ will be identical to the values on $x_0=\var{l0}$. Correspondingly,
337 \var{periodic1}=\False and \var{periodic2}=\False sets periodic boundary conditions
338 in $x_1$-direction and $x_2$-direction, respectively.
339 If \var{optimize}=\True mesh node relabeling will be attempted to reduce the computation and also ParMETIS will be used to improve the mesh partition if running on multiple CPUs with MPI.
340 \end{funcdesc}
341
342 \begin{funcdesc}{GlueFaces}{meshList,safetyFactor=0.2,tolerance=1.e-13}
343 Generates a new \Domain object from the list \var{meshList} of \finley meshes.
344 Nodes in face elements whose difference of coordinates is less then \var{tolerance} times the
345 diameter of the domain are merged. The corresponding face elements are removed from the mesh.
346
347 TODO: explain \var{safetyFactor} and show an example.
348 \end{funcdesc}
349
350 \begin{funcdesc}{JoinFaces}{meshList,safetyFactor=0.2,tolerance=1.e-13}
351 Generates a new \Domain object from the list \var{meshList} of \finley meshes.
352 Face elements whose nodes coordinates have difference is less then \var{tolerance} times the
353 diameter of the domain are combined to form a contact element \index{element!contact}
354 The corresponding face elements are removed from the mesh.
355
356 TODO: explain \var{safetyFactor} and show an example.
357 \end{funcdesc}

Properties

Name Value
svn:eol-style native
svn:keywords Author Date Id Revision

  ViewVC Help
Powered by ViewVC 1.1.26