Belle II Software development
KLMDisplacementGenerator.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9/* Own header. */
10#include <klm/calibration/KLMDisplacementGenerator.h>
11
12/* KLM headers. */
13#include <klm/eklm/geometry/AlignmentChecker.h>
14#include <klm/eklm/geometry/GeometryData.h>
15
16/* Basf2 headers. */
17#include <framework/gearbox/Unit.h>
18
19/* ROOT headers. */
20#include <TFile.h>
21#include <TRandom.h>
22#include <TTree.h>
23
24using namespace Belle2;
25
27 m_GeoDat(&(EKLM::GeometryData::Instance())),
28 m_ElementNumbers(&(KLMElementNumbers::Instance())),
29 m_eklmElementNumbers(&(EKLMElementNumbers::Instance()))
30{
31}
32
34{
35}
36
38 EKLMAlignment* alignment, EKLMSegmentAlignment* segmentAlignment)
39{
40 KLMAlignmentData alignmentData(0, 0, 0, 0, 0, 0);
41 int iSection, iLayer, iSector, iPlane, iSegment, segment;
42 for (iSection = 1; iSection <= m_GeoDat->getNSections(); iSection++) {
43 for (iLayer = 1; iLayer <= m_GeoDat->getNDetectorLayers(iSection);
44 iLayer++) {
45 for (iSector = 1; iSector <= m_GeoDat->getNSectors(); iSector++) {
46 int module = m_ElementNumbers->moduleNumberEKLM(iSection, iSector, iLayer);
47 alignment->setModuleAlignment(module, &alignmentData);
48 for (iPlane = 1; iPlane <= m_GeoDat->getNPlanes(); iPlane++) {
49 for (iSegment = 1; iSegment <= m_GeoDat->getNSegments(); iSegment++) {
51 iSection, iLayer, iSector, iPlane, iSegment);
52 segmentAlignment->setSegmentAlignment(segment, &alignmentData);
53 }
54 }
55 }
56 }
57 }
58}
59
61 double deltaU, double deltaV, double deltaGamma)
62{
63 EKLMAlignment alignment;
64 EKLMSegmentAlignment segmentAlignment;
65 KLMAlignmentData moduleAlignment(deltaU, deltaV, 0, 0, 0, deltaGamma);
66 int iSection, iLayer, iSector;
67 fillZeroDisplacements(&alignment, &segmentAlignment);
68 for (iSection = 1; iSection <= m_GeoDat->getNSections(); iSection++) {
69 for (iLayer = 1; iLayer <= m_GeoDat->getNDetectorLayers(iSection);
70 iLayer++) {
71 for (iSector = 1; iSector <= m_GeoDat->getNSectors(); iSector++) {
72 int module = m_ElementNumbers->moduleNumberEKLM(iSection, iSector, iLayer);
73 alignment.setModuleAlignment(module, &moduleAlignment);
74 }
75 }
76 }
77}
78
80 EKLMAlignment* alignment, EKLMSegmentAlignment* segmentAlignment,
81 bool displaceModule, bool displaceSegment, bool moduleSameDisplacement,
82 bool moduleZeroDeltaU, bool moduleZeroDeltaV, bool moduleZeroDeltaGamma)
83{
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;
96 KLMAlignmentData moduleAlignmentData, segmentAlignmentData, *alignmentData;
97 EKLM::AlignmentChecker alignmentChecker(false);
98 int iSection, iLayer, iSector, iPlane, iSegment, segment;
99 double d;
100 fillZeroDisplacements(alignment, segmentAlignment);
101 for (iSection = 1; iSection <= m_GeoDat->getNSections(); iSection++) {
102 for (iLayer = 1; iLayer <= m_GeoDat->getNDetectorLayers(iSection);
103 iLayer++) {
104 for (iSector = 1; iSector <= m_GeoDat->getNSectors(); iSector++) {
105 if (moduleSameDisplacement) {
106 if (iSection > 1 || iLayer > 1 || iSector > 1) {
107 int module = m_ElementNumbers->moduleNumberEKLM(1, 1, 1);
108 alignmentData = const_cast<KLMAlignmentData*>(
109 alignment->getModuleAlignment(module));
110 module = m_ElementNumbers->moduleNumberEKLM(iSection, iSector, iLayer);
111 alignment->setModuleAlignment(module, alignmentData);
112 for (iPlane = 1; iPlane <= m_GeoDat->getNPlanes(); iPlane++) {
113 for (iSegment = 1; iSegment <= m_GeoDat->getNSegments();
114 iSegment++) {
115 segment = m_eklmElementNumbers->segmentNumber(1, 1, 1, iPlane, iSegment);
116 alignmentData = const_cast<KLMAlignmentData*>(
117 segmentAlignment->getSegmentAlignment(segment));
119 iSection, iLayer, iSector, iPlane, iSegment);
120 segmentAlignment->setSegmentAlignment(segment, alignmentData);
121 }
122 }
123 continue;
124 }
125 }
126module:
127 if (displaceModule) {
128 do {
129 if (moduleZeroDeltaU)
130 d = 0;
131 else
132 d = gRandom->Uniform(moduleMinDeltaU, moduleMaxDeltaU);
133 moduleAlignmentData.setDeltaU(d);
134 if (moduleZeroDeltaV)
135 d = 0;
136 else
137 d = gRandom->Uniform(moduleMinDeltaV, moduleMaxDeltaV);
138 moduleAlignmentData.setDeltaV(d);
139 if (moduleZeroDeltaGamma)
140 d = 0;
141 else
142 d = gRandom->Uniform(moduleMinDeltaGamma, moduleMaxDeltaGamma);
143 moduleAlignmentData.setDeltaGamma(d);
144 } while (!alignmentChecker.checkSectorAlignment(
145 iSection, iLayer, iSector, &moduleAlignmentData));
146 int module = m_ElementNumbers->moduleNumberEKLM(iSection, iSector, iLayer);
147 alignment->setModuleAlignment(module, &moduleAlignmentData);
148 }
149 if (displaceSegment) {
150 for (iPlane = 1; iPlane <= m_GeoDat->getNPlanes(); iPlane++) {
151 for (iSegment = 1; iSegment <= m_GeoDat->getNSegments();
152 iSegment++) {
153 do {
154 segmentAlignmentData.setDeltaU(
155 gRandom->Uniform(segmentMinDeltaU, segmentMaxDeltaU));
156 segmentAlignmentData.setDeltaV(
157 gRandom->Uniform(segmentMinDeltaV, segmentMaxDeltaV));
158 segmentAlignmentData.setDeltaGamma(
159 gRandom->Uniform(segmentMinDeltaGamma, segmentMaxDeltaGamma));
160 } while (!alignmentChecker.checkSegmentAlignment(
161 iSection, iLayer, iSector, iPlane, iSegment,
162 &moduleAlignmentData, &segmentAlignmentData, false));
164 iSection, iLayer, iSector, iPlane, iSegment);
165 segmentAlignment->setSegmentAlignment(
166 segment, &segmentAlignmentData);
167 }
168 }
169 } else {
170 for (iPlane = 1; iPlane <= m_GeoDat->getNPlanes(); iPlane++) {
171 for (iSegment = 1; iSegment <= m_GeoDat->getNSegments();
172 iSegment++) {
173 if (!alignmentChecker.checkSegmentAlignment(
174 iSection, iLayer, iSector, iPlane, iSegment,
175 &moduleAlignmentData, &segmentAlignmentData, false))
176 goto module;
177 }
178 }
179 }
180 }
181 }
182 }
183}
184
186 EKLMAlignment* alignment, EKLMSegmentAlignment* segmentAlignment,
187 const char* inputFile)
188{
189 /* cppcheck-suppress variableScope */
190 int i, n, iSection, iLayer, iSector, iPlane, iSegment, segment, param;
191 float value;
192 TFile* f;
193 TTree* tEKLMModule, *tEKLMSegment;
194 KLMAlignmentData* alignmentData;
195 fillZeroDisplacements(alignment, segmentAlignment);
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", &param);
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", &param);
210 tEKLMSegment->SetBranchAddress("value", &value);
211 n = tEKLMModule->GetEntries();
212 for (i = 0; i < n; i++) {
213 tEKLMModule->GetEntry(i);
214 int module = m_ElementNumbers->moduleNumberEKLM(iSection, iSector, iLayer);
215 alignmentData = const_cast<KLMAlignmentData*>(
216 alignment->getModuleAlignment(module));
217 switch (param) {
218 case 1:
219 alignmentData->setDeltaU(value);
220 break;
221 case 2:
222 alignmentData->setDeltaV(value);
223 break;
224 case 3:
225 alignmentData->setDeltaGamma(value);
226 break;
227 }
228 }
229 n = tEKLMSegment->GetEntries();
230 for (i = 0; i < n; i++) {
231 tEKLMSegment->GetEntry(i);
233 iSection, iLayer, iSector, iPlane, iSegment);
234 alignmentData = const_cast<KLMAlignmentData*>(
235 segmentAlignment->getSegmentAlignment(segment));
236 switch (param) {
237 case 1:
238 alignmentData->setDeltaV(value);
239 break;
240 case 2:
241 alignmentData->setDeltaGamma(value);
242 break;
243 }
244 }
245}
246
248{
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;
258 EKLMAlignment alignment;
259 EKLMSegmentAlignment segmentAlignment;
260 KLMAlignmentData alignmentDataRandom;
261 EKLM::AlignmentChecker alignmentChecker(false);
262 TTree* t;
263 f->cd();
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++) {
270 fillZeroDisplacements(&alignment, &segmentAlignment);
271 deltaU = gRandom->Uniform(minDeltaU, maxDeltaU);
272 deltaV = gRandom->Uniform(minDeltaV, maxDeltaV);
273 deltaGamma = gRandom->Uniform(minDeltaGamma, maxDeltaGamma);
274 alignmentDataRandom.setDeltaU(deltaU);
275 alignmentDataRandom.setDeltaV(deltaV);
276 alignmentDataRandom.setDeltaGamma(deltaGamma);
277 for (iSection = 1; iSection <= m_GeoDat->getNSections(); iSection++) {
278 for (iLayer = 1; iLayer <= m_GeoDat->getNDetectorLayers(iSection);
279 iLayer++) {
280 for (iSector = 1; iSector <= m_GeoDat->getNSectors(); iSector++) {
281 int module = m_ElementNumbers->moduleNumberEKLM(iSection, iSector, iLayer);
282 alignment.setModuleAlignment(module, &alignmentDataRandom);
283 }
284 }
285 }
286 if (alignmentChecker.checkAlignment(&alignment, &segmentAlignment))
287 alignmentStatus = 1;
288 else
289 alignmentStatus = 0;
290 t->Fill();
291 }
292 t->Write();
293}
294
296{
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;
307 EKLMAlignment alignment;
308 EKLMSegmentAlignment segmentAlignment;
309 KLMAlignmentData alignmentDataRandom;
310 EKLM::AlignmentChecker alignmentChecker(false);
311 TTree* t;
312 f->cd();
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");
320 for (jPlane = 1; jPlane <= m_GeoDat->getNPlanes(); jPlane++) {
321 printf("Plane %d\n", jPlane);
322 for (jSegment = 1; jSegment <= m_GeoDat->getNSegments(); jSegment++) {
323 printf("Segment %d\n", jSegment);
324 for (i = 0; i < nPoints; i++) {
325 fillZeroDisplacements(&alignment, &segmentAlignment);
326 deltaU = gRandom->Uniform(minDeltaU, maxDeltaU);
327 deltaV = gRandom->Uniform(minDeltaV, maxDeltaV);
328 deltaGamma = gRandom->Uniform(minDeltaGamma, maxDeltaGamma);
329 alignmentDataRandom.setDeltaU(deltaU);
330 alignmentDataRandom.setDeltaV(deltaV);
331 alignmentDataRandom.setDeltaGamma(deltaGamma);
332 for (iSection = 1; iSection <= m_GeoDat->getNSections(); iSection++) {
333 for (iLayer = 1; iLayer <= m_GeoDat->getNDetectorLayers(iSection);
334 iLayer++) {
335 for (iSector = 1; iSector <= m_GeoDat->getNSectors(); iSector++) {
337 iSection, iLayer, iSector, jPlane, jSegment);
338 segmentAlignment.setSegmentAlignment(
339 segment, &alignmentDataRandom);
340 }
341 }
342 }
343 if (alignmentChecker.checkAlignment(&alignment, &segmentAlignment))
344 alignmentStatus = 1;
345 else
346 alignmentStatus = 0;
347 t->Fill();
348 }
349 }
350 }
351 t->Write();
352}
353
355{
356 TFile* f;
357 f = new TFile(outputFile, "recreate");
360 delete f;
361}
362
364 EKLMAlignment* alignment, EKLMSegmentAlignment* segmentAlignment,
365 const char* outputFile)
366{
367 int iSection, iLayer, iSector, iPlane, iSegment, segment, param;
368 float value;
369 KLMAlignmentData* alignmentData;
370 TFile* f;
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", &param, "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", &param, "param/I");
386 tEKLMSegment->Branch("value", &value, "value/F");
387 for (iSection = 1; iSection <= m_GeoDat->getNSections(); iSection++) {
388 for (iLayer = 1; iLayer <= m_GeoDat->getNDetectorLayers(iSection);
389 iLayer++) {
390 for (iSector = 1; iSector <= m_GeoDat->getNSectors(); iSector++) {
391 int module = m_ElementNumbers->moduleNumberEKLM(iSection, iSector, iLayer);
392 alignmentData = const_cast<KLMAlignmentData*>(
393 alignment->getModuleAlignment(module));
394 param = 1;
395 value = alignmentData->getDeltaU();
396 tEKLMModule->Fill();
397 /* cppcheck-suppress redundantAssignment */
398 param = 2;
399 /* cppcheck-suppress redundantAssignment */
400 value = alignmentData->getDeltaV();
401 tEKLMModule->Fill();
402 /* cppcheck-suppress redundantAssignment */
403 param = 6;
404 /* cppcheck-suppress redundantAssignment */
405 value = alignmentData->getDeltaGamma();
406 tEKLMModule->Fill();
407 for (iPlane = 1; iPlane <= m_GeoDat->getNPlanes(); iPlane++) {
408 for (iSegment = 1; iSegment <= m_GeoDat->getNSegments(); iSegment++) {
410 iSection, iLayer, iSector, iPlane, iSegment);
411 alignmentData = const_cast<KLMAlignmentData*>(
412 segmentAlignment->getSegmentAlignment(segment));
413 param = 2;
414 value = alignmentData->getDeltaV();
415 tEKLMSegment->Fill();
416 /* cppcheck-suppress redundantAssignment */
417 param = 6;
418 /* cppcheck-suppress redundantAssignment */
419 value = alignmentData->getDeltaGamma();
420 tEKLMSegment->Fill();
421 }
422 }
423 }
424 }
425 }
426 f->cd();
427 tEKLMModule->Write();
428 tEKLMSegment->Write();
429 delete tEKLMModule;
430 delete tEKLMSegment;
431 delete f;
432}
Class to store EKLM alignment data in the database.
Definition: EKLMAlignment.h:30
EKLM element numbers.
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.
KLM Alignment data.
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.
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).
KLM element numbers.
KLMModuleNumber moduleNumberEKLM(int section, int sector, int layer) const
Get module number for EKLM.
static const double rad
Standard of [angle].
Definition: Unit.h:50
static const double cm
Standard units with the value = 1.
Definition: Unit.h:47
Abstract base class for different kinds of events.