CMS 3D CMS Logo

CTPPSModifiedOpticalFunctionsESSource.cc
Go to the documentation of this file.
1 // authors: Jan Kaspar (jan.kaspar@gmail.com)
2 
8 
10 
13 
15 
17 public:
20 
21  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
22 
23  std::shared_ptr<LHCInterpolatedOpticalFunctionsSetCollection> produce(const CTPPSInterpolatedOpticsRcd &);
24 
25 private:
27 
29 
30  double factor_;
31 
33 };
34 
35 //----------------------------------------------------------------------------------------------------
36 //----------------------------------------------------------------------------------------------------
37 
39  : scenario_(iConfig.getParameter<std::string>("scenario")),
40  factor_(iConfig.getParameter<double>("factor")),
41  rpDecId_45_N_(iConfig.getParameter<unsigned int>("rpId_45_N")),
42  rpDecId_45_F_(iConfig.getParameter<unsigned int>("rpId_45_F")),
43  rpDecId_56_N_(iConfig.getParameter<unsigned int>("rpId_56_N")),
44  rpDecId_56_F_(iConfig.getParameter<unsigned int>("rpId_56_F"))
45 {
46  auto cc = setWhatProduced(this, iConfig.getParameter<std::string>("outputOpticsLabel"));
47 
49  (edm::ESInputTag("", iConfig.getParameter<std::string>("inputOpticsLabel")));
50 }
51 
52 //----------------------------------------------------------------------------------------------------
53 
56 
57  desc.add<std::string>("inputOpticsLabel", "")->setComment("label of the input optics records");
58  desc.add<std::string>("outputOpticsLabel", "modified")->setComment("label of the output optics records");
59 
60  desc.add<std::string>("scenario", "none")->setComment("name of modification scenario");
61 
62  desc.add<double>("factor", 0.)->setComment("size of modification (number of sigmas)");
63 
64  desc.add<unsigned int>("rpId_45_N", 0)->setComment("decimal RP id for 45 near");
65  desc.add<unsigned int>("rpId_45_F", 0)->setComment("decimal RP id for 45 far");
66  desc.add<unsigned int>("rpId_56_N", 0)->setComment("decimal RP id for 56 near");
67  desc.add<unsigned int>("rpId_56_F", 0)->setComment("decimal RP id for 56 far");
68 
69  descriptions.add("ctppsModifiedOpticalFunctionsESSource", desc);
70 }
71 
72 //----------------------------------------------------------------------------------------------------
73 
74 std::shared_ptr<LHCInterpolatedOpticalFunctionsSetCollection> CTPPSModifiedOpticalFunctionsESSource::produce(
75  const CTPPSInterpolatedOpticsRcd &iRecord) {
76  // get input
78 
79  // prepare output
80  std::shared_ptr<LHCInterpolatedOpticalFunctionsSetCollection> output =
81  std::make_shared<LHCInterpolatedOpticalFunctionsSetCollection>(input);
82 
83  // premare arm/RP id map
84  struct ArmInfo {
85  unsigned int rpId_N = 0, rpId_F = 0;
86  };
87 
88  std::map<unsigned int, ArmInfo> armInfo;
89 
90  for (const auto &fsp : *output) {
91  CTPPSDetId rpId(fsp.first);
92  unsigned int rpDecId = 100 * rpId.arm() + 10 * rpId.station() + rpId.rp();
93 
94  if (rpDecId == rpDecId_45_N_)
95  armInfo[0].rpId_N = fsp.first;
96  if (rpDecId == rpDecId_45_F_)
97  armInfo[0].rpId_F = fsp.first;
98  if (rpDecId == rpDecId_56_N_)
99  armInfo[1].rpId_N = fsp.first;
100  if (rpDecId == rpDecId_56_F_)
101  armInfo[1].rpId_F = fsp.first;
102  }
103 
104  // loop over arms
105  bool applied = false;
106 
107  for (const auto &ap : armInfo) {
108  const auto &arm = ap.first;
109 
110  //printf("* arm %u\n", arm);
111 
112  auto &of_N = output->find(ap.second.rpId_N)->second;
113  auto &of_F = output->find(ap.second.rpId_F)->second;
114 
115  const double z_N = of_N.getScoringPlaneZ();
116  const double z_F = of_F.getScoringPlaneZ();
117  const double de_z = (arm == 0) ? z_N - z_F : z_F - z_N;
118 
119  //printf(" z_N = %.3f m, z_F = %.3f m\n", z_N*1E-2, z_F*1E-2);
120 
121  if (of_N.m_xi_values.size() != of_N.m_xi_values.size())
122  throw cms::Exception("CTPPSModifiedOpticalFunctionsESSource")
123  << "Different xi sampling of optical functions in near and far RP.";
124 
125  // loop over sampling points (in xi)
126  for (unsigned int i = 0; i < of_N.m_xi_values.size(); ++i) {
127  const double xi = of_N.m_xi_values[i];
128 
129  double x_d_N = of_N.m_fcn_values[LHCOpticalFunctionsSet::exd][i];
130  double x_d_F = of_F.m_fcn_values[LHCOpticalFunctionsSet::exd][i];
131 
132  double L_x_N = of_N.m_fcn_values[LHCOpticalFunctionsSet::eLx][i];
133  double L_x_F = of_F.m_fcn_values[LHCOpticalFunctionsSet::eLx][i];
134  double Lp_x = of_N.m_fcn_values[LHCOpticalFunctionsSet::eLpx][i];
135 
136  //printf(" xi = %.3f, Lp_x = %.3f, %.3f\n", xi, Lp_x, (L_x_F - L_x_N) / de_z);
137 
138  // apply modification scenario
139  if (scenario_ == "none")
140  applied = true;
141 
142  if (scenario_ == "Lx") {
143  const double a = 3180., b = 40.; // cm
144  const double de_L_x = factor_ * (a * xi + b);
145  L_x_N += de_L_x;
146  L_x_F += de_L_x;
147  applied = true;
148  }
149 
150  if (scenario_ == "Lpx") {
151  const double a = 0.42, b = 0.015; // dimensionless
152  const double de_Lp_x = factor_ * (a * xi + b) * Lp_x;
153  Lp_x += de_Lp_x;
154  L_x_N -= de_Lp_x * de_z / 2.;
155  L_x_F += de_Lp_x * de_z / 2.;
156  applied = true;
157  }
158 
159  if (scenario_ == "xd") {
160  const double d = 0.08; // dimensionless
161  x_d_N += x_d_N * d * factor_;
162  x_d_F += x_d_F * d * factor_;
163  applied = true;
164  }
165 
166  // TODO: for test only
167  if (scenario_ == "Lx-scale") {
168  L_x_N *= factor_;
169  L_x_F *= factor_;
170  applied = true;
171  }
172 
173  // store updated values
174  of_N.m_fcn_values[LHCOpticalFunctionsSet::exd][i] = x_d_N;
175  of_F.m_fcn_values[LHCOpticalFunctionsSet::exd][i] = x_d_F;
176 
177  of_N.m_fcn_values[LHCOpticalFunctionsSet::eLx][i] = L_x_N;
178  of_F.m_fcn_values[LHCOpticalFunctionsSet::eLx][i] = L_x_F;
179 
180  of_N.m_fcn_values[LHCOpticalFunctionsSet::eLpx][i] = Lp_x;
181  of_F.m_fcn_values[LHCOpticalFunctionsSet::eLpx][i] = Lp_x;
182  }
183 
184  // re-initialise splines
185  of_N.initializeSplines();
186  of_F.initializeSplines();
187  }
188 
189  // modification applied?
190  if (!applied)
191  edm::LogError("CTPPSModifiedOpticalFunctionsESSource") << "Could not apply scenario `" + scenario_ + "'.";
192 
193  // save modified output
194  return output;
195 }
196 
197 //----------------------------------------------------------------------------------------------------
198 
uint32_t station() const
Definition: CTPPSDetId.h:58
T getParameter(std::string const &) const
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:124
void setComment(std::string const &value)
std::shared_ptr< LHCInterpolatedOpticalFunctionsSetCollection > produce(const CTPPSInterpolatedOpticsRcd &)
static std::string const input
Definition: EdmProvDump.cc:48
ProductT const & get(ESGetToken< ProductT, DepRecordT > const &iToken) const
uint32_t arm() const
Definition: CTPPSDetId.h:51
ParameterDescriptionBase * add(U const &iLabel, T const &value)
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:60
double b
Definition: hdecay.h:120
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
double a
Definition: hdecay.h:121
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::ESGetToken< LHCInterpolatedOpticalFunctionsSetCollection, CTPPSInterpolatedOpticsRcd > inputOpticsToken_
uint32_t rp() const
Definition: CTPPSDetId.h:65