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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // Base class for ALICE modules. Both sensitive modules (Modules) and //
21 // non-sensitive ones are described by this base class. This class //
22 // supports the hit and digit trees produced by the simulation and also //
23 // the objects produced by the reconstruction. //
25 // This class is also responsible for building the geometry of the //
30 <img src="picts/AliModuleClass.gif">
34 ///////////////////////////////////////////////////////////////////////////////
37 #include <TObjArray.h>
38 #include <TClonesArray.h>
41 #include <TDirectory.h>
42 #include <TVirtualMC.h>
43 #include <TGeoManager.h>
46 #include "AliConfig.h"
47 #include "AliLoader.h"
49 #include "AliModule.h"
51 #include "AliTrackReference.h"
53 #include "../RAW/AliRawDataHeader.h"
57 //_______________________________________________________________________
58 AliModule::AliModule():
72 fCurrentIterTrackRef(0),
76 // Default constructor for the AliModule class
80 //_______________________________________________________________________
81 AliModule::AliModule(const char* name,const char *title):
85 fIdtmed(new TArrayI(100)),
86 fIdmate(new TArrayI(100)),
90 fHistograms(new TList()),
94 fTrackReferences(new TClonesArray("AliTrackReference", 100)),
96 fCurrentIterTrackRef(0),
100 // Normal constructor invoked by all Modules.
101 // Create the list for Module specific histograms
102 // Add this Module to the global list of Modules in Run.
104 // Get the Module numeric ID
105 Int_t id = gAlice->GetModuleID(name);
107 // Module already added !
108 AliWarning(Form("Module: %s already present at %d",name,id));
112 // Add this Module to the list of Modules
114 gAlice->AddModule(this);
118 // Clear space for tracking media and material indexes
120 for(Int_t i=0;i<100;i++) (*fIdmate)[i]=(*fIdtmed)[i]=0;
123 //_______________________________________________________________________
124 AliModule::AliModule(const AliModule &mod):
142 fCurrentIterTrackRef(0),
151 //_______________________________________________________________________
152 AliModule::~AliModule()
158 // Remove this Module from the list of Modules
160 TObjArray * modules = gAlice->Modules();
161 if (modules) modules->Remove(this);
163 // Delete ROOT geometry
171 fHistograms->Clear();
175 // Delete track references
176 if (fTrackReferences) {
177 fTrackReferences->Delete();
178 delete fTrackReferences;
179 fTrackReferences = 0;
181 // Delete TArray objects
187 //_______________________________________________________________________
188 void AliModule::Copy(TObject & /* mod */) const
191 // Copy *this onto mod, not implemented for AliModule
193 AliFatal("Not implemented!");
196 //_______________________________________________________________________
197 void AliModule::Disable()
200 // Disable Module on viewer
206 // Loop through geometry to disable all
207 // nodes for this Module
208 while((node = dynamic_cast<TNode*>(next()))) {
209 node->SetVisibility(-1);
213 //_______________________________________________________________________
214 void AliModule::Enable()
217 // Enable Module on the viewver
223 // Loop through geometry to enable all
224 // nodes for this Module
225 while((node = dynamic_cast<TNode*>(next()))) {
226 node->SetVisibility(3);
230 //_______________________________________________________________________
231 void AliModule::AliMaterial(Int_t imat, const char* name, Float_t a,
232 Float_t z, Float_t dens, Float_t radl,
233 Float_t absl, Float_t *buf, Int_t nwbuf) const
236 // Store the parameters for a material
238 // imat the material index will be stored in (*fIdmate)[imat]
239 // name material name
243 // radl radiation length
244 // absl absorbtion length
245 // buf adress of an array user words
246 // nwbuf number of user words
249 TString uniquename = GetName();
250 uniquename.Append(name);
251 if(gAlice->IsRootGeometry()){
252 TGeoMaterial *mat = gGeoManager->GetMaterial(uniquename.Data());
253 kmat = mat->GetUniqueID();
254 (*fIdmate)[imat]=kmat;
256 gMC->Material(kmat, uniquename.Data(), a, z, dens, radl, absl, buf, nwbuf);
257 (*fIdmate)[imat]=kmat;
261 //_______________________________________________________________________
262 void AliModule::AliGetMaterial(Int_t imat, char* name, Float_t &a,
263 Float_t &z, Float_t &dens, Float_t &radl,
267 // Store the parameters for a material
269 // imat the material index will be stored in (*fIdmate)[imat]
270 // name material name
274 // radl radiation length
275 // absl absorbtion length
276 // buf adress of an array user words
277 // nwbuf number of user words
282 kmat=(*fIdmate)[imat];
283 gMC->Gfmate(kmat, name, a, z, dens, radl, absl, buf, nwbuf);
286 //_______________________________________________________________________
287 void AliModule::AliMixture(Int_t imat, const char *name, Float_t *a,
288 Float_t *z, Float_t dens, Int_t nlmat,
292 // Defines mixture or compound imat as composed by
293 // nlmat materials defined by arrays a, z and wmat
295 // If nlmat > 0 wmat contains the proportion by
296 // weights of each basic material in the mixture
298 // If nlmat < 0 wmat contains the number of atoms
299 // of eack kind in the molecule of the compound
300 // In this case, wmat is changed on output to the relative weigths.
302 // imat the material index will be stored in (*fIdmate)[imat]
303 // name material name
304 // a array of atomic masses
305 // z array of atomic numbers
307 // nlmat number of components
308 // wmat array of concentrations
311 TString uniquename = GetName();
312 uniquename.Append(name);
313 if(gAlice->IsRootGeometry()){
314 TGeoMaterial *mat = gGeoManager->GetMaterial(uniquename.Data());
315 kmat = mat->GetUniqueID();
316 (*fIdmate)[imat]=kmat;
318 gMC->Mixture(kmat, uniquename.Data(), a, z, dens, nlmat, wmat);
319 (*fIdmate)[imat]=kmat;
323 //_______________________________________________________________________
324 void AliModule::AliMedium(Int_t numed, const char *name, Int_t nmat,
325 Int_t isvol, Int_t ifield, Float_t fieldm,
326 Float_t tmaxfd, Float_t stemax, Float_t deemax,
327 Float_t epsil, Float_t stmin, Float_t *ubuf,
331 // Store the parameters of a tracking medium
333 // numed the medium number is stored into (*fIdtmed)[numed]
335 // nmat the material number is stored into (*fIdmate)[nmat]
336 // isvol sensitive volume if isvol!=0
337 // ifield magnetic field flag (see below)
338 // fieldm maximum magnetic field
339 // tmaxfd maximum deflection angle due to magnetic field
340 // stemax maximum step allowed
341 // deemax maximum fractional energy loss in one step
342 // epsil tracking precision in cm
343 // stmin minimum step due to continuous processes
345 // ifield = 0 no magnetic field
346 // = -1 user decision in guswim
347 // = 1 tracking performed with Runge Kutta
348 // = 2 tracking performed with helix
349 // = 3 constant magnetic field along z
352 TString uniquename = GetName();
353 uniquename.Append(name);
354 if(gAlice->IsRootGeometry()){
355 TList *medialist = gGeoManager->GetListOfMedia();
356 TIter next(medialist);
358 while((med = (TGeoMedium*) next())){
359 if(!strcmp(uniquename.Data(),med->GetName())){
361 (*fIdtmed)[numed]=kmed;
366 gMC->Medium(kmed, uniquename.Data(), (*fIdmate)[nmat], isvol, ifield,
367 fieldm, tmaxfd, stemax, deemax, epsil, stmin, ubuf, nbuf);
368 (*fIdtmed)[numed]=kmed;
372 //_______________________________________________________________________
373 void AliModule::AliMatrix(Int_t &nmat, Float_t theta1, Float_t phi1,
374 Float_t theta2, Float_t phi2, Float_t theta3,
378 // Define a rotation matrix. Angles are in degrees.
380 // nmat on output contains the number assigned to the rotation matrix
381 // theta1 polar angle for axis I
382 // phi1 azimuthal angle for axis I
383 // theta2 polar angle for axis II
384 // phi2 azimuthal angle for axis II
385 // theta3 polar angle for axis III
386 // phi3 azimuthal angle for axis III
388 gMC->Matrix(nmat, theta1, phi1, theta2, phi2, theta3, phi3);
391 //_______________________________________________________________________
392 Float_t AliModule::ZMin() const
397 //_______________________________________________________________________
398 Float_t AliModule::ZMax() const
403 //_______________________________________________________________________
404 void AliModule::SetEuclidFile(char* material, char* geometry)
407 // Sets the name of the Euclid file
409 fEuclidMaterial=material;
411 fEuclidGeometry=geometry;
413 char* name = new char[strlen(material)];
414 strcpy(name,material);
415 strcpy(&name[strlen(name)-4],".euc");
416 fEuclidGeometry=name;
421 //_______________________________________________________________________
422 void AliModule::ReadEuclid(const char* filnam, char* topvol)
425 // read in the geometry of the detector in euclid file format
427 // id_det : the detector identification (2=its,...)
428 // topvol : return parameter describing the name of the top
429 // volume of geometry.
434 // several changes have been made by miroslav helbich
435 // subroutine is rewrited to follow the new established way of memory
436 // booking for tracking medias and rotation matrices.
437 // all used tracking media have to be defined first, for this you can use
438 // subroutine greutmed.
439 // top volume is searched as only volume not positioned into another
442 Int_t i, nvol, iret, itmed, irot, numed, npar, ndiv, iaxe;
443 Int_t ndvmx, nr, flag;
444 char key[5], card[77], natmed[21];
445 char name[5], mother[5], shape[5], konly[5], volst[7000][5];
448 Float_t teta1, phi1, teta2, phi2, teta3, phi3, orig, step;
450 const Int_t kMaxRot=5000;
451 Int_t idrot[kMaxRot],istop[7000];
454 // *** The input filnam name will be with extension '.euc'
455 filtmp=gSystem->ExpandPathName(filnam);
456 lun=fopen(filtmp,"r");
459 AliError(Form("Could not open file %s",filnam));
462 //* --- definition of rotation matrix 0 ---
463 TArrayI &idtmed = *fIdtmed;
464 for(i=1; i<kMaxRot; ++i) idrot[i]=-99;
468 for(i=0;i<77;i++) card[i]=0;
469 iret=fscanf(lun,"%77[^\n]",card);
470 if(iret<=0) goto L20;
475 if (!strcmp(key,"TMED")) {
476 sscanf(&card[5],"%d '%[^']'",&itmed,natmed);
477 if( itmed<0 || itmed>=100 ) {
478 AliError(Form("TMED illegal medium number %d for %s",itmed,natmed));
481 //Pad the string with blanks
484 while(i<20) natmed[i++]=' ';
487 if( idtmed[itmed]<=0 ) {
488 AliError(Form("TMED undefined medium number %d for %s",itmed,natmed));
491 gMC->Gckmat(idtmed[itmed],natmed);
493 } else if (!strcmp(key,"ROTM")) {
494 sscanf(&card[4],"%d %f %f %f %f %f %f",&irot,&teta1,&phi1,&teta2,&phi2,&teta3,&phi3);
495 if( irot<=0 || irot>=kMaxRot ) {
496 AliError(Form("ROTM rotation matrix number %d illegal",irot));
499 AliMatrix(idrot[irot],teta1,phi1,teta2,phi2,teta3,phi3);
501 } else if (!strcmp(key,"VOLU")) {
502 sscanf(&card[5],"'%[^']' '%[^']' %d %d", name, shape, &numed, &npar);
504 for(i=0;i<npar;i++) fscanf(lun,"%f",&par[i]);
507 gMC->Gsvolu( name, shape, idtmed[numed], par, npar);
508 //* save the defined volumes
509 strcpy(volst[++nvol],name);
512 } else if (!strcmp(key,"DIVN")) {
513 sscanf(&card[5],"'%[^']' '%[^']' %d %d", name, mother, &ndiv, &iaxe);
514 gMC->Gsdvn ( name, mother, ndiv, iaxe );
516 } else if (!strcmp(key,"DVN2")) {
517 sscanf(&card[5],"'%[^']' '%[^']' %d %d %f %d",name, mother, &ndiv, &iaxe, &orig, &numed);
518 gMC->Gsdvn2( name, mother, ndiv, iaxe, orig,idtmed[numed]);
520 } else if (!strcmp(key,"DIVT")) {
521 sscanf(&card[5],"'%[^']' '%[^']' %f %d %d %d", name, mother, &step, &iaxe, &numed, &ndvmx);
522 gMC->Gsdvt ( name, mother, step, iaxe, idtmed[numed], ndvmx);
524 } else if (!strcmp(key,"DVT2")) {
525 sscanf(&card[5],"'%[^']' '%[^']' %f %d %f %d %d", name, mother, &step, &iaxe, &orig, &numed, &ndvmx);
526 gMC->Gsdvt2 ( name, mother, step, iaxe, orig, idtmed[numed], ndvmx );
528 } else if (!strcmp(key,"POSI")) {
529 sscanf(&card[5],"'%[^']' %d '%[^']' %f %f %f %d '%[^']'", name, &nr, mother, &xo, &yo, &zo, &irot, konly);
530 if( irot<0 || irot>=kMaxRot ) {
531 Error("ReadEuclid","POSI %s#%d rotation matrix number %d illegal\n",name,nr,irot);
534 if( idrot[irot] == -99) {
535 Error("ReadEuclid","POSI %s#%d undefined matrix number %d\n",name,nr,irot);
538 //*** volume name cannot be the top volume
539 for(i=1;i<=nvol;i++) {
540 if (!strcmp(volst[i],name)) istop[i]=0;
543 gMC->Gspos ( name, nr, mother, xo, yo, zo, idrot[irot], konly );
545 } else if (!strcmp(key,"POSP")) {
546 sscanf(&card[5],"'%[^']' %d '%[^']' %f %f %f %d '%[^']' %d", name, &nr, mother, &xo, &yo, &zo, &irot, konly, &npar);
547 if( irot<0 || irot>=kMaxRot ) {
548 Error("ReadEuclid","POSP %s#%d rotation matrix number %d illegal\n",name,nr,irot);
551 if( idrot[irot] == -99) {
552 Error("ReadEuclid","POSP %s#%d undefined matrix number %d\n",name,nr,irot);
556 for(i=0;i<npar;i++) fscanf(lun,"%f",&par[i]);
559 //*** volume name cannot be the top volume
560 for(i=1;i<=nvol;i++) {
561 if (!strcmp(volst[i],name)) istop[i]=0;
564 gMC->Gsposp ( name, nr, mother, xo,yo,zo, idrot[irot], konly, par, npar);
567 if (strcmp(key,"END")) goto L10;
568 //* find top volume in the geometry
570 for(i=1;i<=nvol;i++) {
571 if (istop[i] && flag) {
572 AliWarning(Form(" %s is another possible top volume",volst[i]));
574 if (istop[i] && !flag) {
575 strcpy(topvol,volst[i]);
576 AliDebug(2, Form("volume %s taken as a top volume",topvol));
581 AliWarning("top volume not found");
585 //* commented out only for the not cernlib version
586 AliDebug(1, Form("file: %s is now read in",filnam));
591 AliError("reading error or premature end of file");
594 //_______________________________________________________________________
595 void AliModule::ReadEuclidMedia(const char* filnam)
598 // read in the materials and tracking media for the detector
599 // in euclid file format
601 // filnam: name of the input file
602 // id_det: id_det is the detector identification (2=its,...)
604 // author : miroslav helbich
606 Float_t sxmgmx = gAlice->Field()->Max();
607 Int_t isxfld = gAlice->Field()->Integ();
608 Int_t end, i, iret, itmed;
609 char key[5], card[130], natmed[21], namate[21];
614 Int_t nwbuf, isvol, ifield, nmat;
615 Float_t a, z, dens, radl, absl, fieldm, tmaxfd, stemax, deemax, epsil, stmin;
618 for(i=0;i<end;i++) if(filnam[i]=='.') {
623 // *** The input filnam name will be with extension '.euc'
624 AliDebug(1, Form("The file name is %s",filnam)); //Debug
625 filtmp=gSystem->ExpandPathName(filnam);
626 lun=fopen(filtmp,"r");
629 AliWarning(Form("Could not open file %s",filnam));
633 // Retrieve Mag Field parameters
634 Int_t globField=gAlice->Field()->Integ();
635 Float_t globMaxField=gAlice->Field()->Max();
636 // TArrayI &idtmed = *fIdtmed;
639 for(i=0;i<130;i++) card[i]=0;
640 iret=fscanf(lun,"%4s %[^\n]",key,card);
641 if(iret<=0) goto L20;
645 if (!strcmp(key,"MATE")) {
646 sscanf(card,"%d '%[^']' %f %f %f %f %f %d",&imate,namate,&a,&z,&dens,&radl,&absl,&nwbuf);
647 if (nwbuf>0) for(i=0;i<nwbuf;i++) fscanf(lun,"%f",&ubuf[i]);
648 //Pad the string with blanks
651 while(i<20) namate[i++]=' ';
654 AliMaterial(imate,namate,a,z,dens,radl,absl,ubuf,nwbuf);
655 //* read tracking medium
656 } else if (!strcmp(key,"TMED")) {
657 sscanf(card,"%d '%[^']' %d %d %d %f %f %f %f %f %f %d",
658 &itmed,natmed,&nmat,&isvol,&ifield,&fieldm,&tmaxfd,
659 &stemax,&deemax,&epsil,&stmin,&nwbuf);
660 if (nwbuf>0) for(i=0;i<nwbuf;i++) fscanf(lun,"%f",&ubuf[i]);
661 if (ifield<0) ifield=isxfld;
662 if (fieldm<0) fieldm=sxmgmx;
663 //Pad the string with blanks
666 while(i<20) natmed[i++]=' ';
669 AliMedium(itmed,natmed,nmat,isvol,globField,globMaxField,tmaxfd,
670 stemax,deemax,epsil,stmin,ubuf,nwbuf);
671 // (*fImedia)[idtmed[itmed]-1]=id_det;
675 if (strcmp(key,"END")) goto L10;
678 //* commented out only for the not cernlib version
679 AliDebug(1, Form("file %s is now read in",filnam));
684 AliWarning("reading error or premature end of file");
687 //_______________________________________________________________________
688 void AliModule::RemapTrackReferencesIDs(Int_t *map)
691 // Remapping track reference
692 // Called at finish primary
694 if (!fTrackReferences) return;
695 for (Int_t i=0;i<fTrackReferences->GetEntries();i++){
696 AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
698 Int_t newID = map[ref->GetTrack()];
699 if (newID>=0) ref->SetTrack(newID);
702 ref->SetBit(kNotDeleted,kFALSE);
703 fTrackReferences->RemoveAt(i);
707 fTrackReferences->Compress();
712 //_______________________________________________________________________
713 AliTrackReference* AliModule::FirstTrackReference(Int_t track)
716 // Initialise the hit iterator
717 // Return the address of the first hit for track
718 // If track>=0 the track is read from disk
719 // while if track<0 the first hit of the current
724 if (fRunLoader == 0x0)
725 AliFatal("AliRunLoader not initialized. Can not proceed");
726 fRunLoader->GetAliRun()->GetMCApp()->ResetTrackReferences();
727 fRunLoader->TreeTR()->GetEvent(track);
730 fMaxIterTrackRef = fTrackReferences->GetEntriesFast();
731 fCurrentIterTrackRef = 0;
732 if(fMaxIterTrackRef) return dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(0));
736 //_______________________________________________________________________
737 AliTrackReference* AliModule::NextTrackReference()
740 // Return the next hit for the current track
742 if(fMaxIterTrackRef) {
743 if(++fCurrentIterTrackRef<fMaxIterTrackRef)
744 return dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(fCurrentIterTrackRef));
748 AliWarning("Iterator called without calling FistTrackReference before");
754 //_______________________________________________________________________
755 void AliModule::ResetTrackReferences()
758 // Reset number of hits and the hits array
760 fMaxIterTrackRef = 0;
761 if (fTrackReferences) fTrackReferences->Clear();
764 //_____________________________________________________________________________
766 AliLoader* AliModule::MakeLoader(const char* /*topfoldername*/)
771 //PH Merged with v3-09-08 |
773 //_____________________________________________________________________________
775 void AliModule::SetTreeAddress()
778 // Set branch address for track reference Tree
783 // Branch address for track reference tree
784 TTree *treeTR = TreeTR();
786 if (treeTR && fTrackReferences) {
787 branch = treeTR->GetBranch(GetName());
790 AliDebug(3, Form("(%s) Setting for TrackRefs",GetName()));
791 branch->SetAddress(&fTrackReferences);
795 //can be called before MakeBranch and than does not make sense to issue the warning
796 AliDebug(1, Form("(%s) Failed for Track References. Can not find branch in tree.",
802 //_____________________________________________________________________________
803 void AliModule::AddTrackReference(Int_t label){
805 // add a trackrefernce to the list
806 if (!fTrackReferences) {
807 AliError("Container trackrefernce not active");
810 Int_t nref = fTrackReferences->GetEntriesFast();
811 TClonesArray &lref = *fTrackReferences;
812 new(lref[nref]) AliTrackReference(label);
816 //_____________________________________________________________________________
817 void AliModule::MakeBranchTR(Option_t */*option*/)
820 // Makes branch in treeTR
822 AliDebug(2,Form("Making Track Refs. Branch for %s",GetName()));
823 TTree * tree = TreeTR();
824 if (fTrackReferences && tree)
826 TBranch *branch = tree->GetBranch(GetName());
829 AliDebug(2,Form("Branch %s is already in tree.",GetName()));
833 branch = tree->Branch(GetName(),&fTrackReferences);
837 AliDebug(2,Form("FAILED for %s: tree=%#x fTrackReferences=%#x",
838 GetName(),tree,fTrackReferences));
842 //_____________________________________________________________________________
843 TTree* AliModule::TreeTR()
846 // Return TR tree pointer
848 if ( fRunLoader == 0x0)
850 AliError("Can not get the run loader");
854 TTree* tree = fRunLoader->TreeTR();
859 //_____________________________________________________________________________
860 void AliModule::Digits2Raw()
862 // This is a dummy version that just copies the digits file contents
863 // to a raw data file.
865 AliWarning(Form("Dummy version called for %s", GetName()));
867 const Int_t kNDetectors = 17;
868 const char* kDetectors[kNDetectors] = {"TPC", "ITSSPD", "ITSSDD", "ITSSSD", "TRD", "TOF", "PHOS", "RICH", "EMCAL", "MUON", "MUTR", "ZDC", "PMD", "START", "VZERO", "CRT", "FMD"};
869 const Int_t kDetectorDDLs[kNDetectors] = {216, 20, 12, 16, 18, 72, 20, 20, 22, 20, 2, 1, 6, 1, 1, 1, 3};
872 for (Int_t i = 0; i < kNDetectors; i++) {
873 if (strcmp(GetName(), kDetectors[i]) == 0) {
874 nDDLs = kDetectorDDLs[i];
875 ddlOffset = 0x100 * i;
879 if (!GetLoader()) return;
880 fstream digitsFile(GetLoader()->GetDigitsFileName(), ios::in);
881 if (!digitsFile) return;
883 digitsFile.seekg(0, ios::end);
884 UInt_t size = digitsFile.tellg();
885 UInt_t ddlSize = 4 * (size / (4*nDDLs));
886 Char_t* buffer = new Char_t[ddlSize+1];
888 for (Int_t iDDL = 0; iDDL < nDDLs; iDDL++) {
890 sprintf(fileName, "%s_%d.ddl", GetName(), iDDL + ddlOffset);
891 fstream rawFile(fileName, ios::out);
892 if (!rawFile) return;
894 AliRawDataHeader header;
895 header.fSize = ddlSize + sizeof(header);
896 rawFile.write((char*) &header, sizeof(header));
898 digitsFile.read(buffer, ddlSize);
899 rawFile.write(buffer, ddlSize);
908 //_____________________________________________________________________________
909 Int_t AliModule::GetDebug() const
911 AliWarning("Don't use this method any more, use AliDebug instead");