1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
17 // AliMuonAccEffSubmitter : a class to help submit Acc x Eff simulations
18 // anchored to real runs for J/psi, upsilon, single muons, etc...
20 // This class is dealing with 3 different directories :
22 // - template directory ($ALICE_ROOT/PWG/muondep/AccEffTemplates) containing the
23 // basic template files to be used for a simuation. A template can contain
24 // some variables that will be replaced during during the copy from template
27 // - local directory, where the files from the template directory, are copied
28 // once the class has been configured properly (i.e. using the various Set, Use,
29 // etc... methods). Some other files (e.g. JDL ones) are generated from
30 // scratch and also copied into this directory.
31 // At this point one could(should) check the files, as they are the ones
32 // to be copied to the remote directory for the production
34 // - remote directory, the alien directory where the files will be copied
35 // (from the local directory) before the actual submission
37 // ==========================================================
41 // AliMuonAccEffSubmitter a;
42 // a.UseOCDBSnapshots(kFALSE);
43 // a.SetRemoteDir("/alice/cern.ch/user/l/laphecet/Analysis/LHC13d/simjpsi/pp503z0");
44 // a.ShouldOverwriteFiles(true);
45 // a.MakeNofEventsPropToTriggerCount("CMUL7-B-NOPF-MUON");
46 // a.SetVar("VAR_GENLIB_PARNAME","\"pp 5.03\"");
47 // a.SetRunList(195682);
49 // a.Run("test"); // will do everything but the submit
50 // a.Submit(false); // actual submission
53 // author: Laurent Aphecetche (Subatech
56 #include "AliMuonAccEffSubmitter.h"
58 #include "AliAnalysisTriggerScalers.h"
62 #include "TGridResult.h"
65 #include "TObjString.h"
77 ClassImp(AliMuonAccEffSubmitter)
79 //______________________________________________________________________________
80 AliMuonAccEffSubmitter::AliMuonAccEffSubmitter(const char* generator)
81 : AliMuonGridSubmitter(AliMuonGridSubmitter::kAccEff),
83 fFixedNofEvents(10000),
84 fMaxEventsPerChunk(5000),
86 fSplitMaxInputFileNumber(20),
89 fUseOCDBSnapshots(kFALSE),
91 fUseAODMerging(kFALSE)
95 SetOCDBPath("raw://");
97 SetLocalDirectory("Snapshot",LocalDir());
99 SetVar("VAR_OCDB_PATH","\"raw://\"");
101 SetVar("VAR_GENPARAM_GENLIB_TYPE","AliGenMUONlib::kJpsi");
102 SetVar("VAR_GENPARAM_GENLIB_PARNAME","\"pPb 5.03\"");
104 SetVar("VAR_GENCORRHF_QUARK","5");
105 SetVar("VAR_GENCORRHF_ENERGY","5");
107 // some default values for J/psi
108 SetVar("VAR_GENPARAMCUSTOM_PDGPARTICLECODE","443");
110 // default values below are from J/psi p+Pb (from muon_calo pass)
111 SetVar("VAR_GENPARAMCUSTOM_Y_P0","4.08E5");
112 SetVar("VAR_GENPARAMCUSTOM_Y_P1","7.1E4");
114 SetVar("VAR_GENPARAMCUSTOM_PT_P0","1.13E9");
115 SetVar("VAR_GENPARAMCUSTOM_PT_P1","18.05");
116 SetVar("VAR_GENPARAMCUSTOM_PT_P2","2.05");
117 SetVar("VAR_GENPARAMCUSTOM_PT_P3","3.34");
119 // some default values for single muons
120 SetVar("VAR_GENPARAMCUSTOMSINGLE_PTMIN","0.35");
122 SetVar("VAR_GENPARAMCUSTOMSINGLE_PT_P0","4.05962");
123 SetVar("VAR_GENPARAMCUSTOMSINGLE_PT_P1","1.0");
124 SetVar("VAR_GENPARAMCUSTOMSINGLE_PT_P2","2.46187");
125 SetVar("VAR_GENPARAMCUSTOMSINGLE_PT_P3","2.08644");
127 SetVar("VAR_GENPARAMCUSTOMSINGLE_Y_P0","0.729545");
128 SetVar("VAR_GENPARAMCUSTOMSINGLE_Y_P1","0.53837");
129 SetVar("VAR_GENPARAMCUSTOMSINGLE_Y_P2","0.141776");
130 SetVar("VAR_GENPARAMCUSTOMSINGLE_Y_P3","0.0130173");
132 UseOCDBSnapshots(fUseOCDBSnapshots);
134 SetGenerator(generator);
136 MakeNofEventsPropToTriggerCount();
138 AddToTemplateFileList("CheckESD.C");
139 AddToTemplateFileList("CheckAOD.C");
140 AddToTemplateFileList("AODtrain.C");
141 AddToTemplateFileList("validation.sh");
143 AddToTemplateFileList("Config.C");
144 AddToTemplateFileList("rec.C");
145 AddToTemplateFileList("sim.C");
146 AddToTemplateFileList("simrun.C");
147 AddToTemplateFileList(RunJDLName().Data());
149 UseExternalConfig(fExternalConfig);
152 //______________________________________________________________________________
153 AliMuonAccEffSubmitter::~AliMuonAccEffSubmitter()
158 ///______________________________________________________________________________
159 Bool_t AliMuonAccEffSubmitter::Generate(const char* jdlname) const
161 if ( TString(jdlname).Contains("merge",TString::kIgnoreCase) )
163 return GenerateMergeJDL(jdlname);
167 return GenerateRunJDL(jdlname);
171 ///______________________________________________________________________________
172 Bool_t AliMuonAccEffSubmitter::GenerateMergeJDL(const char* name) const
174 /// Create the JDL for merging jobs
175 /// FIXME: not checked !
179 std::ostream* os = CreateJDLFile(name);
186 Bool_t final = TString(name).Contains("merge",TString::kIgnoreCase);
188 (*os) << "# Generated merging jdl (production mode)" << std::endl
189 << "# $1 = run number" << std::endl
190 << "# $2 = merging stage" << std::endl
191 << "# Stage_<n>.xml made via: find <OutputDir> *Stage<n-1>/*root_archive.zip" << std::endl;
193 OutputToJDL(*os,"Packages",
194 GetMapValue("AliRoot"),
195 GetMapValue("Geant3"),
199 OutputToJDL(*os,"Executable","AOD_merge.sh");
201 OutputToJDL(*os,"Price","1");
205 OutputToJDL(*os,"Jobtag","comment: AliMuonAccEffSubmitter final merging");
209 OutputToJDL(*os,"Jobtag","comment: AliMuonAccEffSubmitter merging stage $2");
212 OutputToJDL(*os,"Workdirectorysize","5000MB");
214 OutputToJDL(*os,"Validationcommand",Form("%s/validation_merge.sh",RemoteDir().Data()));
216 OutputToJDL(*os,"TTL","7200");
218 OutputToJDL(*os,"OutputArchive",
219 "log_archive.zip:stderr,stdout@disk=1",
220 "root_archive.zip:AliAOD.root,AliAOD.Muons.root,AnalysisResults.root@disk=3"
223 OutputToJDL(*os,"Arguments",(final ? "2":"1")); // for AOD_merge.sh, 1 means intermediate merging stage, 2 means final merging
227 OutputToJDL(*os,"InputFile",Form("LF:%s/AODtrain.C",RemoteDir().Data()));
228 OutputToJDL(*os,"OutputDir",Form("%s/$1/Stage_$2/#alien_counter_03i#",RemoteDir().Data()));
229 OutputToJDL(*os,"InputDataCollection",Form("%s/$1/Stage_$2.xml,nodownload",RemoteDir().Data()));
230 OutputToJDL(*os,"split","se");
231 OutputToJDL(*os,"SplitMaxInputFileNumber",GetSplitMaxInputFileNumber());
232 OutputToJDL(*os,"InputDataListFormat","xml-single");
233 OutputToJDL(*os,"InputDataList","wn.xml");
237 OutputToJDL(*os,"InputFile",Form("LF:%s/AODtrain.C",RemoteDir().Data()),
238 Form("LF:%s/$1/wn.xml",RemoteDir().Data()));
239 OutputToJDL(*os,"OutputDir",Form("%s/$1",RemoteDir().Data()));
245 //______________________________________________________________________________
246 Bool_t AliMuonAccEffSubmitter::GenerateRunJDL(const char* name) const
248 /// Generate (locally) the JDL to perform the simulation+reco+aod filtering
249 /// (to be then copied to the grid and finally submitted)
253 std::ostream* os = CreateJDLFile(name);
260 OutputToJDL(*os,"Packages",
261 GetMapValue("AliRoot"),
262 GetMapValue("Geant3"),
266 OutputToJDL(*os,"Jobtag","comment: AliMuonAccEffSubmitter RUN $1");
268 OutputToJDL(*os,"split","production:1-$2");
270 OutputToJDL(*os,"Price","1");
272 OutputToJDL(*os,"OutputDir",Form("%s/$1/#alien_counter_03i#",RemoteDir().Data()));
274 OutputToJDL(*os,"Executable","/alice/bin/aliroot_new");
277 files.SetOwner(kTRUE);
278 TIter next(TemplateFileList());
281 while ( ( file = static_cast<TObjString*>(next())) )
283 if ( !file->String().Contains(".jdl",TString::kIgnoreCase) ||
284 !file->String().Contains("OCDB_") )
286 files.Add(new TObjString(Form("LF:%s/%s",RemoteDir().Data(),file->String().Data())));
290 if ( fUseOCDBSnapshots )
292 files.Add(new TObjString(Form("LF:%s/OCDB/$1/OCDB_sim.root",RemoteDir().Data())));
293 files.Add(new TObjString(Form("LF:%s/OCDB/$1/OCDB_rec.root",RemoteDir().Data())));
296 OutputToJDL(*os,"InputFile",files);
298 if ( CompactMode() == 0 )
301 OutputToJDL(*os,"OutputArchive", "log_archive.zip:stderr,stdout,aod.log,checkaod.log,checkesd.log,rec.log,sim.log@disk=1",
302 "root_archive.zip:galice*.root,Kinematics*.root,TrackRefs*.root,AliESDs.root,AliAOD.root,AliAOD.Muons.root,Merged.QA.Data.root,Run*.root@disk=2");
304 else if ( CompactMode() == 1 )
306 // keep only muon AODs and QA
307 OutputToJDL(*os,"OutputArchive", "log_archive.zip:stderr,stdout,aod.log,checkaod.log,checkesd.log,rec.log,sim.log@disk=1",
308 "root_archive.zip:galice*.root,AliAOD.Muons.root,Merged.QA.Data.root@disk=2");
312 AliError(Form("Unknown CompactMode %d",CompactMode()));
317 OutputToJDL(*os,"splitarguments","simrun.C --run $1 --chunk #alien_counter# --event $3");
319 OutputToJDL(*os,"Workdirectorysize","5000MB");
321 OutputToJDL(*os,"JDLVariables","Packages","OutputDir");
323 OutputToJDL(*os,"Validationcommand",Form("%s/validation.sh",RemoteDir().Data()));
325 OutputToJDL(*os,"TTL","72000");
330 //______________________________________________________________________________
331 Bool_t AliMuonAccEffSubmitter::MakeOCDBSnapshots()
333 /// Run sim.C and rec.C in a special mode to generate OCDB snapshots
334 /// Can only be done after the templates have been copied locally
336 if (!IsValid()) return kFALSE;
338 if (!fUseOCDBSnapshots) return kTRUE;
340 if (!NofRuns()) return kFALSE;
346 const std::vector<int>& runs = RunList();
348 for ( std::vector<int>::size_type i = 0; i < runs.size(); ++i )
350 Int_t runNumber = runs[i];
352 TString ocdbSim(Form("%s/OCDB/%d/OCDB_sim.root",SnapshotDir().Data(),runNumber));
353 TString ocdbRec(Form("%s/OCDB/%d/OCDB_rec.root",SnapshotDir().Data(),runNumber));
355 if ( !gSystem->AccessPathName(ocdbSim.Data()) &&
356 !gSystem->AccessPathName(ocdbRec.Data()) )
358 AliWarning(Form("Local OCDB snapshots already there for run %d. Will not redo them. If you want to force them, delete them by hand !",runNumber));
363 gSystem->Exec(Form("aliroot -b -q -x simrun.C --run %d --snapshot",runNumber));
365 if ( gSystem->AccessPathName(ocdbSim.Data()) )
367 AliError(Form("Could not create OCDB snapshot for simulation"));
371 if ( gSystem->AccessPathName(ocdbRec.Data()) )
373 AliError(Form("Could not create OCDB snapshot for reconstruction"));
378 LocalFileList()->Add(new TObjString(ocdbSim));
379 LocalFileList()->Add(new TObjString(ocdbRec));
385 //______________________________________________________________________________
386 Bool_t AliMuonAccEffSubmitter::Merge(Int_t stage, Bool_t dryRun)
388 /// Submit multiple merging jobs with the format "submit AOD_merge(_final).jdl run# (stage#)".
389 /// Also produce the xml collection before sending jobs
390 /// Initial AODs will be taken from fRemoteDir/[RUNNUMBER] while the merged
391 /// ones will be put into fMergedDir/AODs/[RUNNUMBER]
394 /// - inDir = "/alice/sim/2012/LHC12a10_bis" (where to find the data to merge)
395 /// = 0x0 --> inDir = homeDir/outDir/resDir
396 /// - outDir = "Sim/LHC11h/embedding/AODs" (where to store merged results)
397 /// - runList.txt must contains the list of run number
398 /// - stage=0 --> final merging / stage>0 --> intermediate merging i
401 if (!RemoteDirectoryExists(MergedDir().Data())) {
402 AliError(Form("directory %s does not exist", MergedDir().Data()));
406 gGrid->Cd(MergedDir().Data());
408 TString jdl = MergeJDLName(stage==0);
410 if (!RemoteFileExists(jdl.Data()))
412 AliError(Form("file %s does not exist in %s\n", jdl.Data(), RemoteDir().Data()));
416 const std::vector<int>& runs = RunList();
420 AliError("No run to work with");
426 gSystem->Exec("rm -f __failed__");
427 Bool_t failedRun = kFALSE;
429 for ( std::vector<int>::size_type i = 0; i < runs.size(); ++i )
432 AliInfo(Form("\n --- processing run %d ---\n", run));
434 TString runDir = Form("%s/%d", MergedDir().Data(), run);
436 if (!RemoteDirectoryExists(runDir.Data()))
438 AliInfo(Form(" - creating output directory %s\n", runDir.Data()));
439 gSystem->Exec(Form("alien_mkdir -p %s", runDir.Data()));
442 if (RemoteFileExists(Form("%s/root_archive.zip", runDir.Data())))
444 AliWarning(" ! final merging already done");
448 Int_t lastStage = GetLastStage(runDir.Data());
450 if (stage > 0 && stage != lastStage+1)
452 AliError(Form(" ! lastest merging stage = %d. Next must be stage %d or final stage\n", lastStage, lastStage+1));
456 TString wn = (stage > 0) ? Form("Stage_%d.xml", stage) : "wn.xml";
457 TString find = (lastStage == 0) ?
458 Form("alien_find -x %s %s/%d *root_archive.zip", wn.Data(), RemoteDir().Data(), run) :
459 Form("alien_find -x %s %s/%d/Stage_%d *root_archive.zip", wn.Data(), RemoteDir().Data(), run, lastStage);
460 gSystem->Exec(Form("%s 1> %s 2>/dev/null", find.Data(), wn.Data()));
461 gSystem->Exec(Form("grep -c /event %s > __nfiles__", wn.Data()));
462 ifstream f2("__nfiles__");
464 nFiles.ReadLine(f2,kTRUE);
466 gSystem->Exec("rm -f __nfiles__");
467 printf(" - number of files to merge = %d\n", nFiles.Atoi());
468 if (nFiles.Atoi() == 0) {
469 printf(" ! collection of files to merge is empty\n");
470 gSystem->Exec(Form("rm -f %s", wn.Data()));
472 } else if (stage > 0 && nFiles.Atoi() <= splitLevel && !reply.BeginsWith("y")) {
473 if (!reply.BeginsWith("n")) {
474 printf(" ! number of files to merge <= split level (%d). Continue? [Y/n] ", splitLevel);
476 reply.Gets(stdin,kTRUE);
479 if (reply.BeginsWith("n")) {
480 gSystem->Exec(Form("rm -f %s", wn.Data()));
487 TString dirwn = Form("%s/%s", runDir.Data(), wn.Data());
488 if (RemoteFileExists(dirwn.Data())) gGrid->Rm(dirwn.Data());
489 gSystem->Exec(Form("alien_cp file:%s alien://%s", wn.Data(), dirwn.Data()));
490 gSystem->Exec(Form("rm -f %s", wn.Data()));
494 if (stage > 0) query = Form("submit %s %d %d", jdl.Data(), run, stage);
495 else query = Form("submit %s %d", jdl.Data(), run);
496 printf(" - %s ...", query.Data());
505 Bool_t done = kFALSE;
506 TGridResult *res = gGrid->Command(query);
509 TString cjobId1 = res->GetKey(0,"jobId");
510 if (!cjobId1.IsDec())
518 AliInfo(Form(" DONE\n --> the job Id is: %s \n", cjobId1.Data()));
530 gSystem->Exec(Form("echo %d >> __failed__", run));
538 AliInfo("\n--------------------\n");
539 AliInfo("list of failed runs:\n");
540 gSystem->Exec("cat __failed__");
541 gSystem->Exec("rm -f __failed__");
548 //______________________________________________________________________________
549 void AliMuonAccEffSubmitter::Print(Option_t* opt) const
553 AliMuonGridSubmitter::Print(opt);
557 std::cout << Form("For each run, will generate %5.2f times the number of real events for trigger %s",
558 fRatio,ReferenceTrigger().Data()) << std::endl;
562 std::cout << Form("For each run, will generate %10d events",fFixedNofEvents) << std::endl;
565 std::cout << "MaxEventsPerChunk = " << fMaxEventsPerChunk << std::endl;
567 std::cout << "Will" << (fUseOCDBSnapshots ? "" : " NOT") << " use OCDB snaphosts" << std::endl;
570 //______________________________________________________________________________
571 Bool_t AliMuonAccEffSubmitter::Run(const char* mode)
573 /// mode can be one of (case insensitive)
575 /// LOCAL : copy the template files from the template directory to the local one
576 /// UPLOAD : copy the local files to the grid (requires LOCAL)
577 /// OCDB : make ocdb snapshots (requires LOCAL)
578 /// SUBMIT : submit the jobs (requires LOCAL + UPLOAD)
579 /// FULL : all of the above (requires all of the above)
581 /// TEST : as SUBMIT, but in dry mode (does not actually submit the jobs)
583 if (!IsValid()) return kFALSE;
588 if ( smode == "FULL")
590 return ( Run("LOCAL") && Run("OCDB") && Run("UPLOAD") && Run("SUBMIT") );
593 if ( smode == "LOCAL")
595 return CopyTemplateFilesToLocal();
598 if ( smode == "UPLOAD" )
600 return (CopyLocalFilesToRemote());
603 if ( smode == "OCDB" )
605 Bool_t ok = Run("LOCAL");
608 ok = MakeOCDBSnapshots();
613 if ( smode == "TEST" )
615 Bool_t ok = Run("LOCAL") && Run("OCDB") && Run("UPLOAD");
618 ok = (Submit(kTRUE)>0);
623 if ( smode == "FULL" )
625 Bool_t ok = Run("LOCAL") && Run("OCDB") && Run("UPLOAD");
628 ok = (Submit(kFALSE)>0);
633 if( smode == "SUBMIT" )
635 return (Submit(kFALSE)>0);
641 //______________________________________________________________________________
642 Bool_t AliMuonAccEffSubmitter::SetGenerator(const char* generator)
644 // set the variable to select the generator macro in Config.C
646 gSystem->Load("libEVGEN");
650 TString generatorFile(Form("%s/%s.C",TemplateDir().Data(),generator));
652 Int_t nofMissingVariables(0);
654 // first check we indeed have such a macro
655 if (!gSystem->AccessPathName(generatorFile.Data()))
657 TObjArray* variables = GetVariables(generatorFile.Data());
659 TIter next(variables);
662 while ( ( var = static_cast<TObjString*>(next())) )
664 if ( !Vars()->GetValue(var->String()) )
666 ++nofMissingVariables;
667 AliError(Form("file %s expect the variable %s to be defined, but we've not defined it !",generatorFile.Data(),var->String().Data()));
673 if ( !nofMissingVariables )
675 if (CheckCompilation(generatorFile.Data()))
678 SetVar("VAR_GENERATOR",Form("%s",generator));
679 AddToTemplateFileList(Form("%s.C",generator));
690 AliError(Form("Can not work with the macro %s",generatorFile.Data()));
695 //______________________________________________________________________________
696 void AliMuonAccEffSubmitter::SetOCDBPath(const char* ocdbPath)
698 /// Sets the OCDB path to be used
700 SetMapKeyValue("OCDBPath",ocdbPath);
704 //______________________________________________________________________________
705 void AliMuonAccEffSubmitter::SetOCDBSnapshotDir(const char* dir)
707 // change the directory used for snapshot
709 if (gSystem->AccessPathName(Form("%s/OCDB",dir)))
711 AliError(Form("Snapshot top directory (%s) should contain an OCDB subdir with runnumbers in there",dir));
715 SetMapKeyValue("OCDBSnapshot",dir);
719 //______________________________________________________________________________
720 void AliMuonAccEffSubmitter::MakeNofEventsPropToTriggerCount(const char* trigger, Float_t ratio)
722 SetMapKeyValue("ReferenceTrigger",trigger);
726 //______________________________________________________________________________
727 void AliMuonAccEffSubmitter::MakeNofEventsFixed(Int_t nevents)
729 fFixedNofEvents = nevents;
731 SetMapKeyValue("ReferenceTrigger","");
735 //______________________________________________________________________________
736 Int_t AliMuonAccEffSubmitter::Submit(Bool_t dryRun)
738 /// Submit multiple production jobs with the format "submit jdl 000run#.xml 000run#".
740 /// Return the number of submitted (master) jobs
743 /// - outputDir = "/alice/cern.ch/user/p/ppillot/Sim/LHC10h/JPsiPbPb276/AlignRawVtxRaw/ESDs"
744 /// - runList must contains the list of run number
745 /// - trigger is the (fully qualified) trigger name used to compute the base number of events
746 /// - mult is the factor to apply to the number of trigger to get the number of events to be generated
747 /// (# generated events = # triggers x mult
749 if (!IsValid()) return 0;
753 gGrid->Cd(RemoteDir());
755 if (!RemoteFileExists(RunJDLName()))
757 AliError(Form("file %s does not exist in %s", RunJDLName().Data(), RemoteDir().Data()));
763 AliError("No run list set. Use SetRunList");
766 const std::vector<int>& runs = RunList();
770 AliError("No run to work with");
774 // cout << "total number of selected MB events = " << totEvt << endl;
775 // cout << "required number of generated events = " << nGenEvents << endl;
776 // cout << "number of generated events per MB event = " << ratio << endl;
779 std::cout << "run\tchunks\tevents" << std::endl;
780 std::cout << "----------------------" << std::endl;
785 AliAnalysisTriggerScalers* ts(0x0);
787 for (std::vector<int>::size_type i=0; i < runs.size(); ++i)
789 Int_t runNumber = runs[i];
791 Int_t nEvtRun(fFixedNofEvents);
797 AliInfo(Form("Creating AliAnalysisTriggerScalers from OCDB=%s",OCDBPath().Data()));
798 ts = new AliAnalysisTriggerScalers(runs,OCDBPath().Data());
801 AliAnalysisTriggerScalerItem* trigger = ts->GetTriggerScaler(runNumber, "L2A", ReferenceTrigger().Data());
805 AliError(Form("Could not get trigger %s for run %09d",ReferenceTrigger().Data(),runNumber));
808 nEvtRun = TMath::Nint(fRatio * trigger->Value());
813 while (nEvtRun/nChunk+0.5 > MaxEventsPerChunk())
818 Int_t nEvtChunk = TMath::Nint(nEvtRun/nChunk + 0.5);
822 nEvts += nChunk*nEvtChunk;
824 std::cout << runNumber << "\t" << nChunk << "\t" << nEvtChunk << std::endl;
826 TString query(Form("submit %s %d %d %d", RunJDLName().Data(), runNumber, nChunk, nEvtChunk));
828 std::cout << query.Data() << " ..." << std::flush;
830 TGridResult* res = 0x0;
834 res = gGrid->Command(query);
839 TString cjobId1 = res->GetKey(0,"jobId");
841 if (!cjobId1.Length())
843 std::cout << " FAILED" << std::endl << std::endl;
849 std::cout << "DONE" << std::endl;
850 std::cout << Form(" --> the job Id is: %s",cjobId1.Data()) << std::endl << std::endl;
855 std::cout << " FAILED" << std::endl << std::endl;
861 std::cout << std::endl
862 << "total number of jobs = " << nJobs << std::endl
863 << "total number of generated events = " << nEvts << std::endl
871 //______________________________________________________________________________
872 void AliMuonAccEffSubmitter::UpdateLocalFileList(Bool_t clearSnapshots)
874 /// Update the list of local files
876 AliMuonGridSubmitter::UpdateLocalFileList();
878 if (!NofRuns()) return;
880 if ( clearSnapshots )
882 TIter next(LocalFileList());
885 while ( ( file = static_cast<TObjString*>(next())) )
887 if ( file->String().Contains("OCDB_") )
889 LocalFileList()->Remove(file);
892 LocalFileList()->Compress();
895 const char* type[] = { "sim","rec" };
897 const std::vector<int>& runs = RunList();
899 for ( std::vector<int>::size_type i = 0; i < runs.size(); ++i )
901 Int_t runNumber = runs[i];
903 for ( Int_t t = 0; t < 2; ++t )
905 TString snapshot(Form("%s/OCDB/%d/OCDB_%s.root",SnapshotDir().Data(),runNumber,type[t]));
907 if ( !gSystem->AccessPathName(snapshot.Data()) )
909 if ( !LocalFileList()->FindObject(snapshot.Data()) )
911 LocalFileList()->Add(new TObjString(snapshot));
918 //______________________________________________________________________________
919 void AliMuonAccEffSubmitter::UseOCDBSnapshots(Bool_t flag)
921 /// Whether or not to use OCDB snapshots
922 /// Using OCDB snapshots will speed-up both the sim and reco initialization
923 /// phases on each worker node, but takes time to produce...
924 /// So using them is not always a win-win...
926 fUseOCDBSnapshots = flag;
929 SetVar("VAR_OCDB_SNAPSHOT","kTRUE");
933 SetVar("VAR_OCDB_SNAPSHOT","kFALSE");
936 UpdateLocalFileList();
939 //______________________________________________________________________________
940 void AliMuonAccEffSubmitter::UseAODMerging(Bool_t flag)
942 /// whether or not we should generate JDL for merging AODs
944 fUseAODMerging = flag;
946 AddToTemplateFileList(MergeJDLName(kFALSE).Data());
947 AddToTemplateFileList(MergeJDLName(kTRUE).Data());
948 AddToTemplateFileList("AOD_merge.sh");
949 AddToTemplateFileList("validation_merge.sh");
952 //______________________________________________________________________________
953 void AliMuonAccEffSubmitter::UseExternalConfig(const char* externalConfigFullFilePath)
955 // use an external config (or the default Config.C if externalConfigFullFilePath="")
957 fExternalConfig = externalConfigFullFilePath;
958 if ( fExternalConfig.Length() > 0 )
960 AddToTemplateFileList(fExternalConfig);
964 AddToTemplateFileList("Config.C");