Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
60 changes: 42 additions & 18 deletions ALICE3/Core/FastTracker.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -17,7 +17,11 @@
#include "TMatrixD.h"
#include "TMatrixDSymEigen.h"
#include "TRandom.h"
#include <TEnv.h>
#include <THashList.h>
#include <TObject.h>

#include <fstream>
#include <string>
#include <vector>

Expand Down Expand Up @@ -220,8 +224,8 @@
for (int i = 0; i < kNPassiveBound; i++) {
AddLayer(Form("tpc_boundary%d", i), rBoundary[i], zLength, radLBoundary[i], xrhoBoundary[i], 0); // dummy errors
}
for (Int_t k = 0; k < tpcRows; k++) {

Check failure on line 227 in ALICE3/Core/FastTracker.cxx

View workflow job for this annotation

GitHub Actions / O2 linter

[root/entity]

Replace ROOT entities with equivalents from standard C++ or from O2.
Float_t rowRadius = 0;

Check failure on line 228 in ALICE3/Core/FastTracker.cxx

View workflow job for this annotation

GitHub Actions / O2 linter

[root/entity]

Replace ROOT entities with equivalents from standard C++ or from O2.
if (k < innerRows)
rowRadius = rowOneRadius + k * tpcInnerRadialPitch;
else if (k >= innerRows && k < (innerRows + middleRows))
Expand Down Expand Up @@ -260,27 +264,47 @@
return;
}

std::ifstream infile(filename);
if (!infile.is_open()) {
LOG(fatal) << "Could not open detector configuration file: " << filename;
return;
}

std::string line;
while (std::getline(infile, line)) {
if (line.empty() || line[0] == '#')
continue; // skip comments and empty lines
std::istringstream iss(line);
std::string name;
float r, z, x0, xrho, resRPhi, resZ, eff;
int type;
if (!(iss >> name >> r >> z >> x0 >> xrho >> resRPhi >> resZ >> eff >> type)) {
LOG(error) << "Malformed line in detector config: " << line;
TEnv env(filename.c_str());
THashList* table = env.GetTable();
std::vector<std::string> layers;
for (int i = 0; i < table->GetEntries(); ++i) {
const std::string key = table->At(i)->GetName();
// key should contain exactly one dot
if (key.find('.') == std::string::npos || key.find('.') != key.rfind('.')) {
LOG(fatal) << "Key " << key << " does not contain exactly one dot";
continue;
}
AddLayer(name.c_str(), r, z, x0, xrho, resRPhi, resZ, eff, type);
const std::string firstPart = key.substr(0, key.find('.'));
if (std::find(layers.begin(), layers.end(), firstPart) == layers.end()) {
layers.push_back(firstPart);
}
}
// env.Print();
// Layers
for (const auto& layer : layers) {
LOG(info) << " Reading layer " << layer;

auto getKey = [&layer, &env](const std::string& name) {
std::string key = layer + "." + name;
if (!env.Defined(key.c_str())) {
LOG(warning) << "Key " << key << " not defined in configuration file, getting the default value";
}
LOG(info) << " Getting key " << key;
return key;
};
const float r = env.GetValue(getKey("r").c_str(), -1.0f);
LOG(info) << " Layer " << layer << " has radius " << r;
const float z = env.GetValue(getKey("z").c_str(), -1.0f);
const float x0 = env.GetValue(getKey("x0").c_str(), 0.0f);
const float xrho = env.GetValue(getKey("xrho").c_str(), 0.0f);
const float resRPhi = env.GetValue(getKey("resRPhi").c_str(), 0.0f);
const float resZ = env.GetValue(getKey("resZ").c_str(), 0.0f);
const float eff = env.GetValue(getKey("eff").c_str(), 0.0f);
const int type = env.GetValue(getKey("type").c_str(), 0);

// void AddLayer(TString name, float r, float z, float x0, float xrho, float resRPhi = 0.0f, float resZ = 0.0f, float eff = 0.0f, int type = 0);
AddLayer(layer.c_str(), r, z, x0, xrho, resRPhi, resZ, eff, type);
}
infile.close();
}

float FastTracker::Dist(float z, float r)
Expand All @@ -298,7 +322,7 @@
index = 1;
z0 = -4 * sigmaD + i * dz0;
dist += index * (dz0 / 3.) * (1 / o2::math_utils::sqrt(o2::constants::math::TwoPI) / sigmaD) * std::exp(-z0 * z0 / 2. / sigmaD / sigmaD) * (1 / o2::math_utils::sqrt((z - z0) * (z - z0) + r * r));
if (index != 4)

Check failure on line 325 in ALICE3/Core/FastTracker.cxx

View workflow job for this annotation

GitHub Actions / O2 linter

[magic-number]

Avoid magic numbers in expressions. Assign the value to a clearly named variable or constant.
index = 4;
else
index = 2;
Expand Down Expand Up @@ -428,7 +452,7 @@
// check if layer is reached
float targetX = 1e+3;
inputTrack.getXatLabR(layers[il].getRadius(), targetX, magneticField);
if (targetX > 999.f) {

Check failure on line 455 in ALICE3/Core/FastTracker.cxx

View workflow job for this annotation

GitHub Actions / O2 linter

[magic-number]

Avoid magic numbers in expressions. Assign the value to a clearly named variable or constant.
LOGF(debug, "Failed to find intercept for layer %d at radius %.2f cm", il, layers[il].getRadius());
break; // failed to find intercept
}
Expand Down Expand Up @@ -505,7 +529,7 @@

float targetX = 1e+3;
inputTrack.getXatLabR(layers[il].getRadius(), targetX, magneticField);
if (targetX > 999)

Check failure on line 532 in ALICE3/Core/FastTracker.cxx

View workflow job for this annotation

GitHub Actions / O2 linter

[magic-number]

Avoid magic numbers in expressions. Assign the value to a clearly named variable or constant.
continue; // failed to find intercept

if (!inputTrack.propagateTo(targetX, magneticField)) {
Expand Down Expand Up @@ -577,7 +601,7 @@
// backpropagate to original radius
float finalX = 1e+3;
bool inPropStatus = inwardTrack.getXatLabR(initialRadius, finalX, magneticField);
if (finalX > 999) {

Check failure on line 604 in ALICE3/Core/FastTracker.cxx

View workflow job for this annotation

GitHub Actions / O2 linter

[magic-number]

Avoid magic numbers in expressions. Assign the value to a clearly named variable or constant.
LOG(debug) << "Failed to find intercept for initial radius " << initialRadius << " cm, x = " << finalX << " and status " << inPropStatus << " and sn = " << inwardTrack.getSnp() << " r = " << inwardTrack.getY() * inwardTrack.getY();
return -3; // failed to find intercept
}
Expand All @@ -587,7 +611,7 @@
}

// only attempt to continue if intercepts are at least four
if (nIntercepts < 4)

Check failure on line 614 in ALICE3/Core/FastTracker.cxx

View workflow job for this annotation

GitHub Actions / O2 linter

[magic-number]

Avoid magic numbers in expressions. Assign the value to a clearly named variable or constant.
return nIntercepts;

// generate efficiency
Expand Down Expand Up @@ -631,7 +655,7 @@
if (fcovm[ii][jj] * fcovm[ii][jj] > std::abs(fcovm[ii][ii] * fcovm[jj][jj])) {
rubenConditional = true;
if (makePositiveDefinite) {
fcovm[ii][jj] = TMath::Sign(1, fcovm[ii][jj]) * covMatFactor * sqrt(std::abs(fcovm[ii][ii] * fcovm[jj][jj]));

Check failure on line 658 in ALICE3/Core/FastTracker.cxx

View workflow job for this annotation

GitHub Actions / O2 linter

[std-prefix]

Use std:: prefix for names from the std namespace.
}
}
}
Expand Down Expand Up @@ -669,7 +693,7 @@
for (int j = 0; j < 5; ++j)
val += eigVec[j][ii] * outputTrack.getParam(j);
// smear parameters according to eigenvalues
params_[ii] = gRandom->Gaus(val, sqrt(eigVal[ii]));

Check failure on line 696 in ALICE3/Core/FastTracker.cxx

View workflow job for this annotation

GitHub Actions / O2 linter

[std-prefix]

Use std:: prefix for names from the std namespace.
}

// invert eigenvector matrix
Expand All @@ -682,7 +706,7 @@
outputTrack.setParam(val, ii);
}
// should make a sanity check that par[2] sin(phi) is in [-1, 1]
if (fabs(outputTrack.getParam(2)) > 1.) {

Check failure on line 709 in ALICE3/Core/FastTracker.cxx

View workflow job for this annotation

GitHub Actions / O2 linter

[std-prefix]

Use std:: prefix for names from the std namespace.
LOG(info) << " --- smearTrack failed sin(phi) sanity check: " << outputTrack.getParam(2);
return -2;
}
Expand Down
30 changes: 30 additions & 0 deletions ALICE3/macros/testFastTracker.cxx
Original file line number Diff line number Diff line change
@@ -0,0 +1,30 @@
// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
// All rights not expressly granted are reserved.
//
// This software is distributed under the terms of the GNU General Public
// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
//
// In applying this license CERN does not waive the privileges and immunities
// granted to it by virtue of its status as an Intergovernmental Organization
// or submit itself to any jurisdiction.

/// \file testFastTracker.C
/// \author Nicolò Jacazio nicolo.jacazio@cern.ch
/// \brief Test the FastTracker functionality

#include "ALICE3/Core/FastTracker.h"

#include "CCDB/BasicCCDBManager.h"
#include "DataFormatsParameters/GRPLHCIFData.h"

void testFastTracker(std::string geometryFile = "a3geo.ini")
{

// auto& ccdb = o2::ccdb::BasicCCDBManager::instance();
// ccdb.setURL("http://alice-ccdb.cern.ch");
o2::fastsim::FastTracker fastTracker;
fastTracker.AddGenericDetector(geometryFile);
// fastTracker.AddGenericDetector(geometryFile, &ccdb);
fastTracker.Print();
}
Loading