Dune35tSimChannelExtractService_service.cc
Go to the documentation of this file.
1 // Dune35tSimChannelExtractService.cxx
2 
5 #include <string>
8 
9 using std::string;
10 using std::endl;
11 
12 //**********************************************************************
13 
14 namespace {
15 
16 std::ostream& operator<<(std::ostream& out, const std::vector<float>& vals) {
17  out << "[";
18  unsigned int nval = 0;
19  for ( float val : vals ) {
20  out << val;
21  if ( ++nval < vals.size() ) out << ", ";
22  }
23  out << "]";
24  return out;
25 }
26 
27 } // end unnamed namespace
28 
29 //**********************************************************************
30 
33 : m_ntick(0),
34  fFirstCollectionChannel(9999999) {
35  fFractUUCollect = pset.get< std::vector<float> >("FractUUCollect");
36  fFractUVCollect = pset.get< std::vector<float> >("FractUVCollect");
37  fFractVUCollect = pset.get< std::vector<float> >("FractVUCollect");
38  fFractVVCollect = pset.get< std::vector<float> >("FractVVCollect");
39  fFractUUMiss = pset.get< std::vector<float> >("FractUUMiss");
40  fFractUVMiss = pset.get< std::vector<float> >("FractUVMiss");
41  fFractVUMiss = pset.get< std::vector<float> >("FractVUMiss");
42  fFractVVMiss = pset.get< std::vector<float> >("FractVVMiss");
43  fFractZUMiss = pset.get< std::vector<float> >("FractZUMiss");
44  fFractZVMiss = pset.get< std::vector<float> >("FractZVMiss");
45  fFractHorizGapUMiss = pset.get< std::vector<float> >("FractHorizGapUMiss");
46  fFractVertGapUMiss = pset.get< std::vector<float> >("FractVertGapUMiss");
47  fFractHorizGapVMiss = pset.get< std::vector<float> >("FractHorizGapVMiss");
48  fFractVertGapVMiss = pset.get< std::vector<float> >("FractVertGapVMiss");
49  fFractHorizGapZMiss = pset.get< std::vector<float> >("FractHorizGapZMiss");
50  fFractVertGapZMiss = pset.get< std::vector<float> >("FractVertGapZMiss");
51  fFractHorizGapUCollect = pset.get< std::vector<float> >("FractHorizGapUCollect");
52  fFractVertGapUCollect = pset.get< std::vector<float> >("FractVertGapUCollect");
53  fFractHorizGapVCollect = pset.get< std::vector<float> >("FractHorizGapVCollect");
54  fFractVertGapVCollect = pset.get< std::vector<float> >("FractVertGapVCollect");
55  init();
56 }
57 
58 //**********************************************************************
59 
62  const sim::SimChannel* psc, AdcSignalVector& fChargeWork) const {
63  int dflag = 0;
64  fChargeWork.clear();
65  fChargeWork.resize(m_ntick, 0.0);
66  if ( psc == nullptr ) return 0;
67  AdcSignalVector fChargeWorkCollInd(m_ntick, 0.0);
68  string fname = "Dune35tSimChannelExtractService::extract";
70  unsigned int chan = psc->Channel();
71  const geo::View_t view = geo->View(chan);
72  fChargeWorkCollInd.clear();
73  fChargeWorkCollInd.resize(fChargeWork.size(), 0.);
74  for ( size_t t=0; t<fChargeWork.size(); ++t ) {
75  const std::vector<sim::IDE> ides = psc->TrackIDsAndEnergies(t,t);
76  for ( auto const &ide : ides ) {
77  GapType_t gaptype = combtest35t(ide.x,ide.y,ide.z);
78  switch ( gaptype ) {
79  case ACTIVE:
80  fChargeWork[t] += ide.numElectrons;
81  break;
82  case UCOMB:
83  {
84  dflag = GapHasDeflector(ide.x,ide.y,ide.z);
85  switch (view) {
86  case geo::kU:
87  fChargeWork[t] += ide.numElectrons * (1.0 - fFractUUCollect[dflag] - fFractUUMiss[dflag]);
88  fChargeWorkCollInd[t] += ide.numElectrons * fFractUUCollect[dflag];
89  break;
90  case geo::kV:
91  fChargeWork[t] += ide.numElectrons * (1.0 - fFractVUCollect[dflag] - fFractUUCollect[dflag] - fFractVUMiss[dflag]);
92  fChargeWorkCollInd[t] += ide.numElectrons * fFractVUCollect[dflag];
93  break;
94  case geo::kZ:
95  fChargeWork[t] += ide.numElectrons * (1.0-fFractVUCollect[dflag]-fFractUUCollect[dflag]-fFractZUMiss[dflag]);
96  break;
97  default:
98  throw cet::exception(fname) << "ILLEGAL VIEW Type: " << view <<"\n";
99  }
100  break;
101  }
102  case VCOMB:
103  {
104  dflag = GapHasDeflector(ide.x,ide.y,ide.z);
105  switch (view) {
106  case geo::kU:
107  fChargeWork[t] += ide.numElectrons * (1.0 - fFractUVCollect[dflag] - fFractUVMiss[dflag]);
108  fChargeWorkCollInd[t] += ide.numElectrons * fFractUVCollect[dflag];
109  break;
110  case geo::kV:
111  fChargeWork[t] += ide.numElectrons * (1.0 - fFractVVCollect[dflag] - fFractUVCollect[dflag] - fFractVVMiss[dflag]);
112  fChargeWorkCollInd[t] += ide.numElectrons * fFractVVCollect[dflag];
113  break;
114  case geo::kZ:
115  fChargeWork[t] += ide.numElectrons * (1.0-fFractVVCollect[dflag]-fFractUVCollect[dflag]-fFractZVMiss[dflag]);
116  break;
117  default:
118  throw cet::exception(fname) << "ILLEGAL VIEW Type: " << view <<"\n";
119  }
120  break;
121  }
122  case HORIZGAP:
123  {
124  dflag = GapHasDeflector(ide.x,ide.y,ide.z);
125  switch (view) {
126  case geo::kU:
127  fChargeWork[t] += ide.numElectrons * (1.0-fFractHorizGapUMiss[dflag]-fFractHorizGapUCollect[dflag]);
128  fChargeWorkCollInd[t] += ide.numElectrons * fFractHorizGapUCollect[dflag];
129  break;
130  case geo::kV:
131  fChargeWork[t] += ide.numElectrons * (1.0-fFractHorizGapVMiss[dflag]-fFractHorizGapUCollect[dflag]-fFractHorizGapVCollect[dflag]);
132  fChargeWorkCollInd[t] += ide.numElectrons * fFractHorizGapVCollect[dflag];
133  break;
134  case geo::kZ:
135  fChargeWork[t] += ide.numElectrons * (1.0-fFractHorizGapZMiss[dflag]-fFractHorizGapUCollect[dflag]-fFractHorizGapVCollect[dflag]);
136  break;
137  default:
138  throw cet::exception(fname) << "ILLEGAL VIEW Type: " << view <<"\n";
139  }
140  break;
141  }
142  case VERTGAP:
143  {
144  dflag = GapHasDeflector(ide.x,ide.y,ide.z);
145  switch (view) {
146  case geo::kU:
147  fChargeWork[t] += ide.numElectrons * (1.0-fFractVertGapUMiss[dflag]-fFractVertGapUCollect[dflag]);
148  fChargeWorkCollInd[t] += ide.numElectrons * fFractVertGapUCollect[dflag];
149  break;
150  case geo::kV:
151  fChargeWork[t] += ide.numElectrons * (1.0-fFractVertGapVMiss[dflag]-fFractVertGapUCollect[dflag]-fFractVertGapVCollect[dflag]);
152  fChargeWorkCollInd[t] += ide.numElectrons * fFractVertGapVCollect[dflag];
153  break;
154  case geo::kZ:
155  fChargeWork[t] += ide.numElectrons * (1.0-fFractVertGapZMiss[dflag]-fFractVertGapUCollect[dflag]-fFractVertGapVCollect[dflag]);
156  break;
157  default:
158  throw cet::exception(fname) << "ILLEGAL VIEW Type: " << view <<"\n";
159  }
160  break;
161  }
162  case NONACTIVE:
163  break;
164  } // end switch gaptype
165  } // end loop over ides
166  } // end loop over ticks
167  // Convolute and combine charges.
168  m_psss->Convolute(clockData, chan, fChargeWork);
169  m_psss->Convolute(clockData, fFirstCollectionChannel, fChargeWorkCollInd);
170  for ( unsigned int itck=0; itck<m_ntick; ++itck ) {
171  fChargeWork[itck] += fChargeWorkCollInd[itck];
172  }
173  return 0;
174 }
175 
176 //**********************************************************************
177 
178 std::ostream& Dune35tSimChannelExtractService::print(std::ostream& out, std::string prefix) const {
179  out << prefix << "Dune35tSimChannelExtractService";
180  out << prefix << " FractUUCollect: " << fFractUUCollect << endl;
181  out << prefix << " FractUVCollect: " << fFractUVCollect << endl;
182  out << prefix << " FractVUCollect: " << fFractVUCollect << endl;
183  out << prefix << " FractVVCollect: " << fFractVVCollect << endl;
184  out << prefix << " FractUUMiss: " << fFractUUMiss << endl;
185  out << prefix << " FractUVMiss: " << fFractUVMiss << endl;
186  out << prefix << " FractVUMiss: " << fFractVUMiss << endl;
187  out << prefix << " FractVVMiss: " << fFractVVMiss << endl;
188  out << prefix << " FractZUMiss: " << fFractZUMiss << endl;
189  out << prefix << " FractZVMiss: " << fFractZVMiss << endl;
190  out << prefix << " FractHorizGapUMiss: " << fFractHorizGapUMiss << endl;
191  out << prefix << " FractVertGapUMiss: " << fFractVertGapUMiss << endl;
192  out << prefix << " FractHorizGapVMiss: " << fFractHorizGapVMiss << endl;
193  out << prefix << " FractVertGapVMiss: " << fFractVertGapVMiss << endl;
194  out << prefix << " FractHorizGapZMiss: " << fFractHorizGapZMiss << endl;
195  out << prefix << " FractVertGapZMiss: " << fFractVertGapZMiss << endl;
196  out << prefix << " FractHorizGapUCollect: " << fFractHorizGapUCollect << endl;
197  out << prefix << " FractVertGapUCollect: " << fFractVertGapUCollect << endl;
198  out << prefix << " FractHorizGapVCollect: " << fFractHorizGapVCollect << endl;
199  out << prefix << " FractVertGapVCollect: " << fFractVertGapVCollect;
200  return out;
201 }
202 
203 //**********************************************************************
204 
206 
207  if ( m_init ) return;
208 
209  // Fetch the number of ticks.
210  m_ntick = m_pfft->FFTSize();
211  if ( m_ntick%2 != 0 )
212  throw cet::exception("SimChannelExtractAllService")
213  << "FFT size is not a power of 2.";
215 
216  // Find the first collection channel.
217  for ( uint32_t ichan=0; ichan<geo->Nchannels(); ++ichan ) {
218  if ( geo->View(ichan) == geo::kZ ) {
219  fFirstCollectionChannel = ichan;
220  break;
221  }
222  }
223 
224  // initialize the comb test positions. This is clumsy here mainly due to the irregular geometry
225  // should write something more systematic for the FD. There is also some duplication as the
226  // vertical positions of APA's 0 and 3 are assumed to be the same. Could think about either adding
227  // an exception if they're not, or defining more y positions to hold different APA positions if we want
228  // them to be different at a later time. Simulation may always be perfect though.
229 
230  // WireEndPoints takes cryostat, tpc, plane, wire, as ints and returns endpoints
231  //geo->WireEndPoints(c,t,p,w,xyzbeg,xyzend);
232 
233  // wire endpoints are at the places where the wire hits the comb that supports it. Bits of
234  // wire running over the comb are beyond the endpoints. So we need to extrapolate.
235 
236  double xyzbeg[3],xyzend[3];
237  int lastwire = 0;
238 
239  // Numbers in comments are from Geometry V3 for debugging purposes.
240 
241  // APA 0
242 
243  geo->WireEndPoints(0,0,0,0,xyzbeg,xyzend); // first U wire in TPC 0.
244  zcomb2 = xyzbeg[2]; // 0.0
245  ycomb5 = xyzend[1]; // 113.142
246 
247  lastwire = geo->Nwires(0,0,0)-1; // 358 in v3
248  geo->WireEndPoints(0,0,0,lastwire,xyzbeg,xyzend); // last U wire in TPC 0.
249  zcomb5 = xyzend[2]; // 50.8929
250  ycomb2 = xyzbeg[1]; // -82.9389
251 
252  geo->WireEndPoints(0,0,1,0,xyzbeg,xyzend); // first V wire in TPC 0.
253  zcomb4 = xyzend[2]; // 50.5774
254  ycomb4 = xyzbeg[1]; // 113.142
255 
256  lastwire = geo->Nwires(1,0,0)-1; // 344 in v3
257  geo->WireEndPoints(0,0,1,lastwire,xyzbeg,xyzend); // last V wire in TPC 0.
258  zcomb3 = xyzbeg[2]; // 0.3155
259  ycomb3 = xyzend[1]; // -82.6234
260 
261  // the collection wires appear to end where they meet their comb.
262  //geo->WireEndPoints(0,0,2,0,xyzbeg,xyzend); // first collection wire in TPC 0
263  //ycomb3 = xyzbeg[2]; // -82.308
264  //ycomb4 = xyzend[2]; // 113.142
265 
266  // need to get zcomb1, zcomb6, ycomb1, and ycomb6 -- extrapolate
267 
268  zcomb1 = zcomb2 - (zcomb3 - zcomb2);
269  zcomb6 = zcomb5 + (zcomb5 - zcomb4);
270  ycomb1 = ycomb2 - (ycomb3 - ycomb2);
271  ycomb6 = ycomb5 + (ycomb5 - ycomb4);
272 
273 
274  // APA 1
275 
276  geo->WireEndPoints(0,2,0,0,xyzbeg,xyzend); // first U wire in TPC 2.
277  zcomb11 = xyzend[2]; // 102.817
278  ycomb8 = xyzbeg[1]; // -85.221
279 
280  lastwire = geo->Nwires(0,2,0)-1; // 194 in v3
281  geo->WireEndPoints(0,2,0,lastwire,xyzbeg,xyzend); // last U wire in TPC 2.
282  zcomb8 = xyzbeg[2]; // 51.924
283  ycomb11 = xyzend[1]; // -0.831
284 
285  geo->WireEndPoints(0,2,1,0,xyzbeg,xyzend); // first V wire in TPC 2.
286  zcomb9 = xyzbeg[2]; // 52.2395
287  ycomb9 = xyzend[1]; // -85.222
288 
289  lastwire = geo->Nwires(1,2,0)-1; // 188 in v3
290  geo->WireEndPoints(0,2,1,lastwire,xyzbeg,xyzend); // last V wire in TPC 2.
291  zcomb10 = xyzend[2]; // 102.501
292  ycomb10 = xyzbeg[1]; // -1.14655
293 
294  // need to get zcomb7, zcomb12, ycomb7, and ycomb12 -- extrapolate
295 
296  zcomb7 = zcomb8 - (zcomb9 - zcomb8);
297  zcomb12 = zcomb11 + (zcomb11 - zcomb10);
298  ycomb7 = ycomb8 - (ycomb9 - ycomb8);
299  ycomb12 = ycomb11 + (ycomb11 - ycomb10);
300 
301  // APA 2
302 
303  geo->WireEndPoints(0,4,0,0,xyzbeg,xyzend); // first U wire in TPC 4.
304  zcomb8 = xyzbeg[2]; // 51.924 -- same as above
305  ycomb17 = xyzend[1]; // 113.142 -- same as above
306 
307  lastwire = geo->Nwires(0,4,0)-1; // 235 in v3
308  geo->WireEndPoints(0,4,0,lastwire,xyzbeg,xyzend); // last U wire in TPC 4.
309  zcomb11 = xyzend[2]; // 102.817 -- same as above
310  ycomb14 = xyzbeg[1]; // 0.83105
311 
312  geo->WireEndPoints(0,4,1,0,xyzbeg,xyzend); // first V wire in TPC 4.
313  zcomb10 = xyzend[2]; // 102.501 -- same as above
314  ycomb16 = xyzbeg[1]; // 113.142 -- everything ends here in y
315 
316  lastwire = geo->Nwires(1,4,0)-1; // 227 in v3
317  geo->WireEndPoints(0,4,1,lastwire,xyzbeg,xyzend); // last V wire in TPC 4.
318  zcomb9 = xyzbeg[2]; // 52.2395 -- same as above
319  ycomb15 = xyzend[1]; // 1.14655
320 
321  zcomb7 = zcomb8 - (zcomb9 - zcomb8);
322  zcomb12 = zcomb11 + (zcomb11 - zcomb10);
323  ycomb13 = ycomb14 - (ycomb15 - ycomb14);
324  ycomb18 = ycomb17 + (ycomb17 - ycomb16);
325 
326  // APA 3 -- a lot like APA 0
327 
328  geo->WireEndPoints(0,6,0,0,xyzbeg,xyzend); // first U wire in TPC 6.
329  zcomb14 = xyzbeg[2]; // 103.84
330  ycomb5 = xyzend[1]; // 113.142 -- same as above
331 
332  lastwire = geo->Nwires(0,6,0)-1; // 358 in v3
333  geo->WireEndPoints(0,6,0,lastwire,xyzbeg,xyzend); // last U wire in TPC 6.
334  zcomb17 = xyzend[2]; // 154.741
335  ycomb2 = xyzbeg[1]; // -82.9389 -- same as above
336 
337  geo->WireEndPoints(0,6,1,0,xyzbeg,xyzend); // first V wire in TPC 6.
338  zcomb16 = xyzend[2]; // 154.425
339  ycomb4 = xyzbeg[1]; // 113.142 -- same as above
340 
341  lastwire = geo->Nwires(1,6,0)-1; // 344 in v3
342  geo->WireEndPoints(0,6,1,lastwire,xyzbeg,xyzend); // last V wire in TPC 6.
343  zcomb15 = xyzbeg[2]; // 104.164
344  ycomb3 = xyzend[1]; // -82.6234 -- same as above
345 
346  // need to get zcomb13, zcomb18, ycomb1, and ycomb6 -- extrapolate
347  // the ycomb1 and ycomb6 are just copies.
348 
349  zcomb13 = zcomb14 - (zcomb15 - zcomb14);
350  zcomb18 = zcomb17 + (zcomb17 - zcomb16);
351  ycomb1 = ycomb2 - (ycomb3 - ycomb2);
352  ycomb6 = ycomb5 + (ycomb5 - ycomb4);
353 
354  m_init = true;
355 
356 }
357 
358 //**********************************************************************
359 
360 // see the ASCII cartoon of APA's at the bottom of this file for a picture of what all the boundaries are
361 
363 Dune35tSimChannelExtractService::combtest35t(double x, double y, double z) const {
364  if (z<zcomb1) return VERTGAP; // off to the side of the first APA -- kind of like being in a vertical gap
365  if (z<zcomb2) return UCOMB; // over U comb
366  if (z<zcomb3) return VCOMB; // over V comb
367  if (z<zcomb4) {
368  if (y<ycomb1) return HORIZGAP; // below the bottom
369  if (y<ycomb2) return UCOMB; // over U comb
370  if (y<ycomb3) return VCOMB; // over V comb
371  if (y<ycomb4) return ACTIVE; // active volume
372  if (y<ycomb5) return VCOMB; // over V comb
373  if (y<ycomb6) return UCOMB; // over U comb
374  return HORIZGAP; // outside top edge
375  }
376  if (z<zcomb5) return VCOMB; // over V comb
377  if (z<zcomb6) return UCOMB; // over U comb
378  if (z<zcomb7) return VERTGAP; // in gap
379  if (z<zcomb8) return UCOMB; // over U comb
380  if (z<zcomb9) return VCOMB; // over V comb
381  if (z<zcomb10) {
382  if (y<ycomb7) return HORIZGAP; // off the bottom
383  if (y<ycomb8) return UCOMB; // over U comb
384  if (y<ycomb9) return VCOMB; // over V comb
385  if (y<ycomb10) return ACTIVE; // active
386  if (y<ycomb11) return VCOMB; // over V comb
387  if (y<ycomb12) return UCOMB; // over U comb
388  if (y<ycomb13) return HORIZGAP; // over gap
389  if (y<ycomb14) return UCOMB; // over U comb
390  if (y<ycomb15) return VCOMB; // over V comb
391  if (y<ycomb16) return ACTIVE; // active volume
392  if (y<ycomb17) return VCOMB; // over V comb
393  if (y<ycomb18) return UCOMB; // over U comb
394  return HORIZGAP; // above the top edge
395  }
396  if (z<zcomb11) return VCOMB; // over V comb
397  if (z<zcomb12) return UCOMB; // over U comb
398 
399  if (z<zcomb13) return VERTGAP; // outside first APA
400  if (z<zcomb14) return UCOMB; // over U comb
401  if (z<zcomb15) return VCOMB; // over V comb
402  if (z<zcomb16) {
403  if (y<ycomb1) return HORIZGAP; // below the bottom
404  if (y<ycomb2) return UCOMB; // over U comb
405  if (y<ycomb3) return VCOMB; // over V comb
406  if (y<ycomb4) return ACTIVE; // active volume
407  if (y<ycomb5) return VCOMB; // over V comb
408  if (y<ycomb6) return UCOMB; // over U comb
409  return HORIZGAP; // outside top edge
410  }
411  if (z<zcomb17) return VCOMB; // over V comb
412  if (z<zcomb18) return UCOMB; // over U comb
413  return VERTGAP; // off the end in Z.
414 }
415 
416  // see the ASCII cartoon of APA's at the bottom of this file for a picture of what all the boundaries are
417  // returns 0 if this gap does not have a deflector as described by Bo Yu, LBNE DocDB 10073. Returns 1 if this gap does.
418  // Also returns 0 if we are not in a gap. This is an int instead of a bool so it can be used as the index into the parameter array
419  // the deflector is between APA's 0 and 1 in the numbering corresponding to the cartoon below
420 
421  int Dune35tSimChannelExtractService::GapHasDeflector(double x, double y, double z) const
422  {
423  if ( y < ycomb12 && y > ycomb7 && x > 0 && z < zcomb9 && z > zcomb4 ) return 1;
424  return 0;
425  }
426 
427 
428 /* -------------------------------------------------
429  APA Cartoons for the combtest35t method
430 
431  z->
432  ^
433  |
434  y
435 
436 
437  zcomb1 zcomb6
438  zcomb2 zcomb5
439  zcomb3 zcomb4
440  ______________________________________________ ycomb6
441  |____________________________________________| ycomb5
442  ||__________________________________________|| ycomb4
443  ||| |||
444  ||| |||
445  ||| |||
446  ||| |||
447  ||| |||
448  ||| |||
449  ||| |||
450  ||| |||
451  ||| |||
452  ||| |||
453  ||| |||
454  ||| |||
455  ||| |||
456  ||| |||
457  ||| |||
458  ||| APA0 |||
459  ||| |||
460  ||| |||
461  ||| |||
462  ||| |||
463  ||| |||
464  ||| |||
465  ||| |||
466  ||| |||
467  ||| |||
468  ||| |||
469  ||| |||
470  ||| |||
471  ||| |||
472  ||| |||
473  ||| |||
474  ||| |||
475  ||| |||
476  ||| |||
477  ||__________________________________________|| ycomb3
478  |____________________________________________| ycomb2
479  ______________________________________________ ycomb1
480 
481 
482  z->
483 
484  ^
485  |
486  y
487 
488 
489  zcomb7 zcomb12
490  zcomb8 zcomb11
491  zcomb9 zcomb10
492  ______________________________________________ ycomb18
493  |____________________________________________| ycomb17
494  ||__________________________________________|| ycomb16
495  ||| |||
496  ||| |||
497  ||| |||
498  ||| |||
499  ||| |||
500  ||| |||
501  ||| |||
502  ||| |||
503  ||| |||
504  ||| |||
505  ||| |||
506  ||| APA2 |||
507  ||| |||
508  ||| |||
509  ||| |||
510  ||| |||
511  ||| |||
512  ||| |||
513  ||| |||
514  ||| |||
515  ||| |||
516  ||| |||
517  ||| |||
518  ||__________________________________________|| ycomb15
519  |____________________________________________| ycomb14
520  ______________________________________________ ycomb13
521 
522  ______________________________________________ ycomb12
523  |____________________________________________| ycomb11
524  ||__________________________________________|| ycomb10
525  ||| |||
526  ||| |||
527  ||| |||
528  ||| |||
529  ||| |||
530  ||| APA1 |||
531  ||| |||
532  ||| |||
533  ||| |||
534  ||| |||
535  ||| |||
536  ||__________________________________________|| ycomb9
537  |____________________________________________| ycomb8
538  ______________________________________________ ycomb7
539 
540 
541  APA 0 Cartoon:
542 
543  z->
544 
545  ^
546  |
547  y
548 
549 
550  zcomb13 zcomb18
551  zcomb14 zcomb17
552  zcomb15 zcomb16
553  ______________________________________________ ycomb6
554  |____________________________________________| ycomb5
555  ||__________________________________________|| ycomb4
556  ||| |||
557  ||| |||
558  ||| |||
559  ||| |||
560  ||| |||
561  ||| |||
562  ||| |||
563  ||| |||
564  ||| |||
565  ||| |||
566  ||| |||
567  ||| |||
568  ||| |||
569  ||| |||
570  ||| |||
571  ||| |||
572  ||| APA3 |||
573  ||| |||
574  ||| |||
575  ||| |||
576  ||| |||
577  ||| |||
578  ||| |||
579  ||| |||
580  ||| |||
581  ||| |||
582  ||| |||
583  ||| |||
584  ||| |||
585  ||| |||
586  ||| |||
587  ||| |||
588  ||| |||
589  ||| |||
590  ||__________________________________________|| ycomb3
591  |____________________________________________| ycomb2
592  ______________________________________________ ycomb1
593 
594 
595 
596 */
597 
598 //**********************************************************************
599 
601 
602 //**********************************************************************
Energy deposited on a readout channel by simulated tracks.
Definition: SimChannel.h:140
Dune35tSimChannelExtractService(fhicl::ParameterSet const &pset, art::ActivityRegistry &)
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
std::string string
Definition: nybbler.cc:12
Planes which measure V.
Definition: geo_types.h:130
int GapHasDeflector(double x, double y, double z) const
Planes which measure Z direction.
Definition: geo_types.h:132
unsigned int Nwires(unsigned int p, unsigned int tpc=0, unsigned int cstat=0) const
Returns the total number of wires in the specified plane.
art framework interface to geometry description
unsigned int Nchannels() const
Returns the number of TPC readout channels in the detector.
art::ServiceHandle< util::SignalShapingServiceDUNE > m_psss
Planes which measure U.
Definition: geo_types.h:129
int FFTSize() const
Definition: LArFFT.h:69
art::ServiceHandle< util::LArFFT > m_pfft
int extract(detinfo::DetectorClocksData const &clockData, const sim::SimChannel *psc, AdcSignalVector &sig) const
T get(std::string const &key) const
Definition: ParameterSet.h:271
GapType_t combtest35t(double x, double y, double z) const
View_t View(geo::PlaneID const &pid) const
Returns the view (wire orientation) on the channels of specified TPC plane.
raw::ChannelID_t Channel() const
Returns the readout channel this object describes.
Definition: SimChannel.h:329
void Convolute(detinfo::DetectorClocksData const &clockData, Channel channel, std::vector< T > &func) const
std::ostream & print(std::ostream &out=std::cout, std::string prefix="") const
Contains all timing reference information for the detector.
std::vector< sim::IDE > TrackIDsAndEnergies(TDC_t startTDC, TDC_t endTDC) const
Return all the recorded energy deposition within a time interval.
Definition: SimChannel.cxx:180
void WireEndPoints(geo::WireID const &wireid, double *xyzStart, double *xyzEnd) const
Fills two arrays with the coordinates of the wire end points.
std::vector< AdcSignal > AdcSignalVector
Definition: AdcTypes.h:22
list x
Definition: train.py:276
LArSoft geometry interface.
Definition: ChannelGeo.h:16
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
QTextStream & endl(QTextStream &s)
#define DEFINE_ART_SERVICE_INTERFACE_IMPL(svc, iface)