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

Parity mock data generator. More...

#include <iostream>
#include <boost/random.hpp>
#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 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 29 of file QwMockDataGenerator.cc.

Referenced by main().

Function Documentation

◆ main()

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

Without anything, print usage

Definition at line 48 of file QwMockDataGenerator.cc.

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

+ Here is the call graph for this function:

◆ stringify()

std::string stringify ( int i)
inline

Definition at line 42 of file QwMockDataGenerator.cc.

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

Variable Documentation

◆ kBeamTrips

const bool kBeamTrips = true
static

Definition at line 36 of file QwMockDataGenerator.cc.

◆ kDebug

const bool kDebug = false
static

Definition at line 39 of file QwMockDataGenerator.cc.

Referenced by main().

◆ kMultiplet

const int kMultiplet = 64
static

Definition at line 33 of file QwMockDataGenerator.cc.

Referenced by main().