25 #include "sll_assert.h"
26 #include "sll_errors.h"
27 #include "sll_memory.h"
28 #include "sll_working_precision.h"
36 sll_f_cubic_spline_2d_get_x1_delta, &
37 sll_f_cubic_spline_2d_get_x1_max, &
38 sll_f_cubic_spline_2d_get_x1_min, &
39 sll_f_cubic_spline_2d_get_x2_delta, &
40 sll_f_cubic_spline_2d_get_x2_max, &
41 sll_f_cubic_spline_2d_get_x2_min, &
81 sll_real64,
dimension(:, :),
pointer :: interpolation_points
136 const_eta1_min_slope, &
137 const_eta1_max_slope, &
138 const_eta2_min_slope, &
139 const_eta2_max_slope, &
148 sll_int32,
intent(in) :: npts1
149 sll_int32,
intent(in) :: npts2
150 sll_real64,
intent(in) :: eta1_min
151 sll_real64,
intent(in) :: eta1_max
152 sll_real64,
intent(in) :: eta2_min
153 sll_real64,
intent(in) :: eta2_max
154 sll_int32,
intent(in) :: eta1_bc_type
155 sll_int32,
intent(in) :: eta2_bc_type
156 sll_real64,
intent(in),
optional :: const_eta1_min_slope
157 sll_real64,
intent(in),
optional :: const_eta1_max_slope
158 sll_real64,
intent(in),
optional :: const_eta2_min_slope
159 sll_real64,
intent(in),
optional :: const_eta2_max_slope
160 sll_real64,
dimension(:),
intent(in),
optional :: eta1_min_slopes
161 sll_real64,
dimension(:),
intent(in),
optional :: eta1_max_slopes
162 sll_real64,
dimension(:),
intent(in),
optional :: eta2_min_slopes
163 sll_real64,
dimension(:),
intent(in),
optional :: eta2_max_slopes
166 sll_allocate(interpolator, ierr)
168 call interpolator%init( &
177 const_eta1_min_slope, &
178 const_eta1_max_slope, &
179 const_eta2_min_slope, &
180 const_eta2_max_slope, &
205 const_eta1_min_slope, &
206 const_eta1_max_slope, &
207 const_eta2_min_slope, &
208 const_eta2_max_slope, &
215 sll_int32,
intent(in) :: npts1
216 sll_int32,
intent(in) :: npts2
217 sll_real64,
intent(in) :: eta1_min
218 sll_real64,
intent(in) :: eta1_max
219 sll_real64,
intent(in) :: eta2_min
220 sll_real64,
intent(in) :: eta2_max
221 sll_int32,
intent(in) :: eta1_bc_type
222 sll_int32,
intent(in) :: eta2_bc_type
223 sll_real64,
intent(in),
optional :: const_eta1_min_slope
224 sll_real64,
intent(in),
optional :: const_eta1_max_slope
225 sll_real64,
intent(in),
optional :: const_eta2_min_slope
226 sll_real64,
intent(in),
optional :: const_eta2_max_slope
227 sll_real64,
dimension(:),
intent(in),
optional :: eta1_min_slopes
228 sll_real64,
dimension(:),
intent(in),
optional :: eta1_max_slopes
229 sll_real64,
dimension(:),
intent(in),
optional :: eta2_min_slopes
230 sll_real64,
dimension(:),
intent(in),
optional :: eta2_max_slopes
232 interpolator%npts1 = npts1
233 interpolator%npts2 = npts2
234 interpolator%bc_type1 = eta1_bc_type
235 interpolator%bc_type2 = eta2_bc_type
238 interpolator%spline, &
247 const_slope_x1_min=const_eta1_min_slope, &
248 const_slope_x1_max=const_eta1_max_slope, &
249 const_slope_x2_min=const_eta2_min_slope, &
250 const_slope_x2_max=const_eta2_max_slope, &
251 x1_min_slopes=eta1_min_slopes, &
252 x1_max_slopes=eta1_max_slopes, &
253 x2_min_slopes=eta2_min_slopes, &
254 x2_max_slopes=eta2_max_slopes)
266 sll_real64,
dimension(:, :),
intent(in) :: data_array
267 sll_real64,
dimension(:),
intent(in),
optional :: eta1_coords
268 sll_real64,
dimension(:),
intent(in),
optional :: eta2_coords
269 sll_int32,
intent(in),
optional :: size_eta1_coords
270 sll_int32,
intent(in),
optional :: size_eta2_coords
272 if (
present(eta1_coords) .or.
present(eta2_coords) &
273 .or.
present(size_eta1_coords) .or.
present(size_eta2_coords))
then
274 sll_error(
'compute_interpolants_cs2d',
'This case is not yet implemented')
284 sll_real64,
intent(in) :: eta1
285 sll_real64,
intent(in) :: eta2
292 sll_real64,
intent(in) :: eta1
293 sll_real64,
intent(in) :: eta2
300 sll_real64,
intent(in) :: eta1
301 sll_real64,
intent(in) :: eta2
308 eta1, eta2, data_out)
311 sll_int32,
intent(in) :: num_points1
312 sll_int32,
intent(in) :: num_points2
313 sll_real64,
dimension(:, :),
intent(in) :: eta1
314 sll_real64,
dimension(:, :),
intent(in) :: eta2
315 sll_real64,
dimension(:, :),
intent(in) :: data_in
316 sll_real64,
intent(out) :: data_out(num_points1, num_points2)
321 do j = 1, num_points2
322 do i = 1, num_points1
323 data_out(i, j) = this%interpolate_from_interpolant_value(eta1(i, j), eta2(i, j))
339 sll_int32,
intent(in) :: num_points1
340 sll_int32,
intent(in) :: num_points2
341 sll_real64,
dimension(:, :),
intent(in) :: alpha1
342 sll_real64,
dimension(:, :),
intent(in) :: alpha2
343 sll_real64,
dimension(:, :),
intent(in) :: data_in
344 sll_real64,
intent(out) :: data_out(num_points1, num_points2)
346 sll_real64 :: eta1_min
347 sll_real64 :: eta1_max
348 sll_real64 :: delta_eta1
350 sll_real64 :: eta2_min
351 sll_real64 :: eta2_max
352 sll_real64 :: delta_eta2
356 eta1_min = sll_f_cubic_spline_2d_get_x1_min(this%spline)
357 eta1_max = sll_f_cubic_spline_2d_get_x1_max(this%spline)
358 eta2_min = sll_f_cubic_spline_2d_get_x2_min(this%spline)
359 eta2_max = sll_f_cubic_spline_2d_get_x2_max(this%spline)
360 delta_eta1 = sll_f_cubic_spline_2d_get_x1_delta(this%spline)
361 delta_eta2 = sll_f_cubic_spline_2d_get_x2_delta(this%spline)
365 if (this%bc_type1 == sll_p_periodic .and. &
366 this%bc_type2 == sll_p_periodic)
then
368 do j = 1, num_points2
369 do i = 1, num_points1
370 eta1 = eta1_min + (i - 1)*delta_eta1
371 eta2 = eta2_min + (j - 1)*delta_eta2
373 modulo(eta1 - eta1_min + alpha1(i, j), eta1_max - eta1_min)
375 modulo(eta2 - eta2_min + alpha2(i, j), eta2_max - eta2_min)
376 data_out(i, j) = this%interpolate_from_interpolant_value(eta1, eta2)
380 else if (this%bc_type1 == sll_p_hermite .and. &
381 this%bc_type2 == sll_p_hermite)
then
383 do j = 1, num_points2
384 do i = 1, num_points1
385 eta1 = eta1_min + (i - 1)*delta_eta1 + alpha1(i, j)
386 eta2 = eta2_min + (j - 1)*delta_eta2 + alpha2(i, j)
387 eta1 = min(eta1, eta1_max)
388 eta2 = min(eta2, eta2_max)
389 eta1 = max(eta1, eta1_min)
390 eta2 = max(eta2, eta2_min)
391 data_out(i, j) = this%interpolate_from_interpolant_value(eta1, eta2)
397 do j = 1, num_points2
398 do i = 1, num_points1
399 eta1 = eta1_min + (i - 1)*delta_eta1 + alpha1(i, j)
400 eta2 = eta2_min + (j - 1)*delta_eta2 + alpha2(i, j)
401 sll_assert(eta1_min <= eta1 .and. eta1 <= eta1_max)
402 sll_assert(eta2_min <= eta2 .and. eta2 <= eta2_max)
403 data_out(i, j) = this%interpolate_from_interpolant_value(eta1, eta2)
420 sll_real64,
dimension(:),
intent(in),
optional :: coeffs_1d
421 sll_real64,
dimension(:, :),
intent(in),
optional :: coeffs_2d
423 sll_int32,
intent(in),
optional :: coeff2d_size1
424 sll_int32,
intent(in),
optional :: coeff2d_size2
425 sll_real64,
dimension(:),
intent(in),
optional :: knots1
426 sll_real64,
dimension(:),
intent(in),
optional :: knots2
427 sll_int32,
intent(in),
optional :: size_knots1
428 sll_int32,
intent(in),
optional :: size_knots2
429 print *,
'set_coefficients_cs2d(): ERROR: This function has not been ', &
431 print *, interpolator%npts1
432 if (
present(coeffs_1d))
then
433 print *,
'coeffs_1d present but not used'
435 if (
present(coeffs_2d))
then
436 print *,
'coeffs_2d present but not used'
438 if (
present(coeff2d_size1))
then
439 print *,
'coeff2d_size1 present but not used'
441 if (
present(coeff2d_size2))
then
442 print *,
'coeff2d_size2 present but not used'
444 if (
present(knots1))
then
445 print *,
'knots1 present but not used'
447 if (
present(knots2))
then
448 print *,
'knots2 present but not used'
450 if (
present(size_knots1))
then
451 print *,
'size_knots1 present but not used'
453 if (
present(size_knots2))
then
454 print *,
'size_knots2 present but not used'
482 print *,
'get_coefficients_cs2d(): ERROR: This function has not been ', &
485 print *, interpolator%npts1
493 print *,
'coefficients_are_set_cs2d(): this function has not been implemented yet.'
494 print *,
'#', interpolator%npts1
Deallocate the interpolator object.
Describe different boundary conditions.
Class for the cubic spline sll_c_interpolator_2d.
subroutine compute_interpolants_cs2d(interpolator, data_array, eta1_coords, size_eta1_coords, eta2_coords, size_eta2_coords)
subroutine spline_interpolate2d(this, num_points1, num_points2, data_in, eta1, eta2, data_out)
subroutine initialize_cs2d_interpolator(interpolator, npts1, npts2, eta1_min, eta1_max, eta2_min, eta2_max, 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)
function interpolate_deriv2_cs2d(interpolator, eta1, eta2)
pointer get_coefficients_cs2d(interpolator)
type(sll_t_cubic_spline_interpolator_2d) function, pointer, public sll_f_new_cubic_spline_interpolator_2d(npts1, npts2, eta1_min, eta1_max, eta2_min, eta2_max, 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)
Function that return a pointer to a cubic spline interpolator 2d object. The result can be the target...
logical function coefficients_are_set_cs2d(interpolator)
subroutine set_coefficients_cs2d(interpolator, coeffs_1d, coeffs_2d, coeff2d_size1, coeff2d_size2, knots1, size_knots1, knots2, size_knots2)
subroutine delete_sll_cubic_spline_interpolator_2d(interpolator)
function interpolate_deriv1_cs2d(interpolator, eta1, eta2)
subroutine spline_interpolate2d_disp(this, num_points1, num_points2, data_in, alpha1, alpha2, data_out)
function interpolate_value_cs2d(interpolator, eta1, eta2)
Provides capabilities for data and derivative interpolation with cubic B-splines and different bounda...
subroutine, public sll_s_cubic_spline_2d_compute_interpolant(data, spline)
Computes the spline coefficients for the given data. The coeffcients are first computed in the second...
real(kind=f64) function, public sll_f_cubic_spline_2d_eval_deriv_x2(x1, x2, spline)
Returns the interpolated value of the derivative.
real(kind=f64) function, public sll_f_cubic_spline_2d_eval_deriv_x1(x1, x2, spline)
Returns the interpolated value of the derivative in the x1 direction at the point (x1,...
real(kind=f64) function, public sll_f_cubic_spline_2d_eval(x1, x2, spline)
Returns the interpolated value of the image of the point (x1,x2) using the spline decomposition store...
subroutine, public sll_s_cubic_spline_2d_init(this, num_pts_x1, num_pts_x2, x1_min, x1_max, x2_min, x2_max, x1_bc_type, x2_bc_type, const_slope_x1_min, const_slope_x1_max, const_slope_x2_min, const_slope_x2_max, x1_min_slopes, x1_max_slopes, x2_min_slopes, x2_max_slopes)
subroutine, public sll_s_cubic_spline_2d_free(spline)
abstract data type for 2d interpolation
Pointer to this interpolator derived type.
The spline-based interpolator is only a wrapper around the capabilities of the cubic splines.
Basic type for one-dimensional cubic spline data.
Base class/basic interface for 2D interpolators.