StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFixedVertexFinder.cxx
1 /*
2  * StFixedVertexFinder.cxx
3  * $Id: StFixedVertexFinder.cxx,v 1.8 2018/03/24 20:10:41 jwebb Exp $
4  *
5  * Author Lee Barnby (University of Birmingham) May 2006.
6  *
7  */
8 
9 #include "StChain/StMaker.h"
10 #include "StGenericVertexMaker/StFixedVertexFinder.h"
11 #include "StMcEvent/StMcEvent.hh"
12 #include "StMcEvent/StMcVertex.hh"
13 #include "St_base/StMessMgr.h"
14 
15 
16 StFixedVertexFinder::StFixedVertexFinder() :
18  mFixedX(0),
19  mFixedY(0),
20  mFixedZ(0),
21  mFixedEx(0),
22  mFixedEy(0),
23  mFixedEz(0)
24 {
25  mMode=0;
26 }
27 
28 int StFixedVertexFinder::fit(StEvent* event)
29 {
30  StPrimaryVertex prim;
31  StVertexFinderId VFId;
32 
33  // All errors are zero... by default
34  // if user has set a vertex error, we use it, regardless of mode
35 
36  float s2x = mFixedEx * mFixedEx;
37  float s2y = mFixedEy * mFixedEy;
38  float s2z = mFixedEz * mFixedEz;
39 
40  float cov[6] = {s2x,
41  0.0,s2y,
42  0.0,0.0,s2z};
43 
44 
45  if (mMode == 0){
46  // Really the default which takes the SetPos() TBI
47  LOG_DEBUG << "StFixedVertexFinder::fit() fixing a vertex" << endm;
48  StThreeVectorD pos(mFixedX,mFixedY,mFixedZ);
49  prim.setPosition(pos);
50  VFId = undefinedVertexFinder;
51 
52  } else {
53 
54  LOG_DEBUG << "StFixedVertexFinder::fit() reading a vertex from StMcEvent" << endm;
55  StMcEvent* mcEvent= (StMcEvent*) StMaker::GetChain()->GetDataSet("StMcEvent");
56 
57  if ( ! mcEvent){
58  // Not sure what to do here, will leave and do nothing (so no vertex)
59  // and it will be super obvious that something went terribly wrong
60  LOG_ERROR << "StFixedVertexFinder::fit() no StMcEvent" << endm;
61  return 0;
62 
63  } else {
64  StMcVertex *pv = mcEvent->primaryVertex();
65  StThreeVectorD pos(pv->position().x(),
66  pv->position().y(),
67  pv->position().z());
68  //cout << "DEBUG :: saving McVertex " << pos << endl;
69  prim.setPosition(pos);
70  }
71  VFId = mcEventVertexFFinder;
72  }
73 
74  prim.setCovariantMatrix(cov);
75  prim.setFlag(1); // So that we know its the primary vertex
76  prim.setRanking(-5); // Have to have something
77  prim.setVertexFinderId(VFId); // Id depends on MC or fixed position used
78  addVertex(prim);
79 
80  return size();
81 }
82 
83 void StFixedVertexFinder::printInfo(ostream& os)const{
84  os << "StFixedVertexFinder - fixed vertex" << endl;
85  os << "Fixed position: x=" << mFixedX << " y=" << mFixedY << " z=" << mFixedZ << endl;
86 }
87 
88 void StFixedVertexFinder::UseVertexConstraint(){
89  LOG_WARN << "StFixedVertexFinder::UseVertexConstraint() - vertex beam constraint NOT implemented in context of fixed vertex finder" << endm;
90 
91 }
92 
93 void StFixedVertexFinder::SetVertexPosition(double x, double y, double z){
94  mFixedX=x;
95  mFixedY=y;
96  mFixedZ=z;
97 }
98 void StFixedVertexFinder::SetVertexError(double x, double y, double z){
99  mFixedEx=x;
100  mFixedEy=y;
101  mFixedEz=z;
102 }
103 
104 /*
105  * $Log: StFixedVertexFinder.cxx,v $
106  * Revision 1.8 2018/03/24 20:10:41 jwebb
107  * Added option for user to specify the uncertainties on the vertex. Useful
108  * in embedding jobs in order to get the track association with primary
109  * vertex correct (especially when tracks are from precision tracking, eg
110  * HFT).
111  *
112  * Revision 1.7 2017/05/12 18:37:23 smirnovd
113  * Cosmetic changes
114  *
115  * Removed log messages from source files
116  * Prefixed included headers with paths to respective modules
117  *
118  * Revision 1.6 2017/01/03 22:17:36 smirnovd
119  * [Stylistic] Changed public addVertex() to accept references
120  *
121  * Avoid unnecessary gymnastics with pointers
122  *
123  * Revision 1.5 2016/08/18 17:46:14 smirnovd
124  * Squashed commit of the following refactoring changes:
125  *
126  * Date: Wed Jul 27 18:31:18 2016 -0400
127  *
128  * Removed unused arguments in UseVertexConstraint()
129  *
130  * In StiPPVertexFinder and StvPPVertexFinder this method does nothing
131  *
132  * Date: Wed Jul 27 16:47:58 2016 -0400
133  *
134  * Make old UseVertexConstraint private virtual and call it from its public replacement in the base class
135  *
136  * also mark methods as private explicitly
137  *
138  * Date: Wed Jul 27 16:52:02 2016 -0400
139  *
140  * Removed unused private data member mWeight
141  *
142  * Date: Wed Jul 27 16:50:42 2016 -0400
143  *
144  * Prefer base class static beamline parameters rather than this class private members
145  *
146  * Date: Wed Jul 27 16:21:49 2016 -0400
147  *
148  * StPPVertexFinder: Got rid of unused private beamline parameters
149  *
150  * The equivalent measurements are available from the base class
151  * StGenericVertexFinder
152  *
153  * Date: Wed Jul 27 16:19:19 2016 -0400
154  *
155  * StPPVertexFinder: For beamline position use equivalent static methods from parent class
156  *
157  * Date: Wed Jul 27 16:05:50 2016 -0400
158  *
159  * StGenericVertexMaker: Assigning once is enough
160  *
161  * Date: Mon Aug 15 10:43:49 2016 -0400
162  *
163  * StGenericVertexFinder: Print out beamline parameters
164  *
165  * Print beamline values as extracted from the database before any modification.
166  *
167  * Date: Wed Jul 6 15:33:02 2016 -0400
168  *
169  * Stylistic changes and minor refactoring
170  *
171  * Whitespace and comments for improved readability
172  * s/track/stiKalmanTrack/
173  *
174  * Date: Wed Jul 6 15:28:16 2016 -0400
175  *
176  * StPPVertexFinder: Switched to cleaner c++11 range loop syntax
177  *
178  * Date: Wed Jul 6 15:22:14 2016 -0400
179  *
180  * StPPVertexFinder: Minor c++ refactoring
181  *
182  * - Removed unused counter
183  * - c-style array to std::array
184  *
185  * Date: Wed Jul 6 15:20:11 2016 -0400
186  *
187  * Deleted commented out code
188  *
189  * Removed unused #include's StMinuitVertexFinder
190  *
191  * Revision 1.4 2006/05/18 19:14:24 lbarnby
192  * Added SetVertexPosition function. Tidied up comments/docs
193  *
194  */
Event data structure to hold all information from a Monte Carlo simulation. This class is the interfa...
Definition: StMcEvent.hh:169