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
3b7efe26
Commit
3b7efe26
authored
7 years ago
by
Stephen A. Wood
Committed by
Mark K Jones
7 years ago
Browse files
Options
Downloads
Patches
Plain Diff
Apply oop correction before computing kinematics
parent
12664e6e
No related branches found
No related tags found
No related merge requests found
Changes
2
Hide whitespace changes
Inline
Side-by-side
Showing
2 changed files
src/THcPrimaryKine.cxx
+24
-7
24 additions, 7 deletions
src/THcPrimaryKine.cxx
src/THcPrimaryKine.h
+3
-2
3 additions, 2 deletions
src/THcPrimaryKine.h
with
27 additions
and
9 deletions
src/THcPrimaryKine.cxx
+
24
−
7
View file @
3b7efe26
...
@@ -6,7 +6,7 @@ These are usually the electron kinematics.
...
@@ -6,7 +6,7 @@ These are usually the electron kinematics.
*/
*/
#include
"THcPrimaryKine.h"
#include
"THcPrimaryKine.h"
#include
"TH
aTrackingModule
.h"
#include
"TH
cHallCSpectrometer
.h"
#include
"THcGlobals.h"
#include
"THcGlobals.h"
#include
"THcParmList.h"
#include
"THcParmList.h"
#include
"THaRunBase.h"
#include
"THaRunBase.h"
...
@@ -106,17 +106,21 @@ THaAnalysisObject::EStatus THcPrimaryKine::Init( const TDatime& run_time )
...
@@ -106,17 +106,21 @@ THaAnalysisObject::EStatus THcPrimaryKine::Init( const TDatime& run_time )
// Locate the spectrometer apparatus named in fSpectroName and save
// Locate the spectrometer apparatus named in fSpectroName and save
// pointer to it.
// pointer to it.
fSpectro
=
dynamic_cast
<
THaTrackingModule
*>
fSpectro
=
dynamic_cast
<
THcHallCSpectrometer
*>
(
FindModule
(
fSpectroName
.
Data
(),
"THaTrackingModule"
));
(
FindModule
(
fSpectroName
.
Data
(),
"THcHallCSpectrometer"
));
if
(
!
fSpectro
)
if
(
!
fSpectro
)
{
fStatus
=
kInitError
;
return
fStatus
;
return
fStatus
;
}
// Optional beam apparatus
// Optional beam apparatus
if
(
fBeamName
.
Length
()
>
0
)
{
if
(
fBeamName
.
Length
()
>
0
)
{
fBeam
=
dynamic_cast
<
THaBeamModule
*>
fBeam
=
dynamic_cast
<
THaBeamModule
*>
(
FindModule
(
fBeamName
.
Data
(),
"THaBeamModule"
)
);
(
FindModule
(
fBeamName
.
Data
(),
"THaBeamModule"
)
);
if
(
!
fBeam
)
if
(
!
fBeam
)
{
fStatus
=
kInitError
;
return
fStatus
;
return
fStatus
;
}
if
(
fM
<=
0.0
)
if
(
fM
<=
0.0
)
fM
=
fBeam
->
GetBeamInfo
()
->
GetM
();
fM
=
fBeam
->
GetBeamInfo
()
->
GetM
();
}
}
...
@@ -144,6 +148,11 @@ Int_t THcPrimaryKine::Process( const THaEvData& )
...
@@ -144,6 +148,11 @@ Int_t THcPrimaryKine::Process( const THaEvData& )
// Determine 4-momentum of incident particle.
// Determine 4-momentum of incident particle.
// If a beam module given, use it to get the beam momentum. This
// If a beam module given, use it to get the beam momentum. This
// module may apply corrections for beam energy loss, variations, etc.
// module may apply corrections for beam energy loss, variations, etc.
Double_t
xptar
=
trkifo
->
GetTheta
()
+
fOopCentralOffset
;
TVector3
pvect
;
fSpectro
->
TransportToLab
(
trkifo
->
GetP
(),
xptar
,
trkifo
->
GetPhi
(),
pvect
);
if
(
fBeam
)
{
if
(
fBeam
)
{
fP0
.
SetVectM
(
fBeam
->
GetBeamInfo
()
->
GetPvect
(),
fM
);
fP0
.
SetVectM
(
fBeam
->
GetBeamInfo
()
->
GetPvect
(),
fM
);
}
else
{
}
else
{
...
@@ -152,7 +161,7 @@ Int_t THcPrimaryKine::Process( const THaEvData& )
...
@@ -152,7 +161,7 @@ Int_t THcPrimaryKine::Process( const THaEvData& )
fP0
.
SetXYZM
(
0.0
,
0.0
,
p_in
,
fM
);
fP0
.
SetXYZM
(
0.0
,
0.0
,
p_in
,
fM
);
}
}
fP1
.
SetVectM
(
trkifo
->
GetP
vect
()
,
fM
);
fP1
.
SetVectM
(
p
vect
,
fM
);
fA
.
SetXYZM
(
0.0
,
0.0
,
0.0
,
fMA
);
// Assume target at rest
fA
.
SetXYZM
(
0.0
,
0.0
,
0.0
,
fMA
);
// Assume target at rest
// proton mass (for x_bj)
// proton mass (for x_bj)
...
@@ -161,10 +170,11 @@ Int_t THcPrimaryKine::Process( const THaEvData& )
...
@@ -161,10 +170,11 @@ Int_t THcPrimaryKine::Process( const THaEvData& )
// Standard electron kinematics
// Standard electron kinematics
fQ
=
fP0
-
fP1
;
fQ
=
fP0
-
fP1
;
// cqx, cqy, cqz, omega
fQ2
=
-
fQ
.
M2
();
fQ2
=
-
fQ
.
M2
();
fQ3mag
=
fQ
.
P
();
fQ3mag
=
fQ
.
P
();
fOmega
=
fQ
.
E
();
fOmega
=
fQ
.
E
();
// cqxzabs = TMath::Sqrt(fQ.X()*fQ.X() + fQ.Y()*fQ.Y());
fA1
=
fA
+
fQ
;
fA1
=
fA
+
fQ
;
// fW2 = fA1.M2();
// fW2 = fA1.M2();
fMp1
=
fMp
+
fQ
;
fMp1
=
fMp
+
fQ
;
...
@@ -189,8 +199,15 @@ Int_t THcPrimaryKine::ReadDatabase( const TDatime& date )
...
@@ -189,8 +199,15 @@ Int_t THcPrimaryKine::ReadDatabase( const TDatime& date )
cout
<<
"In THcPrimaryKine::ReadDatabase()"
<<
endl
;
cout
<<
"In THcPrimaryKine::ReadDatabase()"
<<
endl
;
#endif
#endif
char
prefix
[
2
];
prefix
[
0
]
=
tolower
(
GetName
()[
0
]);
prefix
[
1
]
=
'\0'
;
fOopCentralOffset
=
0.0
;
DBRequest
list
[]
=
{
DBRequest
list
[]
=
{
{
"gtargmass_amu"
,
&
fMA_amu
,
kDouble
,
0
,
1
},
{
"gtargmass_amu"
,
&
fMA_amu
,
kDouble
,
0
,
1
},
{
Form
(
"%s_oopcentral_offset"
,
prefix
),
&
fOopCentralOffset
,
kDouble
,
0
,
1
},
{
0
}
{
0
}
};
};
gHcParms
->
LoadParmValues
((
DBRequest
*
)
&
list
);
gHcParms
->
LoadParmValues
((
DBRequest
*
)
&
list
);
...
...
This diff is collapsed.
Click to expand it.
src/THcPrimaryKine.h
+
3
−
2
View file @
3b7efe26
...
@@ -11,7 +11,7 @@
...
@@ -11,7 +11,7 @@
#include
"TLorentzVector.h"
#include
"TLorentzVector.h"
#include
"TString.h"
#include
"TString.h"
class
TH
aTrackingModule
;
class
TH
cHallCSpectrometer
;
class
THaBeamModule
;
class
THaBeamModule
;
typedef
TLorentzVector
FourVect
;
typedef
TLorentzVector
FourVect
;
...
@@ -80,12 +80,13 @@ protected:
...
@@ -80,12 +80,13 @@ protected:
Double_t
fM
;
// Mass of incident particle (GeV/c^2)
Double_t
fM
;
// Mass of incident particle (GeV/c^2)
Double_t
fMA
;
// Target mass (GeV/c^2)
Double_t
fMA
;
// Target mass (GeV/c^2)
Double_t
fMA_amu
;
// Target mass (amu)
Double_t
fMA_amu
;
// Target mass (amu)
Double_t
fOopCentralOffset
;
// Out plane offset of spectrometer
virtual
Int_t
DefineVariables
(
EMode
mode
=
kDefine
);
virtual
Int_t
DefineVariables
(
EMode
mode
=
kDefine
);
TString
fSpectroName
;
// Name of spectrometer to consider
TString
fSpectroName
;
// Name of spectrometer to consider
TString
fBeamName
;
// Name of beam position apparatus
TString
fBeamName
;
// Name of beam position apparatus
TH
aTrackingModule
*
fSpectro
;
// Pointer to spectrometer object
TH
cHallCSpectrometer
*
fSpectro
;
// Pointer to spectrometer object
THaBeamModule
*
fBeam
;
// Pointer to beam position apparatus
THaBeamModule
*
fBeam
;
// Pointer to beam position apparatus
ClassDef
(
THcPrimaryKine
,
0
)
//Single arm kinematics module
ClassDef
(
THcPrimaryKine
,
0
)
//Single arm kinematics module
...
...
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