10#include <klm/calibration/KLMDisplacementGenerator.h>
13#include <klm/eklm/geometry/AlignmentChecker.h>
14#include <klm/eklm/geometry/GeometryData.h>
17#include <framework/gearbox/Unit.h>
27 m_GeoDat(&(EKLM::GeometryData::Instance())),
41 int iSection, iLayer, iSector, iPlane, iSegment, segment;
47 alignment->setModuleAlignment(module, &alignmentData);
51 iSection, iLayer, iSector, iPlane, iSegment);
61 double deltaU,
double deltaV,
double deltaGamma)
66 int iSection, iLayer, iSector;
73 alignment.setModuleAlignment(module, &moduleAlignment);
81 bool displaceModule,
bool displaceSegment,
bool moduleSameDisplacement,
82 bool moduleZeroDeltaU,
bool moduleZeroDeltaV,
bool moduleZeroDeltaGamma)
84 const double moduleMaxDeltaU = 5. *
Unit::cm;
85 const double moduleMinDeltaU = -5. *
Unit::cm;
86 const double moduleMaxDeltaV = 5. *
Unit::cm;
87 const double moduleMinDeltaV = -5. *
Unit::cm;
88 const double moduleMaxDeltaGamma = 0.02 *
Unit::rad;
89 const double moduleMinDeltaGamma = -0.02 *
Unit::rad;
90 const double segmentMaxDeltaU = 1. *
Unit::cm;
91 const double segmentMinDeltaU = -1. *
Unit::cm;
92 const double segmentMaxDeltaV = 0.2 *
Unit::cm;
93 const double segmentMinDeltaV = -0.2 *
Unit::cm;
94 const double segmentMaxDeltaGamma = 0.003 *
Unit::rad;
95 const double segmentMinDeltaGamma = -0.003 *
Unit::rad;
98 int iSection, iLayer, iSector, iPlane, iSegment, segment;
105 if (moduleSameDisplacement) {
106 if (iSection > 1 || iLayer > 1 || iSector > 1) {
109 alignment->getModuleAlignment(module));
111 alignment->setModuleAlignment(module, alignmentData);
119 iSection, iLayer, iSector, iPlane, iSegment);
127 if (displaceModule) {
129 if (moduleZeroDeltaU)
132 d = gRandom->Uniform(moduleMinDeltaU, moduleMaxDeltaU);
134 if (moduleZeroDeltaV)
137 d = gRandom->Uniform(moduleMinDeltaV, moduleMaxDeltaV);
139 if (moduleZeroDeltaGamma)
142 d = gRandom->Uniform(moduleMinDeltaGamma, moduleMaxDeltaGamma);
145 iSection, iLayer, iSector, &moduleAlignmentData));
147 alignment->setModuleAlignment(module, &moduleAlignmentData);
149 if (displaceSegment) {
155 gRandom->Uniform(segmentMinDeltaU, segmentMaxDeltaU));
157 gRandom->Uniform(segmentMinDeltaV, segmentMaxDeltaV));
159 gRandom->Uniform(segmentMinDeltaGamma, segmentMaxDeltaGamma));
161 iSection, iLayer, iSector, iPlane, iSegment,
162 &moduleAlignmentData, &segmentAlignmentData,
false));
164 iSection, iLayer, iSector, iPlane, iSegment);
166 segment, &segmentAlignmentData);
174 iSection, iLayer, iSector, iPlane, iSegment,
175 &moduleAlignmentData, &segmentAlignmentData,
false))
187 const char* inputFile)
190 int i, n, iSection, iLayer, iSector, iPlane, iSegment, segment, param;
193 TTree* tEKLMModule, *tEKLMSegment;
196 f =
new TFile(inputFile);
197 tEKLMModule = (TTree*)f->Get(
"eklm_module");
198 tEKLMModule->SetBranchAddress(
"section", &iSection);
199 tEKLMModule->SetBranchAddress(
"layer", &iLayer);
200 tEKLMModule->SetBranchAddress(
"sector", &iSector);
201 tEKLMModule->SetBranchAddress(
"param", ¶m);
202 tEKLMModule->SetBranchAddress(
"value", &value);
203 tEKLMSegment = (TTree*)f->Get(
"eklm_segment");
204 tEKLMSegment->SetBranchAddress(
"section", &iSection);
205 tEKLMSegment->SetBranchAddress(
"layer", &iLayer);
206 tEKLMSegment->SetBranchAddress(
"sector", &iSector);
207 tEKLMSegment->SetBranchAddress(
"plane", &iPlane);
208 tEKLMSegment->SetBranchAddress(
"segment", &iSegment);
209 tEKLMSegment->SetBranchAddress(
"param", ¶m);
210 tEKLMSegment->SetBranchAddress(
"value", &value);
211 n = tEKLMModule->GetEntries();
212 for (i = 0; i < n; i++) {
213 tEKLMModule->GetEntry(i);
216 alignment->getModuleAlignment(module));
219 alignmentData->setDeltaU(value);
222 alignmentData->setDeltaV(value);
225 alignmentData->setDeltaGamma(value);
229 n = tEKLMSegment->GetEntries();
230 for (i = 0; i < n; i++) {
231 tEKLMSegment->GetEntry(i);
233 iSection, iLayer, iSector, iPlane, iSegment);
238 alignmentData->setDeltaV(value);
241 alignmentData->setDeltaGamma(value);
249 const int nPoints = 1000;
250 const float maxDeltaU = 5. *
Unit::cm;
251 const float minDeltaU = -5. *
Unit::cm;
252 const float maxDeltaV = 5. *
Unit::cm;
253 const float minDeltaV = -5. *
Unit::cm;
254 const float maxDeltaGamma = 0.02 *
Unit::rad;
255 const float minDeltaGamma = -0.02 *
Unit::rad;
256 float deltaU, deltaV, deltaGamma;
257 int i, alignmentStatus, iSection, iLayer, iSector;
264 t =
new TTree(
"module",
"");
265 t->Branch(
"deltaU", &deltaU,
"deltaU/F");
266 t->Branch(
"deltaV", &deltaV,
"deltaV/F");
267 t->Branch(
"deltaGamma", &deltaGamma,
"deltaGamma/F");
268 t->Branch(
"status", &alignmentStatus,
"status/I");
269 for (i = 0; i < nPoints; i++) {
271 deltaU = gRandom->Uniform(minDeltaU, maxDeltaU);
272 deltaV = gRandom->Uniform(minDeltaV, maxDeltaV);
273 deltaGamma = gRandom->Uniform(minDeltaGamma, maxDeltaGamma);
282 alignment.setModuleAlignment(module, &alignmentDataRandom);
286 if (alignmentChecker.
checkAlignment(&alignment, &segmentAlignment))
297 const int nPoints = 1000;
298 const float maxDeltaU = 9. *
Unit::cm;
299 const float minDeltaU = -4. *
Unit::cm;
300 const float maxDeltaV = 0.2 *
Unit::cm;
301 const float minDeltaV = -0.2 *
Unit::cm;
302 const float maxDeltaGamma = 0.003 *
Unit::rad;
303 const float minDeltaGamma = -0.003 *
Unit::rad;
304 float deltaU, deltaV, deltaGamma;
305 int i, alignmentStatus;
306 int iSection, iLayer, iSector, jPlane, jSegment, segment;
313 t =
new TTree(
"segment",
"");
314 t->Branch(
"plane", &jPlane,
"plane/I");
315 t->Branch(
"segment", &jSegment,
"segment/I");
316 t->Branch(
"deltaU", &deltaU,
"deltaU/F");
317 t->Branch(
"deltaV", &deltaV,
"deltaV/F");
318 t->Branch(
"deltaGamma", &deltaGamma,
"deltaGamma/F");
319 t->Branch(
"status", &alignmentStatus,
"status/I");
321 printf(
"Plane %d\n", jPlane);
323 printf(
"Segment %d\n", jSegment);
324 for (i = 0; i < nPoints; i++) {
326 deltaU = gRandom->Uniform(minDeltaU, maxDeltaU);
327 deltaV = gRandom->Uniform(minDeltaV, maxDeltaV);
328 deltaGamma = gRandom->Uniform(minDeltaGamma, maxDeltaGamma);
337 iSection, iLayer, iSector, jPlane, jSegment);
339 segment, &alignmentDataRandom);
343 if (alignmentChecker.
checkAlignment(&alignment, &segmentAlignment))
357 f =
new TFile(outputFile,
"recreate");
365 const char* outputFile)
367 int iSection, iLayer, iSector, iPlane, iSegment, segment, param;
371 TTree* tEKLMModule, *tEKLMSegment;
372 f =
new TFile(outputFile,
"recreate");
373 tEKLMModule =
new TTree(
"eklm_module",
"");
374 tEKLMModule->Branch(
"section", &iSection,
"section/I");
375 tEKLMModule->Branch(
"layer", &iLayer,
"layer/I");
376 tEKLMModule->Branch(
"sector", &iSector,
"sector/I");
377 tEKLMModule->Branch(
"param", ¶m,
"param/I");
378 tEKLMModule->Branch(
"value", &value,
"value/F");
379 tEKLMSegment =
new TTree(
"eklm_segment",
"");
380 tEKLMSegment->Branch(
"section", &iSection,
"section/I");
381 tEKLMSegment->Branch(
"layer", &iLayer,
"layer/I");
382 tEKLMSegment->Branch(
"sector", &iSector,
"sector/I");
383 tEKLMSegment->Branch(
"plane", &iPlane,
"plane/I");
384 tEKLMSegment->Branch(
"segment", &iSegment,
"segment/I");
385 tEKLMSegment->Branch(
"param", ¶m,
"param/I");
386 tEKLMSegment->Branch(
"value", &value,
"value/F");
393 alignment->getModuleAlignment(module));
395 value = alignmentData->getDeltaU();
400 value = alignmentData->getDeltaV();
405 value = alignmentData->getDeltaGamma();
410 iSection, iLayer, iSector, iPlane, iSegment);
414 value = alignmentData->getDeltaV();
415 tEKLMSegment->Fill();
419 value = alignmentData->getDeltaGamma();
420 tEKLMSegment->Fill();
427 tEKLMModule->Write();
428 tEKLMSegment->Write();
Class to store EKLM alignment data in the database.
int segmentNumber(int section, int layer, int sector, int plane, int segment) const
Get segment number.
int getNPlanes() const
Get number of planes.
int getNSections() const
Get number of sections.
int getNDetectorLayers(int section) const
Get number of detector layers.
int getNSegments() const
Get number of segments.
int getNSectors() const
Get number of sectors.
Class to store EKLM alignment data in the database.
const KLMAlignmentData * getSegmentAlignment(EKLMSegmentNumber segment) const
Get segment alignment data.
void setSegmentAlignment(EKLMSegmentNumber segment, KLMAlignmentData *dat)
Set segment alignment data.
Class for EKLM alignment checking.
bool checkAlignment(const EKLMAlignment *alignment, const EKLMSegmentAlignment *segmentAlignment) const
Check alignment.
bool checkSectorAlignment(int section, int layer, int sector, const KLMAlignmentData *sectorAlignment) const
Check sector alignment.
bool checkSegmentAlignment(int section, int layer, int sector, int plane, int segment, const KLMAlignmentData *sectorAlignment, const KLMAlignmentData *segmentAlignment, bool calledFromSectorCheck) const
Check segment alignment.
void setDeltaGamma(float deltaGamma)
Set rotation in alpha.
void setDeltaU(float deltaU)
Set shift in U.
void setDeltaV(float deltaV)
Set shift in V.
const EKLM::GeometryData * m_GeoDat
Geometry data.
void studyAlignmentLimits(const char *outputFile)
Generate random displacements and check if they are correct (no overlaps).
void readDisplacementFromROOTFile(EKLMAlignment *alignment, EKLMSegmentAlignment *segmentAlignment, const char *inputFile)
Read displacement from ROOT file.
void studyModuleAlignmentLimits(TFile *f)
Generate random module displacements and check if they are correct (no overlaps).
const KLMElementNumbers * m_ElementNumbers
Element numbers.
const EKLMElementNumbers * m_eklmElementNumbers
EKLM element numbers.
~KLMDisplacementGenerator()
Destructor.
void generateRandomDisplacement(EKLMAlignment *alignment, EKLMSegmentAlignment *segmentAlignment, bool displaceModule, bool displaceSegment, bool moduleSameDisplacement=false, bool moduleZeroDeltaU=false, bool moduleZeroDeltaV=false, bool moduleZeroDeltaGamma=false)
Generation of random displacements.
void fillZeroDisplacements(EKLMAlignment *alignment, EKLMSegmentAlignment *segmentAlignment)
Fill EKLMAlignment with zero displacements.
void saveDisplacement(EKLMAlignment *alignment, EKLMSegmentAlignment *segmentAlignment, const char *outputFile)
Save displacements to a ROOT file.
void generateFixedModuleDisplacement(double deltaU, double deltaV, double deltaGamma)
Generation of fixed module displacements.
void studySegmentAlignmentLimits(TFile *f)
Generate random segment displacements and check if they are correct (no overlaps).
KLMDisplacementGenerator()
Constructor.
KLMModuleNumber moduleNumberEKLM(int section, int sector, int layer) const
Get module number for EKLM.
static const double rad
Standard of [angle].
static const double cm
Standard units with the value = 1.
Abstract base class for different kinds of events.