ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_sampleWeight.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
81
82!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
83
85
86 use pm_kind, only: SK, IK, LK
87
88 implicit none
89
90 character(*, SK), parameter :: MODULE_NAME = "@pm_sampleWeight"
91
92!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
93
115 type, abstract :: weight_type
116 end type
117
118!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
119
149 type, extends(weight_type) :: aweight_type
150 end type
151
175 type(aweight_type), parameter :: aweight = aweight_type()
176#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
177 !DIR$ ATTRIBUTES DLLEXPORT :: aweight
178#endif
179
180!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
181
209 type, extends(weight_type) :: fweight_type
210 end type
211
233 type(fweight_type), parameter :: fweight = fweight_type()
234#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
235 !DIR$ ATTRIBUTES DLLEXPORT :: fweight
236#endif
237
238!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
239
267 type, extends(weight_type) :: rweight_type
268 end type
269
290 type(rweight_type), parameter :: rweight = rweight_type()
291#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
292 !DIR$ ATTRIBUTES DLLEXPORT :: rweight
293#endif
294
295!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
296
374 interface getReweight
375
376 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
377
378#if IK5_ENABLED
379 PURE module function getReweight_IK_IK5(weight, skip) result(reweight)
380#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
381 !DEC$ ATTRIBUTES DLLEXPORT :: getReweight_IK_IK5
382#endif
383 use pm_kind, only: IKG => IK5
384 integer(IKG) , intent(in) , contiguous :: weight(:)
385 integer(IKG) , intent(in) :: skip
386 integer(IKG) :: reweight(size(weight, 1, IK))
387 end function
388#endif
389
390#if IK4_ENABLED
391 PURE module function getReweight_IK_IK4(weight, skip) result(reweight)
392#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
393 !DEC$ ATTRIBUTES DLLEXPORT :: getReweight_IK_IK4
394#endif
395 use pm_kind, only: IKG => IK4
396 integer(IKG) , intent(in) , contiguous :: weight(:)
397 integer(IKG) , intent(in) :: skip
398 integer(IKG) :: reweight(size(weight, 1, IK))
399 end function
400#endif
401
402#if IK3_ENABLED
403 PURE module function getReweight_IK_IK3(weight, skip) result(reweight)
404#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
405 !DEC$ ATTRIBUTES DLLEXPORT :: getReweight_IK_IK3
406#endif
407 use pm_kind, only: IKG => IK3
408 integer(IKG) , intent(in) , contiguous :: weight(:)
409 integer(IKG) , intent(in) :: skip
410 integer(IKG) :: reweight(size(weight, 1, IK))
411 end function
412#endif
413
414#if IK2_ENABLED
415 PURE module function getReweight_IK_IK2(weight, skip) result(reweight)
416#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
417 !DEC$ ATTRIBUTES DLLEXPORT :: getReweight_IK_IK2
418#endif
419 use pm_kind, only: IKG => IK2
420 integer(IKG) , intent(in) , contiguous :: weight(:)
421 integer(IKG) , intent(in) :: skip
422 integer(IKG) :: reweight(size(weight, 1, IK))
423 end function
424#endif
425
426#if IK1_ENABLED
427 PURE module function getReweight_IK_IK1(weight, skip) result(reweight)
428#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
429 !DEC$ ATTRIBUTES DLLEXPORT :: getReweight_IK_IK1
430#endif
431 use pm_kind, only: IKG => IK1
432 integer(IKG) , intent(in) , contiguous :: weight(:)
433 integer(IKG) , intent(in) :: skip
434 integer(IKG) :: reweight(size(weight, 1, IK))
435 end function
436#endif
437
438 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
439
440#if RK5_ENABLED
441 PURE module function getReweight_RK_RK5(weight, skip) result(reweight)
442#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
443 !DEC$ ATTRIBUTES DLLEXPORT :: getReweight_RK_RK5
444#endif
445 use pm_kind, only: RKG => RK5
446 real(RKG) , intent(in) , contiguous :: weight(:)
447 real(RKG) , intent(in) :: skip
448 real(RKG) :: reweight(size(weight, 1, IK))
449 end function
450#endif
451
452#if RK4_ENABLED
453 PURE module function getReweight_RK_RK4(weight, skip) result(reweight)
454#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
455 !DEC$ ATTRIBUTES DLLEXPORT :: getReweight_RK_RK4
456#endif
457 use pm_kind, only: RKG => RK4
458 real(RKG) , intent(in) , contiguous :: weight(:)
459 real(RKG) , intent(in) :: skip
460 real(RKG) :: reweight(size(weight, 1, IK))
461 end function
462#endif
463
464#if RK3_ENABLED
465 PURE module function getReweight_RK_RK3(weight, skip) result(reweight)
466#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
467 !DEC$ ATTRIBUTES DLLEXPORT :: getReweight_RK_RK3
468#endif
469 use pm_kind, only: RKG => RK3
470 real(RKG) , intent(in) , contiguous :: weight(:)
471 real(RKG) , intent(in) :: skip
472 real(RKG) :: reweight(size(weight, 1, IK))
473 end function
474#endif
475
476#if RK2_ENABLED
477 PURE module function getReweight_RK_RK2(weight, skip) result(reweight)
478#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
479 !DEC$ ATTRIBUTES DLLEXPORT :: getReweight_RK_RK2
480#endif
481 use pm_kind, only: RKG => RK2
482 real(RKG) , intent(in) , contiguous :: weight(:)
483 real(RKG) , intent(in) :: skip
484 real(RKG) :: reweight(size(weight, 1, IK))
485 end function
486#endif
487
488#if RK1_ENABLED
489 PURE module function getReweight_RK_RK1(weight, skip) result(reweight)
490#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
491 !DEC$ ATTRIBUTES DLLEXPORT :: getReweight_RK_RK1
492#endif
493 use pm_kind, only: RKG => RK1
494 real(RKG) , intent(in) , contiguous :: weight(:)
495 real(RKG) , intent(in) :: skip
496 real(RKG) :: reweight(size(weight, 1, IK))
497 end function
498#endif
499
500 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
501
502#if RK5_ENABLED
503 PURE module function getReweight_IK_RK5(weight, skip) result(reweight)
504#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
505 !DEC$ ATTRIBUTES DLLEXPORT :: getReweight_IK_RK5
506#endif
507 use pm_kind, only: RKG => RK5
508 integer(IK) , intent(in) , contiguous :: weight(:)
509 real(RKG) , intent(in) :: skip
510 integer(IK) :: reweight(size(weight, 1, IK))
511 end function
512#endif
513
514#if RK4_ENABLED
515 PURE module function getReweight_IK_RK4(weight, skip) result(reweight)
516#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
517 !DEC$ ATTRIBUTES DLLEXPORT :: getReweight_IK_RK4
518#endif
519 use pm_kind, only: RKG => RK4
520 integer(IK) , intent(in) , contiguous :: weight(:)
521 real(RKG) , intent(in) :: skip
522 integer(IK) :: reweight(size(weight, 1, IK))
523 end function
524#endif
525
526#if RK3_ENABLED
527 PURE module function getReweight_IK_RK3(weight, skip) result(reweight)
528#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
529 !DEC$ ATTRIBUTES DLLEXPORT :: getReweight_IK_RK3
530#endif
531 use pm_kind, only: RKG => RK3
532 integer(IK) , intent(in) , contiguous :: weight(:)
533 real(RKG) , intent(in) :: skip
534 integer(IK) :: reweight(size(weight, 1, IK))
535 end function
536#endif
537
538#if RK2_ENABLED
539 PURE module function getReweight_IK_RK2(weight, skip) result(reweight)
540#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
541 !DEC$ ATTRIBUTES DLLEXPORT :: getReweight_IK_RK2
542#endif
543 use pm_kind, only: RKG => RK2
544 integer(IK) , intent(in) , contiguous :: weight(:)
545 real(RKG) , intent(in) :: skip
546 integer(IK) :: reweight(size(weight, 1, IK))
547 end function
548#endif
549
550#if RK1_ENABLED
551 PURE module function getReweight_IK_RK1(weight, skip) result(reweight)
552#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
553 !DEC$ ATTRIBUTES DLLEXPORT :: getReweight_IK_RK1
554#endif
555 use pm_kind, only: RKG => RK1
556 integer(IK) , intent(in) , contiguous :: weight(:)
557 real(RKG) , intent(in) :: skip
558 integer(IK) :: reweight(size(weight, 1, IK))
559 end function
560#endif
561
562 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
563
564 end interface
565
566!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
567
640 interface setReweight
641
642 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
643
644#if IK5_ENABLED
645 PURE module subroutine setReweight_IK_IK5(weight, skip)
646#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
647 !DEC$ ATTRIBUTES DLLEXPORT :: setReweight_IK_IK5
648#endif
649 use pm_kind, only: IKG => IK5
650 integer(IKG) , intent(inout) , contiguous :: weight(:)
651 integer(IKG) , intent(in) :: skip
652 integer(IKG) :: reweight(size(weight, 1, IK))
653 end subroutine
654#endif
655
656#if IK4_ENABLED
657 PURE module subroutine setReweight_IK_IK4(weight, skip)
658#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
659 !DEC$ ATTRIBUTES DLLEXPORT :: setReweight_IK_IK4
660#endif
661 use pm_kind, only: IKG => IK4
662 integer(IKG) , intent(inout) , contiguous :: weight(:)
663 integer(IKG) , intent(in) :: skip
664 integer(IKG) :: reweight(size(weight, 1, IK))
665 end subroutine
666#endif
667
668#if IK3_ENABLED
669 PURE module subroutine setReweight_IK_IK3(weight, skip)
670#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
671 !DEC$ ATTRIBUTES DLLEXPORT :: setReweight_IK_IK3
672#endif
673 use pm_kind, only: IKG => IK3
674 integer(IKG) , intent(inout) , contiguous :: weight(:)
675 integer(IKG) , intent(in) :: skip
676 integer(IKG) :: reweight(size(weight, 1, IK))
677 end subroutine
678#endif
679
680#if IK2_ENABLED
681 PURE module subroutine setReweight_IK_IK2(weight, skip)
682#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
683 !DEC$ ATTRIBUTES DLLEXPORT :: setReweight_IK_IK2
684#endif
685 use pm_kind, only: IKG => IK2
686 integer(IKG) , intent(inout) , contiguous :: weight(:)
687 integer(IKG) , intent(in) :: skip
688 integer(IKG) :: reweight(size(weight, 1, IK))
689 end subroutine
690#endif
691
692#if IK1_ENABLED
693 PURE module subroutine setReweight_IK_IK1(weight, skip)
694#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
695 !DEC$ ATTRIBUTES DLLEXPORT :: setReweight_IK_IK1
696#endif
697 use pm_kind, only: IKG => IK1
698 integer(IKG) , intent(inout) , contiguous :: weight(:)
699 integer(IKG) , intent(in) :: skip
700 integer(IKG) :: reweight(size(weight, 1, IK))
701 end subroutine
702#endif
703
704 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
705
706#if RK5_ENABLED
707 PURE module subroutine setReweight_RK_RK5(weight, skip)
708#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
709 !DEC$ ATTRIBUTES DLLEXPORT :: setReweight_RK_RK5
710#endif
711 use pm_kind, only: RKG => RK5
712 real(RKG) , intent(inout) , contiguous :: weight(:)
713 real(RKG) , intent(in) :: skip
714 end subroutine
715#endif
716
717#if RK4_ENABLED
718 PURE module subroutine setReweight_RK_RK4(weight, skip)
719#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
720 !DEC$ ATTRIBUTES DLLEXPORT :: setReweight_RK_RK4
721#endif
722 use pm_kind, only: RKG => RK4
723 real(RKG) , intent(inout) , contiguous :: weight(:)
724 real(RKG) , intent(in) :: skip
725 end subroutine
726#endif
727
728#if RK3_ENABLED
729 PURE module subroutine setReweight_RK_RK3(weight, skip)
730#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
731 !DEC$ ATTRIBUTES DLLEXPORT :: setReweight_RK_RK3
732#endif
733 use pm_kind, only: RKG => RK3
734 real(RKG) , intent(inout) , contiguous :: weight(:)
735 real(RKG) , intent(in) :: skip
736 end subroutine
737#endif
738
739#if RK2_ENABLED
740 PURE module subroutine setReweight_RK_RK2(weight, skip)
741#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
742 !DEC$ ATTRIBUTES DLLEXPORT :: setReweight_RK_RK2
743#endif
744 use pm_kind, only: RKG => RK2
745 real(RKG) , intent(inout) , contiguous :: weight(:)
746 real(RKG) , intent(in) :: skip
747 end subroutine
748#endif
749
750#if RK1_ENABLED
751 PURE module subroutine setReweight_RK_RK1(weight, skip)
752#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
753 !DEC$ ATTRIBUTES DLLEXPORT :: setReweight_RK_RK1
754#endif
755 use pm_kind, only: RKG => RK1
756 real(RKG) , intent(inout) , contiguous :: weight(:)
757 real(RKG) , intent(in) :: skip
758 end subroutine
759#endif
760
761 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
762
763#if RK5_ENABLED
764 PURE module subroutine setReweight_IK_RK5(weight, skip)
765#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
766 !DEC$ ATTRIBUTES DLLEXPORT :: setReweight_IK_RK5
767#endif
768 use pm_kind, only: RKG => RK5
769 integer(IK) , intent(inout) , contiguous :: weight(:)
770 real(RKG) , intent(in) :: skip
771 end subroutine
772#endif
773
774#if RK4_ENABLED
775 PURE module subroutine setReweight_IK_RK4(weight, skip)
776#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
777 !DEC$ ATTRIBUTES DLLEXPORT :: setReweight_IK_RK4
778#endif
779 use pm_kind, only: RKG => RK4
780 integer(IK) , intent(inout) , contiguous :: weight(:)
781 real(RKG) , intent(in) :: skip
782 end subroutine
783#endif
784
785#if RK3_ENABLED
786 PURE module subroutine setReweight_IK_RK3(weight, skip)
787#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
788 !DEC$ ATTRIBUTES DLLEXPORT :: setReweight_IK_RK3
789#endif
790 use pm_kind, only: RKG => RK3
791 integer(IK) , intent(inout) , contiguous :: weight(:)
792 real(RKG) , intent(in) :: skip
793 end subroutine
794#endif
795
796#if RK2_ENABLED
797 PURE module subroutine setReweight_IK_RK2(weight, skip)
798#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
799 !DEC$ ATTRIBUTES DLLEXPORT :: setReweight_IK_RK2
800#endif
801 use pm_kind, only: RKG => RK2
802 integer(IK) , intent(inout) , contiguous :: weight(:)
803 real(RKG) , intent(in) :: skip
804 end subroutine
805#endif
806
807#if RK1_ENABLED
808 PURE module subroutine setReweight_IK_RK1(weight, skip)
809#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
810 !DEC$ ATTRIBUTES DLLEXPORT :: setReweight_IK_RK1
811#endif
812 use pm_kind, only: RKG => RK1
813 integer(IK) , intent(inout) , contiguous :: weight(:)
814 real(RKG) , intent(in) :: skip
815 end subroutine
816#endif
817
818 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
819
820 end interface
821
822!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
823
824end module pm_sampleWeight ! LCOV_EXCL_LINE
Generate and return a reweighting of the input weight vector, such that the sequence represented by o...
Generate and return a reweighting of the input weight vector, such that the sequence represented by o...
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter IK3
Definition: pm_kind.F90:368
integer, parameter RK5
Definition: pm_kind.F90:478
integer, parameter RK4
Definition: pm_kind.F90:489
integer, parameter IK1
Definition: pm_kind.F90:382
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 IK2
Definition: pm_kind.F90:375
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 IK4
Definition: pm_kind.F90:361
integer, parameter IK5
Definition: pm_kind.F90:354
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 the types, classes, and procedures relevant to weights of random samples.
type(rweight_type), parameter rweight
This is a scalar parameter object of type rweight_type that is exclusively used to signify the Symmet...
type(fweight_type), parameter fweight
This is a scalar parameter object of type fweight_type that is exclusively used to signify the Symmet...
character(*, SK), parameter MODULE_NAME
type(aweight_type), parameter aweight
This is a scalar parameter object of type aweight_type that is exclusively used to signify the Symmet...
This is a concrete derived type whose instances are exclusively used to signify the aweight type of s...
This is a concrete derived type whose instances are exclusively used to signify the fweight type of s...
This is a concrete derived type whose instances are exclusively used to signify the rweight type of s...
This is an abstract derived type for constructing concrete derived types to distinguish various proce...