1 /****************************************************************************
2 * Very important, delicate and rather obscure macro. *
4 * Creates list of "trackable" tracks, *
5 * sorts tracks for matching with the ITS, *
6 * calculates efficiency, resolutions etc. *
8 * Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch *
9 * with several nice improvements by: M.Ivanov, GSI, m.ivanov@gsi.de *
10 ****************************************************************************/
12 #if !defined(__CINT__) || defined(__MAKECINT__)
16 #include<TClassTable.h>
17 #include<TClonesArray.h>
27 #include <AliTPCtracker.h>
32 multievent comparison of reconstructed tracks with "exact tracks"
39 Int_t fEventN; //event number
46 enum tagprimary {kPrimaryCharged = 0x4000};
48 Int_t good_tracks(GoodTrackTPC *gt, Int_t max, Int_t firstev, Int_t eventn);
50 Int_t FindPrimaries(Int_t nparticles);
52 Int_t AliTPCComparison2(Int_t firstev=0, Int_t eventn=1) {
54 const Int_t MAX = 20000;
55 const Bool_t kOLD = kFALSE; // if tracking was done with a previous version
56 /***********************************************************************/
58 TFile *inkin = TFile::Open("rfio:galice.root");
59 // if(gAlice)delete gAlice; COMMENTED BECAUSE OF A BUG (IN COMPILED MODE)
60 gAlice = (AliRun*)inkin->Get("gAlice");
61 cout<<"AliRun object found on file "<<gAlice<<endl;
62 AliKalmanTrack::SetConvConst(1000/0.299792458/gAlice->Field()->SolenoidField());
65 delete gAlice; COMMENTED BECAUSE OF A BUG IN COMPILED MODE
68 /***********************************************************************/
69 cerr<<"Doing comparison...\n";
71 gBenchmark->Start("AliTPCComparison2");
76 cf=TFile::Open("AliTPCclusters.root");
77 if (!cf->IsOpen()) {cerr<<"Can't open AliTPCclusters.root !\n"; return 1;}
78 digp= (AliTPCParam*)cf->Get("75x40_100x60_150x60");
79 if (!digp) { cerr<<"TPC parameters have not been found !\n"; return 2; }
82 AliTPCtracker *tracker =0;
83 TObjArray tarray(MAX);
84 AliTPCtrack *iotrack=0;
87 TFile *tf=TFile::Open("AliTPCtracks.root");
91 eventptr[firstev+1]=0;
92 for (Int_t event=firstev;event<eventn; event++){
93 cout<<"================================================"<<endl;
94 cout<<"Processing event "<<event<<endl;
95 cout<<"================================================"<<endl;
98 tracker = new AliTPCtracker(digp,event);
99 tracker->LoadInnerSectors();
100 tracker->LoadOuterSectors();
105 sprintf(tname,"TreeT_TPC");
108 sprintf(tname,"TreeT_TPC_%d",event);
112 if (!tf->IsOpen()) {cerr<<"Can't open AliTPCtracks.root !\n"; return 3;}
114 tracktree=(TTree*)tf->Get(tname);
115 if (!tracktree) {cerr<<"Can't get a tree with TPC tracks !\n"; return 4;}
116 TBranch *tbranch=tracktree->GetBranch("tracks");
117 Int_t nentr0=(Int_t)tracktree->GetEntries();
119 for (i=0; i<nentr0; i++) {
120 iotrack=new AliTPCtrack;
121 tbranch->SetAddress(&iotrack);
122 tracktree->GetEvent(i);
123 if(kOLD)tracker->CookLabel(iotrack,0.1);
124 tarray.AddLast(iotrack);
126 eventptr[event+1] = nentr; //store end of the event
133 GoodTrackTPC gt[MAX];
136 ifstream in("good_tracks_tpc");
138 cerr<<"Reading good tracks...\n";
139 while (in>>gt[ngood].fEventN>>gt[ngood].lab>>gt[ngood].code>>
140 gt[ngood].px>>gt[ngood].py>>gt[ngood].pz>>
141 gt[ngood].x >>gt[ngood].y >>gt[ngood].z) {
146 cerr<<"Too many good tracks !\n";
151 if (!in.eof()) cerr<<"Read error (good_tracks_tpc) !\n";
153 cerr<<"Marking good tracks (this will take a while)...\n";
154 ngood=good_tracks(gt,45000,firstev,eventn);
155 printf("Goood %d\n", ngood);
156 ofstream out("good_tracks_tpc");
158 cout<<"File good_tracks_tpc opened\n";
159 for (Int_t ngd=0; ngd<ngood; ngd++) {
160 out<<gt[ngd].fEventN<<' '<<gt[ngd].lab<<' '<<gt[ngd].code<<' '<<
161 gt[ngd].px<<' '<<gt[ngd].py<<' '<<gt[ngd].pz<<' '<<
162 gt[ngd].x <<' '<<gt[ngd].y <<' '<<gt[ngd].z <<endl;
164 } else cerr<<"Can not open file (good_tracks_tpc) !\n";
168 ofstream out2("good_tracks_tpc_par");
171 //cout<<"File good_tracks_tpc opened\n";
172 for (Int_t ngd=0; ngd<ngood; ngd++) {
173 Float_t pt = TMath::Sqrt(gt[ngd].px*gt[ngd].px+gt[ngd].py*gt[ngd].py);
175 if (TMath::Abs(pt)>0.01) angle = gt[ngd].pz/pt;
176 out2<<gt[ngd].fEventN<<"\t"<<gt[ngd].lab<<"\t"<<gt[ngd].code<<"\t"<<
177 pt<<"\t"<<angle<<"\t"<<endl;
179 } else cerr<<"Can not open file (good_tracks_tpc) !\n";
184 cerr<<"Number of good tracks : "<<ngood<<endl;
185 cout<<"Number of good tracks : "<<ngood<<endl;
186 if(ngood==0)return 5;
187 TH1F *hp=new TH1F("hp","PHI resolution",50,-20.,20.); hp->SetFillColor(4);
188 TH1F *hl=new TH1F("hl","LAMBDA resolution",50,-20,20);hl->SetFillColor(4);
189 TH1F *hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.);
190 hpt->SetFillColor(2);
191 TH1F *hmpt=new TH1F("hmpt","Relative Pt resolution (pt>4GeV/c)",30,-60,60);
192 hmpt->SetFillColor(6);
194 AliTPCtrack *trk=(AliTPCtrack*)tarray.UncheckedAt(0);
196 Double_t pmin=0.1*(100/0.299792458/0.2/trk->GetConvConst());
197 Double_t pmax=6.0+pmin;
199 TH1F *hgood=new TH1F("hgood","Good tracks",30,pmin,pmax);
200 TH1F *hfound=new TH1F("hfound","Found tracks",30,pmin,pmax);
201 TH1F *hfake=new TH1F("hfake","Fake tracks",30,pmin,pmax);
202 TH1F *hg=new TH1F("hg","",30,pmin,pmax); //efficiency for good tracks
203 hg->SetLineColor(4); hg->SetLineWidth(2);
204 TH1F *hf=new TH1F("hf","Efficiency for fake tracks",30,pmin,pmax);
205 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
207 TH1F *he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);
208 TH2F *hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,250.);
209 hep->SetMarkerStyle(8);
210 hep->SetMarkerSize(0.4);
213 Int_t track_notfound[MAX], itrack_notfound=0;
214 Int_t track_fake[MAX], itrack_fake=0;
215 Int_t track_multifound[MAX], track_multifound_n[MAX], itrack_multifound=0;
217 Int_t lab=gt[ngood].lab,tlab=-1;
219 TMath::Sqrt(gt[ngood].px*gt[ngood].px + gt[ngood].py*gt[ngood].py);
221 if (ptg<pmin) continue; // Should be 1e-33 instead pmin?
224 Int_t ievent = gt[ngood].fEventN;
226 AliTPCtrack *track=0;
227 for (i=eventptr[ievent]; i<eventptr[ievent+1]; i++) {
229 track=(AliTPCtrack*)tarray.UncheckedAt(i);
230 tlab=track->GetLabel();
231 if (lab==TMath::Abs(tlab)) break;
233 if (i==eventptr[ievent+1]) {
234 track_notfound[itrack_notfound++]=lab;
240 AliTPCtrack * mitrack;
241 for (mi=eventptr[ievent]; mi<eventptr[ievent+1]; mi++) {
243 mitrack=(AliTPCtrack*)tarray.UncheckedAt(mi);
244 if (lab==TMath::Abs(mitrack->GetLabel())) micount++;
247 track_multifound[itrack_multifound]=lab;
248 track_multifound_n[itrack_multifound]=micount;
253 Double_t xk=gt[ngood].x;
254 if (!track) continue;
255 // printf("Track =%p\n",track);
256 track->PropagateTo(xk);
258 if (lab==tlab) hfound->Fill(ptg);
260 track_fake[itrack_fake++]=lab;
263 Double_t par[5]; track->GetExternalParameters(xk,par);
264 Float_t phi=TMath::ASin(par[2]) + track->GetAlpha();
265 if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
266 if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
267 Float_t lam=TMath::ATan(par[3]);
268 Float_t pt_1=TMath::Abs(par[4]);
270 if (TMath::Abs(gt[ngood].code)==11 && ptg>4.) {
271 hmpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
275 Float_t phig=TMath::ATan2(gt[ngood].py,gt[ngood].px);
276 hp->Fill((phi - phig)*1000.);
278 Float_t lamg=TMath::ATan2(gt[ngood].pz,ptg);
279 hl->Fill((lam - lamg)*1000.);
281 hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
284 Float_t mom=1./(pt_1*TMath::Cos(lam));
285 Float_t dedx=track->GetdEdx();
286 hep->Fill(mom,dedx,1.);
287 if (TMath::Abs(gt[ngood].code)==211)
288 if (mom>0.4 && mom<0.5) {
294 cout<<"\nList of Not found tracks :\n";
295 for ( i = 0; i< itrack_notfound; i++){
296 cout<<track_notfound[i]<<"\t";
297 if ((i%5)==4) cout<<"\n";
299 cout<<"\nList of fake tracks :\n";
300 for ( i = 0; i< itrack_fake; i++){
301 cout<<track_fake[i]<<"\t";
302 if ((i%5)==4) cout<<"\n";
304 cout<<"\nList of multiple found tracks :\n";
305 for ( i=0; i<itrack_multifound; i++) {
306 cout<<"id. "<<track_multifound[i]
307 <<" found - "<<track_multifound_n[i]<<"times\n";
310 Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
311 if (ng!=0) cerr<<"\n\nIntegral efficiency is about "<<nf/ng*100.<<" %\n";
312 if (ng!=0) cout<<"\n\nIntegral efficiency is about "<<nf/ng*100.<<" %\n";
313 cout<<"Total number of found tracks ="<<nentr<<endl;
314 cout<<"Total number of \"good\" tracks ="
315 <<mingood<<" (selected for comparison: "<<ng<<')'<<endl<<endl;
318 gStyle->SetOptStat(111110);
319 gStyle->SetOptFit(1);
321 TCanvas *c1=new TCanvas("c1","",0,0,700,850);
323 TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw();
324 p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10);
325 hp->SetFillColor(4); hp->SetXTitle("(mrad)"); hp->Fit("gaus"); c1->cd();
327 TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw();
328 p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10);
329 hl->SetXTitle("(mrad)"); hl->Fit("gaus"); c1->cd();
331 TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw();
332 p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10);
333 hpt->SetXTitle("(%)"); hpt->Fit("gaus"); c1->cd();
335 TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw();
336 p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10);
337 hmpt->SetXTitle("(%)"); hmpt->Fit("gaus"); c1->cd();
339 TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd();
340 p5->SetFillColor(41); p5->SetFrameFillColor(10);
341 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
342 hg->Divide(hfound,hgood,1,1.,"b");
343 hf->Divide(hfake,hgood,1,1.,"b");
345 hg->SetYTitle("Tracking efficiency");
346 hg->SetXTitle("Pt (GeV/c)");
349 TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
351 TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
355 hf->SetFillStyle(3013);
358 hf->Draw("histsame");
359 TText *text = new TText(0.461176,0.248448,"Fake tracks");
360 text->SetTextSize(0.05);
362 text = new TText(0.453919,1.11408,"Good tracks");
363 text->SetTextSize(0.05);
368 TCanvas *c2=new TCanvas("c2","",320,32,530,590);
370 TPad *p6=new TPad("p6","",0.,0.,1.,.5); p6->Draw();
371 p6->cd(); p6->SetFillColor(42); p6->SetFrameFillColor(10);
372 he->SetFillColor(2); he->SetFillStyle(3005);
373 he->SetXTitle("Arbitrary Units");
374 he->Fit("gaus"); c2->cd();
376 TPad *p7=new TPad("p7","",0.,0.5,1.,1.); p7->Draw();
377 p7->cd(); p7->SetFillColor(42); p7->SetFrameFillColor(10);
378 hep->SetXTitle("p (Gev/c)"); hep->SetYTitle("dE/dX (Arb. Units)");
379 hep->Draw(); c1->cd();
381 gBenchmark->Stop("AliTPCComparison2");
382 gBenchmark->Show("AliTPCComparison2");
389 Int_t good_tracks(GoodTrackTPC *gt, Int_t max, Int_t firstev, Int_t eventn) {
390 //eventn - number of events in file
392 TFile *file=TFile::Open("galice.root");
393 if (!file->IsOpen()) {cerr<<"Can't open galice.root !\n"; exit(4);}
394 // delete gAlice; gAlice = 0;
395 if (!(gAlice=(AliRun*)file->Get("gAlice"))) {
396 cerr<<"gAlice have not been found on galice.root !\n";
400 TFile *fdigit = TFile::Open("digits.root");
403 AliTPC *TPC=(AliTPC*)gAlice->GetDetector("TPC");
404 Int_t ver = TPC->IsVersion();
405 cerr<<"TPC version "<<ver<<" has been found !\n";
407 AliTPCParam *digp=(AliTPCParam*)file->Get("75x40_100x60");
409 cerr<<"2 pad-lenght geom hits with 3 pad-length geom digits...\n";
411 digp = new AliTPCParamSR();
415 digp =(AliTPCParam*)gDirectory->Get("75x40_100x60_150x60");
418 if (!digp) { cerr<<"TPC parameters have not been found !\n"; exit(6); }
421 Int_t nrow_up=digp->GetNRowUp();
422 Int_t nrows=digp->GetNRowLow()+nrow_up;
423 Int_t zero=digp->GetZeroSup();
424 Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap);
425 Int_t good_number=Int_t(0.4*nrows);
427 Int_t nt=0; //reset counter
428 char treeName[100]; //declare event identifier
430 for (Int_t event=firstev;event<eventn;event++){
431 Int_t np=gAlice->GetEvent(event);
432 Int_t nopr = FindPrimaries(np);
433 cout<<"function good_tracks -- event "<<event<<", Charged primaries:"<<nopr<<"\n";
435 Int_t *good=new Int_t[np];
436 for (Int_t ii=0; ii<np; ii++) good[ii]=0;
439 Int_t sectors_by_rows=0;
441 AliSimDigits da, *digits=&da;
446 TFile *cf=TFile::Open("AliTPCclusters.root");
447 if (!cf->IsOpen()){cerr<<"Can't open AliTPCclusters.root !\n";exit(5);}
448 AliTPCClustersArray *ca=new AliTPCClustersArray;
450 ca->SetClusterType("AliTPCcluster");
451 ca->ConnectTree("TreeC_TPC_0");
452 Int_t nrows=Int_t(ca->GetTree()->GetEntries());
453 for (Int_t n=0; n<nrows; n++) {
454 AliSegmentID *s=ca->LoadEntry(n);
456 digp->AdjustSectorRow(s->GetID(),sec,row);
457 AliTPCClustersRow &clrow = *ca->GetRow(sec,row);
458 Int_t ncl=clrow.GetArray()->GetEntriesFast();
460 AliTPCcluster *c=(AliTPCcluster*)clrow[ncl];
461 Int_t lab=c->GetLabel(0);
462 if (lab<0) continue; //noise cluster
464 if (sec>=digp->GetNInnerSector())
465 if (row==nrow_up-1 ) good[lab]|=0x1000;
466 if (sec>=digp->GetNInnerSector())
467 if (row==nrow_up-1-gap) good[lab]|=0x800;
470 ca->ClearRow(sec,row);
478 sprintf(treeName,"TreeD_75x40_100x60_150x60_%d",event);
479 // TD=(TTree*)gDirectory->Get(treeName);
480 TD=(TTree*)fdigit->Get(treeName);
481 TD->GetBranch("Segment")->SetAddress(&digits);
482 count = new Int_t[np];
484 for (i=0; i<np; i++) count[i]=0;
485 sectors_by_rows=(Int_t)TD->GetEntries();
486 for (i=0; i<sectors_by_rows; i++) {
487 if (!TD->GetEvent(i)) continue;
489 digp->AdjustSectorRow(digits->GetID(),sec,row);
491 digits->ExpandTrackBuffer();
492 do { //Many thanks to J.Chudoba who noticed this
493 Int_t it=digits->CurrentRow(), ip=digits->CurrentColumn();
494 // Short_t dig = digits->GetDigit(it,ip);
495 Short_t dig = digits->CurrentDigit();
497 Int_t idx0=digits->GetTrackIDFast(it,ip,0)-2;
498 Int_t idx1=digits->GetTrackIDFast(it,ip,1)-2;
499 Int_t idx2=digits->GetTrackIDFast(it,ip,2)-2;
500 if (idx0>=0 && dig>=zero && idx0<np ) count[idx0]+=1; //background events - higher track ID's
501 if (idx1>=0 && dig>=zero && idx1<np ) count[idx1]+=1;
502 if (idx2>=0 && dig>=zero && idx2<np ) count[idx2]+=1;
503 } while (digits->Next());
504 for (Int_t j=0; j<np; j++) {
506 if (sec>=digp->GetNInnerSector())
507 if (row==nrow_up-1 ) good[j]|=0x4000;
508 if (sec>=digp->GetNInnerSector())
509 if (row==nrow_up-1-gap) good[j]|=0x1000;
511 if (sec>=digp->GetNInnerSector())
512 if (row==nrow_up-1-shift) good[j]|=0x2000;
513 if (sec>=digp->GetNInnerSector())
514 if (row==nrow_up-1-gap-shift) good[j]|=0x800;
522 delete TD; //Thanks to Mariana Bondila
526 cerr<<"Invalid TPC version !\n";
531 /** select tracks which are "good" enough **/
532 //printf("\t %d \n",np);
533 for (Int_t i=0; i<np; i++) {
534 if ((good[i]&0x5000) != 0x5000)
535 if ((good[i]&0x2800) != 0x2800) continue;
536 if ((good[i]&0x7FF ) < good_number) continue;
538 TParticle *p = (TParticle*)gAlice->Particle(i);
540 if (p->Pt()<0.100) continue;
541 if (TMath::Abs(p->Pz()/p->Pt())>0.999) continue;
542 // THIS IS GOOD IF YOU WANT SECONDARIES
543 Int_t j=p->GetFirstMother();
545 TParticle *pp = (TParticle*)gAlice->Particle(j);
546 //if (!(pp->TestBit(kPrimaryCharged))) continue; //only one decay is allowed
549 if(!(p->TestBit(kPrimaryCharged)))continue; // only primaries
552 gt[nt].fEventN=event;
554 gt[nt].code=p->GetPdgCode();
555 gt[nt].px=0.; gt[nt].py=0.; gt[nt].pz=0.;
556 gt[nt].x=0.; gt[nt].y=0.; gt[nt].z=0.;
558 if (nt==max) {cerr<<"Too many good tracks !\n"; break;}
562 /** check if there is also information at the entrance of the TPC **/
563 TTree *TH=gAlice->TreeH(); np=(Int_t)TH->GetEntries();
565 for (Int_t i=0; i<np; i++) {
568 AliTPChit *phit = (AliTPChit*)TPC->FirstHit(-1);
569 for ( ; phit; phit=(AliTPChit*)TPC->NextHit() ) {
570 if (phit->fQ !=0. ) continue;
572 Double_t px=phit->X(), py=phit->Y(), pz=phit->Z();
574 if ((phit=(AliTPChit*)TPC->NextHit())==0) break;
575 if (phit->fQ != 0.) continue;
577 Double_t x=phit->X(), y=phit->Y(), z=phit->Z();
578 if (TMath::Sqrt(x*x+y*y)>90.) continue;
580 Int_t j, lab=phit->Track();
581 for (j=0; j<nt; j++) {if (gt[j].fEventN==event && gt[j].lab==lab) break;}
584 // (px,py,pz) - in global coordinate system, (x,y,z) - in local !
585 gt[j].px=px; gt[j].py=py; gt[j].pz=pz;
586 Float_t cs,sn; digp->AdjustCosSin(phit->fSector,cs,sn);
587 gt[j].x = x*cs + y*sn;
588 gt[j].y =-x*sn + y*cs;
593 //printf("\n%d\n",nt);
596 } // /// loop on events
597 delete gAlice; gAlice=0;
600 gBenchmark->Stop("AliTPCComparison2");
601 gBenchmark->Show("AliTPCComparison2");
606 Int_t FindPrimaries(Int_t nparticles){
609 Double_t vertcut = 0.001;
610 Double_t decacut = 3.;
611 Double_t timecut = 0.;
613 TParticle * part = gAlice->Particle(0);
614 Double_t xori = part->Vx();
615 Double_t yori = part->Vy();
616 Double_t zori = part->Vz();
617 for(Int_t iprim = 0; iprim<nparticles; iprim++){ //loop on tracks
619 part = gAlice->Particle(iprim);
620 Double_t ptot=TMath::Sqrt(part->Px()*part->Px()+part->Py()*part->Py()+part->Pz()*part->Pz());
621 if (ptot<0.01) continue;
622 char *xxx=strstr(part->GetName(),"XXX");
625 TParticlePDG *ppdg = part->GetPDG();
626 if(TMath::Abs(ppdg->Charge())!=3)continue; // only charged (no quarks)
628 Double_t dist=TMath::Sqrt((part->Vx()-xori)*(part->Vx()-xori)+(part->Vy()-yori)*(part->Vy()-yori)+(part->Vz()-zori)*(part->Vz()-zori));
629 if(dist>vertcut)continue; // cut on the vertex
631 if(part->T()>timecut)continue;
633 //Double_t ptot=TMath::Sqrt(part->Px()*part->Px()+part->Py()*part->Py()+part->Pz()*part->Pz());
634 if(ptot==(TMath::Abs(part->Pz())))continue; // no beam particles
636 Bool_t prmch = kTRUE; // candidate primary track
637 Int_t fidau=part->GetFirstDaughter(); // cut on daughters
641 lasdau=part->GetLastDaughter();
645 for(Int_t j=fidau;j<=lasdau;j++){
646 TParticle *dau=gAlice->Particle(j);
647 Double_t distd=TMath::Sqrt((dau->Vx()-xori)*(dau->Vx()-xori)+(dau->Vy()-yori)*(dau->Vy()-yori)+(dau->Vz()-zori)*(dau->Vz()-zori));
648 if(distd<decacut)prmch=kFALSE; // eliminate if the decay is near the vertex
654 part->SetBit(kPrimaryCharged);