JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
QwMockDataGenerator.cc File Reference

Parity mock data generator. More...

#include <iostream>
#include <random>
#include "QwLog.h"
#include "QwBeamLine.h"
#include "QwOptionsParity.h"
#include "QwEventBuffer.h"
#include "QwHelicity.h"
#include "QwHelicityPattern.h"
#include "QwBlindDetectorArray.h"
#include "QwSubsystemArrayParity.h"
#include "QwDetectorArray.h"
#include "TStopwatch.h"
+ Include dependency graph for QwMockDataGenerator.cc:

Go to the source code of this file.

Macros

#define NVARS   3
 

Functions

std::string stringify (int i)
 
int main (int argc, char *argv[])
 

Variables

static const int kMultiplet = 64
 
static const bool kBeamTrips = true
 
static const bool kDebug = false
 

Detailed Description

Parity mock data generator.

Definition in file QwMockDataGenerator.cc.

Macro Definition Documentation

◆ NVARS

#define NVARS   3

Definition at line 28 of file QwMockDataGenerator.cc.

Referenced by main().

Function Documentation

◆ main()

int main ( int argc,
char * argv[] )

Without anything, print usage

Definition at line 47 of file QwMockDataGenerator.cc.

48{
49 // Define the command line options
51
52 /// Without anything, print usage
53 if (argc == 1) {
54 gQwOptions.Usage();
55 exit(0);
56 }
57
58 // Fill the search paths for the parameter files; this sets a static
59 // variable within the QwParameterFile class which will be used by
60 // all instances.
61 // The "scratch" directory should be first.
63 QwParameterFile::AppendToSearchPath(getenv_safe_string("QWANALYSIS") + "/Analysis/prminput");
64 QwParameterFile::AppendToSearchPath(getenv_safe_string("QWANALYSIS") + "/Parity/prminput");
65
66 // Set the command line arguments and the configuration filename,
67 // and we define the options that can be used in them (using QwOptions).
68 gQwOptions.SetCommandLine(argc, argv);
69 gQwOptions.SetConfigFile("qwmockdataanalysis.conf");
70
71 // Event buffer
72 QwEventBuffer eventbuffer;
73 eventbuffer.ProcessOptions(gQwOptions);
74
75 // Detector array
77 detectors.ProcessOptions(gQwOptions);
78
79 // Get the helicity
80 QwHelicity* helicity = dynamic_cast<QwHelicity*>(detectors.GetSubsystemByName("Helicity Info"));
81 if (! helicity) QwWarning << "No helicity subsystem defined!" << QwLog::endl;
82
83 // Possible scenarios:
84 // - everything is random, no correlations at all, no asymmetries at all
85 // - variations in the mean of position, current, yield over the course of
86 // a run (linearly with run number, change mean/sigma as function of event
87 // number)
88 // - one parameter has a helicity-correlated asymmetry, every other parameter
89 // is random and independently distributed
90 // - two parameters have independent helicity-correlated asymmetries
91 // - two parameters have correlated helicity-correlated asymmetries
92 // - beam modulation
93
94 // Get the beamline channels we want to correlate
95 detectors.LoadMockDataParameters("mock_parameters_list.map");
96
97//-----------------------------------------------------------------------------------------------
98 // Get the main detector channels we want to correlate
99// QwDetectorArray* maindetector =
100// dynamic_cast<QwDetectorArray*>(detectors.GetSubsystemByName("Main Detector"));
101// if (! maindetector) QwWarning << "No main detector subsystem defined!" << QwLog::endl;
102
103 // new vectors for GetSubsystemByType
104 std::vector <QwDetectorArray*> detchannels;
105 std::vector <VQwSubsystem*> tempvector = detectors.GetSubsystemByType("QwDetectorArray");
106 // detectors.GetSubsystemByType;
107
108 for (std::size_t i = 0; i < tempvector.size(); i++){
109 detchannels.push_back(dynamic_cast<QwDetectorArray*>(tempvector[i]));
110
111 //return detchannels;
112 }
113
114/*
115if(1==2){
116 Double_t bar_mean = 2.0e4;
117 Double_t bar_sigma = 3.0e2;
118 Double_t bar_asym = 4.0e-4;
119 maindetector->SetRandomEventParameters(bar_mean, bar_sigma);
120 maindetector->SetRandomEventAsymmetry(bar_asym);
121 // Specific values
122 // disabled, wdc, 2010-07-23
123 maindetector->GetIntegrationPMT("MD2Neg")->SetRandomEventAsymmetry(2.0e-2);
124 maindetector->GetIntegrationPMT("MD2Pos")->SetRandomEventAsymmetry(2.0e-2);
125 maindetector->GetIntegrationPMT("MD3Neg")->SetRandomEventAsymmetry(3.0e-3);
126 maindetector->GetIntegrationPMT("MD3Pos")->SetRandomEventAsymmetry(3.0e-3);
127 maindetector->GetIntegrationPMT("MD4Neg")->SetRandomEventAsymmetry(4.0e-4);
128 maindetector->GetIntegrationPMT("MD4Pos")->SetRandomEventAsymmetry(4.0e-4);
129 maindetector->GetIntegrationPMT("MD5Neg")->SetRandomEventAsymmetry(5.0e-5);
130 maindetector->GetIntegrationPMT("MD5Pos")->SetRandomEventAsymmetry(5.0e-5);
131 maindetector->GetIntegrationPMT("MD6Neg")->SetRandomEventAsymmetry(6.0e-6);
132 maindetector->GetIntegrationPMT("MD6Pos")->SetRandomEventAsymmetry(6.0e-6);
133 maindetector->GetIntegrationPMT("MD7Neg")->SetRandomEventAsymmetry(7.0e-7);
134 maindetector->GetIntegrationPMT("MD7Pos")->SetRandomEventAsymmetry(7.0e-7);
135 maindetector->GetIntegrationPMT("MD8Neg")->SetRandomEventAsymmetry(8.0e-8);
136 maindetector->GetIntegrationPMT("MD8Pos")->SetRandomEventAsymmetry(8.0e-8);
137
138 // Set a asymmetric helicity asymmetry on one of the bars
139 maindetector->GetIntegrationPMT("MD1Neg")->SetRandomEventAsymmetry(5.0e-5);
140 maindetector->GetIntegrationPMT("MD1Pos")->SetRandomEventAsymmetry(-5.0e-5);
141
142 // Set a drift component (amplitude, phase, frequency)
143 maindetector->GetIntegrationPMT("MD3Neg")->AddRandomEventDriftParameters(3.0e6, 0, 60*Qw::Hz);
144 maindetector->GetIntegrationPMT("MD3Neg")->AddRandomEventDriftParameters(6.0e5, 0, 120*Qw::Hz);
145 maindetector->GetIntegrationPMT("MD3Neg")->AddRandomEventDriftParameters(4.5e5, 0, 240*Qw::Hz);
146
147} //end if(1==2)
148*/
149
150
151 // Initialize randomness provider and distribution
152 std::mt19937 randomnessGenerator(999); // Mersenne twister with seed (see below)
153 std::normal_distribution<double> normalDistribution;
154 auto normal = [&]() -> double { return normalDistribution(randomnessGenerator); };
155 // WARNING: This variate_generator will return the SAME random values as the
156 // variate_generator in QwVQWK_Channel when used with the same default seed!
157 // Therefore you really should give an explicitly different seed for the
158 // mt19937 randomness generator.
159
160 // Initialize the stopwatch
161 TStopwatch stopwatch;
162
163 // Loop over all runs
164 UInt_t runnumber_min = (UInt_t) gQwOptions.GetIntValuePairFirst("run");
165 UInt_t runnumber_max = (UInt_t) gQwOptions.GetIntValuePairLast("run");
166 for (UInt_t run = runnumber_min;
167 run <= runnumber_max;
168 run++) {
169
170 // Set the random seed for this run
171 randomnessGenerator.seed(run);
172 QwCombinedBCM<QwVQWK_Channel>::SetTripSeed(0x56781234 ^ (run*run));
173
174 // Open new output file
175 // (giving run number as argument to OpenDataFile confuses the segment search)
176
177
178 TString filename = Form("%sQwMock_%u.log", eventbuffer.GetDataDirectory().Data(), run);
179 if (eventbuffer.IsOnline()) {
180 if (eventbuffer.ReOpenStream() != CODA_OK) {
181 std::cout << "Error: could not open ET stream!" << std::endl;
182 return 0;
183 }
184 } else {
185 if (eventbuffer.OpenDataFile(filename,"W") != CODA_OK) {
186 std::cout << "Error: could not open file!" << std::endl;
187 return 0;
188 }
189 }
190 eventbuffer.ResetControlParameters();
191 eventbuffer.EncodePrestartEvent(run, 0); // prestart: runnumber, runtype
192 eventbuffer.EncodeGoEvent();
193
194
195 // Helicity initialization loop
196 helicity->SetEventPatternPhase(-1, -1, -1);
197 // 24-bit seed, should be larger than 0x1, 0x55 = 0101 0101
198 // Consecutive runs should have no trivially related seeds:
199 // e.g. with 0x2 * run, the first two files will be just 1 MPS offset...
200 unsigned int seed = 0x1234 ^ run;
201 helicity->SetFirstBits(24, seed & 0xFFFFFF);
202
203
204 // Retrieve the requested range of event numbers
205 if (kDebug) std::cout << "Starting event loop..." << std::endl;
206 Int_t eventnumber_min = gQwOptions.GetIntValuePairFirst("event");
207 Int_t eventnumber_max = gQwOptions.GetIntValuePairLast("event");
208
209 // Warn when only few events are requested, probably a problem in the input
210 if (abs(eventnumber_max - eventnumber_min) < 10)
211 QwWarning << "Only " << abs(eventnumber_max - eventnumber_min)
212 << " events will be generated." << QwLog::endl;
213
214 // Event generation loop
215 for (Int_t event = eventnumber_min; event <= eventnumber_max; event++) {
216
217 // First clear the event
218 detectors.ClearEventData();
219
220 // Set the event, pattern and phase number
221 // - event number increments for every event
222 // - pattern number increments for every multiplet
223 // - phase number gives position in multiplet
224 helicity->SetEventPatternPhase(event, event / kMultiplet, event % kMultiplet + 1);
225
226 // Run the helicity predictor
227 helicity->RunPredictor();
228 // Concise helicity printout
229 if (kDebug) {
230 // - actual helicity
231 if (helicity->GetHelicityActual() == 0) std::cout << "-";
232 else if (helicity->GetHelicityActual() == 1) std::cout << "+";
233 else std::cout << "?";
234 // - delayed helicity
235 if (helicity->GetHelicityDelayed() == 0) std::cout << "(-) ";
236 else if (helicity->GetHelicityDelayed() == 1) std::cout << "(+) ";
237 else std::cout << "(?) ";
238 if (event % kMultiplet + 1 == 4) {
239 std::cout << std::hex << helicity->GetRandomSeedActual() << std::dec << ", \t";
240 std::cout << std::hex << helicity->GetRandomSeedDelayed() << std::dec << std::endl;
241 }
242 }
243
244 // Calculate the time assuming one ms for every helicity window
245 double time = event * detectors.GetWindowPeriod();
246
247 // Fill the detectors with randomized data
248
249 int myhelicity = helicity->GetHelicityActual() ? +1 : -1;
250 //std::cout << myhelicity << std::endl;
251
252 // Secondly introduce correlations between variables
253 //
254 // N-dimensional correlated normal random variables:
255 // X = C' * Z
256 // with
257 // X correlated and normally distributed,
258 // Z independent and normally distributed,
259 // C the Cholesky decomposition of the positive-definite covariance matrix
260 // (C should probably be calculated offline)
261 //
262 /* Sigma =
263 1.00000 0.50000 0.50000
264 0.50000 2.00000 0.30000
265 0.50000 0.30000 1.50000
266
267 C =
268 1.00000 0.50000 0.50000
269 0.00000 1.32288 0.03780
270 0.00000 0.00000 1.11739
271
272 Sigma = C' * C
273 */
274 double z[NVARS], x[NVARS];
275 double C[NVARS][NVARS];
276 for (int var = 0; var < NVARS; var++) {
277 x[var] = 0.0;
278 z[var] = normal();
279 C[var][var] = 1.0;
280 }
281 C[0][0] = 1.0; C[0][1] = 0.5; C[0][2] = 0.5;
282 C[1][0] = 0.0; C[1][1] = 1.32288; C[1][2] = 0.03780;
283 C[2][0] = 0.0; C[2][1] = 0.0; C[2][2] = 1.11739;
284 for (int i = 0; i < NVARS; i++)
285 for (int j = 0; j < NVARS; j++)
286 x[i] += C[j][i] * z[j];
287
288 // Assign to data elements
289 //maindetector->GetChannel("MD2Neg")->SetExternalRandomVariable(x[0]);
290 //lumidetector->GetChannel("dlumi1")->SetExternalRandomVariable(x[1]);
291 //beamline->GetBCM("qwk_bcm0l07")->SetExternalRandomVariable(x[2]);
292
293
294 // Randomize data for this event
295 detectors.RandomizeEventData(myhelicity, time);
296// detectors.ProcessEvent();
297// beamline-> ProcessEvent(); //Do we need to keep this line now? Check the maindetector correlation with beamline devices with and without it.
298
299 for (std::size_t i = 0; i < detchannels.size(); i++){
300 detchannels[i]->ExchangeProcessedData();
301 detchannels[i]->RandomizeMollerEvent(myhelicity);
302 }
303
304 // Write this event to file
305 Int_t status = eventbuffer.EncodeSubsystemData(detectors);
306 if (status != CODA_OK) {
307 QwError << "Error: could not write event " << event << QwLog::endl;
308 break;
309 }
310
311 // Periodically print event number
312 constexpr int nevents = kDebug ? 1000 : 10000;
313 if (event % nevents == 0) {
314 QwMessage << "Generated " << event << " events ";
315 stopwatch.Stop();
316 QwMessage << "(" << stopwatch.RealTime()*1e3/nevents << " ms per event)";
317 stopwatch.Reset();
318 stopwatch.Start();
320 }
321
322 } // end of event loop
323
324
325 eventbuffer.EncodeEndEvent();
326 eventbuffer.CloseDataFile();
327 eventbuffer.ReportRunSummary();
328
329 if (eventbuffer.IsOnline()) {
330 QwMessage << "Wrote mock data run to ET stream successfully." << QwLog::endl;
331 } else {
332 QwMessage << "Wrote mock data run " << filename << " successfully." << QwLog::endl;
333 }
334
335 } // end of run loop
336
337} // end of main
const std::string getenv_safe_string(const char *name)
Definition QwOptions.h:43
#define gQwOptions
Definition QwOptions.h:31
#define QwError
Predefined log drain for errors.
Definition QwLog.h:39
#define QwWarning
Predefined log drain for warnings.
Definition QwLog.h:44
#define QwMessage
Predefined log drain for regular messages.
Definition QwLog.h:49
void DefineOptionsParity(QwOptions &options)
#define NVARS
static const bool kDebug
static const int kMultiplet
Event buffer management for reading and processing CODA data.
Int_t EncodeSubsystemData(QwSubsystemArray &subsystems)
Int_t EncodeGoEvent()
const TString & GetDataDirectory() const
void ReportRunSummary()
Bool_t IsOnline()
Int_t EncodeEndEvent()
Int_t EncodePrestartEvent(int runnumber, int runtype=0)
void ProcessOptions(QwOptions &options)
Sets internal flags based on the QwOptions.
void ResetControlParameters()
Int_t OpenDataFile(UInt_t current_run, Short_t seg)
static std::ostream & endl(std::ostream &)
End of the line.
Definition QwLog.cc:297
static void AppendToSearchPath(const TString &searchdir)
Add a directory to the search path.
static void SetTripSeed(uint seedval)
Subsystem for managing arrays of PMT detectors with integration and combination.
Subsystem for helicity state management and pattern recognition.
Definition QwHelicity.h:35
UInt_t GetRandomSeedDelayed()
Definition QwHelicity.h:78
UInt_t GetRandomSeedActual()
Definition QwHelicity.h:77
Int_t GetHelicityActual()
void SetEventPatternPhase(Int_t event, Int_t pattern, Int_t phase)
void SetFirstBits(UInt_t nbits, UInt_t firstbits)
Int_t GetHelicityDelayed()
Subsystem array container specialized for parity analysis with asymmetry calculations.

References QwParameterFile::AppendToSearchPath(), QwSubsystemArray::ClearEventData(), QwEventBuffer::CloseDataFile(), DefineOptionsParity(), QwEventBuffer::EncodeEndEvent(), QwEventBuffer::EncodeGoEvent(), QwEventBuffer::EncodePrestartEvent(), QwEventBuffer::EncodeSubsystemData(), QwLog::endl(), QwEventBuffer::GetDataDirectory(), getenv_safe_string(), QwHelicityBase::GetHelicityActual(), QwHelicityBase::GetHelicityDelayed(), QwHelicity::GetRandomSeedActual(), QwHelicity::GetRandomSeedDelayed(), QwSubsystemArrayParity::GetSubsystemByName(), QwSubsystemArray::GetSubsystemByType(), QwSubsystemArray::GetWindowPeriod(), gQwOptions, QwEventBuffer::IsOnline(), kDebug, kMultiplet, QwSubsystemArrayParity::LoadMockDataParameters(), NVARS, QwEventBuffer::OpenDataFile(), QwEventBuffer::ProcessOptions(), QwSubsystemArray::ProcessOptions(), QwError, QwMessage, QwWarning, QwSubsystemArray::RandomizeEventData(), QwEventBuffer::ReOpenStream(), QwEventBuffer::ReportRunSummary(), QwEventBuffer::ResetControlParameters(), QwHelicityBase::RunPredictor(), QwHelicityBase::SetEventPatternPhase(), QwHelicityBase::SetFirstBits(), and QwCombinedBCM< T >::SetTripSeed().

+ Here is the call graph for this function:

◆ stringify()

std::string stringify ( int i)
inline

Definition at line 41 of file QwMockDataGenerator.cc.

41 {
42 std::ostringstream stream;
43 stream << i;
44 return stream.str();
45}

Variable Documentation

◆ kBeamTrips

const bool kBeamTrips = true
static

Definition at line 35 of file QwMockDataGenerator.cc.

◆ kDebug

const bool kDebug = false
static

Definition at line 38 of file QwMockDataGenerator.cc.

Referenced by main().

◆ kMultiplet

const int kMultiplet = 64
static

Definition at line 32 of file QwMockDataGenerator.cc.

Referenced by main().