Neko  0.8.1
A portable framework for high-order spectral element flow simulations
time_interpolator.f90
Go to the documentation of this file.
1 ! Copyright (c) 2022, 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
38  use device_math, only : device_add3s2
39  use math, only : add3s2
40  use utils, only : neko_error
41  use, intrinsic :: iso_c_binding
42  implicit none
43  private
44 
46  type, public :: time_interpolator_t
48  integer :: order
49  contains
51  procedure, pass(this) :: init => time_interpolator_init
53  procedure, pass(this) :: free => time_interpolator_free
55  procedure, pass(this) :: interpolate => time_interpolator_interpolate
56 
57  end type time_interpolator_t
58 
59 contains
60 
63  subroutine time_interpolator_init(this, order)
64  class(time_interpolator_t), intent(inout) :: this
65  integer, intent(in), target :: order
66 
68  call this%free()
69 
71  this%order = order
72 
73  end subroutine time_interpolator_init
74 
75 
77  subroutine time_interpolator_free(this)
78  class(time_interpolator_t), intent(inout) :: this
79 
80  end subroutine time_interpolator_free
81 
89  subroutine time_interpolator_interpolate(this, t, f, t_past, f_past, t_future, f_future)
90  class(time_interpolator_t), intent(inout) :: this
91  real(kind=rp), intent(inout) :: t, t_past, t_future
92  type(field_t), intent(inout) :: f, f_past, f_future
93  real(kind=rp) :: w_past, w_future !Weights for the interpolation
94  integer :: n
95 
96  if (this%order .eq. 2) then
97 
98  n = f%dof%size()
99  w_past = ( t_future - t ) / ( t_future - t_past )
100  w_future = ( t - t_past ) / ( t_future - t_past )
101 
102  if (neko_bcknd_device .eq. 1) then
103  call device_add3s2(f%x_d, f_past%x_d, f_future%x_d, &
104  w_past, w_future, n)
105  else
106  call add3s2(f%x, f_past%x, f_future%x, w_past, w_future, n)
107  end if
108 
109  else
110  call neko_error("Time interpolation of required order is not implemented")
111  end if
112 
113  end subroutine time_interpolator_interpolate
114 
115 end module time_interpolator
subroutine, public device_add3s2(a_d, b_d, c_d, c1, c2, n)
Defines a field.
Definition: field.f90:34
Definition: math.f90:60
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
Definition: math.f90:686
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.
Utilities.
Definition: utils.f90:35
Provides a tool to perform interpolation in time.