Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_hermite_interpolator_2d.F90
Go to the documentation of this file.
1 !**************************************************************
2 ! Copyright INRIA
3 ! Authors :
4 ! CALVI project team
5 !
6 ! This code SeLaLib (for Semi-Lagrangian-Library)
7 ! is a parallel library for simulating the plasma turbulence
8 ! in a tokamak.
9 !
10 ! This software is governed by the CeCILL-B license
11 ! under French law and abiding by the rules of distribution
12 ! of free software. You can use, modify and redistribute
13 ! the software under the terms of the CeCILL-B license as
14 ! circulated by CEA, CNRS and INRIA at the following URL
15 ! "http://www.cecill.info".
16 !**************************************************************
17 
32 
33 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
34 #include "sll_assert.h"
35 #include "sll_errors.h"
36 #include "sll_memory.h"
37 #include "sll_working_precision.h"
38 
44 
45  use sll_m_interpolators_2d_base, only: &
47 
48  implicit none
49 
50  public :: &
52 
53  private
54 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
55 
66  type(sll_t_hermite_interpolation_2d), pointer :: hermite
68  sll_int32 :: npts1
70  sll_int32 :: npts2
71  contains
73  procedure, pass(interpolator) :: init => initialize_hermite_interpolator_2d
75  procedure :: compute_interpolants => wrap_compute_interpolants_hermite_2d
77  procedure :: interpolate_from_interpolant_value => wrap_interpolate_value_hermite_2d
79  procedure :: interpolate_from_interpolant_derivative_eta1 => wrap_interpolate_deriv1_hermite_2d
81  procedure :: interpolate_from_interpolant_derivative_eta2 => wrap_interpolate_deriv2_hermite_2d
83  procedure, pass :: interpolate_array => wrap_interpolate_array_hermite_2d
85  procedure, pass :: interpolate_array_disp => wrap_interpolate2d_disp_hermite_2d
87  procedure, pass :: set_coefficients => wrap_set_coefficients_hermite_2d
89  procedure, pass :: get_coefficients => wrap_get_coefficients_hermite_2d
91  procedure, pass :: coefficients_are_set => wrap_coefficients_are_set_hermite_2d
93  procedure, pass :: delete => delete_sll_hermite_interpolator_2d
95 
96 contains
97 
100  npts1, &
101  npts2, &
102  eta1_min, &
103  eta1_max, &
104  eta2_min, &
105  eta2_max, &
106  degree1, &
107  degree2, &
108  eta1_hermite_continuity, &
109  eta2_hermite_continuity, &
110  eta1_bc_type, &
111  eta2_bc_type, &
112  const_eta1_min_slope, &
113  const_eta1_max_slope, &
114  const_eta2_min_slope, &
115  const_eta2_max_slope, &
116  eta1_min_slopes, &
117  eta1_max_slopes, &
118  eta2_min_slopes, &
119  eta2_max_slopes) &
120  result(interpolator)
121 
122  type(sll_hermite_interpolator_2d), pointer :: interpolator
123  sll_int32, intent(in) :: npts1
124  sll_int32, intent(in) :: npts2
125  sll_real64, intent(in) :: eta1_min
126  sll_real64, intent(in) :: eta1_max
127  sll_real64, intent(in) :: eta2_min
128  sll_real64, intent(in) :: eta2_max
129  sll_int32, intent(in) :: degree1
130  sll_int32, intent(in) :: degree2
131  sll_int32, intent(in) :: eta1_hermite_continuity
132  sll_int32, intent(in) :: eta2_hermite_continuity
133  sll_int32, intent(in) :: eta1_bc_type
134  sll_int32, intent(in) :: eta2_bc_type
135  sll_real64, intent(in), optional :: const_eta1_min_slope
136  sll_real64, intent(in), optional :: const_eta1_max_slope
137  sll_real64, intent(in), optional :: const_eta2_min_slope
138  sll_real64, intent(in), optional :: const_eta2_max_slope
139  sll_real64, dimension(:), intent(in), optional :: eta1_min_slopes
140  sll_real64, dimension(:), intent(in), optional :: eta1_max_slopes
141  sll_real64, dimension(:), intent(in), optional :: eta2_min_slopes
142  sll_real64, dimension(:), intent(in), optional :: eta2_max_slopes
143  sll_int32 :: ierr
144 
145  sll_allocate(interpolator, ierr)
146 
147  interpolator%npts1 = npts1
148  interpolator%npts2 = npts2
149 
150  call interpolator%init( &
151  npts1, &
152  npts2, &
153  eta1_min, &
154  eta1_max, &
155  eta2_min, &
156  eta2_max, &
157  degree1, &
158  degree2, &
159  eta1_hermite_continuity, &
160  eta2_hermite_continuity, &
161  eta1_bc_type, &
162  eta2_bc_type, &
163  const_eta1_min_slope, &
164  const_eta1_max_slope, &
165  const_eta2_min_slope, &
166  const_eta2_max_slope, &
167  eta1_min_slopes, &
168  eta1_max_slopes, &
169  eta2_min_slopes, &
170  eta2_max_slopes)
171 
173 
175  interpolator, &
176  npts1, &
177  npts2, &
178  eta1_min, &
179  eta1_max, &
180  eta2_min, &
181  eta2_max, &
182  degree1, &
183  degree2, &
184  eta1_hermite_continuity, &
185  eta2_hermite_continuity, &
186  eta1_bc_type, &
187  eta2_bc_type, &
188  const_eta1_min_slope, &
189  const_eta1_max_slope, &
190  const_eta2_min_slope, &
191  const_eta2_max_slope, &
192  eta1_min_slopes, &
193  eta1_max_slopes, &
194  eta2_min_slopes, &
195  eta2_max_slopes)
196 
197  class(sll_hermite_interpolator_2d), intent(inout) :: interpolator
198  sll_int32, intent(in) :: npts1
199  sll_int32, intent(in) :: npts2
200  sll_real64, intent(in) :: eta1_min
201  sll_real64, intent(in) :: eta1_max
202  sll_real64, intent(in) :: eta2_min
203  sll_real64, intent(in) :: eta2_max
204  sll_int32, intent(in) :: degree1
205  sll_int32, intent(in) :: degree2
206  sll_int32, intent(in) :: eta1_hermite_continuity
207  sll_int32, intent(in) :: eta2_hermite_continuity
208  sll_int32, intent(in) :: eta1_bc_type
209  sll_int32, intent(in) :: eta2_bc_type
210  sll_real64, intent(in), optional :: const_eta1_min_slope
211  sll_real64, intent(in), optional :: const_eta1_max_slope
212  sll_real64, intent(in), optional :: const_eta2_min_slope
213  sll_real64, intent(in), optional :: const_eta2_max_slope
214  sll_real64, dimension(:), intent(in), optional :: eta1_min_slopes
215  sll_real64, dimension(:), intent(in), optional :: eta1_max_slopes
216  sll_real64, dimension(:), intent(in), optional :: eta2_min_slopes
217  sll_real64, dimension(:), intent(in), optional :: eta2_max_slopes
218 
219  interpolator%hermite => sll_f_new_hermite_interpolation_2d( &
220  npts1, &
221  npts2, &
222  eta1_min, &
223  eta1_max, &
224  eta2_min, &
225  eta2_max, &
226  degree1, &
227  degree2, &
228  eta1_hermite_continuity, &
229  eta2_hermite_continuity, &
230  eta1_bc_type, &
231  eta2_bc_type, &
232  const_eta1_min_slope, &
233  const_eta1_max_slope, &
234  const_eta2_min_slope, &
235  const_eta2_max_slope, &
236  eta1_min_slopes, &
237  eta1_max_slopes, &
238  eta2_min_slopes, &
239  eta2_max_slopes)
240 
242 
244  interpolator, &
245  data_array, &
246  eta1_coords, &
247  size_eta1_coords, &
248  eta2_coords, &
249  size_eta2_coords)
250  class(sll_hermite_interpolator_2d), intent(inout) :: interpolator
251  sll_real64, dimension(:, :), intent(in) :: data_array
252  sll_real64, dimension(:), intent(in), optional :: eta1_coords
253  sll_real64, dimension(:), intent(in), optional :: eta2_coords
254  sll_int32, intent(in), optional :: size_eta1_coords
255  sll_int32, intent(in), optional :: size_eta2_coords
256 
257  if (present(eta1_coords) .or. present(eta2_coords) &
258  .or. present(size_eta1_coords) .or. present(size_eta2_coords)) then
259  sll_error('wrap_compute_interpolants_hermite_2d', 'This case is not yet implemented')
260  end if
261 
262  call sll_s_compute_interpolants_hermite_2d(interpolator%hermite, data_array)
263 
265 
266  function wrap_interpolate_value_hermite_2d(interpolator, eta1, eta2) result(val)
267  class(sll_hermite_interpolator_2d), intent(in) :: interpolator
268  sll_real64 :: val
269  sll_real64, intent(in) :: eta1
270  sll_real64, intent(in) :: eta2
271  val = sll_f_interpolate_value_hermite_2d(eta1, eta2, interpolator%hermite)
272 
274 
275  function wrap_interpolate_deriv1_hermite_2d(interpolator, eta1, eta2) result(val)
276  class(sll_hermite_interpolator_2d), intent(in) :: interpolator
277  sll_real64 :: val
278  sll_real64, intent(in) :: eta1
279  sll_real64, intent(in) :: eta2
280  val = 0._f64 + eta1 + eta2
281  print *, '#wrap_interpolate_deriv1_hermite_2d'
282  print *, '#not implemented for the moment'
283  stop
284  sll_assert(interpolator%npts1 > 0)
285  !interpolate_x1_derivative_2D(eta1,eta2,interpolator%spline)
287 
288  function wrap_interpolate_deriv2_hermite_2d(interpolator, eta1, eta2) result(val)
289  class(sll_hermite_interpolator_2d), intent(in) :: interpolator
290  sll_real64 :: val
291  sll_real64, intent(in) :: eta1
292  sll_real64, intent(in) :: eta2
293  val = 0._f64 + eta1 + eta2
294  print *, '#wrap_interpolate_deriv1_hermite_2d'
295  print *, '#not implemented for the moment'
296  stop
297  sll_assert(interpolator%npts1 > 0)
298  !interpolate_x1_derivative_2D(eta1,eta2,interpolator%spline)
300 
302  this, &
303  num_points1, &
304  num_points2, &
305  data_in, &
306  eta1, &
307  eta2, &
308  data_out)
309  class(sll_hermite_interpolator_2d), intent(inout) :: this
310  sll_int32, intent(in) :: num_points1
311  sll_int32, intent(in) :: num_points2
312  sll_real64, dimension(:, :), intent(in) :: eta1
313  sll_real64, dimension(:, :), intent(in) :: eta2
314  sll_real64, dimension(:, :), intent(in) :: data_in
315  sll_real64, intent(out) :: data_out(num_points1, num_points2)
316  sll_int32 :: i
317  sll_int32 :: j
318  call sll_s_compute_interpolants_hermite_2d(this%hermite, data_in)
319  do j = 1, num_points2
320  do i = 1, num_points1
321  data_out(i, j) = this%interpolate_from_interpolant_value(eta1(i, j), eta2(i, j))
322  end do
323  end do
324  !print *,'#wrap_interpolate_array_hermite_2d'
325  !print *,'#not implemented for the moment'
326  !stop
327  end subroutine wrap_interpolate_array_hermite_2d
328 
330  this, &
331  num_points1, &
332  num_points2, &
333  data_in, &
334  alpha1, &
335  alpha2, &
336  data_out)
337  class(sll_hermite_interpolator_2d), intent(inout) :: this
338  sll_int32, intent(in) :: num_points1
339  sll_int32, intent(in) :: num_points2
340  sll_real64, dimension(:, :), intent(in) :: alpha1
341  sll_real64, dimension(:, :), intent(in) :: alpha2
342  sll_real64, dimension(:, :), intent(in) :: data_in
343  sll_real64, intent(out) :: data_out(num_points1, num_points2)
344  print *, '#wrap_interpolate2d_disp_hermite_2d'
345  print *, '#not implemented for the moment'
346  data_out = 0.0_f64*data_in
347  stop
348  sll_assert(size(alpha1, 1) == num_points1)
349  sll_assert(size(alpha2, 1) == num_points1)
350  sll_assert(this%npts1 == num_points1)
352 
354  interpolator, &
355  coeffs_1d, &
356  coeffs_2d, &
357  coeff2d_size1, &
358  coeff2d_size2, &
359  knots1, &
360  size_knots1, &
361  knots2, &
362  size_knots2)
363  class(sll_hermite_interpolator_2d), intent(inout) :: interpolator
364  sll_real64, dimension(:), intent(in), optional :: coeffs_1d
365  sll_real64, dimension(:, :), intent(in), optional :: coeffs_2d
366  sll_int32, intent(in), optional :: coeff2d_size1
367  sll_int32, intent(in), optional :: coeff2d_size2
368  sll_real64, dimension(:), intent(in), optional :: knots1
369  sll_real64, dimension(:), intent(in), optional :: knots2
370  sll_int32, intent(in), optional :: size_knots1
371  sll_int32, intent(in), optional :: size_knots2
372  print *, 'wrap_set_coefficients_hermite_2d(): ERROR: This function has not been ', &
373  'implemented yet.'
374  print *, interpolator%npts1
375  if (present(coeffs_1d)) then
376  print *, 'coeffs_1d present but not used'
377  end if
378  if (present(coeffs_2d)) then
379  print *, 'coeffs_2d present but not used'
380  end if
381  if (present(coeff2d_size1)) then
382  print *, 'coeff2d_size1 present but not used'
383  end if
384  if (present(coeff2d_size2)) then
385  print *, 'coeff2d_size2 present but not used'
386  end if
387  if (present(knots1)) then
388  print *, 'knots1 present but not used'
389  end if
390  if (present(knots2)) then
391  print *, 'knots2 present but not used'
392  end if
393  if (present(size_knots1)) then
394  print *, 'size_knots1 present but not used'
395  end if
396  if (present(size_knots2)) then
397  print *, 'size_knots2 present but not used'
398  end if
399  stop
400  end subroutine wrap_set_coefficients_hermite_2d
401 
402  function wrap_get_coefficients_hermite_2d(interpolator) result(res)
403  class(sll_hermite_interpolator_2d), intent(in) :: interpolator
404  sll_real64, dimension(:, :), pointer :: res
405 
406  print *, 'wrap_get_coefficients_hermite_2d: ERROR: This function has not been ', &
407  'implemented yet.'
408  res => null()
409  print *, interpolator%npts1
410  stop
412 
413  function wrap_coefficients_are_set_hermite_2d(interpolator) result(res)
414  class(sll_hermite_interpolator_2d), intent(in) :: interpolator
415  logical :: res
416  res = .false.
417  print *, 'wrap_coefficients_are_set_hermite_2d(): '
418  print *, 'this function has not been implemented yet.'
419  print *, '#', interpolator%npts1
420  !stop
422 
423  subroutine delete_sll_hermite_interpolator_2d(interpolator)
424  class(sll_hermite_interpolator_2d), intent(inout) :: interpolator
425  print *, '#warning delete_sll_hermite_interpolator_2d'
426  print *, '#not implemented for the moment'
427  sll_assert(interpolator%npts1 > 0)
429 
subroutine, public sll_s_compute_interpolants_hermite_2d(interp, f)
type(sll_t_hermite_interpolation_2d) function, pointer, public sll_f_new_hermite_interpolation_2d(npts1, npts2, eta1_min, eta1_max, eta2_min, eta2_max, degree1, degree2, eta1_hermite_continuity, eta2_hermite_continuity, eta1_bc_type, eta2_bc_type, const_eta1_min_slope, const_eta1_max_slope, const_eta2_min_slope, const_eta2_max_slope, eta1_min_slopes, eta1_max_slopes, eta2_min_slopes, eta2_max_slopes)
real(kind=f64) function, public sll_f_interpolate_value_hermite_2d(eta1, eta2, interp)
subroutine delete_sll_hermite_interpolator_2d(interpolator)
type(sll_hermite_interpolator_2d) function, pointer, public sll_f_new_hermite_interpolator_2d(npts1, npts2, eta1_min, eta1_max, eta2_min, eta2_max, degree1, degree2, eta1_hermite_continuity, eta2_hermite_continuity, eta1_bc_type, eta2_bc_type, const_eta1_min_slope, const_eta1_max_slope, const_eta2_min_slope, const_eta2_max_slope, eta1_min_slopes, eta1_max_slopes, eta2_min_slopes, eta2_max_slopes)
real(kind=f64) function, dimension(:, :), pointer wrap_get_coefficients_hermite_2d(interpolator)
logical function wrap_coefficients_are_set_hermite_2d(interpolator)
subroutine wrap_set_coefficients_hermite_2d(interpolator, coeffs_1d, coeffs_2d, coeff2d_size1, coeff2d_size2, knots1, size_knots1, knots2, size_knots2)
subroutine initialize_hermite_interpolator_2d(interpolator, npts1, npts2, eta1_min, eta1_max, eta2_min, eta2_max, degree1, degree2, eta1_hermite_continuity, eta2_hermite_continuity, eta1_bc_type, eta2_bc_type, const_eta1_min_slope, const_eta1_max_slope, const_eta2_min_slope, const_eta2_max_slope, eta1_min_slopes, eta1_max_slopes, eta2_min_slopes, eta2_max_slopes)
subroutine wrap_compute_interpolants_hermite_2d(interpolator, data_array, eta1_coords, size_eta1_coords, eta2_coords, size_eta2_coords)
real(kind=f64) function wrap_interpolate_deriv2_hermite_2d(interpolator, eta1, eta2)
real(kind=f64) function wrap_interpolate_deriv1_hermite_2d(interpolator, eta1, eta2)
subroutine wrap_interpolate_array_hermite_2d(this, num_points1, num_points2, data_in, eta1, eta2, data_out)
real(kind=f64) function wrap_interpolate_value_hermite_2d(interpolator, eta1, eta2)
subroutine wrap_interpolate2d_disp_hermite_2d(this, num_points1, num_points2, data_in, alpha1, alpha2, data_out)
abstract data type for 2d interpolation
The hermite-based interpolator is only a wrapper around the capabilities of the hermite interpolation...
Base class/basic interface for 2D interpolators.
    Report Typos and Errors