1 |
|
2 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
3 |
% |
4 |
% Copyright (c) 2003-2009 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 |
|
15 |
\chapter{ The Module \finley} |
16 |
\label{CHAPTER ON FINLEY} |
17 |
|
18 |
\begin{figure} |
19 |
\centerline{\includegraphics[width=\figwidth]{figures/FinleyMesh}} |
20 |
\caption{Subdivision of an Ellipse into triangles order 1 (\finleyelement{Tri3})} |
21 |
\label{FINLEY FIG 0} |
22 |
\end{figure} |
23 |
|
24 |
\begin{figure} |
25 |
\centerline{\includegraphics[width=\figwidth]{figures/FinleyContact}} |
26 |
\caption{Mesh around a contact region (\finleyelement{Rec4})} |
27 |
\label{FINLEY FIG 01} |
28 |
\end{figure} |
29 |
|
30 |
\declaremodule{extension}{finley} \modulesynopsis{Solving linear, steady partial differential equations using |
31 |
finite elements} |
32 |
|
33 |
{\it finley} is a library of C functions solving linear, steady partial differential equations |
34 |
\index{partial differential equations} (PDEs) or systems of PDEs using isoparametrical finite |
35 |
elements \index{FEM!isoparametrical}. |
36 |
It supports unstructured, 1D, 2D and 3D meshes. The module \finley provides an access to the |
37 |
library through the \LinearPDE class of \escript supporting its full functionality. {\it finley} |
38 |
is parallelized using the OpenMP \index{OpenMP} paradigm. |
39 |
|
40 |
\section{Formulation} |
41 |
|
42 |
For a single PDE with a solution with a single component the linear PDE is defined in the |
43 |
following form: |
44 |
\begin{equation}\label{FINLEY.SINGLE.1} |
45 |
\begin{array}{cl} & |
46 |
\displaystyle{ |
47 |
\int\hackscore{\Omega} |
48 |
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 } \\ |
49 |
+ & \displaystyle{\int\hackscore{\Gamma} d \cdot vu \; d{\Gamma} } |
50 |
+ \displaystyle{\int\hackscore{\Gamma^{contact}} d^{contact} \cdot [v][u] \; d{\Gamma} } \\ |
51 |
= & \displaystyle{\int\hackscore{\Omega} X\hackscore{j} \cdot v\hackscore{,j}+ Y \cdot v \; d\Omega }\\ |
52 |
+ & \displaystyle{\int\hackscore{\Gamma} y \cdot v \; d{\Gamma}} + |
53 |
\displaystyle{\int\hackscore{\Gamma^{contact}} y^{contact}\cdot [v] \; d{\Gamma}} \\ |
54 |
\end{array} |
55 |
\end{equation} |
56 |
|
57 |
\section{Meshes} |
58 |
To understand the usage of \finley one needs to have an understanding of how the finite element meshes |
59 |
\index{FEM!mesh} are defined. \fig{FINLEY FIG 0} shows an example of the |
60 |
subdivision of an ellipse into so called elements \index{FEM!elements} \index{element}. |
61 |
In this case, triangles have been used but other forms of subdivisions |
62 |
can be constructed, e.g. into quadrilaterals or, in the three dimensional case, into tetrahedrons |
63 |
and hexahedrons. The idea of the finite element method is to approximate the solution by a function |
64 |
which is a polynomial of a certain order and is continuous across it boundary to neighbor elements. |
65 |
In the example of \fig{FINLEY FIG 0} a linear polynomial is used on each triangle. As one can see, the triangulation |
66 |
is quite a poor approximation of the ellipse. It can be improved by introducing a midpoint on each element edge then |
67 |
positioning those nodes located on an edge expected to describe the boundary, onto the boundary. |
68 |
In this case the triangle gets a curved edge which requires a parametrization of the triangle using a |
69 |
quadratic polynomial. For this case, the solution is also approximated by a piecewise quadratic polynomial |
70 |
(which explains the name isoparametrical elements), see \Ref{Zienc,NumHand} for more details. |
71 |
|
72 |
The union of all elements defines the domain of the PDE. |
73 |
Each element is defined by the nodes used to describe its shape. In \fig{FINLEY FIG 0} the element, |
74 |
which has type \finleyelement{Tri3}, |
75 |
with element reference number $19$ \index{element!reference number} is defined by the nodes |
76 |
with reference numbers $9$, $11$ and $0$ \index{node!reference number}. Notice that the order is counterclockwise. |
77 |
The coefficients of the PDE are evaluated at integration nodes with each individual element. |
78 |
For quadrilateral elements a Gauss quadrature scheme is used. In the case of triangular elements a |
79 |
modified form is applied. The boundary of the domain is also subdivided into elements. \index{element!face} In \fig{FINLEY FIG 0} |
80 |
line elements with two nodes are used. The elements are also defined by their describing nodes, e.g. |
81 |
the face element reference number $20$ which has type \finleyelement{Line2} is defined by the nodes |
82 |
with the reference numbers $11$ and $0$. Again the order is crucial, if moving from the first |
83 |
to second node the domain has to lie on the left hand side (in the case of a two dimension surface element |
84 |
the domain has to lie on the left hand side when moving counterclockwise). If the gradient on the |
85 |
surface of the domain is to be calculated rich face elements face to be used. Rich elements on a face |
86 |
are identical to interior elements but with a modified order of nodes such that the 'first' face of the element aligns |
87 |
with the surface of the domain. In \fig{FINLEY FIG 0} |
88 |
elements of the type \finleyelement{Tri3Face} are used. |
89 |
The face element reference number $20$ as a rich face element is defined by the nodes |
90 |
with reference numbers $11$, $0$ and $9$. Notice that the face element $20$ is identical to the |
91 |
interior element $19$ except that, in this case, the order of the node is different to align the first |
92 |
edge of the triangle (which is the edge starting with the first node) with the boundary of the domain. |
93 |
|
94 |
Be aware that face elements and elements in the interior of the domain must match, i.e. a face element must be the face |
95 |
of an interior element or, in case of a rich face element, it must be identical to an interior element. |
96 |
If no face elements are specified |
97 |
\finley implicitly assumes homogeneous natural boundary conditions \index{natural boundary conditions!homogeneous}, |
98 |
i.e. \var{d}=$0$ and \var{y}=$0$, on the entire boundary of the domain. For |
99 |
inhomogeneous natural boundary conditions \index{natural boundary conditions!inhomogeneous}, |
100 |
the boundary must be described by face elements. |
101 |
|
102 |
If discontinuities of the PDE solution are considered contact elements |
103 |
\index{element!contact}\index{contact conditions} are introduced to describe the contact region $\Gamma^{contact}$ |
104 |
even if $d^{contact}$ and $y^{contact}$ are zero. \fig{FINLEY FIG 01} shows a simple example of a mesh |
105 |
of rectangular elements around a contact region $\Gamma^{contact}$ \index{element!contact}. |
106 |
The contact region is described by the |
107 |
elements $4$, $3$ and $6$. Their element type is \finleyelement{Line2_Contact}. |
108 |
The nodes $9$, $12$, $6$, $5$ define contact element $4$, where the coordinates of nodes $12$ and $5$ and |
109 |
nodes $4$ and $6$ are identical with the idea that nodes $12$ and $9$ are located above and |
110 |
nodes $5$ and $6$ below the contact region. |
111 |
Again, the order of the nodes within an element is crucial. There is also the option of using rich elements |
112 |
if the gradient is to be calculated on the contact region. Similarly to the rich face elements |
113 |
these are constructed from two interior elements by reordering the nodes such that |
114 |
the 'first' face of the element above and the 'first' face of the element below the |
115 |
contact regions line up. The rich version of element |
116 |
$4$ is of type \finleyelement{Rec4Face_Contact} and is defined by the nodes $9$, $12$, $16$, $18$, $6$, $5$, $0$ and |
117 |
$2$. |
118 |
|
119 |
\tab{FINLEY TAB 1} shows the interior element types and the corresponding element types to be used |
120 |
on the face and contacts. \fig{FINLEY.FIG:1}, \fig{FINLEY.FIG:2} and \fig{FINLEY.FIG:4} show the ordering of |
121 |
the nodes within an element. |
122 |
|
123 |
\begin{table} |
124 |
\begin{tablev}{l|llll}{textrm}{interior}{face}{rich face}{contact}{rich contact} |
125 |
\linev{\finleyelement{Line2}}{\finleyelement{Point1}}{\finleyelement{Line2Face}}{\finleyelement{Point1_Contact}}{\finleyelement{Line2Face_Contact}} |
126 |
\linev{\finleyelement{Line3}}{\finleyelement{Point1}}{\finleyelement{Line3Face}}{\finleyelement{Point1_Contact}}{\finleyelement{Line3Face_Contact}} |
127 |
\linev{\finleyelement{Tri3}}{\finleyelement{Line2}}{\finleyelement{Tri3Face}}{\finleyelement{Line2_Contact}}{\finleyelement{Tri3Face_Contact}} |
128 |
\linev{\finleyelement{Tri6}}{\finleyelement{Line3}}{\finleyelement{Tri6Face}}{\finleyelement{Line3_Contact}}{\finleyelement{Tri6Face_Contact}} |
129 |
\linev{\finleyelement{Rec4}}{\finleyelement{Line2}}{\finleyelement{Rec4Face}}{\finleyelement{Line2_Contact}}{\finleyelement{Rec4Face_Contact}} |
130 |
\linev{\finleyelement{Rec8}}{\finleyelement{Line3}}{\finleyelement{Rec8Face}}{\finleyelement{Line3_Contact}}{\finleyelement{Rec8Face_Contact}} |
131 |
\linev{\finleyelement{Rec9}}{\finleyelement{Line3}}{\finleyelement{Rec9Face}}{\finleyelement{Line3_Contact}}{\finleyelement{Rec9Face_Contact}} |
132 |
\linev{\finleyelement{Tet4}}{\finleyelement{Tri6}}{\finleyelement{Tet4Face}}{\finleyelement{Tri6_Contact}}{\finleyelement{Tet4Face_Contact}} |
133 |
\linev{\finleyelement{Tet10}}{\finleyelement{Tri9}}{\finleyelement{Tet10Face}}{\finleyelement{Tri9_Contact}}{\finleyelement{Tet10Face_Contact}} |
134 |
\linev{\finleyelement{Hex8}}{\finleyelement{Rec4}}{\finleyelement{Hex8Face}}{\finleyelement{Rec4_Contact}}{\finleyelement{Hex8Face_Contact}} |
135 |
\linev{\finleyelement{Hex20}}{\finleyelement{Rec8}}{\finleyelement{Hex20Face}}{\finleyelement{Rec8_Contact}}{\finleyelement{Hex20Face_Contact}} |
136 |
\end{tablev} |
137 |
\caption{Finley elements and corresponding elements to be used on domain faces and contacts. |
138 |
The rich types have to be used if the gradient of function is to be calculated on faces and contacts, respectively.} |
139 |
\label{FINLEY TAB 1} |
140 |
\end{table} |
141 |
|
142 |
The native \finley file format is defined as follows. |
143 |
Each node \var{i} has \var{dim} spatial coordinates \var{Node[i]}, a reference number |
144 |
\var{Node_ref[i]}, a degree of freedom \var{Node_DOF[i]} and tag \var{Node_tag[i]}. |
145 |
In most cases \var{Node_DOF[i]}=\var{Node_ref[i]} however, for periodic boundary conditions, |
146 |
\var{Node_DOF[i]} is chosen differently, see example below. The tag can be used to mark nodes sharing |
147 |
the same properties. Element \var{i} is defined by the \var{Element_numNodes} nodes \var{Element_Nodes[i]} |
148 |
which is a list of node reference numbers. The order is crucial. |
149 |
It has a reference number \var{Element_ref[i]} and a tag \var{Element_tag[i]}. The tag |
150 |
can be used to mark elements sharing the same properties. For instance elements above |
151 |
a contact region are marked with $2$ and elements below a contact region are marked with $1$. |
152 |
\var{Element_Type} and \var{Element_Num} give the element type and the number of elements in the mesh. |
153 |
Analogue notations are used for face and contact elements. The following Python script |
154 |
prints the mesh definition in the \finley file format: |
155 |
\begin{python} |
156 |
print "%s\n"%mesh_name |
157 |
# node coordinates: |
158 |
print "%dD-nodes %d\n"%(dim,numNodes) |
159 |
for i in range(numNodes): |
160 |
print "%d %d %d"%(Node_ref[i],Node_DOF[i],Node_tag[i]) |
161 |
for j in range(dim): print " %e"%Node[i][j] |
162 |
print "\n" |
163 |
# interior elements |
164 |
print "%s %d\n"%(Element_Type,Element_Num) |
165 |
for i in range(Element_Num): |
166 |
print "%d %d"%(Element_ref[i],Element_tag[i]) |
167 |
for j in range(Element_numNodes): print " %d"%Element_Nodes[i][j] |
168 |
print "\n" |
169 |
# face elements |
170 |
print "%s %d\n"%(FaceElement_Type,FaceElement_Num) |
171 |
for i in range(FaceElement_Num): |
172 |
print "%d %d"%(FaceElement_ref[i],FaceElement_tag[i]) |
173 |
for j in range(FaceElement_numNodes): print " %d"%FaceElement_Nodes[i][j] |
174 |
print "\n" |
175 |
# contact elements |
176 |
print "%s %d\n"%(ContactElement_Type,ContactElement_Num) |
177 |
for i in range(ContactElement_Num): |
178 |
print "%d %d"%(ContactElement_ref[i],ContactElement_tag[i]) |
179 |
for j in range(ContactElement_numNodes): print " %d"%ContactElement_Nodes[i][j] |
180 |
print "\n" |
181 |
# point sources (not supported yet) |
182 |
write("Point1 0",face_element_type,numFaceElements) |
183 |
\end{python} |
184 |
|
185 |
The following example of a mesh file defines the mesh shown in \fig{FINLEY FIG 01}: |
186 |
\begin{verbatim} |
187 |
Example 1 |
188 |
2D Nodes 16 |
189 |
0 0 0 0. 0. |
190 |
2 2 0 0.33 0. |
191 |
3 3 0 0.66 0. |
192 |
7 4 0 1. 0. |
193 |
5 5 0 0. 0.5 |
194 |
6 6 0 0.33 0.5 |
195 |
8 8 0 0.66 0.5 |
196 |
10 10 0 1.0 0.5 |
197 |
12 12 0 0. 0.5 |
198 |
9 9 0 0.33 0.5 |
199 |
13 13 0 0.66 0.5 |
200 |
15 15 0 1.0 0.5 |
201 |
16 16 0 0. 1.0 |
202 |
18 18 0 0.33 1.0 |
203 |
19 19 0 0.66 1.0 |
204 |
20 20 0 1.0 1.0 |
205 |
Rec4 6 |
206 |
0 1 0 2 6 5 |
207 |
1 1 2 3 8 6 |
208 |
2 1 3 7 10 8 |
209 |
5 2 12 9 18 16 |
210 |
7 2 13 19 18 9 |
211 |
10 2 20 19 13 15 |
212 |
Line2 0 |
213 |
Line2_Contact 3 |
214 |
4 0 9 12 6 5 |
215 |
3 0 13 9 8 6 |
216 |
6 0 15 13 10 8 |
217 |
Point1 0 |
218 |
\end{verbatim} |
219 |
Notice that the order in which the nodes and elements are given is arbitrary. |
220 |
In the case that rich contact elements are used the contact element section gets |
221 |
the form |
222 |
\begin{verbatim} |
223 |
Rec4Face_Contact 3 |
224 |
4 0 9 12 16 18 6 5 0 2 |
225 |
3 0 13 9 18 19 8 6 2 3 |
226 |
6 0 15 13 19 20 10 8 3 7 |
227 |
\end{verbatim} |
228 |
Periodic boundary condition \index{boundary conditions!periodic} can be introduced by altering \var{Node_DOF}. |
229 |
It allows identification of nodes even if they have different physical locations. For instance, to |
230 |
enforce periodic boundary conditions at the face $x_0=0$ and $x_0=1$ one identifies |
231 |
the degrees of freedom for nodes $0$, $5$, $12$ and $16$ with the degrees of freedom for |
232 |
$7$, $10$, $15$ and $20$, respectively. The node section of the \finley mesh gets now the form: |
233 |
\begin{verbatim} |
234 |
2D Nodes 16 |
235 |
0 0 0 0. 0. |
236 |
2 2 0 0.33 0. |
237 |
3 3 0 0.66 0. |
238 |
7 0 0 1. 0. |
239 |
5 5 0 0. 0.5 |
240 |
6 6 0 0.33 0.5 |
241 |
8 8 0 0.66 0.5 |
242 |
10 5 0 1.0 0.5 |
243 |
12 12 0 0. 0.5 |
244 |
9 9 0 0.33 0.5 |
245 |
13 13 0 0.66 0.5 |
246 |
15 12 0 1.0 0.5 |
247 |
16 16 0 0. 1.0 |
248 |
18 18 0 0.33 1.0 |
249 |
19 19 0 0.66 1.0 |
250 |
20 16 0 1.0 1.0 |
251 |
\end{verbatim} |
252 |
|
253 |
\clearpage |
254 |
\input{finleyelements} |
255 |
\clearpage |
256 |
|
257 |
|
258 |
\begin{table} |
259 |
{\scriptsize |
260 |
\begin{tabular}{l||c|c|c|c|c|c|c|c} |
261 |
\member{setSolverMethod} & \member{DIRECT}& \member{PCG} & \member{GMRES} & \member{TFQMR} & \member{MINRES} & \member{PRES20} & \member{BICGSTAB} & \member{LUMPING} \\ |
262 |
\hline |
263 |
\hline |
264 |
\member{setReordering} & $\checkmark$ & & & & & &\\ |
265 |
\hline \member{setRestart} & & & $\checkmark$ & & & $20$ & \\ |
266 |
\hline\member{setTruncation} & & & $\checkmark$ & & & $5$ & \\ |
267 |
\hline\member{setIterMax} & & $\checkmark$& $\checkmark$ & $\checkmark$& $\checkmark$& $\checkmark$ & $\checkmark$ \\ |
268 |
\hline\member{setTolerance} & & $\checkmark$& $\checkmark$ & $\checkmark$& $\checkmark$& $\checkmark$ & $\checkmark$ \\ |
269 |
\hline\member{setAbsoluteTolerance} & & $\checkmark$& $\checkmark$ & $\checkmark$& $\checkmark$& $\checkmark$ & $\checkmark$ \\ |
270 |
\hline\member{setReordering} & $\checkmark$ & & & & & & & \\ |
271 |
\end{tabular} |
272 |
} |
273 |
\caption{Solvers available for |
274 |
\finley |
275 |
and the \PASO package and the relevant options in \class{SolverOptions}. |
276 |
\MKL supports |
277 |
\MINIMUMFILLIN |
278 |
and |
279 |
\NESTEDDESCTION |
280 |
reordering. |
281 |
Currently the \UMFPACK interface does not support any reordering. |
282 |
\label{TAB FINLEY SOLVER OPTIONS 1} } |
283 |
\end{table} |
284 |
|
285 |
\begin{table} |
286 |
{\scriptsize |
287 |
\begin{tabular}{l||c|c|c|c|c|c|c|c} |
288 |
\member{setPreconditioner} & |
289 |
\member{NO_PRECONDITIONER} & |
290 |
\member{AMG} & |
291 |
\member{JACOBI} & |
292 |
\member{GAUSS_SEIDEL}& |
293 |
\member{REC_ILU}& |
294 |
\member{RILU} & |
295 |
\member{ILU0} & |
296 |
\member{DIRECT} \\ |
297 |
\hline |
298 |
status: & |
299 |
later & |
300 |
later & |
301 |
$\checkmark$ & |
302 |
$\checkmark$& |
303 |
$\checkmark$ & |
304 |
later & |
305 |
$\checkmark$ & |
306 |
later \\ |
307 |
\hline |
308 |
\hline |
309 |
\member{setCoarsening}& |
310 |
& |
311 |
$\checkmark$ & |
312 |
& |
313 |
& |
314 |
& |
315 |
& |
316 |
& |
317 |
\\ |
318 |
|
319 |
|
320 |
\hline\member{setLevelMax}& |
321 |
& |
322 |
$\checkmark$ & |
323 |
& |
324 |
& |
325 |
& |
326 |
& |
327 |
& |
328 |
\\ |
329 |
|
330 |
\hline\member{setCoarseningThreshold}& |
331 |
& |
332 |
$\checkmark$ & |
333 |
& |
334 |
& |
335 |
& |
336 |
& |
337 |
& |
338 |
\\ |
339 |
|
340 |
\hline\member{setMinCoarseMatrixSize} & |
341 |
& |
342 |
$\checkmark$ & |
343 |
& |
344 |
& |
345 |
& |
346 |
& |
347 |
& |
348 |
\\ |
349 |
|
350 |
\hline\member{setNumSweeps} & |
351 |
& |
352 |
& |
353 |
$\checkmark$ & |
354 |
$\checkmark$ & |
355 |
& |
356 |
& |
357 |
& |
358 |
\\ |
359 |
|
360 |
\hline\member{setNumPreSweeps}& |
361 |
& |
362 |
$\checkmark$ & |
363 |
& |
364 |
& |
365 |
& |
366 |
& |
367 |
& |
368 |
\\ |
369 |
|
370 |
\hline\member{setNumPostSweeps} & |
371 |
& |
372 |
$\checkmark$ & |
373 |
& |
374 |
& |
375 |
& |
376 |
& |
377 |
& |
378 |
\\ |
379 |
|
380 |
\hline\member{setInnerTolerance}& |
381 |
& |
382 |
& |
383 |
& |
384 |
& |
385 |
& |
386 |
& |
387 |
& |
388 |
\\ |
389 |
|
390 |
\hline\member{setDropTolerance}& |
391 |
& |
392 |
& |
393 |
& |
394 |
& |
395 |
& |
396 |
& |
397 |
& |
398 |
\\ |
399 |
|
400 |
\hline\member{setDropStorage}& |
401 |
& |
402 |
& |
403 |
& |
404 |
& |
405 |
& |
406 |
& |
407 |
& |
408 |
\\ |
409 |
|
410 |
\hline\member{setRelaxationFactor}& |
411 |
& |
412 |
& |
413 |
& |
414 |
& |
415 |
& |
416 |
$\checkmark$ & |
417 |
& |
418 |
\\ |
419 |
|
420 |
\hline\member{adaptInnerTolerance}& |
421 |
& |
422 |
& |
423 |
& |
424 |
& |
425 |
& |
426 |
& |
427 |
& |
428 |
\\ |
429 |
|
430 |
\hline\member{setInnerIterMax}& |
431 |
& |
432 |
& |
433 |
& |
434 |
& |
435 |
& |
436 |
& |
437 |
& |
438 |
\\ |
439 |
\end{tabular} |
440 |
} |
441 |
\caption{Preconditioners available for \finley and the \PASO package and the relevant options in \class{SolverOptions}. \label{TAB FINLEY SOLVER OPTIONS 2}} |
442 |
\end{table} |
443 |
|
444 |
\subsection{Linear Solvers in \SolverOptions} |
445 |
Table~\ref{TAB FINLEY SOLVER OPTIONS 1} and |
446 |
Table~\ref{TAB FINLEY SOLVER OPTIONS 2} show the solvers and preconditioners supported by |
447 |
\finley through the \PASO library. Currently direct solvers are not supported under MPI. |
448 |
By default, \finley is using the iterative solvers \PCG for symmetric and \BiCGStab for non-symmetric problems. |
449 |
If the direct solver is selected which can be useful when solving very ill-posedequations |
450 |
\finley uses the \MKL solver package. If \MKL is not available \UMFPACK is used. If \UMFPACK is not available |
451 |
a suitable iterative solver from the \PASO is used. |
452 |
|
453 |
\subsection{Functions} |
454 |
\begin{funcdesc}{ReadMesh}{fileName,integrationOrder=-1} |
455 |
creates a \Domain object form the FEM mesh defined in |
456 |
file \var{fileName}. The file must be given the \finley file format. |
457 |
If \var{integrationOrder} is positive, a numerical integration scheme |
458 |
chosen which is accurate on each element up to a polynomial of |
459 |
degree \var{integrationOrder} \index{integration order}. Otherwise |
460 |
an appropriate integration order is chosen independently. |
461 |
\end{funcdesc} |
462 |
|
463 |
\begin{funcdesc}{load}{fileName} |
464 |
recovers a \Domain object from a dump file created by the \ |
465 |
eateseates a \Domain object form the FEM mesh defined in |
466 |
file \var{fileName}. The file must be given the \finley file format. |
467 |
If \var{integrationOrder} is positive, a numerical integration scheme |
468 |
chosen which is accurate on each element up to a polynomial of |
469 |
degree \var{integrationOrder} \index{integration order}. Otherwise |
470 |
an appropriate integration order is chosen independently. |
471 |
\end{funcdesc} |
472 |
|
473 |
\begin{funcdesc}{Rectangle}{n0,n1,order=1,l0=1.,l1=1., integrationOrder=-1, \\ |
474 |
periodic0=\False,periodic1=\False,useElementsOnFace=\False,optimize=\False} |
475 |
Generates a \Domain object representing a two dimensional rectangle between |
476 |
$(0,0)$ and $(l0,l1)$ with orthogonal edges. The rectangle is filled with |
477 |
\var{n0} elements along the $x_0$-axis and |
478 |
\var{n1} elements along the $x_1$-axis. |
479 |
For \var{order}=1 and \var{order}=2 |
480 |
\finleyelement{Rec4} and |
481 |
\finleyelement{Rec8} are used, respectively. |
482 |
In the case of \var{useElementsOnFace}=\False, |
483 |
\finleyelement{Line2} and |
484 |
\finleyelement{Line3} are used to subdivide the edges of the rectangle, respectively. |
485 |
In the case of \var{useElementsOnFace}=\True (this option should be used if gradients |
486 |
are calculated on domain faces), |
487 |
\finleyelement{Rec4Face} and |
488 |
\finleyelement{Rec8Face} are used on the edges, respectively. |
489 |
If \var{integrationOrder} is positive, a numerical integration scheme |
490 |
chosen which is accurate on each element up to a polynomial of |
491 |
degree \var{integrationOrder} \index{integration order}. Otherwise |
492 |
an appropriate integration order is chosen independently. If |
493 |
\var{periodic0}=\True, periodic boundary conditions \index{periodic boundary conditions} |
494 |
along the $x_0$-directions are enforced. That means when for any solution of a PDE solved by \finley |
495 |
the value on the line $x_0=0$ will be identical to the values on $x_0=\var{l0}$. |
496 |
Correspondingly, |
497 |
\var{periodic1}=\False sets periodic boundary conditions |
498 |
in $x_1$-direction. |
499 |
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. |
500 |
\end{funcdesc} |
501 |
|
502 |
\begin{funcdesc}{Brick}{n0,n1,n2,order=1,l0=1.,l1=1.,l2=1., integrationOrder=-1, \\ |
503 |
periodic0=\False,periodic1=\False,periodic2=\False,useElementsOnFace=\False,optimize=\False} |
504 |
Generates a \Domain object representing a three dimensional brick between |
505 |
$(0,0,0)$ and $(l0,l1,l2)$ with orthogonal faces. The brick is filled with |
506 |
\var{n0} elements along the $x_0$-axis, |
507 |
\var{n1} elements along the $x_1$-axis and |
508 |
\var{n2} elements along the $x_2$-axis. |
509 |
For \var{order}=1 and \var{order}=2 |
510 |
\finleyelement{Hex8} and |
511 |
\finleyelement{Hex20} are used, respectively. |
512 |
In the case of \var{useElementsOnFace}=\False, |
513 |
\finleyelement{Rec4} and |
514 |
\finleyelement{Rec8} are used to subdivide the faces of the brick, respectively. |
515 |
In the case of \var{useElementsOnFace}=\True (this option should be used if gradients |
516 |
are calculated on domain faces), |
517 |
\finleyelement{Hex8Face} and |
518 |
\finleyelement{Hex20Face} are used on the brick faces, respectively. |
519 |
If \var{integrationOrder} is positive, a numerical integration scheme |
520 |
chosen which is accurate on each element up to a polynomial of |
521 |
degree \var{integrationOrder} \index{integration order}. Otherwise |
522 |
an appropriate integration order is chosen independently. If |
523 |
\var{periodic0}=\True, periodic boundary conditions \index{periodic boundary conditions} |
524 |
along the $x_0$-directions are enforced. That means when for any solution of a PDE solved by \finley |
525 |
the value on the plane $x_0=0$ will be identical to the values on $x_0=\var{l0}$. Correspondingly, |
526 |
\var{periodic1}=\False and \var{periodic2}=\False sets periodic boundary conditions |
527 |
in $x_1$-direction and $x_2$-direction, respectively. |
528 |
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. |
529 |
\end{funcdesc} |
530 |
|
531 |
\begin{funcdesc}{GlueFaces}{meshList,safetyFactor=0.2,tolerance=1.e-13} |
532 |
Generates a new \Domain object from the list \var{meshList} of \finley meshes. |
533 |
Nodes in face elements whose difference of coordinates is less then \var{tolerance} times the |
534 |
diameter of the domain are merged. The corresponding face elements are removed from the mesh. |
535 |
|
536 |
TODO: explain \var{safetyFactor} and show an example. |
537 |
\end{funcdesc} |
538 |
|
539 |
\begin{funcdesc}{JoinFaces}{meshList,safetyFactor=0.2,tolerance=1.e-13} |
540 |
Generates a new \Domain object from the list \var{meshList} of \finley meshes. |
541 |
Face elements whose nodes coordinates have difference is less then \var{tolerance} times the |
542 |
diameter of the domain are combined to form a contact element \index{element!contact} |
543 |
The corresponding face elements are removed from the mesh. |
544 |
|
545 |
TODO: explain \var{safetyFactor} and show an example. |
546 |
\end{funcdesc} |