DDASToys for NSCLDAQ  6.2-000
Macros | Functions
cudafit_analytic.cu File Reference

Provide trace fitting using the libucdafit library. More...

#include "cudafit_analytic.cuh"
#include <limits>
#include <ctime>
#include <iostream>
#include <stdexcept>
#include <string>
#include <float.h>
#include <DE_Optimizer.h>
#include <PSO_Optimizer.h>
#include "fit_extensions.h"
#include "functions_analytic.h"
#include "reductions.cu"
Include dependency graph for cudafit_analytic.cu:

Functions

__host__ __device__ float logistic (float A, float k, float x1, float x)
 Evaluate a logistic function for the specified parameters and point. More...
 
__host__ __device__ float decay (float A, float k, float x1, float x)
 Exponential decay function. More...
 
__host__ __device__ float singlePulse (float A1, float k1, float k2, float x1, float C, float x)
 Evaluate the value of a single pulse in accordance with our canonical functional form.
More...
 
__host__ __device__ float doublePulse (float A1, float k1, float k2, float x1, float A2, float k3, float k4, float x2, float C, float x)
 Evaluate the canonical form of a double pulse.
More...
 
__host__ __device__ float chiFitness1 (const float *pParams, float x, float y, float wt)
 Computes chi-square fitness for one point in one solution given that d_fitness has pulled out what we need. This fitness is for a single-pulse fit. More...
 
__host__ __device__ float chiFitness2 (const float *pParams, float x, float y, float wt)
 Computes chi-square fitness contribution for one point in one solution given that our caller has pulled out what we need. More...
 
__global__ void d_fitness1 (const float *pSolutions, float *pFitnesses, int nParams, int nSol, unsigned short *pXcoords, unsigned short *pYcoords, float *pWeights, int nPoints)
 This function lives in the GPU and: More...
 
__global__ void d_fitness2 (const float *pSolutions, float *pFitnesses, int nParams, int nSol, unsigned short *pXcoords, unsigned short *pYcoords, float *pWeights, int nPoints)
 This function lives in the GPU and: More...
 
void h_fitSingle (const CudaOptimize::SolutionSet *solutions, CudaOptimize::FitnessSet *fitnesses, dim3 grid, dim3 block)
 Invokes the kernel that produces the fitness measure. The fitness is computed in the GPU and is the chi square. More...
 
void h_fitDouble (const CudaOptimize::SolutionSet *solutions, CudaOptimize::FitnessSet *fitnesses, dim3 grid, dim3 block)
 Host part to setup computation of the fitnesses across the swarm for our fits for a double pulse. Really this just sets up the kernel call for fitness2 which does the rest. More...
 

Detailed Description

Provide trace fitting using the libucdafit library.

Author
Ron Foxfox@n.nosp@m.scl..nosp@m.msu.e.nosp@m.du
Note
We provide call compatible interfaces with lmfit1 and lmfit2.
This fit will not thread due to libcudaoptimize's need for us to have global data for the device pointers to the trace.

Function Documentation

◆ chiFitness1()

__host__ __device__ float chiFitness1 ( const float *  pParams,
float  x,
float  y,
float  wt 
)

Computes chi-square fitness for one point in one solution given that d_fitness has pulled out what we need. This fitness is for a single-pulse fit.

The choice of weight determines whether the returned value is Neyman's (weight by variance estimated from data wt = 1/y) or Pearson's (weight by variance estimated from fit wt = 1/fit).

Parameters
pParamsPointer to this solutions parameters.
xx-coordinate.
yy-coordinate.
wtWeight for this coordinate.
Returns
Chi-square between fit function and trace data.

◆ chiFitness2()

__host__ __device__ float chiFitness2 ( const float *  pParams,
float  x,
float  y,
float  wt 
)

Computes chi-square fitness contribution for one point in one solution given that our caller has pulled out what we need.

The choice of weight determines whether the returned value is Neyman's (weight by variance estimated from data wt = 1/y) or Pearson's (weight by variance estimated from fit wt = 1/fit).

Parameters
pParamsPointer to this solutions parameters.
xx-coordinate.
yy-coordinate.
wtWeight for this coordinate.
Returns
Chi-square between fit function and trace data.

◆ d_fitness1()

__global__ void d_fitness1 ( const float *  pSolutions,
float *  pFitnesses,
int  nParams,
int  nSol,
unsigned short *  pXcoords,
unsigned short *  pYcoords,
float *  pWeights,
int  nPoints 
)

This function lives in the GPU and:

  • Computes the chi-square contribution for a single point for a single solution in the swarm for a single pulse with an offset.
  • Uses reduceToSum to sum the chisquare contributions over the entire trace. The result is put into the fitness value for our solution.
Parameters
pSolutionsPointer to solutions array in the GPU.
pFitnessesPointer to the array of fitnesse for all solutions in the swarm.
nParamsNumber of parameters in the fit (should be 5).
nSolNumber of solutions in the swarm.
pXcoordsTrace x-coordinates array.
pYcoordsTrace y-coordinates array.
pWeightsy weights to apply.
nPointsNumber of points in the trace.

◆ d_fitness2()

__global__ void d_fitness2 ( const float *  pSolutions,
float *  pFitnesses,
int  nParams,
int  nSol,
unsigned short *  pXcoords,
unsigned short *  pYcoords,
float *  pWeights,
int  nPoints 
)

This function lives in the GPU and:

  • Computes the chi-square contribution for a single point for a single solution in the swarm for a single pulse with an offset.
  • Uses reduceToSum to sum the chisquare contributions over the entire trace. The result is put into the fitness value for our solution.
Parameters
pSolutionsPointer to solutions array in the GPU.
pFitnessesPointer to the array of fitnesse for all solutions in the swarm.
nParamsNumber of parameters in the fit (should be 5).
nSolNumber of solutions in the swarm.
pXcoordsTrace x-coordinates array.
pYcoordsTrace y-coordinates array.
pWeightsy weights to apply.
nPointsNumber of points in the trace.

◆ decay()

__host__ __device__ float decay ( float  A,
float  k,
float  x1,
float  x 
)

Exponential decay function.

Signals from detectors usually have a falling shape that approximates an exponential. This function evaluates this decay at some point.

Parameters
AAmplitude of the signal
kDecay time factor f the signal.
x1Position of the pulse.
xWhere to evaluate the signal.
Returns
Exponential decay evaluated at x.

◆ doublePulse()

__host__ __device__ float doublePulse ( float  A1,
float  k1,
float  k2,
float  x1,
float  A2,
float  k3,
float  k4,
float  x2,
float  C,
float  x 
)

Evaluate the canonical form of a double pulse.

This is done by summing two single pulses. The constant term is thrown into the first pulse. The second pulse gets a constant term of 0.

Parameters
A1Amplitude of the first pulse.
k1Steepness of first pulse rise.
k2Decay time of the first pulse.
x1Position of the first pulse.
A2Amplitude of the second pulse.
k3Steepness of second pulse rise.
k4Decay time of second pulse.
x2Position of second pulse.
CConstant offset the pulses sit on.
xPosition at which to evaluate the pulse.
Returns
Value of the double-pulse function evaluated at x.

◆ h_fitDouble()

void h_fitDouble ( const CudaOptimize::SolutionSet *  solutions,
CudaOptimize::FitnessSet *  fitnesses,
dim3  grid,
dim3  block 
)

Host part to setup computation of the fitnesses across the swarm for our fits for a double pulse. Really this just sets up the kernel call for fitness2 which does the rest.

Parameters
solutionsPointer to the current Solution set.
fitnessesPointer to the current fitness set
gridComputaional grid geometry.
blockShapes of the blocks within the grid.

◆ h_fitSingle()

void h_fitSingle ( const CudaOptimize::SolutionSet *  solutions,
CudaOptimize::FitnessSet *  fitnesses,
dim3  grid,
dim3  block 
)

Invokes the kernel that produces the fitness measure. The fitness is computed in the GPU and is the chi square.

Parameters
solutionsPointer to the Cuda solution set.
fitnessesPointer to the current fitness set
gridComputational grid being used.
blockShapes of blocks within the grid.

◆ logistic()

__host__ __device__ float logistic ( float  A,
float  k,
float  x1,
float  x 
)

Evaluate a logistic function for the specified parameters and point.

A logistic function is a function with a sigmoidal shape. We use it to fit the rising edge of signals DDAS digitizes from detectors. See e.g. https://en.wikipedia.org/wiki/Logistic_function for a discussion of this function.

Parameters
AAmplitude of the signal.
kSteepness of the signal (related to the rise time).
x1Mid point of the rise of the sigmoid.
xLocation at which to evaluate the function.
Returns
Logistic function evaluated at x.

◆ singlePulse()

__host__ __device__ float singlePulse ( float  A1,
float  k1,
float  k2,
float  x1,
float  C,
float  x 
)

Evaluate the value of a single pulse in accordance with our canonical functional form.

The form is a sigmoid rise with an exponential decay that sits on top of a constant offset. The exponential decay is turned on with switchOn() above when x > the rise point of the sigmoid.

Parameters
A1Pulse amplitiude
k1Sigmoid rise steepness.
k2Exponential decay time constant.
x1Sigmoid position.
CConstant offset.
xPosition at which to evaluate this function.
Returns
Single pulse evaluated at x.