forked from johannesgerer/jburkardt-f
-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathfem2d_stokes.html
613 lines (555 loc) · 18.4 KB
/
fem2d_stokes.html
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
<html>
<head>
<title>
FEM2D_STOKES - Finite Element Solution of the 2D Stokes Equations
</title>
</head>
<body bgcolor="#EEEEEE" link="#CC0000" alink="#FF3300" vlink="#000055">
<h1 align = "center">
FEM2D_STOKES <br>
Steady Incompressible Stokes Equations in 2D <br>
Finite Element Solution <br>
Banded Storage
</h1>
<hr>
<p>
<b>FEM2D_STOKES</b>
is a FORTRAN90 program which
applies the finite element method to solve
a form of the steady incompressible Stokes's equations over an arbitrary
triangulated region in 2D.
</p>
<p>
The geometry is entirely external to the program. The user
specifies one file of nodal coordinates, and one file that
describes the triangles in terms of six node coordinates.
</p>
<p>
The program makes a default assumption that all boundary velocities
correspond to Dirichlet boundary conditions, and that one pressure
is specified (for uniqueness of the pressure system). The user
can adjust these boundary conditions (and even specify Dirichlet
constraints on any variable at any node) by setting the appropriate
data in certain user routines.
</p>
<p>
<i>At the moment, Neumann conditions, if specified, must have a
zero right hand side. The machinery to integrate a nonzero
Neumann condition has not been set up yet.</i>
</p>
<p>
The linear system is created, stored, and solved using a
form of the LINPACK/LAPACK "general band" storage and versions
of LINPACK's DGBFA and DGBSL factorization and solution routines.
</p>
<h2 align = "center">
Computational Region
</h2>
<p>
The computational region is initially unknown by the program. The user
specifies it by preparing a file containing the coordinates of
the nodes, and a file containing the indices of nodes that make
up triangles that form a triangulation of the region. For the
following ridiculously small example:
<pre>
10-11-12
|\ |\
| \ | \
6 7 8 9
| \ | \
| \| \
1--2--3--4--5
</pre>
the node file could be:
<pre>
0.0 0.0
1.0 0.0
2.0 0.0
3.0 0.0
0.0 1.0
1.0 1.0
2.0 1.0
3.0 1.0
0.0 2.0
1.0 2.0
2.0 2.0
</pre>
and the triangle file would be
<pre>
1 3 10 2 7 6
12 10 3 11 7 8
3 5 12 4 9 8
</pre>
</p>
<h2 align = "center">
The Stokes Equations
</h2>
<p>
The state variables are a velocity vector (U,V)(X,Y) and a scalar
pressure P(X,Y). The state variables are constrained by the
momentum and continuity equations, which apply inside the region:
<pre><b>
- nu * ( Uxx + Uyy ) + dP/dx = U_RHS(x,y)
- nu * ( Vxx + Vyy ) + dP/dy = V_RHS(x,y)
dU/dx + dV/dy = P_RHS(x,y)
</b></pre>
where, typically, the right hand side functions are zero. However,
the user is free to assign nonzero values to these functions
through a user routine.
</p>
<h2 align = "center">
Boundary Conditions
</h2>
<p>
At every point on the boundary of the region, the program assumes
that both components of the velocity are specified.
<pre><b>
U(node) = U_BC(node)
V(node) = V_BC(node)
</b></pre>
This is known as a "Dirichlet boundary condition". The right hand side
of the boundary condition is left unspecified until the user
routine is called. If a wall is intended, then the appropriate
condition has U_BC and V_BC zero. An inlet might have a
particular flow profile function used to assign nonzero values.
</p>
<p>
At one point in the region, the program assumes
that the value of the pressure is specified.
<pre><b>
P(node) = P_BC(node)
</b></pre>
Such a condition must be supplied; otherwise the pressure cannot
be uniquely determined, since it is essentially a potential function,
and so is unique only "up to a constant". Note that the program allows
the user to specify pressure conditions anyway, and these can be
of Dirichlet or Neumann type. In general, however, this is not
a physically or mathematically or computationally good thing to do!
</p>
<p>
The user routine <b>boundary_type</b> can be used to modify the
type of the boundary conditions associated with a degree of freedom
at a boundary node - or even to add constraints to variables
associated with nodes in the interior.
</p>
<p>
One choice that the user can make is to reset certain boundary
conditions to be of Neumann type:
<pre><b>
dU/dn(node) = U_BC(node)
dV/dn(node) = V_BC(node)
</b></pre>
The right hand side of the boundary condition is left unspecified
until the user routine is called. As mentioned earlier, the program
cannot currently handle Neumann conditions with nonzero right hand side.
(A nonzero value is simply ignored, but won't actually cause the
program to fail.)
</p>
<h2 align = "center">
Computational Procedure
</h2>
<p>
We use linear finite elements for the pressure function, and to
generate these, we only need the nodes that are the vertices of
the triangles. (We will call these vertices "pressure nodes.")
Because quadratic basis functions are to be used
for the velocity, however, each triangle will also have three extra
midside nodes available for that.
</p>
<p>
We now assume that the unknown velocity component functions U(x,y)
and V(x,y) can be represented as linear combinations of the quadratic
basis functions associated with each node, and that the scalar
pressure P(x,y) can be represented as a linear combination of the linear
basis functions associated with each pressure node.
</p>
<p>
For every node, we can determine a quadratic velocity basis function
PSI(I)(x,y). For every pressure node I, we can determine a linear
basis function PHI(I)(x,y). If we assume that our solutions are linear
combinations of these basis functions, then we need to solve for
the coefficients.
</p>
<p>
To do so, we apply the Galerkin-Petrov method. For each pressure
node, and its corresponding basis function PHI(I), we generate a
copy of the continuity equation, multiplied by that basis function,
and integrated over the region:
<blockquote>
Integral ( Ux(x,y) + Vy(x,y) ) * PHI(I)(x,y) dx dy =
Integral ( P_RHS(x,y) * PHI(I)(x,y) dx dy )
</blockquote>
</p>
<p>
Similarly, for each node and its corresponding velocity basis function
PSI(I), we generate two copies of the momentum equation, one for
each component. We multiply by PSI(I), integrate over the region,
and use integration by parts to lower the order of differentiation.
This gives us:
<blockquote>
Integral nu * ( Ux(x,y) * PSIx(I)(x,y) + Uy(x,y) * PSIy(I)(x,y) )
+ Px(x,y) * PSI(I)(x,y) dx dy =
Integral ( U_RHS(x,y) * PSI(I)(x,y) dx dy )
<br>
Integral nu * ( Vx(x,y) * PSIx(I)(x,y) + Vy(x,y) * PSIy(I)(x,y) )
+ Py(x,y) * PSI(I)(x,y) dx dy =
Integral ( V_RHS(x,y) * PSI(I)(x,y) dx dy )
</blockquote>
</p>
<p>
After adjusting for the boundary conditions, the set of all such
equations yields a linear system for the coefficients of the
finite element representation of the solution.
</p>
<h2 align = "center">
User Input Routines
</h2>
<p>
The program requires the user to supply the following routines:
</p>
<p>
The default boundary condition types are passed to the user, along
with other information. The user modifies any data as necessary, and
returns it. This is done by a user-supplied routine of the form
<blockquote><b>
subroutine boundary_type ( node_num, node_xy, node_boundary, node_type,
node_u_condition, node_v_condition, node_p_condition )
</b></blockquote>
</p>
<p>
The right hand sides of any Dirichlet boundary conditions are
determined by a user-supplied routine of the form
<blockquote><b>
subroutine dirichlet_condition ( node_num, node_xy, u_bc, v_bc, p_bc )
</b></blockquote>
</p>
<p>
The right hand sides of any Neumann boundary conditions are
determined by a user-supplied routine of the form
<blockquote><b>
subroutine neumann_condition ( node_num, node_xy, u_bc, v_bc, p_bc )
</b></blockquote>
</p>
<p>
The right hand sides (or "source terms") of the state equations
are determined by a user-supplied routine of the form
<blockquote><b>
subroutine rhs ( node_num, node_xy, u_rhs, v_rhs, p_rhs )
</b></blockquote>
</p>
<h2 align = "center">
Program Output
</h2>
<p>
The program writes out various node, triangle, pressure and velocity
data files that can be used to create plots of the geometry and
the solution.
</p>
<p>
Graphics files created include:
<ul>
<li>
<b>nodes6.eps</b>, an image of the nodes;
</li>
<li>
<b>triangles6.eps</b>, an image of the elements;
</li>
</ul>
</p>
<p>
Data files created include:
<ul>
<li>
<b>nodes3.txt</b>, the nodes associated with pressure;
</li>
<li>
<b>triangles3.txt</b>, the linear triangles associated with pressure;
</li>
<li>
<b>pressure3.txt</b>, the pressure at the pressure nodes;
</li>
<li>
<b>velocity6.txt</b>, the velocity at the velocity nodes.
</li>
</ul>
</p>
<h3 align = "center">
Usage:
</h3>
<p>
To run the program, the user must write a file, called, perhaps,
<i>myprog.f90</i>, containing routines defining certain data,
compile this file with <b>FEM2D_STOKES.F90</b>,
and run the executable, supplying the node and triangle files
on the command line.
</p>
<p>
<dl>
<dt>
f90 fem2d_stokes.f90 <i>myprog.f90</i>
</dt>
<dd>
compiles <b>fem2d_stokes.f90</b> and your program,
and creates an executable called <i>a.out</i>.
</dd>
<dt>
a.out <i>nodes.txt</i> <i>triangles.txt</i>
</dt>
<dd>
runs the program with the geometry defined in
<i>nodes.txt</i> and <i>triangles.txt</i>.
</dd>
</dl>
</p>
<h3 align = "center">
Licensing:
</h3>
<p>
The computer code and data files described and made available on this web page
are distributed under
<a href = "../../txt/gnu_lgpl.txt">the GNU LGPL license.</a>
</p>
<h3 align = "center">
Languages:
</h3>
<p>
<b>FEM2D_STOKES</b> is available in
<a href = "../../cpp_src/fem2d_stokes/fem2d_stokes.html">a C++ version</a> and
<a href = "../../f_src/fem2d_stokes/fem2d_stokes.html">a FORTRAN90 version</a> and
<a href = "../../m_src/fem2d_stokes/fem2d_stokes.html">a MATLAB version</a>.
</p>
<h3 align = "center">
Related Data and Programs:
</h3>
<p>
<a href = "../../f_src/fem2d_stokes_cavity/fem2d_stokes_cavity.html">
FEM2D_STOKES_CAVITY</a>,
a FORTRAN90 library which
contains the user-supplied routines necessary to run <b>fem2d_stokes</b>
on the "cavity" problem.
</p>
<p>
<a href = "../../f_src/fem2d_stokes_channel/fem2d_stokes_channel.html">
FEM2D_STOKES_CHANNEL</a>,
a FORTRAN90 library which
contains the user-supplied routines necessary to run <b>fem2d_stokes</b>
on the "channel" problem.
</p>
<p>
<a href = "../../f_src/fem2d_stokes_inout/fem2d_stokes_inout.html">
FEM2D_STOKES_INOUT</a>,
a FORTRAN90 library which
contains the user-supplied routines necessary to run <b>fem2d_stokes</b>
on the "inout" problem.
</p>
<h3 align = "center">
Reference:
</h3>
<p>
<ol>
<li>
Max Gunzburger,<br>
Finite Element Methods for Viscous Incompressible Flows,<br>
A Guide to Theory, Practice, and Algorithms,<br>
Academic Press, 1989,<br>
ISBN: 0-12-307350-2,<br>
LC: TA357.G86.
</li>
<li>
Hans Rudolf Schwarz,<br>
Finite Element Methods,<br>
Academic Press, 1988,<br>
ISBN: 0126330107,<br>
LC: TA347.F5.S3313.
</li>
<li>
Gilbert Strang, George Fix,<br>
An Analysis of the Finite Element Method,<br>
Cambridge, 1973,<br>
ISBN: 096140888X,<br>
LC: TA335.S77.
</li>
<li>
Olgierd Zienkiewicz,<br>
The Finite Element Method,<br>
Sixth Edition,<br>
Butterworth-Heinemann, 2005,<br>
ISBN: 0750663200,<br>
LC: TA640.2.Z54
</li>
</ol>
</p>
<h3 align = "center">
Source code:
</h3>
<p>
<ul>
<li>
<a href = "fem2d_stokes.f90">fem2d_stokes.f90</a>,
the source code;
</li>
<li>
<a href = "fem2d_stokes.sh">fem2d_stokes.sh</a>,
commands to compile the partial program;
</li>
</ul>
</p>
<h3 align = "center">
List of Routines:
</h3>
<p>
<ul>
<li>
<b>MAIN</b> is the main program for FEM2D_STOKES.
</li>
<li>
<b>ASSEMBLE_STOKES</b> assembles the finite element Stokes equations.
</li>
<li>
<b>BANDWIDTH</b> determines the bandwidth of the coefficient matrix.
</li>
<li>
<b>BASIS_MN_T3:</b> all bases at N points for a T3 element.
</li>
<li>
<b>BASIS_MN_T6:</b> all bases at N points for a T6 element.
</li>
<li>
<b>CH_CAP</b> capitalizes a single character.
</li>
<li>
<b>CH_EQI</b> is a case insensitive comparison of two characters for equality.
</li>
<li>
<b>CH_TO_DIGIT</b> returns the integer value of a base 10 digit.
</li>
<li>
<b>DGB_FA</b> performs a LINPACK-style PLU factorization of an DGB matrix.
</li>
<li>
<b>DGB_MXV</b> multiplies a DGB matrix times a vector.
</li>
<li>
<b>DGB_PRINT_SOME</b> prints some of a DGB matrix.
</li>
<li>
<b>DGB_SL</b> solves a system factored by DGB_FA.
</li>
<li>
<b>DIRICHLET_APPLY</b> accounts for Dirichlet boundary conditions.
</li>
<li>
<b>FILE_COLUMN_COUNT</b> counts the number of columns in the first line of a file.
</li>
<li>
<b>FILE_NAME_SPECIFICATION</b> determines the names of the input files.
</li>
<li>
<b>FILE_ROW_COUNT</b> counts the number of row records in a file.
</li>
<li>
<b>GET_UNIT</b> returns a free FORTRAN unit number.
</li>
<li>
<b>I4_HUGE</b> returns a "huge" I4.
</li>
<li>
<b>I4COL_COMPARE</b> compares columns I and J of a integer array.
</li>
<li>
<b>I4COL_SORT_A</b> ascending sorts an I4COL.
</li>
<li>
<b>I4COL_SWAP</b> swaps columns I and J of an I4COL.
</li>
<li>
<b>I4MAT_DATA_READ</b> reads data from an I4MAT file.
</li>
<li>
<b>I4MAT_HEADER_READ</b> reads the header from an I4MAT.
</li>
<li>
<b>I4MAT_TRANSPOSE_PRINT_SOME</b> prints some of the transpose of an I4MAT.
</li>
<li>
<b>LVEC_PRINT</b> prints a logical vector.
</li>
<li>
<b>NEUMANN_APPLY</b> accounts for Neumann boundary conditions.
</li>
<li>
<b>NODES3_WRITE</b> writes the pressure nodes to a file.
</li>
<li>
<b>POINTS_PLOT</b> plots a pointset.
</li>
<li>
<b>PRESSURE3_WRITE</b> writes the pressures to a file.
</li>
<li>
<b>QUAD_RULE</b> sets the quadrature rule for assembly.
</li>
<li>
<b>R8MAT_DATA_READ</b> reads data from an R8MAT file.
</li>
<li>
<b>R8MAT_HEADER_READ</b> reads the header from an R8MAT file.
</li>
<li>
<b>R8MAT_TRANSPOSE_PRINT_SOME</b> prints some of an R8MAT, transposed.
</li>
<li>
<b>R8VEC_PRINT_SOME</b> prints "some" of an R8VEC.
</li>
<li>
<b>REFERENCE_TO_PHYSICAL_T6</b> maps T6 reference points to physical points.
</li>
<li>
<b>S_TO_I4</b> reads an I4 from a string.
</li>
<li>
<b>S_TO_I4VEC</b> reads an I4VEC from a string.
</li>
<li>
<b>S_TO_R8</b> reads an R8 from a string.
</li>
<li>
<b>S_TO_R8VEC</b> reads an R8VEC from a string.
</li>
<li>
<b>S_WORD_COUNT</b> counts the number of "words" in a string.
</li>
<li>
<b>SORT_HEAP_EXTERNAL</b> externally sorts a list of items into ascending order.
</li>
<li>
<b>TIMESTAMP</b> prints the current YMDHMS date as a time stamp.
</li>
<li>
<b>TRIANGLE_AREA_2D</b> computes the area of a triangle in 2D.
</li>
<li>
<b>TRIANGLES3_WRITE</b> writes the pressure triangles to a file.
</li>
<li>
<b>TRIANGULATION_ORDER6_BOUND_NODE</b> indicates which nodes are on the boundary.
</li>
<li>
<b>TRIANGULATION_ORDER6_PLOT</b> plots a 6-node triangulation of a set of nodes.
</li>
<li>
<b>VELOCITY6_WRITE</b> writes the velocities to a file.
</li>
</ul>
</p>
<p>
You can go up one level to <a href = "../f_src.html">
the FORTRAN90 source codes</a>.
</p>
<hr>
<i>
Last revised on 04 January 2011.
</i>
<!-- John Burkardt -->
</body>
</html>