Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_hermite_interpolator_1d.F90
Go to the documentation of this file.
1 
7 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
8 #include "sll_assert.h"
9 #include "sll_errors.h"
10 #include "sll_memory.h"
11 #include "sll_working_precision.h"
12 
18 
19  use sll_m_interpolators_1d_base, only: &
21 
22  implicit none
23 
24  public :: &
26 
27  private
28 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
29 
37 
41  type(sll_t_hermite_interpolation_1d), pointer :: hermite
43  sll_int32 :: npts
44  contains
46  procedure, pass(interpolator) :: init => initialize_hermite_interpolator_1d
48  procedure :: compute_interpolants => wrap_compute_interpolants_hermite_1d
50  procedure :: interpolate_from_interpolant_derivatives_eta1 => interpolate_array_derivatives_hi1d
52  procedure :: interpolate_array => wrap_interpolate_array_hermite_1d
53 
55  procedure :: interpolate_array_disp => interpolate_array_disp_hi1d
57  procedure :: interpolate_array_disp_inplace => interpolate_array_disp_inplace_hi1d
59  procedure :: interpolate_from_interpolant_derivative_eta1 => interpolate_derivative_eta1_hi1d
61  procedure :: interpolate_from_interpolant_array => interpolate_array_values_hi1d
63  procedure :: interpolate_from_interpolant_value => wrap_interpolate_value_hermite_1d
65  procedure, pass :: set_coefficients => set_coefficients_hi1d
67  procedure, pass :: get_coefficients => get_coefficients_hi1d
69 
70 !PN DEFINED BUT NOT USED
71 ! !> Deallocate the class interpolator
72 ! interface sll_o_delete
73 ! module procedure delete_hi1d
74 ! end interface
75 
76 contains !**********************************************************
77 
80  npts, &
81  eta_min, &
82  eta_max, &
83  degree, &
84  eta_hermite_continuity, &
85  eta_bc_type, &
86  const_eta_min_slope, &
87  const_eta_max_slope, &
88  eta_min_slopes, &
89  eta_max_slopes) &
90  result(interpolator)
91 
92  type(sll_hermite_interpolator_1d), pointer :: interpolator
93  sll_int32, intent(in) :: npts
94  sll_real64, intent(in) :: eta_min
95  sll_real64, intent(in) :: eta_max
96  sll_int32, intent(in) :: degree
97  sll_int32, intent(in) :: eta_hermite_continuity
98  sll_int32, intent(in) :: eta_bc_type
99  sll_real64, intent(in), optional :: const_eta_min_slope
100  sll_real64, intent(in), optional :: const_eta_max_slope
101  sll_real64, dimension(:), intent(in), optional :: eta_min_slopes
102  sll_real64, dimension(:), intent(in), optional :: eta_max_slopes
103  sll_int32 :: ierr
104 
105  sll_allocate(interpolator, ierr)
106 
107  interpolator%npts = npts
108 
109  call interpolator%init( &
110  npts, &
111  eta_min, &
112  eta_max, &
113  degree, &
114  eta_hermite_continuity, &
115  eta_bc_type, &
116  const_eta_min_slope, &
117  const_eta_max_slope, &
118  eta_min_slopes, &
119  eta_max_slopes)
120 
122 
124  interpolator, &
125  npts, &
126  eta_min, &
127  eta_max, &
128  degree, &
129  eta_hermite_continuity, &
130  eta_bc_type, &
131  const_eta_min_slope, &
132  const_eta_max_slope, &
133  eta_min_slopes, &
134  eta_max_slopes)
135 
136  class(sll_hermite_interpolator_1d), intent(inout) :: interpolator
137  sll_int32, intent(in) :: npts
138  sll_real64, intent(in) :: eta_min
139  sll_real64, intent(in) :: eta_max
140  sll_int32, intent(in) :: degree
141  sll_int32, intent(in) :: eta_hermite_continuity
142  sll_int32, intent(in) :: eta_bc_type
143  sll_real64, intent(in), optional :: const_eta_min_slope
144  sll_real64, intent(in), optional :: const_eta_max_slope
145  sll_real64, dimension(:), intent(in), optional :: eta_min_slopes
146  sll_real64, dimension(:), intent(in), optional :: eta_max_slopes
147 
148  interpolator%hermite => sll_f_new_hermite_interpolation_1d( &
149  npts, &
150  eta_min, &
151  eta_max, &
152  degree, &
153  eta_hermite_continuity, &
154  eta_bc_type, &
155  const_eta_min_slope, &
156  const_eta_max_slope, &
157  eta_min_slopes, &
158  eta_max_slopes)
159 
161 
163  interpolator, &
164  data_array, &
165  eta_coords, &
166  size_eta_coords)
167  class(sll_hermite_interpolator_1d), intent(inout) :: interpolator
168  sll_real64, dimension(:), intent(in) :: data_array
169  sll_real64, dimension(:), intent(in), optional :: eta_coords
170  sll_int32, intent(in), optional :: size_eta_coords
171 
172  if (present(eta_coords) .or. present(size_eta_coords)) then
173  sll_error('wrap_compute_interpolants_hermite_1d', 'This case is not yet implemented')
174  end if
175 
176  call sll_s_compute_interpolants_hermite_1d(interpolator%hermite, data_array)
177 
179 
180  function wrap_interpolate_value_hermite_1d(interpolator, eta1) result(val)
181  class(sll_hermite_interpolator_1d), intent(in) :: interpolator
182  sll_real64 :: val
183  sll_real64, intent(in) :: eta1
184  val = sll_f_interpolate_value_hermite_1d(eta1, interpolator%hermite)
185 
187 
189  this, &
190  num_pts, &
191  data, &
192  coordinates, &
193  output_array)
194  class(sll_hermite_interpolator_1d), intent(inout) :: this
195  sll_int32, intent(in) :: num_pts
196  sll_real64, dimension(num_pts), intent(in) :: coordinates
197  sll_real64, dimension(:), intent(in) :: data
198  sll_real64, dimension(num_pts), intent(out) :: output_array
199  sll_int32 :: i
200  call sll_s_compute_interpolants_hermite_1d(this%hermite, data)
201  do i = 1, num_pts
202  output_array(i) = this%interpolate_from_interpolant_value(coordinates(i))
203  end do
204  end subroutine wrap_interpolate_array_hermite_1d
205 
206  subroutine interpolate_array_disp_hi1d(this, num_pts, data, alpha, output_array)
207  class(sll_hermite_interpolator_1d), intent(inout) :: this
208  sll_real64, intent(in) :: alpha
209  sll_int32, intent(in) :: num_pts ! size of output array
210  sll_real64, dimension(:), intent(in) :: data ! data to be interpolated points where output is desired
211  sll_real64, dimension(1:num_pts), intent(out) :: output_array
212 
213 !call interpolate_from_interpolant_array(data,alpha,this%hermite)
214 !data_out=this%hermite%data_out
215  print *, 'interpolate_array_disp_hi1d:', &
216  ' not implemented for hermite interpolation'
217  sll_assert(this%npts > 0)
218  output_array = 0.0_f64*alpha + data
219 
220  end subroutine interpolate_array_disp_hi1d
221 
222  subroutine interpolate_array_disp_inplace_hi1d(this, num_pts, data, alpha)
223  class(sll_hermite_interpolator_1d), intent(inout) :: this
224  sll_real64, intent(in) :: alpha
225  sll_int32, intent(in) :: num_pts ! size of output array
226  sll_real64, dimension(num_pts), intent(inout) :: data ! data to be interpolated points where output is desired
227 
228 !call interpolate_from_interpolant_array(data,alpha,this%hermite)
229 !data_out=this%hermite%data_out
230  print *, 'interpolate_array_disp_inplace_hi1d:', &
231  ' not implemented for hermite interpolation'
232  sll_assert(this%npts > 0)
233 
235 
236 !PN DEFINED BUT NOT USED
237 !subroutine delete_hi1d (obj)
238 ! class(sll_hermite_interpolator_1d) :: obj
239 ! !should be fixed; for the moment just commented
240 ! !call delete(obj%hermite)
241 ! SLL_ASSERT(obj%npts>0)
242 !end subroutine delete_hi1d
243 
245  interpolator, &
246  num_pts, &
247  vals_to_interpolate, &
248  output_array)
249  class(sll_hermite_interpolator_1d), intent(inout) :: interpolator
250  sll_int32, intent(in) :: num_pts
251  sll_real64, dimension(num_pts), intent(in) :: vals_to_interpolate
252  sll_real64, dimension(num_pts), intent(out) :: output_array
253  !sll_int32 :: ierr
254  output_array = 0.0_f64
255  print *, 'interpolate_from_interpolant_array:', &
256  ' not implemented for hermite interpolation'
257  print *, num_pts
258  print *, maxval(vals_to_interpolate)
259  output_array = 0._f64
260  !print *,interpolator%bc_type
261  stop
262  sll_assert(interpolator%npts > 0)
263  end subroutine interpolate_array_values_hi1d
264 
266  interpolator, &
267  num_pts, &
268  vals_to_interpolate, &
269  output_array)
270  class(sll_hermite_interpolator_1d), intent(inout) :: interpolator
271  sll_int32, intent(in) :: num_pts
272  sll_real64, dimension(:), intent(in) :: vals_to_interpolate
273  sll_real64, dimension(:), intent(out) :: output_array
274  !sll_int32 :: ierr
275  output_array = 0.0_f64
276  print *, 'interpolate_from_interpolant_derivatives_eta1: ', &
277  'not implemented for hermite interpolation'
278  print *, num_pts
279  print *, maxval(vals_to_interpolate)
280  print *, maxval(output_array)
281  !print *,interpolator%bc_type
282  stop
283  sll_assert(interpolator%npts > 0)
285 
286  function interpolate_derivative_eta1_hi1d(interpolator, eta1) result(val)
287  class(sll_hermite_interpolator_1d), intent(in) :: interpolator
288  sll_real64 :: val
289  sll_real64, intent(in) :: eta1
290  print *, 'interpolate_derivative_eta1_hi1d: ', &
291  'not implemented for hermite interpolation'
292  !print *,interpolator%bc_type
293  print *, eta1
294  val = 0._f64
295  stop
296  sll_assert(interpolator%npts > 0)
297  end function
298 
299 !PN DEFINED BUT NOT USED
300 ! function interpolate_value_hi1d( interpolator, eta1 ) result(val)
301 ! class(sll_hermite_interpolator_1d), intent(in) :: interpolator
302 ! sll_real64 :: val
303 ! sll_real64, intent(in) :: eta1
304 ! print*, 'interpolate_value_hi1d: ', &
305 ! 'not implemented for hermite interpolation'
306 ! val = 0._f64
307 ! print *,eta1
308 ! !print *,interpolator%bc_type
309 ! stop
310 ! end function
311 
312 !PN DEFINED BUT NOT USED
313 ! function interpolate_array_hi1d(this, num_points, data, coordinates) &
314 ! result(data_out)
315 ! class(sll_hermite_interpolator_1d), intent(in) :: this
316 ! !class(sll_spline_1D), intent(in) :: this
317 ! sll_int32, intent(in) :: num_points
318 ! sll_real64, dimension(:), intent(in) :: coordinates
319 ! sll_real64, dimension(:), intent(in) :: data
320 ! sll_real64, dimension(num_points) :: data_out
321 ! ! local variables
322 ! !sll_int32 :: ierr
323 ! print*, 'interpolate_array_hi1d: ', &
324 ! 'not implemented for hermite interpolation'
325 ! print *,maxval(coordinates)
326 ! print *,maxval(data)
327 ! !print *,this%bc_type
328 ! data_out = 0._f64
329 ! stop
330 ! end function
331 
332 !PN DEFINED BUT NOT USED
333 ! subroutine compute_interpolants_hi1d( interpolator, data_array,&
334 ! eta_coords, &
335 ! size_eta_coords)
336 ! class(sll_hermite_interpolator_1d), intent(inout) :: interpolator
337 ! sll_real64, dimension(:), intent(in) :: data_array
338 ! sll_real64, dimension(:), intent(in),optional :: eta_coords
339 ! sll_int32, intent(in),optional :: size_eta_coords
340 ! print*, 'compute_interpolants_hi1d:', &
341 ! ' not implemented for hermite interpolation'
342 ! if(present(eta_coords))then
343 ! print *,'eta_coords present but not used'
344 ! endif
345 ! if(present(size_eta_coords))then
346 ! print *,'size_eta_coords present but not used'
347 ! endif
348 ! print *,maxval(data_array)
349 ! !print *,interpolator%bc_type
350 ! stop
351 ! end subroutine
352 
353  subroutine set_coefficients_hi1d(interpolator, coeffs)
354  class(sll_hermite_interpolator_1d), intent(inout) :: interpolator
355  sll_real64, dimension(:), intent(in), optional :: coeffs
356  print *, 'set_coefficients_hi1d(): ERROR: This function has not been ', &
357  'implemented yet.'
358  if (present(coeffs)) then
359  print *, '#coeffs present but not used'
360  end if
361  !print *,interpolator%bc_type
362  stop
363  sll_assert(interpolator%npts > 0)
364  end subroutine set_coefficients_hi1d
365 
366  function get_coefficients_hi1d(interpolator)
367  class(sll_hermite_interpolator_1d), intent(in) :: interpolator
368  sll_real64, dimension(:), pointer :: get_coefficients_hi1d
369 
370  print *, 'get_coefficients_hi1d(): ERROR: This function has not been ', &
371  'implemented yet.'
372  sll_assert(interpolator%npts > 0)
373  get_coefficients_hi1d => null()
374  stop
375  end function get_coefficients_hi1d
376 
real(kind=f64) function, public sll_f_interpolate_value_hermite_1d(eta, interp)
type(sll_t_hermite_interpolation_1d) function, pointer, public sll_f_new_hermite_interpolation_1d(npts, eta_min, eta_max, degree, eta_hermite_continuity, eta_bc_type, const_eta_min_slope, const_eta_max_slope, eta_min_slopes, eta_max_slopes)
subroutine, public sll_s_compute_interpolants_hermite_1d(interp, f)
Interpolator class and methods of hermite 1D interpolator.
real(kind=f64) function, dimension(:), pointer get_coefficients_hi1d(interpolator)
real(kind=f64) function wrap_interpolate_value_hermite_1d(interpolator, eta1)
subroutine set_coefficients_hi1d(interpolator, coeffs)
subroutine wrap_compute_interpolants_hermite_1d(interpolator, data_array, eta_coords, size_eta_coords)
subroutine interpolate_array_derivatives_hi1d(interpolator, num_pts, vals_to_interpolate, output_array)
subroutine interpolate_array_disp_hi1d(this, num_pts, data, alpha, output_array)
subroutine interpolate_array_disp_inplace_hi1d(this, num_pts, data, alpha)
subroutine interpolate_array_values_hi1d(interpolator, num_pts, vals_to_interpolate, output_array)
real(kind=f64) function interpolate_derivative_eta1_hi1d(interpolator, eta1)
subroutine wrap_interpolate_array_hermite_1d(this, num_pts, data, coordinates, output_array)
type(sll_hermite_interpolator_1d) function, pointer, public sll_f_new_hermite_interpolator_1d(npts, eta_min, eta_max, degree, eta_hermite_continuity, eta_bc_type, const_eta_min_slope, const_eta_max_slope, eta_min_slopes, eta_max_slopes)
subroutine initialize_hermite_interpolator_1d(interpolator, npts, eta_min, eta_max, degree, eta_hermite_continuity, eta_bc_type, const_eta_min_slope, const_eta_max_slope, eta_min_slopes, eta_max_slopes)
Module for 1D interpolation and reconstruction.
The hermite-based interpolator is only a wrapper around the capabilities of the hermite interpolation...
Abstract class for 1D interpolation and reconstruction.
    Report Typos and Errors