00001
00022 #include "TFile.h"
00023 #include "StIOMaker/StIOMaker.h"
00024 #include "StMessMgr.h"
00025 #include "StPrepEmbedMaker.h"
00026 #include "StEvtHddr.h"
00027 #include "TTree.h"
00028 #include "StGenericVertexMaker/StGenericVertexMaker.h"
00029 #include "StGenericVertexMaker/StFixedVertexFinder.h"
00030
00031 #include "tables/St_vertexSeed_Table.h"
00032 #include "TString.h"
00033 #include "TSystem.h"
00034
00035 #include <unistd.h>
00036
00037 ClassImp(StPrepEmbedMaker)
00038 struct embedSettings{
00039 Double_t mult;
00040 Int_t pid;
00041 Double_t ptlow;
00042 Double_t pthigh;
00043 Double_t etalow;
00044 Double_t etahigh;
00045 Double_t philow;
00046 Double_t phihigh;
00047 Double_t temperature;
00048 Int_t rnd1;
00049 Int_t rnd2;
00050 Double_t vzlow;
00051 Double_t vzhigh;
00052 Double_t vr;
00053 Double_t vpdvz;
00054 Int_t NReqTrg;
00055 static const Int_t nTriggerId = 32 ;
00056 Int_t ReqTrgId[nTriggerId];
00057 TString mode;
00058 };
00059
00060 static embedSettings *mSettings = 0;
00061 const Double_t StPrepEmbedMaker::mRapidityMaximumCut = 10.0 ;
00062
00063
00064 StPrepEmbedMaker::StPrepEmbedMaker(const Char_t *name) : StMaker(name)
00065 {
00066 mGeant3=0;
00067 mTagFile = "" ;
00068 mMoreTagsFile = "" ;
00069 mFzFile = "temp.fz";
00070 mEventCounter = 0;
00071
00072 if( !mSettings ){
00073 mSettings = new embedSettings();
00074
00076 mSettings->vzlow = -200.0 ;
00077 mSettings->vzhigh = 200.0 ;
00078
00080 mSettings->mult=0.;
00081
00083 mSettings->temperature=0.3;
00084
00086 mSettings->NReqTrg = 0 ;
00087 for(Int_t itrg=0; itrg<mSettings->nTriggerId; itrg++){
00088 mSettings->ReqTrgId[itrg] = 0 ;
00089 }
00090
00092 mSettings->mode = "flatpt";
00093 }
00094 mFile = 0;
00095 mMoreFile = 0 ;
00096 mTree = 0;
00097 mSkipMode = kFALSE;
00098 mSpreadMode = kFALSE;
00099 mOpenFzFile = kFALSE;
00100 mPrimeMode = kFALSE;
00101 mPrimed = kFALSE;
00102 mVpdVzCutMode = kFALSE;
00103 mRapidityMode = kTRUE;
00104 }
00105
00106 StPrepEmbedMaker::~StPrepEmbedMaker() {
00107 SafeDelete(mFile);
00108 }
00109
00110
00111 Int_t StPrepEmbedMaker::Init()
00112 {
00113 srand((unsigned)time(0));
00114 mSettings->rnd1 = abs(int(rand()*10000)+getpid());
00115 mSettings->rnd2 = abs(int(rand()*10000)+getpid());
00116
00117 if (mSettings->mode.CompareTo("strange", TString::kIgnoreCase) == 0)
00118 {
00119 mSpreadMode= kTRUE;
00120 LOG_INFO <<"StPrepEmbedMaker::Init Setting spreader mode for embedding mode "<<mSettings->mode <<endm;
00121 }
00122
00123 return StMaker::Init();
00124 }
00125
00126
00127 Int_t StPrepEmbedMaker::InitRun(const int runnum)
00128 {
00129
00130 if (! mGeant3) {
00131 mGeant3 = TGiant3::Geant3();
00132 if( ! mGeant3)
00133 {
00134 LOG_ERROR << "StPrepEmbedMaker::InitRun Geant3 pointer not found. exiting."<<endm;
00135 return kStErr;
00136 }
00137 }
00138
00139
00140 if (mTagFile.IsWhitespace()){
00141 LOG_ERROR << "StPrepEmbedMaker::InitRun TagFile has not been defined" << endm;
00142 return kStErr;
00143 }
00144
00145
00146 mFile = TFile::Open(mTagFile);
00147 if (! mFile ) {
00148 LOG_ERROR << "StPrepEmbedMaker::Init TagFile : " << mTagFile << " cannot be opened" << endm;
00149 return kStErr;
00150 }
00151
00152
00153 mTree = (TTree *) mFile->Get("Tag");
00154 if (! mTree ) {
00155 LOG_ERROR << "StPrepEmbedMaker::Init In TagFile : " << mTagFile << " cannot find TTree \"Tag\"" << endm;
00156 return kStErr;
00157 }
00158
00159
00160 if (mSpreadMode){
00161 LOG_INFO << "StPrepEmbedMaker::Init Spreader mode set. Looking for MoreTags file ..."<<endm;
00162 mMoreTagsFile = mTagFile;
00163 const int indx1 = mMoreTagsFile.Index(".tags",0);
00164 const int indx2 = mMoreTagsFile.Last('.');
00165 if (indx1!=indx2) mMoreTagsFile.Remove(indx1+1,(indx2-indx1));
00166 mMoreTagsFile.Insert(indx1+1,"moretags.");
00167 mMoreFile = TFile::Open(mMoreTagsFile);
00168
00169 if (mMoreFile ) {
00170 mMoreTree = (TTree *) mMoreFile->Get("MoreTags");
00171 if (! mMoreTree ) {
00172 LOG_ERROR << "StPrepEmbedMaker::Init In MoreTagsFile : " << mMoreTagsFile << " cannot find TTree \"MoreTags\"" << endm;
00173 return kStErr;
00174 }
00175 }
00176 else {
00177 LOG_INFO << "StPrepEmbedMaker::Init File moretags.root not found. If this embedding is for years 2007 through 2009, this will most likely cause the embedding to quit prematurely."<<endm;
00178 return kStErr ;
00179 }
00180 }
00181
00182
00183
00184 #if 0
00185 if(mSettings->mode.CompareTo("Spectrum", TString::kIgnoreCase) == 0)
00186 {
00187
00188
00189
00190
00191
00192 TString cmd;
00193 cmd = Form("root4star -q \'~starofl/embedding/getVerticiesFromTags.C\(1000,\"./\",\"%s\")\'",
00194 mTagFile.Data());
00195 cmd = Form("~starofl/embedding/GENTX/gentx %s %s %i %i %f %f %f %i %f %f %f %i",
00196 mTagFile.Data(),"temp.fz",
00197 1000,mSettings->mult,mSettings->etalow, mSettings->etahigh,
00198 0.0,0,mSettings->ptlow, mSettings->pthigh,mSettings->temperature,
00199 mSettings->rnd1);
00200 gSystem->Exec(cmd.Data());
00201
00202 Do("gfile p temp.fz");
00203 }
00204 #endif
00205
00206 if( mOpenFzFile ) {
00207
00208
00209
00210 if (StIOMaker* ioMaker = (StIOMaker*)GetMakerInheritsFrom("StIOMaker")) {
00211 if (const Char_t* daqfile = ioMaker->GetFile()) {
00212 LOG_DEBUG << "StPrepEmbedMaker::InitRun daq file: " << daqfile << endm;
00213 mFzFile = gSystem->BaseName(daqfile);
00214 mFzFile.ReplaceAll(".daq",".fz");
00215 }
00216 }
00217 LOG_INFO << "StPrepEmbedMaker::InitRun Open FZ file: " << mFzFile << endm;
00218 Do("user/output o " + mFzFile);
00219 }
00220
00221
00222 Do("make gstar");
00223 gSystem->Load("libgstar");
00224
00225 Do("detp hadr_on");
00226 TString cmd("rndm ");
00227 cmd+=mSettings->rnd1; cmd+=" "; cmd+=mSettings->rnd2;
00228 Do(cmd.Data());
00229
00230 return 0;
00231 }
00232
00233
00234 Int_t StPrepEmbedMaker::Make()
00235 {
00236 mEventCounter++;
00237 StEvtHddr* EvtHddr = (StEvtHddr*) GetDataSet("EvtHddr");
00238 if (! EvtHddr) {
00239 LOG_ERROR << "StPrepEmbedMaker::Make EvtHddr has not been found" << endm;
00240 return kStErr;
00241 }
00242 Int_t nFound = mTree->Draw("uncorrectedNumberOfPrimaries:primaryVertexFlag",
00243 Form("mRunNumber==%i&&mEventNumber==%i",EvtHddr->GetRunNumber(),EvtHddr->GetEventNumber()),
00244 "goff");
00245 if (nFound != 1) {
00246 LOG_ERROR << "StPrepEmbedMaker::Make Run/Event = " << EvtHddr->GetRunNumber() << "/" << EvtHddr->GetEventNumber()
00247 << " has been found in tag file" << nFound << " times" << endm;
00248 return kStErr;
00249 }
00250 LOG_INFO << "StPrepEmbedMaker::Make Run/Event = " << EvtHddr->GetRunNumber()
00251 << "/" << EvtHddr->GetEventNumber()
00252 << " has been found with uncorrectedNumberOfPrimaries = " << mTree->GetV1()[0]
00253 << " and primaryVertexFlag = " << mTree->GetV2()[0] << endm;
00254
00255 if (mTree->GetV1()[0] <= 0 || mTree->GetV2()[0] )
00256 {
00257 LOG_ERROR << "StPrepEmbedMaker::Make reject this event" << endm;
00258 return kStErr;
00259 }
00260
00261
00262 const Int_t numberOfPrimaryTracks = (Int_t) mTree->GetV1()[0];
00263 const Int_t npart = getMultiplicity( *EvtHddr, numberOfPrimaryTracks ) ;
00264
00265 nFound = (Int_t) mTree->Draw("primaryVertexX:primaryVertexY:primaryVertexZ:TriggerId",
00266 Form("mRunNumber==%i&&mEventNumber==%i",
00267 EvtHddr->GetRunNumber(),
00268 EvtHddr->GetEventNumber()),
00269 "goff");
00270 const Double_t xyz[3] = {mTree->GetV1()[0],mTree->GetV2()[0],mTree->GetV3()[0]};
00271
00272
00273 if (fabs(xyz[0])<1e-7 && fabs(xyz[1])<1e-7 && fabs(xyz[2])<1e-7 ){
00274 LOG_INFO << "StPrepEmbedMaker::Make Event " << EvtHddr->GetEventNumber()
00275 << " has tags with vertex approx at (0,0,0) - probably no PV, skipping."
00276 << endm;
00277 return kStSKIP;
00278 }
00279
00280
00281
00282
00283 if (mSkipMode == kTRUE){
00284 const Double_t vr = sqrt(xyz[0]*xyz[0] + xyz[1]*xyz[1]);
00285 if (xyz[2]<mSettings->vzlow || xyz[2]>mSettings->vzhigh || vr>=mSettings->vr ){
00286 LOG_INFO << "StPrepEmbedMaker::Make Event " << EvtHddr->GetEventNumber()
00287 << " has tags with vertex at (" << xyz[0] << "," << xyz[1] << "," << xyz[2]
00288 << "), vr = " << vr
00289 << " - out of Vz or Vr range, skipping." << endm;
00290 return kStSKIP;
00291 }
00292
00293 LOG_INFO << "StPrepEmbedMaker::Make Event " << EvtHddr->GetEventNumber()
00294 << " has tags with vertex at (" << xyz[0] << "," << xyz[1] << "," << xyz[2]
00295 << "), vr = " << vr
00296 << " - within requested Vz and Vr range !" << endm;
00297 }
00298
00299
00300 if (mSkipMode == kTRUE){
00301 LOG_INFO << "StPrepEmbedMaker::Event " << EvtHddr->GetEventNumber()
00302 << " has Triggers: " << endm;
00303 for (Int_t iTrg=0 ; iTrg<mSettings->nTriggerId ; iTrg++){
00304 LOG_INFO << mTree->GetV4()[iTrg] << " ";
00305 }
00306 LOG_INFO << endm;
00307
00308 Bool_t fired = kFALSE;
00309 for (Int_t iTrg=0 ; iTrg<mSettings->nTriggerId ; iTrg++){
00310 for (Int_t iReqTrg=0; iReqTrg<mSettings->NReqTrg ; iReqTrg++) {
00311 if (mTree->GetV4()[iTrg] == mSettings->ReqTrgId[iReqTrg]){
00312 LOG_INFO << "StPrepEmbedMaker::Requested trigger " << mSettings->ReqTrgId[iReqTrg] << " is fired!" << endm;
00313 fired = kTRUE;
00314 }
00315 }
00316 }
00317
00318 if (!fired && mSettings->NReqTrg>0) {
00319 LOG_INFO << "StPrepEmbedMaker::No requested triggers are fired in this event, skipping." << endm;
00320 return kStSKIP;
00321 }
00322 }
00323
00324
00325 Float_t vpdvz;
00326 if(mSkipMode == kTRUE && mVpdVzCutMode == kTRUE)
00327 {
00328
00329 nFound=0;
00330 nFound = (Int_t) mTree->Draw("VpdVz",
00331 Form("mRunNumber==%i&&mEventNumber==%i",
00332 EvtHddr->GetRunNumber(),
00333 EvtHddr->GetEventNumber()),"goff");
00334
00335
00336 if(nFound == -1 && mMoreTree) {
00337 nFound = (Int_t) mMoreTree->Draw("VpdVz",
00338 Form("RunId==%i&&EvtId==%i",
00339 EvtHddr->GetRunNumber(),
00340 EvtHddr->GetEventNumber()),
00341 "goff");
00342
00343 LOG_INFO << "StPrepEmbedMaker::Make Use moretags file to extract VpdVz, nFound =" << nFound << endm ;
00344 }
00345
00346 if (nFound != 1) {
00347 LOG_ERROR << "StPrepEmbedMaker::Make Run/Event = " << EvtHddr->GetRunNumber() << "/" << EvtHddr->GetEventNumber()
00348 << " has been found in moretags file " << nFound << " times" << endm;
00349 return kStErr;
00350 }
00351 vpdvz = mMoreTree->GetV1()[0];
00352 LOG_INFO << vpdvz << endm;
00353
00354
00355 if( fabs(vpdvz) < 1e-7 ) {
00356 LOG_INFO << "StPrepEmbedMaker::Make Event " << EvtHddr->GetEventNumber()
00357 << " has tags with vertex at (" << xyz[0] << "," << xyz[1] << "," << xyz[2]
00358 << "), VpdVz = " << vpdvz
00359 << " - VpdVz is too small (i.e. no BTOF in this run), skipping." << endm;
00360 return kStSKIP;
00361 }
00362 if( fabs(xyz[2]-vpdvz) > mSettings->vpdvz ) {
00363 LOG_INFO << "StPrepEmbedMaker::Make Event " << EvtHddr->GetEventNumber()
00364 << " has tags with vertex at (" << xyz[0] << "," << xyz[1] << "," << xyz[2]
00365 << "), VpdVz = " << vpdvz
00366 << " - out of |Vz-VpdVz| range, skipping." << endm;
00367 return kStSKIP;
00368 }
00369 }
00370
00371
00372
00373
00374
00375
00376
00377 Double_t xyzerr[3] = {0.,0.,0.};
00378
00379
00380
00381 if(mSettings->mode.CompareTo("strange", TString::kIgnoreCase) == 0)
00382 {
00383
00384
00385
00386 nFound=0;
00387 nFound = (Int_t) mTree->Draw("sigmaPVX:sigmaPVY:sigmaPVZ",
00388 Form("mRunNumber==%i&&mEventNumber==%i",
00389 EvtHddr->GetRunNumber(),
00390 EvtHddr->GetEventNumber()),"goff");
00391
00392
00393 if(nFound == -1 && mMoreTree) {
00394 nFound = (Int_t) mMoreTree->Draw("VXERR:VYERR:VZERR",
00395 Form("RunId==%i&&EvtId==%i",
00396 EvtHddr->GetRunNumber(),
00397 EvtHddr->GetEventNumber()),
00398 "goff");
00399
00400 LOG_INFO << "StPrepEmbedMaker::Make Use moretags file to extract vertex errors, nFound =" << nFound << endm ;
00401 }
00402
00403 if (nFound != 1) {
00404 LOG_ERROR << "StPrepEmbedMaker::Make Run/Event = " << EvtHddr->GetRunNumber() << "/" << EvtHddr->GetEventNumber()
00405 << " has been found in moretags file " << nFound << " times" << endm;
00406 return kStErr;
00407 }
00408 xyzerr[0] = mMoreTree->GetV1()[0];
00409 xyzerr[1] = mMoreTree->GetV2()[0];
00410 xyzerr[2] = mMoreTree->GetV3()[0];
00411 LOG_INFO << xyzerr[0] << " " << xyzerr[1] << " " << xyzerr[2] << endm;
00412
00413
00414
00415
00416 StGenericVertexMaker * vmaker = (StGenericVertexMaker*) GetMaker("GenericVertex");
00417 StFixedVertexFinder * vfinder = (StFixedVertexFinder *) vmaker->GetGenericFinder();
00418 vfinder->SetVertexPosition(xyz[0],xyz[1],xyz[2]);
00419 }
00420
00421 if( mPrimeMode && !mPrimed ) {
00422 mSavePid = mSettings->pid;
00423 mSettings->pid = 45;
00424 }
00425
00426
00427 gkine(npart, xyz[2], xyz[2]);
00428
00429
00430 if( mRapidityMode ) phasespace(npart);
00431
00432 Do(Form("gvertex %f %f %f",xyz[0],xyz[1],xyz[2]));
00433 if( mSettings->mode.CompareTo("strange", TString::kIgnoreCase) == 0 )
00434 {
00435 Do(Form("gspread %f %f %f", xyzerr[0],xyzerr[1],xyzerr[2]));
00436 }
00437 else
00438 {
00439 Do("vsig 0 0;");
00440 }
00441
00442
00443 if( mSettings->mode.CompareTo("Spectrum", TString::kIgnoreCase) == 0 ){
00444
00445 if ( mSettings->temperature > 0.0 ){
00446
00447
00448
00449
00450
00451 LOG_INFO << "StPrepEmbedMaker::Make Generate sloped momentum distribution with T="
00452 << mSettings->temperature << " GeV !"
00453 << endm;
00454
00455 Do("subevent 0");
00456 Do(Form("detp mick miky.np=%d code=%d mult=%d slope=%f dy=%f",
00457 1, mSettings->pid, -npart, 1.0/mSettings->temperature, -mSettings->etahigh));
00458
00459 Do(Form("user/input please my.micky"));
00460
00461 }
00462 else{
00463 LOG_ERROR << "StPrepEmbedMaker::Make Input temperature <= 0, T=" << mSettings->temperature
00464 << ", skip to generate sloped momentum distribution"
00465 << endm ;
00466 }
00467 }
00468
00469 Do("trig 1");
00470
00471 if( mPrimeMode && !mPrimed ){
00472 mSettings->pid = mSavePid;
00473 mPrimed = kTRUE;
00474 }
00475
00476 return kStOK;
00477 }
00478
00479
00480 Int_t StPrepEmbedMaker::Finish()
00481 {
00482 if( mOpenFzFile ) {
00484 LOG_INFO << "StPrepEmbedMaker::Finish Write and close fz file: " << mFzFile << endm;
00485 Do("user/output c " + mFzFile);
00486 }
00487 return 0;
00488 }
00489
00490
00491 void StPrepEmbedMaker::Do(const Char_t *job)
00492 {
00493 Int_t l=strlen(job);
00494 if (l) {
00495 LOG_INFO << "StPrepEmbedMaker::Do(" << job << ");" << endm;
00496 mGeant3->Kuexel(job);
00497 }
00498 }
00499
00500
00501 void StPrepEmbedMaker::SetPartOpt(const Int_t pid, const Double_t mult)
00502 {
00503 mSettings->mult=mult; mSettings->pid=pid;
00504 LOG_INFO << "StPrepEmbedMaker::SetPartOpt mult = " << mSettings->mult
00505 << " pid = " << mSettings->pid << endm;
00506 }
00507
00508
00509 void StPrepEmbedMaker::SetOpt(const Double_t ptlow, const Double_t pthigh,
00510 const Double_t etalow, const Double_t etahigh, const Double_t philow,
00511 const Double_t phihigh, const TString type)
00512 {
00513 mSettings->ptlow=ptlow; mSettings->pthigh=pthigh;
00514 mSettings->etalow=etalow; mSettings->etahigh=etahigh;
00515 mSettings->philow=philow; mSettings->phihigh=phihigh;
00516 mSettings->mode=type;
00517 LOG_INFO << "StPrepEmbedMaker::SetOpt ptlow = " << mSettings->ptlow << " pthigh = " << mSettings->pthigh
00518 << " etalow = " << mSettings->etalow << " etahigh = " << mSettings->etahigh
00519 << " philow = " << mSettings->philow << " phihigh = " << mSettings->phihigh
00520 <<" Mode: "<< type.Data() << endm;
00521 }
00522
00523
00524 void StPrepEmbedMaker::SetTemp(const double t)
00525 {
00526 mSettings->temperature=t;
00527 LOG_INFO << "StPrepEmbedMaker::SetTemp set temperature= " << mSettings->temperature << endm;
00528 }
00529
00530
00531 void StPrepEmbedMaker::SetTagFile(const Char_t *file)
00532 {
00533 mTagFile = file;
00534 LOG_INFO << "StPrepEmbedMaker::SetTagFile set tags file= " << mTagFile << endm;
00535 }
00536
00537
00538 void StPrepEmbedMaker::SetSkipMode(const Bool_t flag)
00539 {
00540 LOG_INFO << "StPrepEmbedMaker::SetSkipMode set skip mode= ";
00541 mSkipMode = flag;
00542
00543 if( mSkipMode ){
00544 LOG_INFO << " ON" << endm ;
00545 }
00546 else{
00547 LOG_INFO << " OFF" << endm ;
00548 }
00549 }
00550
00551
00552 void StPrepEmbedMaker::SetSpreadMode(const Bool_t flag)
00553 {
00554 mSpreadMode=flag;
00555
00556 LOG_INFO << "StPrepEmbedMaker::SetSpreadMode set spread mode= ";
00557
00558 if( mSpreadMode ){
00559 LOG_INFO << " ON" << endm ;
00560 }
00561 else{
00562 LOG_INFO << " OFF" << endm ;
00563 }
00564 }
00565
00566
00567 void StPrepEmbedMaker::SetPrimeMode(const Bool_t flag)
00568 {
00569 mPrimeMode=flag;
00570
00571 LOG_INFO << "StPrepEmbedMaker::SetPrimeMode set prime mode= ";
00572
00573 if( mPrimeMode ){
00574 LOG_INFO << " ON" << endm ;
00575 }
00576 else{
00577 LOG_INFO << " OFF" << endm ;
00578 }
00579 }
00580
00581
00582 void StPrepEmbedMaker::SetRapidityMode(const Bool_t flag)
00583 {
00584 mRapidityMode=flag;
00585
00586 LOG_INFO << "StPrepEmbedMaker::SetRapidityMode set rapidity mode= ";
00587
00588 if( mRapidityMode ){
00589 LOG_INFO << " Rapidity" << endm ;
00590 }
00591 else{
00592 LOG_INFO << " Pseudo-rapidity" << endm ;
00593 }
00594 }
00595
00596
00597 void StPrepEmbedMaker::SetVpdVzCutMode(const Bool_t flag)
00598 {
00599
00600 mVpdVzCutMode=flag;
00601
00602 LOG_INFO << "StPrepEmbedMaker::SetVpdVzCutMode set VpdVz cut mode= ";
00603
00604 if( mVpdVzCutMode ){
00605 LOG_INFO << " ON" << endm ;
00606 }
00607 else{
00608 LOG_INFO << " OFF" << endm ;
00609 }
00610
00611 if(flag){
00612
00613 SetSpreadMode(flag);
00614 }
00615 }
00616
00617
00618 Int_t StPrepEmbedMaker::getMultiplicity(const StEvtHddr& EvtHddr, const Int_t nprimarytracks) const
00619 {
00621
00622 Int_t npart = 0;
00623 if(mSettings->mult < 1.)
00624 {
00625 npart=int(mSettings->mult * nprimarytracks);
00626 if (! npart)
00627 {
00628 LOG_INFO << "StPrepEmbedMaker::Event " << EvtHddr.GetEventNumber()
00629 << " has too small numberOfPrimaryTracks " << nprimarytracks << " for the mult fraction requested. Forcing npart to 1." << endm;
00630 npart=1;
00631 }
00632
00633 }
00634 else
00635 {
00636 npart = int (mSettings->mult);
00637 }
00638
00639 return npart ;
00640 }
00641
00642
00643 void StPrepEmbedMaker::SetTrgOpt(const Int_t TrgId)
00644 {
00645
00646 if(TrgId == 0){
00647 LOG_ERROR << "StPrepEmbedMaker::SetTrgOpt Input trigger id = 0. Skip" << endm;
00648 return;
00649 }
00650
00651 if(mSettings->NReqTrg >= mSettings->nTriggerId) {
00652 LOG_ERROR << "StPrepEmbedMaker::SetTrgOpt too many triggers are requested!" <<endm;
00653 return;
00654 }
00655
00656 mSettings->ReqTrgId[mSettings->NReqTrg] = TrgId ;
00657 mSettings->NReqTrg ++ ;
00658 LOG_INFO << "StPrepEmbedMaker::SetTrgOpt trigger " << mSettings->ReqTrgId[mSettings->NReqTrg-1] << " requested" << endm;
00659 }
00660
00661
00662 void StPrepEmbedMaker::SetZVertexCut(const Double_t vzlow, const Double_t vzhigh)
00663 {
00664
00665 if( (vzlow > vzhigh) || (vzlow == 0.0 && vzhigh == 0.0) ){
00666 LOG_ERROR << "StPrepEmbedMaker::SetZVertexCut input vzlow > vzhigh or vzlow = vzhigh = 0" << endm;
00667 return;
00668 }
00669
00670 mSettings->vzlow = vzlow ;
00671 mSettings->vzhigh = vzhigh ;
00672 LOG_INFO << "StPrepEmbedMaker::SetZVertexCut Cut z-vertex in " << mSettings->vzlow
00673 << " < vz < "
00674 << mSettings->vzhigh
00675 << " (cm)" << endm;
00676 }
00677
00678
00679 void StPrepEmbedMaker::SetVrCut(const Double_t vr)
00680 {
00681
00682 if( vr == 0.0 ){
00683 LOG_ERROR << "StPrepEmbedMaker::SetVrCut input vr = 0" << endm;
00684 return;
00685 }
00686
00687 mSettings->vr = vr ;
00688 LOG_INFO << "StPrepEmbedMaker::SetVrCut Cut vr in " << mSettings->vr
00689 << " (cm)" << endm;
00690 }
00691
00692
00693 void StPrepEmbedMaker::SetVpdVzCut(const Double_t vpdvz)
00694 {
00695
00696 if( vpdvz <= 0.0 ){
00697 LOG_ERROR << "StPrepEmbedMaker::SetVpdVzCut input |vpdvz-vz| <= 0" << endm;
00698 return;
00699 }
00700
00701 mSettings->vpdvz = vpdvz ;
00702 LOG_INFO << "StPrepEmbedMaker::SetVpdVzCut Cut |vpdvz-vz| in " << mSettings->vpdvz
00703 << " (cm)" << endm;
00704 }
00705
00706
00707 void StPrepEmbedMaker::OpenFzFile()
00708 {
00709
00710 mOpenFzFile = kTRUE ;
00711 LOG_INFO << "StPrepEmbedMaker::OpenFzFile Write .fz file. File basename will be taken "
00712 << "from daq file basename" << endm;
00713 }
00714
00715
00716
00717 void StPrepEmbedMaker::phasespace(const Int_t mult)
00718 {
00719 Double_t rapidityMin = mSettings->etalow ;
00720 Double_t rapidityMax = mSettings->etahigh ;
00721
00723 if( mSettings->mode.CompareTo("Spectrum", TString::kIgnoreCase) == 0 ){
00724 rapidityMin = -mRapidityMaximumCut ;
00725 rapidityMax = mRapidityMaximumCut ;
00726 }
00727
00728 Do( Form("phasespace %i %i %f %f %f %f;",
00729 mult, mSettings->pid,
00730 mSettings->ptlow, mSettings->pthigh,
00731 rapidityMin, rapidityMax)
00732 );
00733 }
00734
00735
00736 void StPrepEmbedMaker::gkine(const Int_t mult, const Double_t vzmin, const Double_t vzmax)
00737 {
00738 Double_t rapidityMin = mSettings->etalow ;
00739 Double_t rapidityMax = mSettings->etahigh ;
00740
00742 if( mSettings->mode.CompareTo("Spectrum", TString::kIgnoreCase) == 0 ){
00743 rapidityMin = -mRapidityMaximumCut ;
00744 rapidityMax = mRapidityMaximumCut ;
00745 }
00746
00747 Do( Form("gkine %i %i %f %f %f %f %f %f %f %f;",
00748 mult, mSettings->pid,
00749 mSettings->ptlow, mSettings->pthigh,
00750 rapidityMin, rapidityMax,
00751 mSettings->philow, mSettings->phihigh,
00752 vzmin, vzmax)
00753 );
00754 }
00755
00756
00757
00758
00759
00760
00761
00762
00763
00764
00765
00766
00767
00768
00769
00770
00771
00772
00773
00774
00775
00776
00777
00778
00779
00780
00781
00782
00783
00784
00785
00786
00787
00788
00789
00790
00791
00792
00793
00794
00795
00796
00797
00798
00799
00800
00801
00802
00803
00804
00805
00806
00807
00808
00809
00810
00811
00812
00813
00814
00815
00816
00817
00818
00819
00820
00821
00822
00823
00824
00825
00826
00827