ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_ellipsoid.F90
Go to the documentation of this file.
1!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3!!!! !!!!
4!!!! ParaMonte: Parallel Monte Carlo and Machine Learning Library. !!!!
5!!!! !!!!
6!!!! Copyright (C) 2012-present, The Computational Data Science Lab !!!!
7!!!! !!!!
8!!!! This file is part of the ParaMonte library. !!!!
9!!!! !!!!
10!!!! LICENSE !!!!
11!!!! !!!!
12!!!! https://github.com/cdslaborg/paramonte/blob/main/LICENSE.md !!!!
13!!!! !!!!
14!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
15!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
16
102
103!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
104
106
107 use pm_kind, only: SK, IK, LK
108
109 implicit none
110
111 character(*, SK), parameter :: MODULE_NAME = "@pm_ellipsoid"
112
113!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
114
170
171 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
172
173#if RK5_ENABLED
174 PURE elemental module function getVolUnitBallIter_RK5(ndim) result(volUnitBall)
175#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
176 !DEC$ ATTRIBUTES DLLEXPORT :: getVolUnitBallIter_RK5
177#endif
178 use pm_kind, only: RKG => RK5
179 real(RKG) , intent(in) :: ndim
180 real(RKG) :: volUnitBall
181 end function
182#endif
183
184#if RK4_ENABLED
185 PURE elemental module function getVolUnitBallIter_RK4(ndim) result(volUnitBall)
186#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
187 !DEC$ ATTRIBUTES DLLEXPORT :: getVolUnitBallIter_RK4
188#endif
189 use pm_kind, only: RKG => RK4
190 real(RKG) , intent(in) :: ndim
191 real(RKG) :: volUnitBall
192 end function
193#endif
194
195#if RK3_ENABLED
196 PURE elemental module function getVolUnitBallIter_RK3(ndim) result(volUnitBall)
197#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
198 !DEC$ ATTRIBUTES DLLEXPORT :: getVolUnitBallIter_RK3
199#endif
200 use pm_kind, only: RKG => RK3
201 real(RKG) , intent(in) :: ndim
202 real(RKG) :: volUnitBall
203 end function
204#endif
205
206#if RK2_ENABLED
207 PURE elemental module function getVolUnitBallIter_RK2(ndim) result(volUnitBall)
208#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
209 !DEC$ ATTRIBUTES DLLEXPORT :: getVolUnitBallIter_RK2
210#endif
211 use pm_kind, only: RKG => RK2
212 real(RKG) , intent(in) :: ndim
213 real(RKG) :: volUnitBall
214 end function
215#endif
216
217#if RK1_ENABLED
218 PURE elemental module function getVolUnitBallIter_RK1(ndim) result(volUnitBall)
219#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
220 !DEC$ ATTRIBUTES DLLEXPORT :: getVolUnitBallIter_RK1
221#endif
222 use pm_kind, only: RKG => RK1
223 real(RKG) , intent(in) :: ndim
224 real(RKG) :: volUnitBall
225 end function
226#endif
227
228 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
229
230 end interface
231
292
293 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
294
295#if RK5_ENABLED
296 PURE elemental module subroutine setVolUnitBallIter_RK5(volUnitBall, ndim)
297#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
298 !DEC$ ATTRIBUTES DLLEXPORT :: setVolUnitBallIter_RK5
299#endif
300 use pm_kind, only: RKG => RK5
301 real(RKG) , intent(out) :: volUnitBall
302 integer(IK) , intent(in) :: ndim
303 end subroutine
304#endif
305
306#if RK4_ENABLED
307 PURE elemental module subroutine setVolUnitBallIter_RK4(volUnitBall, ndim)
308#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
309 !DEC$ ATTRIBUTES DLLEXPORT :: setVolUnitBallIter_RK4
310#endif
311 use pm_kind, only: RKG => RK4
312 real(RKG) , intent(out) :: volUnitBall
313 integer(IK) , intent(in) :: ndim
314 end subroutine
315#endif
316
317#if RK3_ENABLED
318 PURE elemental module subroutine setVolUnitBallIter_RK3(volUnitBall, ndim)
319#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
320 !DEC$ ATTRIBUTES DLLEXPORT :: setVolUnitBallIter_RK3
321#endif
322 use pm_kind, only: RKG => RK3
323 real(RKG) , intent(out) :: volUnitBall
324 integer(IK) , intent(in) :: ndim
325 end subroutine
326#endif
327
328#if RK2_ENABLED
329 PURE elemental module subroutine setVolUnitBallIter_RK2(volUnitBall, ndim)
330#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
331 !DEC$ ATTRIBUTES DLLEXPORT :: setVolUnitBallIter_RK2
332#endif
333 use pm_kind, only: RKG => RK2
334 real(RKG) , intent(out) :: volUnitBall
335 integer(IK) , intent(in) :: ndim
336 end subroutine
337#endif
338
339#if RK1_ENABLED
340 PURE elemental module subroutine setVolUnitBallIter_RK1(volUnitBall, ndim)
341#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
342 !DEC$ ATTRIBUTES DLLEXPORT :: setVolUnitBallIter_RK1
343#endif
344 use pm_kind, only: RKG => RK1
345 real(RKG) , intent(out) :: volUnitBall
346 integer(IK) , intent(in) :: ndim
347 end subroutine
348#endif
349
350 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
351
352 end interface
353
354!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
355
411
412 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
413
414#if RK5_ENABLED
415 PURE elemental module function getLogVolUnitBallGamm_RK5(ndim) result(logVolUnitBall)
416#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
417 !DEC$ ATTRIBUTES DLLEXPORT :: getLogVolUnitBallGamm_RK5
418#endif
419 use pm_kind, only: RKG => RK5
420 real(RKG) , intent(in) :: ndim
421 real(RKG) :: logVolUnitBall
422 end function
423#endif
424
425#if RK4_ENABLED
426 PURE elemental module function getLogVolUnitBallGamm_RK4(ndim) result(logVolUnitBall)
427#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
428 !DEC$ ATTRIBUTES DLLEXPORT :: getLogVolUnitBallGamm_RK4
429#endif
430 use pm_kind, only: RKG => RK4
431 real(RKG) , intent(in) :: ndim
432 real(RKG) :: logVolUnitBall
433 end function
434#endif
435
436#if RK3_ENABLED
437 PURE elemental module function getLogVolUnitBallGamm_RK3(ndim) result(logVolUnitBall)
438#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
439 !DEC$ ATTRIBUTES DLLEXPORT :: getLogVolUnitBallGamm_RK3
440#endif
441 use pm_kind, only: RKG => RK3
442 real(RKG) , intent(in) :: ndim
443 real(RKG) :: logVolUnitBall
444 end function
445#endif
446
447#if RK2_ENABLED
448 PURE elemental module function getLogVolUnitBallGamm_RK2(ndim) result(logVolUnitBall)
449#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
450 !DEC$ ATTRIBUTES DLLEXPORT :: getLogVolUnitBallGamm_RK2
451#endif
452 use pm_kind, only: RKG => RK2
453 real(RKG) , intent(in) :: ndim
454 real(RKG) :: logVolUnitBall
455 end function
456#endif
457
458#if RK1_ENABLED
459 PURE elemental module function getLogVolUnitBallGamm_RK1(ndim) result(logVolUnitBall)
460#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
461 !DEC$ ATTRIBUTES DLLEXPORT :: getLogVolUnitBallGamm_RK1
462#endif
463 use pm_kind, only: RKG => RK1
464 real(RKG) , intent(in) :: ndim
465 real(RKG) :: logVolUnitBall
466 end function
467#endif
468
469 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
470
471 end interface
472
473!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
474
559
560 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
561
562#if RK5_ENABLED
563 PURE elemental module subroutine setLogVolUnitBallGamm_RK5(logVolUnitBall, ndim)
564#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
565 !DEC$ ATTRIBUTES DLLEXPORT :: setLogVolUnitBallGamm_RK5
566#endif
567 use pm_kind, only: RKG => RK5
568 real(RKG) , intent(out) :: logVolUnitBall
569 real(RKG) , intent(in) :: ndim
570 end subroutine
571#endif
572
573#if RK4_ENABLED
574 PURE elemental module subroutine setLogVolUnitBallGamm_RK4(logVolUnitBall, ndim)
575#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
576 !DEC$ ATTRIBUTES DLLEXPORT :: setLogVolUnitBallGamm_RK4
577#endif
578 use pm_kind, only: RKG => RK4
579 real(RKG) , intent(out) :: logVolUnitBall
580 real(RKG) , intent(in) :: ndim
581 end subroutine
582#endif
583
584#if RK3_ENABLED
585 PURE elemental module subroutine setLogVolUnitBallGamm_RK3(logVolUnitBall, ndim)
586#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
587 !DEC$ ATTRIBUTES DLLEXPORT :: setLogVolUnitBallGamm_RK3
588#endif
589 use pm_kind, only: RKG => RK3
590 real(RKG) , intent(out) :: logVolUnitBall
591 real(RKG) , intent(in) :: ndim
592 end subroutine
593#endif
594
595#if RK2_ENABLED
596 PURE elemental module subroutine setLogVolUnitBallGamm_RK2(logVolUnitBall, ndim)
597#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
598 !DEC$ ATTRIBUTES DLLEXPORT :: setLogVolUnitBallGamm_RK2
599#endif
600 use pm_kind, only: RKG => RK2
601 real(RKG) , intent(out) :: logVolUnitBall
602 real(RKG) , intent(in) :: ndim
603 end subroutine
604#endif
605
606#if RK1_ENABLED
607 PURE elemental module subroutine setLogVolUnitBallGamm_RK1(logVolUnitBall, ndim)
608#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
609 !DEC$ ATTRIBUTES DLLEXPORT :: setLogVolUnitBallGamm_RK1
610#endif
611 use pm_kind, only: RKG => RK1
612 real(RKG) , intent(out) :: logVolUnitBall
613 real(RKG) , intent(in) :: ndim
614 end subroutine
615#endif
616
617 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
618
619#if RK5_ENABLED
620 PURE elemental module subroutine setLogVolUnitBallIter_RK5(logVolUnitBall, ndim)
621#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
622 !DEC$ ATTRIBUTES DLLEXPORT :: setLogVolUnitBallIter_RK5
623#endif
624 use pm_kind, only: RKG => RK5
625 real(RKG) , intent(out) :: logVolUnitBall
626 integer(IK) , intent(in) :: ndim
627 end subroutine
628#endif
629
630#if RK4_ENABLED
631 PURE elemental module subroutine setLogVolUnitBallIter_RK4(logVolUnitBall, ndim)
632#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
633 !DEC$ ATTRIBUTES DLLEXPORT :: setLogVolUnitBallIter_RK4
634#endif
635 use pm_kind, only: RKG => RK4
636 real(RKG) , intent(out) :: logVolUnitBall
637 integer(IK) , intent(in) :: ndim
638 end subroutine
639#endif
640
641#if RK3_ENABLED
642 PURE elemental module subroutine setLogVolUnitBallIter_RK3(logVolUnitBall, ndim)
643#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
644 !DEC$ ATTRIBUTES DLLEXPORT :: setLogVolUnitBallIter_RK3
645#endif
646 use pm_kind, only: RKG => RK3
647 real(RKG) , intent(out) :: logVolUnitBall
648 integer(IK) , intent(in) :: ndim
649 end subroutine
650#endif
651
652#if RK2_ENABLED
653 PURE elemental module subroutine setLogVolUnitBallIter_RK2(logVolUnitBall, ndim)
654#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
655 !DEC$ ATTRIBUTES DLLEXPORT :: setLogVolUnitBallIter_RK2
656#endif
657 use pm_kind, only: RKG => RK2
658 real(RKG) , intent(out) :: logVolUnitBall
659 integer(IK) , intent(in) :: ndim
660 end subroutine
661#endif
662
663#if RK1_ENABLED
664 PURE elemental module subroutine setLogVolUnitBallIter_RK1(logVolUnitBall, ndim)
665#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
666 !DEC$ ATTRIBUTES DLLEXPORT :: setLogVolUnitBallIter_RK1
667#endif
668 use pm_kind, only: RKG => RK1
669 real(RKG) , intent(out) :: logVolUnitBall
670 integer(IK) , intent(in) :: ndim
671 end subroutine
672#endif
673
674 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
675
676 end interface
677
678!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
679
746 interface getLogVolEll
747
748 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
749
750#if RK5_ENABLED
751 PURE module function getLogVolEll_RK5(gramian) result(logVolEll)
752#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
753 !DEC$ ATTRIBUTES DLLEXPORT :: getLogVolEll_RK5
754#endif
755 use pm_kind, only: RKG => RK5
756 real(RKG) , intent(in) , contiguous :: gramian(:,:)
757 real(RKG) :: logVolEll
758 end function
759#endif
760
761#if RK4_ENABLED
762 PURE module function getLogVolEll_RK4(gramian) result(logVolEll)
763#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
764 !DEC$ ATTRIBUTES DLLEXPORT :: getLogVolEll_RK4
765#endif
766 use pm_kind, only: RKG => RK4
767 real(RKG) , intent(in) , contiguous :: gramian(:,:)
768 real(RKG) :: logVolEll
769 end function
770#endif
771
772#if RK3_ENABLED
773 PURE module function getLogVolEll_RK3(gramian) result(logVolEll)
774#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
775 !DEC$ ATTRIBUTES DLLEXPORT :: getLogVolEll_RK3
776#endif
777 use pm_kind, only: RKG => RK3
778 real(RKG) , intent(in) , contiguous :: gramian(:,:)
779 real(RKG) :: logVolEll
780 end function
781#endif
782
783#if RK2_ENABLED
784 PURE module function getLogVolEll_RK2(gramian) result(logVolEll)
785#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
786 !DEC$ ATTRIBUTES DLLEXPORT :: getLogVolEll_RK2
787#endif
788 use pm_kind, only: RKG => RK2
789 real(RKG) , intent(in) , contiguous :: gramian(:,:)
790 real(RKG) :: logVolEll
791 end function
792#endif
793
794#if RK1_ENABLED
795 PURE module function getLogVolEll_RK1(gramian) result(logVolEll)
796#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
797 !DEC$ ATTRIBUTES DLLEXPORT :: getLogVolEll_RK1
798#endif
799 use pm_kind, only: RKG => RK1
800 real(RKG) , intent(in) , contiguous :: gramian(:,:)
801 real(RKG) :: logVolEll
802 end function
803#endif
804
805 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
806
807 end interface
808
809!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
810
889 interface isMemberEll
890
891 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
892 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
893 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
894
895 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
896
897#if RK5_ENABLED
898 PURE module function isMemberSphOrg_RK5(radius, point) result(isMember)
899#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
900 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberSphOrg_RK5
901#endif
902 use pm_kind, only: RKG => RK5
903 real(RKG) , intent(in) , contiguous :: point(:)
904 real(RKG) , intent(in) :: radius
905 logical(LK) :: isMember
906 end function
907#endif
908
909#if RK4_ENABLED
910 PURE module function isMemberSphOrg_RK4(radius, point) result(isMember)
911#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
912 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberSphOrg_RK4
913#endif
914 use pm_kind, only: RKG => RK4
915 real(RKG) , intent(in) , contiguous :: point(:)
916 real(RKG) , intent(in) :: radius
917 logical(LK) :: isMember
918 end function
919#endif
920
921#if RK3_ENABLED
922 PURE module function isMemberSphOrg_RK3(radius, point) result(isMember)
923#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
924 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberSphOrg_RK3
925#endif
926 use pm_kind, only: RKG => RK3
927 real(RKG) , intent(in) , contiguous :: point(:)
928 real(RKG) , intent(in) :: radius
929 logical(LK) :: isMember
930 end function
931#endif
932
933#if RK2_ENABLED
934 PURE module function isMemberSphOrg_RK2(radius, point) result(isMember)
935#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
936 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberSphOrg_RK2
937#endif
938 use pm_kind, only: RKG => RK2
939 real(RKG) , intent(in) , contiguous :: point(:)
940 real(RKG) , intent(in) :: radius
941 logical(LK) :: isMember
942 end function
943#endif
944
945#if RK1_ENABLED
946 PURE module function isMemberSphOrg_RK1(radius, point) result(isMember)
947#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
948 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberSphOrg_RK1
949#endif
950 use pm_kind, only: RKG => RK1
951 real(RKG) , intent(in) , contiguous :: point(:)
952 real(RKG) , intent(in) :: radius
953 logical(LK) :: isMember
954 end function
955#endif
956
957 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
958
959#if RK5_ENABLED
960 PURE module function isMemberSphCen_RK5(radius, center, point) result(isMember)
961#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
962 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberSphCen_RK5
963#endif
964 use pm_kind, only: RKG => RK5
965 real(RKG) , intent(in) , contiguous :: point(:)
966 real(RKG) , intent(in) , contiguous :: center(:)
967 real(RKG) , intent(in) :: radius
968 logical(LK) :: isMember
969 end function
970#endif
971
972#if RK4_ENABLED
973 PURE module function isMemberSphCen_RK4(radius, center, point) result(isMember)
974#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
975 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberSphCen_RK4
976#endif
977 use pm_kind, only: RKG => RK4
978 real(RKG) , intent(in) , contiguous :: point(:)
979 real(RKG) , intent(in) , contiguous :: center(:)
980 real(RKG) , intent(in) :: radius
981 logical(LK) :: isMember
982 end function
983#endif
984
985#if RK3_ENABLED
986 PURE module function isMemberSphCen_RK3(radius, center, point) result(isMember)
987#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
988 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberSphCen_RK3
989#endif
990 use pm_kind, only: RKG => RK3
991 real(RKG) , intent(in) , contiguous :: point(:)
992 real(RKG) , intent(in) , contiguous :: center(:)
993 real(RKG) , intent(in) :: radius
994 logical(LK) :: isMember
995 end function
996#endif
997
998#if RK2_ENABLED
999 PURE module function isMemberSphCen_RK2(radius, center, point) result(isMember)
1000#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1001 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberSphCen_RK2
1002#endif
1003 use pm_kind, only: RKG => RK2
1004 real(RKG) , intent(in) , contiguous :: point(:)
1005 real(RKG) , intent(in) , contiguous :: center(:)
1006 real(RKG) , intent(in) :: radius
1007 logical(LK) :: isMember
1008 end function
1009#endif
1010
1011#if RK1_ENABLED
1012 PURE module function isMemberSphCen_RK1(radius, center, point) result(isMember)
1013#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1014 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberSphCen_RK1
1015#endif
1016 use pm_kind, only: RKG => RK1
1017 real(RKG) , intent(in) , contiguous :: point(:)
1018 real(RKG) , intent(in) , contiguous :: center(:)
1019 real(RKG) , intent(in) :: radius
1020 logical(LK) :: isMember
1021 end function
1022#endif
1023
1024 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1025
1026 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1027 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1028 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1029
1030 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1031
1032#if RK5_ENABLED
1033 PURE module function isMemberEllOrg_RK5(invGram, point) result(isMember)
1034#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1035 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberEllOrg_RK5
1036#endif
1037 use pm_kind, only: RKG => RK5
1038 real(RKG) , intent(in) , contiguous :: point(:)
1039 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1040 logical(LK) :: isMember
1041 end function
1042#endif
1043
1044#if RK4_ENABLED
1045 PURE module function isMemberEllOrg_RK4(invGram, point) result(isMember)
1046#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1047 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberEllOrg_RK4
1048#endif
1049 use pm_kind, only: RKG => RK4
1050 real(RKG) , intent(in) , contiguous :: point(:)
1051 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1052 logical(LK) :: isMember
1053 end function
1054#endif
1055
1056#if RK3_ENABLED
1057 PURE module function isMemberEllOrg_RK3(invGram, point) result(isMember)
1058#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1059 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberEllOrg_RK3
1060#endif
1061 use pm_kind, only: RKG => RK3
1062 real(RKG) , intent(in) , contiguous :: point(:)
1063 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1064 logical(LK) :: isMember
1065 end function
1066#endif
1067
1068#if RK2_ENABLED
1069 PURE module function isMemberEllOrg_RK2(invGram, point) result(isMember)
1070#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1071 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberEllOrg_RK2
1072#endif
1073 use pm_kind, only: RKG => RK2
1074 real(RKG) , intent(in) , contiguous :: point(:)
1075 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1076 logical(LK) :: isMember
1077 end function
1078#endif
1079
1080#if RK1_ENABLED
1081 PURE module function isMemberEllOrg_RK1(invGram, point) result(isMember)
1082#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1083 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberEllOrg_RK1
1084#endif
1085 use pm_kind, only: RKG => RK1
1086 real(RKG) , intent(in) , contiguous :: point(:)
1087 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1088 logical(LK) :: isMember
1089 end function
1090#endif
1091
1092 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1093
1094#if RK5_ENABLED
1095 PURE module function isMemberEllCen_RK5(invGram, center, point) result(isMember)
1096#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1097 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberEllCen_RK5
1098#endif
1099 use pm_kind, only: RKG => RK5
1100 real(RKG) , intent(in) , contiguous :: point(:)
1101 real(RKG) , intent(in) , contiguous :: center(:)
1102 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1103 logical(LK) :: isMember
1104 end function
1105#endif
1106
1107#if RK4_ENABLED
1108 PURE module function isMemberEllCen_RK4(invGram, center, point) result(isMember)
1109#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1110 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberEllCen_RK4
1111#endif
1112 use pm_kind, only: RKG => RK4
1113 real(RKG) , intent(in) , contiguous :: point(:)
1114 real(RKG) , intent(in) , contiguous :: center(:)
1115 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1116 logical(LK) :: isMember
1117 end function
1118#endif
1119
1120#if RK3_ENABLED
1121 PURE module function isMemberEllCen_RK3(invGram, center, point) result(isMember)
1122#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1123 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberEllCen_RK3
1124#endif
1125 use pm_kind, only: RKG => RK3
1126 real(RKG) , intent(in) , contiguous :: point(:)
1127 real(RKG) , intent(in) , contiguous :: center(:)
1128 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1129 logical(LK) :: isMember
1130 end function
1131#endif
1132
1133#if RK2_ENABLED
1134 PURE module function isMemberEllCen_RK2(invGram, center, point) result(isMember)
1135#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1136 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberEllCen_RK2
1137#endif
1138 use pm_kind, only: RKG => RK2
1139 real(RKG) , intent(in) , contiguous :: point(:)
1140 real(RKG) , intent(in) , contiguous :: center(:)
1141 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1142 logical(LK) :: isMember
1143 end function
1144#endif
1145
1146#if RK1_ENABLED
1147 PURE module function isMemberEllCen_RK1(invGram, center, point) result(isMember)
1148#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1149 !DEC$ ATTRIBUTES DLLEXPORT :: isMemberEllCen_RK1
1150#endif
1151 use pm_kind, only: RKG => RK1
1152 real(RKG) , intent(in) , contiguous :: point(:)
1153 real(RKG) , intent(in) , contiguous :: center(:)
1154 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1155 logical(LK) :: isMember
1156 end function
1157#endif
1158
1159 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1160
1161 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1162 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1163 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1164
1165 end interface
1166
1167!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1168
1237
1238 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1239 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1240 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1241
1242 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1243
1244#if RK5_ENABLED
1245 PURE module function getCountMemberSphOrg_RK5(radius, point) result(countMemberEll)
1246#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1247 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberSphOrg_RK5
1248#endif
1249 use pm_kind, only: RKG => RK5
1250 real(RKG) , intent(in) , contiguous :: point(:,:)
1251 real(RKG) , intent(in) :: radius
1252 integer(IK) :: countMemberEll
1253 end function
1254#endif
1255
1256#if RK4_ENABLED
1257 PURE module function getCountMemberSphOrg_RK4(radius, point) result(countMemberEll)
1258#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1259 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberSphOrg_RK4
1260#endif
1261 use pm_kind, only: RKG => RK4
1262 real(RKG) , intent(in) , contiguous :: point(:,:)
1263 real(RKG) , intent(in) :: radius
1264 integer(IK) :: countMemberEll
1265 end function
1266#endif
1267
1268#if RK3_ENABLED
1269 PURE module function getCountMemberSphOrg_RK3(radius, point) result(countMemberEll)
1270#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1271 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberSphOrg_RK3
1272#endif
1273 use pm_kind, only: RKG => RK3
1274 real(RKG) , intent(in) , contiguous :: point(:,:)
1275 real(RKG) , intent(in) :: radius
1276 integer(IK) :: countMemberEll
1277 end function
1278#endif
1279
1280#if RK2_ENABLED
1281 PURE module function getCountMemberSphOrg_RK2(radius, point) result(countMemberEll)
1282#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1283 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberSphOrg_RK2
1284#endif
1285 use pm_kind, only: RKG => RK2
1286 real(RKG) , intent(in) , contiguous :: point(:,:)
1287 real(RKG) , intent(in) :: radius
1288 integer(IK) :: countMemberEll
1289 end function
1290#endif
1291
1292#if RK1_ENABLED
1293 PURE module function getCountMemberSphOrg_RK1(radius, point) result(countMemberEll)
1294#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1295 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberSphOrg_RK1
1296#endif
1297 use pm_kind, only: RKG => RK1
1298 real(RKG) , intent(in) , contiguous :: point(:,:)
1299 real(RKG) , intent(in) :: radius
1300 integer(IK) :: countMemberEll
1301 end function
1302#endif
1303
1304 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1305
1306#if RK5_ENABLED
1307 PURE module function getCountMemberSphCen_RK5(radius, center, point) result(countMemberEll)
1308#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1309 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberSphCen_RK5
1310#endif
1311 use pm_kind, only: RKG => RK5
1312 real(RKG) , intent(in) , contiguous :: point(:,:)
1313 real(RKG) , intent(in) , contiguous :: center(:)
1314 real(RKG) , intent(in) :: radius
1315 integer(IK) :: countMemberEll
1316 end function
1317#endif
1318
1319#if RK4_ENABLED
1320 PURE module function getCountMemberSphCen_RK4(radius, center, point) result(countMemberEll)
1321#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1322 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberSphCen_RK4
1323#endif
1324 use pm_kind, only: RKG => RK4
1325 real(RKG) , intent(in) , contiguous :: point(:,:)
1326 real(RKG) , intent(in) , contiguous :: center(:)
1327 real(RKG) , intent(in) :: radius
1328 integer(IK) :: countMemberEll
1329 end function
1330#endif
1331
1332#if RK3_ENABLED
1333 PURE module function getCountMemberSphCen_RK3(radius, center, point) result(countMemberEll)
1334#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1335 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberSphCen_RK3
1336#endif
1337 use pm_kind, only: RKG => RK3
1338 real(RKG) , intent(in) , contiguous :: point(:,:)
1339 real(RKG) , intent(in) , contiguous :: center(:)
1340 real(RKG) , intent(in) :: radius
1341 integer(IK) :: countMemberEll
1342 end function
1343#endif
1344
1345#if RK2_ENABLED
1346 PURE module function getCountMemberSphCen_RK2(radius, center, point) result(countMemberEll)
1347#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1348 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberSphCen_RK2
1349#endif
1350 use pm_kind, only: RKG => RK2
1351 real(RKG) , intent(in) , contiguous :: point(:,:)
1352 real(RKG) , intent(in) , contiguous :: center(:)
1353 real(RKG) , intent(in) :: radius
1354 integer(IK) :: countMemberEll
1355 end function
1356#endif
1357
1358#if RK1_ENABLED
1359 PURE module function getCountMemberSphCen_RK1(radius, center, point) result(countMemberEll)
1360#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1361 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberSphCen_RK1
1362#endif
1363 use pm_kind, only: RKG => RK1
1364 real(RKG) , intent(in) , contiguous :: point(:,:)
1365 real(RKG) , intent(in) , contiguous :: center(:)
1366 real(RKG) , intent(in) :: radius
1367 integer(IK) :: countMemberEll
1368 end function
1369#endif
1370
1371 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1372
1373 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1374 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1375 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1376
1377 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1378
1379#if RK5_ENABLED
1380 PURE module function getCountMemberEllOrg_RK5(invGram, point) result(countMemberEll)
1381#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1382 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberEllOrg_RK5
1383#endif
1384 use pm_kind, only: RKG => RK5
1385 real(RKG) , intent(in) , contiguous :: point(:,:)
1386 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1387 integer(IK) :: countMemberEll
1388 end function
1389#endif
1390
1391#if RK4_ENABLED
1392 PURE module function getCountMemberEllOrg_RK4(invGram, point) result(countMemberEll)
1393#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1394 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberEllOrg_RK4
1395#endif
1396 use pm_kind, only: RKG => RK4
1397 real(RKG) , intent(in) , contiguous :: point(:,:)
1398 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1399 integer(IK) :: countMemberEll
1400 end function
1401#endif
1402
1403#if RK3_ENABLED
1404 PURE module function getCountMemberEllOrg_RK3(invGram, point) result(countMemberEll)
1405#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1406 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberEllOrg_RK3
1407#endif
1408 use pm_kind, only: RKG => RK3
1409 real(RKG) , intent(in) , contiguous :: point(:,:)
1410 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1411 integer(IK) :: countMemberEll
1412 end function
1413#endif
1414
1415#if RK2_ENABLED
1416 PURE module function getCountMemberEllOrg_RK2(invGram, point) result(countMemberEll)
1417#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1418 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberEllOrg_RK2
1419#endif
1420 use pm_kind, only: RKG => RK2
1421 real(RKG) , intent(in) , contiguous :: point(:,:)
1422 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1423 integer(IK) :: countMemberEll
1424 end function
1425#endif
1426
1427#if RK1_ENABLED
1428 PURE module function getCountMemberEllOrg_RK1(invGram, point) result(countMemberEll)
1429#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1430 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberEllOrg_RK1
1431#endif
1432 use pm_kind, only: RKG => RK1
1433 real(RKG) , intent(in) , contiguous :: point(:,:)
1434 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1435 integer(IK) :: countMemberEll
1436 end function
1437#endif
1438
1439 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1440
1441#if RK5_ENABLED
1442 PURE module function getCountMemberEllCen_RK5(invGram, center, point) result(countMemberEll)
1443#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1444 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberEllCen_RK5
1445#endif
1446 use pm_kind, only: RKG => RK5
1447 real(RKG) , intent(in) , contiguous :: point(:,:)
1448 real(RKG) , intent(in) , contiguous :: center(:)
1449 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1450 integer(IK) :: countMemberEll
1451 end function
1452#endif
1453
1454#if RK4_ENABLED
1455 PURE module function getCountMemberEllCen_RK4(invGram, center, point) result(countMemberEll)
1456#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1457 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberEllCen_RK4
1458#endif
1459 use pm_kind, only: RKG => RK4
1460 real(RKG) , intent(in) , contiguous :: point(:,:)
1461 real(RKG) , intent(in) , contiguous :: center(:)
1462 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1463 integer(IK) :: countMemberEll
1464 end function
1465#endif
1466
1467#if RK3_ENABLED
1468 PURE module function getCountMemberEllCen_RK3(invGram, center, point) result(countMemberEll)
1469#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1470 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberEllCen_RK3
1471#endif
1472 use pm_kind, only: RKG => RK3
1473 real(RKG) , intent(in) , contiguous :: point(:,:)
1474 real(RKG) , intent(in) , contiguous :: center(:)
1475 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1476 integer(IK) :: countMemberEll
1477 end function
1478#endif
1479
1480#if RK2_ENABLED
1481 PURE module function getCountMemberEllCen_RK2(invGram, center, point) result(countMemberEll)
1482#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1483 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberEllCen_RK2
1484#endif
1485 use pm_kind, only: RKG => RK2
1486 real(RKG) , intent(in) , contiguous :: point(:,:)
1487 real(RKG) , intent(in) , contiguous :: center(:)
1488 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1489 integer(IK) :: countMemberEll
1490 end function
1491#endif
1492
1493#if RK1_ENABLED
1494 PURE module function getCountMemberEllCen_RK1(invGram, center, point) result(countMemberEll)
1495#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1496 !DEC$ ATTRIBUTES DLLEXPORT :: getCountMemberEllCen_RK1
1497#endif
1498 use pm_kind, only: RKG => RK1
1499 real(RKG) , intent(in) , contiguous :: point(:,:)
1500 real(RKG) , intent(in) , contiguous :: center(:)
1501 real(RKG) , intent(in) , contiguous :: invGram(:,:)
1502 integer(IK) :: countMemberEll
1503 end function
1504#endif
1505
1506 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1507
1508 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1509 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1510 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1511
1512 end interface
1513
1514!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1515
1516contains
1517
1518!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1519
1520!#if 0
1521! !> \brief
1522! !> Generate Compute the *effective* total volume of all bounded ellipsoidal regions together while considering the possible overlaps.
1523! !> This computation is done aggressively via Monte Carlo simulation.
1524! !>
1525! !> \param[in] nd : The number of dimensions.
1526! !> \param[in] nc : The number of input clusters.
1527! !> \param[in] reltol : The accuracy tolerance of the computed effective sum of volumes (NOT the log-transformed of the sum).
1528! !> \param[in] logVol : An array of size `(nc)` representing the logarithms of the volumes of the ellipsoids normalized to the volume of an nd-ball.
1529! !> \param[in] center : An array of size `(nd,nc)` representing the centers of the ellipsoids.
1530! !> \param[in] choDia : An array of size `(nd,nc)` representing the diagonal elements of the Cholesky lower triangle of the covariance matrices of the ellipsoids.
1531! !> \param[in] logSumVol : `` = sum(logVol)``.
1532! !> \param[in] invCov : An array of size `(nd,nd,nc)` representing the inverse covariance matrices of the ellipsoids.
1533! !> \param[in] choLow : An array of size `(nd,nd,nc)` whose lower-triangle represents the Cholesky lower triangle of the covariance matrices of the ellipsoids.
1534! !> \param[inout] logSumVolEff : The logarithm of the sum of `LogVolNormed`. On output, it will be rewritten with the effective log-sum of the (normalized) volumes.
1535! !> \param[in] ScaleFactorSq : An array of size `(nc)` representing the maximum of Mahalanobis-distance-squared of the members of each cluster from its center.
1536! !> This is needed if the ellipsoid properties have not been scaled, i.e., `ScaleFactorSq > 1` (**optional**, the default is `ScaleFactorSq(1:nc) = 1`).
1537! !
1538! ! scratch:
1539! ! \param[out] Overlap : An array of size `(nc,nc)` representing the diagonal elements of the Cholesky lower triangle of the covariance matrices of the ellipsoids.
1540! subroutine getLogVolUnion ( reltol & ! LCOV_EXCL_LINE
1541! , logVol & ! LCOV_EXCL_LINE
1542! , center & ! LCOV_EXCL_LINE
1543! , choLow & ! LCOV_EXCL_LINE
1544! , invCov & ! LCOV_EXCL_LINE
1545! , logSumVol & ! LCOV_EXCL_LINE
1546! , logSumVolEff & ! LCOV_EXCL_LINE
1547! , ScaleFactorSq & ! LCOV_EXCL_LINE
1548! !, Overlap & ! LCOV_EXCL_LINE
1549! )
1550!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1551! !DEC$ ATTRIBUTES DLLEXPORT :: getLogVolUnion
1552!#endif
1553! use pm_distUnifEll, only: setUnifEllRand, lowDia
1554! use pm_mathCumSum, only: setCumSum
1555! use pm_kind, only: RK
1556! implicit none
1557!
1558! integer(IK) , intent(in) :: reltol
1559! real(RK) , intent(in) :: logSumVol
1560! real(RK) , intent(in) :: logVol(nc)
1561! real(RK) , intent(in) :: center(nd,nc)
1562! real(RK) , intent(in) :: invCov(nd,nd,nc)
1563! real(RK) , intent(in) :: choLow(nd,nd,nc)
1564! real(RK) , intent(out) :: logSumVolEff
1565! real(RK) , intent(in) , optional :: ScaleFactorSq(nc)
1566! !integer(IK) , intent(out) :: Overlap(nc,nc)
1567!
1568! integer(IK) :: ic,jc
1569! integer(IK) :: nsim,isim
1570! integer(IK) :: membershipCount
1571! real(RK) :: CumSumVolNormed(nc)
1572! real(RK) :: ScaleFacSq(nc)
1573! real(RK) :: VolNormed(nc)
1574! real(RK) :: isimAccepted
1575! real(RK) :: simsam(nd)
1576! real(RK) :: unifrnd
1577! real(RK) :: mahalSq
1578!
1579! if (present(ScaleFactorSq)) then
1580! ScaleFacSq(1:nc) = ScaleFactorSq
1581! else
1582! ScaleFacSq(1:nc) = 1._RK
1583! end if
1584!
1585! ! Compute the expected number of simulated points in each bounded region.
1586!
1587! nsim = nint(1._RK / reltol**2)
1588! VolNormed = exp(logVol - logSumVol)
1589! call setCumSum(CumSumVolNormed, VolNormed)
1590! CumSumVolNormed(1:nc) = CumSumVolNormed(1:nc) / CumSumVolNormed(nc)
1591!
1592!#if CHECK_ENABLED
1593! if (abs(CumSumVolNormed(nc) - 1._RK) > 1.e-6) then
1594! write(*,*) "abs(CumSumVolNormed(nc) - 1._RK) > 1.e-6", CumSumVolNormed(nc)
1595! error stop
1596! end if
1597!#endif
1598!
1599! ! Simulate points within each cluster and compute the effective sum of volumes by excluding overlaps.
1600!
1601! isim = 0._RK
1602! isimAccepted = 0._RK
1603! loopGeneratePoint: do
1604!
1605! call random_number(unifrnd)
1606! ic = minloc(CumSumVolNormed, mask = CumSumVolNormed >= unifrnd, dim = 1)
1607!
1608! call setUnifEllRand( rand = simsam & ! LCOV_EXCL_LINE
1609! , mean = center(1:nd,ic) & ! LCOV_EXCL_LINE
1610! , chol = choLow(1:nd,1:nd, ic) & ! LCOV_EXCL_LINE
1611! , subset = lowDia & ! LCOV_EXCL_LINE
1612! )
1613!
1614!
1615! membershipCount = 1_IK
1616! do jc = 1, nc
1617! if (jc/=ic) then
1618! simsam(1:nd) = simsam(1:nd) - center(1:nd,jc)
1619! mahalSq = dot_product( simsam(1:nd) , matmul(invCov(1:nd,1:nd,jc), simsam(1:nd)) )
1620! simsam(1:nd) = simsam(1:nd) + center(1:nd,jc)
1621! if (mahalSq < ScaleFacSq(jc)) membershipCount = membershipCount + 1
1622! end if
1623! end do
1624!
1625! isim = isim + 1
1626! isimAccepted = isimAccepted + 1._RK / membershipCount
1627! if (isim < nsim) cycle loopGeneratePoint
1628! exit loopGeneratePoint
1629!
1630! end do loopGeneratePoint
1631!
1632! if (isimAccepted < isim) then
1633! logSumVolEff = log(isimAccepted / isim) + logSumVol
1634! else
1635! logSumVolEff = logSumVol
1636! end if
1637!
1638! end subroutine getLogVolUnion
1639!#endif
1640!
1641!
1642! !> \brief
1643! !> Return the (bounding) ellipsoid of the input set of `Point`s.
1644! !>
1645! !> \param[in] ndim : The number of dimensions.
1646! !> \param[in] np : The number of input points.
1647! !> \param[in] Point : The input array of points of shape `(ndim,np)`.
1648! !> \param[in] isBounding : A logical input value that if `.true.` will yield the bounding ellipsoid of the input `Point` dataset.
1649! !> \param[out] failed : A logical output value that is `.true.` only if the Cholesky factorization fails.
1650! !>
1651! !> \return
1652! !> `ellipsoid` : An object of type [ellipsoid_type](@ref ellipsoid_type).
1653! !>
1654! !> \warning
1655! !> If any error occurs during the ellipsoid construction, the first element of the `choDia`
1656! !> component of the output `ellipsoid` object will be set to a negative number on return.
1657! !> This is to preserve the purity of the function and to avoid the need for error handling.
1658! !>
1659! !>
1660! !>
1661! !> \remark
1662! !> When `isBounding = .true.`, the `scale`, `scaleSq`, `logVolNormed` are computed as scaled.
1663! !>
1664! !> \todo
1665! !> \phigh The `isBounding = .false.` scenario still needs work. Perhaps this argument should be changed to `boundingStatus` with three
1666! !> possible values, corresponding to three different degrees of computations: 1) no bounding, 2) computing scalefactor and logVolNormed,
1667! !> 3) computing (2) as well as scaling the properties.
1668! function ellipsoid_typer(ndim, np, Point, isBounding, failed) result(ellipsoid)
1669!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1670! !DEC$ ATTRIBUTES DLLEXPORT :: ellipsoid_typer
1671!#endif
1672! use pm_kind, only: IK, LK, RK ! LCOV_EXCL_LINE
1673! use pm_matrixInv, only: getMatInvFromChoLow
1674! use pm_matrixChol, only: setChoLow
1675! implicit none
1676! integer(IK) , intent(in) :: ndim, np
1677! real(RKG) , intent(in) :: Point(ndim,np)
1678! logical(LK) , intent(in) :: isBounding
1679! logical(LK) , intent(out) :: failed
1680! type(ellipsoid_type) :: ellipsoid
1681!
1682! real(RKG) :: scaleSqInverse
1683! real(RKG) :: NormedPoint(ndim,np)
1684! integer(IK) :: id, ip
1685!
1686! ! Compute the mean.
1687!
1688! ellipsoid%center = sum(Point, dim = 2) / np
1689!
1690! do concurrent(ip = 1:np)
1691! NormedPoint(1:ndim,ip) = Point(1:ndim,ip) - ellipsoid%center
1692! end do
1693!
1694! ! Compute the upper representative matrix of the cluster representative matrices.
1695!
1696! allocate(ellipsoid%choLowCovUpp(1:ndim,1:ndim), ellipsoid%choDia(1:ndim), source = 0._RK)
1697! do ip = 1, np
1698! do id = 1, ndim
1699! ellipsoid%choLowCovUpp(1:id,id) = ellipsoid%choLowCovUpp(1:id,id) + NormedPoint(1:id,ip) * NormedPoint(id,ip)
1700! end do
1701! end do
1702! ellipsoid%choLowCovUpp(1:ndim,1:ndim) = ellipsoid%choLowCovUpp(1:ndim,1:ndim) / real(np - 1, RK)
1703!
1704! ! Compute the Cholesky Factor of the cluster representative matrices.
1705!
1706! call setChoLow(ellipsoid%choLowCovUpp(1:ndim,1:ndim), ellipsoid%choDia(1:ndim), ndim)
1707! failed = ellipsoid%choDia(1) < 0._RK
1708! if (failed) return
1709!
1710! ! Compute the inverse of the cluster representative matrices.
1711!
1712! ellipsoid%invCov = getMatInvFromChoLow(ndim, ellipsoid%choLowCovUpp, ellipsoid%choDia)
1713!
1714! ! Compute the mahalSq of the points.
1715!
1716! allocate(ellipsoid%mahalSq(np))
1717!
1718!! if (isBounding) then
1719!
1720! ! Compute the scaleFcator of the bounding region.
1721!
1722! ellipsoid%scaleSq = NEGBIG_RK
1723! do ip = 1, np
1724! ellipsoid%mahalSq(ip) = dot_product( NormedPoint(1:ndim,ip) , matmul(ellipsoid%invCov, NormedPoint(1:ndim,ip)) )
1725! if (ellipsoid%scaleSq < ellipsoid%mahalSq(ip)) ellipsoid%scaleSq = ellipsoid%mahalSq(ip)
1726! end do
1727!
1728! ellipsoid%scale = sqrt(ellipsoid%scaleSq)
1729! scaleSqInverse = 1._RK / ellipsoid%scaleSq
1730!
1731! ellipsoid%choDia(1:ndim) = ellipsoid%choDia * ellipsoid%scale
1732! ellipsoid%invCov(1:ndim,1:ndim) = ellipsoid%invCov * scaleSqInverse
1733! do concurrent(id = 1:ndim); ellipsoid%choLowCovUpp(id+1:ndim,id) = ellipsoid%choLowCovUpp(id+1:ndim,id) * ellipsoid%scale; end do
1734!
1735!! else
1736!!
1737!! do concurrent(ip = 1:np)
1738!! ellipsoid%mahalSq(ip) = dot_product( NormedPoint(1:ndim,ip) , matmul(ellipsoid%invCov, NormedPoint(1:ndim,ip)) )
1739!! end do
1740!!
1741!! ellipsoid%scale = 1._RK
1742!! ellipsoid%scaleSq = 1._RK
1743!! ellipsoid%logVolNormed = sum(log(ellipsoid%choDia(1:ndim)))
1744!!
1745!! end if
1746!
1747! ellipsoid%logVolNormed = sum(log(ellipsoid%choDia(1:ndim)))
1748!
1749! end function ellipsoid_typer
1750
1751!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1752!
1753! pure function getLogLikeFitness(count, logVolRatio) result(logLikeFitness)
1754!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1755! !DEC$ ATTRIBUTES DLLEXPORT :: getLogLikeFitness
1756!#endif
1757! use pm_kind, only: RKG
1758! implicit none
1759! integer(IK) , intent(in) :: count
1760! real(RKG) , intent(in) :: logVolRatio
1761! real(RKG) :: logLikeFitness
1762! logLikeFitness = count * (1._RK + logVolRatio - exp(logVolRatio))
1763! end function getLogLikeFitness
1764!
1765!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1766
1767end module pm_ellipsoid ! LCOV_EXCL_LINE
Generate and return the number of points that are members (i.e., inside) of the specified -dimensiona...
Generate and return the natural logarithm of the volume of an -dimensional ellipsoid.
Generate and return the natural logarithm of the volume of an -dimensional ball of unit-radius.
Generate and return the natural logarithm of the volume of an -dimensional ball of unit-radius.
Generate and return .true. if and only if the input point is a member (i.e., inside) of the specified...
Return the natural logarithm of the volume of an -dimensional ball of unit-radius.
Return the volume of an -dimensional ball of unit-radius.
This module contains classes and procedures for setting up and computing the properties of the hyper-...
character(*, SK), parameter MODULE_NAME
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter RK5
Definition: pm_kind.F90:478
integer, parameter RK4
Definition: pm_kind.F90:489
integer, parameter RK2
Definition: pm_kind.F90:511
integer, parameter RK3
Definition: pm_kind.F90:500
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
integer, parameter RK1
Definition: pm_kind.F90:522