Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
164 changes: 164 additions & 0 deletions src/bspline_oo_module.f90
Original file line number Diff line number Diff line change
Expand Up @@ -21,6 +21,10 @@ module bspline_oo_module

private

integer,parameter :: int_size = storage_size(1) !! size of a default integer [bits]
integer,parameter :: logical_size = storage_size(.true.) !! size of a default logical [bits]
integer,parameter :: real_size = storage_size(1.0_wp) !! size of a real(wp) [bits]

type,public,abstract :: bspline_class
!! Base class for the b-spline types
private
Expand All @@ -31,17 +35,27 @@ module bspline_oo_module
private
procedure,non_overridable :: destroy_base !! destructor for the abstract type
procedure(destroy_func),deferred,public :: destroy !! destructor
procedure(size_func),deferred,public :: size_of !! size of the structure in bits
procedure,public,non_overridable :: status_ok !! returns true if the last `iflag` status code was `=0`.
procedure,public,non_overridable :: status_message => get_bspline_status_message !! retrieve the last status message
procedure,public,non_overridable :: clear_flag => clear_bspline_flag !! to reset the `iflag` saved in the class.
end type bspline_class

abstract interface

pure subroutine destroy_func(me) !! interface for bspline destructor routines
import :: bspline_class
implicit none
class(bspline_class),intent(inout) :: me
end subroutine destroy_func

pure function size_func(me) result(s) !! interface for size routines
import :: bspline_class
implicit none
class(bspline_class),intent(in) :: me
integer :: s !! size of the structure in bits
end function size_func

end interface

type,extends(bspline_class),public :: bspline_1d
Expand All @@ -58,6 +72,7 @@ end subroutine destroy_func
procedure :: initialize_1d_specify_knots
procedure,public :: evaluate => evaluate_1d
procedure,public :: destroy => destroy_1d
procedure,public :: size_of => size_1d
final :: finalize_1d
end type bspline_1d

Expand All @@ -80,6 +95,7 @@ end subroutine destroy_func
procedure :: initialize_2d_specify_knots
procedure,public :: evaluate => evaluate_2d
procedure,public :: destroy => destroy_2d
procedure,public :: size_of => size_2d
final :: finalize_2d
end type bspline_2d

Expand Down Expand Up @@ -107,6 +123,7 @@ end subroutine destroy_func
procedure :: initialize_3d_specify_knots
procedure,public :: evaluate => evaluate_3d
procedure,public :: destroy => destroy_3d
procedure,public :: size_of => size_3d
final :: finalize_3d
end type bspline_3d

Expand Down Expand Up @@ -139,6 +156,7 @@ end subroutine destroy_func
procedure :: initialize_4d_specify_knots
procedure,public :: evaluate => evaluate_4d
procedure,public :: destroy => destroy_4d
procedure,public :: size_of => size_4d
final :: finalize_4d
end type bspline_4d

Expand Down Expand Up @@ -176,6 +194,7 @@ end subroutine destroy_func
procedure :: initialize_5d_specify_knots
procedure,public :: evaluate => evaluate_5d
procedure,public :: destroy => destroy_5d
procedure,public :: size_of => size_5d
final :: finalize_5d
end type bspline_5d

Expand Down Expand Up @@ -218,6 +237,7 @@ end subroutine destroy_func
procedure :: initialize_6d_specify_knots
procedure,public :: evaluate => evaluate_6d
procedure,public :: destroy => destroy_6d
procedure,public :: size_of => size_6d
final :: finalize_6d
end type bspline_6d

Expand Down Expand Up @@ -329,6 +349,150 @@ pure function get_bspline_status_message(me,iflag) result(msg)
end function get_bspline_status_message
!*****************************************************************************************

!*****************************************************************************************
!>
! Actual size of a [[bspline_1d]] structure in bits.

pure function size_1d(me) result(s)

implicit none

class(bspline_1d),intent(in) :: me
integer :: s !! size of the structure in bits

s = 2*int_size + logical_size + 2*int_size

if (allocated(me%bcoef)) s = s + real_size*size(me%bcoef)
if (allocated(me%tx)) s = s + real_size*size(me%tx)

end function size_1d
!*****************************************************************************************

!*****************************************************************************************
!>
! Actual size of a [[bspline_2d]] structure in bits.

pure function size_2d(me) result(s)

implicit none

class(bspline_2d),intent(in) :: me
integer :: s !! size of the structure in bits

s = 2*int_size + logical_size + 6*int_size

if (allocated(me%bcoef)) s = s + real_size*size(me%bcoef,1)*&
size(me%bcoef,2)
if (allocated(me%tx)) s = s + real_size*size(me%tx)
if (allocated(me%ty)) s = s + real_size*size(me%ty)

end function size_2d
!*****************************************************************************************

!*****************************************************************************************
!>
! Actual size of a [[bspline_3d]] structure in bits.

pure function size_3d(me) result(s)

implicit none

class(bspline_3d),intent(in) :: me
integer :: s !! size of the structure in bits

s = 2*int_size + logical_size + 10*int_size

if (allocated(me%bcoef)) s = s + real_size*size(me%bcoef,1)*&
size(me%bcoef,2)*&
size(me%bcoef,3)
if (allocated(me%tx)) s = s + real_size*size(me%tx)
if (allocated(me%ty)) s = s + real_size*size(me%ty)
if (allocated(me%tz)) s = s + real_size*size(me%tz)

end function size_3d
!*****************************************************************************************

!*****************************************************************************************
!>
! Actual size of a [[bspline_4d]] structure in bits.

pure function size_4d(me) result(s)

implicit none

class(bspline_4d),intent(in) :: me
integer :: s !! size of the structure in bits

s = 2*int_size + logical_size + 14*int_size

if (allocated(me%bcoef)) s = s + real_size*size(me%bcoef,1)*&
size(me%bcoef,2)*&
size(me%bcoef,3)*&
size(me%bcoef,4)
if (allocated(me%tx)) s = s + real_size*size(me%tx)
if (allocated(me%ty)) s = s + real_size*size(me%ty)
if (allocated(me%tz)) s = s + real_size*size(me%tz)
if (allocated(me%tq)) s = s + real_size*size(me%tq)

end function size_4d
!*****************************************************************************************

!*****************************************************************************************
!>
! Actual size of a [[bspline_5d]] structure in bits.

pure function size_5d(me) result(s)

implicit none

class(bspline_5d),intent(in) :: me
integer :: s !! size of the structure in bits

s = 2*int_size + logical_size + 18*int_size

if (allocated(me%bcoef)) s = s + real_size*size(me%bcoef,1)*&
size(me%bcoef,2)*&
size(me%bcoef,3)*&
size(me%bcoef,4)*&
size(me%bcoef,5)
if (allocated(me%tx)) s = s + real_size*size(me%tx)
if (allocated(me%ty)) s = s + real_size*size(me%ty)
if (allocated(me%tz)) s = s + real_size*size(me%tz)
if (allocated(me%tq)) s = s + real_size*size(me%tq)
if (allocated(me%tr)) s = s + real_size*size(me%tr)

end function size_5d
!*****************************************************************************************

!*****************************************************************************************
!>
! Actual size of a [[bspline_6d]] structure in bits.

pure function size_6d(me) result(s)

implicit none

class(bspline_6d),intent(in) :: me
integer :: s !! size of the structure in bits

s = 2*int_size + logical_size + 22*int_size

if (allocated(me%bcoef)) s = s + real_size*size(me%bcoef,1)*&
size(me%bcoef,2)*&
size(me%bcoef,3)*&
size(me%bcoef,4)*&
size(me%bcoef,5)*&
size(me%bcoef,6)
if (allocated(me%tx)) s = s + real_size*size(me%tx)
if (allocated(me%ty)) s = s + real_size*size(me%ty)
if (allocated(me%tz)) s = s + real_size*size(me%tz)
if (allocated(me%tq)) s = s + real_size*size(me%tq)
if (allocated(me%tr)) s = s + real_size*size(me%tr)
if (allocated(me%ts)) s = s + real_size*size(me%ts)

end function size_6d
!*****************************************************************************************

!*****************************************************************************************
!>
! Destructor for contents of the base [[bspline_class]] class.
Expand Down
9 changes: 9 additions & 0 deletions src/tests/test_oo.f90
Original file line number Diff line number Diff line change
Expand Up @@ -107,6 +107,15 @@ program bspline_oo_test
end do
end if

write(*,*) ''
write(*,*) 'size of 1d structure: ', s1%size_of()*8, 'bytes'
write(*,*) 'size of 2d structure: ', s2%size_of()*8, 'bytes'
write(*,*) 'size of 3d structure: ', s3%size_of()*8, 'bytes'
write(*,*) 'size of 4d structure: ', s4%size_of()*8, 'bytes'
write(*,*) 'size of 5d structure: ', s5%size_of()*8, 'bytes'
write(*,*) 'size of 6d structure: ', s6%size_of()*8, 'bytes'
write(*,*) ''

! compute max error at interpolation points

errmax = 0.0_wp
Expand Down