CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
BeamSpotWrite2DB.cc
Go to the documentation of this file.
1 
8 // C++ standard
9 #include <string>
10 #include <fstream>
11 
12 // CMS
26 
27 // ROOT
28 #include "TFile.h"
29 #include "TTree.h"
30 
31 class BeamSpotWrite2DB : public edm::one::EDAnalyzer<edm::one::SharedResources> {
32 public:
33  explicit BeamSpotWrite2DB(const edm::ParameterSet&);
34  ~BeamSpotWrite2DB() override;
36 
37 private:
38  void analyze(const edm::Event&, const edm::EventSetup&) override;
39  void endJob() override;
40 
41  std::ifstream fasciiFile;
43 };
44 
46  usesResource("PoolDBOutputService");
47  fasciiFileName = iConfig.getUntrackedParameter<std::string>("OutputFileName");
48  if (!fasciiFileName.empty()) {
49  fasciiFile.open(fasciiFileName.c_str());
50  } else {
51  throw cms::Exception("Inconsistent Data") << " expected input file name is null\n";
52  }
53 }
54 
56 
58 
60  edm::LogPrint("BeamSpotWrite2DB") << " Read beam spot data from text file: " << fasciiFileName;
61  edm::LogPrint("BeamSpotWrite2DB") << " please see plugins/BeamSpotWrite2DB.cc for format of text file.";
62 
63  edm::LogInfo("BeamSpotWrite2DB")
64  << " Content of the file is expected to have this format with the first column as a keyword:";
65  edm::LogInfo("BeamSpotWrite2DB") << " x\n y\n z\n sigmaZ\n dxdz\n dydz\n beamWidthX\n beamWidthY";
66  for (int i = 0; i < 7; i++) {
67  for (int j = 0; j < 7; j++) {
68  edm::LogInfo("BeamSpotWrite2DB") << " cov[" << i << "][" << j << "] cov[" << i << "][" << j << "] cov[" << i
69  << "][" << j << "] cov[" << i << "][" << j << "] cov[" << i << "][" << j
70  << "] cov[" << j << "][" << j << "] cov[" << i << "][" << j << "]";
71  }
72  }
73 
74  // extract from file
75  double x, y, z, sigmaZ, dxdz, dydz, beamWidthX, beamWidthY, emittanceX, emittanceY, betastar;
77  double cov[7][7];
78  int type;
79 
80  fasciiFile >> tag >> type;
81  fasciiFile >> tag >> x;
82  fasciiFile >> tag >> y;
83  fasciiFile >> tag >> z;
84  fasciiFile >> tag >> sigmaZ;
85  fasciiFile >> tag >> dxdz;
86  fasciiFile >> tag >> dydz;
87  fasciiFile >> tag >> beamWidthX;
88  fasciiFile >> tag >> beamWidthY;
89  fasciiFile >> tag >> cov[0][0] >> cov[0][1] >> cov[0][2] >> cov[0][3] >> cov[0][4] >> cov[0][5] >> cov[0][6];
90  fasciiFile >> tag >> cov[1][0] >> cov[1][1] >> cov[1][2] >> cov[1][3] >> cov[1][4] >> cov[1][5] >> cov[1][6];
91  fasciiFile >> tag >> cov[2][0] >> cov[2][1] >> cov[2][2] >> cov[2][3] >> cov[2][4] >> cov[2][5] >> cov[2][6];
92  fasciiFile >> tag >> cov[3][0] >> cov[3][1] >> cov[3][2] >> cov[3][3] >> cov[3][4] >> cov[3][5] >> cov[3][6];
93  fasciiFile >> tag >> cov[4][0] >> cov[4][1] >> cov[4][2] >> cov[4][3] >> cov[4][4] >> cov[4][5] >> cov[4][6];
94  fasciiFile >> tag >> cov[5][0] >> cov[5][1] >> cov[5][2] >> cov[5][3] >> cov[5][4] >> cov[5][5] >> cov[5][6];
95  fasciiFile >> tag >> cov[6][0] >> cov[6][1] >> cov[6][2] >> cov[6][3] >> cov[6][4] >> cov[6][5] >> cov[6][6];
96  fasciiFile >> tag >> emittanceX;
97  fasciiFile >> tag >> emittanceY;
98  fasciiFile >> tag >> betastar;
99 
100  BeamSpotObjects abeam;
101 
102  abeam.setType(type);
103  abeam.setPosition(x, y, z);
104  abeam.setSigmaZ(sigmaZ);
105  abeam.setdxdz(dxdz);
106  abeam.setdydz(dydz);
107  abeam.setBeamWidthX(beamWidthX);
108  abeam.setBeamWidthY(beamWidthY);
109  abeam.setEmittanceX(emittanceX);
110  abeam.setEmittanceY(emittanceY);
111  abeam.setBetaStar(betastar);
112 
113  for (int i = 0; i < 7; ++i) {
114  for (int j = 0; j < 7; ++j) {
115  abeam.setCovariance(i, j, cov[i][j]);
116  }
117  }
118 
119  edm::LogPrint("BeamSpotWrite2DB") << " write results to DB...";
120 
122  if (poolDbService.isAvailable()) {
123  edm::LogPrint("BeamSpotWrite2DB") << "poolDBService available";
124  if (poolDbService->isNewTagRequest("BeamSpotObjectsRcd")) {
125  edm::LogPrint("BeamSpotWrite2DB") << "new tag requested";
126  poolDbService->createOneIOV<BeamSpotObjects>(abeam, poolDbService->beginOfTime(), "BeamSpotObjectsRcd");
127  } else {
128  edm::LogPrint("BeamSpotWrite2DB") << "no new tag requested";
129  poolDbService->appendOneIOV<BeamSpotObjects>(abeam, poolDbService->currentTime(), "BeamSpotObjectsRcd");
130  }
131  }
132  edm::LogPrint("BeamSpotWrite2DB") << "[BeamSpotWrite2DB] endJob done \n";
133 }
134 
137  desc.setComment(
138  "Writes out a DB file containing a BeamSpotObjects payload, according to parameters defined in ASCII file");
139  desc.addUntracked<std::string>("OutputFileName", {});
140  descriptions.addWithDefaultLabel(desc);
141 }
142 
143 //define this as a plug-in
T getUntrackedParameter(std::string const &, T const &) const
float dydz
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
void setCovariance(int i, int j, double val)
set i,j element of the full covariance matrix 7x7
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
void setEmittanceY(double val)
set emittance
float dxdz
std::ifstream fasciiFile
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::string fasciiFileName
void setBetaStar(double val)
set beta star
void endJob() override
~BeamSpotWrite2DB() override
void setType(int type)
set beam type
void createOneIOV(const T &payload, cond::Time_t firstSinceTime, const std::string &recordName)
void appendOneIOV(const T &payload, cond::Time_t sinceTime, const std::string &recordName)
void setComment(std::string const &value)
void setEmittanceX(double val)
set emittance
int iEvent
Definition: GenABIO.cc:224
bool isNewTagRequest(const std::string &recordName)
void setdydz(double val)
set dydz slope, crossing angle in XZ
bool isAvailable() const
Definition: Service.h:40
Log< level::Warning, true > LogPrint
void setdxdz(double val)
set dxdz slope, crossing angle
Log< level::Info, false > LogInfo
void setSigmaZ(double val)
set sigma Z, RMS bunch length
static void fillDescriptions(edm::ConfigurationDescriptions &)
void analyze(const edm::Event &, const edm::EventSetup &) override
void setPosition(double x, double y, double z)
set XYZ position
void setBeamWidthX(double val)
set average transverse beam width X
BeamSpotWrite2DB(const edm::ParameterSet &)
void setBeamWidthY(double val)
set average transverse beam width Y