Rivet Analyses Reference

ATLAS_2014_I1282447

W + charm production at 7 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1282447
Status: VALIDATED
Authors:
  • Kristin Lohwasser
References:Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • W+c and inclusive double production in pp collisions

This routine implements the measurement of $W$ boson production in association with a single charm quark with $4.6 \text{fb}^{-1}$ of $pp$ collision data at $\sqrt{s} = 7$ TeV collected with the ATLAS detector at the Large Hadron Collider. Results are quoted for the $W$ boson production decaying into either an electron or muon and the charm quark being identified as any charmed hadron above 5 GeV inside $\Delta R = 0.4$ of a jet with more than 25 GeV. Alternatively the presence of the charm quark is indicated by the presence of a $D$ or a $D^*$ meson with a $p_\text{T}$ above 8 GeV. The cross sections are quoted for the number of opposite sign minus same sign events, where the signs under consideration are the charge of the $W$ boson and the charmed hadrons tagging the event. Given are the integrated and differential cross sections as a function of the pseudorapidity of the lepton from the $W$-boson decay are measured. Additionally, the cross-section ratios $\sigma(W^+ + bar{c}) / \sigma(W^- + c)$ as well as the $p_\text{T}$ dependent cross sections of the $D$/$D^*$ mesons normalized to the $W$ inclusive cross sections are published. The measured data is unfolded to the Born level. One should therefore take care to run on samples without QED radiation off of the electrons. IMPORTANT NOTICE --- For the MC predictions in the paper, the branching fractions to $D$ and $D^*$ mesons have been adapted to be 0.2256 ($D$) and 0.2287 ($D^*$) (LEP/HERA combination). Some suggestion on how to do post-processing -- in case separate samples for $W^+ + c$ and $W^- + c$ and inclusive $W$ are used -- are included in the cc-file. This post-processing is needed to properly display the histograms for the cross section ratio plots $\sigma(W^+ + bar{c}) / \sigma(W^- + c)$ as well as for the cross section ratios of $W + D^{(*)}$ production over inclusive $W$ production ($\sigma(W^{+/-} D^{(*)}) / \sigma(W^{+/-})$) as a function of the $D^{(*)}$ meson transverse momentum.

Source code: ATLAS_2014_I1282447.cc
  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
// -*- C++ -*-
// ATLAS W+c analysis

//////////////////////////////////////////////////////////////////////////
/*
  Description of rivet analysis ATLAS_2014_I1282447 W+c production

  This rivet routine implements the ATLAS W+c analysis.
  Apart from those histograms, described and published on HEP Data, here
  are some helper histograms defined, these are:

  d02-x01-y01, d02-x01-y02 and d08-x01-y01 are ratios, the nominator ("_plus")
  and denominator ("_minus") histograms are also given, so that the ratios can
  be reconstructed if need be (e.g. when running on separate samples).

  d05 and d06 are ratios over inclusive W production.
  The routine has to be run on a sample for inclusive W production in order to
  make sure the denominator ("_winc") is correctly filled.

  The ratios can be constructed using the following sample code:
  python divideWCharm.py

  import yoda
  hists_wc   = yoda.read("Rivet_Wc.yoda")
  hists_winc = yoda.read("Rivet_Winc.yoda")

  ## division histograms --> ONLY for different plus minus runs
  # (merge before using yodamerge Rivet_plus.yoda Rivet_minus.yoda > Rivet_Wc.yoda)

  d02y01_plus = hists_wc["/ATLAS_2014_I1282447/d02-x01-y01_plus"]
  d02y01_minus = hists_wc["/ATLAS_2014_I1282447/d02-x01-y01_minus"]
  ratio_d02y01 =  d02y01_plus.divide(d02y01_minus)
  ratio_d02y01.path = "/ATLAS_2014_I1282447/d02-x01-y01"

  d02y02_plus = hists_wc["/ATLAS_2014_I1282447/d02-x01-y02_plus"]
  d02y02_minus = hists_wc["/ATLAS_2014_I1282447/d02-x01-y02_minus"]
  ratio_d02y02=  d02y02_plus.divide(d02y02_minus)
  ratio_d02y02.path = "/ATLAS_2014_I1282447/d02-x01-y02"

  d08y01_plus = hists_wc["/ATLAS_2014_I1282447/d08-x01-y01_plus"]
  d08y01_minus = hists_wc["/ATLAS_2014_I1282447/d08-x01-y01_minus"]
  ratio_d08y01=  d08y01_plus.divide(d08y01_minus)
  ratio_d08y01.path = "/ATLAS_2014_I1282447/d08-x01-y01"

  # inclusive cross section
  h_winc = hists_winc["/ATLAS_2014_I1282447/d05-x01-y01"]
  h_d    = hists_wc["/ATLAS_2014_I1282447/d01-x01-y02"]
  h_dstar= hists_wc["/ATLAS_2014_I1282447/d01-x01-y03"]

  ratio_wd      =  h_d.divide(h_winc)
  ratio_wd.path = "/ATLAS_2014_I1282447/d05-x01-y02"

  ratio_wdstar      =  h_d.divide(h_winc)
  ratio_wdstar.path = "/ATLAS_2014_I1282447/d05-x01-y03"

  # pT differential
  h_winc_plus  = hists_winc["/ATLAS_2014_I1282447/d06-x01-y01_winc"]
  h_winc_minus = hists_winc["/ATLAS_2014_I1282447/d06-x01-y02_winc"]

  h_wd_plus      = hists_wc["/ATLAS_2014_I1282447/d06-x01-y01_wplus"]
  h_wd_minus     = hists_wc["/ATLAS_2014_I1282447/d06-x01-y02_wminus"]
  h_wdstar_plus  = hists_wc["/ATLAS_2014_I1282447/d06-x01-y03_wplus"]
  h_wdstar_minus = hists_wc["/ATLAS_2014_I1282447/d06-x01-y04_wminus"]

  ratio_wd_plus       =  h_wd_plus.divide(h_winc_plus)
  ratio_wd_plus.path  = "/ATLAS_2014_I1282447/d06-x01-y01"
  ratio_wd_minus      =  h_wd_plus.divide(h_winc_minus)
  ratio_wd_minus.path = "/ATLAS_2014_I1282447/d06-x01-y02"

  ratio_wdstar_plus       =  h_wdstar_plus.divide(h_winc_plus)
  ratio_wdstar_plus.path  = "/ATLAS_2014_I1282447/d06-x01-y03"
  ratio_wdstar_minus      =  h_wdstar_plus.divide(h_winc_minus)
  ratio_wdstar_minus.path = "/ATLAS_2014_I1282447/d06-x01-y04"

  ratio_wd_plus =  h_wd_plus.divide(h_winc_plus)
  ratio_wd_plus.path = "/ATLAS_2014_I1282447/d06-x01-y01"
  ratio_wd_minus =  h_wd_plus.divide(h_winc_minus)
  ratio_wd_minus.path = "/ATLAS_2014_I1282447/d06-x01-y02"

  h_winc_plus= hists_winc["/ATLAS_2014_I1282447/d06-x01-y01_winc"]
  h_winc_minus= hists_winc["/ATLAS_2014_I1282447/d06-x01-y02_winc"]

  ## copy other histograms for plotting

  d01x01y01= hists_wc["/ATLAS_2014_I1282447/d01-x01-y01"]
  d01x01y01.path = "/ATLAS_2014_I1282447/d01-x01-y01"

  d01x01y02= hists_wc["/ATLAS_2014_I1282447/d01-x01-y02"]
  d01x01y02.path = "/ATLAS_2014_I1282447/d01-x01-y02"

  d01x01y03= hists_wc["/ATLAS_2014_I1282447/d01-x01-y03"]
  d01x01y03.path = "/ATLAS_2014_I1282447/d01-x01-y03"

  d03x01y01= hists_wc["/ATLAS_2014_I1282447/d03-x01-y01"]
  d03x01y01.path = "/ATLAS_2014_I1282447/d03-x01-y01"

  d03x01y02= hists_wc["/ATLAS_2014_I1282447/d03-x01-y02"]
  d03x01y02.path = "/ATLAS_2014_I1282447/d03-x01-y02"

  d04x01y01= hists_wc["/ATLAS_2014_I1282447/d04-x01-y01"]
  d04x01y01.path = "/ATLAS_2014_I1282447/d04-x01-y01"

  d04x01y02= hists_wc["/ATLAS_2014_I1282447/d04-x01-y02"]
  d04x01y02.path = "/ATLAS_2014_I1282447/d04-x01-y02"

  d04x01y03= hists_wc["/ATLAS_2014_I1282447/d04-x01-y03"]
  d04x01y03.path = "/ATLAS_2014_I1282447/d04-x01-y03"

  d04x01y04= hists_wc["/ATLAS_2014_I1282447/d04-x01-y04"]
  d04x01y04.path = "/ATLAS_2014_I1282447/d04-x01-y04"

  d07x01y01= hists_wc["/ATLAS_2014_I1282447/d07-x01-y01"]
  d07x01y01.path = "/ATLAS_2014_I1282447/d07-x01-y01"

  yoda.write([ratio_d02y01,ratio_d02y02,ratio_d08y01, ratio_wd ,ratio_wdstar,ratio_wd_plus,ratio_wd_minus ,ratio_wdstar_plus,ratio_wdstar_minus,d01x01y01,d01x01y02,d01x01y03,d03x01y01,d03x01y02,d04x01y01,d04x01y02,d04x01y03,d04x01y04,d07x01y01],"validation.yoda")

*/
//////////////////////////////////////////////////////////////////////////

#include "Rivet/Analysis.hh"
#include "Rivet/Projections/UnstableParticles.hh"
#include "Rivet/Projections/WFinder.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/ChargedFinalState.hh"
#include "Rivet/Projections/VetoedFinalState.hh"

namespace Rivet {



  class ATLAS_2014_I1282447 : public Analysis {
  public:

    /// Constructor
    ATLAS_2014_I1282447() : Analysis("ATLAS_2014_I1282447")
    {

    }


    /// @name Analysis methods
    //@{

    /// Book histograms and initialise projections before the run
    void init() {

      /// @todo Initialise and register projections here
      UnstableParticles fs;

      Cut cuts = Cuts::etaIn(-2.5, 2.5) & (Cuts::pT > 20*GeV);

      /// should use sample WITHOUT QED radiation off the electron
      WFinder wfinder_born_el(fs, cuts, PID::ELECTRON, 25*GeV, 8000*GeV, 15*GeV, 0.1, WFinder::ChargedLeptons::PROMPT, WFinder::ClusterPhotons::ALL, WFinder::AddPhotons::YES);
      declare(wfinder_born_el, "WFinder_born_el");

      WFinder wfinder_born_mu(fs, cuts, PID::MUON    , 25*GeV, 8000*GeV, 15*GeV, 0.1, WFinder::ChargedLeptons::PROMPT, WFinder::ClusterPhotons::ALL, WFinder::AddPhotons::YES);
      declare(wfinder_born_mu, "WFinder_born_mu");

      // all hadrons that could be coming from a charm decay --
      // -- for safety, use region -3.5 - 3.5
      declare(UnstableParticles(Cuts::abseta <3.5), "hadrons");

      // Input for the jets: no neutrinos, no muons, and no electron which passed the electron cuts
      // also: NO electron, muon or tau (needed due to ATLAS jet truth reconstruction feature)
      VetoedFinalState veto;

      veto.addVetoOnThisFinalState(wfinder_born_el);
      veto.addVetoOnThisFinalState(wfinder_born_mu);
      veto.addVetoPairId(PID::ELECTRON);
      veto.addVetoPairId(PID::MUON);
      veto.addVetoPairId(PID::TAU);

      FastJets jets(veto, FastJets::ANTIKT, 0.4);
      declare(jets, "jets");

      // Book histograms

      // charge separated integrated cross sections
      book(_hist_wcjet_charge  ,"d01-x01-y01");
      book(_hist_wd_charge     ,"d01-x01-y02");
      book(_hist_wdstar_charge ,"d01-x01-y03");

      // charge integrated total cross sections
      book(_hist_wcjet_ratio,"d02-x01-y01");
      book(_hist_wd_ratio   ,"d02-x01-y02");

      book(_hist_wcjet_minus ,"d02-x01-y01_minus");
      book(_hist_wd_minus    ,"d02-x01-y02_minus");

      book(_hist_wcjet_plus  ,"d02-x01-y01_plus");
      book(_hist_wd_plus     ,"d02-x01-y02_plus");

      // eta distributions
      book(_hist_wplus_wcjet_eta_lep   ,"d03-x01-y01");
      book(_hist_wminus_wcjet_eta_lep  ,"d03-x01-y02");

      book(_hist_wplus_wdminus_eta_lep ,"d04-x01-y01");
      book(_hist_wminus_wdplus_eta_lep ,"d04-x01-y02");
      book(_hist_wplus_wdstar_eta_lep  ,"d04-x01-y03");
      book(_hist_wminus_wdstar_eta_lep ,"d04-x01-y04");

      // ratio of cross section (WD over W inclusive) // postprocess!
      book(_hist_w_inc             ,"d05-x01-y01");
      book(_hist_wd_winc_ratio    ,"d05-x01-y02");
      book(_hist_wdstar_winc_ratio,"d05-x01-y03");

      // ratio of cross section (WD over W inclusive -- function of pT of D meson)
      book(_hist_wplusd_wplusinc_pt_ratio      ,"d06-x01-y01");
      book(_hist_wminusd_wminusinc_pt_ratio    ,"d06-x01-y02");
      book(_hist_wplusdstar_wplusinc_pt_ratio  ,"d06-x01-y03");
      book(_hist_wminusdstar_wminusinc_pt_ratio,"d06-x01-y04");

      // could use for postprocessing!
      book(_hist_wplusd_wplusinc_pt       ,"d06-x01-y01_wplus");
      book(_hist_wminusd_wminusinc_pt     ,"d06-x01-y02_wminus");
      book(_hist_wplusdstar_wplusinc_pt   ,"d06-x01-y03_wplus");
      book(_hist_wminusdstar_wminusinc_pt ,"d06-x01-y04_wminus");

      book(_hist_wplus_winc  ,"d06-x01-y01_winc");
      book(_hist_wminus_winc ,"d06-x01-y02_winc");

      // jet multiplicity of charge integrated W+cjet cross section (+0 or +1 jet in addition to the charm jet)
      book(_hist_wcjet_jets  ,"d07-x01-y01");

      // jet multiplicity of W+cjet cross section ratio (+0 or +1 jet in addition to the charm jet)
      book(_hist_wcjet_jets_ratio ,"d08-x01-y01");
      book(_hist_wcjet_jets_plus   ,"d08-x01-y01_plus");
      book(_hist_wcjet_jets_minus  ,"d08-x01-y01_minus");

    }


    /// Perform the per-event analysis
    void analyze(const Event& event) {

      double charge_weight = 0; // account for OS/SS events

      int    lepton_charge = 0;
      double lepton_eta    = 0.;

      /// Find leptons
      const WFinder& wfinder_born_el = apply<WFinder>(event, "WFinder_born_el");
      const WFinder& wfinder_born_mu = apply<WFinder>(event, "WFinder_born_mu");

      if (wfinder_born_el.empty() && wfinder_born_mu.empty()) {
        MSG_DEBUG("No W bosons found");
        vetoEvent;
      }

      bool keepevent = false;

      //check electrons
      if (!wfinder_born_el.empty()) {
        const FourMomentum nu = wfinder_born_el.constituentNeutrinos()[0];
        if (wfinder_born_el.mT() > 40*GeV && nu.pT() > 25*GeV) {
          keepevent = true;
          lepton_charge = wfinder_born_el.constituentLeptons()[0].charge();
          lepton_eta = fabs(wfinder_born_el.constituentLeptons()[0].pseudorapidity());
        }
      }

      //check muons
      if (!wfinder_born_mu.empty()) {
        const FourMomentum nu = wfinder_born_mu.constituentNeutrinos()[0];
        if (wfinder_born_mu.mT() > 40*GeV && nu.pT() > 25*GeV) {
          keepevent = true;
          lepton_charge = wfinder_born_mu.constituentLeptons()[0].charge();
          lepton_eta = fabs(wfinder_born_mu.constituentLeptons()[0].pseudorapidity());
        }
      }

      if (!keepevent) {
        MSG_DEBUG("Event does not pass mT and MET cuts");
        vetoEvent;
      }

      if (lepton_charge > 0) {
        _hist_wplus_winc->fill(10.);
        _hist_wplus_winc->fill(16.);
        _hist_wplus_winc->fill(30.);
        _hist_wplus_winc->fill(60.);
        _hist_w_inc->fill(+1);
      }
      else if (lepton_charge < 0) {
        _hist_wminus_winc->fill(10.);
        _hist_wminus_winc->fill(16.);
        _hist_wminus_winc->fill(30.);
        _hist_wminus_winc->fill(60.);
        _hist_w_inc->fill(-1);
      }

      // Find hadrons in the event
      const UnstableParticles& fs = apply<UnstableParticles>(event, "hadrons");

      /// FIND Different channels
      // 1: wcjet
      // get jets
      const Jets& jets = apply<FastJets>(event, "jets").jetsByPt(Cuts::pT>25.0*GeV && Cuts::abseta<2.5);
      // loop over jets to select jets used to match to charm
      Jets js;
      int    matched_charmHadron = 0;
      double charm_charge = 0.;
      int    njets = 0;
      int    nj = 0;
      bool   mat_jet = false;

      double ptcharm = 0;
      if (matched_charmHadron > -1) {
        for (const Jet& j : jets) {
          mat_jet = false;

          njets += 1;
          for (const Particle& p : fs.particles()) {
            /// @todo Avoid touching HepMC!
            ConstGenParticlePtr part = p.genParticle();
            if (p.hasCharm()) {
              //if(isFromBDecay(p)) continue;
              if (p.fromBottom()) continue;
              if (p.pT() < 5*GeV ) continue;
              if (hasCharmedChildren(part)) continue;
              if (deltaR(p, j) < 0.3) {
                mat_jet = true;
                if (p.pT() > ptcharm) {
                  charm_charge = part->pdg_id();
                  ptcharm = p.pT();
                }
              }
            }
          }
          if (mat_jet) nj++;
        }

        if (charm_charge * lepton_charge > 0)  charge_weight = -1;
        else charge_weight = +1;

        if (nj == 1)  {
          if (lepton_charge > 0) {
            _hist_wcjet_charge        ->fill(         1, charge_weight);
            _hist_wcjet_plus          ->fill(         0, charge_weight);
            _hist_wplus_wcjet_eta_lep ->fill(lepton_eta, charge_weight);
            _hist_wcjet_jets_plus     ->fill(njets-1   , charge_weight);
          }
          else if (lepton_charge < 0) {
            _hist_wcjet_charge        ->fill(        -1, charge_weight);
            _hist_wcjet_minus         ->fill(         0, charge_weight);
            _hist_wminus_wcjet_eta_lep->fill(lepton_eta, charge_weight);
            _hist_wcjet_jets_minus    ->fill(njets-1   , charge_weight);
          }

          _hist_wcjet_jets->fill(njets-1, charge_weight);
        }
      }

      // // 1/2: w+d(*) meson

      for (const Particle& p : fs.particles()) {

        /// @todo Avoid touching HepMC!
        ConstGenParticlePtr part = p.genParticle();
        if (p.pT() < 8*GeV)       continue;
        if (fabs(p.eta()) > 2.2)  continue;

        // W+D
        if (abs(part->pdg_id()) == 411) {
          if (lepton_charge * part->pdg_id() > 0)  charge_weight = -1;
          else                                     charge_weight = +1;

          // fill histos
          if (lepton_charge > 0) {
            _hist_wd_charge            ->fill(         1, charge_weight);
            _hist_wd_plus              ->fill(         0, charge_weight);
            _hist_wplus_wdminus_eta_lep->fill(lepton_eta, charge_weight);
            _hist_wplusd_wplusinc_pt   ->fill(    p.pT(), charge_weight);
          }
          else if (lepton_charge < 0) {
            _hist_wd_charge            ->fill(        -1, charge_weight);
            _hist_wd_minus             ->fill(         0, charge_weight);
            _hist_wminus_wdplus_eta_lep->fill(lepton_eta, charge_weight);
            _hist_wminusd_wminusinc_pt ->fill(p.pT()    , charge_weight);
          }
        }

        // W+Dstar
        if ( abs(part->pdg_id()) == 413 ) {
          if (lepton_charge*part->pdg_id() > 0) charge_weight = -1;
          else charge_weight = +1;

          if (lepton_charge > 0) {
            _hist_wdstar_charge->fill(+1, charge_weight);
            _hist_wd_plus->fill( 0, charge_weight);
            _hist_wplus_wdstar_eta_lep->fill( lepton_eta, charge_weight);
            _hist_wplusdstar_wplusinc_pt->fill(  p.pT(), charge_weight);
          }
          else if (lepton_charge < 0) {
            _hist_wdstar_charge->fill(-1, charge_weight);
            _hist_wd_minus->fill(0, charge_weight);
            _hist_wminus_wdstar_eta_lep->fill(lepton_eta, charge_weight);
            _hist_wminusdstar_wminusinc_pt->fill(p.pT(), charge_weight);
          }
        }
      }

    }


    /// Normalise histograms etc., after the run
    void finalize() {

      const double sf = crossSection() / sumOfWeights();

      // norm to cross section
      // d01
      scale(_hist_wcjet_charge,  sf);
      scale(_hist_wd_charge,     sf);
      scale(_hist_wdstar_charge, sf);

      //d02
      scale(_hist_wcjet_plus,  sf);
      scale(_hist_wcjet_minus, sf);
      scale(_hist_wd_plus,     sf);
      scale(_hist_wd_minus,    sf);

      divide(_hist_wcjet_plus, _hist_wcjet_minus, _hist_wcjet_ratio);
      divide(_hist_wd_plus,    _hist_wd_minus,    _hist_wd_ratio   );

      //d03
      scale(_hist_wplus_wcjet_eta_lep,  sf);
      scale(_hist_wminus_wcjet_eta_lep, sf);

      //d04
      scale(_hist_wplus_wdminus_eta_lep, crossSection()/sumOfWeights());
      scale(_hist_wminus_wdplus_eta_lep, crossSection()/sumOfWeights());

      scale(_hist_wplus_wdstar_eta_lep , crossSection()/sumOfWeights());
      scale(_hist_wminus_wdstar_eta_lep, crossSection()/sumOfWeights());

      //d05
      scale(_hist_w_inc, 0.01 * sf); // in percent --> /100
      divide(_hist_wd_charge,     _hist_w_inc, _hist_wd_winc_ratio    );
      divide(_hist_wdstar_charge, _hist_w_inc, _hist_wdstar_winc_ratio);

      //d06, in percentage!
      scale(_hist_wplusd_wplusinc_pt,       sf);
      scale(_hist_wminusd_wminusinc_pt,     sf);
      scale(_hist_wplusdstar_wplusinc_pt,   sf);
      scale(_hist_wminusdstar_wminusinc_pt, sf);

      scale(_hist_wplus_winc,  0.01 * sf); // in percent --> /100
      scale(_hist_wminus_winc, 0.01 * sf); // in percent --> /100

      divide(_hist_wplusd_wplusinc_pt,       _hist_wplus_winc , _hist_wplusd_wplusinc_pt_ratio      );
      divide(_hist_wminusd_wminusinc_pt,     _hist_wminus_winc, _hist_wminusd_wminusinc_pt_ratio    );
      divide(_hist_wplusdstar_wplusinc_pt,   _hist_wplus_winc , _hist_wplusdstar_wplusinc_pt_ratio  );
      divide(_hist_wminusdstar_wminusinc_pt, _hist_wminus_winc, _hist_wminusdstar_wminusinc_pt_ratio);


      //d07
      scale(_hist_wcjet_jets, sf);

      //d08
      scale(_hist_wcjet_jets_minus, sf);
      scale(_hist_wcjet_jets_plus,  sf);
      divide(_hist_wcjet_jets_plus, _hist_wcjet_jets_minus , _hist_wcjet_jets_ratio);
    }

    //@}


  private:

    // Data members like post-cuts event weight counters go here

    // Check whether particle comes from b-decay
    bool isFromBDecay(const Particle& p) {
      
      /// @todo I think we can just replicated the original behaviour with this call
      /// Note slight difference to Rivet's native Particle::fromBottom method!
      return p.hasAncestorWith([](const Particle &p)->bool{return p.hasBottom();});
      /*
      bool isfromB = false;

      if (p.genParticle() == nullptr)  return false;

      ConstGenParticlePtr part = p.genParticle();
      ConstGenVertexPtr ivtx = part->production_vertex();
      while (ivtx) {
        if (ivtx->particles_in().size() < 1) {
          isfromB = false;
          break;
        }
        const HepMC::GenVertex::particles_in_const_iterator iPart_invtx = ivtx->particles_in_const_begin();
        part = (*iPart_invtx);
        if (!part) {
          isfromB = false;
          break;
        }
        isfromB = PID::hasBottom(part->pdg_id());

        if (isfromB == true)  break;
        ivtx = part->production_vertex();
        if ( part->pdg_id() == 2212 || !ivtx )  break; // reached beam
      }
      return isfromB;
       */
    }


    // Check whether particle has charmed children
    /// @todo Use built-in method and avoid HepMC!
    bool hasCharmedChildren(ConstGenParticlePtr part) {

      bool hasCharmedChild = false;
      if (part == nullptr)  return false;

      ConstGenVertexPtr ivtx = part->end_vertex();
      if (ivtx == nullptr)  return false;

      // if (ivtx->particles_out_size() < 2) return false;
      //HepMC::GenVertex::particles_out_const_iterator iPart_invtx = ivtx->particles_out_const_begin();
      //HepMC::GenVertex::particles_out_const_iterator end_invtx = ivtx->particles_out_const_end();

      for(ConstGenParticlePtr p2: HepMCUtils::particles(ivtx, Relatives::CHILDREN)){
        if (p2 == part)  continue;
        hasCharmedChild = PID::hasCharm(p2->pdg_id());
        if (hasCharmedChild == true)  break;
        hasCharmedChild = hasCharmedChildren(p2);
        if (hasCharmedChild == true)  break;
      }
      return hasCharmedChild;
    }


  private:

    /// @name Histograms
    //@{

    //d01-x01-
    Histo1DPtr   _hist_wcjet_charge;
    Histo1DPtr   _hist_wd_charge;
    Histo1DPtr   _hist_wdstar_charge;

    //d02-x01-
    Scatter2DPtr _hist_wcjet_ratio;
    Scatter2DPtr _hist_wd_ratio;
    Histo1DPtr _hist_wcjet_plus;
    Histo1DPtr _hist_wd_plus;

    Histo1DPtr _hist_wcjet_minus;
    Histo1DPtr _hist_wd_minus;

    //d03-x01-
    Histo1DPtr _hist_wplus_wcjet_eta_lep;
    Histo1DPtr _hist_wminus_wcjet_eta_lep;

    //d04-x01-
    Histo1DPtr _hist_wplus_wdminus_eta_lep;
    Histo1DPtr _hist_wminus_wdplus_eta_lep;

    //d05-x01-
    Histo1DPtr _hist_wplus_wdstar_eta_lep;
    Histo1DPtr _hist_wminus_wdstar_eta_lep;


    // postprocessing histos
    //d05-x01
    Histo1DPtr _hist_w_inc;
    Scatter2DPtr _hist_wd_winc_ratio;
    Scatter2DPtr _hist_wdstar_winc_ratio;

    //d06-x01
    Histo1DPtr _hist_wplus_winc;
    Histo1DPtr _hist_wminus_winc;

    Scatter2DPtr _hist_wplusd_wplusinc_pt_ratio;
    Scatter2DPtr _hist_wminusd_wminusinc_pt_ratio;
    Scatter2DPtr _hist_wplusdstar_wplusinc_pt_ratio;
    Scatter2DPtr _hist_wminusdstar_wminusinc_pt_ratio;

    Histo1DPtr _hist_wplusd_wplusinc_pt ;
    Histo1DPtr _hist_wminusd_wminusinc_pt;
    Histo1DPtr _hist_wplusdstar_wplusinc_pt;
    Histo1DPtr _hist_wminusdstar_wminusinc_pt;

    // d07-x01
    Histo1DPtr _hist_wcjet_jets ;

    //d08-x01
    Scatter2DPtr  _hist_wcjet_jets_ratio ;
    Histo1DPtr    _hist_wcjet_jets_plus ;
    Histo1DPtr    _hist_wcjet_jets_minus;
    //@}

  };


  // The hook for the plugin system
  RIVET_DECLARE_PLUGIN(ATLAS_2014_I1282447);

}