Skip to content

Commit e971904

Browse files
authored
Expose spline BC to C/C++ (#31)
2 parents c7a0f0a + 2d47aa2 commit e971904

File tree

3 files changed

+49
-1
lines changed

3 files changed

+49
-1
lines changed

include/fitpack_curves_c.h

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -59,6 +59,8 @@ FP_FLAG fitpack_curve_c_all_derivatives(fitpack_curve_c *self, FP_REAL x, FP_RE
5959
FP_REAL fitpack_curve_c_smoothing (fitpack_curve_c *self);
6060
FP_REAL fitpack_curve_c_mse (fitpack_curve_c *self);
6161
FP_SIZE fitpack_curve_c_degree (fitpack_curve_c *self);
62+
void fitpack_curve_c_set_bc (fitpack_curve_c *self, FP_FLAG bc);
63+
FP_FLAG fitpack_curve_c_get_bc (fitpack_curve_c *self);
6264

6365
#ifdef __cplusplus
6466
}

include/fpCurve.hpp

Lines changed: 12 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -133,6 +133,18 @@ class fpCurve
133133
return fitpack_curve_c_integral(&cptr, from, to);
134134
}
135135

136+
// Set spline behavior outside the support
137+
void set_bc(FP_FLAG bc)
138+
{
139+
fitpack_curve_c_set_bc(&cptr, bc);
140+
}
141+
142+
// Get spline behavior outside the support
143+
FP_FLAG get_bc()
144+
{
145+
return fitpack_curve_c_get_bc(&cptr);
146+
}
147+
136148
// Get fourier coefficients
137149
FP_FLAG fourier(const vector<FP_REAL> &alpha,
138150
vector<FP_REAL> &A, vector<FP_REAL> &B)

src/fitpack_curves_c.f90

Lines changed: 35 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -24,7 +24,8 @@
2424
! **************************************************************************************************
2525
module fitpack_curves_c
2626
use fitpack_curves, only: fitpack_curve
27-
use fitpack_core, only: FP_SIZE, FP_REAL, FP_FLAG, FP_BOOL
27+
use fitpack_core, only: FP_SIZE, FP_REAL, FP_FLAG, FP_BOOL, &
28+
OUTSIDE_EXTRAPOLATE, OUTSIDE_ZERO, OUTSIDE_NOT_ALLOWED, OUTSIDE_NEAREST_BND
2829
use, intrinsic :: iso_c_binding
2930
implicit none
3031
private
@@ -48,6 +49,8 @@ module fitpack_curves_c
4849
public :: fitpack_curve_c_smoothing
4950
public :: fitpack_curve_c_mse
5051
public :: fitpack_curve_c_degree
52+
public :: fitpack_curve_c_set_bc
53+
public :: fitpack_curve_c_get_bc
5154

5255
!> Opaque-pointer C derived type
5356
type, public, bind(C) :: fitpack_curve_c
@@ -358,5 +361,36 @@ integer(FP_SIZE) function fitpack_curve_c_degree(this) bind(c,name='fitpack_curv
358361

359362
end function fitpack_curve_c_degree
360363

364+
!> Set spline behavior outside the support
365+
subroutine fitpack_curve_c_set_bc(this, bc) bind(c,name='fitpack_curve_c_set_bc')
366+
type(fitpack_curve_c), intent(inout) :: this
367+
integer(FP_FLAG), intent(in), value :: bc
368+
369+
type(fitpack_curve), pointer :: fcurve
370+
371+
!> Get object; allocate it in case
372+
call fitpack_curve_c_pointer(this,fcurve)
373+
374+
!> Validate extrapolation behavior value
375+
if (bc /= OUTSIDE_EXTRAPOLATE .and. bc /= OUTSIDE_ZERO .and. &
376+
bc /= OUTSIDE_NOT_ALLOWED .and. bc /= OUTSIDE_NEAREST_BND) return
377+
378+
fcurve%bc = bc
379+
380+
end subroutine fitpack_curve_c_set_bc
381+
382+
!> Get spline behavior outside the support
383+
integer(FP_FLAG) function fitpack_curve_c_get_bc(this) bind(c,name='fitpack_curve_c_get_bc')
384+
type(fitpack_curve_c), intent(inout) :: this
385+
386+
type(fitpack_curve), pointer :: fcurve
387+
388+
!> Get object; allocate it in case
389+
call fitpack_curve_c_pointer(this,fcurve)
390+
391+
fitpack_curve_c_get_bc = fcurve%bc
392+
393+
end function fitpack_curve_c_get_bc
394+
361395

362396
end module fitpack_curves_c

0 commit comments

Comments
 (0)