Hall A ROOT/C++ Analyzer (podd)
Loading...
Searching...
No Matches
THaReacPointFoil.cxx
Go to the documentation of this file.
1//*-- Author : BR
2
4//
5// THaReacPointFoil
6//
7// Calculate vertex for all the tracks from the beam
8// apparatus, assuming a single foil target at z = 0
9//
10//
12
13#include "THaReacPointFoil.h"
14#include "THaSpectrometer.h"
15#include "THaTrack.h"
16#include "THaBeam.h"
17#include "TMath.h"
18
20
21//_____________________________________________________________________________
22THaReacPointFoil::THaReacPointFoil( const char* name, const char* description,
23 const char* spectro, const char* beam ) :
24 THaPhysicsModule(name,description), fSpectroName(spectro),
25 fBeamName(beam), fSpectro(nullptr), fBeam(nullptr)
26{
27 // Normal constructor.
28
29}
30
31//_____________________________________________________________________________
33{
34 // Destructor
35
37}
38
39//_____________________________________________________________________________
41{
42 // Clear all event-by-event variables.
43
46}
47
48//_____________________________________________________________________________
50{
51 // Define/delete analysis variables
52
54}
55
56//_____________________________________________________________________________
58{
59 // Initialize the module.
60 // Locate the spectrometer apparatus named in fSpectroName and save
61 // pointer to it.
62
63 // Standard initialization. Calls this object's DefineVariables().
64 if( THaPhysicsModule::Init( run_time ) != kOK )
65 return fStatus;
66
67 fSpectro = static_cast<THaSpectrometer*>
68 ( FindModule( fSpectroName.Data(), "THaSpectrometer"));
69 if( !fSpectro )
70 return fStatus;
71
72 if( fBeamName.Length() > 0 )
73 fBeam = static_cast<THaBeam*>( FindModule( fBeamName.Data(), "THaBeam") );
74
75 return fStatus;
76}
77
78//_____________________________________________________________________________
80{
81 // Calculate the vertex coordinates.
82
83 if( !IsOK() ) return -1;
84
85 Int_t ntracks = fSpectro->GetNTracks();
86 if( ntracks == 0 ) return 0;
87
89 if( !tracks ) return -2;
90
91 TVector3 beam_org, beam_ray( 0.0, 0.0, 1.0 );
92 if( fBeam ) {
93 beam_org = fBeam->GetPosition();
94 beam_ray = fBeam->GetDirection();
95 }
96 static const TVector3 yax( 0.0, 1.0, 0.0 );
97 static const TVector3 xax( 1.0, 0.0, 0.0 );
98 TVector3 org, v;
99
100 for( Int_t i = 0; i<ntracks; i++ ) {
101 auto* theTrack = static_cast<THaTrack*>( tracks->At(i) );
102 // Ignore junk tracks
103 if( !theTrack || !theTrack->HasTarget() )
104 continue;
105 org.SetX( 0. );
106 org.SetZ( 0. );
107 org.SetY( 0. );
108 Double_t t = 0; // dummy
109 if( !IntersectPlaneWithRay( xax, yax, org,
110 beam_org, beam_ray, t, v ))
111 continue; // Oops, track and beam parallel?
112 theTrack->SetVertex(v);
113
114 // FIXME: preliminary
115 if( theTrack == fSpectro->GetGoldenTrack() ) {
116 fVertex = theTrack->GetVertex();
117 fVertexOK = true;
118 }
119 // FIXME: calculate vertex coordinate errors here (need beam errors)
120
121
122 }
123 return 0;
124}
125
int Int_t
double Double_t
const char Option_t
Option_t Option_t TPoint TPoint const char mode
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t org
char name[80]
static Int_t DefineVarsFromList(const void *list, EType type, EMode mode, const char *def_prefix, const TObject *obj, const char *prefix, const char *here, const char *comment_subst="")
THaAnalysisObject * FindModule(const char *name, const char *classname, bool do_error=true)
static Bool_t IntersectPlaneWithRay(const TVector3 &xax, const TVector3 &yax, const TVector3 &org, const TVector3 &ray_start, const TVector3 &ray_vect, Double_t &length, TVector3 &intersect)
Bool_t IsOK() const
virtual const TVector3 & GetDirection() const
Definition THaBeam.h:24
virtual const TVector3 & GetPosition() const
Definition THaBeam.h:23
virtual void Clear(Option_t *opt="")
THaReacPointFoil(const char *name, const char *description, const char *spectro="", const char *beam="")
virtual Int_t Process(const THaEvData &)
virtual Int_t DefineVariables(EMode mode=kDefine)
THaSpectrometer * fSpectro
virtual void Clear(Option_t *opt="")
THaTrack * GetGoldenTrack() const
Int_t GetNTracks() const
TClonesArray * GetTracks() const
virtual void VertexClear()
static const RVarDef * GetRVarDef()
void SetX(SCoord_t x)
void SetY(SCoord_t y)
Ssiz_t Length() const
const char * Data() const
v
ClassImp(TPyArg)
void tracks()