Neko  0.8.1
A portable framework for high-order spectral element flow simulations
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
36  implicit none
37 
39  abstract interface
40  function blasius_profile(y, delta, u)
41  import rp
42  real(kind=rp), intent(in) :: y, delta, u
43  real(kind=rp) :: blasius_profile
44  end function blasius_profile
45  end interface
46 
47 contains
48 
51  function blasius_linear(y, delta, u)
52  real(kind=rp), intent(in) :: y, delta, u
53  real(kind=rp) :: blasius_linear
54  real(kind=rp) :: arg
55 
56  arg = y /delta
57  if (arg .gt. 1.0_rp) then
58  blasius_linear = u
59  else
60  blasius_linear = u * (y / delta)
61  end if
62 
63  end function blasius_linear
64 
67  function blasius_quadratic(y, delta, u)
68  real(kind=rp), intent(in) :: y, delta, u
69  real(kind=rp) :: blasius_quadratic
70  real(kind=rp) :: arg
71 
72  arg = ( 2.0_rp * (y / delta) - (y / delta)**2 )
73 
74  if (arg .gt. 1.0_rp) then
76  else
77  blasius_quadratic = u * arg
78  end if
79 
80  end function blasius_quadratic
81 
84  function blasius_cubic(y, delta, u)
85  real(kind=rp), intent(in) :: y, delta, u
86  real(kind=rp) :: blasius_cubic
87  real(kind=rp) :: arg
88 
89  arg = ( 3.0_rp / 2.0_rp * (y / delta) - 0.5_rp * (y / delta)**3 )
90 
91  if (arg .gt. 1.0_rp) then
92  blasius_cubic = 1.0_rp
93  else
94  blasius_cubic = u * arg
95  end if
96 
97  end function blasius_cubic
98 
102  function blasius_quartic(y, delta, u)
103  real(kind=rp), intent(in) :: y, delta, u
104  real(kind=rp) :: blasius_quartic
105  real(kind=rp) :: arg
106 
107  arg = 2.0_rp * (y / delta) - 2.0_rp * (y / delta)**3 + (y / delta)**4
108 
109  if (arg .gt. 1.0_rp) then
110  blasius_quartic = u
111  else
112  blasius_quartic = u * arg
113  end if
114 
115  end function blasius_quartic
116 
119  function blasius_sin(y, delta, u)
120  real(kind=rp), intent(in) :: y, delta, u
121  real(kind=rp) :: blasius_sin
122  real(kind=rp), parameter :: pi = 4.0_rp * atan(1.0_rp)
123  real(kind=rp) :: arg
124 
125  arg = (pi / 2.0_rp) * (y/delta)
126 
127  if (arg .gt. 0.5_rp * pi) then
128  blasius_sin = u
129  else
130  blasius_sin = u * sin(arg)
131  end if
132 
133  end function blasius_sin
134 
135 end module flow_profile
Abstract interface for computing a Blasius flow profile.
Defines a flow profile.
real(kind=rp) function blasius_quadratic(y, delta, u)
Quadratic approximate Blasius Profile .
real(kind=rp) function blasius_sin(y, delta, u)
Sinusoidal approximate Blasius Profile .
real(kind=rp) function blasius_cubic(y, delta, u)
Cubic approximate Blasius Profile .
real(kind=rp) function blasius_linear(y, delta, u)
Linear approximate Blasius profile .
real(kind=rp) function blasius_quartic(y, delta, u)
Quartic approximate Blasius Profile .
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12