Module helixfit_fitpar_mod

module helixfit_fitpar_mod

        ! Uses
    use precision_mod
    use namelist_mod
    use physicalconstants_mod
    use track_mod
    use helixfit_common_mod

        ! Types
    public type fitpar_type

        ! Variables
    integer, public, PARAMETER :: kParTypeNone = 0
    integer, public, PARAMETER :: kParTypeQ0 = 1
    integer, public, PARAMETER :: kParTypeABCFP = 2
    integer, public, PARAMETER :: kParTypeCRFP = 3
    integer, public, PARAMETER :: kParTypePtAf = 4
    integer, public, PARAMETER :: kParTypePT = 5
    integer, public, PARAMETER :: kMaxKinks = 20
    integer, public, PARAMETER :: kMaxPar = 5 + 2 * kMaxKinks + 1
    logical, public, PARAMETER :: HelixFitParTref = .true.
    real (kind=R8), private, dimension (5) :: parstepWC_ABCFP = (/ 0.4000, 0.4000, 1.6000, 0.01, 0.4 /)
    real (kind=R8), private, dimension (5) :: parstepDT_ABCFP = (/ 0.0050, 0.0050, 0.0050, 0.0005, 0.005 /)
    real (kind=R8), private, dimension (5) :: parstepWC_CRFP = (/ 0.4000, 0.4000, 0.4000, 0.01, 0.4 /)
    real (kind=R8), private, dimension (5) :: parstepDT_CRFP = (/ 0.2000, 0.2000, 0.2000, 0.0025, 0.05 /)
    real (kind=R8), private, dimension (5) :: parstepWC_PTAF = (/ 0.0000, 0.0000, 0.000, 0.0000, 0.0000 /)
    real (kind=R8), private, dimension (5) :: parstepWC_PT = (/ 0.0000, 0.0000, 0.000, 0.0000, 0.0000 /)
    real (kind=R8), private, dimension (5) :: parstepDT_PTAF = (/ 0.0100, 0.0100, 0.050, 0.0005, 0.1000 /)
    real (kind=R8), private, dimension (5) :: parstepDT_PT = (/ 0.0100, 0.0100, 0.100, 0.0020, 0.0500 /)
    real (kind=R8), public, dimension (kMaxPar) :: parstepWC
    real (kind=R8), public, dimension (kMaxPar) :: parstepDT

        ! Subroutines and functions
    public subroutine clearFitpar (fitpar)
    public function setFitpar (fitpar, partype, V, P, Q, mass, Tref, Zmean) result (success)
    public subroutine scaleFitpar (fitpar, offset, scale)
    public subroutine unscaleFitpar (fitpar)
    public subroutine getFitparVPQ (fitpar, V, P, Q)
    public function tryGetFitparVPQ (fitpar, V, P, Q) result (ierror)
    public subroutine updateFitpar (fitpar, spar)
    public subroutine addFitparTref (fitpar)
    public subroutine copyFitpar (fitpar, partype, old)
    public subroutine VPQtoPPCP (v, p, q, pmag, pt, cost, phi)
    public subroutine getFitparPPCP (fitpar, pmag, pt, cost, phi)
    private function jsign (a) result (ival)
    private function v3mag (a) result (s)
    private function sqr (a) result (b)

end module helixfit_fitpar_mod

Description of Types

fitpar_type

public type fitpar_type
    integer :: partype
            
 Helix fit parameters. See setFitpar() for definitions
 of different parametrizations

 parametrization type kParTypeXXX
    integer :: numpar
             number of parameters
    real (kind=R8), dimension (kMaxPar) :: spar
             (scaled) parameters
    real (kind=R8), dimension (kMaxPar) :: parOffset
            
 Parameter scaling:
   realPar   = (scaledPar * parScale) + parOffset
   scaledPar = (realPar - parOffset) / parScale

 parameter offset
    real (kind=R8), dimension (kMaxPar) :: parScale
             parameter scale
    integer :: numkinks
            
 Kink information

 number of kinks
    integer :: trefoffset
            
 Tref fit information

 offset of tref in the par() array
    integer :: Q
            
 Helix constant and cached parameters-
 kept consistent with par() by the
 getpar and setpar subroutines

 charge
    real (kind=R8) :: mass
             mass in MeV/c^2
    real (kind=R8) :: z
             z coordinate
    real (kind=R8) :: zmean
             mean z coordinate
    real (kind=R8) :: cost
             cos(angle between P and Z)
    real (kind=R8) :: tref
             time reference (decay time)
    real (kind=R8) :: elossHe
             energy loss
    real (kind=R8) :: elossDC
             energy loss
    real (kind=R8) :: elossFoil
             energy loss
    real (kind=R8), dimension (3) :: none_v
             used only by kParTypeNone
    real (kind=R8), dimension (3) :: none_p
             used only by kParTypeNone
    real (kind=R8) :: fldoffpmag
             used only by kParTypeQ0
end type fitpar_type

Description of Variables

kParTypeNone

integer, public, PARAMETER :: kParTypeNone = 0
 Fit parametrizations

kParTypeQ0

integer, public, PARAMETER :: kParTypeQ0 = 1
 B=0, Q=0, straight track params

kParTypeABCFP

integer, public, PARAMETER :: kParTypeABCFP = 2
 fit circle (A-B-C param) and angular frequency

kParTypeCRFP

integer, public, PARAMETER :: kParTypeCRFP = 3
 fit center, radius and angular frequency

kParTypePtAf

integer, public, PARAMETER :: kParTypePtAf = 4
 fit pt and angular frequency

kParTypePT

integer, public, PARAMETER :: kParTypePT = 5
 fit Pmag, tan(theta), starting point, phase

kMaxKinks

integer, public, PARAMETER :: kMaxKinks = 20
 maximum number of kinks

kMaxPar

integer, public, PARAMETER :: kMaxPar = 5 + 2 * kMaxKinks + 1
 maximum number of fit parameters: 5 + 2*kinks + tref

HelixFitParTref

logical, public, PARAMETER :: HelixFitParTref = .true.
 fit tref using main fitter

parstepWC_ABCFP

real (kind=R8), private, dimension (5) :: parstepWC_ABCFP = (/ 0.4000, 0.4000, 1.6000, 0.01, 0.4 /)
 step in A
 step in B
 step in C
 step in angular frequency
 step in phase

parstepDT_ABCFP

real (kind=R8), private, dimension (5) :: parstepDT_ABCFP = (/ 0.0050, 0.0050, 0.0050, 0.0005, 0.005 /)
 step in A
 step in B
 step in C
 step in angular frequency
 step in phase

parstepWC_CRFP

real (kind=R8), private, dimension (5) :: parstepWC_CRFP = (/ 0.4000, 0.4000, 0.4000, 0.01, 0.4 /)
 step in U
 step in V
 step in R
 step in angular frequency
 step in phase

parstepDT_CRFP

real (kind=R8), private, dimension (5) :: parstepDT_CRFP = (/ 0.2000, 0.2000, 0.2000, 0.0025, 0.05 /)
 step in U
 step in V
 step in R
 step in angular frequency
 step in phase

parstepWC_PTAF

real (kind=R8), private, dimension (5) :: parstepWC_PTAF = (/ 0.0000, 0.0000, 0.000, 0.0000, 0.0000 /)
 step in U
 step in V
 step in Pt, MeV/c
 step in angular frequency
 step in phase

parstepWC_PT

real (kind=R8), private, dimension (5) :: parstepWC_PT = (/ 0.0000, 0.0000, 0.000, 0.0000, 0.0000 /)
 step in U
 step in V
 step in Pmag, MeV/c
 step in tan(theta)
 step in phase

parstepDT_PTAF

real (kind=R8), private, dimension (5) :: parstepDT_PTAF = (/ 0.0100, 0.0100, 0.050, 0.0005, 0.1000 /)
 step in U
 step in V
 step in Pt, MeV/c
 step in angular frequency
 step in phase

parstepDT_PT

real (kind=R8), private, dimension (5) :: parstepDT_PT = (/ 0.0100, 0.0100, 0.100, 0.0020, 0.0500 /)
 step in U
 step in V
 step in Pmag, MeV/c
 step in tan(theta)
 step in phase

parstepWC

real (kind=R8), public, dimension (kMaxPar) :: parstepWC

parstepDT

real (kind=R8), public, dimension (kMaxPar) :: parstepDT

Description of Subroutines and Functions

clearFitpar

public subroutine clearFitpar (fitpar)
    type (fitpar_type), INTENT(out) :: fitpar
end subroutine clearFitpar

setFitpar

public function setFitpar (fitpar, partype, V, P, Q, mass, Tref, Zmean) result (success)
    type (fitpar_type), INTENT(out) :: fitpar
            
 Description: initialize fitpar_type

 Arguments:
    integer, INTENT(in) :: partype
             paramertization type
    real (kind=R8), INTENT(in), dimension (3) :: v
             track starting point
    real (kind=R8), INTENT(in), dimension (3) :: p
             track momentum
    integer, INTENT(in) :: q
             charge
    real (kind=R8), INTENT(in) :: mass
             mass in MeV/c^2
    real (kind=R8), INTENT(in) :: Tref
             time reference
    real (kind=R8), INTENT(in) :: Zmean
             mean z coordinate
    logical :: success
    ! Calls: die
end function setFitpar

scaleFitpar

public subroutine scaleFitpar (fitpar, offset, scale)
    type (fitpar_type), INTENT(inout) :: fitpar
            
 Description: scale this parametrization
    real (kind=R8), INTENT(in), dimension (:) :: offset
    real (kind=R8), INTENT(in), dimension (:) :: scale
end subroutine scaleFitpar

unscaleFitpar

public subroutine unscaleFitpar (fitpar)
    type (fitpar_type), INTENT(inout) :: fitpar
            
 Description: undo scaling of parametrization
    ! Calls: scaleFitpar
end subroutine unscaleFitpar

getFitparVPQ

public subroutine getFitparVPQ (fitpar, V, P, Q)
    type (fitpar_type), INTENT(in) :: fitpar
            
 Extract fit parameters, see tryGetFitparVPQ()
    real (kind=R8), INTENT(out), dimension (3) :: V
    real (kind=R8), INTENT(out), dimension (3) :: P
    integer, INTENT(out) :: Q
    ! Calls: die
end subroutine getFitparVPQ

tryGetFitparVPQ

public function tryGetFitparVPQ (fitpar, V, P, Q) result (ierror)
    type (fitpar_type), INTENT(in) :: fitpar
            
 Description: extract V,P and Q from fitpar_type

 Arguments:
    real (kind=R8), INTENT(out), dimension (3) :: V
    real (kind=R8), INTENT(out), dimension (3) :: P
    integer, INTENT(out) :: Q
    integer :: ierror
             return value
    ! Calls: die
end function tryGetFitparVPQ

updateFitpar

public subroutine updateFitpar (fitpar, spar)
    type (fitpar_type), INTENT(inout) :: fitpar
            
 Description: update fit parameters

 Arguments:
    real (kind=R8), INTENT(in), dimension (kMaxPar) :: spar
    ! Calls: die
end subroutine updateFitpar

addFitparTref

public subroutine addFitparTref (fitpar)
    type (fitpar_type), INTENT(inout) :: fitpar
            
 Description: add tref to list of fit parameters

 Arguments:
    ! Calls: die
end subroutine addFitparTref

copyFitpar

public subroutine copyFitpar (fitpar, partype, old)
    type (fitpar_type), INTENT(out) :: fitpar
            
 Description: copy "old" to "fitpar"

 Arguments:
    integer, INTENT(in) :: partype
    type (fitpar_type), INTENT(in) :: old
    ! Calls: GetFitparVPQ, clearFitpar
end subroutine copyFitpar

VPQtoPPCP

public subroutine VPQtoPPCP (v, p, q, pmag, pt, cost, phi)
    real (kind=R8), INTENT(in), dimension (3) :: v
            
 Description: extract P, Pt, cos(theta) and phi from v, p and q

 Arguments:

 position
    real (kind=R8), INTENT(in), dimension (3) :: p
             momentum
    integer, INTENT(in) :: q
             charge
    real (kind=R8), INTENT(out) :: pmag
    real (kind=R8), INTENT(out) :: pt
    real (kind=R8), INTENT(out) :: cost
    real (kind=R8), INTENT(out) :: phi
end subroutine VPQtoPPCP

getFitparPPCP

public subroutine getFitparPPCP (fitpar, pmag, pt, cost, phi)
    type (fitpar_type), INTENT(in) :: fitpar
            
 Description: extract P, Pt, cos(theta) and phi from fitpar_type

 Arguments:
    real (kind=R8), INTENT(out) :: pmag
    real (kind=R8), INTENT(out) :: pt
    real (kind=R8), INTENT(out) :: cost
    real (kind=R8), INTENT(out) :: phi
    ! Calls: VPQtoPPCP, getFitparVPQ
end subroutine getFitparPPCP

jsign

private function jsign (a) result (ival)
    real (kind=R8), INTENT(in) :: a
            
 Description: return sign of a real(r8) number
 Returns: +1 if positive
          -1 if negative
           0 if zero
    integer :: ival
end function jsign

v3mag

private function v3mag (a) result (s)
    real (kind=R8), INTENT(in), dimension (3) :: a
            
 Description: compute vector length
    real (kind=R8) :: s
end function v3mag

sqr

private function sqr (a) result (b)
    real (kind=R8), INTENT(in) :: a
             Description: return square of a real(r8) number
    real (kind=R8) :: b
end function sqr