CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
AlignPCLThresholdsWriter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: CondFormats/PCLConfig
4 // Class: AlignPCLThresholdsWriter
5 //
11 //
12 // Original Author: Marco Musich
13 // Created: Wed, 22 Feb 2017 12:04:36 GMT
14 //
15 //
16 
17 // system include files
18 #include <memory>
19 
20 // user include files
23 
26 
31 
32 //
33 // class declaration
34 //
35 
36 namespace DOFs {
37  enum dof { X, Y, Z, thetaX, thetaY, thetaZ, extraDOF };
38 }
39 
41 public:
43  ~AlignPCLThresholdsWriter() override = default;
44 
45  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
46 
47 private:
48  void analyze(const edm::Event&, const edm::EventSetup&) override;
50 
51  // ----------member data ---------------------------
53  const unsigned int m_minNrecords;
54  const std::vector<edm::ParameterSet> m_parameters;
55 };
56 
57 //
58 // constructors and destructor
59 //
61  : m_record(iConfig.getParameter<std::string>("record")),
62  m_minNrecords(iConfig.getParameter<unsigned int>("minNRecords")),
63  m_parameters(iConfig.getParameter<std::vector<edm::ParameterSet> >("thresholds")) {}
64 
65 //
66 // member functions
67 //
68 
69 // ------------ method called for each event ------------
71  using namespace edm;
72 
73  // output object
74  AlignPCLThresholds myThresholds{};
75 
76  edm::LogInfo("AlignPCLThresholdsWriter") << "Size of AlignPCLThresholds object " << myThresholds.size() << std::endl;
77 
78  // loop on the PSet and insert the conditions
79  std::array<std::string, 6> mandatories = {{"X", "Y", "Z", "thetaX", "thetaY", "thetaZ"}};
80  std::vector<std::string> alignables;
81 
82  // fill the list of alignables
83  for (auto& thePSet : m_parameters) {
84  const std::string alignableId(thePSet.getParameter<std::string>("alignableId"));
85  // only if it is not yet in the list
86  if (std::find(alignables.begin(), alignables.end(), alignableId) == alignables.end()) {
87  alignables.push_back(alignableId);
88  }
89  }
90 
91  for (auto& alignable : alignables) {
98 
99  std::vector<std::string> presentDOF;
100 
101  // extra degrees of freedom
102  std::vector<AlignPCLThreshold::coordThresholds> extraDOFs = std::vector<AlignPCLThreshold::coordThresholds>();
103 
104  for (auto& thePSet : m_parameters) {
105  const std::string alignableId(thePSet.getParameter<std::string>("alignableId"));
106  const std::string DOF(thePSet.getParameter<std::string>("DOF"));
107 
108  const double cutoff(thePSet.getParameter<double>("cut"));
109  const double sigCut(thePSet.getParameter<double>("sigCut"));
110  const double maxMoveCut(thePSet.getParameter<double>("maxMoveCut"));
111  const double maxErrorCut(thePSet.getParameter<double>("maxErrorCut"));
112 
113  if (alignableId == alignable) {
114  presentDOF.push_back(DOF);
115  // create the objects
116 
117  switch (mapOntoEnum(DOF)) {
118  case DOFs::X:
119  my_X.setThresholds(cutoff, sigCut, maxErrorCut, maxMoveCut, DOF);
120  break;
121  case DOFs::Y:
122  my_Y.setThresholds(cutoff, sigCut, maxErrorCut, maxMoveCut, DOF);
123  break;
124  case DOFs::Z:
125  my_Z.setThresholds(cutoff, sigCut, maxErrorCut, maxMoveCut, DOF);
126  break;
127  case DOFs::thetaX:
128  my_tX.setThresholds(cutoff, sigCut, maxErrorCut, maxMoveCut, DOF);
129  break;
130  case DOFs::thetaY:
131  my_tY.setThresholds(cutoff, sigCut, maxErrorCut, maxMoveCut, DOF);
132  break;
133  case DOFs::thetaZ:
134  my_tZ.setThresholds(cutoff, sigCut, maxErrorCut, maxMoveCut, DOF);
135  break;
136  default:
137  edm::LogInfo("AlignPCLThresholdsWriter")
138  << "Appending Extra degree of freeedom: " << DOF << " " << mapOntoEnum(DOF) << std::endl;
140  ExtraDOF.setThresholds(cutoff, sigCut, maxErrorCut, maxMoveCut, DOF);
141  extraDOFs.push_back(ExtraDOF);
142  }
143 
144  AlignPCLThreshold a(my_X, my_tX, my_Y, my_tY, my_Z, my_tZ, extraDOFs);
145  myThresholds.setAlignPCLThreshold(alignableId, a);
146 
147  } // if alignable is found in the PSet
148  } // loop on the PSets
149 
150  // checks if all mandatories are present
151  edm::LogInfo("AlignPCLThresholdsWriter")
152  << "Size of AlignPCLThresholds object " << myThresholds.size() << std::endl;
153  for (auto& mandatory : mandatories) {
154  if (std::find(presentDOF.begin(), presentDOF.end(), mandatory) == presentDOF.end()) {
155  edm::LogWarning("AlignPCLThresholdsWriter")
156  << "Configuration for DOF: " << mandatory << " for alignable " << alignable << "is not present \n"
157  << "Will build object with defaults!" << std::endl;
158  }
159  }
160 
161  } // ends loop on the alignable units
162 
163  // set the minimum number of records to be used in pede
164  myThresholds.setNRecords(m_minNrecords);
165  edm::LogInfo("AlignPCLThresholdsWriter") << "Content of AlignPCLThresholds " << std::endl;
166 
167  // use buil-in method in the CondFormat
168  myThresholds.printAll();
169 
170  // Form the data here
172  if (poolDbService.isAvailable()) {
173  cond::Time_t valid_time = poolDbService->currentTime();
174  // this writes the payload to begin in current run defined in cfg
175  poolDbService->writeOneIOV(myThresholds, valid_time, m_record);
176  }
177 }
178 
180  if (coord == "X") {
181  return DOFs::X;
182  } else if (coord == "Y") {
183  return DOFs::Y;
184  } else if (coord == "Z") {
185  return DOFs::Z;
186  } else if (coord == "thetaX") {
187  return DOFs::thetaX;
188  } else if (coord == "thetaY") {
189  return DOFs::thetaY;
190  } else if (coord == "thetaZ") {
191  return DOFs::thetaZ;
192  } else {
193  return DOFs::extraDOF;
194  }
195 }
196 
197 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
200  desc.setComment("Plugin to write payloads of type AlignPCLThresholds");
201  desc.add<std::string>("record", "AlignPCLThresholdsRcd");
202  desc.add<unsigned int>("minNRecords", 25000);
203  edm::ParameterSetDescription desc_thresholds;
204 
205  desc_thresholds.add<std::string>("alignableId");
206  desc_thresholds.add<std::string>("DOF");
207  desc_thresholds.add<double>("cut");
208  desc_thresholds.add<double>("sigCut");
209  desc_thresholds.add<double>("maxMoveCut");
210  desc_thresholds.add<double>("maxErrorCut");
211 
212  std::vector<edm::ParameterSet> default_thresholds(1);
213  desc.addVPSet("thresholds", desc_thresholds, default_thresholds);
214  descriptions.addWithDefaultLabel(desc);
215 }
216 
217 //define this as a plug-in
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ParameterDescriptionBase * addVPSet(U const &iLabel, ParameterSetDescription const &validator, std::vector< ParameterSet > const &defaults)
const std::vector< edm::ParameterSet > m_parameters
DOFs::dof mapOntoEnum(std::string coord)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
AlignPCLThresholdsWriter(const edm::ParameterSet &)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
void analyze(const edm::Event &, const edm::EventSetup &) override
void setComment(std::string const &value)
int iEvent
Definition: GenABIO.cc:224
unsigned long long Time_t
Definition: Time.h:14
~AlignPCLThresholdsWriter() override=default
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Log< level::Info, false > LogInfo
double a
Definition: hdecay.h:119
Log< level::Warning, false > LogWarning
void setThresholds(float theCut, float theSigCut, float theErrorCut, float theMaxMoveCut, const std::string &theLabel)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)