//******************************************************** // DFDCPseudo_factory.cc - factory producing first-order // reconstructed points for the FDC // Author: Craig Bookwalter (craigb at jlab.org) // Date: March 2006 // UVX cathode-wire matching revised by Simon Taylor, Aug 2006 //******************************************************** #include "DFDCPseudo_factory.h" #include "HDGEOMETRY/DGeometry.h" #include "DFDCGeometry.h" #define HALF_CELL 0.5 #define MAX_DEFLECTION 0.15 #define X0 0 #define QA 1 #define K2 2 #define ITER_MAX 100 #define TOLX 1e-4 #define TOLF 1e-4 #define A_OVER_H 0.475 #define ALPHA 1e-4 // rate parameter for Newton step backtracking algorithm #define CHARGE_TO_ENERGY 0.01 //place holder /// /// DFDCCathodeCluster_gLayer_cmp(): /// non-member function passed to std::sort() to sort DFDCCathodeCluster pointers /// by their gLayer attributes. /// bool DFDCCathodeCluster_gLayer_cmp(const DFDCCathodeCluster* a, const DFDCCathodeCluster* b) { return a->gLayer < b->gLayer; } /// /// DFDCAnode_gLayer_cmp(): /// non-member function passed to std::sort() to sort DFDCHit pointers /// for the anode wires by their gLayer attributes. /// bool DFDCAnode_gLayer_cmp(const DFDCHit* a, const DFDCHit* b) { return a->gLayer < b->gLayer; } /// /// DFDCPseudo_factory::DFDCPseudo_factory(): /// default constructor -- initializes log file /// DFDCPseudo_factory::DFDCPseudo_factory() { logFile = new ofstream("DFDCPseudo_factory.log"); _log = new JStreamLog(*logFile, "PSEUDO"); *_log << "File initialized." << endMsg; } /// /// DFDCPseudo_factory::~DFDCPseudo_factory(): /// default destructor -- closes log file /// DFDCPseudo_factory::~DFDCPseudo_factory() { logFile->close(); delete logFile; delete _log; } //------------------ // brun //------------------ jerror_t DFDCPseudo_factory::brun(JEventLoop *loop, int runnumber) { // Get pointer to DGeometry object DApplication* dapp=dynamic_cast(eventLoop->GetJApplication()); const DGeometry *dgeom = dapp->GetDGeometry(runnumber); dgeom->GetFDCWires(fdcwires); return NOERROR; } /// /// DFDCPseudo_factory::evnt(): /// this is the place that anode hits and DFDCCathodeClusters are organized into pseudopoints. /// jerror_t DFDCPseudo_factory::evnt(JEventLoop* eventLoop, int eventNo) { vector fdcHits; vector xHits; vector cathClus; vector uClus; vector oneLayerU; vector vClus; vector oneLayerV; vector oneLayerX; // Get all FDC hits (anode and cathode) eventLoop->Get(fdcHits); // For events with a very large number of hits, assume // we can't reconstruct them so bail early // Feb. 8, 2008 D.L. if(fdcHits.size()>(5.0+5.0+1.0)*25.0*6.0){ _DBG_<<"Too many hits in FDC! Psuedopoint reconstruction in FDC bypassed for event "<GetJEvent().GetEventNumber()<Get(cathClus); // Ensure clusters are in order of ascending Z position. std::sort(cathClus.begin(), cathClus.end(), DFDCCathodeCluster_gLayer_cmp); // Sift through hits and select out anode hits. for (unsigned int i=0; i < fdcHits.size(); i++) if (fdcHits[i]->type == 0) xHits.push_back(fdcHits[i]); // Make sure the wires are also in order of ascending z position std::sort(xHits.begin(), xHits.end(), DFDCAnode_gLayer_cmp); // Sift through clusters and put U and V clusters into respective vectors. for (unsigned int i=0; i < cathClus.size(); i++) { if (cathClus[i]->plane == 1) vClus.push_back(cathClus[i]); else uClus.push_back(cathClus[i]); } vector::iterator uIt = uClus.begin(); vector::iterator vIt = vClus.begin(); vector::iterator xIt = xHits.begin(); // For each layer, get its sets of V, X, and U hits, and then pass them to the geometrical // organization routine, DFDCPseudo_factory::makePseudo() for (int iLayer=1; iLayer <= 24; iLayer++) { for (; ((uIt != uClus.end() && (*uIt)->gLayer == iLayer)); uIt++) oneLayerU.push_back(*uIt); for (; ((vIt != vClus.end() && (*vIt)->gLayer == iLayer)); vIt++) oneLayerV.push_back(*vIt); for (; ((xIt != xHits.end() && (*xIt)->gLayer == iLayer)); xIt++) oneLayerX.push_back(*xIt); if (oneLayerU.size()>0 && oneLayerV.size()>0 && oneLayerX.size()>0) makePseudo(oneLayerX, oneLayerU, oneLayerV,iLayer); oneLayerU.clear(); oneLayerV.clear(); oneLayerX.clear(); } return NOERROR; } /// /// DFDCPseudo_factory::makePseudo(): /// performs UV+X matching to create pseudopoints /// void DFDCPseudo_factory::makePseudo(vector& x, vector& u, vector& v, int layer) { vector::iterator xIt; vector::iterator uIt; vector::iterator vIt; centroid_t temp; float E1,E2,pos1,pos2; // Loop over all U and V clusters looking for peaks upeaks.clear(); for (uIt = u.begin(); uIt != u.end(); uIt++){ vector::const_iterator strip=(*uIt)->members.begin(); switch((*uIt)->members.size()){ case 0: // Make sure we have data!! break; case 1: // One isolated hit in the cluster: use element number itself temp.pos=(*strip)->element; temp.q=2.*((*strip)->q); // Each cathode view should see half the // anode charge temp.numstrips=1; CalcMeanTime((*uIt)->members, temp.t, temp.t_rms); upeaks.push_back(temp); break; case 2: //Two adjacent hits: use average for the centroid pos1=(*strip)->element; pos2=(*(strip+1))->element; E1=(*strip)->q; E2=(*(strip+1))->q; temp.pos=(pos1*E1+pos2*E2)/(E1+E2); temp.q=2.*(E1+E2); temp.numstrips=2; CalcMeanTime((*uIt)->members, temp.t, temp.t_rms); upeaks.push_back(temp); break; default: for (strip=(*uIt)->members.begin(); strip!=(*uIt)->members.end();strip++){ FindCentroid((*uIt)->members,strip,upeaks); } break; } } vpeaks.clear(); for (vIt = v.begin(); vIt != v.end(); vIt++){ vector::const_iterator strip=(*vIt)->members.begin(); switch((*vIt)->members.size()){ case 0: // Make sure we have data!! break; case 1: // One isolated hit in the cluster: use element number itself temp.pos=(*strip)->element; temp.q=2.*((*strip)->q); temp.numstrips=1; CalcMeanTime((*vIt)->members, temp.t, temp.t_rms); vpeaks.push_back(temp); break; case 2: //Two adjacent hits: use average for the centroid pos1=(*strip)->element; pos2=(*(strip+1))->element; E1=(*strip)->q; E2=(*(strip+1))->q; temp.pos=(pos1*E1+pos2*E2)/(E1+E2); temp.q=2.*(E1+E2); temp.numstrips=2; CalcMeanTime((*vIt)->members, temp.t, temp.t_rms); vpeaks.push_back(temp); break; default: for (strip=(*vIt)->members.begin(); strip!=(*vIt)->members.end();strip++){ FindCentroid((*vIt)->members,strip,vpeaks); } break; } } if (upeaks.size()*vpeaks.size()>0){ //Loop over all u and v centroids looking for matches with wires for (unsigned int i=0;ielement<=WIRES_PER_PLANE && (*xIt)->element>0){ // First check if mean times of strips and wire are close. // There are 3 values to compare so we look at the RMS // of the 3 differences. (I'm just making this up!) // 1/2/2008 D. L. float dt1 = (*xIt)->t - upeaks[i].t; float dt2 = (*xIt)->t - vpeaks[j].t; float dt3 = upeaks[i].t - vpeaks[j].t; float trms = sqrt((dt1*dt1 + dt2*dt2 + dt3*dt3)/3.0); if(trms>20.0)continue; float x_from_wire=DFDCGeometry::getWireR(*xIt); if (fabs(x_from_wire-x_from_strips)w = x_from_wire; newPseu->dw = xres; newPseu->s = y_from_strips; newPseu->ds = yres; newPseu->wire = fdcwires[layer-1][(*xIt)->element-1]; newPseu->time = (*xIt)->t; newPseu->dist = newPseu->time*DRIFT_SPEED; newPseu->status = status; newPseu->dE = CHARGE_TO_ENERGY*(upeaks[i].q+vpeaks[j].q)/2.; // It can occur (although rarely) that newPseu->wire is NULL // which causes us to crash below. In these cases, we can't really // make a psuedo point so we delete the current object // and just go on to the next one. if(newPseu->wire==NULL){ _DBG_<<"newPseu->wire=NULL! This shouldn't happen. Complain to staylor@jlab.org"<wire->udir(0); cosangle=newPseu->wire->udir(1); newPseu->x=(newPseu->w)*cosangle+(newPseu->s)*sinangle; newPseu->y=-(newPseu->w)*sinangle+(newPseu->s)*cosangle; double sigx2=HALF_CELL*HALF_CELL/3.; double sigy2=MAX_DEFLECTION*MAX_DEFLECTION/3.; newPseu->covxx=sigx2*cosangle*cosangle+sigy2*sinangle*sinangle; newPseu->covyy=sigx2*sinangle*sinangle+sigy2*cosangle*cosangle; newPseu->covxy=(sigy2-sigx2)*sinangle*cosangle; _data.push_back(newPseu); } // match in x } else _DBG_ << "Bad wire " << (*xIt)->element <& H, float &t, float &t_rms) { // Calculate mean t=0.0; for(unsigned int i=0; it; if(H.size()>0)t/=(float)H.size(); // Calculate RMS t_rms=0.0; for(unsigned int i=0; it-t),2.0); if(H.size()>0)t_rms = sqrt(t_rms/(float)H.size()); } // /// DFDCPseudo_factory::FindCentroid() /// Uses the Newton-Raphson method for solving the set of non-linear /// equations describing the charge distribution over 3 strips for the peak /// position x0, the anode charge qa, and the "width" parameter k2. /// See Numerical Recipes in C p.379-383. /// Updates list of centroids. /// jerror_t DFDCPseudo_factory::FindCentroid(const vector& H, vector::const_iterator peak, vector¢roids){ centroid_t temp; // Make sure we do not exceed the range of the vector if (peak>H.begin() && peak+1!=H.end()){ float err_diff1=0.,err_diff2=0.; // Fill in time info in temp 1/2/2008 D.L. CalcMeanTime(H, temp.t, temp.t_rms); // Some code for checking for significance of fluctuations. // Currently disabled. //float dq1=(*(peak-1))->dq; //float dq2=(*peak)->dq; //float dq3=(*(peak+1))->dq; //err_diff1=sqrt(dq1*dq1+dq2*dq2); //err_diff2=sqrt(dq2*dq2+dq3*dq3); // Check for a peak in three adjacent strips if ((*peak)->q-(*(peak-1))->q > err_diff1 && (*peak)->q-(*(peak+1))->q > err_diff2){ // Define some matrices for use in the Newton-Raphson iteration DMatrix J(3,3); //Jacobean matrix DMatrix F(3,1),N(3,1),X(3,1),par(3,1); DMatrix newpar(3,1); int i=0; double sum=0.; // Initialize the matrices to some suitable starting values par(X0,0)=double((*peak)->element); par(K2,0)=1.; for (vector::const_iterator j=peak-1;j<=peak+1;j++){ X(i,0)=double((*j)->element); N(i++,0)=double((*j)->q); sum+=double((*j)->q); } par(QA,0)=2.*sum; newpar=par; // Newton-Raphson procedure double errf=0.,errx=0; for (int iter=1;iter<=ITER_MAX;iter++){ errf=0.; errx=0.; // Compute Jacobian matrix: J_ij = dF_i/dx_j. for (i=0;i<3;i++){ double argp=par(K2,0)*(par(X0,0)-X(i,0)+A_OVER_H); double argm=par(K2,0)*(par(X0,0)-X(i,0)-A_OVER_H); J(i,QA)=-(tanh(argp)-tanh(argm))/4.; J(i,K2)=-par(QA,0)/4.*(argp/par(K2,0)*(1.-tanh(argp)*tanh(argp)) -argm/par(K2,0)*(1.-tanh(argm)*tanh(argm))); J(i,X0)=-par(QA,0)*par(K2,0)/4. *(tanh(argm)*tanh(argm)-tanh(argp)*tanh(argp)); F(i,0)=N(i,0)-par(QA,0)/4.*(tanh(argp)-tanh(argm)); if (fabs(F(i,0))>errf) errf=fabs(F(i,0)); } // Check for convergence if (errf