/** Executes the algorithm
*
*  @throw runtime_error Thrown if algorithm cannot execute
*/
void DiffractionEventCalibrateDetectors::exec() {
  // Try to retrieve optional properties
  const int maxIterations = getProperty("MaxIterations");
  const double peakOpt = getProperty("LocationOfPeakToOptimize");

  // Get the input workspace
  EventWorkspace_sptr inputW = getProperty("InputWorkspace");

  // retrieve the properties
  const std::string rb_params = getProperty("Params");

  // Get some stuff from the input workspace
  Instrument_const_sptr inst = inputW->getInstrument();

  // Build a list of Rectangular Detectors
  std::vector<boost::shared_ptr<RectangularDetector>> detList;
  // --------- Loading only one bank ----------------------------------
  std::string onebank = getProperty("BankName");
  bool doOneBank = (onebank != "");
  for (int i = 0; i < inst->nelements(); i++) {
    boost::shared_ptr<RectangularDetector> det;
    boost::shared_ptr<ICompAssembly> assem;
    boost::shared_ptr<ICompAssembly> assem2;

    det = boost::dynamic_pointer_cast<RectangularDetector>((*inst)[i]);
    if (det) {
      if (det->getName().compare(onebank) == 0)
        detList.push_back(det);
      if (!doOneBank)
        detList.push_back(det);
    } else {
      // Also, look in the first sub-level for RectangularDetectors (e.g. PG3).
      // We are not doing a full recursive search since that will be very long
      // for lots of pixels.
      assem = boost::dynamic_pointer_cast<ICompAssembly>((*inst)[i]);
      if (assem) {
        for (int j = 0; j < assem->nelements(); j++) {
          det = boost::dynamic_pointer_cast<RectangularDetector>((*assem)[j]);
          if (det) {
            if (det->getName().compare(onebank) == 0)
              detList.push_back(det);
            if (!doOneBank)
              detList.push_back(det);

          } else {
            // Also, look in the second sub-level for RectangularDetectors (e.g.
            // PG3).
            // We are not doing a full recursive search since that will be very
            // long for lots of pixels.
            assem2 = boost::dynamic_pointer_cast<ICompAssembly>((*assem)[j]);
            if (assem2) {
              for (int k = 0; k < assem2->nelements(); k++) {
                det = boost::dynamic_pointer_cast<RectangularDetector>(
                    (*assem2)[k]);
                if (det) {
                  if (det->getName().compare(onebank) == 0)
                    detList.push_back(det);
                  if (!doOneBank)
                    detList.push_back(det);
                }
              }
            }
          }
        }
      }
    }
  }

  // set-up minimizer

  std::string inname = getProperty("InputWorkspace");
  std::string outname = inname + "2"; // getProperty("OutputWorkspace");

  IAlgorithm_sptr algS = createChildAlgorithm("SortEvents");
  algS->setProperty("InputWorkspace", inputW);
  algS->setPropertyValue("SortBy", "X Value");
  algS->executeAsChildAlg();

  // Write DetCal File
  std::string filename = getProperty("DetCalFilename");
  std::fstream outfile;
  outfile.open(filename.c_str(), std::ios::out);

  if (detList.size() > 1) {
    outfile << "#\n";
    outfile << "#  Mantid Optimized .DetCal file for SNAP with TWO detector "
               "panels\n";
    outfile << "#  Old Panel, nominal size and distance at -90 degrees.\n";
    outfile << "#  New Panel, nominal size and distance at +90 degrees.\n";
    outfile << "#\n";
    outfile << "# Lengths are in centimeters.\n";
    outfile << "# Base and up give directions of unit vectors for a local\n";
    outfile << "# x,y coordinate system on the face of the detector.\n";
    outfile << "#\n";
    outfile << "# " << DateAndTime::getCurrentTime().toFormattedString("%c")
            << "\n";
    outfile << "#\n";
    outfile << "6         L1     T0_SHIFT\n";
    IComponent_const_sptr source = inst->getSource();
    IComponent_const_sptr sample = inst->getSample();
    outfile << "7  " << source->getDistance(*sample) * 100 << "            0\n";
    outfile << "4 DETNUM  NROWS  NCOLS  WIDTH   HEIGHT   DEPTH   DETD   "
               "CenterX   CenterY   CenterZ    BaseX    BaseY    BaseZ      "
               "UpX      UpY      UpZ\n";
  }

  Progress prog(this, 0.0, 1.0, detList.size());
  for (int det = 0; det < static_cast<int>(detList.size()); det++) {
    std::string par[6];
    par[0] = detList[det]->getName();
    par[1] = inname;
    par[2] = outname;
    std::ostringstream strpeakOpt;
    strpeakOpt << peakOpt;
    par[3] = strpeakOpt.str();
    par[4] = rb_params;

    // --- Create a GroupingWorkspace for this detector name ------
    CPUTimer tim;
    IAlgorithm_sptr alg2 =
        AlgorithmFactory::Instance().create("CreateGroupingWorkspace", 1);
    alg2->initialize();
    alg2->setProperty("InputWorkspace", inputW);
    alg2->setPropertyValue("GroupNames", detList[det]->getName());
    std::string groupWSName = "group_" + detList[det]->getName();
    alg2->setPropertyValue("OutputWorkspace", groupWSName);
    alg2->executeAsChildAlg();
    par[5] = groupWSName;
    std::cout << tim << " to CreateGroupingWorkspace\n";

    const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex;
    gsl_multimin_fminimizer *s = nullptr;
    gsl_vector *ss, *x;
    gsl_multimin_function minex_func;

    // finally do the fitting

    int nopt = 6;
    int iter = 0;
    int status = 0;

    /* Starting point */
    x = gsl_vector_alloc(nopt);
    gsl_vector_set(x, 0, 0.0);
    gsl_vector_set(x, 1, 0.0);
    gsl_vector_set(x, 2, 0.0);
    gsl_vector_set(x, 3, 0.0);
    gsl_vector_set(x, 4, 0.0);
    gsl_vector_set(x, 5, 0.0);

    /* Set initial step sizes to 0.1 */
    ss = gsl_vector_alloc(nopt);
    gsl_vector_set_all(ss, 0.1);

    /* Initialize method and iterate */
    minex_func.n = nopt;
    minex_func.f = &Mantid::Algorithms::gsl_costFunction;
    minex_func.params = &par;

    s = gsl_multimin_fminimizer_alloc(T, nopt);
    gsl_multimin_fminimizer_set(s, &minex_func, x, ss);

    do {
      iter++;
      status = gsl_multimin_fminimizer_iterate(s);

      if (status)
        break;

      double size = gsl_multimin_fminimizer_size(s);
      status = gsl_multimin_test_size(size, 1e-2);

    } while (status == GSL_CONTINUE && iter < maxIterations &&
             s->fval != -0.000);

    // Output summary to log file
    if (s->fval != -0.000)
      movedetector(gsl_vector_get(s->x, 0), gsl_vector_get(s->x, 1),
                   gsl_vector_get(s->x, 2), gsl_vector_get(s->x, 3),
                   gsl_vector_get(s->x, 4), gsl_vector_get(s->x, 5), par[0],
                   getProperty("InputWorkspace"));
    else {
      gsl_vector_set(s->x, 0, 0.0);
      gsl_vector_set(s->x, 1, 0.0);
      gsl_vector_set(s->x, 2, 0.0);
      gsl_vector_set(s->x, 3, 0.0);
      gsl_vector_set(s->x, 4, 0.0);
      gsl_vector_set(s->x, 5, 0.0);
    }

    std::string reportOfDiffractionEventCalibrateDetectors =
        gsl_strerror(status);
    if (s->fval == -0.000)
      reportOfDiffractionEventCalibrateDetectors = "No events";

    g_log.information() << "Detector = " << det << "\n"
                        << "Method used = "
                        << "Simplex"
                        << "\n"
                        << "Iteration = " << iter << "\n"
                        << "Status = "
                        << reportOfDiffractionEventCalibrateDetectors << "\n"
                        << "Minimize PeakLoc-" << peakOpt << " = " << s->fval
                        << "\n";
    // Move in cm for small shifts
    g_log.information() << "Move (X)   = " << gsl_vector_get(s->x, 0) * 0.01
                        << "  \n";
    g_log.information() << "Move (Y)   = " << gsl_vector_get(s->x, 1) * 0.01
                        << "  \n";
    g_log.information() << "Move (Z)   = " << gsl_vector_get(s->x, 2) * 0.01
                        << "  \n";
    g_log.information() << "Rotate (X) = " << gsl_vector_get(s->x, 3) << "  \n";
    g_log.information() << "Rotate (Y) = " << gsl_vector_get(s->x, 4) << "  \n";
    g_log.information() << "Rotate (Z) = " << gsl_vector_get(s->x, 5) << "  \n";

    Kernel::V3D CalCenter =
        V3D(gsl_vector_get(s->x, 0) * 0.01, gsl_vector_get(s->x, 1) * 0.01,
            gsl_vector_get(s->x, 2) * 0.01);
    Kernel::V3D Center = detList[det]->getPos() + CalCenter;
    int pixmax = detList[det]->xpixels() - 1;
    int pixmid = (detList[det]->ypixels() - 1) / 2;
    BoundingBox box;
    detList[det]->getAtXY(pixmax, pixmid)->getBoundingBox(box);
    double baseX = box.xMax();
    double baseY = box.yMax();
    double baseZ = box.zMax();
    Kernel::V3D Base = V3D(baseX, baseY, baseZ) + CalCenter;
    pixmid = (detList[det]->xpixels() - 1) / 2;
    pixmax = detList[det]->ypixels() - 1;
    detList[det]->getAtXY(pixmid, pixmax)->getBoundingBox(box);
    double upX = box.xMax();
    double upY = box.yMax();
    double upZ = box.zMax();
    Kernel::V3D Up = V3D(upX, upY, upZ) + CalCenter;
    Base -= Center;
    Up -= Center;
    // Rotate around x
    baseX = Base[0];
    baseY = Base[1];
    baseZ = Base[2];
    double deg2rad = M_PI / 180.0;
    double angle = gsl_vector_get(s->x, 3) * deg2rad;
    Base = V3D(baseX, baseY * cos(angle) - baseZ * sin(angle),
               baseY * sin(angle) + baseZ * cos(angle));
    upX = Up[0];
    upY = Up[1];
    upZ = Up[2];
    Up = V3D(upX, upY * cos(angle) - upZ * sin(angle),
             upY * sin(angle) + upZ * cos(angle));
    // Rotate around y
    baseX = Base[0];
    baseY = Base[1];
    baseZ = Base[2];
    angle = gsl_vector_get(s->x, 4) * deg2rad;
    Base = V3D(baseZ * sin(angle) + baseX * cos(angle), baseY,
               baseZ * cos(angle) - baseX * sin(angle));
    upX = Up[0];
    upY = Up[1];
    upZ = Up[2];
    Up = V3D(upZ * cos(angle) - upX * sin(angle), upY,
             upZ * sin(angle) + upX * cos(angle));
    // Rotate around z
    baseX = Base[0];
    baseY = Base[1];
    baseZ = Base[2];
    angle = gsl_vector_get(s->x, 5) * deg2rad;
    Base = V3D(baseX * cos(angle) - baseY * sin(angle),
               baseX * sin(angle) + baseY * cos(angle), baseZ);
    upX = Up[0];
    upY = Up[1];
    upZ = Up[2];
    Up = V3D(upX * cos(angle) - upY * sin(angle),
             upX * sin(angle) + upY * cos(angle), upZ);
    Base.normalize();
    Up.normalize();
    Center *= 100.0;
    // << det+1  << "  "
    outfile << "5  " << detList[det]->getName().substr(4) << "  "
            << detList[det]->xpixels() << "  " << detList[det]->ypixels()
            << "  " << 100.0 * detList[det]->xsize() << "  "
            << 100.0 * detList[det]->ysize() << "  "
            << "0.2000"
            << "  " << Center.norm() << "  ";
    Center.write(outfile);
    outfile << "  ";
    Base.write(outfile);
    outfile << "  ";
    Up.write(outfile);
    outfile << "\n";

    // clean up dynamically allocated gsl stuff
    gsl_vector_free(x);
    gsl_vector_free(ss);
    gsl_multimin_fminimizer_free(s);

    // Remove the now-unneeded grouping workspace
    AnalysisDataService::Instance().remove(groupWSName);
    prog.report(detList[det]->getName());
  }

  // Closing
  outfile.close();
}
示例#2
0
void SofQWCentre::exec() {
  using namespace Geometry;

  MatrixWorkspace_const_sptr inputWorkspace = getProperty("InputWorkspace");

  // Do the full check for common binning
  if (!WorkspaceHelpers::commonBoundaries(*inputWorkspace)) {
    g_log.error(
        "The input workspace must have common binning across all spectra");
    throw std::invalid_argument(
        "The input workspace must have common binning across all spectra");
  }

  std::vector<double> verticalAxis;
  MatrixWorkspace_sptr outputWorkspace = setUpOutputWorkspace(
      inputWorkspace, getProperty("QAxisBinning"), verticalAxis);
  setProperty("OutputWorkspace", outputWorkspace);

  // Holds the spectrum-detector mapping
  std::vector<specnum_t> specNumberMapping;
  std::vector<detid_t> detIDMapping;

  m_EmodeProperties.initCachedValues(*inputWorkspace, this);
  int emode = m_EmodeProperties.m_emode;

  // Get a pointer to the instrument contained in the workspace
  Instrument_const_sptr instrument = inputWorkspace->getInstrument();

  // Get the distance between the source and the sample (assume in metres)
  IComponent_const_sptr source = instrument->getSource();
  IComponent_const_sptr sample = instrument->getSample();
  V3D beamDir = sample->getPos() - source->getPos();
  beamDir.normalize();

  try {
    double l1 = source->getDistance(*sample);
    g_log.debug() << "Source-sample distance: " << l1 << '\n';
  } catch (Exception::NotFoundError &) {
    g_log.error("Unable to calculate source-sample distance");
    throw Exception::InstrumentDefinitionError(
        "Unable to calculate source-sample distance",
        inputWorkspace->getTitle());
  }

  // Conversion constant for E->k. k(A^-1) = sqrt(energyToK*E(meV))
  const double energyToK = 8.0 * M_PI * M_PI * PhysicalConstants::NeutronMass *
                           PhysicalConstants::meV * 1e-20 /
                           (PhysicalConstants::h * PhysicalConstants::h);

  // Loop over input workspace bins, reassigning data to correct bin in output
  // qw workspace
  const size_t numHists = inputWorkspace->getNumberHistograms();
  const size_t numBins = inputWorkspace->blocksize();
  Progress prog(this, 0.0, 1.0, numHists);
  for (int64_t i = 0; i < int64_t(numHists); ++i) {
    try {
      // Now get the detector object for this histogram
      IDetector_const_sptr spectrumDet = inputWorkspace->getDetector(i);
      if (spectrumDet->isMonitor())
        continue;

      const double efixed = m_EmodeProperties.getEFixed(*spectrumDet);

      // For inelastic scattering the simple relationship q=4*pi*sinTheta/lambda
      // does not hold. In order to
      // be completely general we must calculate the momentum transfer by
      // calculating the incident and final
      // wave vectors and then use |q| = sqrt[(ki - kf)*(ki - kf)]
      DetectorGroup_const_sptr detGroup =
          boost::dynamic_pointer_cast<const DetectorGroup>(spectrumDet);
      std::vector<IDetector_const_sptr> detectors;
      if (detGroup) {
        detectors = detGroup->getDetectors();
      } else {
        detectors.push_back(spectrumDet);
      }

      const size_t numDets = detectors.size();
      // cache to reduce number of static casts
      const double numDets_d = static_cast<double>(numDets);
      const auto &Y = inputWorkspace->y(i);
      const auto &E = inputWorkspace->e(i);
      const auto &X = inputWorkspace->x(i);

      // Loop over the detectors and for each bin calculate Q
      for (size_t idet = 0; idet < numDets; ++idet) {
        IDetector_const_sptr det = detectors[idet];
        // Calculate kf vector direction and then Q for each energy bin
        V3D scatterDir = (det->getPos() - sample->getPos());
        scatterDir.normalize();
        for (size_t j = 0; j < numBins; ++j) {
          const double deltaE = 0.5 * (X[j] + X[j + 1]);
          // Compute ki and kf wave vectors and therefore q = ki - kf
          double ei(0.0), ef(0.0);
          if (emode == 1) {
            ei = efixed;
            ef = efixed - deltaE;
            if (ef < 0) {
              std::string mess =
                  "Energy transfer requested in Direct mode exceeds incident "
                  "energy.\n Found for det ID: " +
                  std::to_string(idet) + " bin No " + std::to_string(j) +
                  " with Ei=" + boost::lexical_cast<std::string>(efixed) +
                  " and energy transfer: " +
                  boost::lexical_cast<std::string>(deltaE);
              throw std::runtime_error(mess);
            }
          } else {
            ei = efixed + deltaE;
            ef = efixed;
            if (ef < 0) {
              std::string mess =
                  "Incident energy of a neutron is negative. Are you trying to "
                  "process Direct data in Indirect mode?\n Found for det ID: " +
                  std::to_string(idet) + " bin No " + std::to_string(j) +
                  " with efied=" + boost::lexical_cast<std::string>(efixed) +
                  " and energy transfer: " +
                  boost::lexical_cast<std::string>(deltaE);
              throw std::runtime_error(mess);
            }
          }

          if (ei < 0)
            throw std::runtime_error(
                "Negative incident energy. Check binning.");

          const V3D ki = beamDir * sqrt(energyToK * ei);
          const V3D kf = scatterDir * (sqrt(energyToK * (ef)));
          const double q = (ki - kf).norm();

          // Test whether it's in range of the Q axis
          if (q < verticalAxis.front() || q > verticalAxis.back())
            continue;
          // Find which q bin this point lies in
          const MantidVec::difference_type qIndex =
              std::upper_bound(verticalAxis.begin(), verticalAxis.end(), q) -
              verticalAxis.begin() - 1;

          // Add this spectra-detector pair to the mapping
          specNumberMapping.push_back(
              outputWorkspace->getSpectrum(qIndex).getSpectrumNo());
          detIDMapping.push_back(det->getID());

          // And add the data and it's error to that bin, taking into account
          // the number of detectors contributing to this bin
          outputWorkspace->mutableY(qIndex)[j] += Y[j] / numDets_d;
          // Standard error on the average
          outputWorkspace->mutableE(qIndex)[j] =
              sqrt((pow(outputWorkspace->e(qIndex)[j], 2) + pow(E[j], 2)) /
                   numDets_d);
        }
      }

    } catch (Exception::NotFoundError &) {
      // Get to here if exception thrown when calculating distance to detector
      // Presumably, if we get to here the spectrum will be all zeroes anyway
      // (from conversion to E)
      continue;
    }
    prog.report();
  }

  // If the input workspace was a distribution, need to divide by q bin width
  if (inputWorkspace->isDistribution())
    this->makeDistribution(outputWorkspace, verticalAxis);

  // Set the output spectrum-detector mapping
  SpectrumDetectorMapping outputDetectorMap(specNumberMapping, detIDMapping);
  outputWorkspace->updateSpectraUsing(outputDetectorMap);

  // Replace any NaNs in outputWorkspace with zeroes
  if (this->getProperty("ReplaceNaNs")) {
    auto replaceNans = this->createChildAlgorithm("ReplaceSpecialValues");
    replaceNans->setChild(true);
    replaceNans->initialize();
    replaceNans->setProperty("InputWorkspace", outputWorkspace);
    replaceNans->setProperty("OutputWorkspace", outputWorkspace);
    replaceNans->setProperty("NaNValue", 0.0);
    replaceNans->setProperty("InfinityValue", 0.0);
    replaceNans->setProperty("BigNumberThreshold", DBL_MAX);
    replaceNans->execute();
  }
}
/**
 * Cache frequently accessed values
 * @param instrument : The instrument for this run
 * @param detID : The det ID for this observation
 */
void CachedExperimentInfo::initCaches(
    const Geometry::Instrument_const_sptr &instrument, const detid_t detID) {
  // Throws if detector does not exist
  // Takes into account possible detector mapping
  IDetector_const_sptr det = m_exptInfo.getDetectorByID(detID);

  // Instrument distances
  boost::shared_ptr<const ReferenceFrame> refFrame =
      instrument->getReferenceFrame();
  m_beam = refFrame->pointingAlongBeam();
  m_up = refFrame->pointingUp();
  m_horiz = refFrame->pointingHorizontal();

  IComponent_const_sptr source = instrument->getSource();
  IComponent_const_sptr sample = instrument->getSample();
  IComponent_const_sptr aperture =
      instrument->getComponentByName("aperture", 1);
  if (!aperture) {
    throw std::invalid_argument(
        "No component named \"aperture\" found in instrument.");
  }
  IObjComponent_const_sptr firstChopper = instrument->getChopperPoint(0);
  const Kernel::V3D samplePos = sample->getPos();
  const Kernel::V3D beamDir = samplePos - source->getPos();

  // Cache
  m_twoTheta = det->getTwoTheta(samplePos, beamDir);
  m_phi = det->getPhi();
  m_modToChop = firstChopper->getDistance(*source);
  m_apertureToChop = firstChopper->getDistance(*aperture);
  m_chopToSample = sample->getDistance(*firstChopper);
  m_sampleToDet = det->getDistance(*sample);

  // Aperture
  Geometry::BoundingBox apertureBox;
  aperture->getBoundingBox(apertureBox);
  if (apertureBox.isNull()) {
    throw std::invalid_argument("CachedExperimentInfo::initCaches - Aperture "
                                "has no bounding box, cannot sample from it");
  }
  m_apertureSize.first = apertureBox.maxPoint()[0] - apertureBox.minPoint()[0];
  m_apertureSize.second = apertureBox.maxPoint()[1] - apertureBox.minPoint()[1];

  // Sample volume
  const API::Sample &sampleDescription = m_exptInfo.sample();
  const Geometry::Object &shape = sampleDescription.getShape();
  m_sampleWidths = shape.getBoundingBox().width();

  // Detector volume
  // Make sure it encompasses all possible detectors
  det->getBoundingBox(m_detBox);
  if (m_detBox.isNull()) {
    throw std::invalid_argument("CachedExperimentInfo::initCaches - Detector "
                                "has no bounding box, cannot sample from it. "
                                "ID:" +
                                boost::lexical_cast<std::string>(det->getID()));
  }

  const double rad2deg = 180. / M_PI;
  const double thetaInDegs = twoTheta() * rad2deg;
  const double phiInDegs = phi() * rad2deg;

  m_gonimeter = new Goniometer;
  m_gonimeter->makeUniversalGoniometer();
  m_gonimeter->setRotationAngle("phi", thetaInDegs);
  m_gonimeter->setRotationAngle("chi", phiInDegs);
  m_sampleToDetMatrix =
      m_exptInfo.sample().getOrientedLattice().getU() * m_gonimeter->getR();

  // EFixed
  m_efixed = m_exptInfo.getEFixed(det);
}
void ModeratorTzero::execEvent(const std::string &emode) {
  g_log.information("Processing event workspace");

  const MatrixWorkspace_const_sptr matrixInputWS =
      getProperty("InputWorkspace");
  EventWorkspace_const_sptr inputWS =
      boost::dynamic_pointer_cast<const EventWorkspace>(matrixInputWS);

  // generate the output workspace pointer
  const size_t numHists = static_cast<size_t>(inputWS->getNumberHistograms());
  Mantid::API::MatrixWorkspace_sptr matrixOutputWS =
      getProperty("OutputWorkspace");
  EventWorkspace_sptr outputWS;
  if (matrixOutputWS == matrixInputWS) {
    outputWS = boost::dynamic_pointer_cast<EventWorkspace>(matrixOutputWS);
  } else {
    // Make a brand new EventWorkspace
    outputWS = boost::dynamic_pointer_cast<EventWorkspace>(
        WorkspaceFactory::Instance().create("EventWorkspace", numHists, 2, 1));
    // Copy geometry over.
    WorkspaceFactory::Instance().initializeFromParent(inputWS, outputWS, false);
    // You need to copy over the data as well.
    outputWS->copyDataFrom((*inputWS));
    // Cast to the matrixOutputWS and save it
    matrixOutputWS = boost::dynamic_pointer_cast<MatrixWorkspace>(outputWS);
    setProperty("OutputWorkspace", matrixOutputWS);
  }

  // Get pointers to sample and source
  IComponent_const_sptr source = m_instrument->getSource();
  IComponent_const_sptr sample = m_instrument->getSample();
  double Lss = source->getDistance(*sample); // distance from source to sample

  // calculate tof shift once for all neutrons if emode==Direct
  double t0_direct(-1);
  if (emode == "Direct") {
    Kernel::Property *eiprop = inputWS->run().getProperty("Ei");
    double Ei = boost::lexical_cast<double>(eiprop->value());
    mu::Parser parser;
    parser.DefineVar("incidentEnergy", &Ei); // associate E1 to this parser
    parser.SetExpr(m_formula);
    t0_direct = parser.Eval();
  }

  // Loop over the spectra
  Progress prog(this, 0.0, 1.0, numHists); // report progress of algorithm
  PARALLEL_FOR1(outputWS)
  for (int i = 0; i < static_cast<int>(numHists); ++i) {
    PARALLEL_START_INTERUPT_REGION
    size_t wsIndex = static_cast<size_t>(i);
    EventList &evlist = outputWS->getEventList(wsIndex);
    if (evlist.getNumberEvents() > 0) // don't bother with empty lists
    {
      IDetector_const_sptr det;
      double L1(Lss); // distance from source to sample
      double L2(-1);  // distance from sample to detector

      try {
        det = inputWS->getDetector(i);
        if (det->isMonitor()) {
          // redefine the sample as the monitor
          L1 = source->getDistance(*det);
          L2 = 0;
        } else {
          L2 = sample->getDistance(*det);
        }
      } catch (Exception::NotFoundError &) {
        g_log.error() << "Unable to calculate distances to/from detector" << i
                      << std::endl;
      }

      if (L2 >= 0) {
        // One parser for each parallel processor needed (except Edirect mode)
        double E1;
        mu::Parser parser;
        parser.DefineVar("incidentEnergy", &E1); // associate E1 to this parser
        parser.SetExpr(m_formula);

        // fast neutrons are shifted by min_t0_next, irrespective of tof
        double v1_max = L1 / m_t1min;
        E1 = m_convfactor * v1_max * v1_max;
        double min_t0_next = parser.Eval();

        if (emode == "Indirect") {
          double t2(-1.0); // time from sample to detector. (-1) signals error
          if (det->isMonitor()) {
            t2 = 0.0;
          } else {
            static const double convFact =
                1.0e-6 * sqrt(2 * PhysicalConstants::meV /
                              PhysicalConstants::NeutronMass);
            std::vector<double> wsProp = det->getNumberParameter("Efixed");
            if (!wsProp.empty()) {
              double E2 = wsProp.at(0);        //[E2]=meV
              double v2 = convFact * sqrt(E2); //[v2]=meter/microsec
              t2 = L2 / v2;
            } else {
              // t2 is kept to -1 if no Efixed is found
              g_log.debug() << "Efixed not found for detector " << i
                            << std::endl;
            }
          }
          if (t2 >= 0) // t2 < 0 when no detector info is available
          {
            // fix the histogram bins
            MantidVec &x = evlist.dataX();
            for (double &tof : x) {
              if (tof < m_t1min + t2)
                tof -= min_t0_next;
              else
                tof -= CalculateT0indirect(tof, L1, t2, E1, parser);
            }

            MantidVec tofs = evlist.getTofs();
            for (double &tof : tofs) {
              if (tof < m_t1min + t2)
                tof -= min_t0_next;
              else
                tof -= CalculateT0indirect(tof, L1, t2, E1, parser);
            }
            evlist.setTofs(tofs);
            evlist.setSortOrder(Mantid::DataObjects::EventSortType::UNSORTED);
          } // end of if( t2>= 0)
        }   // end of if(emode=="Indirect")
        else if (emode == "Elastic") {
          // Apply t0 correction to histogram bins
          MantidVec &x = evlist.dataX();
          for (double &tof : x) {
            if (tof < m_t1min * (L1 + L2) / L1)
              tof -= min_t0_next;
            else
              tof -= CalculateT0elastic(tof, L1 + L2, E1, parser);
          }

          MantidVec tofs = evlist.getTofs();
          for (double &tof : tofs) {
            // add a [-0.1,0.1] microsecond noise to avoid artifacts
            // resulting from original tof data
            if (tof < m_t1min * (L1 + L2) / L1)
              tof -= min_t0_next;
            else
              tof -= CalculateT0elastic(tof, L1 + L2, E1, parser);
          }
          evlist.setTofs(tofs);
          evlist.setSortOrder(Mantid::DataObjects::EventSortType::UNSORTED);

          MantidVec tofs_b = evlist.getTofs();
          MantidVec xarray = evlist.readX();
        } // end of else if(emode=="Elastic")
        else if (emode == "Direct") {
          // fix the histogram bins
          MantidVec &x = evlist.dataX();
          for (double &tof : x) {
            tof -= t0_direct;
          }

          MantidVec tofs = evlist.getTofs();
          for (double &tof : tofs) {
            tof -= t0_direct;
          }
          evlist.setTofs(tofs);
          evlist.setSortOrder(Mantid::DataObjects::EventSortType::UNSORTED);
        } // end of else if(emode=="Direct")
      }   // end of if(L2 >= 0)
    }     // end of if (evlist.getNumberEvents() > 0)
    prog.report();
    PARALLEL_END_INTERUPT_REGION
  } // end of for (int i = 0; i < static_cast<int>(numHists); ++i)
  PARALLEL_CHECK_INTERUPT_REGION
  outputWS->clearMRU(); // Clears the Most Recent Used lists */
} // end of void ModeratorTzero::execEvent()