ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_distanceHellinger.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
149
150!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
151
153
154 use pm_kind, only: SK, IK, LK
155 use pm_mathConst, only: ninf_type
156 use pm_mathConst, only: pinf_type
157 use pm_except, only: getInfNeg
158 use pm_except, only: getInfPos
159
160 implicit none
161
162 character(*, SK), parameter :: MODULE_NAME = "@pm_distanceHellinger"
163
164!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
165
271
272 ! PMF
273
274 interface getDisHellSq
275
276 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
277
278#if RK5_ENABLED
279 PURE module function getDisHellSqPMF_RK5(p, q) result(hellsq)
280#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
281 !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPMF_RK5
282#endif
283 use pm_kind, only: RKG => RK5
284 real(RKG) , intent(in), contiguous :: p(:), q(:)
285 real(RKG) :: hellsq
286 end function
287#endif
288
289#if RK4_ENABLED
290 PURE module function getDisHellSqPMF_RK4(p, q) result(hellsq)
291#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
292 !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPMF_RK4
293#endif
294 use pm_kind, only: RKG => RK4
295 real(RKG) , intent(in), contiguous :: p(:), q(:)
296 real(RKG) :: hellsq
297 end function
298#endif
299
300#if RK3_ENABLED
301 PURE module function getDisHellSqPMF_RK3(p, q) result(hellsq)
302#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
303 !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPMF_RK3
304#endif
305 use pm_kind, only: RKG => RK3
306 real(RKG) , intent(in), contiguous :: p(:), q(:)
307 real(RKG) :: hellsq
308 end function
309#endif
310
311#if RK2_ENABLED
312 PURE module function getDisHellSqPMF_RK2(p, q) result(hellsq)
313#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
314 !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPMF_RK2
315#endif
316 use pm_kind, only: RKG => RK2
317 real(RKG) , intent(in), contiguous :: p(:), q(:)
318 real(RKG) :: hellsq
319 end function
320#endif
321
322#if RK1_ENABLED
323 PURE module function getDisHellSqPMF_RK1(p, q) result(hellsq)
324#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
325 !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPMF_RK1
326#endif
327 use pm_kind, only: RKG => RK1
328 real(RKG) , intent(in), contiguous :: p(:), q(:)
329 real(RKG) :: hellsq
330 end function
331#endif
332
333 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
334
335 end interface
336
337 ! Proc
338
339 interface getDisHellSq
340
341 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
342 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
343 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
344
345 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
346
347#if RK5_ENABLED
348 module function getDisHellSqPRC_RK5(p, q, lb, ub, failed, msg) result(hellsq)
349#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
350 !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_RK5
351#endif
352 use pm_kind, only: RKG => RK5
353 procedure(real(RKG)) :: p, q
354 real(RKG) , intent(in) , optional :: lb
355 real(RKG) , intent(in) , optional :: ub
356 character(*, SK), intent(out) , optional :: msg
357 logical(LK) , intent(out) , optional :: failed
358 real(RKG) :: hellsq
359 end function
360#endif
361
362#if RK4_ENABLED
363 module function getDisHellSqPRC_RK4(p, q, lb, ub, failed, msg) result(hellsq)
364#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
365 !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_RK4
366#endif
367 use pm_kind, only: RKG => RK4
368 procedure(real(RKG)) :: p, q
369 real(RKG) , intent(in) , optional :: lb
370 real(RKG) , intent(in) , optional :: ub
371 character(*, SK), intent(out) , optional :: msg
372 logical(LK) , intent(out) , optional :: failed
373 real(RKG) :: hellsq
374 end function
375#endif
376
377#if RK3_ENABLED
378 module function getDisHellSqPRC_RK3(p, q, lb, ub, failed, msg) result(hellsq)
379#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
380 !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_RK3
381#endif
382 use pm_kind, only: RKG => RK3
383 procedure(real(RKG)) :: p, q
384 real(RKG) , intent(in) , optional :: lb
385 real(RKG) , intent(in) , optional :: ub
386 character(*, SK), intent(out) , optional :: msg
387 logical(LK) , intent(out) , optional :: failed
388 real(RKG) :: hellsq
389 end function
390#endif
391
392#if RK2_ENABLED
393 module function getDisHellSqPRC_RK2(p, q, lb, ub, failed, msg) result(hellsq)
394#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
395 !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_RK2
396#endif
397 use pm_kind, only: RKG => RK2
398 procedure(real(RKG)) :: p, q
399 real(RKG) , intent(in) , optional :: lb
400 real(RKG) , intent(in) , optional :: ub
401 character(*, SK), intent(out) , optional :: msg
402 logical(LK) , intent(out) , optional :: failed
403 real(RKG) :: hellsq
404 end function
405#endif
406
407#if RK1_ENABLED
408 module function getDisHellSqPRC_RK1(p, q, lb, ub, failed, msg) result(hellsq)
409#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
410 !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_RK1
411#endif
412 use pm_kind, only: RKG => RK1
413 procedure(real(RKG)) :: p, q
414 real(RKG) , intent(in) , optional :: lb
415 real(RKG) , intent(in) , optional :: ub
416 character(*, SK), intent(out) , optional :: msg
417 logical(LK) , intent(out) , optional :: failed
418 real(RKG) :: hellsq
419 end function
420#endif
421
422 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
423
424 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
425 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
426 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
427!
428! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
429!
430!#if RK5_ENABLED
431! module function getDisHellSqPRC_FF_RK5(p, q, lb, ub) result(hellsq)
432!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
433! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_FF_RK5
434!#endif
435! use pm_kind, only: RKG => RK5
436! procedure(real(RKG)) :: p, q
437! real(RKG) , intent(in) :: lb
438! real(RKG) , intent(in) :: ub
439! real(RKG) :: hellsq
440! end function
441!#endif
442!
443!#if RK4_ENABLED
444! module function getDisHellSqPRC_FF_RK4(p, q, lb, ub) result(hellsq)
445!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
446! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_FF_RK4
447!#endif
448! use pm_kind, only: RKG => RK4
449! procedure(real(RKG)) :: p, q
450! real(RKG) , intent(in) :: lb
451! real(RKG) , intent(in) :: ub
452! real(RKG) :: hellsq
453! end function
454!#endif
455!
456!#if RK3_ENABLED
457! module function getDisHellSqPRC_FF_RK3(p, q, lb, ub) result(hellsq)
458!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
459! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_FF_RK3
460!#endif
461! use pm_kind, only: RKG => RK3
462! procedure(real(RKG)) :: p, q
463! real(RKG) , intent(in) :: lb
464! real(RKG) , intent(in) :: ub
465! real(RKG) :: hellsq
466! end function
467!#endif
468!
469!#if RK2_ENABLED
470! module function getDisHellSqPRC_FF_RK2(p, q, lb, ub) result(hellsq)
471!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
472! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_FF_RK2
473!#endif
474! use pm_kind, only: RKG => RK2
475! procedure(real(RKG)) :: p, q
476! real(RKG) , intent(in) :: lb
477! real(RKG) , intent(in) :: ub
478! real(RKG) :: hellsq
479! end function
480!#endif
481!
482!#if RK1_ENABLED
483! module function getDisHellSqPRC_FF_RK1(p, q, lb, ub) result(hellsq)
484!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
485! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_FF_RK1
486!#endif
487! use pm_kind, only: RKG => RK1
488! procedure(real(RKG)) :: p, q
489! real(RKG) , intent(in) :: lb
490! real(RKG) , intent(in) :: ub
491! real(RKG) :: hellsq
492! end function
493!#endif
494!
495! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
496!
497! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
498! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
499! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
500!
501! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
502!
503!#if RK5_ENABLED
504! module function getDisHellSqPRC_FI_RK5(p, q, lb, ub) result(hellsq)
505!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
506! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_FI_RK5
507!#endif
508! use pm_kind, only: RKG => RK5
509! procedure(real(RKG)) :: p, q
510! real(RKG) , intent(in) :: lb
511! type(pinf_type) , intent(in) :: ub
512! real(RKG) :: hellsq
513! end function
514!#endif
515!
516!#if RK4_ENABLED
517! module function getDisHellSqPRC_FI_RK4(p, q, lb, ub) result(hellsq)
518!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
519! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_FI_RK4
520!#endif
521! use pm_kind, only: RKG => RK4
522! procedure(real(RKG)) :: p, q
523! real(RKG) , intent(in) :: lb
524! type(pinf_type) , intent(in) :: ub
525! real(RKG) :: hellsq
526! end function
527!#endif
528!
529!#if RK3_ENABLED
530! module function getDisHellSqPRC_FI_RK3(p, q, lb, ub) result(hellsq)
531!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
532! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_FI_RK3
533!#endif
534! use pm_kind, only: RKG => RK3
535! procedure(real(RKG)) :: p, q
536! real(RKG) , intent(in) :: lb
537! type(pinf_type) , intent(in) :: ub
538! real(RKG) :: hellsq
539! end function
540!#endif
541!
542!#if RK2_ENABLED
543! module function getDisHellSqPRC_FI_RK2(p, q, lb, ub) result(hellsq)
544!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
545! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_FI_RK2
546!#endif
547! use pm_kind, only: RKG => RK2
548! procedure(real(RKG)) :: p, q
549! real(RKG) , intent(in) :: lb
550! type(pinf_type) , intent(in) :: ub
551! real(RKG) :: hellsq
552! end function
553!#endif
554!
555!#if RK1_ENABLED
556! module function getDisHellSqPRC_FI_RK1(p, q, lb, ub) result(hellsq)
557!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
558! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_FI_RK1
559!#endif
560! use pm_kind, only: RKG => RK1
561! procedure(real(RKG)) :: p, q
562! real(RKG) , intent(in) :: lb
563! type(pinf_type) , intent(in) :: ub
564! real(RKG) :: hellsq
565! end function
566!#endif
567!
568! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
569!
570! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
571! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
572! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
573!
574! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
575!
576!#if RK5_ENABLED
577! module function getDisHellSqPRC_IF_RK5(p, q, lb, ub) result(hellsq)
578!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
579! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_IF_RK5
580!#endif
581! use pm_kind, only: RKG => RK5
582! procedure(real(RKG)) :: p, q
583! type(ninf_type) , intent(in) :: lb
584! real(RKG) , intent(in) :: ub
585! real(RKG) :: hellsq
586! end function
587!#endif
588!
589!#if RK4_ENABLED
590! module function getDisHellSqPRC_IF_RK4(p, q, lb, ub) result(hellsq)
591!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
592! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_IF_RK4
593!#endif
594! use pm_kind, only: RKG => RK4
595! procedure(real(RKG)) :: p, q
596! type(ninf_type) , intent(in) :: lb
597! real(RKG) , intent(in) :: ub
598! real(RKG) :: hellsq
599! end function
600!#endif
601!
602!#if RK3_ENABLED
603! module function getDisHellSqPRC_IF_RK3(p, q, lb, ub) result(hellsq)
604!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
605! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_IF_RK3
606!#endif
607! use pm_kind, only: RKG => RK3
608! procedure(real(RKG)) :: p, q
609! type(ninf_type) , intent(in) :: lb
610! real(RKG) , intent(in) :: ub
611! real(RKG) :: hellsq
612! end function
613!#endif
614!
615!#if RK2_ENABLED
616! module function getDisHellSqPRC_IF_RK2(p, q, lb, ub) result(hellsq)
617!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
618! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_IF_RK2
619!#endif
620! use pm_kind, only: RKG => RK2
621! procedure(real(RKG)) :: p, q
622! type(ninf_type) , intent(in) :: lb
623! real(RKG) , intent(in) :: ub
624! real(RKG) :: hellsq
625! end function
626!#endif
627!
628!#if RK1_ENABLED
629! module function getDisHellSqPRC_IF_RK1(p, q, lb, ub) result(hellsq)
630!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
631! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_IF_RK1
632!#endif
633! use pm_kind, only: RKG => RK1
634! procedure(real(RKG)) :: p, q
635! type(ninf_type) , intent(in) :: lb
636! real(RKG) , intent(in) :: ub
637! real(RKG) :: hellsq
638! end function
639!#endif
640!
641! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
642!
643! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
644! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
645! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
646!
647! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
648!
649!#if RK5_ENABLED
650! module function getDisHellSqPRC_II_RK5(p, q, lb, ub) result(hellsq)
651!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
652! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_II_RK5
653!#endif
654! use pm_kind, only: RKG => RK5
655! procedure(real(RKG)) :: p, q
656! type(ninf_type) , intent(in) :: ninf_type
657! type(pinf_type) , intent(in) :: pinf_type
658! real(RKG) :: hellsq
659! end function
660!#endif
661!
662!#if RK4_ENABLED
663! module function getDisHellSqPRC_II_RK4(p, q, lb, ub) result(hellsq)
664!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
665! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_II_RK4
666!#endif
667! use pm_kind, only: RKG => RK4
668! procedure(real(RKG)) :: p, q
669! type(ninf_type) , intent(in) :: ninf_type
670! type(pinf_type) , intent(in) :: pinf_type
671! real(RKG) :: hellsq
672! end function
673!#endif
674!
675!#if RK3_ENABLED
676! module function getDisHellSqPRC_II_RK3(p, q, lb, ub) result(hellsq)
677!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
678! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_II_RK3
679!#endif
680! use pm_kind, only: RKG => RK3
681! procedure(real(RKG)) :: p, q
682! type(ninf_type) , intent(in) :: ninf_type
683! type(pinf_type) , intent(in) :: pinf_type
684! real(RKG) :: hellsq
685! end function
686!#endif
687!
688!#if RK2_ENABLED
689! module function getDisHellSqPRC_II_RK2(p, q, lb, ub) result(hellsq)
690!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
691! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_II_RK2
692!#endif
693! use pm_kind, only: RKG => RK2
694! procedure(real(RKG)) :: p, q
695! type(ninf_type) , intent(in) :: ninf_type
696! type(pinf_type) , intent(in) :: pinf_type
697! real(RKG) :: hellsq
698! end function
699!#endif
700!
701!#if RK1_ENABLED
702! module function getDisHellSqPRC_II_RK1(p, q, lb, ub) result(hellsq)
703!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
704! !DEC$ ATTRIBUTES DLLEXPORT :: getDisHellSqPRC_II_RK1
705!#endif
706! use pm_kind, only: RKG => RK1
707! procedure(real(RKG)) :: p, q
708! type(ninf_type) , intent(in) :: ninf_type
709! type(pinf_type) , intent(in) :: pinf_type
710! real(RKG) :: hellsq
711! end function
712!#endif
713!
714! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
715!
716! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
717! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
718! !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
719 end interface
720
721!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
722
723end module pm_distanceHellinger ! LCOV_EXCL_LINE
Generate and return the square of the Hellinger distance of two univariate (discrete or continuous) d...
Generate and return an IEEE-compliant negative infinity.
Definition: pm_except.F90:1719
Generate and return an IEEE-compliant positive infinity.
Definition: pm_except.F90:1189
This module contains classes and procedures for computing the Hellinger statistical distance between ...
character(*, SK), parameter MODULE_NAME
This module contains procedures and generic interfaces for testing for exceptional cases at runtime.
Definition: pm_except.F90:46
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
This module contains relevant mathematical constants.
This is the indicator type for generating instances of objects that indicate the use of the negative ...
This is the indicator type for generating instances of objects that indicate the use of the positive ...