Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
flow_profile.f90
Go to the documentation of this file.
1! Copyright (c) 2021, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
35 use num_types, only : rp
36 implicit none
37 private
38
40 abstract interface
41 function blasius_profile(y, delta, u)
42 import rp
43 real(kind=rp), intent(in) :: y, delta, u
44 real(kind=rp) :: blasius_profile
45 end function blasius_profile
46 end interface
47
50
51contains
52
55 function blasius_linear(y, delta, u)
56 real(kind=rp), intent(in) :: y, delta, u
57 real(kind=rp) :: blasius_linear
58 real(kind=rp) :: arg
59
60 arg = y /delta
61 if (arg .gt. 1.0_rp) then
63 else
64 blasius_linear = u * (y / delta)
65 end if
66
67 end function blasius_linear
68
71 function blasius_quadratic(y, delta, u)
72 real(kind=rp), intent(in) :: y, delta, u
73 real(kind=rp) :: blasius_quadratic
74 real(kind=rp) :: arg
75
76 arg = ( 2.0_rp * (y / delta) - (y / delta)**2 )
77
78 if (arg .gt. 1.0_rp) then
80 else
81 blasius_quadratic = u * arg
82 end if
83
84 end function blasius_quadratic
85
88 function blasius_cubic(y, delta, u)
89 real(kind=rp), intent(in) :: y, delta, u
90 real(kind=rp) :: blasius_cubic
91 real(kind=rp) :: arg
92
93 arg = ( 3.0_rp / 2.0_rp * (y / delta) - 0.5_rp * (y / delta)**3 )
94
95 if (arg .gt. 1.0_rp) then
96 blasius_cubic = 1.0_rp
97 else
98 blasius_cubic = u * arg
99 end if
100
101 end function blasius_cubic
102
106 function blasius_quartic(y, delta, u)
107 real(kind=rp), intent(in) :: y, delta, u
108 real(kind=rp) :: blasius_quartic
109 real(kind=rp) :: arg
110
111 arg = 2.0_rp * (y / delta) - 2.0_rp * (y / delta)**3 + (y / delta)**4
112
113 if (arg .gt. 1.0_rp) then
115 else
116 blasius_quartic = u * arg
117 end if
118
119 end function blasius_quartic
120
123 function blasius_sin(y, delta, u)
124 real(kind=rp), intent(in) :: y, delta, u
125 real(kind=rp) :: blasius_sin
126 real(kind=rp), parameter :: pi = 4.0_rp * atan(1.0_rp)
127 real(kind=rp) :: arg
128
129 arg = (pi / 2.0_rp) * (y/delta)
130
131 if (arg .gt. 0.5_rp * pi) then
132 blasius_sin = u
133 else
134 blasius_sin = u * sin(arg)
135 end if
136
137 end function blasius_sin
138
139end module flow_profile
Abstract interface for computing a Blasius flow profile.
Defines a flow profile.
real(kind=rp) function, public blasius_quadratic(y, delta, u)
Quadratic approximate Blasius Profile .
real(kind=rp) function, public blasius_quartic(y, delta, u)
Quartic approximate Blasius Profile .
real(kind=rp) function, public blasius_sin(y, delta, u)
Sinusoidal approximate Blasius Profile .
real(kind=rp) function, public blasius_cubic(y, delta, u)
Cubic approximate Blasius Profile .
real(kind=rp) function, public blasius_linear(y, delta, u)
Linear approximate Blasius profile .
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12