source: trunk/Mars/fact/analysis/produce_drs_time_fits_file.C@ 19344

Last change on this file since 19344 was 17991, checked in by dneise, 10 years ago
write NCHIPS and NCELLS into the header. was forgotten
File size: 6.4 KB
Line 
1// The FACT collaboration
2// Dominik Neise and Sebastian Mueller October 2014
3// A callisto like CERN-Root script to process the FACT DRS
4// (Domino-Ring-Sampler) time calibration files.
5//
6// The basic script is taken from the October 2014 analysis framework running on
7// La Palma, maintained by Daniela Dorner and Thomas Bretz. There it was called
8// "callisto_drs_time.C"
9
10// Dominik and Sebastian:
11// Our first approach was to use the MWriteFitsFile class but we have not been
12// able to use it.
13
14// MWriteFitsFile fitsfile(outfile,
15// MWriteFitsFile::kSingleFile,
16// "RECREATE",
17// "fitsfile",
18// "atitle");
19// fitsfile.SetBytesPerSample("Data", 2);
20// fitsfile.AddContainer("MDrsCalibrationTime","DANDELION");
21
22#include "MLogManip.h"
23
24int produce_drs_time_fits_file(
25 const char *drs_time_file, //* path to raw drs-time-calib-file
26 const char *drs_file, //* path to drs.fits(.gz) file
27 const char *outfile //* output path like yyyymmdd_rrr.drs.time.fits
28){
29 gLog.Separator("produce_drs_time_fits_file");
30 gLog << all;
31 gLog << "DRS Timing " << drs_time_file << '\n';
32 gLog << "DRS 1024 " << drs_file << '\n';
33 gLog << endl;
34
35 MParList plist0;
36
37 MTaskList tlist0;
38 plist0.AddToList(&tlist0);
39
40 MDrsCalibration drscalib1024;
41 if (!drscalib1024.ReadFits(drs_file)){
42 gLog << "Error while opening drs amplitude calibration file:"
43 << drs_file << "Aborting" << endl;
44 return 52; // DN: Why 52 I don't know.
45 }
46 plist0.AddToList(&drscalib1024);
47
48 MDrsCalibrationTime timecam;
49 plist0.AddToList(&timecam);
50
51 MEvtLoop loop0("DetermineTimeCal");
52 loop0.SetParList(&plist0);
53
54 factfits drstimeRawFile(drs_time_file);
55 const int NumberOfChips = drstimeRawFile.GetInt("NPIX")/9; // --> 160
56 const int NumberOfCells = drstimeRawFile.GetInt("NROI"); // --> 1024
57 drstimeRawFile.close();
58
59 MRawFitsRead read0(drs_time_file);
60 tlist0.AddToList(&read0);
61 MGeomApply apply0;
62 tlist0.AddToList(&apply0);
63 MDrsCalibApply drsapply0;
64 tlist0.AddToList(&drsapply0);
65 MContinue cont0("MRawEvtHeader.GetTriggerID!=33792", "SelectTim");
66 tlist0.AddToList(&cont0);
67 MFillH fill0("MHDrsCalibrationTime");
68 fill0.SetNameTab("DeltaT");
69 tlist0.AddToList(&fill0);
70
71 if (!loop0.Eventloop(0)){
72 gLog << "Error performing the loop over the drs-time-run:"
73 << drs_time_file << "Aborting" << endl;
74 return 8; // DN: why 8 I don't know.
75 }
76
77 ofits drstimeFile(outfile);
78 drstimeFile.SetDefaultKeys();
79 drstimeFile.AddColumnDouble(
80 NumberOfChips*NumberOfCells,
81 "SamplingTimeDeviation", "nominal slices",
82 "see following COMMENT."
83 );
84 drstimeFile.SetInt("NCHIPS", NumberOfChips, "number of drs chips");
85 drstimeFile.SetInt("NCELLS", NumberOfCells, "number of cells of each chip");
86 drstimeFile.AddComment(
87"The SamplingTimeDeviation specifies the deviation of the actual ");
88 drstimeFile.AddComment(
89"sampling time of a certain cell, from its nominal sampling time.");
90 drstimeFile.AddComment(
91"It is measured in nominal slices.");
92 drstimeFile.AddComment(
93"In order to convert the sample index into the actual sampling time");
94 drstimeFile.AddComment(
95"just do: ");
96 drstimeFile.AddComment(
97" cell_index + td(cell_index) - td(stop_cell_index)");
98 drstimeFile.AddComment(
99"where td(i) refers to the SamplingTimeDeviation of cell i of the ");
100 drstimeFile.AddComment(
101"current DRS4 chip.");
102 drstimeFile.AddComment(
103"This gives you the actual sampling time of the cell measured in ");
104 drstimeFile.AddComment(
105"nominal slices [nsl] with respect to the stop cells sampling time.");
106 drstimeFile.AddComment(
107"Convert it to ns, simply by multiplication ");
108 drstimeFile.AddComment(
109"with 0.5 ns/nsl.");
110 drstimeFile.AddComment(
111"");
112
113 drstimeFile.AddComment(
114"In FACT it became common practice to measure un-time-calibrated ");
115 drstimeFile.AddComment(
116" DRS4 data in 'time slices' or simply 'samples', while ");
117 drstimeFile.AddComment(
118" time-calibrated DRS4 data usually measured in ns. Since the method,");
119 drstimeFile.AddComment(
120" that is emplyed to measure the DRS4 time calibration constants has ");
121 drstimeFile.AddComment(
122" no possibility to assert the actual length of a slice is really half");
123 drstimeFile.AddComment(
124" a nanosecond, this practice is not advisable.");
125 drstimeFile.AddComment(
126" I propose instead to call un-time-calibrated data 'samples',");
127 drstimeFile.AddComment(
128" since this is what they are. If one wants to stress the fact,");
129 drstimeFile.AddComment(
130" that no drs time calibration has been applied one should refer to ");
131 drstimeFile.AddComment(
132"'uncalibrated slices'. Since it *should* be common practice to apply");
133 drstimeFile.AddComment(
134" the drs4 time calibration in all cases, I also propose to use the");
135 drstimeFile.AddComment(
136" short term of 'slices' or 'sl'. If one wants to stress, that the");
137 drstimeFile.AddComment(
138" drs4 time calibration has actually been applied to the data,");
139 drstimeFile.AddComment(
140" the term 'calibrated slices' or 'nominal slices' or short 'nsl'.");
141
142 drstimeFile.WriteTableHeader("DRS_CELL_TIMES");
143
144 // By drs_sampling_time_deviations we refer to the deviations (measured in nominal slices) of
145 // the actual sampling time compared to the nominal sampling time of every cell of every chip.
146 double *drs_sampling_time_deviations = new double[ NumberOfChips * NumberOfCells ];
147 for (int chip = 0; chip < NumberOfChips; chip++){
148 for (int cell = 0; cell < NumberOfCells; cell++){
149 // Dominik and Sebastian:
150 // We ended with using DrsCalibrateTime.Sum() in order to retrieve
151 // the contents DrsCalibrateTime.fStat.
152 // First we wanted to access the member fStat of class DrsCalibrateTime
153 // which is declared public but it did not work out of the box.
154 // drs_sampling_time_deviations[chip*NumberOfCells+cell] =
155 // timecam.fStat[chip*NumberOfCells+cell].first;
156 drs_sampling_time_deviations[ chip * NumberOfCells + cell ] =
157 timecam.Sum( chip * NumberOfCells + cell );
158 }
159 }
160
161 drstimeFile.WriteRow(drs_sampling_time_deviations, sizeof(double) * NumberOfChips * NumberOfCells);
162 drstimeFile.close();
163
164 return 0;
165}
Note: See TracBrowser for help on using the repository browser.