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  const double factor_;
31 
33 };
34 
35 //----------------------------------------------------------------------------------------------------
36 //----------------------------------------------------------------------------------------------------
37 
39  : inputOpticsToken_(setWhatProduced(this, iConfig.getParameter<std::string>("outputOpticsLabel"))
40  .consumes(edm::ESInputTag("", iConfig.getParameter<std::string>("inputOpticsLabel")))),
41  scenario_(iConfig.getParameter<std::string>("scenario")),
42  factor_(iConfig.getParameter<double>("factor")),
43  rpDecId_45_N_(iConfig.getParameter<unsigned int>("rpId_45_N")),
44  rpDecId_45_F_(iConfig.getParameter<unsigned int>("rpId_45_F")),
45  rpDecId_56_N_(iConfig.getParameter<unsigned int>("rpId_56_N")),
46  rpDecId_56_F_(iConfig.getParameter<unsigned int>("rpId_56_F")) {}
47 
48 //----------------------------------------------------------------------------------------------------
49 
52 
53  desc.add<std::string>("inputOpticsLabel", "")->setComment("label of the input optics records");
54  desc.add<std::string>("outputOpticsLabel", "modified")->setComment("label of the output optics records");
55 
56  desc.add<std::string>("scenario", "none")->setComment("name of modification scenario");
57 
58  desc.add<double>("factor", 0.)->setComment("size of modification (number of sigmas)");
59 
60  desc.add<unsigned int>("rpId_45_N", 0)->setComment("decimal RP id for 45 near");
61  desc.add<unsigned int>("rpId_45_F", 0)->setComment("decimal RP id for 45 far");
62  desc.add<unsigned int>("rpId_56_N", 0)->setComment("decimal RP id for 56 near");
63  desc.add<unsigned int>("rpId_56_F", 0)->setComment("decimal RP id for 56 far");
64 
65  descriptions.add("ctppsModifiedOpticalFunctionsESSource", desc);
66 }
67 
68 //----------------------------------------------------------------------------------------------------
69 
70 std::shared_ptr<LHCInterpolatedOpticalFunctionsSetCollection> CTPPSModifiedOpticalFunctionsESSource::produce(
71  const CTPPSInterpolatedOpticsRcd &iRecord) {
72  // get input
74 
75  // prepare output
76  std::shared_ptr<LHCInterpolatedOpticalFunctionsSetCollection> output =
77  std::make_shared<LHCInterpolatedOpticalFunctionsSetCollection>(input);
78 
79  // premare arm/RP id map
80  struct ArmInfo {
81  unsigned int rpId_N = 0, rpId_F = 0;
82  };
83 
84  std::map<unsigned int, ArmInfo> armInfo;
85 
86  for (const auto &fsp : *output) {
87  CTPPSDetId rpId(fsp.first);
88  unsigned int rpDecId = 100 * rpId.arm() + 10 * rpId.station() + rpId.rp();
89 
90  if (rpDecId == rpDecId_45_N_)
91  armInfo[0].rpId_N = fsp.first;
92  if (rpDecId == rpDecId_45_F_)
93  armInfo[0].rpId_F = fsp.first;
94  if (rpDecId == rpDecId_56_N_)
95  armInfo[1].rpId_N = fsp.first;
96  if (rpDecId == rpDecId_56_F_)
97  armInfo[1].rpId_F = fsp.first;
98  }
99 
100  // loop over arms
101  bool applied = false;
102 
103  for (const auto &ap : armInfo) {
104  const auto &arm = ap.first;
105 
106  //printf("* arm %u\n", arm);
107 
108  auto &of_N = output->find(ap.second.rpId_N)->second;
109  auto &of_F = output->find(ap.second.rpId_F)->second;
110 
111  const double z_N = of_N.getScoringPlaneZ();
112  const double z_F = of_F.getScoringPlaneZ();
113  const double de_z = (arm == 0) ? z_N - z_F : z_F - z_N;
114 
115  //printf(" z_N = %.3f m, z_F = %.3f m\n", z_N*1E-2, z_F*1E-2);
116 
117  if (of_N.m_xi_values.size() != of_N.m_xi_values.size())
118  throw cms::Exception("CTPPSModifiedOpticalFunctionsESSource")
119  << "Different xi sampling of optical functions in near and far RP.";
120 
121  // loop over sampling points (in xi)
122  for (unsigned int i = 0; i < of_N.m_xi_values.size(); ++i) {
123  const double xi = of_N.m_xi_values[i];
124 
125  double x_d_N = of_N.m_fcn_values[LHCOpticalFunctionsSet::exd][i];
126  double x_d_F = of_F.m_fcn_values[LHCOpticalFunctionsSet::exd][i];
127 
128  double L_x_N = of_N.m_fcn_values[LHCOpticalFunctionsSet::eLx][i];
129  double L_x_F = of_F.m_fcn_values[LHCOpticalFunctionsSet::eLx][i];
130  double Lp_x = of_N.m_fcn_values[LHCOpticalFunctionsSet::eLpx][i];
131 
132  double L_y_N = of_N.m_fcn_values[LHCOpticalFunctionsSet::eLy][i];
133  double L_y_F = of_F.m_fcn_values[LHCOpticalFunctionsSet::eLy][i];
134  double Lp_y = of_N.m_fcn_values[LHCOpticalFunctionsSet::eLpy][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  if (scenario_ == "Lpy") {
174  const double a = 2.66, b = 0.015; // dimensionless
175  const double de_Lp_y = factor_ * (a * xi + b) * Lp_y;
176  Lp_y += de_Lp_y;
177  L_y_N -= de_Lp_y * de_z / 2.;
178  L_y_F += de_Lp_y * de_z / 2.;
179  applied = true;
180  }
181 
182  // store updated values
183  of_N.m_fcn_values[LHCOpticalFunctionsSet::exd][i] = x_d_N;
184  of_F.m_fcn_values[LHCOpticalFunctionsSet::exd][i] = x_d_F;
185 
186  of_N.m_fcn_values[LHCOpticalFunctionsSet::eLx][i] = L_x_N;
187  of_F.m_fcn_values[LHCOpticalFunctionsSet::eLx][i] = L_x_F;
188 
189  of_N.m_fcn_values[LHCOpticalFunctionsSet::eLpx][i] = Lp_x;
190  of_F.m_fcn_values[LHCOpticalFunctionsSet::eLpx][i] = Lp_x;
191 
192  of_N.m_fcn_values[LHCOpticalFunctionsSet::eLy][i] = L_y_N;
193  of_F.m_fcn_values[LHCOpticalFunctionsSet::eLy][i] = L_y_F;
194 
195  of_N.m_fcn_values[LHCOpticalFunctionsSet::eLpy][i] = Lp_y;
196  of_F.m_fcn_values[LHCOpticalFunctionsSet::eLpy][i] = Lp_y;
197  }
198 
199  // re-initialise splines
200  of_N.initializeSplines();
201  of_F.initializeSplines();
202  }
203 
204  // modification applied?
205  if (!applied)
206  edm::LogError("CTPPSModifiedOpticalFunctionsESSource") << "Could not apply scenario `" + scenario_ + "'.";
207 
208  // save modified output
209  return output;
210 }
211 
212 //----------------------------------------------------------------------------------------------------
213 
const edm::ESGetToken< LHCInterpolatedOpticalFunctionsSetCollection, CTPPSInterpolatedOpticsRcd > inputOpticsToken_
std::shared_ptr< LHCInterpolatedOpticalFunctionsSetCollection > produce(const CTPPSInterpolatedOpticsRcd &)
uint32_t arm() const
Definition: CTPPSDetId.h:57
Log< level::Error, false > LogError
static std::string const input
Definition: EdmProvDump.cc:50
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:61
d
Definition: ztail.py:151
uint32_t rp() const
Definition: CTPPSDetId.h:71
uint32_t station() const
Definition: CTPPSDetId.h:64
double b
Definition: hdecay.h:118
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
HLT enums.
double a
Definition: hdecay.h:119
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ProductT const & get(ESGetToken< ProductT, DepRecordT > const &iToken) const