Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Class List | Directories | File List | Namespace Members | Class Members | File Members | Related Pages

SpillTimeFinder.cxx

Go to the documentation of this file.
00001 #include "SpillTimeFinder.h"
00002 #include "Conventions/Munits.h"
00003 #include "MessageService/MsgService.h"
00004 #include "DatabaseInterface/DbiSqlContext.h"
00005 #include <cmath>
00006 
00007 CVSID( "$Id: SpillTimeFinder.cxx,v 1.35 2008/05/30 18:35:43 tagg Exp $" );
00008 
00009 SpillTimeFinder* SpillTimeFinder::fgInstance = 0;
00010 
00011 static SpillTimeND kSpillTime_VeryOld(VldTimeStamp(0,0));             // 1970
00012 static SpillTimeND kSpillTime_Future(VldTimeStamp(2147480000,0));    // 2038
00013 
00014 SpillTimeFinder& SpillTimeFinder::Instance()
00015 {
00018 
00019   // Cleaner destructor calls DigitSpillTimeFinder dtor
00020   static Cleaner cleaner;
00021   
00022   // SpillTimeFinder mode "1" the default.  This is hardwired for now.
00023   if (!fgInstance) {
00024     cleaner.UseMe();
00025     fgInstance = new SpillTimeFinder();
00026     if (!fgInstance) {
00027       MSG("SpillTime", Msg::kError)
00028         << "No DigitSpillTimeFinder Instance - fatal." << endl;
00029       assert(fgInstance);               // Kill job is there is no instance
00030     }
00031   }
00032   return *fgInstance;
00033 }
00034 
00035 
00036 SpillTimeFinder::SpillTimeFinder() :
00037   fPrevTable(0),
00038   fNextTable(0),
00039   fLastQueryPrev(kSpillTime_Future),
00040   fLastQueryNext(kSpillTime_VeryOld)
00041 {
00042   Registry r;
00043   r.Set("Task",3);
00044   InitializeConfig(r);
00045 }
00046 
00047 SpillTimeFinder::~SpillTimeFinder()
00048 {
00049   if(fPrevTable) delete fPrevTable;
00050   if(fNextTable) delete fNextTable;
00051 }
00052 
00053 void SpillTimeFinder::ResetCache()
00054 {
00055   fLastQueryPrev=kSpillTime_Future;
00056   fLastQueryNext=kSpillTime_VeryOld;
00057 }
00058 
00059 void SpillTimeFinder::ConfigModified()
00060 {
00061   bool ok = true;
00062   ok = ok && GetConfig().Get("Task",fTask);
00063   if(!ok) MSG("SpillTimeFinder",Msg::kWarning) << "Problem configuring SpillTimeFinder." << endl;  
00064   ResetCache();
00065 }
00066 
00067 Bool_t SpillTimeFinder::DataIsAvailable( const VldContext& context )
00068 {
00073 
00074   // This is the relevant time at the kicker fire $74 MIBS
00075   VldTimeStamp kickerTime = context.GetTimeStamp();
00076   Double_t offset = GetOffset(context,fTask);
00077   kickerTime.Add(-offset); // Move requested time back to get kicker time.
00078 
00079   VldContext kickerContext(Detector::kNear,
00080                            context.GetSimFlag(),
00081                            kickerTime);
00082 
00083   // Query on the current context.
00084   fCurrentTable.NewQuery(kickerContext,fTask);
00085   
00086   // <Refresh cache here if any.>
00087 
00088 
00089   const DbiValidityRec* myVrec = fCurrentTable.GetValidityRec();
00090   if ( myVrec->IsGap () ) return false; // No data available
00091   
00092   return true;  // Data available.
00093 }
00094 
00095 
00096 
00097 void SpillTimeFinder::FindClosestEntries(const VldContext& context,
00098                                          const SpillTimeND* &outPrev,
00099                                          const SpillTimeND* &outNext )
00100 {
00101   
00106 
00107   if(context.GetSimFlag()==SimFlag::kMC) 
00108       MSG("SpillTimeFinder",Msg::kError) << "MC context in SpillTimeFinder request. This will all end in tears, I know it..." << endl;
00109 
00110   outPrev = outNext = 0;
00111 
00112   // This is the relevant time at the kicker fire $74 MIBS
00113   VldTimeStamp kickerTime = context.GetTimeStamp();
00114   Double_t offset = GetOffset(context,fTask);
00115   kickerTime.Add(-offset); // Move requested time back to get kicker time.
00116 
00117   VldContext kickerContext(Detector::kNear,
00118                            context.GetSimFlag(),
00119                            kickerTime );
00120 
00121   MSG("SpillTimeFinder",Msg::kVerbose) 
00122     << "STF Query: " << context.AsString() << endl;
00123   MSG("SpillTimeFinder",Msg::kVerbose) 
00124     << "Kicker Context: " << kickerContext.AsString() << endl;
00125   
00126   
00127   // Optimization:
00128   // If this request is in (prev,next) of the last call,
00129   // the same prev,next results should still apply.
00130   if(kickerTime < fLastQueryNext.GetTimeStamp()) {
00131     if(kickerTime > fLastQueryPrev.GetTimeStamp()) {
00132       // Match! Just give the last results.
00133       MSG("SpillTimeFinder",Msg::kDebug) << "Using last call's cached result." << endl;
00134       outPrev = &fLastQueryPrev;
00135       outNext = &fLastQueryNext;
00136       MSG("SpillTimeFinder",Msg::kVerbose) << "Result: " << outPrev->GetTimeStamp().AsString() 
00137                                            << "  " << outNext->GetTimeStamp().AsString() << endl;
00138 
00139       return;
00140     }     
00141   }
00142 
00143 
00144   // Query on the current context.
00145   fCurrentTable.NewQuery(kickerContext,fTask);
00146 
00147   // <Refresh cache here if any.>
00148 
00149   MSG("SpillTimeFinder",Msg::kDebug) << "Requested context, got " << fCurrentTable.GetNumRows() << " rows" << endl;  
00150   // See what rows we can get.
00151   FindBestRows(fCurrentTable,kickerTime, outPrev, outNext);
00152 
00153   MSG("SpillTimeFinder",Msg::kVerbose) << "First query: Prev = "
00154                                        << ((outPrev==0)?"<null>":outPrev->GetTimeStamp().AsString())
00155                                        << "   " 
00156                                        << ((outNext==0)?"<null>":outNext->GetTimeStamp().AsString())
00157                                        << endl;
00158 
00159 
00160 
00161   // Are we done?
00162   if((outPrev)&&(outNext)) return;
00163 
00164 
00166   // Ok, we need to hunt for more info.
00167 
00168   // Default: if we didn't find a table at all, set the search param
00169   // as the requested param.
00170   VldTimeStamp curStartVldTime;
00171   VldTimeStamp curStopVldTime;
00172 
00173   // Get the time limits of the current query.
00174   const DbiValidityRec* curValidity = fCurrentTable.GetValidityRec();
00175   if(curValidity) {
00176     // Use the start and stop times here.
00177     curStartVldTime = curValidity->GetVldRange().GetTimeStart();
00178     curStopVldTime = curValidity->GetVldRange().GetTimeEnd();
00179   } else {
00180     curStartVldTime = kickerTime;
00181     curStopVldTime = kickerTime;
00182   }
00183 
00184   fLastQueryRangeBeg = curStartVldTime; 
00185   fLastQueryRangeBeg.Add(GetOffset(context,fTask));
00186   fLastQueryRangeEnd = curStopVldTime;  
00187   fLastQueryRangeEnd.Add(GetOffset(context,fTask));
00188 
00189 
00190   // Did we miss the early side?
00191   while(outPrev ==0) {
00192     if(curStartVldTime < VldTimeStamp(1975,0,0,0,0,0)){
00193       MSG("SpillTimeFinder",Msg::kDebug) << "Hit start of vld range. Hit starting bookend (<1975 AD)" << endl;
00194       outPrev = &kSpillTime_VeryOld;
00195       break;
00196     }
00197     
00198     MSG("SpillTimeFinder",Msg::kDebug) << "Hit start of vld range.. looking earlier..." << endl;
00199     // Look for an earlier table.
00200     const char* sqltxt = Form("(TIMEEND<='%s') and (TASK=%d) and (DETECTORMASK & %d) and (SIMMASK & %d) order by TIMEEND desc limit 1",
00201                               curStartVldTime.AsString("s"),
00202                               fTask,
00203                               kickerContext.GetDetector(),
00204                               kickerContext.GetSimFlag() );
00205     MSG("SpillTimeFinder",Msg::kDebug) << "Request: " << sqltxt << endl;
00206     
00207     DbiSqlContext myContext(sqltxt);
00208     if(fPrevTable) {
00209       fPrevTable->NewQuery(myContext,Dbi::kAnyTask);
00210     } else {
00211       fPrevTable = new DbiResultPtr<SpillTimeND>("SPILLTIMEND",myContext);
00212     }
00213     
00214     if(fPrevTable->GetNumRows()>0) {
00215       // Found a table with data! Look for the best row..
00216       MSG("SpillTimeFinder",Msg::kDebug) << "...earlier table gotten. Rows: " << fPrevTable->GetNumRows() << endl;
00217       const SpillTimeND* n;    
00218       FindBestRows(*fPrevTable,kickerTime,outPrev,n);
00219       if(n>0) MSG("SpillTimeFinder",Msg::kError) << "Conflict! Found 'next' time in 'previous' table: " << n->GetTimeStamp().AsString() << endl;    
00220       
00221     } else {
00222       VldTimeStamp nexttry = fPrevTable->GetValidityRec()->GetVldRange().GetTimeStart();
00223       if(nexttry >= curStartVldTime) nexttry = VldTimeStamp(curStartVldTime.GetSec()-1,0); // protect against infinite loop.
00224       curStartVldTime = nexttry;
00225     }
00226   }
00227  
00228 
00229   // Did we miss the late side?
00230   while(outNext ==0) {
00231     MSG("SpillTimeFinder",Msg::kDebug) << "Hit end of vld range.. looking later..." << endl;
00232     const char* sqltxt = Form("(TIMESTART>='%s') and (TASK=%d) and (DETECTORMASK & %d) and (SIMMASK & %d) order by TIMESTART asc limit 1",
00233                               curStopVldTime.AsString("s"),
00234                               fTask,
00235                               kickerContext.GetDetector(),
00236                               kickerContext.GetSimFlag()
00237                               );
00238     
00239     MSG("SpillTimeFinder",Msg::kDebug) << "Request: " << sqltxt << endl;
00240     
00241     DbiSqlContext myContext(sqltxt);
00242     if(fNextTable)
00243       fNextTable->NewQuery(myContext,Dbi::kAnyTask);
00244     else 
00245       fNextTable = new DbiResultPtr<SpillTimeND>("SPILLTIMEND",myContext);
00246     
00247     if(fNextTable->GetNumRows()>0) {
00248       // Got a table with data! Look for the best row:
00249       MSG("SpillTimeFinder",Msg::kDebug) << "...later table gotten. Rows: " << fNextTable->GetNumRows() << endl;
00250       const SpillTimeND* p;
00251       
00252       FindBestRows(*fNextTable,kickerTime,p,outNext);
00253       if(p>0) MSG("SpillTimeFinder",Msg::kError) << "Conflict! Found 'previous' time in 'next' table: " << p->GetTimeStamp().AsString() << endl;
00254       
00255     } else {
00256      VldTimeStamp nexttry = fNextTable->GetValidityRec()->GetVldRange().GetTimeStart();
00257      if(nexttry <= curStopVldTime) nexttry = VldTimeStamp(curStopVldTime.GetSec()+1,0); // protect against infinite loop.
00258      curStopVldTime = nexttry;
00259     }
00260   }
00261 
00262   // Sanity check.
00263   if(outPrev->GetTimeStamp() > kickerTime) MSG("SpillTimeFinder",Msg::kError) << "Sanity check failed. Prev time is late!";
00264   if(outNext->GetTimeStamp() < kickerTime) MSG("SpillTimeFinder",Msg::kError) << "Sanity check failed. Next time is early!";
00265 
00266   // Load the cache.
00267   fLastQueryPrev = *outPrev;
00268   fLastQueryNext = *outNext;
00269 
00270   MSG("SpillTimeFinder",Msg::kVerbose) << "Result: " << outPrev->GetTimeStamp().AsString() 
00271                                        << "  " << outNext->GetTimeStamp().AsString() << endl;
00272 
00273 }
00274 
00275 
00276 
00277 Double_t SpillTimeFinder::GetOffset(const VldContext& context, Int_t task) 
00278 {
00285 
00286   Double_t offset = 0;
00287   
00288   switch(context.GetDetector()) {
00289   case Detector::kFar:
00290     offset += GetOffsetNDNuToFDNu(context); //fallthrough
00291   case Detector::kNear:
00292     if(task == SpillTimeND::kTask_Vtm) {
00293       offset += GetOffsetSgateToNDNu(context);
00294     } else {
00295       offset += GetOffsetKickerToNDNu(context); // Get offset from kicker to ND
00296     }
00297     break;
00298   default:
00299     break;    
00300   }
00301   return offset;
00302 }
00303 
00304 Double_t SpillTimeFinder::GetOffsetNDNuToFDNu(const VldContext& cx) 
00305 {
00311   // paper, DocDB 1870.
00312   
00313   // A note on signs:
00314   // A positive offset indicates the event will appear to be later at the FD
00315   // A negative offset indicates the event will appear to be earlier
00316   // Things that make the ND clock run slow will give a positive offset
00317   // Things that make the FD clock run slow will give a negative offset
00318 
00319   // New version - May 2008. Get offset calibration from DB, if available
00320   
00321   double offset = 0;
00322   fCalibrationNDNuToFDNu.NewQuery(cx,0);
00323   if(fCalibrationNDNuToFDNu.GetNumRows()>0) {
00324     static int sCall = 0;
00325     if(sCall==0) MSG("SpillTimeFinder",Msg::kInfo) << "Initial Spill Timing Calibration constants (NDNuToFDNu):" << endl;
00326     for(unsigned int i = 0; i<fCalibrationNDNuToFDNu.GetNumRows(); i++) {
00327       if(sCall==0) MSG("SpillTimeFinder",Msg::kInfo) 
00328         << Form("%10.5e",fCalibrationNDNuToFDNu.GetRow(i)->GetOffset())
00329         << " +- "
00330         << fCalibrationNDNuToFDNu.GetRow(i)->GetError() 
00331         << "\t"
00332         << fCalibrationNDNuToFDNu.GetRow(i)->GetDescription() << endl;
00333 
00334       offset += fCalibrationNDNuToFDNu.GetRow(i)->GetOffset();
00335     }
00336     sCall++;
00337     return offset;
00338   }
00339   
00340   MAXMSG("SpillTimeFinder",Msg::kWarning,1) << "Could not find SPILLTIMECALIBRATION table." << endl;
00341   MAXMSG("SpillTimeFinder",Msg::kWarning,1) << "Resorting to hard-coded calibration." << endl;
00342 
00343   const VldTimeStamp ndSwitch(2006,7,27,14,40,0,0); // Time Alfons changed ND cable.
00344   
00346   // ND-FD time of flight, 734,298.62 meters (Wes Smart), at c
00347   // Error is 0.7m, or 2.3ns.  
00348   offset += 2.44935653818215832*Munits::millisecond;
00349 
00350   //--------------------------------------
00351   // ND latencies
00352 
00354   // ND antenna:
00355   // OLD: total guess is 300ft, with same prop. speed as FD antenna.
00356   // This delays counter reset at ND
00357   // Error is +-60%, or 0.3us. (It's unlikely to be longer than 500ft, or shorter than 200)
00358   //+0.5*Munits::microsecond 
00359   // NEW April 5,2006. measurement by Chuck Andrews, reported by Andrew Rader
00360   // Fibre is in fact 765ft +/- 1% between the upstairs tranciever and U/G DAQ rack.
00361   // Andy gives the refrativce index as 1.4677 at 1300nm. 
00362   // Error is is 1%, or 12ns.
00363   // offset += 1.1415*Munits::microsecond;
00364   //
00365   // NEW Dec 6 2006:
00366   // Measured again with the BOLT. Got a new length
00367   // of 227.5m. Multiply by 1.470/3e8 s/m (the BOLT conversion constant).
00368   // So offset =  1.1155us +- 28ns
00369   // Take difference between measurements as the error: 
00370   //   1.1415-1.1155 = 0.027us
00371   offset += 1.1147*Munits::microsecond;
00372   
00374   // ND cable run from DAQ racks to GPS reciever:
00375   if(cx.GetTimeStamp() < ndSwitch) {
00376     // Cable run from ND tranciever to reciever (DAQ rack to timing rack)
00377     // Cat J sez: 5 segments of 32ns RG59 unioned together. +/-10% or 16 ns.
00378     offset += 0.160*Munits::microsecond;
00379   } else {
00380     // Alfons: we actually replaced 2m optical patch cable in the beams division crate
00381     // with a 50m optical patch cable, which leads now to the timing rack.
00382     // The chain of N x 10' RG58 is still in the cable tray, but not used.
00383     // Instead, the optical->elect converter is now in the timing rack and
00384     // connected with a (2+/-0.5) m RG58 to the GPS receiver    
00385     offset += 6.7*Munits::nanosecond // 2m electrical cable
00386       +0.2446*Munits::microsecond; // 50m optical cable
00387     // Another ~5ns of error here.
00388   }
00389   
00391   // ND timing system latency, MCC to MTM
00392   // But, it takes finite time to propagate CNTRST to the front end.
00393   // Tom Fitzpatrick measured this too (email 23/05/05) as being 450ns
00394   // Take reasonable error as +-50ns
00395   offset += 0.450*Munits::microsecond;
00396 
00398   // ND timing system latency, MINDER
00399   // The timing signals from the PMT to the timestamper is
00400   // finely timed in, so as to be less than 1 bucket. --Dave Reyna
00401   // Take offset as 10+-10ns
00402   offset += 10*Munits::nanosecond;
00403 
00405   // CandStrip reconstruction bias.
00406   // The ND measures the mean time of each strip, not the time of the first digit
00407   // The time of the first digit is a better estimate (compared to the FD), so correct for this
00408   // The bias varies as a function of pulse height, but mean bias is 7ns, +- 9ns rms
00409   offset += -7.0*Munits::nanosecond;
00410     
00413   // FD antenna: re-measured by Dave Saranen. Fibre length ~3415ft Email 23/05/05
00414   //
00415   // This delays counter reset at FD.
00416   // Error is on precision of the ping measurment: +-0.1 us
00417   // offset += -5.1*Munits::microsecond;
00418   
00419   // Much better measurement by Alec Habig, using our custom-bought BOLT unit.
00420   // Email Nov 8, 2006
00421   // Bolt gave a measurement of 1049m for the antenna run, and uses a group velocity of 
00422   // n_eff=1.470 and c=3e8
00423   // Disagreement between Orlando (FNAL) and BOLT measurement in for the
00424   // multimode fibre between VAX room and rescue room is 4m.  Error on BOLT is claimed to
00425   // be 2.5m.  Error on Orlando's measurement was 6m. 
00426   // Use differences in multimode fibre to see the error.
00427   offset += -5140*Munits::nanosecond;
00428   
00430   // FD timing system latency. Measured by using the RawVaTimingMonitorBlock
00431   // which gives a stable value of 639980381 TDC counts for the PPS fiducial tick.
00432   // Precision is +-1.5ns, but varies from run to run within ~100 ns.
00433   // Negative: the FD clock is even slower because of this latency.
00434   offset += -30.6547*Munits::microsecond;
00435   
00437   // New, Sept 2006:FD timing system fiducial cable length.
00438   // Makes the above measurement just a wee 
00439   // bit wrong. +-6ns (the rise time of the test pulse).
00440   // Negative.. it means the fiducial round-trip was even longer than we thought.
00441   offset += -153.0 * Munits::nanosecond; 
00442 
00444   // New Dec 2006: Average timing constants at the FD are not exactly
00445   // zero, but have an average offset of -27ns. This makes the events appear to be too early.
00446   // Error is half the offset, 14ns.
00447   offset += -27 * Munits::nanosecond;
00448   
00449   
00451   // Differences between ND and FD transceivers and recivers.
00452   // Difference = 0 +/- 20ns.
00453 
00454 
00456   // Detector effects. Total strip half-lenth plus WLS pigtail plus clear fibre readout
00457   // is 5m@ND and 8m@FD. (Very roughly; depends on the actual centroid of
00458   // accepted events vertices.) Thus, the events will appear at the FD about
00459   // 16ns later than events in the ND, due to the optical transit time.
00460   // Error on this is about 1m, or 5ns.
00461   offset += 16.0*Munits::nanosecond;
00462 
00463 
00465   // PMT transit time difference.
00466   // Both M16 and M64 have a 10.9ns transit time. Error is 1ns.
00467   
00468   
00469   return offset;
00470 }
00471 
00472 Double_t SpillTimeFinder::GetOffsetKickerToNDNu(const VldContext& cx) 
00473 {
00478 
00479   double offset = 0;
00480   fCalibrationKickerToNDNu.NewQuery(cx,1);
00481   if(fCalibrationKickerToNDNu.GetNumRows()>0) {
00482     static int sCall = 0;
00483     if(sCall==0) MSG("SpillTimeFinder",Msg::kInfo) << "Initial Spill Timing Calibration constants (KickerToNDNu):" << endl;
00484     for(unsigned int i = 0; i<fCalibrationKickerToNDNu.GetNumRows(); i++) {
00485       if(sCall==0) MSG("SpillTimeFinder",Msg::kInfo) 
00486         << Form("%10.5e",fCalibrationKickerToNDNu.GetRow(i)->GetOffset())
00487         << " +- "
00488         << fCalibrationKickerToNDNu.GetRow(i)->GetError() 
00489         << "\t"
00490         << fCalibrationKickerToNDNu.GetRow(i)->GetDescription() << endl;
00491 
00492       offset += fCalibrationKickerToNDNu.GetRow(i)->GetOffset();
00493     }
00494     sCall++;
00495     return offset;
00496   }
00497 
00498 
00500   const double kOffset = 
00501     // The following two values are taken directly from the 
00502     // TimeMaster control settings: SgateDelayRevs and SgateDelayFine:
00503     +19*11.064*Munits::microsecond // $74 MIBS until actual kicker fire, in MI turns
00504     +6.2*Munits::microsecond       // Fine tuning of SGATE window to $74
00505 
00506     // This tuning gets the first neutrinos at 1.5us into the SGATE
00507     // so add this on.  Error is +-0.1 us
00508     //+1.5*Munits::microsecond
00509     // Fixed 7/6/2006: retune based upon TOF study; move by -2.4 us
00510     -0.9*Munits::microsecond
00511 
00512     // This delay is not precise, because it does not truely wait an integer
00513     // number of turns, but rather waits for the NEXT n turns. Therefore, this 
00514     // delay is actually shorter. Tom Fitzpatrick measured this as being ~1us. (email 23/05/05)
00515     // Error not known. Estimate worst-case as 100%, or 1.0us
00516     - 1.0*Munits::microsecond
00517 
00518     // Now we have time between $74 and SGATE at the MCC. 
00519     // But, it takes finite time to propagate SGATE to the front end.
00520     // Tom Fitzpatrick measured this too (email 23/05/05) as being 450ns
00521     // Error is negligable compared to other offsets
00522     - 0.450*Munits::microsecond
00523 
00524   ;
00525 
00526   return kOffset;
00527 }
00528 
00529 Double_t SpillTimeFinder::GetOffsetSgateToNDNu(const VldContext& cx) 
00530 {
00535 
00536   double offset = 0;
00537   fCalibrationSgateToNDNu.NewQuery(cx,2);
00538   if(fCalibrationSgateToNDNu.GetNumRows()>0) {
00539     static int sCall = 0;
00540     if(sCall==0) MSG("SpillTimeFinder",Msg::kInfo) << "Initial Spill Timing Calibration constants (SgateToNDNu):" << endl;
00541     for(unsigned int i = 0; i<fCalibrationSgateToNDNu.GetNumRows(); i++) {
00542       if(sCall==0) MSG("SpillTimeFinder",Msg::kInfo) 
00543         << Form("%10.5e",fCalibrationSgateToNDNu.GetRow(i)->GetOffset())
00544         << " +- "
00545         << fCalibrationSgateToNDNu.GetRow(i)->GetError() 
00546         << "\t"
00547         << fCalibrationSgateToNDNu.GetRow(i)->GetDescription() << endl;
00548 
00549       offset += fCalibrationSgateToNDNu.GetRow(i)->GetOffset();
00550     }
00551     sCall++;
00552     return offset;
00553   }
00554 
00555 
00556   // This changes due to changing the fine SGATE delay at the ND.
00557   // Tuning this doesn' affect the offset above, but only changes
00558   // when the neutrinos are seen relative to the SGATE window.
00559    UInt_t sec = cx.GetTimeStamp().GetSec();
00560   if(sec < 1115450000) {
00561     return 0.4 * Munits::microsecond; // 400ns into the window, as reported by Peter Shannahan
00562   } else if(sec < 1115977300) {
00563     return -0.7 * Munits::microsecond; // A mistuning introduced May 6 by me
00564   } else {
00565     return 1.5 * Munits::microsecond; // The final tuning introduced May 12 to give a 1.5us window before the spill.
00566   }
00567 
00568   
00569 }
00570 
00571 
00572  void SpillTimeFinder::FindBestRows(DbiResultPtr<SpillTimeND> &table,
00573                                    const VldTimeStamp& time, 
00574                                    const SpillTimeND* &outPrev,
00575                                    const SpillTimeND* &outNext)
00576 {
00582   
00583   //
00584   outPrev = outNext = 0;
00585 
00586   UInt_t low = 0;
00587   UInt_t high= table.GetNumRows();
00588 
00589   if(high <= 0) return; // No rows available. (Query must have failed.)
00590 
00591   // Can't depend on table being sorted, so I'll just brute-force the mother.
00592 
00593   double dt;
00594   double minneg=1e99;
00595   double minpos=1e99;
00596 
00597 
00598   for(UInt_t i=low;i<high;i++) {
00599     const SpillTimeND* row = table.GetRow(i);
00600     dt = GetTimeDifference(row->GetTimeStamp(),time);
00601     if(dt <=0) {
00602       if( (-dt)<minneg ) {
00603         minneg = -dt;
00604         outPrev = row;
00605       }
00606     } else {
00607       if( (dt)<minpos ) {
00608         minpos = -dt;
00609         outNext = row;
00610       }
00611     }
00612   }
00613 
00614   // Throw out any solutions that are just plain stupid.
00615   const VldTimeStamp tooEarly(1975,0,0,0,0,0,0);
00616 
00617   if(outPrev) 
00618         if(outPrev->GetTimeStamp() < tooEarly) outPrev = 0;
00619   if(outNext)
00620         if(outNext->GetTimeStamp() < tooEarly) outNext = 0;
00621 }
00622 
00623 Double_t SpillTimeFinder::GetTimeDifference(const VldTimeStamp& a, 
00624                                             const VldTimeStamp& b)
00625 {
00626   // Return a-b
00627   double d1 = a.GetSec() - b.GetSec();
00628   double d2 = (a.GetNanoSec() - b.GetNanoSec())*1e-9;
00629   return d1+d2;  
00630 }
00631 
00632 
00633 // User calls.
00634 
00635 
00636 
00637 const SpillTimeND& SpillTimeFinder::GetRecentSpill(const VldContext& context)
00638 {
00639   // Return reference to most recent spill record.
00640 
00641   const SpillTimeND* prev;
00642   const SpillTimeND* next;
00643   FindClosestEntries(context,prev,next);
00644   if(prev) return *prev;
00645   else     return kSpillTime_VeryOld;
00646 }
00647 
00648 const SpillTimeND& SpillTimeFinder::GetNextSpill(const VldContext& context)
00649 {
00650   // Return reference to next spill record.
00651 
00652   const SpillTimeND* prev;
00653   const SpillTimeND* next;
00654   FindClosestEntries(context,prev,next);
00655   if(next) return *next;
00656   else     return kSpillTime_Future;
00657 }
00658 
00659 const SpillTimeND& SpillTimeFinder::GetNearestSpill(const VldContext& context)
00660 {
00661   // return referenct to prev or next spill record, whichever is closer.
00662 
00663   const SpillTimeND* prev;
00664   const SpillTimeND* next;
00665   FindClosestEntries(context,prev,next);
00666 
00667   if((next==0)  && (prev==0)) return kSpillTime_VeryOld;
00668 
00669   if(!prev) return *next;
00670   if(!next) return *prev;
00671 
00672   double offset = GetOffset(context,fTask);
00673   VldTimeStamp prevts = prev->GetTimeStamp();
00674   prevts.Add(offset);
00675   VldTimeStamp nextts = next->GetTimeStamp();
00676   nextts.Add(offset);
00677 
00678   double dt1 = GetTimeDifference(prevts,context.GetTimeStamp());
00679   double dt2 = GetTimeDifference(nextts,context.GetTimeStamp());
00680 
00681   if(fabs(dt1)<=fabs(dt2)) {
00682     return *prev;
00683   }
00684   // else
00685   return *next; 
00686 }
00687 
00688 VldTimeStamp SpillTimeFinder::GetTimeOfRecentSpill(const VldContext& cx)
00689 {
00690   VldTimeStamp tspill = GetRecentSpill(cx).GetTimeStamp();
00691   tspill.Add(GetOffset(cx,fTask));
00692   return tspill;
00693 }
00694 
00695 VldTimeStamp SpillTimeFinder::GetTimeOfNextSpill(const VldContext& cx)
00696 { 
00697   VldTimeStamp tspill = GetNextSpill(cx).GetTimeStamp();
00698   tspill.Add(GetOffset(cx,fTask));
00699   return tspill;
00700 }
00701 
00702 VldTimeStamp SpillTimeFinder::GetTimeOfNearestSpill(const VldContext& cx)
00703 {
00704   VldTimeStamp tspill = GetNearestSpill(cx).GetTimeStamp();
00705   tspill.Add(GetOffset(cx,fTask));
00706   return tspill;
00707 }
00708 
00709 

Generated on Mon Feb 15 11:07:38 2010 for loon by  doxygen 1.3.9.1