GooFit  v2.1.3
MappedPdf.cu
Go to the documentation of this file.
1 #include <goofit/PDFs/combine/MappedPdf.h>
2 
3 namespace GooFit {
4 
5 __device__ fptype device_Mapped(fptype *evt, fptype *p, unsigned int *indices) {
6  // Structure : nP mapFunctionIndex mapParamIndex functionIndex1 parameterIndex1 functionIndex2 parameterIndex2 ...
7 
8  // Find mapping between event variables and function to evaluate
9  unsigned int mapFunction = RO_CACHE(indices[1]);
10  // This is an index into the MappedPdf's list of functions
11  // int targetFunction = (int) floor(0.5 +
12  // (*(reinterpret_cast<device_function_ptr>(device_function_table[mapFunction])))(evt, p, paramIndices +
13  // indices[2]));
14  auto targetFunction = static_cast<int>(floor(0.5 + callFunction(evt, mapFunction, RO_CACHE(indices[2]))));
15 
16  targetFunction *= 2; // Because there are two pieces of information about each function
17  targetFunction += 3; // Because first function information begins at index 3
18 
19  // fptype ret = (*(reinterpret_cast<device_function_ptr>(device_function_table[indices[targetFunction]])))(evt, p,
20  // paramIndices + indices[targetFunction + 1]);
21  fptype ret = callFunction(evt, RO_CACHE(indices[targetFunction]), RO_CACHE(indices[targetFunction + 1]));
22  ret *= normalisationFactors[RO_CACHE(indices[targetFunction + 1])];
23  // if (gpuDebug & 1)
24  // if ((gpuDebug & 1) && (0 == BLOCKIDX) && (0 == THREADIDX))
25  // printf("[%i, %i] Mapped: %i (%f %f %f %f) %f\n", BLOCKIDX, THREADIDX, targetFunction, evt[0], evt[1], evt[2],
26  // evt[3], ret);
27  return ret;
28 }
29 
30 __device__ device_function_ptr ptr_to_Mapped = device_Mapped;
31 
32 __host__ MappedPdf::MappedPdf(std::string n, GooPdf *m, std::vector<GooPdf *> &t)
33  : GooPdf(n) {
34  components.push_back(m);
35  std::vector<unsigned int> pindices;
36  pindices.push_back(m->getFunctionIndex());
37  pindices.push_back(m->getParameterIndex());
38 
39  std::set<int> functionIndicesUsed;
40 
41  for(GooPdf *f : t) {
42  components.push_back(f);
43  pindices.push_back(f->getFunctionIndex());
44  pindices.push_back(f->getParameterIndex());
45  functionIndicesUsed.insert(f->getFunctionIndex());
46  }
47 
48  if(functionIndicesUsed.size() > 1) {
49  std::cout << "Warning: More than one function type given to MappedPdf " << getName()
50  << " constructor. This may slow execution by causing sequential evaluations.\n";
51  }
52 
53  observables = getObservables();
54  GET_FUNCTION_ADDR(ptr_to_Mapped);
55  initialize(pindices);
56 }
57 
58 __host__ fptype MappedPdf::normalize() const {
59  // std::cout << "Normalising MappedPdf " << getName() << std::endl;
60  fptype ret = 0;
61 
62  for(unsigned int i = 1; i < components.size(); ++i) { // No need to normalize mapping function.
63  fptype curr = components[i]->normalize();
64  ret += curr;
65  }
66 
67  host_normalisation[parameters] = 1.0;
68  return ret;
69 }
70 } // namespace GooFit