CMS 3D CMS Logo

Functions
HcalPulseShapes.cc File Reference
#include "CalibCalorimetry/HcalAlgos/interface/HcalPulseShapes.h"
#include "FWCore/Utilities/interface/Exception.h"
#include "FWCore/Framework/interface/ConsumesCollector.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "CLHEP/Random/RandFlat.h"
#include <cmath>
#include <iostream>
#include <fstream>
#include "TMath.h"

Go to the source code of this file.

Functions

double gexp (double t, double A, double c, double t0, double s)
 
double onePulse (double t, double A, double sigma, double theta, double m)
 

Function Documentation

◆ gexp()

double gexp ( double  t,
double  A,
double  c,
double  t0,
double  s 
)
inline

Definition at line 559 of file HcalPulseShapes.cc.

References A, c, JetChargeProducer_cfi::exp, alignCSCRings::s, mathSSE::sqrt(), submitPVValidationJobs::t, and FrontierCondition_GT_autoExpress_cfi::t0.

Referenced by HcalPulseShapes::analyticPulseShapeSiPMHO().

559  {
560  static double const root2(sqrt(2));
561  return -A * 0.5 * exp(c * t + 0.5 * c * c * s * s - c * s) * (erf(-0.5 * root2 / s * (t - t0 + c * s * s)) - 1);
562 }
T sqrt(T t)
Definition: SSEVec.h:19
Definition: APVGainStruct.h:7

◆ onePulse()

double onePulse ( double  t,
double  A,
double  sigma,
double  theta,
double  m 
)
inline

Definition at line 564 of file HcalPulseShapes.cc.

References A, visualization-live-secondInstance_cfg::m, submitPVValidationJobs::t, and theta().

Referenced by HcalPulseShapes::analyticPulseShapeSiPMHE().

564  {
565  return (t < theta) ? 0 : A * TMath::LogNormal(t, sigma, theta, m);
566 }
Definition: APVGainStruct.h:7
Geom::Theta< T > theta() const