ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_batse.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
66
67!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
68
70
71 use pm_kind, only: SK, IK, RKD, RKB
72
73 implicit none
74
75 character(*, SK), parameter :: MODULE_NAME = "@pm_batse"
76
77 real(RKB) , parameter :: LOG_TEN = log(10._RKB)
78 real(RKB) , parameter :: MIN_LOG10PH53_4_LOGPBOLZERO = 4.92_RKB
79 real(RKB) , parameter :: MAX_LOG10PH53_4_LOGPBOLZERO = 6.318167895318538_RKB
83
84 ! Parameters for the effective peak flux of SGRBs (P64ms conversion to P1024):
85
86! ! The scale of the change in BATSE efficiency for different GRB durations.
87! real(RKB) , parameter :: THRESH_ERFC_BASE = +0.146314238936889_RKB
88!
89! ! The scale of the change in BATSE efficiency for different GRB durations.
90! real(RKB) , parameter :: THRESH_ERFC_AMP = +0.570285374263156_RKB
91!
92! ! Mean duration in the Error function used to model the connection between the peak fluxes in 64 and 1024 ms.
93! real(RKB) , parameter :: THRESH_ERFC_AVG = -0.480811530417719_RKB
94!
95! ! scale of the duration in the Error function used to model the connection between the peak fluxes in 64 and 1024 ms.
96! real(RKB) , parameter :: THRESH_ERFC_STD = +1.292443569094922_RKB
97
98 real(RKB) , parameter :: LOGPF53_MINUS_LOGPBOL = 11.328718657530706_RKB
100
102 real(RKB) , parameter :: THRESH_ERFC_BASE = +0.146314238936889_RKB * LOG_TEN
103
105 real(RKB) , parameter :: THRESH_ERFC_AMP = +0.282313526464596_RKB * LOG_TEN ! Serfc
106
108 real(RKB) , parameter :: THRESH_ERFC_AVG = -0.483553339256463_RKB * LOG_TEN ! meandur
109
111 real(RKB) , parameter :: THRESH_ERFC_STD = 1.0514698984694800_RKB * LOG_TEN ! scaledur
112
114 real(RKB) , parameter :: THRESH_ERFC_STD_INV = 1._RKB / THRESH_ERFC_STD ! inverse scaledur
115
116 ! The height of the ERFC function.
117 real(RKB) , parameter :: THRESH_ERFC_HEIGHT = 2 * THRESH_ERFC_AMP
118
120 ! Effective LogPbol limit above which trigger efficiency is 100%, for any Log(Epk) and Log(dur). It is equivalent to maximum Log(Pbol) at very long durations.
122
123!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
124
135 type :: grb_type
136 real(RKD) :: pbol = -huge(0._RKD) ! < \public The scalar `real` of kind \RKD, containing a measure of the bolometric energy flux of a BATSE GRB event.
137 real(RKD) :: epeak = -huge(0._RKD) ! < \public The scalar `real` of kind \RKD, containing a measure of the spectral peak energy of a BATSE GRB event.
138 real(RKD) :: sbol = -huge(0._RKD) ! < \public The scalar `real` of kind \RKD, containing a measure of the bolometric energy fluence of a BATSE GRB event.
139 real(RKD) :: t90 = -huge(0._RKD) ! < \public The scalar `real` of kind \RKD, containing a measure of the duration of a BATSE GRB event.
140 real(RKD) :: pph1024 = -huge(0._RKD) ! < \public The scalar `real` of kind \RKD, containing a measure of the peak 1024ms photon flux of a BATSE GRB event in the BATSE nominal detection energy window \f$[50, 300]\f$.
141 end type
142
143!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
144
192
193 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
194
195#if RK5_ENABLED
196 pure elemental module function getCorrectionLogEffPPF_RK5(logT90) result(correction)
197#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
198 !DEC$ ATTRIBUTES DLLEXPORT :: getCorrectionLogEffPPF_RK5
199#endif
200 use pm_kind, only: RKG => RK5
201 real(RKG) , intent(in) :: logT90
202 real(RKG) :: correction
203 end function
204#endif
205
206#if RK4_ENABLED
207 pure elemental module function getCorrectionLogEffPPF_RK4(logT90) result(correction)
208#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
209 !DEC$ ATTRIBUTES DLLEXPORT :: getCorrectionLogEffPPF_RK4
210#endif
211 use pm_kind, only: RKG => RK4
212 real(RKG) , intent(in) :: logT90
213 real(RKG) :: correction
214 end function
215#endif
216
217#if RK3_ENABLED
218 pure elemental module function getCorrectionLogEffPPF_RK3(logT90) result(correction)
219#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
220 !DEC$ ATTRIBUTES DLLEXPORT :: getCorrectionLogEffPPF_RK3
221#endif
222 use pm_kind, only: RKG => RK3
223 real(RKG) , intent(in) :: logT90
224 real(RKG) :: correction
225 end function
226#endif
227
228#if RK2_ENABLED
229 pure elemental module function getCorrectionLogEffPPF_RK2(logT90) result(correction)
230#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
231 !DEC$ ATTRIBUTES DLLEXPORT :: getCorrectionLogEffPPF_RK2
232#endif
233 use pm_kind, only: RKG => RK2
234 real(RKG) , intent(in) :: logT90
235 real(RKG) :: correction
236 end function
237#endif
238
239#if RK1_ENABLED
240 pure elemental module function getCorrectionLogEffPPF_RK1(logT90) result(correction)
241#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
242 !DEC$ ATTRIBUTES DLLEXPORT :: getCorrectionLogEffPPF_RK1
243#endif
244 use pm_kind, only: RKG => RK1
245 real(RKG) , intent(in) :: logT90
246 real(RKG) :: correction
247 end function
248#endif
249
250 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
251
252 end interface getCorrectionLogEffPPF
253
254!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
255
306 interface getLogEffPPF
307
308 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
309
310#if RK5_ENABLED
311 pure elemental module function getLogEffPPF_RK5(logPPF64, logT90) result(logEffPPF)
312#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
313 !DEC$ ATTRIBUTES DLLEXPORT :: getLogEffPPF_RK5
314#endif
315 use pm_kind, only: RKG => RK5
316 real(RKG) , intent(in) :: logPPF64, logT90
317 real(RKG) :: logEffPPF
318 end function
319#endif
320
321#if RK4_ENABLED
322 pure elemental module function getLogEffPPF_RK4(logPPF64, logT90) result(logEffPPF)
323#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
324 !DEC$ ATTRIBUTES DLLEXPORT :: getLogEffPPF_RK4
325#endif
326 use pm_kind, only: RKG => RK4
327 real(RKG) , intent(in) :: logPPF64, logT90
328 real(RKG) :: logEffPPF
329 end function
330#endif
331
332#if RK3_ENABLED
333 pure elemental module function getLogEffPPF_RK3(logPPF64, logT90) result(logEffPPF)
334#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
335 !DEC$ ATTRIBUTES DLLEXPORT :: getLogEffPPF_RK3
336#endif
337 use pm_kind, only: RKG => RK3
338 real(RKG) , intent(in) :: logPPF64, logT90
339 real(RKG) :: logEffPPF
340 end function
341#endif
342
343#if RK2_ENABLED
344 pure elemental module function getLogEffPPF_RK2(logPPF64, logT90) result(logEffPPF)
345#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
346 !DEC$ ATTRIBUTES DLLEXPORT :: getLogEffPPF_RK2
347#endif
348 use pm_kind, only: RKG => RK2
349 real(RKG) , intent(in) :: logPPF64, logT90
350 real(RKG) :: logEffPPF
351 end function
352#endif
353
354#if RK1_ENABLED
355 pure elemental module function getLogEffPPF_RK1(logPPF64, logT90) result(logEffPPF)
356#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
357 !DEC$ ATTRIBUTES DLLEXPORT :: getLogEffPPF_RK1
358#endif
359 use pm_kind, only: RKG => RK1
360 real(RKG) , intent(in) :: logPPF64, logT90
361 real(RKG) :: logEffPPF
362 end function
363#endif
364
365 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
366
367 end interface getLogEffPPF
368
369!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
370
427 interface getLogPbol
428
429 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
430
431#if RK5_ENABLED
432 pure elemental module function getLogPbol_RK5(logEpk, logPF53) result(logPbol)
433#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
434 !DEC$ ATTRIBUTES DLLEXPORT :: getLogPbol_RK5
435#endif
436 use pm_kind, only: RKG => RK5
437 real(RKG) , intent(in) :: logEpk, logPF53
438 real(RKG) :: logPbol
439 end function
440#endif
441
442#if RK4_ENABLED
443 pure elemental module function getLogPbol_RK4(logEpk, logPF53) result(logPbol)
444#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
445 !DEC$ ATTRIBUTES DLLEXPORT :: getLogPbol_RK4
446#endif
447 use pm_kind, only: RKG => RK4
448 real(RKG) , intent(in) :: logEpk, logPF53
449 real(RKG) :: logPbol
450 end function
451#endif
452
453#if RK3_ENABLED
454 pure elemental module function getLogPbol_RK3(logEpk, logPF53) result(logPbol)
455#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
456 !DEC$ ATTRIBUTES DLLEXPORT :: getLogPbol_RK3
457#endif
458 use pm_kind, only: RKG => RK3
459 real(RKG) , intent(in) :: logEpk, logPF53
460 real(RKG) :: logPbol
461 end function
462#endif
463
464#if RK2_ENABLED
465 pure elemental module function getLogPbol_RK2(logEpk, logPF53) result(logPbol)
466#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
467 !DEC$ ATTRIBUTES DLLEXPORT :: getLogPbol_RK2
468#endif
469 use pm_kind, only: RKG => RK2
470 real(RKG) , intent(in) :: logEpk, logPF53
471 real(RKG) :: logPbol
472 end function
473#endif
474
475#if RK1_ENABLED
476 pure elemental module function getLogPbol_RK1(logEpk, logPF53) result(logPbol)
477#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
478 !DEC$ ATTRIBUTES DLLEXPORT :: getLogPbol_RK1
479#endif
480 use pm_kind, only: RKG => RK1
481 real(RKG) , intent(in) :: logEpk, logPF53
482 real(RKG) :: logPbol
483 end function
484#endif
485
486 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
487
488 end interface getLogPbol
489
490!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
491
548 interface getLogPF53
549
550 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
551
552#if RK5_ENABLED
553 pure elemental module function getLogPF53_RK5(logEpk, logPbol) result(logPF53)
554#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
555 !DEC$ ATTRIBUTES DLLEXPORT :: getLogPF53_RK5
556#endif
557 use pm_kind, only: RKG => RK5
558 real(RKG) , intent(in) :: logEpk, logPbol
559 real(RKG) :: logPF53
560 end function
561#endif
562
563#if RK4_ENABLED
564 pure elemental module function getLogPF53_RK4(logEpk, logPbol) result(logPF53)
565#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
566 !DEC$ ATTRIBUTES DLLEXPORT :: getLogPF53_RK4
567#endif
568 use pm_kind, only: RKG => RK4
569 real(RKG) , intent(in) :: logEpk, logPbol
570 real(RKG) :: logPF53
571 end function
572#endif
573
574#if RK3_ENABLED
575 pure elemental module function getLogPF53_RK3(logEpk, logPbol) result(logPF53)
576#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
577 !DEC$ ATTRIBUTES DLLEXPORT :: getLogPF53_RK3
578#endif
579 use pm_kind, only: RKG => RK3
580 real(RKG) , intent(in) :: logEpk, logPbol
581 real(RKG) :: logPF53
582 end function
583#endif
584
585#if RK2_ENABLED
586 pure elemental module function getLogPF53_RK2(logEpk, logPbol) result(logPF53)
587#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
588 !DEC$ ATTRIBUTES DLLEXPORT :: getLogPF53_RK2
589#endif
590 use pm_kind, only: RKG => RK2
591 real(RKG) , intent(in) :: logEpk, logPbol
592 real(RKG) :: logPF53
593 end function
594#endif
595
596#if RK1_ENABLED
597 pure elemental module function getLogPF53_RK1(logEpk, logPbol) result(logPF53)
598#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
599 !DEC$ ATTRIBUTES DLLEXPORT :: getLogPF53_RK1
600#endif
601 use pm_kind, only: RKG => RK1
602 real(RKG) , intent(in) :: logEpk, logPbol
603 real(RKG) :: logPF53
604 end function
605#endif
606
607 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
608
609 end interface getLogPF53
610
611!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
612
660 interface getLog10PF53
661
662 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
663
664#if RK5_ENABLED
665 pure elemental module function getLog10PF53_RK5(log10Epk, log10Pbol) result(log10PF53)
666#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
667 !DEC$ ATTRIBUTES DLLEXPORT :: getLog10PF53_RK5
668#endif
669 use pm_kind, only: RKG => RK5
670 real(RKG) , intent(in) :: log10Epk, log10Pbol
671 real(RKG) :: log10PF53
672 end function
673#endif
674
675#if RK4_ENABLED
676 pure elemental module function getLog10PF53_RK4(log10Epk, log10Pbol) result(log10PF53)
677#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
678 !DEC$ ATTRIBUTES DLLEXPORT :: getLog10PF53_RK4
679#endif
680 use pm_kind, only: RKG => RK4
681 real(RKG) , intent(in) :: log10Epk, log10Pbol
682 real(RKG) :: log10PF53
683 end function
684#endif
685
686#if RK3_ENABLED
687 pure elemental module function getLog10PF53_RK3(log10Epk, log10Pbol) result(log10PF53)
688#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
689 !DEC$ ATTRIBUTES DLLEXPORT :: getLog10PF53_RK3
690#endif
691 use pm_kind, only: RKG => RK3
692 real(RKG) , intent(in) :: log10Epk, log10Pbol
693 real(RKG) :: log10PF53
694 end function
695#endif
696
697#if RK2_ENABLED
698 pure elemental module function getLog10PF53_RK2(log10Epk, log10Pbol) result(log10PF53)
699#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
700 !DEC$ ATTRIBUTES DLLEXPORT :: getLog10PF53_RK2
701#endif
702 use pm_kind, only: RKG => RK2
703 real(RKG) , intent(in) :: log10Epk, log10Pbol
704 real(RKG) :: log10PF53
705 end function
706#endif
707
708#if RK1_ENABLED
709 pure elemental module function getLog10PF53_RK1(log10Epk, log10Pbol) result(log10PF53)
710#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
711 !DEC$ ATTRIBUTES DLLEXPORT :: getLog10PF53_RK1
712#endif
713 use pm_kind, only: RKG => RK1
714 real(RKG) , intent(in) :: log10Epk, log10Pbol
715 real(RKG) :: log10PF53
716 end function
717#endif
718
719 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
720
721 end interface getLog10PF53
722
723!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
724
725contains
726
727!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
728
729! !> Return all log of data in natural (Neper) base.
730! !>
731! !> \param[in] inFilePath : The path to the input BATSE file.
732! !> \param[in] outFilePath : The path to the output BATSE file.
733! !> \param[in] isLgrb : A logical flag indicating what type of input file is being processed.
734! ! integer(IK) , parameter :: NLGRB = 1366_IK, NSGRB = 565_IK, NVAR = 4_IK
735!
736! ! GRB attributes
737! type :: Event_type
738! real(RK) :: logPbol, logEpk, logSbol, logT90, logPF53
739! end type Event_type
740!
741! type :: GRB_type
742! integer(IK) :: count
743! type(Event_type), allocatable :: Event(:)
744! !type(Event_type) :: Event(NLGRB)
745! end type GRB_type
746!
747! !> \cond excluded
748!#if CAF_ENABLED
749! type(GRB_type) :: GRB[*]
750!#else
751! type(GRB_type) :: GRB
752!#endif
753! !> \endcond excluded
754!
755! integer(IK), allocatable :: Trigger(:)
756! !integer(IK) :: Trigger(NLGRB)
757! !integer(IK) :: TriggerSGRB(NSGRB)
758!
759! subroutine readDataGRB(inFilePath,outFilePath,isLgrb)
760!#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
761! !DEC$ ATTRIBUTES DLLEXPORT :: readDataGRB
762!#endif
763!
764! use pm_parallelism, only: image_type
765! use pm_kind, only: IK, LK, RK
766! implicit none
767! character(*, SK), intent(in) :: inFilePath, outFilePath
768! integer(IK) :: inFileUnit, outFileUnit, igrb
769! logical(LK) , intent(in) :: isLgrb
770! type(image_type) :: image
771!
772! image = image_type()
773!
774! if (isLgrb) then
775! GRB%count = NLGRB
776! else
777! GRB%count = NSGRB
778! end if
779!
780! if (allocated(GRB%Event)) deallocate(GRB%Event); allocate(GRB%Event(GRB%count))
781! if (allocated(Trigger)) deallocate(Trigger); allocate(Trigger(GRB%count))
782!
783! open( newunit = inFileUnit & ! LCOV_EXCL_LINE
784! , file = inFilePath & ! LCOV_EXCL_LINE
785! , status = "old" & ! LCOV_EXCL_LINE
786!#if __INTEL_COMPILER && WINDOWS_ENABLED
787! , SHARED & ! LCOV_EXCL_LINE
788!#endif
789! )
790!
791! if (image%is%first) then
792! ! LCOV_EXCL_START
793! open(newunit=outFileUnit,file=outFilePath,status="replace")
794! write(outFileUnit,"(9a30)" ) "trigger" &
795! , "logPbol_1eV_20MeV" &
796! , "logSbol_1eV_20MeV" &
797! , "logEpk" &
798! , "logEPR1024" &
799! , "logEFR" &
800! , "logFPR1024" &
801! , "logT90" &
802! , "logEffPF53"
803! end if
804! ! LCOV_EXCL_STOP
805!
806! ! skip the header row in the input file
807! read(inFileUnit,*)
808!
809! ! read BATSE GRB data
810! do igrb = 1, GRB%count
811!
812! read(inFileUnit,* ) Trigger(igrb) & ! LCOV_EXCL_LINE
813! , GRB%Event(igrb)%logPF53 &
814! , GRB%Event(igrb)%logEpk &
815! , GRB%Event(igrb)%logSbol &
816! , GRB%Event(igrb)%logT90
817!
818! ! convert all values to logarithm in base Neper
819!
820! GRB%Event(igrb)%logPF53 = LOG_TEN * GRB%Event(igrb)%logPF53
821! GRB%Event(igrb)%logEpk = LOG_TEN * GRB%Event(igrb)%logEpk
822! GRB%Event(igrb)%logSbol = LOG_TEN * GRB%Event(igrb)%logSbol
823! GRB%Event(igrb)%logT90 = LOG_TEN * GRB%Event(igrb)%logT90
824!
825! ! convert photon count data to energy in units of ergs
826!
827! GRB%Event(igrb)%logPbol = getLogPbol( GRB%Event(igrb)%logEpk, GRB%Event(igrb)%logPF53 )
828! if (isLgrb) then
829! GRB%Event(igrb)%logSbol = getLogPbol( GRB%Event(igrb)%logEpk, GRB%Event(igrb)%logSbol )
830! else
831! GRB%Event(igrb)%logPF53 = GRB%Event(igrb)%logPF53 - THRESH_ERFC_AMP * erfc( (GRB%Event(igrb)%logT90-THRESH_ERFC_AVG) * THRESH_ERFC_STD_INV )
832! end if
833!
834! ! write the converted data to output file
835!
836! if (image%is%first) then
837! ! LCOV_EXCL_START
838! write(outFileUnit,"(I30,8E30.6)") Trigger(igrb) &
839! , GRB%Event(igrb)%logPbol &
840! , GRB%Event(igrb)%logSbol &
841! , GRB%Event(igrb)%logEpk &
842! , GRB%Event(igrb)%logEpk-GRB%Event(igrb)%logPbol &
843! , GRB%Event(igrb)%logEpk-GRB%Event(igrb)%logSbol &
844! , GRB%Event(igrb)%logSbol-GRB%Event(igrb)%logPbol &
845! , GRB%Event(igrb)%logT90 &
846! , GRB%Event(igrb)%logPF53
847!
848! end if
849! ! LCOV_EXCL_STOP
850!
851! end do
852!
853! if (image%is%first) close(outFileUnit)
854!
855! close(inFileUnit)
856!
857!!#if CAF_ENABLED
858!! sync images(*)
859!!
860!! else
861!!
862!! sync images(1)
863!! do igrb = 1, GRB%count
864!! GRB%Event(igrb)%logPbol = GRB[1]%Event(igrb)%logPbol
865!! GRB%Event(igrb)%logSbol = GRB[1]%Event(igrb)%logSbol
866!! GRB%Event(igrb)%logPF53 = GRB[1]%Event(igrb)%logPF53
867!! GRB%Event(igrb)%logEpk = GRB[1]%Event(igrb)%logEpk
868!! GRB%Event(igrb)%logT90 = GRB[1]%Event(igrb)%logT90
869!! end do
870!!
871!! end if
872!!#endif
873!
874! end subroutine readDataGRB
875!
876!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
877
878end module pm_batse
Generate and return the correction required for converting an input natural-log peak photon flux in 6...
Definition: pm_batse.F90:191
[LEGACY code] This generic interface is identical to the generic interface getLogPF53 with the only ...
Definition: pm_batse.F90:660
Generate and return the conversion of an input natural logarithm of peak photon flux in 64ms timescal...
Definition: pm_batse.F90:306
Generate and return the conversion of an input natural logarithm of a bolometric ( ) energy flux/flue...
Definition: pm_batse.F90:548
Generate and return the conversion of an input natural logarithm of photon flux/fluence over a given ...
Definition: pm_batse.F90:427
This module contains procedures and generic interfaces for modeling data and detectors of the BATSE G...
Definition: pm_batse.F90:69
real(RKB), parameter MAX_LOG10PH53_4_LOGPBOLZERO
Definition: pm_batse.F90:79
real(RKB), parameter DIF_LOGPH53_4_LOGPBOLZERO
Definition: pm_batse.F90:82
real(RKB), parameter THRESH_ERFC_AMP
The scale of the change in BATSE efficiency for different GRB durations.
Definition: pm_batse.F90:105
real(RKB), parameter THRESH_ERFC_STD
Scale of the duration in the Error function used to model the connection between the peak fluxes in 6...
Definition: pm_batse.F90:111
real(RKB), parameter THRESH_LOGPBOL64_CORRECTION
Correction that must be added to logPbol64ms to convert it to effective peak flux.
Definition: pm_batse.F90:121
real(RKB), parameter THRESH_ERFC_AVG
Mean duration in the Error function used to model the connection between the peak fluxes in 64 and 10...
Definition: pm_batse.F90:108
real(RKB), parameter LOG_TEN
Definition: pm_batse.F90:77
real(RKB), parameter MAX_LOGPH53_4_LOGPBOLZERO
Definition: pm_batse.F90:80
real(RKB), parameter MIN_LOG10PH53_4_LOGPBOLZERO
Definition: pm_batse.F90:78
real(RKB), parameter MIN_LOGPH53_4_LOGPBOLZERO
Definition: pm_batse.F90:81
real(RKB), parameter LOG10PF53_MINUS_LOG10PBOL
Definition: pm_batse.F90:99
real(RKB), parameter THRESH_ERFC_STD_INV
Inverse scale of the duration in the Error function used to model the connection between the peak flu...
Definition: pm_batse.F90:114
real(RKB), parameter THRESH_ERFC_HEIGHT
Definition: pm_batse.F90:117
real(RKB), parameter LOGPF53_MINUS_LOGPBOL
Definition: pm_batse.F90:98
real(RKB), parameter THRESH_ERFC_BASE
The scale of the change in BATSE efficiency for different SGRB durations.
Definition: pm_batse.F90:102
character(*, SK), parameter MODULE_NAME
Definition: pm_batse.F90:75
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 RKB
The scalar integer constant of intrinsic default kind, representing the Best-precision real kind supp...
Definition: pm_kind.F90:1371
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 RKD
The double precision real kind in Fortran mode. On most platforms, this is an 64-bit real kind.
Definition: pm_kind.F90:568
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 is the derived type for generating objects that contain attributes of BATSE catalog GRBs.
Definition: pm_batse.F90:135