a Code for the Combination of Indirect and Direct Constraints on High Energy Physics Models Logo
GenerateEvent.cpp
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2013 HEPfit Collaboration
3  *
4  *
5  * For the licensing terms see doc/COPYING.
6  */
7 
8 #include "GenerateEvent.h"
9 #include <TSystem.h>
10 #include <TRandom3.h>
11 #include <iomanip>
12 
14  const std::string& ModelConf_i,
15  const std::string& OutDirName_i,
16  const std::string& JobTag_i
17  /*const bool noMC_i*/)
18 : myInputParser(ModelF, ThObsF)
19 {
20  outputTerm = 0;
21  ModelConf = ModelConf_i;
22  JobTag = JobTag_i;
23  if (OutDirName_i != ""){
24  OutDirName = "GeneratedEvents/" + OutDirName_i + JobTag;
25  OldOutDirName = ("GeneratedEvents/DELETED/" + OutDirName_i + JobTag);
26  ObsDirName = OutDirName + "/Observables";
27  ParsDirName = OutDirName + "/Parameters";
28  CGODirName = OutDirName + "/CGO";
29  SMDebugDirName = OutDirName + "/SMDebug";
30  SUSYDebugDirName = OutDirName + "/SUSYDebug";
31  outputTerm = 1;
32  }
33  //noMC = noMC_i;
34 
35 #ifdef _MPI
36  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
37  MPI_Comm_size(MPI_COMM_WORLD, &procnum);
38  MPI_Get_processor_name(processorName, &nameLen);
39 #else
40  rank = 0;
41  procnum = 1;
42 #endif
43 }
44 
46 {
47  Obs.clear();
48  boost::ptr_vector<Observable>().swap(Obs);
49  Mod = NULL;
50 }
51 
52 void GenerateEvent::generate(int unsigned nIteration_i, int seed, bool weight_i)
53 {
54  weight = weight_i;
55  nIteration = nIteration_i;
56  if (nIteration == 0) outputTerm = 0;
57  int it_comp = 0;
58  int rem_iteration = nIteration + 1;
59  int itno = 0;
60 #ifdef _MPI
61  gRandom->SetSeed(seed * (/*(unsigned)time(NULL)*/ seed + rank * procnum + nameLen));
62 #else
63  gRandom->SetSeed(seed);
64 #endif
65 
66  try {
67  /* initialize model and set model parameters */
68  initModel();
69 
70  /* define parameter distributions */
72 
73  /* create directory structure and necessary files */
74  if (rank == 0) createDirectories();
75 
76  /* memory allocation for MPI */
77  sendbuff = new double[buffersize];
78  sendbuff_w = new double[buffersize];
79  sendbuff_int = new int[1];
80  if( rank == 0 ){
81  buff = new double*[procnum];
82  buff[0]=new double[procnum * buffersize];
83  for(int i = 1; i < procnum; i++) buff[i] = buff[i - 1] + buffersize;
84 
85  buff_w = new double*[procnum];
86  buff_w[0]=new double[procnum * buffersize];
87  for(int i = 1; i < procnum; i++) buff_w[i] = buff_w[i - 1] + buffersize;
88  } else {
89  buff = new double*[1];
90  buff[0] = new double[1];
91  buff_w = new double*[1];
92  buff_w[0] = new double[1];
93  }
94  buff_int = new int*[procnum];
95  buff_int[0]=new int[procnum];
96  for(int i = 1; i < procnum; i++) buff_int[i] = buff_int[i - 1] + 1;
97 
98  /* computation of observables */
99  while (rem_iteration > 0) {
100  itno = rank + nIteration - rem_iteration + 1;
101  if (rank < rem_iteration){
102  generateRandomEvent(itno);
103  for (boost::ptr_vector<Observable>::iterator it = Obs.begin(); it < Obs.end(); it++) {
104  if (it->getObsType().compare("HiggsObservable") != 0) {
105  sendbuff[positionID] = it->computeTheoryValue();
106  if (weight && it->getDistr().compare("noweight") != 0) sendbuff_w[positionID] = it->computeWeight();
107  else sendbuff_w[positionID] = 0.;
108  positionID++;
109  } else {
110  double hw = it->computeWeight();
111  std::vector<double> theoryValues;
112  it->getTheoryValues(theoryValues);
113  for (int j = 0; j < it->getNTheoryValues() + 1 + it->getNChannels(); j++) {
114  sendbuff[positionID] = theoryValues.at(j);
115  if (weight) sendbuff_w[positionID] = hw;
116  positionID++;
117  }
118  }
119  }
120  for (std::vector<CorrelatedGaussianObservables>::iterator it1 = CGO.begin(); it1 < CGO.end(); it1++) {
121  std::vector<Observable> ObsInCGO = it1->getObs();
122  for (std::vector<Observable>::iterator it2 = ObsInCGO.begin(); it2 < ObsInCGO.end(); it2++) {
123  sendbuff[positionID] = it2->computeTheoryValue();
124  if (weight && it2->getDistr().compare("noweight") != 0) sendbuff_w[positionID] = it1->computeWeight();
125  else sendbuff_w[positionID] = 0.;
126  positionID++;
127  }
128  }
129  sendbuff_int[0] = rank + 1;
130  } else {
131  for (int i = 0; i < buffersize; i++){
132  sendbuff[i] = 0.;
133  sendbuff_w[i] = 0.;
134  }
135  sendbuff_int[0] = 0;
136  }
137 
138 #ifdef _MPI
139  MPI_Gather(sendbuff, buffersize, MPI_DOUBLE, buff[0], buffersize, MPI_DOUBLE, 0, MPI_COMM_WORLD);
140  if (weight) MPI_Gather(sendbuff_w, buffersize, MPI_DOUBLE, buff_w[0], buffersize, MPI_DOUBLE, 0, MPI_COMM_WORLD);
141  MPI_Allgather(sendbuff_int, 1, MPI_INT, buff_int[0], 1, MPI_INT, MPI_COMM_WORLD);
142 #else
143  buff[0] = sendbuff;
144  buff_w[0] = sendbuff_w;
145  buff_int[0] = sendbuff_int;
146 #endif
147 
148  /* output generation */
149  if (rank == 0) {
150  itno = nIteration - rem_iteration + 1;
151  for (int iproc = 0; iproc < procnum; iproc++) {
152  positionID = 0;
153  if (buff_int[iproc][0] != 0) {
154  for (std::vector<ModelParameter>::iterator it = ModParsVar.begin(); it < ModParsVar.end(); it++) {
155  if (outputTerm == 1) (*ParsOut[it->getname()]) << buff[iproc][positionID] << std::endl;
156  positionID++;
157  }
158  if (outputTerm == 0 && itno != 0) std::cout << "\nEvent No.: " << itno++ << std::endl;
159  else if (outputTerm == 0 && itno++ == 0) std::cout << "\nThe central event is: " << std::endl;
160  if (outputTerm == 0 && Obs.size() > 0) std::cout << "\nObservables: \n" << std::endl;
161  for (boost::ptr_vector<Observable>::iterator it = Obs.begin(); it < Obs.end(); it++) {
162  if (it->getObsType().compare("HiggsObservable") != 0) {
163  if (outputTerm == 0) {
164  if (weight && it->getDistr().compare("noweight") != 0) std::cout << it->getName() << " = " << buff[iproc][positionID] << " (weight: " << buff_w[iproc][positionID] << ")" << std::endl;
165  else std::cout << it->getName() << " = " << buff[iproc][positionID] << std::endl;
166  positionID++;
167  } else {
168  if (weight && it->getDistr().compare("noweight") != 0) (*ObsOut[it->getName()]) << buff[iproc][positionID] << "\t" << buff_w[iproc][positionID] << std::endl;
169  else (*ObsOut[it->getName()]) << buff[iproc][positionID] << std::endl;
170  positionID++;
171  }
172  } else {
173  if (outputTerm == 0) {
174  if (weight) std::cout << "\n" << it->getName() << " :" << " (weight: " << buff_w[iproc][positionID] << ")" << std::endl;
175  else std::cout << "\n" << it->getName() << " :" << std::endl;
176  for (int i = 0; i < it->getNTheoryValues(); i++) {
177  std::cout << " sigma(" << i + 1 << ") = " << buff[iproc][positionID] << std::endl;
178  positionID++;
179  }
180  std::cout << " BR = " << buff[iproc][positionID] << std::endl;
181  positionID++;
182  for (int i = it->getNTheoryValues() + 1; i < it->getNTheoryValues() + 1 + it->getNChannels(); i++) {
183  std::cout << " mu(" << i - (it->getNTheoryValues() + 1) + 1 << ") = " << buff[iproc][positionID] << std::endl;
184  positionID++;
185  }
186  } else {
187  for (int i = 0; i < it->getNTheoryValues(); i++) {
188  (*ObsOut[it->getName()]) << std::setw(10) << std::left << buff[iproc][positionID];
189  positionID++;
190  }
191  (*ObsOut[it->getName()]) << std::setw(10) << std::left << buff[iproc][positionID];
192  positionID++;
193  for (int i = it->getNTheoryValues() + 1; i < it->getNTheoryValues() + 1 + it->getNChannels(); i++) {
194  (*ObsOut[it->getName()]) << std::setw(10) << std::left << buff[iproc][positionID];
195  positionID++;
196  }
197  if (weight) (*ObsOut[it->getName()]) << std::setw(10) << std::left << buff_w[iproc][positionID - 1];
198  (*ObsOut[it->getName()]) << std::endl;
199  }
200  }
201  }
202  if (outputTerm == 0 && CGO.size() > 0) std::cout << "\nCorrelated Gaussian Observables: \n" << std::endl;
203  for (std::vector<CorrelatedGaussianObservables>::iterator it1 = CGO.begin(); it1 < CGO.end(); it1++) {
204  if (outputTerm == 0) {
205  if (weight && !(it1->isPrediction())) std::cout << it1->getName() << ": (weight: " << buff_w[iproc][positionID] << ")" << std::endl;
206  else std::cout << it1->getName() << std::endl;
207  }
208  std::vector<Observable> ObsInCGO = it1->getObs();
209  for (std::vector<Observable>::iterator it2 = ObsInCGO.begin(); it2 < ObsInCGO.end(); it2++) {
210  if (outputTerm == 0) {
211  if (weight && it2->getDistr().compare("noweight") != 0) std::cout << it2->getName() << " = " << buff[iproc][positionID] << std::endl;
212  else std::cout << it2->getName() << " = " << buff[iproc][positionID] << std::endl;
213  positionID++;
214  } else {
215  if (weight && it2->getDistr().compare("noweight") != 0) (*CGOOut[it1->getName()]) << buff[iproc][positionID] << "\t";
216  else (*CGOOut[it1->getName()]) << buff[iproc][positionID] << "\t";
217  positionID++;
218  }
219  }
220  if (outputTerm == 0) std::cout << std::endl;
221  else (*CGOOut[it1->getName()]) << buff_w[iproc][positionID - 1] << std::endl;
222  }
223  }
224  }
225  }
226 
227  /* computation of remaining iterations */
228  it_comp = 0;
229  for(int iproc = 0; iproc < procnum; iproc++){
230  it_comp = std::max(it_comp, buff_int[iproc][0]);
231  }
232  rem_iteration -= it_comp;
233  }
234 
235  /* closing of output streams */
236  if (outputTerm == 1 && rank == 0){
237  for (std::vector<ModelParameter>::iterator it = ModParsVar.begin(); it < ModParsVar.end(); it++) {
238  ParsOut[it->getname()]->close();
239  }
240  for (boost::ptr_vector<Observable>::iterator it = Obs.begin(); it < Obs.end(); it++){
241  ObsOut[it->getName()]->close();
242  }
243  for (std::vector<CorrelatedGaussianObservables>::iterator it1 = CGO.begin(); it1 < CGO.end(); it1++){
244  CGOOut[it1->getName()]->close();
245  }
246  summary << "Events\t" << nIteration << "\n";
247  summary.close();
248  std::cout << "Output written to: " << OutDirName << "\n" << std::endl;
249  } else {
250  if (rank == 0) std::cout << std::endl;
251  }
252  return;
253  } catch (std::string message) {
254  std::cerr << message << std::endl;
255  exit(EXIT_FAILURE);
256  }
257  delete buff[0];
258  delete [] buff;
259  delete [] sendbuff;
260 
261  delete buff_w[0];
262  delete [] buff_w;
263  delete [] sendbuff_w;
264 
265  delete buff_int[0];
266  delete [] buff_int;
267  delete [] sendbuff_int;
268 }
269 
271 {
272  if (nIteration > 0 && outputTerm == 1){
273  FileStat_t info;
274  if (gSystem->GetPathInfo("GeneratedEvents", info) != 0) {
275  gSystem->MakeDirectory("GeneratedEvents");
276  }
277  if (gSystem->GetPathInfo("GeneratedEvents/DELETED", info) != 0) {
278  gSystem->MakeDirectory("GeneratedEvents/DELETED");
279  }
280  if (gSystem->GetPathInfo(OutDirName.c_str(), info) == 0){
281  if (gSystem->GetPathInfo(OldOutDirName.c_str(), info) == 0){
282  gSystem->Exec(("rm -rf " + OldOutDirName).c_str());
283  std::cout << "\nWARNING: Removed " << OldOutDirName << std::endl;
284  }
285  gSystem->Exec(("mv -f " + OutDirName + " " + "GeneratedEvents/DELETED/").c_str());
286  std::cout << "WARNING: " << OutDirName << " exists!\nWARNING: " << OutDirName <<" moved to " << OldOutDirName << "\n"<< std::endl;
287  }
288  gSystem->MakeDirectory(OutDirName.c_str());
289  gSystem->MakeDirectory(ObsDirName.c_str());
290  gSystem->MakeDirectory(CGODirName.c_str());
291  gSystem->MakeDirectory(ParsDirName.c_str());
292 
293  if (outputTerm == 1){
294  summary.open((OutDirName + "/Summary.txt").c_str(), std::ios::out);
295  summary << "Model\t" << ModelName << "\n";
296  }
297 
298  for (std::vector<ModelParameter>::iterator it = ModParsVar.begin(); it < ModParsVar.end(); it++) {
299  ParsOut[it->getname()] = boost::make_shared<std::ofstream>((ParsDirName + "/" + it->getname() + ".txt").c_str(), std::ios::out);
300  summary << "Parameter\t" << it->getname() << "\n";
301  }
302 
303  for (boost::ptr_vector<Observable>::iterator it = Obs.begin(); it < Obs.end(); it++) {
304  ObsOut[it->getName()] = boost::make_shared<std::ofstream>((ObsDirName + "/" + it->getName() + ".txt").c_str(), std::ios::out);
305  if (outputTerm != 0 && it->getObsType().compare("HiggsObservable") == 0) {
306  for (int i = 0; i < it->getNTheoryValues(); i++) (*ObsOut[it->getName()]) << std::setw(6) << std::left << "sigma(" << std::setw(1) << i << std::setw(3) << std::left <<")";
307  (*ObsOut[it->getName()]) << std::setw(10) << std::left << "BR";
308  for (int i = it->getNTheoryValues() + 1; i < it->getNTheoryValues() + 1 + it->getNChannels(); i++) (*ObsOut[it->getName()]) << std::setw(3) << std::left << "mu(" << std::setw(1) << i << std::setw(6) << std::left << ")";
309  if (weight) (*ObsOut[it->getName()]) << std::setw(10) << std::left << "weight";
310  (*ObsOut[it->getName()]) << std::endl;
311  }
312  summary << "Observable\t" << it->getName() << "\n";
313  }
314  for (std::vector<CorrelatedGaussianObservables>::iterator it = CGO.begin(); it < CGO.end(); it++){
315  CGOOut[it->getName()] = boost::make_shared<std::ofstream>((CGODirName + "/" + it->getName() + ".txt").c_str(), std::ios::out);
316  summary << "CGO\t" << it->getName() << "\n";
317  }
318 
319  std::cout << "\nRunning in Event Generation mode..." /*\nWARNING: The output should not be used for any statistical analysis. \n Neither randomness or completness of the sample is gauranteed!!\n"*/ << std::endl;
320  } else if(nIteration > 0 && outputTerm == 0) {
321  std::cout.precision(16);
322  std::cout << "\nRunning in Event Generation mode... \nWARNING: Output being sent to terminal, no data written to disk!!" /* \nWARNING: The output should not be used for any statistical analysis. \n Neither randomness or completness of the sample is gauranteed!!\n "*/ << std::endl;
323  } else {
324  std::cout.precision(16);
325  std::cout << std::endl << "Running in Single Event mode...\nNo data will be written to disk!\n" << std::endl;
326  }
327 }
328 
332  if (Obs.size() == 0 && CGO.size() == 0) {
333  if (rank == 0) throw std::runtime_error("\nGenerateEvent::initModel(): No observables or correlated Gaussian observables defined in " + ModelConf + " file\n");
334  else sleep(2);
335  }
336  std::map<std::string, double> DP;
337  for (std::vector<ModelParameter>::iterator it = ModPars.begin(); it < ModPars.end(); it++) {
338  if (it->IsCorrelated()) {
339  for (unsigned int i = 0; i < CGP.size(); i++) {
340  if (CGP[i].getName().compare(it->getCgp_name()) == 0) {
341  std::string index = it->getname().substr(CGP[i].getName().size());
342  long int lindex = strtol(index.c_str(), NULL, 10);
343  if (lindex > 0) DP[CGP[i].getPar(lindex - 1).getname()] = CGP[i].getPar(lindex - 1).getave();
344  else {
345  std::stringstream out;
346  out << it->getname();
347  throw std::runtime_error("GenerateEvents::initModel(): " + out.str() + "seems to be part of a CorrelatedGaussianParameters object, but I couldn't find the corresponding object");
348  }
349  }
350  }
351  } else DP[it->getname()] = it->getave();
352  }
353  if (!myInputParser.getModel()->Init(DP)) {
354  if (myInputParser.getModel()->getmissingModelParameters().size() > 0) {
355  if (rank == 0) std::cout << "\nPlease set the following parameters in the model configuration files:\n" << std::endl;
356  for (std::vector<std::string>::iterator it = myInputParser.getModel()->getmissingModelParameters().begin(); it != myInputParser.getModel()->getmissingModelParameters().end(); it++) {
357  if (rank == 0) std::cout << "ModelParameter\t" << *it << std::endl;
358  }
359  std::cout << std::endl;
361  }
362  if (rank == 0) throw std::runtime_error("\nERROR: " + ModelName + " cannot be initialized");
363  else sleep(2);
364  }
365  std::vector<std::string> unknownParameters = Mod->getUnknownParameters();
366  for (std::vector<ModelParameter>::iterator it = ModPars.begin(); it < ModPars.end(); it++) {
367  if (it->geterrg() > 0. || it->geterrf() > 0.) {
368  if (std::find(unknownParameters.begin(), unknownParameters.end(), it->getname()) == unknownParameters.end()) ModParsVar.push_back(*it);
369  }
370  }
371 
372  if (unknownParameters.size() > 0 && rank == 0) {
373  std::cout << "\n" << std::endl;
374  for (std::vector<std::string>::iterator it = unknownParameters.begin(); it != unknownParameters.end(); it++)
375  std::cout << "WARNING: unknown parameter " << *it << " in model initialization not added the Event Generation." << std::endl;
376  }
377 
378  buffersize = ModParsVar.size();
379  for (std::vector<CorrelatedGaussianObservables>::iterator it = CGO.begin(); it < CGO.end(); it++) {
380  buffersize += it->getObs().size();
381  }
382  for (boost::ptr_vector<Observable>::iterator it = Obs.begin(); it < Obs.end(); it++) {
383  if (it->getObsType().compare("HiggsObservable") == 0) {
384  buffersize += it->getNTheoryValues() + 1 + it->getNChannels();
385  } else {
386  buffersize += 1;
387  }
388  }
389 }
390 
392  if (outputTerm == 0 && rank == 0) std::cout << "\nParameters varied in Event Generation:" << std::endl;
393  for (std::vector<ModelParameter>::const_iterator it = ModParsVar.begin(); it < ModParsVar.end(); it++) {
394  if (it->IsCorrelated()) {
395  for (unsigned int i = 0; i < CGP.size(); i++) {
396  if (CGP[i].getName().compare(it->getCgp_name()) == 0) {
397  std::string index = it->getname().substr(CGP[i].getName().size());
398  long int lindex = strtol(index.c_str(), NULL, 10);
399  if (lindex > 0) {
400  if (outputTerm == 0 && rank == 0) std::cout << CGP[i].getPar(lindex - 1).getname() << ", ";
401  } else {
402  std::stringstream out;
403  out << it->getname();
404  throw std::runtime_error("MonteCarlo::Run(): " + out.str() + "seems to be part of a CorrelatedGaussianParameters object, but I couldn't find the corresponding object");
405  }
406  }
407  }
408  } else if (outputTerm == 0 && rank == 0) std::cout << it->getname() << ", "; //<< k << std::endl;
409  if (it->geterrf() == 0. && it->geterrg() != 0.) {
410  DDist[it->getname()] = new TF1(it->getname().c_str(),
411  "1./sqrt(2.*TMath::Pi())/[1] * exp(-(x-[0])*(x-[0])/2./[1]/[1])",
412  it->getmin(), it->getmax());
413  DDist[it->getname()]->SetParameter(0, it->getave());
414  DDist[it->getname()]->SetParameter(1, it->geterrg());
415  } else if (it->geterrg() == 0. && it->geterrf() != 0.) {
416  DDist[it->getname()] = new TF1(it->getname().c_str(),
417  "1/([1]-[0])",
418  it->getmin(), it->getmax());
419  DDist[it->getname()]->SetParameter(0, it->getmin());
420  DDist[it->getname()]->SetParameter(1, it->getmax());
421  } else {
422  DDist[it->getname()] = new TF1(it->getname().c_str(),
423  "(TMath::Erf((x-[0]+[2])/sqrt(2.)/[1])-TMath::Erf((x-[0]-[2])/sqrt(2.)/[1]))/4./[2]",
424  it->getmin(), it->getmax());
425  DDist[it->getname()]->SetParameter(0, it->getave());
426  DDist[it->getname()]->SetParameter(1, it->geterrg());
427  DDist[it->getname()]->SetParameter(2, it->geterrf());
428  }
429  }
430  if (outputTerm == 0 && rank == 0) std::cout << std::endl;
431 }
432 
433 void GenerateEvent::generateRandomEvent(int iterationNo) {
434  positionID = 0;
435 
436  std::vector<double> vec(ModParsVar.size(), 0.);
437  for (unsigned int i = 0; i < vec.size(); i++) {
438  if (iterationNo == 0) vec[i] = ModParsVar[i].getave();
439  else vec[i] = DDist[ModParsVar[i].getname()]->GetRandom();
440  }
441 
443  for (std::vector<ModelParameter>::const_iterator it = ModParsVar.begin(); it < ModParsVar.end(); it++) {
444  if (it->IsCorrelated()) {
445  for (unsigned int i = 0; i < CGP.size(); i++) {
446  if (CGP[i].getName().compare(it->getCgp_name()) == 0) {
447  std::string index = it->getname().substr(CGP[i].getName().size());
448  long int lindex = strtol(index.c_str(), NULL, 10);
449  if (lindex > 0) {
450  sendbuff[positionID] = DPars[CGP[i].getPar(lindex - 1).getname()];
451  sendbuff_w[positionID] = DPars[CGP[i].getPar(lindex - 1).getname()];
452  positionID++;
453  } else {
454  std::stringstream out;
455  out << it->getname();
456  throw std::runtime_error("MonteCarlo::Run(): " + out.str() + "seems to be part of a CorrelatedGaussianParameters object, but I couldn't find the corresponding object");
457  }
458  }
459  }
460  } else {
461  sendbuff[positionID] = DPars[it->getname()];
462  sendbuff_w[positionID] = DPars[it->getname()];
463  positionID++;
464  }
465  }
466  Mod->Update(DPars);
467 }
468 
469 void GenerateEvent::setDParsFromParameters(const std::vector<double>& parameters, std::map<std::string, double>& DPars_i) {
470  std::map<std::string, std::vector<double> > cgpmap;
471 
472  unsigned int k = 0;
473  for (std::vector<ModelParameter>::const_iterator it = ModParsVar.begin(); it < ModParsVar.end(); it++) {
474 
475  if (it->IsCorrelated()) {
476  std::string index = it->getname().substr(it->getCgp_name().size());
477  unsigned long int lindex = strtol(index.c_str(), NULL, 10);
478  if (lindex - 1 == cgpmap[it->getCgp_name()].size()) cgpmap[it->getCgp_name()].push_back(parameters[k]);
479  else {
480  std::stringstream out;
481  out << it->getname() << " " << lindex;
482  throw std::runtime_error("GenerateEvent::setDParsFromParameters(): " + out.str() + "seems to be a CorrelatedGaussianParameters object but the corresponding parameters are missing or not in the right order");
483  }
484 
485  } else DPars_i[it->getname()] = parameters[k];
486  k++;
487  }
488 
489  for (unsigned int j = 0; j < CGP.size(); j++) {
490  std::vector<double> current = cgpmap.at(CGP[j].getName());
491  if (current.size() != CGP[j].getPars().size()) {
492  std::stringstream out;
493  out << CGP[j].getName();
494  throw std::runtime_error("GenerateEvent::setDParsFromParameters(): " + out.str() + " appears to be represented in cgpmap with a wrong size");
495  }
496 
497  std::vector<double> porig = CGP[j].getOrigParsValue(current);
498 
499  for (unsigned int l = 0; l < porig.size(); l++) {
500  DPars_i[CGP[j].getPar(l).getname()] = porig[l];
501  }
502  }
503 }
504 
505 void GenerateEvent::addCustomObservableType(const std::string name, boost::function<Observable*() > funct){
507 }
ModelFactory
A class for.
Definition: ModelFactory.h:25
GenerateEvent::GenerateEvent
GenerateEvent(ModelFactory &ModelF, ThObsFactory &ThObsF, const std::string &ModelConf_i, const std::string &OutDirName_i, const std::string &JobTag_i)
Constructor.
Definition: GenerateEvent.cpp:13
GenerateEvent::CGODirName
std::string CGODirName
String for the name of the Correlated Gaussian Observables output directory.
Definition: GenerateEvent.h:131
GenerateEvent::buff_int
int ** buff_int
Definition: GenerateEvent.h:153
GenerateEvent::Obs2D
std::vector< Observable2D > Obs2D
Vector for the Observables2D defined in SomeModel.conf.
Definition: GenerateEvent.h:124
GenerateEvent::CGO
std::vector< CorrelatedGaussianObservables > CGO
vector for the Correlated Gaussian Observables defined in SomeModel.conf.
Definition: GenerateEvent.h:125
GenerateEvent::initModel
void initModel()
Definition: GenerateEvent.cpp:329
GenerateEvent::ObsDirName
std::string ObsDirName
String for the name of the observables output directory.
Definition: GenerateEvent.h:130
GenerateEvent::DPars
std::map< std::string, double > DPars
Map of parameters to be passed to Model().
Definition: GenerateEvent.h:115
GenerateEvent::weight
bool weight
Definition: GenerateEvent.h:155
GenerateEvent::outputTerm
bool outputTerm
Flag to specify output stream storage.
Definition: GenerateEvent.h:137
GenerateEvent::nIteration
int nIteration
Definition: GenerateEvent.h:145
GenerateEvent::sendbuff
double * sendbuff
Definition: GenerateEvent.h:148
Model::getmissingModelParameters
std::vector< std::string > getmissingModelParameters()
Definition: Model.h:237
GenerateEvent::addCustomObservableType
void addCustomObservableType(const std::string name, boost::function< Observable *() > funct)
Definition: GenerateEvent.cpp:505
GenerateEvent::Mod
StandardModel * Mod
Name of the model as defined in SomeModel.conf.
Definition: GenerateEvent.h:120
GenerateEvent::buffersize
int buffersize
Definition: GenerateEvent.h:147
GenerateEvent::OldOutDirName
std::string OldOutDirName
String for the name of the backup output directory.
Definition: GenerateEvent.h:129
GenerateEvent::processorName
char processorName[MPI_MAX_PROCESSOR_NAME]
Definition: GenerateEvent.h:142
GenerateEvent::ModelConf
std::string ModelConf
String for the name of the SomeModel.conf file.
Definition: GenerateEvent.h:127
GenerateEvent::ModParsVar
std::vector< ModelParameter > ModParsVar
Vector for the model parameters varied in SomeModel.conf.
Definition: GenerateEvent.h:122
GenerateEvent::rank
int rank
Definition: GenerateEvent.h:138
GenerateEvent::JobTag
std::string JobTag
String for the optional JobTag argument to be passes to the executable.
Definition: GenerateEvent.h:135
GenerateEvent::CGP
std::vector< CorrelatedGaussianParameters > CGP
vector for the Correlated Gaussian Parameters defined in SomeModel.conf.
Definition: GenerateEvent.h:126
GenerateEvent::ParsOut
std::map< std::string, boost::shared_ptr< std::ofstream > > ParsOut
Map of output stream for parameters.
Definition: GenerateEvent.h:119
InputParser::ReadParameters
std::string ReadParameters(const std::string filename_i, const int rank, std::vector< ModelParameter > &ModelPars, boost::ptr_vector< Observable > &Observables, std::vector< Observable2D > &Observables2D, std::vector< CorrelatedGaussianObservables > &CGO, std::vector< CorrelatedGaussianParameters > &CGP)
The member that parses the Observable2D directives from SomeModel.conf file.
Definition: InputParser.cpp:35
ThObsFactory
A class for.
Definition: ThObsFactory.h:26
GenerateEvent::Obs
boost::ptr_vector< Observable > Obs
Vector for the observables defined in SomeModel.conf.
Definition: GenerateEvent.h:123
GenerateEvent::defineParameterDistributions
void defineParameterDistributions()
The parameter distributions as specified in SomeModel.conf is defined here.
Definition: GenerateEvent.cpp:391
StandardModel::Init
virtual bool Init(const std::map< std::string, double > &DPars)
A method to initialize the model parameters.
Definition: StandardModel.cpp:159
GenerateEvent::buff
double ** buff
Definition: GenerateEvent.h:150
InputParser::getModel
StandardModel * getModel() const
A get method to access the pointer to the object of the StandardModel class.
Definition: InputParser.h:109
GenerateEvent::SUSYDebugDirName
std::string SUSYDebugDirName
String for the name of the Correlated Gaussian Observables output directory.
Definition: GenerateEvent.h:133
GenerateEvent::generate
void generate(int unsigned nIteration_i, int seed=0, bool weight_i=false)
The method used to generate events and format output.
Definition: GenerateEvent.cpp:52
Observable
A class for observables.
Definition: Observable.h:28
GenerateEvent::positionID
int positionID
Definition: GenerateEvent.h:146
GenerateEvent::myInputParser
InputParser myInputParser
An oject of the InputParser() class.
Definition: GenerateEvent.h:114
GenerateEvent::SMDebugDirName
std::string SMDebugDirName
String for the name of the observables output directory.
Definition: GenerateEvent.h:132
GenerateEvent::sendbuff_w
double * sendbuff_w
Definition: GenerateEvent.h:149
GenerateEvent::CGOOut
std::map< std::string, boost::shared_ptr< std::ofstream > > CGOOut
Map of output stream for corellated Gaussian observables.
Definition: GenerateEvent.h:118
GenerateEvent::OutDirName
std::string OutDirName
String for the name of the output directory.
Definition: GenerateEvent.h:128
GenerateEvent::ObsOut
std::map< std::string, boost::shared_ptr< std::ofstream > > ObsOut
Map of output stream for observables.
Definition: GenerateEvent.h:117
InputParser::addCustomObservableType
void addCustomObservableType(const std::string name, boost::function< Observable *() > funct)
Definition: InputParser.cpp:250
GenerateEvent.h
StandardModel::Update
virtual bool Update(const std::map< std::string, double > &DPars)
The update method for StandardModel.
Definition: StandardModel.cpp:183
GenerateEvent::nameLen
int nameLen
Definition: GenerateEvent.h:141
GenerateEvent::procnum
int procnum
Definition: GenerateEvent.h:139
GenerateEvent::sendbuff_int
int * sendbuff_int
Definition: GenerateEvent.h:152
GenerateEvent::buff_w
double ** buff_w
Definition: GenerateEvent.h:151
GenerateEvent::createDirectories
void createDirectories()
Definition: GenerateEvent.cpp:270
GenerateEvent::DDist
std::map< std::string, TF1 * > DDist
Map of parameter distributions.
Definition: GenerateEvent.h:116
GenerateEvent::generateRandomEvent
void generateRandomEvent(int iterationNo)
This member generates random numbers for the parameters being varied in the model....
Definition: GenerateEvent.cpp:433
GenerateEvent::ModPars
std::vector< ModelParameter > ModPars
Vector for the model parameters defined in SomeModel.conf.
Definition: GenerateEvent.h:121
GenerateEvent::~GenerateEvent
virtual ~GenerateEvent()
The default destructor.
Definition: GenerateEvent.cpp:45
GenerateEvent::ParsDirName
std::string ParsDirName
String for the name of the parameters output directory.
Definition: GenerateEvent.h:134
GenerateEvent::ModelName
std::string ModelName
Definition: GenerateEvent.h:154
QCD::getUnknownParameters
std::vector< std::string > getUnknownParameters()
A method to get the vector of the parameters that have been specified in the configuration file but n...
Definition: QCD.h:467
GenerateEvent::summary
std::ofstream summary
Definition: GenerateEvent.h:144
GenerateEvent::setDParsFromParameters
void setDParsFromParameters(const std::vector< double > &parameters, std::map< std::string, double > &DPars_i)
Definition: GenerateEvent.cpp:469