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.

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
84 // Get the beamline channels we want to correlate
85 detectors.LoadMockDataParameters("mock_parameters_list.map");
86
87 // new vectors for GetSubsystemByType
88 std::vector <QwDetectorArray*> detchannels;
89 std::vector <VQwSubsystem*> tempvector = detectors.GetSubsystemByType("QwDetectorArray");
90 // detectors.GetSubsystemByType;
91
92 for (std::size_t i = 0; i < tempvector.size(); i++){
93 detchannels.push_back(dynamic_cast<QwDetectorArray*>(tempvector[i]));
94
95 //return detchannels;
96 }
97
98 // Initialize the stopwatch
99 TStopwatch stopwatch;
100
101 // Loop over all runs
102 UInt_t runnumber_min = (UInt_t) gQwOptions.GetIntValuePairFirst("run");
103 UInt_t runnumber_max = (UInt_t) gQwOptions.GetIntValuePairLast("run");
104 for (UInt_t run = runnumber_min;
105 run <= runnumber_max;
106 run++) {
107
108 // Set the random seed for this run
110 QwCombinedBCM<QwVQWK_Channel>::SetTripSeed(0x56781234 ^ (run*run));
111
112 // Open new output file
113 // (giving run number as argument to OpenDataFile confuses the segment search)
114
115
116 TString filename = Form("%sQwMock_%u.log", eventbuffer.GetDataDirectory().Data(), run);
117 if (eventbuffer.IsOnline()) {
118 if (eventbuffer.ReOpenStream() != CODA_OK) {
119 std::cout << "Error: could not open ET stream!" << std::endl;
120 return 0;
121 }
122 } else {
123 if (eventbuffer.OpenDataFile(filename,"W") != CODA_OK) {
124 std::cout << "Error: could not open file!" << std::endl;
125 return 0;
126 }
127 }
128 eventbuffer.ResetControlParameters();
129 eventbuffer.EncodePrestartEvent(run, 0); // prestart: runnumber, runtype
130 eventbuffer.EncodeGoEvent();
131
132
133 // Helicity initialization loop
134 helicity->SetEventPatternPhase(-1, -1, -1);
135 // 24-bit seed, should be larger than 0x1, 0x55 = 0101 0101
136 // Consecutive runs should have no trivially related seeds:
137 // e.g. with 0x2 * run, the first two files will be just 1 MPS offset...
138 unsigned int seed = 0x1234 ^ run;
139 helicity->SetFirstBits(24, seed & 0xFFFFFF);
140
141
142 // Retrieve the requested range of event numbers
143 if (kDebug) std::cout << "Starting event loop..." << std::endl;
144 Int_t eventnumber_min = gQwOptions.GetIntValuePairFirst("event");
145 Int_t eventnumber_max = gQwOptions.GetIntValuePairLast("event");
146
147 // Warn when only few events are requested, probably a problem in the input
148 if (abs(eventnumber_max - eventnumber_min) < 10)
149 QwWarning << "Only " << abs(eventnumber_max - eventnumber_min)
150 << " events will be generated." << QwLog::endl;
151
152 // Event generation loop
153 for (Int_t event = eventnumber_min; event <= eventnumber_max; event++) {
154
155 // First clear the event
156 detectors.ClearEventData();
157
158 // Set the event, pattern and phase number
159 // - event number increments for every event
160 // - pattern number increments for every multiplet
161 // - phase number gives position in multiplet
162 helicity->SetEventPatternPhase(event, event / kMultiplet, event % kMultiplet + 1);
163
164 // Run the helicity predictor
165 helicity->RunPredictor();
166 // Concise helicity printout
167 if (kDebug) {
168 // - actual helicity
169 if (helicity->GetHelicityActual() == 0) std::cout << "-";
170 else if (helicity->GetHelicityActual() == 1) std::cout << "+";
171 else std::cout << "?";
172 // - delayed helicity
173 if (helicity->GetHelicityDelayed() == 0) std::cout << "(-) ";
174 else if (helicity->GetHelicityDelayed() == 1) std::cout << "(+) ";
175 else std::cout << "(?) ";
176 if (event % kMultiplet + 1 == 4) {
177 std::cout << std::hex << helicity->GetRandomSeedActual() << std::dec << ", \t";
178 std::cout << std::hex << helicity->GetRandomSeedDelayed() << std::dec << std::endl;
179 }
180 }
181
182 // Calculate the time assuming one ms for every helicity window
183 double time = event * detectors.GetWindowPeriod();
184
185 // Fill the detectors with randomized data
186
187 int myhelicity = helicity->GetHelicityActual() ? +1 : -1;
188 //std::cout << myhelicity << std::endl;
189
190 // Randomize data for this event
191 detectors.RandomizeEventData(myhelicity, time);
192// detectors.ProcessEvent();
193// beamline-> ProcessEvent(); //Do we need to keep this line now? Check the maindetector correlation with beamline devices with and without it.
194
195 for (std::size_t i = 0; i < detchannels.size(); i++){
196 detchannels[i]->ExchangeProcessedData();
197 detchannels[i]->RandomizeMollerEvent(myhelicity);
198 }
199
200 // Write this event to file
201 Int_t status = eventbuffer.EncodeSubsystemData(detectors);
202 if (status != CODA_OK) {
203 QwError << "Error: could not write event " << event << QwLog::endl;
204 break;
205 }
206
207 // Periodically print event number
208 constexpr int nevents = kDebug ? 1000 : 10000;
209 if (event % nevents == 0) {
210 QwMessage << "Generated " << event << " events ";
211 stopwatch.Stop();
212 QwMessage << "(" << stopwatch.RealTime()*1e3/nevents << " ms per event)";
213 stopwatch.Reset();
214 stopwatch.Start();
216 }
217
218 } // end of event loop
219
220
221 eventbuffer.EncodeEndEvent();
222 eventbuffer.CloseDataFile();
223 eventbuffer.ReportRunSummary();
224
225 if (eventbuffer.IsOnline()) {
226 QwMessage << "Wrote mock data run to ET stream successfully." << QwLog::endl;
227 } else {
228 QwMessage << "Wrote mock data run " << filename << " successfully." << QwLog::endl;
229 }
230
231 } // end of run loop
232
233} // 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
static const bool kDebug
static const int kMultiplet
void DefineOptionsParity(QwOptions &options)
static void Seed(uint seedval)
Seed the internal Random Variable.
Definition MQwMockable.h:72
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(), QwEventBuffer::OpenDataFile(), QwEventBuffer::ProcessOptions(), QwSubsystemArray::ProcessOptions(), QwError, QwMessage, QwWarning, QwSubsystemArray::RandomizeEventData(), QwEventBuffer::ReOpenStream(), QwEventBuffer::ReportRunSummary(), QwEventBuffer::ResetControlParameters(), QwHelicityBase::RunPredictor(), MQwMockable::Seed(), 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().