GoodRunsListSelectionToolというCP toolを導入してみましょう。 GoodRunsListSelectionToolはxAODEventInfoを入力にします。 GoodRunListというのは物理解析にとって良いeventであるかどうかがxmlで記述されたものです。 つまりATLASで検出器の状態が良くなかったと判断されたevent (正確にはlumi block単位) をこのツールで除いていくことが出来ます。
MyAnalysisAlg.hでGoodRunsListSelectionTool パッケージをインクルードします。
// GRL
#include "GoodRunsLists/GoodRunsListSelectionTool.h"
新しいメンバー変数を追加します。
GoodRunsListSelectionTool *m_grl; //!
この//!は呪文だと思ってつけておいて下さい。
MyAnalysisAlg.cxxに移動し、一番上に:
#include <TSystem.h>
initialize()に以下を追加します。
// GRL
m_grl = new GoodRunsListSelectionTool("GoodRunsListSelectionTool");
const char* GRLFilePath =
"$ROOTCOREBIN/data/MyNewPackage/data16_13TeV.periodAllYear_DetStatus-v84-pro20-16_DQDefects-00-02-04_PHYS_StandardGRL_All_Good_25ns.xml";
const char* fullGRLFilePath = gSystem->ExpandPathName (GRLFilePath);
std::vector<std::string> vecStringGRL;
vecStringGRL.push_back(fullGRLFilePath);
ANA_CHECK(m_grl->setProperty( "GoodRunsListVec", vecStringGRL));
ANA_CHECK(m_grl->setProperty("PassThrough", false));
ANA_CHECK(m_grl->initialize());
data16_13TeV.periodAllYear_DetStatus-v84-pro20-16_DQDefects-00-02-04_PHYS_StandardGRL_All_Good_25ns.xmlというのがGoodRunsListです。 GoodRunsListの情報は以下で確認できます:
https://twiki.cern.ch/twiki/bin/view/AtlasProtected/GoodRunListsForAnalysisRun2
login.icepp.jpの以下のディレクトリにコピーしておいたので、shareというディレクトリを作ってそこにコピーします。
$ mkdir MyNewPackage/share
$ cp /gpfs/fs2001/nobe/data2086b/sample_forTutorial2016/GRL/data16_13TeV.periodAllYear_DetStatus-v84-pro20-16_DQDefects-00-02-04_PHYS_StandardGRL_All_Good_25ns.xml MyNewPackage/share/
share以下にファイルを置いておけば、コンパイルしたときに$ROOTCOREBIN/data/MyNewPackage/にリンクが張られます。
さて、引き続きMyxAODAnalysis.cxxを編集していきます。execute()のEventInfoを取得した直後に以下を追加します:
bool isMC = false;
if(eventInfo->eventType( xAOD::EventInfo::IS_SIMULATION)) {
isMC = true;
}
if(!isMC) {
if(!m_grl->passRunLB(*eventInfo)) {
Info("execute()", "Run = %lu : Event = %lu is discarded by GRL.", runNumber, eventNumber);
return EL::StatusCode::SUCCESS;
}
}
finalize()に以下を追加します:
if (m_grl) {
delete m_grl;
m_grl = 0;
}
最後にRootCoreにGoodRunsListを追加します:
PACKAGE_DEP = EventLoop xAODRootAccess AsgTools xAODEventInfo xAODMuon xAODEgamma GoodRunsLists
コンパイルして実行してみて下さい。
先程練習問題で触ってみたelectronに専用のCPツールを使ってエネルギー補正をかけてみましょう。
MyAnalysisAlg.hの先頭に以下を追加してツールをincludeします。
#include <ElectronPhotonFourMomentumCorrection/EgammaCalibrationAndSmearingTool.h>
新しいメンバ関数を追加します。
CP::IEgammaCalibrationAndSmearingTool *m_electronCalibrationAndSmearingTool; //!
MyAnalysisAlg.cxxに移動して、先頭に以下を追加します。
#include "PATInterfaces/CorrectionCode.h"
#include "xAODCore/ShallowAuxContainer.h"
#include "xAODCore/ShallowCopy.h"
次にinitialize()に、
m_electronCalibrationAndSmearingTool = new CP::EgammaCalibrationAndSmearingTool("myEgammaCalibrationAndSmearingTool");
ANA_CHECK(asg::setProperty(m_electronCalibrationAndSmearingTool, "ESModel", "es2016PRE"));
ANA_CHECK(asg::setProperty(m_electronCalibrationAndSmearingTool, "decorrelationModel", "1NPCOR_PLUS_UNCOR"));
ANA_CHECK(m_electronCalibrationAndSmearingTool->initialize());
finalize()に
if(m_electronCalibrationAndSmearingTool){
delete m_electronCalibrationAndSmearingTool;
m_electronCalibrationAndSmearingTool = 0;
}
を追加します。そして、execute()の中でさっき自分で書いたelectron loopの下に以下を追加です
std::pair< xAOD::ElectronContainer*, xAOD::ShallowAuxContainer* > electrons_shallowCopy = xAOD::shallowCopyContainer( *electrons );
xAOD::ElectronContainer::iterator electronSC_itr = (electrons_shallowCopy.first)->begin();
xAOD::ElectronContainer::iterator electronSC_end = (electrons_shallowCopy.first)->end();
for( ; electronSC_itr != electronSC_end; ++electronSC_itr ) {
if(m_electronCalibrationAndSmearingTool->applyCorrection(**electronSC_itr) == CP::CorrectionCode::Error){
Error("execute()", "ElectronCalibrationAndSmearingTool returns Error CorrectionCode");
}
Info("execute()", " corrected electron pt = %.2f GeV", ((*electronSC_itr)->pt() * 0.001));
}
delete electrons_shallowCopy.first;
delete electrons_shallowCopy.second;
元々のコンテナの値を変更することはできないので、shallow copyと呼ばれるコピーを作成して、 そこにツールをあてています。
最後にMakefile.RootCoreに使用するパッケージを追加します。
PACKAGE_DEP = EventLoop xAODRootAccess AsgTools xAODEventInfo [............] xAODCore ElectronPhotonFourMomentumCorrection
コンパイルして実行してください。電子のpTが補正されたのが確認できましたか?
このままではジェットを電子と誤認した事象をたくさん含んでいます。カロリメータのシャワー形状や、内部飛跡検出器の情報を組み合わせて本物の電子を選別していきます。たくさんの変数を組み合わせて定義したLikelihood変数を用いてカットします。専用のCPツールが用意されています。信号に対するefficiencyに応じてloose, medium, tightといったworking pointが定義されています。ここではmedium working pointを要求してみましょう。
MyAnalysisAlg.hで
#include "ElectronPhotonSelectorTools/AsgElectronLikelihoodTool.h"
メンバ関数に
AsgElectronLikelihoodTool* m_MediumLH; //!
MyAnalysisAlg.cxxのinitialize()で
m_MediumLH = new AsgElectronLikelihoodTool ("MediumLH");
// select the working point
ANA_CHECK(m_MediumLH->setProperty("WorkingPoint", "MediumLHElectron"));
ANA_CHECK(m_MediumLH->initialize());
finalize()で
if(m_MediumLH){
delete m_MediumLH;
m_MediumLH = 0;
}
先程追加したshallowCopyのループに以下の太字の2行を追加します。
std::pair< xAOD::ElectronContainer*, xAOD::ShallowAuxContainer* > electrons_shallowCopy = xAOD::shallowCopyContainer( *electrons );
xAOD::ElectronContainer::iterator electronSC_itr = (electrons_shallowCopy.first)->begin();
xAOD::ElectronContainer::iterator electronSC_end = (electrons_shallowCopy.first)->end();
for( ; electronSC_itr != electronSC_end; ++electronSC_itr ) {
if(m_electronCalibrationAndSmearingTool->applyCorrection(**electronSC_itr) == CP::CorrectionCode::Error){
Error("execute()", "ElectronCalibrationAndSmearingTool returns Error CorrectionCode");
}
if ( !LHMediumSel ) continue;
Info("execute()", " corrected electron pt = %.2f GeV", ((*electronSC_itr)->pt() * 0.001));
}
delete electrons_shallowCopy.first;
delete electrons_shallowCopy.second;
Makefile.RootCoreにElectronPhotonSelectorToolsを追加して、コンパイルして、実行です。
バックグラウンドが落ち、さっきと比べて電子の数が大分減ったと思います。