Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
H
hcana
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Snippets
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Container Registry
Model registry
Operate
Environments
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
jlab
hallc
analyzer_software
hcana
Commits
1dce3be1
Commit
1dce3be1
authored
7 years ago
by
bduran
Committed by
Mark K Jones
7 years ago
Browse files
Options
Downloads
Patches
Plain Diff
added new variables for the 2nd raster, edited the code to be able to read both pairs of rasters
parent
239cd15b
No related branches found
No related tags found
No related merge requests found
Changes
4
Hide whitespace changes
Inline
Side-by-side
Showing
4 changed files
src/THcRaster.cxx
+189
-94
189 additions, 94 deletions
src/THcRaster.cxx
src/THcRaster.h
+49
-19
49 additions, 19 deletions
src/THcRaster.h
src/THcRasterRawHit.cxx
+1
-68
1 addition, 68 deletions
src/THcRasterRawHit.cxx
src/THcRasterRawHit.h
+5
-21
5 additions, 21 deletions
src/THcRasterRawHit.h
with
244 additions
and
202 deletions
src/THcRaster.cxx
+
189
−
94
View file @
1dce3be1
/** \class THcRaster
/** \class THcRaster
\ingroup DetSupport
\ingroup DetSupport
A class to decode the fast raster signals.
A class to
the fast raster signals.
Measures the two magnet currents which are proportional to horizontal and
Measures the two magnet currents which are proportional to horizontal and
vertical beam position
vertical beam position
\author Buddhini Waidyawansa
\author Buddhini Waidyawansa
\author Burcu Duran - Melanie Rehfuss (2017)
*/
*/
#include
"TMath.h"
#include
"TMath.h"
...
@@ -18,6 +20,9 @@
...
@@ -18,6 +20,9 @@
#include
"THcParmList.h"
#include
"THcParmList.h"
#include
"THcGlobals.h"
#include
"THcGlobals.h"
#include
"THaGlobals.h"
#include
"THaGlobals.h"
#include
"THaApparatus.h"
#include
"THcRawAdcHit.h"
#include
"THcSignalHit.h"
//#include "THcHitList.h"
//#include "THcHitList.h"
...
@@ -38,37 +43,95 @@ THcRaster::THcRaster( const char* name, const char* description,
...
@@ -38,37 +43,95 @@ THcRaster::THcRaster( const char* name, const char* description,
fAnalyzePedestals
=
0
;
fAnalyzePedestals
=
0
;
fNPedestalEvents
=
0
;
fNPedestalEvents
=
0
;
fRawXADC
=
0
;
FRXA_rawadc
=
0
;
fRawYADC
=
0
;
FRXB_rawadc
=
0
;
fXADC
=
0
;
FRYA_rawadc
=
0
;
fYADC
=
0
;
FRYB_rawadc
=
0
;
fXpos
=
0
;
fXA_ADC
=
0
;
fYpos
=
0
;
fYA_ADC
=
0
;
fXB_ADC
=
0
;
fYB_ADC
=
0
;
fXA_pos
=
0
;
fYA_pos
=
0
;
fXB_pos
=
0
;
fYB_pos
=
0
;
fFrCalMom
=
0
;
fFrCalMom
=
0
;
fFrXADCperCM
=
0
;
fFrXA_ADCperCM
=
1.0
;
fFrXADCperCM
=
0
;
fFrYA_ADCperCM
=
1.0
;
fFrXB_ADCperCM
=
1.0
;
fFrYB_ADCperCM
=
1.0
;
fFrXA_ADC_zero_offset
=
0
;
fFrXA_ADC_zero_offset
=
0
;
fFrXA_ADC_zero_offset
=
0
;
fFrXA_ADC_zero_offset
=
0
;
for
(
Int_t
i
=
0
;
i
<
2
;
i
++
){
frPosAdcPulseIntRaw
=
NULL
;
fPedADC
[
i
]
=
0
;
fAvgPedADC
[
i
]
=
0
;
for
(
Int_t
i
=
0
;
i
<
4
;
i
++
){
}
fPedADC
[
i
]
=
0
;
}
InitArrays
();
}
}
//_____________________________________________________________________________
THcRaster
::
THcRaster
(
)
:
THaBeamDet
(
"THcRaster"
)
// no default constructor available
{
frPosAdcPulseIntRaw
=
NULL
;
InitArrays
();
}
//_____________________________________________________________________________
//_____________________________________________________________________________
THcRaster
::~
THcRaster
()
THcRaster
::~
THcRaster
()
{
{
// Destructor
delete
frPosAdcPulseIntRaw
;
frPosAdcPulseIntRaw
=
NULL
;
DeleteArrays
();
}
}
//_____________________________________________________________________________
THaAnalysisObject
::
EStatus
THcRaster
::
Init
(
const
TDatime
&
date
)
{
//cout << "THcRaster::Init()" << endl;
char
EngineDID
[]
=
"xRASTER"
;
EngineDID
[
0
]
=
toupper
(
GetApparatus
()
->
GetName
()[
0
]);
// Fill detector map with RASTER type channels
if
(
gHcDetectorMap
->
FillMap
(
fDetMap
,
EngineDID
)
<
0
)
{
static
const
char
*
const
here
=
"Init()"
;
Error
(
Here
(
here
),
"Error filling detectormap for %s."
,
EngineDID
);
return
kInitError
;
}
InitHitList
(
fDetMap
,
"THcRasterRawHit"
,
fDetMap
->
GetTotNumChan
()
+
1
);
EStatus
status
;
if
(
(
status
=
THaBeamDet
::
Init
(
date
))
)
return
fStatus
=
status
;
return
fStatus
=
kOK
;
}
//_____________________________________________________________________________
//_____________________________________________________________________________
Int_t
THcRaster
::
ReadDatabase
(
const
TDatime
&
date
)
Int_t
THcRaster
::
ReadDatabase
(
const
TDatime
&
date
)
{
{
// Read parameters such as calibration factor, of this detector from the database.
// Read parameters such as calibration factor, of this detector from the database.
cout
<<
"THcRaster::ReadDatabase()"
<<
endl
;
//
cout << "THcRaster::ReadDatabase()" << endl;
char
prefix
[
2
];
char
prefix
[
2
];
...
@@ -78,12 +141,19 @@ Int_t THcRaster::ReadDatabase( const TDatime& date )
...
@@ -78,12 +141,19 @@ Int_t THcRaster::ReadDatabase( const TDatime& date )
prefix
[
0
]
=
'g'
;
prefix
[
0
]
=
'g'
;
prefix
[
1
]
=
'\0'
;
prefix
[
1
]
=
'\0'
;
// string names;
// string names;
DBRequest
list
[]
=
{
DBRequest
list
[]
=
{
{
"fr_cal_mom"
,
&
fFrCalMom
,
kDouble
},
{
"fr_cal_mom"
,
&
fFrCalMom
,
kDouble
},
{
"frx_adcpercm"
,
&
fFrXADCperCM
,
kDouble
},
{
"frxa_adcpercm"
,
&
fFrXA_ADCperCM
,
kDouble
},
{
"fry_adcpercm"
,
&
fFrYADCperCM
,
kDouble
},
{
"frya_adcpercm"
,
&
fFrYA_ADCperCM
,
kDouble
},
{
"frxb_adcpercm"
,
&
fFrXB_ADCperCM
,
kDouble
},
{
"fryb_adcpercm"
,
&
fFrYB_ADCperCM
,
kDouble
},
{
"frxa_adc_zero_offset"
,
&
fFrXA_ADC_zero_offset
,
kDouble
},
{
"frya_adc_zero_offset"
,
&
fFrYA_ADC_zero_offset
,
kDouble
},
{
"frxb_adc_zero_offset"
,
&
fFrXB_ADC_zero_offset
,
kDouble
},
{
"fryb_adc_zero_offset"
,
&
fFrYB_ADC_zero_offset
,
kDouble
},
{
"pbeam"
,
&
fgpbeam
,
kDouble
},
{
"pbeam"
,
&
fgpbeam
,
kDouble
},
{
"frx_dist"
,
&
fgfrx_dist
,
kDouble
},
{
"frx_dist"
,
&
fgfrx_dist
,
kDouble
},
{
"fry_dist"
,
&
fgfry_dist
,
kDouble
},
{
"fry_dist"
,
&
fgfry_dist
,
kDouble
},
...
@@ -101,9 +171,13 @@ Int_t THcRaster::ReadDatabase( const TDatime& date )
...
@@ -101,9 +171,13 @@ Int_t THcRaster::ReadDatabase( const TDatime& date )
fgbeam_yoff
=
0.0
;
fgbeam_yoff
=
0.0
;
fgbeam_ypoff
=
0.0
;
fgbeam_ypoff
=
0.0
;
fgusefr
=
0
;
fgusefr
=
0
;
// get the calibration factors from gbeam.param file
// get the calibration factors from gbeam.param file
gHcParms
->
LoadParmValues
((
DBRequest
*
)
&
list
,
prefix
);
gHcParms
->
LoadParmValues
((
DBRequest
*
)
&
list
,
prefix
);
frPosAdcPulseIntRaw
=
new
TClonesArray
(
"THcSignalHit"
,
4
);
return
kOK
;
return
kOK
;
}
}
...
@@ -114,7 +188,7 @@ Int_t THcRaster::DefineVariables( EMode mode )
...
@@ -114,7 +188,7 @@ Int_t THcRaster::DefineVariables( EMode mode )
{
{
// Initialize global variables for histogramming and tree
// Initialize global variables for histogramming and tree
cout
<<
"THcRaster::DefineVariables called "
<<
GetName
()
<<
endl
;
//
cout << "THcRaster::DefineVariables called " << GetName() << endl;
if
(
mode
==
kDefine
&&
fIsSetup
)
return
kOK
;
if
(
mode
==
kDefine
&&
fIsSetup
)
return
kOK
;
fIsSetup
=
(
mode
==
kDefine
);
fIsSetup
=
(
mode
==
kDefine
);
...
@@ -122,12 +196,18 @@ Int_t THcRaster::DefineVariables( EMode mode )
...
@@ -122,12 +196,18 @@ Int_t THcRaster::DefineVariables( EMode mode )
// Register variables in global list
// Register variables in global list
RVarDef
vars
[]
=
{
RVarDef
vars
[]
=
{
{
"frx_raw_adc"
,
"Raster X raw ADC"
,
"fRawXADC"
},
{
"frxaRawAdc"
,
"Raster XA raw ADC"
,
"FRXA_rawadc"
},
{
"fry_raw_adc"
,
"Raster Y raw ADC"
,
"fRawYADC"
},
{
"fryaRawAdc"
,
"Raster YA raw ADC"
,
"FRXB_rawadc"
},
{
"frx_adc"
,
"Raster X ADC"
,
"fXADC"
},
{
"frxbRawAdc"
,
"Raster XB raw ADC"
,
"FRYA_rawadc"
},
{
"fry_adc"
,
"Raster Y ADC"
,
"fYADC"
},
{
"frybRawAdc"
,
"Raster YB raw ADC"
,
"FRYB_rawadc"
},
{
"frx"
,
"Raster X position"
,
"fXpos"
},
{
"frxa_adc"
,
"Raster XA ADC"
,
"fXA_ADC"
},
{
"fry"
,
"Raster Y position"
,
"fYpos"
},
{
"frya_adc"
,
"Raster YA ADC"
,
"fYA_ADC"
},
{
"frxb_adc"
,
"Raster XB ADC"
,
"fXB_ADC"
},
{
"fryb_adc"
,
"Raster YB ADC"
,
"fYB_ADC"
},
{
"fr_xa"
,
"Raster XA position"
,
"fXA_pos"
},
{
"fr_ya"
,
"Raster YA position"
,
"fYA_pos"
},
{
"fr_xb"
,
"Raster XB position"
,
"fXB_pos"
},
{
"fr_yb"
,
"Raster YB position"
,
"fYB_pos"
},
{
0
}
{
0
}
};
};
...
@@ -135,29 +215,15 @@ Int_t THcRaster::DefineVariables( EMode mode )
...
@@ -135,29 +215,15 @@ Int_t THcRaster::DefineVariables( EMode mode )
}
}
//_____________________________________________________________________________
//_____________________________________________________________________________
THaAnalysisObject
::
EStatus
THcRaster
::
Init
(
const
TDatime
&
date
)
{
cout
<<
"THcRaster::Init()"
<<
endl
;
// Fill detector map with RASTER type channels
if
(
gHcDetectorMap
->
FillMap
(
fDetMap
,
"RASTER"
)
<
0
)
{
static
const
char
*
const
here
=
"Init()"
;
Error
(
Here
(
here
),
"Error filling detectormap for %s."
,
"RASTER"
);
return
kInitError
;
}
THcHitList
::
InitHitList
(
fDetMap
,
"THcRasterRawHit"
,
fDetMap
->
GetTotNumChan
()
+
1
);
inline
void
THcRaster
::
Clear
(
Option_t
*
opt
)
EStatus
status
;
{
if
(
(
status
=
THaBeamDet
::
Init
(
date
))
)
return
fStatus
=
status
;
return
fStatu
s
=
kOK
;
fNhit
s
=
0
;
frPosAdcPulseIntRaw
->
Clear
();
}
}
//_____________________________________________________________________________
//_____________________________________________________________________________
void
THcRaster
::
AccumulatePedestals
(
TClonesArray
*
rawhits
)
void
THcRaster
::
AccumulatePedestals
(
TClonesArray
*
rawhits
)
{
{
...
@@ -182,21 +248,32 @@ void THcRaster::AccumulatePedestals(TClonesArray* rawhits)
...
@@ -182,21 +248,32 @@ void THcRaster::AccumulatePedestals(TClonesArray* rawhits)
Int_t
nrawhits
=
rawhits
->
GetLast
()
+
1
;
Int_t
nrawhits
=
rawhits
->
GetLast
()
+
1
;
Int_t
ihit
=
0
;
Int_t
ihit
=
0
;
UInt_t
nrPosAdcHits
=
0
;
while
(
ihit
<
nrawhits
)
{
THcRasterRawHit
*
hit
=
(
THcRasterRawHit
*
)
fRawHitList
->
At
(
ihit
);
THcRawAdcHit
&
rawPosAdcHit
=
hit
->
GetRawAdcHitPos
();
Int_t
nsig
=
hit
->
fCounter
;
for
(
UInt_t
thit
=
0
;
thit
<
rawPosAdcHit
.
GetNPulses
();
++
thit
)
{
((
THcSignalHit
*
)
frPosAdcPulseIntRaw
->
ConstructedAt
(
nrPosAdcHits
))
->
Set
(
nsig
,
rawPosAdcHit
.
GetPulseIntRaw
(
thit
));
++
nrPosAdcHits
;
}
ihit
++
;
}
while
(
ihit
<
nrawhits
)
{
for
(
Int_t
ielem
=
0
;
ielem
<
frPosAdcPulseIntRaw
->
GetEntries
();
ielem
++
)
{
THcRasterRawHit
*
hit
=
(
THcRasterRawHit
*
)
fRawHitList
->
At
(
ihit
);
Int_t
nraster
=
((
THcSignalHit
*
)
frPosAdcPulseIntRaw
->
ConstructedAt
(
ielem
))
->
GetPaddleNumber
()
-
1
;
if
(
hit
->
fADC_xsig
>
0
)
{
Double_t
pulseIntRaw
=
((
THcSignalHit
*
)
frPosAdcPulseIntRaw
->
ConstructedAt
(
ielem
))
->
GetData
();
fPedADC
[
0
]
+=
hit
->
fADC_xsig
;
if
(
nraster
==
0
)
fPedADC
[
0
]
=
pulseIntRaw
;
//std::cout<<" raster x pedestal collect "<<fPedADC[0]<<std::endl;
if
(
nraster
==
1
)
fPedADC
[
1
]
=
pulseIntRaw
;
}
if
(
nraster
==
2
)
fPedADC
[
2
]
=
pulseIntRaw
;
if
(
hit
->
fADC_ysig
>
0
)
{
if
(
nraster
==
3
)
fPedADC
[
3
]
=
pulseIntRaw
;
fPedADC
[
1
]
+=
hit
->
fADC_ysig
;
}
//std::cout<<" raster y pedestal collect "<<fPedADC[1]<<std::endl;
}
ihit
++
;
}
}
}
...
@@ -221,11 +298,13 @@ void THcRaster::CalculatePedestals( )
...
@@ -221,11 +298,13 @@ void THcRaster::CalculatePedestals( )
endif
endif
endif
endif
*/
*/
for
(
Int_t
i
=
0
;
i
<
2
;
i
++
){
fAvgPedADC
[
i
]
=
fPedADC
[
i
]
/
fNPedestalEvents
;
fFrXA_ADC_zero_offset
=
fPedADC
[
0
]
/
fNPedestalEvents
;
// std::cout<<" raster pedestal "<<fAvgPedADC[i]<<std::endl;
fFrYA_ADC_zero_offset
=
fPedADC
[
1
]
/
fNPedestalEvents
;
}
fFrXB_ADC_zero_offset
=
fPedADC
[
2
]
/
fNPedestalEvents
;
fFrYB_ADC_zero_offset
=
fPedADC
[
3
]
/
fNPedestalEvents
;
}
}
...
@@ -233,44 +312,56 @@ void THcRaster::CalculatePedestals( )
...
@@ -233,44 +312,56 @@ void THcRaster::CalculatePedestals( )
Int_t
THcRaster
::
Decode
(
const
THaEvData
&
evdata
)
Int_t
THcRaster
::
Decode
(
const
THaEvData
&
evdata
)
{
{
// loops over all channels defined in the detector map
// copies raw data into local variables
// Get the Hall C style hitlist (fRawHitList) for this event
// performs pedestal subtraction
fNhits
=
DecodeToHitList
(
evdata
);
// Get the Hall C style hitlist (fRawHitList) for this event
Int_t
fNhits
=
THcHitList
::
DecodeToHitList
(
evdata
);
// Get the pedestals from the first 1000 events
// Get the pedestals from the first 1000 events
//if(fNPedestalEvents < 10)
//if(fNPedestalEvents < 10)
if
(
(
gHaCuts
->
Result
(
"Pedestal_event"
)
)
&
(
fNPedestalEvents
<
1000
)){
if
(
gHaCuts
->
Result
(
"Pedestal_event"
)
&
(
fNPedestalEvents
<
1000
))
{
AccumulatePedestals
(
fRawHitList
);
AccumulatePedestals
(
fRawHitList
);
fAnalyzePedestals
=
1
;
// Analyze pedestals first normal events
fAnalyzePedestals
=
1
;
// Analyze pedestals first normal events
fNPedestalEvents
++
;
fNPedestalEvents
++
;
return
(
0
);
return
(
0
);
}
}
if
(
fAnalyzePedestals
)
{
CalculatePedestals
();
if
(
fAnalyzePedestals
)
{
fAnalyzePedestals
=
0
;
// Don't analyze pedestals next event
CalculatePedestals
();
fAnalyzePedestals
=
0
;
// Don't analyze pedestals next event
}
Int_t
ihit
=
0
;
UInt_t
nrPosAdcHits
=
0
;
while
(
ihit
<
fNhits
)
{
THcRasterRawHit
*
hit
=
(
THcRasterRawHit
*
)
fRawHitList
->
At
(
ihit
);
THcRawAdcHit
&
rawPosAdcHit
=
hit
->
GetRawAdcHitPos
();
Int_t
nsig
=
hit
->
fCounter
;
for
(
UInt_t
thit
=
0
;
thit
<
rawPosAdcHit
.
GetNPulses
();
++
thit
)
{
((
THcSignalHit
*
)
frPosAdcPulseIntRaw
->
ConstructedAt
(
nrPosAdcHits
))
->
Set
(
nsig
,
rawPosAdcHit
.
GetPulseIntRaw
(
thit
));
++
nrPosAdcHits
;
}
ihit
++
;
}
}
Int_t
ihit
=
0
;
while
(
ihit
<
fNhits
)
{
for
(
Int_t
ielem
=
0
;
ielem
<
frPosAdcPulseIntRaw
->
GetEntries
();
ielem
++
)
{
THcRasterRawHit
*
hit
=
(
THcRasterRawHit
*
)
fRawHitList
->
At
(
ihit
);
Int_t
nraster
=
((
THcSignalHit
*
)
frPosAdcPulseIntRaw
->
ConstructedAt
(
ielem
))
->
GetPaddleNumber
()
-
1
;
Double_t
pulseIntRaw
=
((
THcSignalHit
*
)
frPosAdcPulseIntRaw
->
ConstructedAt
(
ielem
))
->
GetData
();
if
(
nraster
==
0
)
FRXA_rawadc
=
pulseIntRaw
;
if
(
nraster
==
1
)
FRXB_rawadc
=
pulseIntRaw
;
if
(
nraster
==
2
)
FRYA_rawadc
=
pulseIntRaw
;
if
(
nraster
==
3
)
FRYB_rawadc
=
pulseIntRaw
;
}
if
(
hit
->
fADC_xsig
>
0
)
{
fRawXADC
=
hit
->
fADC_xsig
;
//std::cout<<" Raw X ADC = "<<fRawXADC<<std::endl;
}
if
(
hit
->
fADC_ysig
>
0
)
{
fRawYADC
=
hit
->
fADC_ysig
;
//std::cout<<" Raw Y ADC = "<<fRawYADC<<std::endl;
}
ihit
++
;
}
return
0
;
return
0
;
...
@@ -283,13 +374,12 @@ Int_t THcRaster::Decode( const THaEvData& evdata )
...
@@ -283,13 +374,12 @@ Int_t THcRaster::Decode( const THaEvData& evdata )
//_____________________________________________________________________________
//_____________________________________________________________________________
Int_t
THcRaster
::
Process
(
){
Int_t
THcRaster
::
Process
(
){
Double_t
e
Beam
=
0.001
;
Double_t
fgp
Beam
=
0.001
;
DBRequest
list
[]
=
{
DBRequest
list
[]
=
{
{
"gpbeam"
,
&
e
Beam
,
kDouble
,
0
,
1
},
{
"gpbeam"
,
&
fgp
Beam
,
kDouble
,
0
,
1
},
{
0
}
{
0
}
};
};
gHcParms
->
LoadParmValues
(
list
);
gHcParms
->
LoadParmValues
(
list
);
/*
/*
calculate raster position from ADC value.
calculate raster position from ADC value.
From ENGINE/g_analyze_misc.f -
From ENGINE/g_analyze_misc.f -
...
@@ -299,8 +389,11 @@ Int_t THcRaster::Process( ){
...
@@ -299,8 +389,11 @@ Int_t THcRaster::Process( ){
*/
*/
// calculate the raster currents
// calculate the raster currents
fXADC
=
fRawXADC
-
fAvgPedADC
[
0
];
fXA_ADC
=
FRXA_rawadc
-
fFrXA_ADC_zero_offset
;
fYADC
=
fRawYADC
-
fAvgPedADC
[
1
];
fYA_ADC
=
FRXB_rawadc
-
fFrYA_ADC_zero_offset
;
fXB_ADC
=
FRYA_rawadc
-
fFrXB_ADC_zero_offset
;
fYB_ADC
=
FRYB_rawadc
-
fFrYB_ADC_zero_offset
;
//std::cout<<" Raw X ADC = "<<fXADC<<" Raw Y ADC = "<<fYADC<<std::endl;
//std::cout<<" Raw X ADC = "<<fXADC<<" Raw Y ADC = "<<fYADC<<std::endl;
/*
/*
...
@@ -310,8 +403,10 @@ Int_t THcRaster::Process( ){
...
@@ -310,8 +403,10 @@ Int_t THcRaster::Process( ){
gfry = (gfry_adc/gfry_adcpercm)*(gfr_cal_mom/ebeam)
gfry = (gfry_adc/gfry_adcpercm)*(gfr_cal_mom/ebeam)
*/
*/
fXpos
=
(
fXADC
/
fFrXADCperCM
)
*
(
fFrCalMom
/
eBeam
);
fXA_pos
=
(
fXA_ADC
/
fFrXA_ADCperCM
)
*
(
fFrCalMom
/
fgpBeam
);
fYpos
=
(
fYADC
/
fFrYADCperCM
)
*
(
fFrCalMom
/
eBeam
);
fYA_pos
=
(
fYA_ADC
/
fFrYA_ADCperCM
)
*
(
fFrCalMom
/
fgpBeam
);
fXB_pos
=
(
fXB_ADC
/
fFrXB_ADCperCM
)
*
(
fFrCalMom
/
fgpBeam
);
fYB_pos
=
(
fYB_ADC
/
fFrYB_ADCperCM
)
*
(
fFrCalMom
/
fgpBeam
);
// std::cout<<" X = "<<fXpos<<" Y = "<<fYpos<<std::endl;
// std::cout<<" X = "<<fXpos<<" Y = "<<fYpos<<std::endl;
...
@@ -319,9 +414,9 @@ Int_t THcRaster::Process( ){
...
@@ -319,9 +414,9 @@ Int_t THcRaster::Process( ){
Double_t
tt
;
Double_t
tt
;
Double_t
tp
;
Double_t
tp
;
if
(
fgusefr
!=
0
)
{
if
(
fgusefr
!=
0
)
{
fPosition
[
1
].
SetXYZ
(
fXpos
+
fgbeam_xoff
,
fYpos
+
fgbeam_yoff
,
0.0
);
fPosition
[
1
].
SetXYZ
(
fX
A_
pos
+
fgbeam_xoff
,
fY
A_
pos
+
fgbeam_yoff
,
0.0
);
tt
=
fXpos
/
fgfrx_dist
+
fgbeam_xpoff
;
tt
=
fX
A_
pos
/
fgfrx_dist
+
fgbeam_xpoff
;
tp
=
fYpos
/
fgfry_dist
+
fgbeam_ypoff
;
tp
=
fY
A_
pos
/
fgfry_dist
+
fgbeam_ypoff
;
}
else
{
// Just use fixed beam position and angle
}
else
{
// Just use fixed beam position and angle
fPosition
[
0
].
SetXYZ
(
fgbeam_xoff
,
fgbeam_yoff
,
0.0
);
fPosition
[
0
].
SetXYZ
(
fgbeam_xoff
,
fgbeam_yoff
,
0.0
);
tt
=
fgbeam_xpoff
;
tt
=
fgbeam_xpoff
;
...
...
This diff is collapsed.
Click to expand it.
src/THcRaster.h
+
49
−
19
View file @
1dce3be1
//author Burcu Duran - Melanie Rehfuss (2017)
#ifndef ROOT_THcRaster
#ifndef ROOT_THcRaster
#define ROOT_THcRaster
#define ROOT_THcRaster
...
@@ -15,6 +18,7 @@
...
@@ -15,6 +18,7 @@
#include
"THcRasterRawHit.h"
#include
"THcRasterRawHit.h"
#include
"THaCutList.h"
#include
"THaCutList.h"
class
THcRaster
:
public
THaBeamDet
,
public
THcHitList
{
class
THcRaster
:
public
THaBeamDet
,
public
THcHitList
{
public:
public:
...
@@ -22,9 +26,20 @@ class THcRaster : public THaBeamDet, public THcHitList {
...
@@ -22,9 +26,20 @@ class THcRaster : public THaBeamDet, public THcHitList {
THcRaster
(
const
char
*
name
,
const
char
*
description
=
""
,
THaApparatus
*
a
=
NULL
);
THcRaster
(
const
char
*
name
,
const
char
*
description
=
""
,
THaApparatus
*
a
=
NULL
);
~
THcRaster
();
~
THcRaster
();
void
Clear
(
Option_t
*
opt
=
""
);
void
AccumulatePedestals
(
TClonesArray
*
rawhits
);
void
CalculatePedestals
();
Int_t
Decode
(
const
THaEvData
&
);
Int_t
ReadDatabase
(
const
TDatime
&
date
);
Int_t
DefineVariables
(
EMode
mode
);
EStatus
Init
(
const
TDatime
&
run_time
);
EStatus
Init
(
const
TDatime
&
run_time
);
Int_t
Decode
(
const
THaEvData
&
);
void
InitArrays
()
{
/* do nothing */
;}
void
DeleteArrays
()
{
/* do nothing */
;}
Int_t
Process
();
Int_t
Process
();
// TVector3 GetPosition() const { return fPosition[2]; }
// TVector3 GetPosition() const { return fPosition[2]; }
...
@@ -34,12 +49,16 @@ class THcRaster : public THaBeamDet, public THcHitList {
...
@@ -34,12 +49,16 @@ class THcRaster : public THaBeamDet, public THcHitList {
Double_t
GetCurrentX
()
{
return
fRawPos
[
0
];
}
Double_t
GetCurrentX
()
{
return
fRawPos
[
0
];
}
Double_t
GetCurrentY
()
{
return
fRawPos
[
1
];
}
Double_t
GetCurrentY
()
{
return
fRawPos
[
1
];
}
THcRaster
();
protected
:
protected
:
/* void InitializeReconstruction(); */
//event information
Int_t
ReadDatabase
(
const
TDatime
&
date
);
Int_t
DefineVariables
(
EMode
mode
);
Int_t
fNhits
;
/* void InitializeReconstruction(); */
Double_t
fgpbeam
;
//beam momentum
Double_t
fgpbeam
;
//beam momentum
Double_t
fgfrx_dist
;
//Distance of raster to target
Double_t
fgfrx_dist
;
//Distance of raster to target
Double_t
fgfry_dist
;
Double_t
fgfry_dist
;
...
@@ -49,31 +68,42 @@ class THcRaster : public THaBeamDet, public THcHitList {
...
@@ -49,31 +68,42 @@ class THcRaster : public THaBeamDet, public THcHitList {
Double_t
fgbeam_ypoff
;
// Beam offsets
Double_t
fgbeam_ypoff
;
// Beam offsets
Int_t
fgusefr
;
/* Use Raster for beam position */
Int_t
fgusefr
;
/* Use Raster for beam position */
Double_t
fRawXADC
;
// X raw ADC
Double_t
FRXA_rawadc
;
// XA raw ADC
Double_t
fRawYADC
;
// Y raw ADC
Double_t
FRYA_rawadc
;
// YA raw ADC
Double_t
fXADC
;
// X ADC
Double_t
FRXB_rawadc
;
// XB raw ADC
Double_t
fYADC
;
// Y ADC
Double_t
FRYB_rawadc
;
// YB raw ADC
Double_t
fXpos
;
// X position
Double_t
fXA_ADC
;
// XA ADC
Double_t
fYpos
;
// Y position
Double_t
fYA_ADC
;
// YA ADC
Double_t
fXB_ADC
;
// XB ADC
Double_t
fYB_ADC
;
// YB ADC
Double_t
fXA_pos
;
// XA position
Double_t
fYA_pos
;
// YA position
Double_t
fXB_pos
;
// XB position
Double_t
fYB_pos
;
// YB position
Double_t
fFrXA_ADC_zero_offset
;
Double_t
fFrYA_ADC_zero_offset
;
Double_t
fFrXB_ADC_zero_offset
;
Double_t
fFrYB_ADC_zero_offset
;
Double_t
fPedADC
[
2
];
// ADC poedestals
Double_t
fAvgPedADC
[
2
];
// Avergage ADC poedestals
Double_t
fPedADC
[
4
];
// ADC poedestals
//Double_t fAvgPedADC[4]; // Avergage ADC poedestals
Double_t
fRawPos
[
2
];
// current in Raster ADCs for position
Double_t
fRawPos
[
2
];
// current in Raster ADCs for position
TVector3
fPosition
[
2
];
// Beam position at 1st, 2nd BPM or at the target (meters)
TVector3
fPosition
[
4
];
// Beam position at 1st, 2nd BPM or at the target (meters)
TVector3
fDirection
;
TVector3
fDirection
;
TClonesArray
*
frPosAdcPulseIntRaw
;
private
:
private
:
Bool_t
fAnalyzePedestals
;
Bool_t
fAnalyzePedestals
;
Int_t
fNPedestalEvents
;
Int_t
fNPedestalEvents
;
Double_t
fFrCalMom
;
Double_t
fFrCalMom
;
Double_t
fFrXADCperCM
;
Double_t
fFrXA_ADCperCM
;
Double_t
fFrYADCperCM
;
Double_t
fFrYA_ADCperCM
;
Double_t
fFrXB_ADCperCM
;
Double_t
fFrYB_ADCperCM
;
void
CalculatePedestals
();
void
AccumulatePedestals
(
TClonesArray
*
rawhits
);
ClassDef
(
THcRaster
,
0
);
// add THcRaster to ROOT library
ClassDef
(
THcRaster
,
0
);
// add THcRaster to ROOT library
};
};
...
...
This diff is collapsed.
Click to expand it.
src/THcRasterRawHit.cxx
+
1
−
68
View file @
1dce3be1
// Author : Buddhini Waidyawansa
// Author : Buddhini Waidyawansa
// Date : 23-01-2014
// Date : 23-01-2014
//Author : Burcu Duran - Melanie Rehfuss (2017)
///////////////////////////////////////////////////////////////////////////////
///////////////////////////////////////////////////////////////////////////////
...
@@ -19,74 +20,6 @@
...
@@ -19,74 +20,6 @@
#include
"THcRasterRawHit.h"
#include
"THcRasterRawHit.h"
using
namespace
std
;
void
THcRasterRawHit
::
SetData
(
Int_t
signal
,
Int_t
data
)
{
if
(
signal
==
0
)
{
fADC_xsync
=
data
;
}
else
if
(
signal
==
1
)
{
fADC_xsig
=
data
;
}
else
if
(
signal
==
2
)
{
fADC_ysync
=
data
;
}
else
if
(
signal
==
3
)
{
fADC_ysig
=
data
;
}
// std::cout<<" xsync = "<<fADC_xsync
// <<" xsig = "<<fADC_xsig
// <<" ysync = "<<fADC_ysync
// <<" ysig = "<<fADC_ysig << std::endl;
}
Int_t
THcRasterRawHit
::
GetData
(
Int_t
signal
)
{
if
(
signal
==
1
)
{
return
(
fADC_xsync
);
}
else
if
(
signal
==
2
)
{
return
(
fADC_xsig
);
}
else
if
(
signal
==
3
)
{
return
(
fADC_ysync
);
}
else
if
(
signal
==
4
)
{
return
(
fADC_ysig
);
}
return
(
-
1
);
}
// Int_t THcRasterRawHit::Compare(const TObject* obj) const
// {
// // Compare to sort by the plane
// // There is only one raster so no need for an additional check on the counter
// const THcRasterRawHit* hit = dynamic_cast<const THcRasterRawHit*>(obj);
// if(!hit) return -1;
// Int_t p1 = fPlane;
// Int_t p2 = hit->fPlane;
// if(p1 < p2) return -1;
// else if(p1 > p2) return 1;
// }
//_____________________________________________________________________________
THcRasterRawHit
&
THcRasterRawHit
::
operator
=
(
const
THcRasterRawHit
&
rhs
)
{
// Assignment operator.
THcRawHit
::
operator
=
(
rhs
);
if
(
this
!=
&
rhs
)
{
fPlane
=
rhs
.
fPlane
;
fCounter
=
rhs
.
fCounter
;
fADC_xsync
=
rhs
.
fADC_xsync
;
fADC_xsig
=
rhs
.
fADC_xsig
;
fADC_ysync
=
rhs
.
fADC_ysync
;
fADC_ysig
=
rhs
.
fADC_ysig
;
}
return
*
this
;
}
//////////////////////////////////////////////////////////////////////////
//////////////////////////////////////////////////////////////////////////
ClassImp
(
THcRasterRawHit
)
ClassImp
(
THcRasterRawHit
)
This diff is collapsed.
Click to expand it.
src/THcRasterRawHit.h
+
5
−
21
View file @
1dce3be1
//Author : Burcu Duran - Melanie Rehfuss (2017)
#ifndef ROOT_THcRasterRawHit
#ifndef ROOT_THcRasterRawHit
#define ROOT_THcRasterRawHit
#define ROOT_THcRasterRawHit
...
@@ -8,31 +10,13 @@
...
@@ -8,31 +10,13 @@
// //
// //
///////////////////////////////////////////////////////////////////////////////
///////////////////////////////////////////////////////////////////////////////
#include
"THcRawHit.h"
#include
"THcRawH
odoH
it.h"
class
THcRasterRawHit
:
public
THcRawHit
{
class
THcRasterRawHit
:
public
THcRawH
odoH
it
{
public:
public:
THcRasterRawHit
(
Int_t
plane
=
0
,
Int_t
counter
=
0
)
:
THcRawHit
(
plane
,
counter
),
friend
class
THcRaster
;
fADC_xsig
(
-
1
),
fADC_ysig
(
-
1
),
fADC_xsync
(
-
1
),
fADC_ysync
(
-
1
)
{
}
THcRasterRawHit
&
operator
=
(
const
THcRasterRawHit
&
);
~
THcRasterRawHit
()
{}
void
Clear
(
Option_t
*
opt
=
""
)
{
fADC_xsig
=
-
1
;
fADC_ysig
=
-
1
;
fADC_xsync
=
-
1
;
fADC_ysync
=
-
1
;
}
void
SetData
(
Int_t
signal
,
Int_t
data
);
Int_t
GetData
(
Int_t
signal
);
// signals
Int_t
fADC_xsig
;
Int_t
fADC_ysig
;
Int_t
fADC_xsync
;
Int_t
fADC_ysync
;
protected:
protected:
...
...
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment