CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
NtupleManager Class Reference

#include <NtupleManager.h>

Public Member Functions

void BookNtuple ()
 
void FillChi2 ()
 
void FillFitParameters (MatrixMeschach *AtWAMatrix)
 
void FillMeasurements ()
 
void FillNtupleTree ()
 
void FillOptObjects (MatrixMeschach *AtWAMatrix)
 
void InitNtuple ()
 
 NtupleManager ()
 
void WriteNtuple ()
 
 ~NtupleManager ()
 

Static Public Member Functions

static NtupleManagergetInstance ()
 

Private Member Functions

void GetGlobalAngles (const CLHEP::HepRotation &rmGlob, double *theta)
 

Private Attributes

double Chi2CalibratedParameters
 
double Chi2Measurements
 
TClonesArray * CloneCopsMeas
 
TClonesArray * CloneDistancemeter1DimMeas
 
TClonesArray * CloneDistancemeterMeas
 
TClonesArray * CloneFitParam
 
TClonesArray * CloneOptObject
 
TClonesArray * CloneSensor2DMeas
 
TClonesArray * CloneTiltmeterMeas
 
TTree * CocoaTree
 
CopsMeasCopsMeasA
 
Distancemeter1DimMeasDistancemeter1DimMeasA
 
DistancemeterMeasDistancemeterMeasA
 
FitParamFitParamA
 
int NCops
 
int NDegreesOfFreedom
 
int NDistancemeter
 
int NDistancemeter1Dim
 
int NFitParameters
 
int NOptObjects
 
int NSensor2D
 
int NTiltmeter
 
OptObjectOptObjectA
 
Sensor2DMeasSensor2DMeasA
 
TFile * theRootFile
 
TiltmeterMeasTiltmeterMeasA
 

Static Private Attributes

static NtupleManagerinstance = nullptr
 

Detailed Description

Definition at line 20 of file NtupleManager.h.

Constructor & Destructor Documentation

◆ NtupleManager()

NtupleManager::NtupleManager ( )
inline

Definition at line 23 of file NtupleManager.h.

Referenced by getInstance().

23 {};

◆ ~NtupleManager()

NtupleManager::~NtupleManager ( )
inline

Definition at line 24 of file NtupleManager.h.

24 {};

Member Function Documentation

◆ BookNtuple()

void NtupleManager::BookNtuple ( )

Definition at line 39 of file NtupleManager.cc.

References Chi2CalibratedParameters, Chi2Measurements, CloneCopsMeas, CloneDistancemeter1DimMeas, CloneDistancemeterMeas, CloneFitParam, CloneOptObject, CloneSensor2DMeas, CloneTiltmeterMeas, CocoaTree, NCops, NDegreesOfFreedom, NDistancemeter, NDistancemeter1Dim, NFitParameters, NOptObjects, NSensor2D, NTiltmeter, and theRootFile.

Referenced by Fit::startFit().

39  {
40  theRootFile = new TFile("report.root", "RECREATE", "Simple ROOT Ntuple");
41 
42  CocoaTree = new TTree("CocoaTree", "CocoaTree");
43 
44  CocoaTree->Branch("Chi2Measurements", &Chi2Measurements, "Chi2Measurements/D");
45  CocoaTree->Branch("Chi2CalibratedParameters", &Chi2CalibratedParameters, "Chi2CalibratedParameters/D");
46  CocoaTree->Branch("NDegreesOfFreedom", &NDegreesOfFreedom, "NDegreesOfFreedom/I");
47 
48  CloneFitParam = new TClonesArray("FitParam");
49  CocoaTree->Branch("FitParameters", &CloneFitParam, 32000, 2);
50  CocoaTree->Branch("NFitParameters", &NFitParameters, "NFitParameters/I");
51 
52  CloneOptObject = new TClonesArray("OptObject");
53  CocoaTree->Branch("OptObjects", &CloneOptObject, 32000, 2);
54  CocoaTree->Branch("NOptObjects", &NOptObjects, "NOptObjects/I");
55 
56  CloneSensor2DMeas = new TClonesArray("Sensor2DMeas");
57  CocoaTree->Branch("Sensor2DMeasurements", &CloneSensor2DMeas, 32000, 2);
58  CocoaTree->Branch("NSensor2D", &NSensor2D, "NSensor2D/I");
59 
60  CloneDistancemeterMeas = new TClonesArray("DistancemeterMeas");
61  CocoaTree->Branch("DistancemeterMeasurements", &CloneDistancemeterMeas, 32000, 2);
62  CocoaTree->Branch("NDistancemeter", &NDistancemeter, "NDistancemeter/I");
63 
64  CloneDistancemeter1DimMeas = new TClonesArray("Distancemeter1DimMeas");
65  CocoaTree->Branch("Distancemeter1DimMeasurements", &CloneDistancemeter1DimMeas, 32000, 2);
66  CocoaTree->Branch("NDistancemeter1Dim", &NDistancemeter1Dim, "NDistancemeter1Dim/I");
67 
68  CloneTiltmeterMeas = new TClonesArray("TiltmeterMeas");
69  CocoaTree->Branch("TiltmeterMeasurements", &CloneTiltmeterMeas, 32000, 2);
70  CocoaTree->Branch("NTiltmeter", &NTiltmeter, "NTiltmeter/I");
71 
72  CloneCopsMeas = new TClonesArray("CopsMeas");
73  CocoaTree->Branch("CopsMeasurements", &CloneCopsMeas, 32000, 2);
74  CocoaTree->Branch("NCops", &NCops, "NCops/I");
75 
76  theRootFile->Add(CocoaTree);
77 
78  // FitParametersTree = new TTree("FitParametersTree","FitParametersTree");
79  // FitParametersTree->Branch("NFitParameters",&NFitParameters,"NFitParameters/I");
80  // BookFitParameters = false;
81  // theRootFile->Add(FitParametersTree);
82 
83  // MeasurementsTree = new TTree("MeasurementsTree","MeasurementsTree");
84  // MeasurementsTree->Branch("NMeasurements",&NMeasurements,"NMeasurements/I");
85  // BookMeasurements = false;
86  // theRootFile->Add(MeasurementsTree);
87 }
int NDistancemeter1Dim
Definition: NtupleManager.h:70
TClonesArray * CloneOptObject
Definition: NtupleManager.h:48
TClonesArray * CloneTiltmeterMeas
Definition: NtupleManager.h:56
TClonesArray * CloneFitParam
Definition: NtupleManager.h:46
TClonesArray * CloneDistancemeterMeas
Definition: NtupleManager.h:52
TClonesArray * CloneSensor2DMeas
Definition: NtupleManager.h:50
TClonesArray * CloneCopsMeas
Definition: NtupleManager.h:58
TTree * CocoaTree
Definition: NtupleManager.h:42
double Chi2Measurements
Definition: NtupleManager.h:64
TFile * theRootFile
Definition: NtupleManager.h:40
double Chi2CalibratedParameters
Definition: NtupleManager.h:64
TClonesArray * CloneDistancemeter1DimMeas
Definition: NtupleManager.h:54

◆ FillChi2()

void NtupleManager::FillChi2 ( )

Definition at line 119 of file NtupleManager.cc.

References Chi2CalibratedParameters, Chi2Measurements, Model::EntryList(), cuy::ii, Model::MeasurementList(), and NDegreesOfFreedom.

Referenced by Fit::fitNextEvent().

119  {
120  double chi2meas = 0;
121  double chi2cal = 0;
122  ALIint nMeas = 0, nUnk = 0;
123 
124  //----- Calculate the chi2 of measurements
125  std::vector<Measurement*>::const_iterator vmcite;
126  for (vmcite = Model::MeasurementList().begin(); vmcite != Model::MeasurementList().end(); ++vmcite) {
127  for (ALIuint ii = 0; ii < ALIuint((*vmcite)->dim()); ii++) {
128  nMeas++;
129  double c2 = ((*vmcite)->value(ii) - (*vmcite)->valueSimulated(ii)) / (*vmcite)->sigma(ii);
130  chi2meas += c2 * c2;
131  }
132  }
133 
134  //----- Calculate the chi2 of calibrated parameters
135  std::vector<Entry*>::iterator veite;
136  for (veite = Model::EntryList().begin(); veite != Model::EntryList().end(); ++veite) {
137  if ((*veite)->quality() == 2)
138  nUnk++;
139  if ((*veite)->quality() == 1) {
140  // std::cout << " " << (*veite)->valueDisplacementByFitting() << " "
141  // << (*veite)->value << " " << (*veite)->sigma() << std::endl;
142  double c2 = (*veite)->valueDisplacementByFitting() / (*veite)->sigma();
143  chi2cal += c2 * c2;
144  }
145  }
146 
147  Chi2Measurements = chi2meas;
148  Chi2CalibratedParameters = chi2cal;
149  NDegreesOfFreedom = nMeas - nUnk;
150 }
int ALIint
Definition: CocoaGlobals.h:15
ii
Definition: cuy.py:589
double Chi2Measurements
Definition: NtupleManager.h:64
double Chi2CalibratedParameters
Definition: NtupleManager.h:64
static std::vector< Measurement * > & MeasurementList()
Definition: Model.h:88
static std::vector< Entry * > & EntryList()
Definition: Model.h:86
unsigned int ALIuint
Definition: CocoaGlobals.h:17

◆ FillFitParameters()

void NtupleManager::FillFitParameters ( MatrixMeschach AtWAMatrix)

Definition at line 154 of file NtupleManager.cc.

References CloneOptObject, gather_cfg::cout, Model::EntryList(), FitParamA, FitParam::FittedSigma, FitParam::FittedValue, cms::cuda::for(), FittedEntry::getEntryName(), FittedEntry::getName(), FittedEntry::getOptOName(), FittedEntry::getOrder(), FittedEntry::getOrigSigma(), FittedEntry::getOrigValue(), FittedEntry::getQuality(), FittedEntry::getSigma(), FittedEntry::getValue(), if(), cuy::ii, FitParam::InitialSigma, FitParam::InitialValue, MatrixMeschach::Mat(), FitParam::Name, OptObject::Name, NFitParameters, NOptObjects, FitParam::OptObjectIndex, FitParam::Quality, and mathSSE::sqrt().

Referenced by Fit::fitNextEvent().

154  {
155  // double ParValue[1000], ParError[1000];
156  int theMinEntryQuality = 1;
157  int ii = 0;
158  std::vector<Entry*>::const_iterator vecite;
159  for (vecite = Model::EntryList().begin(); vecite != Model::EntryList().end(); ++vecite) {
160  //--- Only for good quality parameters (='unk')
161  if ((*vecite)->quality() >= theMinEntryQuality) {
162  ALIint ipos = (*vecite)->fitPos();
163  FittedEntry* fe = new FittedEntry((*vecite), ipos, sqrt(AtWAMatrix->Mat()->me[ipos][ipos]));
164  // if (!BookFitParameters) {
165  // CocoaTree->Branch("NFitParameters",&NFitParameters,"NFitParameters/I:");
166  // ALIstring partype = fe->getName() + "/D";
167  // FitParametersTree->Branch(fe->getName().c_str(), &ParValue[ii], partype.c_str());
168  // ALIstring parerrname = fe->getName() + "_err";
169  // ALIstring parerrtype = parerrname + "/D";
170  // FitParametersTree->Branch(parerrname.c_str(), &ParError[ii], parerrtype.c_str());
171  // }
172  // ParValue[ii] = fe->getValue();
173  // ParError[ii] = fe->getSigma();
174  std::cout << "EEE " << (*vecite)->ValueDimensionFactor() << " " << (*vecite)->SigmaDimensionFactor() << " "
175  << fe->getOptOName() << " " << fe->getEntryName() << " " << fe->getName() << " " << fe->getOrder()
176  << " " << fe->getQuality() << " " << (*vecite)->type() << " " << std::endl;
177  FitParamA = new ((*CloneFitParam)[ii]) FitParam();
178  FitParamA->Name = fe->getName();
179  if (fe->getQuality() == 1)
180  FitParamA->Quality = "Calibrated";
181  else if (fe->getQuality() == 2)
182  FitParamA->Quality = "Unknown";
183  for (int no = 0; no < NOptObjects; no++) {
184  OptObject* optobj = (OptObject*)CloneOptObject->At(no);
185  if (optobj->Name == fe->getOptOName())
187  }
188  float DF = 1.;
189  if ((*vecite)->type() == "centre" || (*vecite)->type() == "length")
190  DF = 1000.;
191  FitParamA->InitialValue = DF * fe->getOrigValue() * (*vecite)->ValueDimensionFactor();
192  FitParamA->InitialSigma = DF * fe->getOrigSigma() * (*vecite)->SigmaDimensionFactor();
193  FitParamA->FittedValue = DF * fe->getValue() * (*vecite)->ValueDimensionFactor();
194  FitParamA->FittedSigma = DF * fe->getSigma() * (*vecite)->SigmaDimensionFactor();
195  ii++;
196  }
197  }
198  // BookFitParameters = true;
199  NFitParameters = ii;
200  // FitParametersTree->Fill();
201 
202  /*
203  //---------- Loop sets of entries
204  std::vector< FittedEntriesSet* > theFittedEntriesSets;
205  std::vector< FittedEntriesSet* >::const_iterator vfescite;
206  std::vector< FittedEntry* >::const_iterator vfecite;
207  ALIint jj = 1;
208  for( vfescite = theFittedEntriesSets.begin(); vfescite != theFittedEntriesSets.end(); vfescite++) {
209  //---------- Loop entries
210  if( vfescite == theFittedEntriesSets.begin() ) {
211  //----- dump entries names if first set
212  ALIint ii = 0;
213  for( vfecite = ((*vfescite)->FittedEntries()).begin(); vfecite != ((*vfescite)->FittedEntries()).end(); vfecite++) {
214  ALIstring partype = (*vfecite)->getName() + "/D";
215  FitParametersTree->Branch((*vfecite)->getName().c_str(), &ParValue[ii], partype.c_str());
216  ALIstring parerrname = (*vfecite)->getName() + "_err";
217  ALIstring parerrtype = parerrname + "/D";
218  FitParametersTree->Branch(parerrname.c_str(), &ParError[ii], parerrtype.c_str());
219  ii++;
220  }
221  }
222  ALIint ii = 0;
223  for( vfecite = ((*vfescite)->FittedEntries()).begin(); vfecite != ((*vfescite)->FittedEntries()).end(); vfecite++) {
224  ParValue[ii] = (*vfecite)->getValue();
225  ParError[ii] = (*vfecite)->getSigma();
226  ii++;
227  }
228  NFitParameters = ii;
229  FitParametersTree->Fill();
230  jj++;
231  }
232  */
233 }
ALIstring getOptOName() const
Definition: FittedEntry.h:28
ALIint getOrder() const
Definition: FittedEntry.h:35
for(int i=first, nt=offsets[nh];i< nt;i+=gridDim.x *blockDim.x)
TString Name
Definition: NtupleObjects.h:22
int ALIint
Definition: CocoaGlobals.h:15
double FittedValue
Definition: NtupleObjects.h:19
int OptObjectIndex
Definition: NtupleObjects.h:24
TClonesArray * CloneOptObject
Definition: NtupleManager.h:48
double InitialSigma
Definition: NtupleObjects.h:20
ALIstring getName() const
Definition: FittedEntry.h:30
T sqrt(T t)
Definition: SSEVec.h:19
FitParam * FitParamA
Definition: NtupleManager.h:47
TString Name
Definition: NtupleObjects.h:37
ii
Definition: cuy.py:589
const MAT * Mat() const
ALIstring getEntryName() const
Definition: FittedEntry.h:29
ALIint getQuality() const
Definition: FittedEntry.h:36
double InitialValue
Definition: NtupleObjects.h:18
TString Quality
Definition: NtupleObjects.h:23
ALIdouble getOrigSigma() const
Definition: FittedEntry.h:34
ALIdouble getOrigValue() const
Definition: FittedEntry.h:33
static std::vector< Entry * > & EntryList()
Definition: Model.h:86
ALIdouble getValue() const
Definition: FittedEntry.h:31
double FittedSigma
Definition: NtupleObjects.h:21
ALIdouble getSigma() const
Definition: FittedEntry.h:32

◆ FillMeasurements()

void NtupleManager::FillMeasurements ( )

Definition at line 287 of file NtupleManager.cc.

References TiltmeterMeas::AngError, TiltmeterMeas::Angle, CloneOptObject, CopsMeasA, d1, createTree::dd, DistancemeterMeas::DisError, Distancemeter1DimMeas::DisError, DistancemeterMeas::Distance, Distancemeter1DimMeas::Distance, Distancemeter1DimMeasA, DistancemeterMeasA, mps_fire::i, dqmdumpme::last, Model::MeasurementList(), OptObject::Name, Sensor2DMeas::Name, DistancemeterMeas::Name, Distancemeter1DimMeas::Name, TiltmeterMeas::Name, CopsMeas::Name, NCops, NDistancemeter, NDistancemeter1Dim, NOptObjects, NSensor2D, NTiltmeter, Sensor2DMeas::OptObjectIndex, DistancemeterMeas::OptObjectIndex, Distancemeter1DimMeas::OptObjectIndex, TiltmeterMeas::OptObjectIndex, CopsMeas::OptObjectIndex, Sensor2DMeas::PosError, CopsMeas::PosError, Sensor2DMeas::Position, CopsMeas::Position, Sensor2DMeasA, TiltmeterMeas::SimulatedAngle, DistancemeterMeas::SimulatedDistance, Distancemeter1DimMeas::SimulatedDistance, Sensor2DMeas::SimulatedPosition, CopsMeas::SimulatedPosition, contentValuesCheck::ss, TiltmeterMeasA, and groupFilesInBlocks::tt.

Referenced by Fit::fitNextEvent().

287  {
288  //---------- Loop Measurements
289  int ss = 0, dd = 0, d1 = 0, tt = 0, cc = 0;
290  std::vector<Measurement*>::const_iterator vmcite;
291  for (vmcite = Model::MeasurementList().begin(); vmcite != Model::MeasurementList().end(); ++vmcite) {
292  std::vector<ALIstring> optonamelist = (*vmcite)->OptONameList();
293  int last = optonamelist.size() - 1;
294  ALIstring LastOptOName = optonamelist[last];
295  int optoind = -999;
296  for (int no = 0; no < NOptObjects; no++) {
297  OptObject* optobj = (OptObject*)CloneOptObject->At(no);
298  if (optobj->Name == LastOptOName)
299  optoind = no;
300  }
301  //std::cout << "DimSens " << (*vmcite)->type() << " " << (*vmcite)->sigma(0) << " " << LastOptOName << " " << optoind << std::endl;
302  if ((*vmcite)->type() == "SENSOR2D") {
303  Sensor2DMeasA = new ((*CloneSensor2DMeas)[ss]) Sensor2DMeas();
304  Sensor2DMeasA->Name = (*vmcite)->name();
305  Sensor2DMeasA->OptObjectIndex = optoind;
306  for (ALIuint i = 0; i < (*vmcite)->dim(); i++) {
307  Sensor2DMeasA->Position[i] = 1000. * (*vmcite)->value()[i];
308  Sensor2DMeasA->PosError[i] = 1000. * (*vmcite)->sigma()[i];
309  Sensor2DMeasA->SimulatedPosition[i] = 1000. * (*vmcite)->valueSimulated(i);
310  }
311  ss++;
312  }
313  if ((*vmcite)->type() == "DISTANCEMETER") {
314  DistancemeterMeasA = new ((*CloneDistancemeterMeas)[dd]) DistancemeterMeas();
315  DistancemeterMeasA->Name = (*vmcite)->name();
317  DistancemeterMeasA->Distance = 1000. * (*vmcite)->value()[0];
318  DistancemeterMeasA->DisError = 1000. * (*vmcite)->sigma()[0];
319  DistancemeterMeasA->SimulatedDistance = 1000. * (*vmcite)->valueSimulated(0);
320  dd++;
321  }
322  if ((*vmcite)->type() == "DISTANCEMETER1DIM") {
323  Distancemeter1DimMeasA = new ((*CloneDistancemeter1DimMeas)[d1]) Distancemeter1DimMeas();
324  Distancemeter1DimMeasA->Name = (*vmcite)->name();
326  Distancemeter1DimMeasA->Distance = 1000. * (*vmcite)->value()[0];
327  Distancemeter1DimMeasA->DisError = 1000. * (*vmcite)->sigma()[0];
328  Distancemeter1DimMeasA->SimulatedDistance = 1000. * (*vmcite)->valueSimulated(0);
329  d1++;
330  }
331  if ((*vmcite)->type() == "TILTMETER") {
332  TiltmeterMeasA = new ((*CloneTiltmeterMeas)[tt]) TiltmeterMeas();
333  TiltmeterMeasA->Name = (*vmcite)->name();
334  TiltmeterMeasA->OptObjectIndex = optoind;
335  TiltmeterMeasA->Angle = (*vmcite)->value()[0];
336  TiltmeterMeasA->AngError = (*vmcite)->sigma()[0];
337  TiltmeterMeasA->SimulatedAngle = (*vmcite)->valueSimulated(0);
338  tt++;
339  }
340  if ((*vmcite)->type() == "COPS") {
341  CopsMeasA = new ((*CloneCopsMeas)[cc]) CopsMeas();
342  CopsMeasA->Name = (*vmcite)->name();
343  CopsMeasA->OptObjectIndex = optoind;
344  for (ALIuint i = 0; i < (*vmcite)->dim(); i++) {
345  CopsMeasA->Position[i] = 1000. * (*vmcite)->value()[i];
346  CopsMeasA->PosError[i] = 1000. * (*vmcite)->sigma()[i];
347  CopsMeasA->SimulatedPosition[i] = 1000. * (*vmcite)->valueSimulated(i);
348  }
349  cc++;
350  }
351  }
352  NSensor2D = ss;
353  NDistancemeter = dd;
355  NTiltmeter = tt;
356  NCops = cc;
357  // MeasurementsTree->Fill();
358 }
DistancemeterMeas * DistancemeterMeasA
Definition: NtupleManager.h:53
double SimulatedDistance
Definition: NtupleObjects.h:63
int NDistancemeter1Dim
Definition: NtupleManager.h:70
double Position[4]
double SimulatedAngle
Definition: NtupleObjects.h:89
string dd
Definition: createTree.py:154
TClonesArray * CloneOptObject
Definition: NtupleManager.h:48
CopsMeas * CopsMeasA
Definition: NtupleManager.h:59
TString Name
Definition: NtupleObjects.h:51
double PosError[4]
double SimulatedPosition[2]
Definition: NtupleObjects.h:50
TString Name
Definition: NtupleObjects.h:37
Sensor2DMeas * Sensor2DMeasA
Definition: NtupleManager.h:51
std::string ALIstring
Definition: CocoaGlobals.h:9
double PosError[2]
Definition: NtupleObjects.h:49
TString Name
double Position[2]
Definition: NtupleObjects.h:48
static std::vector< Measurement * > & MeasurementList()
Definition: Model.h:88
TiltmeterMeas * TiltmeterMeasA
Definition: NtupleManager.h:57
static constexpr float d1
int OptObjectIndex
double SimulatedPosition[4]
Distancemeter1DimMeas * Distancemeter1DimMeasA
Definition: NtupleManager.h:55
unsigned int ALIuint
Definition: CocoaGlobals.h:17

◆ FillNtupleTree()

void NtupleManager::FillNtupleTree ( )

Definition at line 108 of file NtupleManager.cc.

References CocoaTree.

Referenced by Fit::fitNextEvent().

108 { CocoaTree->Fill(); }
TTree * CocoaTree
Definition: NtupleManager.h:42

◆ FillOptObjects()

void NtupleManager::FillOptObjects ( MatrixMeschach AtWAMatrix)

Definition at line 237 of file NtupleManager.cc.

References OptObject::AnglesGlobal, OptObject::AnglesLocal, OptObject::CentreGlobal, OptObject::CentreLocal, GetGlobalAngles(), mps_fire::i, cuy::ii, OptObject::Name, NOptObjects, OptObjectA, Model::OptOList(), OptObject::Parent, createTree::pp, theta(), OptObject::Type, XCoor, YCoor, and ZCoor.

Referenced by Fit::fitNextEvent().

237  {
238  int ii = 0;
239  std::vector<OpticalObject*>::const_iterator vecobj;
240  for (vecobj = Model::OptOList().begin(); vecobj != Model::OptOList().end(); ++vecobj) {
241  OptObjectA = new ((*CloneOptObject)[ii]) OptObject();
242 
243  OptObjectA->Name = (*vecobj)->name();
244  OptObjectA->Type = (*vecobj)->type();
245 
246  if (!(*vecobj)->parent()) {
247  OptObjectA->Parent = ii;
248  ii++;
249  continue;
250  }
251 
252  int pp = 0;
253  std::vector<OpticalObject*>::const_iterator vecobj2;
254  for (vecobj2 = Model::OptOList().begin(); vecobj2 != Model::OptOList().end(); ++vecobj2) {
255  if ((*vecobj2)->name() == (*vecobj)->parent()->name()) {
256  OptObjectA->Parent = pp;
257  continue;
258  }
259  pp++;
260  }
261 
262  OptObjectA->CentreGlobal[0] = 1000. * (*vecobj)->centreGlobal().x();
263  OptObjectA->CentreGlobal[1] = 1000. * (*vecobj)->centreGlobal().y();
264  OptObjectA->CentreGlobal[2] = 1000. * (*vecobj)->centreGlobal().z();
265 
266  OptObjectA->CentreLocal[0] = 1000. * (*vecobj)->centreLocal().x();
267  OptObjectA->CentreLocal[1] = 1000. * (*vecobj)->centreLocal().y();
268  OptObjectA->CentreLocal[2] = 1000. * (*vecobj)->centreLocal().z();
269 
270  OptObjectA->AnglesLocal[0] = (*vecobj)->getEntryRMangle(XCoor);
271  OptObjectA->AnglesLocal[1] = (*vecobj)->getEntryRMangle(YCoor);
272  OptObjectA->AnglesLocal[2] = (*vecobj)->getEntryRMangle(ZCoor);
273 
274  double theta[3];
275  GetGlobalAngles((*vecobj)->rmGlob(), theta);
276  for (int i = 0; i < 3; i++)
278 
279  ii++;
280  }
281 
282  NOptObjects = ii;
283 }
double AnglesLocal[3]
Definition: NtupleObjects.h:36
double AnglesGlobal[3]
Definition: NtupleObjects.h:34
void GetGlobalAngles(const CLHEP::HepRotation &rmGlob, double *theta)
static std::vector< OpticalObject * > & OptOList()
Definition: Model.h:84
double CentreLocal[3]
Definition: NtupleObjects.h:35
TString Name
Definition: NtupleObjects.h:37
OptObject * OptObjectA
Definition: NtupleManager.h:49
ii
Definition: cuy.py:589
Geom::Theta< T > theta() const
double CentreGlobal[3]
Definition: NtupleObjects.h:33
TString Type
Definition: NtupleObjects.h:38

◆ GetGlobalAngles()

void NtupleManager::GetGlobalAngles ( const CLHEP::HepRotation &  rmGlob,
double *  theta 
)
private

Definition at line 362 of file NtupleManager.cc.

References alpha, HLT_2022v12_cff::beta, funct::cos(), gather_cfg::cout, MillePedeFileConverter_cfg::e, CustomPhysics_cfi::gamma, M_PI, funct::sin(), theta(), geometryCSVtoXML::xx, geometryCSVtoXML::xy, geometryCSVtoXML::xz, geometryCSVtoXML::yy, geometryCSVtoXML::yz, and geometryCSVtoXML::zz.

Referenced by FillOptObjects().

362  {
363  double xx = rmGlob.xx();
364  if (fabs(xx) < 1.e-08)
365  xx = 0.;
366  double xy = rmGlob.xy();
367  if (fabs(xy) < 1.e-08)
368  xy = 0.;
369  double xz = rmGlob.xz();
370  if (fabs(xz) < 1.e-08)
371  xz = 0.;
372  double yx = rmGlob.yx();
373  if (fabs(yx) < 1.e-08)
374  yx = 0.;
375  double yy = rmGlob.yy();
376  if (fabs(yy) < 1.e-08)
377  yy = 0.;
378  double yz = rmGlob.yz();
379  if (fabs(yz) < 1.e-08)
380  yz = 0.;
381  double zx = rmGlob.zx();
382  if (fabs(zx) < 1.e-08)
383  zx = 0.;
384  double zy = rmGlob.zy();
385  if (fabs(zy) < 1.e-08)
386  zy = 0.;
387  double zz = rmGlob.zz();
388  if (fabs(zz) < 1.e-08)
389  zz = 0.;
390 
391  double beta = asin(-zx);
392 
393  double alpha, gamma;
394  if (fabs(zx) != 1.) {
395  double sinalpha = zy / cos(beta);
396  double cosalpha = zz / cos(beta);
397  if (cosalpha >= 0)
398  alpha = asin(sinalpha);
399  else
400  alpha = M_PI - asin(sinalpha);
401  if (alpha > M_PI)
402  alpha -= 2 * M_PI;
403 
404  double singamma = yx / cos(beta);
405  double cosgamma = xx / cos(beta);
406  if (cosgamma >= 0)
407  gamma = asin(singamma);
408  else
409  gamma = M_PI - asin(singamma);
410  if (gamma > M_PI)
411  gamma -= 2 * M_PI;
412 
413  } else {
414  alpha = 0.;
415 
416  double singamma = yz / sin(beta);
417  double cosgamma = yy;
418  if (cosgamma >= 0)
419  gamma = asin(singamma);
420  else
421  gamma = M_PI - asin(singamma);
422  if (gamma > M_PI)
423  gamma -= 2 * M_PI;
424  }
425 
426  int GotGlobalAngles = 0;
427  if (fabs(xy - (sin(alpha) * sin(beta) * cos(gamma) - sin(gamma) * cos(alpha))) > 1.e-08)
428  GotGlobalAngles += 1;
429  if (fabs(xz - (cos(alpha) * sin(beta) * cos(gamma) + sin(gamma) * sin(alpha))) > 1.e-08)
430  GotGlobalAngles += 10;
431  if (fabs(yy - (sin(alpha) * sin(beta) * sin(gamma) + cos(gamma) * cos(alpha))) > 1.e-08)
432  GotGlobalAngles += 100;
433  if (fabs(yz - (cos(alpha) * sin(beta) * sin(gamma) - cos(gamma) * sin(alpha))) > 1.e-08)
434  GotGlobalAngles += 1000;
435  if (GotGlobalAngles > 0)
436  std::cout << "NtupleManager Warning: cannot get global rotation: " << GotGlobalAngles << std::endl;
437 
438  theta[0] = alpha;
439  theta[1] = beta;
440  theta[2] = gamma;
441 }
float alpha
Definition: AMPTWrapper.h:105
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
#define M_PI
Geom::Theta< T > theta() const

◆ getInstance()

NtupleManager * NtupleManager::getInstance ( )
static

Definition at line 29 of file NtupleManager.cc.

References instance, and NtupleManager().

Referenced by Fit::fitNextEvent(), and Fit::startFit().

29  {
30  if (!instance) {
31  instance = new NtupleManager;
32  }
33  return instance;
34 }
static NtupleManager * instance
Definition: NtupleManager.h:36

◆ InitNtuple()

void NtupleManager::InitNtuple ( )

Definition at line 91 of file NtupleManager.cc.

References Chi2CalibratedParameters, Chi2Measurements, CloneFitParam, NCops, NDegreesOfFreedom, NDistancemeter, NDistancemeter1Dim, NFitParameters, NOptObjects, NSensor2D, and NTiltmeter.

Referenced by Fit::fitNextEvent().

91  {
92  CloneFitParam->Clear();
93 
94  Chi2Measurements = 0.;
97  NFitParameters = 0;
98  NOptObjects = 0;
99  NSensor2D = 0;
100  NDistancemeter = 0;
101  NDistancemeter1Dim = 0;
102  NTiltmeter = 0;
103  NCops = 0;
104 }
int NDistancemeter1Dim
Definition: NtupleManager.h:70
TClonesArray * CloneFitParam
Definition: NtupleManager.h:46
double Chi2Measurements
Definition: NtupleManager.h:64
double Chi2CalibratedParameters
Definition: NtupleManager.h:64

◆ WriteNtuple()

void NtupleManager::WriteNtuple ( )

Definition at line 112 of file NtupleManager.cc.

References theRootFile.

Referenced by Fit::startFit().

112  {
113  theRootFile->Write();
114  theRootFile->Close();
115 }
TFile * theRootFile
Definition: NtupleManager.h:40

Member Data Documentation

◆ Chi2CalibratedParameters

double NtupleManager::Chi2CalibratedParameters
private

Definition at line 64 of file NtupleManager.h.

Referenced by BookNtuple(), FillChi2(), and InitNtuple().

◆ Chi2Measurements

double NtupleManager::Chi2Measurements
private

Definition at line 64 of file NtupleManager.h.

Referenced by BookNtuple(), FillChi2(), and InitNtuple().

◆ CloneCopsMeas

TClonesArray* NtupleManager::CloneCopsMeas
private

Definition at line 58 of file NtupleManager.h.

Referenced by BookNtuple().

◆ CloneDistancemeter1DimMeas

TClonesArray* NtupleManager::CloneDistancemeter1DimMeas
private

Definition at line 54 of file NtupleManager.h.

Referenced by BookNtuple().

◆ CloneDistancemeterMeas

TClonesArray* NtupleManager::CloneDistancemeterMeas
private

Definition at line 52 of file NtupleManager.h.

Referenced by BookNtuple().

◆ CloneFitParam

TClonesArray* NtupleManager::CloneFitParam
private

Definition at line 46 of file NtupleManager.h.

Referenced by BookNtuple(), and InitNtuple().

◆ CloneOptObject

TClonesArray* NtupleManager::CloneOptObject
private

Definition at line 48 of file NtupleManager.h.

Referenced by BookNtuple(), FillFitParameters(), and FillMeasurements().

◆ CloneSensor2DMeas

TClonesArray* NtupleManager::CloneSensor2DMeas
private

Definition at line 50 of file NtupleManager.h.

Referenced by BookNtuple().

◆ CloneTiltmeterMeas

TClonesArray* NtupleManager::CloneTiltmeterMeas
private

Definition at line 56 of file NtupleManager.h.

Referenced by BookNtuple().

◆ CocoaTree

TTree* NtupleManager::CocoaTree
private

Definition at line 42 of file NtupleManager.h.

Referenced by BookNtuple(), and FillNtupleTree().

◆ CopsMeasA

CopsMeas* NtupleManager::CopsMeasA
private

Definition at line 59 of file NtupleManager.h.

Referenced by FillMeasurements().

◆ Distancemeter1DimMeasA

Distancemeter1DimMeas* NtupleManager::Distancemeter1DimMeasA
private

Definition at line 55 of file NtupleManager.h.

Referenced by FillMeasurements().

◆ DistancemeterMeasA

DistancemeterMeas* NtupleManager::DistancemeterMeasA
private

Definition at line 53 of file NtupleManager.h.

Referenced by FillMeasurements().

◆ FitParamA

FitParam* NtupleManager::FitParamA
private

Definition at line 47 of file NtupleManager.h.

Referenced by FillFitParameters().

◆ instance

NtupleManager * NtupleManager::instance = nullptr
staticprivate

Definition at line 36 of file NtupleManager.h.

Referenced by getInstance(), and production_tasks.Task::getname().

◆ NCops

int NtupleManager::NCops
private

Definition at line 72 of file NtupleManager.h.

Referenced by BookNtuple(), FillMeasurements(), and InitNtuple().

◆ NDegreesOfFreedom

int NtupleManager::NDegreesOfFreedom
private

Definition at line 65 of file NtupleManager.h.

Referenced by BookNtuple(), FillChi2(), and InitNtuple().

◆ NDistancemeter

int NtupleManager::NDistancemeter
private

Definition at line 69 of file NtupleManager.h.

Referenced by BookNtuple(), FillMeasurements(), and InitNtuple().

◆ NDistancemeter1Dim

int NtupleManager::NDistancemeter1Dim
private

Definition at line 70 of file NtupleManager.h.

Referenced by BookNtuple(), FillMeasurements(), and InitNtuple().

◆ NFitParameters

int NtupleManager::NFitParameters
private

Definition at line 66 of file NtupleManager.h.

Referenced by BookNtuple(), FillFitParameters(), and InitNtuple().

◆ NOptObjects

int NtupleManager::NOptObjects
private

◆ NSensor2D

int NtupleManager::NSensor2D
private

Definition at line 68 of file NtupleManager.h.

Referenced by BookNtuple(), FillMeasurements(), and InitNtuple().

◆ NTiltmeter

int NtupleManager::NTiltmeter
private

Definition at line 71 of file NtupleManager.h.

Referenced by BookNtuple(), FillMeasurements(), and InitNtuple().

◆ OptObjectA

OptObject* NtupleManager::OptObjectA
private

Definition at line 49 of file NtupleManager.h.

Referenced by FillOptObjects().

◆ Sensor2DMeasA

Sensor2DMeas* NtupleManager::Sensor2DMeasA
private

Definition at line 51 of file NtupleManager.h.

Referenced by FillMeasurements().

◆ theRootFile

TFile* NtupleManager::theRootFile
private

Definition at line 40 of file NtupleManager.h.

Referenced by BookNtuple(), and WriteNtuple().

◆ TiltmeterMeasA

TiltmeterMeas* NtupleManager::TiltmeterMeasA
private

Definition at line 57 of file NtupleManager.h.

Referenced by FillMeasurements().