Neko  0.8.99
A portable framework for high-order spectral element flow simulations
time_interpolator.f90
Go to the documentation of this file.
1 ! Copyright (c) 2022-2024, 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  use field, only : field_t
37  use neko_config, only : neko_bcknd_device
38  use device_math, only : device_add3s2
39  use math, only : add3s2, rzero
40  use utils, only : neko_error
41  use, intrinsic :: iso_c_binding
42  use fast3d, only: fd_weights_full
43  implicit none
44  private
45 
47  type, public :: time_interpolator_t
49  integer :: order
50  contains
52  procedure, pass(this) :: init => time_interpolator_init
54  procedure, pass(this) :: free => time_interpolator_free
56  procedure, pass(this) :: interpolate => time_interpolator_interpolate
58  procedure, pass(this) :: interpolate_scalar => time_interpolator_scalar
59 
60  end type time_interpolator_t
61 
62 contains
63 
66  subroutine time_interpolator_init(this, order)
67  class(time_interpolator_t), intent(inout) :: this
68  integer, intent(in), target :: order
69 
71  call this%free()
72 
74  this%order = order
75 
76  end subroutine time_interpolator_init
77 
78 
80  subroutine time_interpolator_free(this)
81  class(time_interpolator_t), intent(inout) :: this
82 
83  end subroutine time_interpolator_free
84 
92  subroutine time_interpolator_interpolate(this, t, f, t_past, f_past, &
93  t_future, f_future)
94  class(time_interpolator_t), intent(inout) :: this
95  real(kind=rp), intent(inout) :: t, t_past, t_future
96  type(field_t), intent(inout) :: f, f_past, f_future
97  real(kind=rp) :: w_past, w_future !Weights for the interpolation
98  integer :: n
99 
100  if (this%order .eq. 2) then
101 
102  n = f%dof%size()
103  w_past = ( t_future - t ) / ( t_future - t_past )
104  w_future = ( t - t_past ) / ( t_future - t_past )
105 
106  if (neko_bcknd_device .eq. 1) then
107  call device_add3s2(f%x_d, f_past%x_d, f_future%x_d, &
108  w_past, w_future, n)
109  else
110  call add3s2(f%x, f_past%x, f_future%x, w_past, w_future, n)
111  end if
112 
113  else
114  call neko_error("Time interpolation of required order &
115  &is not implemented")
116  end if
117 
118  end subroutine time_interpolator_interpolate
119 
127  subroutine time_interpolator_scalar(this, t, f_interpolated, f_n, tlag, n)
128  class(time_interpolator_t), intent(inout) :: this
129  real(kind=rp), intent(in) :: t
130  integer, intent(in) :: n
131  real(kind=rp), dimension(n, 0:this%order - 1), intent(in) :: f_n
132  real(kind=rp), dimension(n), intent(inout) :: f_interpolated
133  real(kind=rp), dimension(0:this%order), intent(in) :: tlag
134  integer :: no, i, l
135 
136 
137  integer, parameter :: lwtmax = 10
138  real(kind=rp) :: wt(0:lwtmax)
139  wt = 0
140 
141  if (this%order .gt. lwtmax) then
142  call neko_error("lwtmax is smaller than the number &
143  &of stored convecting fields")
144  end if
145 
146  no = this%order - 1
147  call fd_weights_full(t, tlag, no, 0, wt) ! interpolation weights
148  call rzero(f_interpolated, n)
149 
150  do i = 1, n
151  do l = 0, no
152  f_interpolated(i) = f_interpolated(i) + wt(l) * f_n(i,l)
153  end do
154  end do
155 
156  end subroutine time_interpolator_scalar
157 
158 end module time_interpolator
subroutine, public device_add3s2(a_d, b_d, c_d, c1, c2, n)
Returns .
Fast diagonalization methods from NEKTON.
Definition: fast3d.f90:61
subroutine, public fd_weights_full(xi, x, n, m, c)
Compute finite-difference stencil weights for evaluating derivatives up to order at a point.
Definition: fast3d.f90:105
Defines a field.
Definition: field.f90:34
Definition: math.f90:60
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
Definition: math.f90:730
subroutine, public rzero(a, n)
Zero a real vector.
Definition: math.f90:184
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Implements type time_interpolator_t.
subroutine time_interpolator_interpolate(this, t, f, t_past, f_past, t_future, f_future)
Interpolate a field at time t from fields at time t-dt and t+dt.
subroutine time_interpolator_free(this)
Destructor.
subroutine time_interpolator_init(this, order)
Constructor.
subroutine time_interpolator_scalar(this, t, f_interpolated, f_n, tlag, n)
Interpolate a scalar field at time t from known scalar fields at different time steps.
Utilities.
Definition: utils.f90:35
Provides a tool to perform interpolation in time.