forked from ECP-Astro/chimera2castro
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathbspline_oo_module.f90
2212 lines (1748 loc) · 73.6 KB
/
bspline_oo_module.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!*****************************************************************************************
!>
! author: Jacob Williams
! license: BSD
! date: 12/6/2015
!
!# Description
!
! Object-oriented style wrappers to [[bspline_sub_module]].
! This module provides classes ([[bspline_1d]], [[bspline_2d]],
! [[bspline_3d]], [[bspline_4d]], [[bspline_5d]], and [[bspline_6d]])
! which can be used instead of the main subroutine interface.
module bspline_oo_module
use,intrinsic :: iso_fortran_env, only: wp => real64
use,intrinsic :: iso_fortran_env, only: error_unit
use bspline_sub_module
implicit none
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
integer :: inbvx = 1 !! internal variable used by dbvalu for efficient processing
integer :: iflag = 1 !! saved `iflag` from the list routine call.
logical :: initialized = .false. !! true if the class is initialized and ready to use
contains
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
!! Class for 1d b-spline interpolation.
private
integer :: nx = 0
integer :: kx = 0
real(wp),dimension(:),allocatable :: bcoef
real(wp),dimension(:),allocatable :: tx
contains
private
generic,public :: initialize => initialize_1d_auto_knots,initialize_1d_specify_knots
procedure :: initialize_1d_auto_knots
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
type,extends(bspline_class),public :: bspline_2d
!! Class for 2d b-spline interpolation.
private
integer :: nx = 0
integer :: ny = 0
integer :: kx = 0
integer :: ky = 0
real(wp),dimension(:,:),allocatable :: bcoef
real(wp),dimension(:),allocatable :: tx
real(wp),dimension(:),allocatable :: ty
integer :: inbvy = 1
integer :: iloy = 1
contains
private
generic,public :: initialize => initialize_2d_auto_knots,initialize_2d_specify_knots
procedure :: initialize_2d_auto_knots
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
type,extends(bspline_class),public :: bspline_3d
!! Class for 3d b-spline interpolation.
private
integer :: nx = 0
integer :: ny = 0
integer :: nz = 0
integer :: kx = 0
integer :: ky = 0
integer :: kz = 0
real(wp),dimension(:,:,:),allocatable :: bcoef
real(wp),dimension(:),allocatable :: tx
real(wp),dimension(:),allocatable :: ty
real(wp),dimension(:),allocatable :: tz
integer :: inbvy = 1
integer :: inbvz = 1
integer :: iloy = 1
integer :: iloz = 1
contains
private
generic,public :: initialize => initialize_3d_auto_knots,initialize_3d_specify_knots
procedure :: initialize_3d_auto_knots
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
type,extends(bspline_class),public :: bspline_4d
!! Class for 4d b-spline interpolation.
private
integer :: nx = 0
integer :: ny = 0
integer :: nz = 0
integer :: nq = 0
integer :: kx = 0
integer :: ky = 0
integer :: kz = 0
integer :: kq = 0
real(wp),dimension(:,:,:,:),allocatable :: bcoef
real(wp),dimension(:),allocatable :: tx
real(wp),dimension(:),allocatable :: ty
real(wp),dimension(:),allocatable :: tz
real(wp),dimension(:),allocatable :: tq
integer :: inbvy = 1
integer :: inbvz = 1
integer :: inbvq = 1
integer :: iloy = 1
integer :: iloz = 1
integer :: iloq = 1
contains
private
generic,public :: initialize => initialize_4d_auto_knots,initialize_4d_specify_knots
procedure :: initialize_4d_auto_knots
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
type,extends(bspline_class),public :: bspline_5d
!! Class for 5d b-spline interpolation.
private
integer :: nx = 0
integer :: ny = 0
integer :: nz = 0
integer :: nq = 0
integer :: nr = 0
integer :: kx = 0
integer :: ky = 0
integer :: kz = 0
integer :: kq = 0
integer :: kr = 0
real(wp),dimension(:,:,:,:,:),allocatable :: bcoef
real(wp),dimension(:),allocatable :: tx
real(wp),dimension(:),allocatable :: ty
real(wp),dimension(:),allocatable :: tz
real(wp),dimension(:),allocatable :: tq
real(wp),dimension(:),allocatable :: tr
integer :: inbvy = 1
integer :: inbvz = 1
integer :: inbvq = 1
integer :: inbvr = 1
integer :: iloy = 1
integer :: iloz = 1
integer :: iloq = 1
integer :: ilor = 1
contains
private
generic,public :: initialize => initialize_5d_auto_knots,initialize_5d_specify_knots
procedure :: initialize_5d_auto_knots
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
type,extends(bspline_class),public :: bspline_6d
!! Class for 6d b-spline interpolation.
private
integer :: nx = 0
integer :: ny = 0
integer :: nz = 0
integer :: nq = 0
integer :: nr = 0
integer :: ns = 0
integer :: kx = 0
integer :: ky = 0
integer :: kz = 0
integer :: kq = 0
integer :: kr = 0
integer :: ks = 0
real(wp),dimension(:,:,:,:,:,:),allocatable :: bcoef
real(wp),dimension(:),allocatable :: tx
real(wp),dimension(:),allocatable :: ty
real(wp),dimension(:),allocatable :: tz
real(wp),dimension(:),allocatable :: tq
real(wp),dimension(:),allocatable :: tr
real(wp),dimension(:),allocatable :: ts
integer :: inbvy = 1
integer :: inbvz = 1
integer :: inbvq = 1
integer :: inbvr = 1
integer :: inbvs = 1
integer :: iloy = 1
integer :: iloz = 1
integer :: iloq = 1
integer :: ilor = 1
integer :: ilos = 1
contains
private
generic,public :: initialize => initialize_6d_auto_knots,initialize_6d_specify_knots
procedure :: initialize_6d_auto_knots
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
interface bspline_1d
procedure :: bspline_1d_constructor_empty,&
bspline_1d_constructor_auto_knots,&
bspline_1d_constructor_specify_knots
end interface
interface bspline_2d
procedure :: bspline_2d_constructor_empty,&
bspline_2d_constructor_auto_knots,&
bspline_2d_constructor_specify_knots
end interface
interface bspline_3d
procedure :: bspline_3d_constructor_empty,&
bspline_3d_constructor_auto_knots,&
bspline_3d_constructor_specify_knots
end interface
interface bspline_4d
procedure :: bspline_4d_constructor_empty,&
bspline_4d_constructor_auto_knots,&
bspline_4d_constructor_specify_knots
end interface
interface bspline_5d
procedure :: bspline_5d_constructor_empty,&
bspline_5d_constructor_auto_knots,&
bspline_5d_constructor_specify_knots
end interface
interface bspline_6d
procedure :: bspline_6d_constructor_empty,&
bspline_6d_constructor_auto_knots,&
bspline_6d_constructor_specify_knots
end interface
contains
!*****************************************************************************************
!*****************************************************************************************
!>
! This routines returns true if the `iflag` code from the last
! routine called was `=0`. Maybe of the routines have output `iflag`
! variables, so they can be checked explicitly, or this routine
! can be used.
!
! If the class is initialized using a function constructor, then
! this is the only way to know if it was properly initialized,
! since those are pure functions with not output `iflag` arguments.
!
! If `status_ok=.false.`, then the error message can be
! obtained from the [[get_bspline_status_message]] routine.
!
! Note: after an error condition, the [[clear_bspline_flag]] routine
! can be called to reset the `iflag` to 0.
elemental function status_ok(me) result(ok)
implicit none
class(bspline_class),intent(in) :: me
logical :: ok
ok = ( me%iflag == 0 )
end function status_ok
!*****************************************************************************************
!*****************************************************************************************
!>
! This sets the `iflag` variable in the class to `0`
! (which indicates that everything is OK). It can be used
! after an error is encountered.
elemental subroutine clear_bspline_flag(me)
implicit none
class(bspline_class),intent(inout) :: me
me%iflag = 0
end subroutine clear_bspline_flag
!*****************************************************************************************
!*****************************************************************************************
!>
! Get the status message from a [[bspline_class]] routine call.
!
! If `iflag` is not included, then the one in the class is used (which
! corresponds to the last routine called.)
! Otherwise, it will convert the
! input `iflag` argument into the appropriate message.
!
! This is a wrapper for [[get_status_message]].
pure function get_bspline_status_message(me,iflag) result(msg)
implicit none
class(bspline_class),intent(in) :: me
character(len=:),allocatable :: msg !! status message associated with the flag
integer,intent(in),optional :: iflag !! the corresponding status code
if (present(iflag)) then
msg = get_status_message(iflag)
else
msg = get_status_message(me%iflag)
end if
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.
! (this routine is called by the extended classes).
pure subroutine destroy_base(me)
implicit none
class(bspline_class),intent(inout) :: me
me%inbvx = 1
me%iflag = 1
me%initialized = .false.
end subroutine destroy_base
!*****************************************************************************************
!*****************************************************************************************
!>
! Destructor for [[bspline_1d]] class.
pure subroutine destroy_1d(me)
implicit none
class(bspline_1d),intent(inout) :: me
call me%destroy_base()
me%nx = 0
me%kx = 0
if (allocated(me%bcoef)) deallocate(me%bcoef)
if (allocated(me%tx)) deallocate(me%tx)
end subroutine destroy_1d
!*****************************************************************************************
!*****************************************************************************************
!>
! Destructor for [[bspline_2d]] class.
pure subroutine destroy_2d(me)
implicit none
class(bspline_2d),intent(inout) :: me
call me%destroy_base()
me%nx = 0
me%ny = 0
me%kx = 0
me%ky = 0
me%inbvy = 1
me%iloy = 1
if (allocated(me%bcoef)) deallocate(me%bcoef)
if (allocated(me%tx)) deallocate(me%tx)
if (allocated(me%ty)) deallocate(me%ty)
end subroutine destroy_2d
!*****************************************************************************************
!*****************************************************************************************
!>
! Destructor for [[bspline_3d]] class.
pure subroutine destroy_3d(me)
implicit none
class(bspline_3d),intent(inout) :: me
call me%destroy_base()
me%nx = 0
me%ny = 0
me%nz = 0
me%kx = 0
me%ky = 0
me%kz = 0
me%inbvy = 1
me%inbvz = 1
me%iloy = 1
me%iloz = 1
if (allocated(me%bcoef)) deallocate(me%bcoef)
if (allocated(me%tx)) deallocate(me%tx)
if (allocated(me%ty)) deallocate(me%ty)
if (allocated(me%tz)) deallocate(me%tz)
end subroutine destroy_3d
!*****************************************************************************************
!*****************************************************************************************
!>
! Destructor for [[bspline_4d]] class.
pure subroutine destroy_4d(me)
implicit none
class(bspline_4d),intent(inout) :: me
me%nx = 0
me%ny = 0
me%nz = 0
me%nq = 0
me%kx = 0
me%ky = 0
me%kz = 0
me%kq = 0
me%inbvy = 1
me%inbvz = 1
me%inbvq = 1
me%iloy = 1
me%iloz = 1
me%iloq = 1
if (allocated(me%bcoef)) deallocate(me%bcoef)
if (allocated(me%tx)) deallocate(me%tx)
if (allocated(me%ty)) deallocate(me%ty)
if (allocated(me%tz)) deallocate(me%tz)
if (allocated(me%tq)) deallocate(me%tq)
end subroutine destroy_4d
!*****************************************************************************************
!*****************************************************************************************
!>
! Destructor for [[bspline_5d]] class.
pure subroutine destroy_5d(me)
implicit none
class(bspline_5d),intent(inout) :: me
me%nx = 0
me%ny = 0
me%nz = 0
me%nq = 0
me%nr = 0
me%kx = 0
me%ky = 0
me%kz = 0
me%kq = 0
me%kr = 0
me%inbvy = 1
me%inbvz = 1
me%inbvq = 1
me%inbvr = 1
me%iloy = 1
me%iloz = 1
me%iloq = 1
me%ilor = 1
if (allocated(me%bcoef)) deallocate(me%bcoef)
if (allocated(me%tx)) deallocate(me%tx)
if (allocated(me%ty)) deallocate(me%ty)
if (allocated(me%tz)) deallocate(me%tz)
if (allocated(me%tq)) deallocate(me%tq)
if (allocated(me%tr)) deallocate(me%tr)
end subroutine destroy_5d
!*****************************************************************************************
!*****************************************************************************************
!>
! Destructor for [[bspline_6d]] class.
pure subroutine destroy_6d(me)
implicit none
class(bspline_6d),intent(inout) :: me
me%nx = 0
me%ny = 0
me%nz = 0
me%nq = 0
me%nr = 0
me%ns = 0
me%kx = 0
me%ky = 0
me%kz = 0
me%kq = 0
me%kr = 0
me%ks = 0
me%inbvy = 1
me%inbvz = 1
me%inbvq = 1
me%inbvr = 1
me%inbvs = 1
me%iloy = 1
me%iloz = 1
me%iloq = 1
me%ilor = 1
me%ilos = 1
if (allocated(me%bcoef)) deallocate(me%bcoef)
if (allocated(me%tx)) deallocate(me%tx)
if (allocated(me%ty)) deallocate(me%ty)
if (allocated(me%tz)) deallocate(me%tz)
if (allocated(me%tq)) deallocate(me%tq)
if (allocated(me%tr)) deallocate(me%tr)
if (allocated(me%ts)) deallocate(me%ts)
end subroutine destroy_6d
!*****************************************************************************************
!*****************************************************************************************
!>
! Finalizer for [[bspline_1d]] class. Just a wrapper for [[destroy_1d]].
pure elemental subroutine finalize_1d(me)
type(bspline_1d),intent(inout) :: me; call me%destroy()
end subroutine finalize_1d
!*****************************************************************************************
!*****************************************************************************************
!>
! Finalizer for [[bspline_2d]] class. Just a wrapper for [[destroy_2d]].
pure elemental subroutine finalize_2d(me)
type(bspline_2d),intent(inout) :: me; call me%destroy()
end subroutine finalize_2d
!*****************************************************************************************
!*****************************************************************************************
!>
! Finalizer for [[bspline_3d]] class. Just a wrapper for [[destroy_3d]].
pure elemental subroutine finalize_3d(me)
type(bspline_3d),intent(inout) :: me; call me%destroy()
end subroutine finalize_3d
!*****************************************************************************************
!*****************************************************************************************
!>
! Finalizer for [[bspline_4d]] class. Just a wrapper for [[destroy_4d]].
pure elemental subroutine finalize_4d(me)
type(bspline_4d),intent(inout) :: me; call me%destroy()
end subroutine finalize_4d
!*****************************************************************************************
!*****************************************************************************************
!>
! Finalizer for [[bspline_5d]] class. Just a wrapper for [[destroy_5d]].
pure elemental subroutine finalize_5d(me)
type(bspline_5d),intent(inout) :: me; call me%destroy()
end subroutine finalize_5d
!*****************************************************************************************
!*****************************************************************************************
!>
! Finalizer for [[bspline_6d]] class. Just a wrapper for [[destroy_6d]].
pure elemental subroutine finalize_6d(me)
type(bspline_6d),intent(inout) :: me; call me%destroy()
end subroutine finalize_6d
!*****************************************************************************************
!*****************************************************************************************
!>
! It returns an empty [[bspline_1d]] type. Note that INITIALIZE still
! needs to be called before it can be used.
! Not really that useful except perhaps in some OpenMP applications.
pure elemental function bspline_1d_constructor_empty() result(me)
implicit none
type(bspline_1d) :: me
end function bspline_1d_constructor_empty
!*****************************************************************************************
!*****************************************************************************************
!>
! Constructor for a [[bspline_1d]] type (auto knots).
! This is a wrapper for [[initialize_1d_auto_knots]].
pure function bspline_1d_constructor_auto_knots(x,fcn,kx) result(me)
implicit none
type(bspline_1d) :: me
real(wp),dimension(:),intent(in) :: x
real(wp),dimension(:),intent(in) :: fcn
integer,intent(in) :: kx
call initialize_1d_auto_knots(me,x,fcn,kx,me%iflag)
end function bspline_1d_constructor_auto_knots
!*****************************************************************************************
!*****************************************************************************************
!>
! Constructor for a [[bspline_1d]] type (user-specified knots).
! This is a wrapper for [[initialize_1d_specify_knots]].
pure function bspline_1d_constructor_specify_knots(x,fcn,kx,tx) result(me)
implicit none
type(bspline_1d) :: me
real(wp),dimension(:),intent(in) :: x
real(wp),dimension(:),intent(in) :: fcn
integer,intent(in) :: kx
real(wp),dimension(:),intent(in) :: tx
call initialize_1d_specify_knots(me,x,fcn,kx,tx,me%iflag)
end function bspline_1d_constructor_specify_knots
!*****************************************************************************************
!*****************************************************************************************
!>
! Initialize a [[bspline_1d]] type (with automatically-computed knots).
! This is a wrapper for [[db1ink]].
pure subroutine initialize_1d_auto_knots(me,x,fcn,kx,iflag)
implicit none
class(bspline_1d),intent(inout) :: me
real(wp),dimension(:),intent(in) :: x
real(wp),dimension(:),intent(in) :: fcn
integer,intent(in) :: kx
integer,intent(out) :: iflag
integer :: iknot
integer :: nx
call me%destroy()
nx = size(x)
me%nx = nx
me%kx = kx
allocate(me%tx(nx+kx))
allocate(me%bcoef(nx))
iknot = 0 !knot sequence chosen by db1ink
call db1ink(x,nx,fcn,kx,iknot,me%tx,me%bcoef,iflag)
me%initialized = iflag==0
me%iflag = iflag
end subroutine initialize_1d_auto_knots
!*****************************************************************************************
!*****************************************************************************************
!>
! Initialize a [[bspline_1d]] type (with user-specified knots).
! This is a wrapper for [[db1ink]].
pure subroutine initialize_1d_specify_knots(me,x,fcn,kx,tx,iflag)
implicit none
class(bspline_1d),intent(inout) :: me
real(wp),dimension(:),intent(in) :: x
real(wp),dimension(:),intent(in) :: fcn
integer,intent(in) :: kx
real(wp),dimension(:),intent(in) :: tx
integer,intent(out) :: iflag
integer :: nx
call me%destroy()
nx = size(x)
call check_knot_vectors_sizes('initialize_1d_specify_knots',nx=nx,kx=kx,tx=tx,iflag=iflag)
if (iflag == 0) then
me%nx = nx
me%kx = kx
allocate(me%tx(nx+kx))
allocate(me%bcoef(nx))
me%tx = tx
call db1ink(x,nx,fcn,kx,1,me%tx,me%bcoef,iflag)
end if
me%initialized = iflag==0
me%iflag = iflag
end subroutine initialize_1d_specify_knots
!*****************************************************************************************
!*****************************************************************************************
!>
! Evaluate a [[bspline_1d]] interpolate. This is a wrapper for [[db1val]].
pure subroutine evaluate_1d(me,xval,idx,f,iflag)
implicit none
class(bspline_1d),intent(inout) :: me
real(wp),intent(in) :: xval
integer,intent(in) :: idx
real(wp),intent(out) :: f
integer,intent(out) :: iflag
if (me%initialized) then
call db1val(xval,idx,me%tx,me%nx,me%kx,me%bcoef,f,iflag,me%inbvx)
else
iflag = 1
end if
me%iflag = iflag
end subroutine evaluate_1d
!*****************************************************************************************
!*****************************************************************************************
!>
! It returns an empty [[bspline_2d]] type. Note that INITIALIZE still
! needs to be called before it can be used.
! Not really that useful except perhaps in some OpenMP applications.
elemental function bspline_2d_constructor_empty() result(me)
implicit none
type(bspline_2d) :: me
end function bspline_2d_constructor_empty
!*****************************************************************************************
!*****************************************************************************************
!>
! Constructor for a [[bspline_2d]] type (auto knots).
! This is a wrapper for [[initialize_2d_auto_knots]].
pure function bspline_2d_constructor_auto_knots(x,y,fcn,kx,ky) result(me)
implicit none
type(bspline_2d) :: me
real(wp),dimension(:),intent(in) :: x
real(wp),dimension(:),intent(in) :: y
real(wp),dimension(:,:),intent(in) :: fcn
integer,intent(in) :: kx
integer,intent(in) :: ky
call initialize_2d_auto_knots(me,x,y,fcn,kx,ky,me%iflag)
end function bspline_2d_constructor_auto_knots
!*****************************************************************************************
!*****************************************************************************************
!>
! Constructor for a [[bspline_2d]] type (user-specified knots).
! This is a wrapper for [[initialize_2d_specify_knots]].
pure function bspline_2d_constructor_specify_knots(x,y,fcn,kx,ky,tx,ty) result(me)
implicit none
type(bspline_2d) :: me
real(wp),dimension(:),intent(in) :: x
real(wp),dimension(:),intent(in) :: y
real(wp),dimension(:,:),intent(in) :: fcn
integer,intent(in) :: kx
integer,intent(in) :: ky
real(wp),dimension(:),intent(in) :: tx
real(wp),dimension(:),intent(in) :: ty
call initialize_2d_specify_knots(me,x,y,fcn,kx,ky,tx,ty,me%iflag)
end function bspline_2d_constructor_specify_knots
!*****************************************************************************************
!*****************************************************************************************
!>
! Initialize a [[bspline_2d]] type (with automatically-computed knots).
! This is a wrapper for [[db2ink]].
pure subroutine initialize_2d_auto_knots(me,x,y,fcn,kx,ky,iflag)
implicit none
class(bspline_2d),intent(inout) :: me
real(wp),dimension(:),intent(in) :: x
real(wp),dimension(:),intent(in) :: y
real(wp),dimension(:,:),intent(in) :: fcn
integer,intent(in) :: kx
integer,intent(in) :: ky
integer,intent(out) :: iflag
integer :: iknot
integer :: nx,ny
call me%destroy()
nx = size(x)
ny = size(y)
me%nx = nx
me%ny = ny
me%kx = kx
me%ky = ky
allocate(me%tx(nx+kx))
allocate(me%ty(ny+ky))
allocate(me%bcoef(nx,ny))
iknot = 0 !knot sequence chosen by db2ink