13#define TRG_SHORT_NAMES
14#define TRGCDC_SHORT_NAMES
16#define ENV_PATH "BELLE2_LOCAL_DIR"
22#include "cdc/geometry/CDCGeometryPar.h"
23#include "trg/trg/Debug.h"
24#include "trg/trg/Utilities.h"
25#include "trg/trg/Time.h"
26#include "trg/trg/Signal.h"
27#include "trg/trg/Constants.h"
28#include "trg/cdc/TRGCDC.h"
29#include "trg/cdc/Layer.h"
30#include "trg/cdc/Cell.h"
31#include "trg/cdc/Wire.h"
32#include "trg/cdc/WireHit.h"
33#include "trg/cdc/HoughFinder.h"
34#include "trg/cdc/Segment.h"
35#include "trg/cdc/SegmentHit.h"
36#include "trg/cdc/Circle.h"
37#include "trg/cdc/TRGCDCTrack.h"
38#include "trg/cdc/Link.h"
39#include "trg/cdc/Relation.h"
40#include "trg/cdc/Fitter3DUtility.h"
41#include "trg/cdc/Fitter3D.h"
42#include "trg/cdc/Helix.h"
43#include "trg/cdc/JLUT.h"
44#include "trg/cdc/JSignal.h"
45#include "trg/cdc/JSignalData.h"
46#include "trg/cdc/WireHitMC.h"
47#include "trg/cdc/TrackMC.h"
48#include "trg/cdc/FrontEnd.h"
49#include "trg/cdc/Merger.h"
50#include "trg/cdc/LUT.h"
51#include "trg/cdc/EventTime.h"
53#ifdef TRGCDC_DISPLAY_HOUGH
54#include "trg/cdc/DisplayRphi.h"
55#include "trg/cdc/DisplayHough.h"
56namespace Belle2_TRGCDC {
57 Belle2::TRGCDCDisplayHough* H0 = 0;
58 Belle2::TRGCDCDisplayHough* H1 = 0;
63#ifdef TRGCDC_DISPLAY_HOUGH
64using namespace Belle2_TRGCDC;
76 return string(
"TRGCDCHoughFinder 5.24");
82 const string& mappingFilePlus,
83 const string& mappingFileMinus,
87 _circleH(
"CircleHough"),
89 _peakFinder(
"PeakFinder"),
96 const string fMap[2] = {mappingFilePlus, mappingFileMinus};
97 const string fName[2] = {string(
"circle hough plus"), string(
"circle hough minus")};
99 for (
unsigned f = 0; f < 2; f++) {
100 const string fn = fMap[f];
101 ifstream infile(fn.c_str(), ios::in);
103 B2FATAL(
"Cannot open Hough mapping file " << fn);
109 while (!isdigit(infile.peek())) {
110 getline(infile, ignores);
119 infile >> nX >> nY >> Ymin >> Ymax;
132#ifdef TRGCDC_DISPLAY_HOUGH
134 H0 =
new TCDisplayHough(
"Plus");
140 H1 =
new TCDisplayHough(
"Minus");
151 unsigned axialSuperLayerId = 0;
154 const unsigned nWires = l->
nCells();
156 if (! nWires)
continue;
157 if ((* l)[0]->stereo())
continue;
161 for (
unsigned j = 0; j < nWires; j++) {
162 const TCCell& w = * (* l)[j];
163 const float x = w.xyPosition().x();
164 const float y = w.xyPosition().y();
167 _plane[0]->
vote(x, y, +1, axialSuperLayerId, 1);
171 _plane[1]->
vote(x, y, -1, axialSuperLayerId, 1);
186#ifdef TRGCDC_DISPLAY_HOUGH
191 cout <<
"TRGCDCHoughFinder ... Hough displays deleted" << endl;
197 vector<TCTrack*>& trackList2D)
const
200 const string sn =
"Hough Finder Finding (trasan version)";
209 for (
unsigned i = 0; i < nLayers; i++) {
211 for (
unsigned j = 0; j < hits.size(); j++) {
212 _plane[0]->
vote(i, hits[j]->cell().localId());
213 _plane[1]->
vote(i, hits[j]->cell().localId());
219#ifdef TRGCDC_DISPLAY_HOUGH
220 string stg =
"2D : Hough : Results of Peak Finding";
223 H0->information(inf);
225 H0->area().append(
_plane[0]);
228 H1->information(inf);
230 H1->area().append(
_plane[1]);
240 for (
unsigned pm = 0; pm < 2; pm++) {
241 for (
unsigned i = 0; i < peaks[pm].size(); i++) {
242 const unsigned peakId = peaks[pm][i];
246 trackList2D.push_back(t);
248#ifdef TRGCDC_DISPLAY_HOUGH
249 vector<const TCTrack*> cc;
251 const string stgNow =
"doFinding : track made";
252 const string infNow =
" ";
255 D->information(infNow);
256 D->area().append(cc, Gdk::Color(
"#FF0066009900"));
272 vector<TCLink*> bests;
273 vector<TCLink*> layers[9];
274 TCLink::separate(links, 9, layers);
277 for (
unsigned i = 0; i < 9; i++) {
284 for (
unsigned i = 0; i < 9; i++) {
285 if (layers[i].size() == 0)
continue;
286 if (layers[i].size() == 1) {
287 bests.push_back(layers[i][0]);
291 TCLink* best = layers[i][0];
292 int timeMin = (layers[i][0]->cell()->signal())[0]->time();
293 for (
unsigned j = 1; j < layers[i].size(); j++) {
294 const TRGTime& t = * (layers[i][j]->cell()->signal())[0];
295 if (t.time() < timeMin) {
301 bests.push_back(best);
308 vector<TRGCDCTrack*>& trackList2DFitted)
const
311 const string sn =
"Hough Finder Fitting (trasan version)";
315 unsigned nCircles = 0;
316 for (
unsigned pm = 0; pm < 2; pm++) {
317 for (
unsigned i = 0; i < peaks[pm].size(); i++) {
318 const unsigned peakId = peaks[pm][i];
321 vector<TCLink*> links;
322 vector<const TCSegment*> segments;
324 for (
unsigned j = 0; j < nLayers; j++) {
325 const vector<unsigned>& ptn =
327 for (
unsigned k = 0; k < ptn.size(); k++) {
329 segments.push_back(& s);
331 TCLink* l =
new TCLink(0, s.hit());
343 c.name(
"CircleFitted_" + TRGUtil::itostring(nCircles));
348 <<
"plane" << pm <<
",serialId=" << peakId << endl;
351 for (
unsigned j = 0; j < segments.size(); j++) {
352 cout << segments[j]->name();
353 if (j != (segments.size() - 1))
363 TCTrack& t = *
new TCTrack(c);
364 t.name(
"Track_" + TRGUtil::itostring(i));
365 trackList2DFitted.push_back(& t);
372#ifdef TRGCDC_DISPLAY_HOUGH
373 vector<const TCCircle*> cc;
375 const string stg =
"doFitting : trasan method";
376 const string inf =
" ";
382 D->area().append(cc, Gdk::Color(
"blue"));
413 std::map<std::string, Belle2::TRGCDCJLUT* >& mLutStorage)
419 mSignalStorage[
"invPhiAxMin"] =
Belle2::TRGCDCJSignal(M_PI, mSignalStorage[
"phi_4"].getToReal(), commonData);
420 mSignalStorage[
"invPhiAxMax"] =
Belle2::TRGCDCJSignal(2 * M_PI, mSignalStorage[
"phi_4"].getToReal(), commonData);
423 for (
unsigned iSt = 0; iSt < 5; iSt++) {
424 string t_inputName =
"phi_" + to_string(iSt);
425 string t_outputName =
"cosPhi_" + to_string(iSt);
427 if (!mLutStorage.count(t_outputName)) {
429 mLutStorage[t_outputName]->setFloatFunction(
430 [ = ](
double aValue) ->
double{
return cos(aValue);},
431 mSignalStorage[t_inputName],
432 mSignalStorage[
"invPhiAxMin"], mSignalStorage[
"invPhiAxMax"], mSignalStorage[t_inputName].getToReal(),
436 mLutStorage[t_outputName]->operate(mSignalStorage[t_inputName], mSignalStorage[t_outputName]);
441 std::map<std::string, Belle2::TRGCDCJLUT* >& mLutStorage)
445 mSignalStorage[
"invPhiAxMin"] =
Belle2::TRGCDCJSignal(-M_PI / 2, mSignalStorage[
"phi_4"].getToReal(), commonData);
446 mSignalStorage[
"invPhiAxMax"] =
Belle2::TRGCDCJSignal(M_PI / 2, mSignalStorage[
"phi_4"].getToReal(), commonData);
449 for (
unsigned iSt = 0; iSt < 5; iSt++) {
450 string t_sininputName =
"phi_" + to_string(iSt);
451 string t_sinoutputName =
"sinPhi_" + to_string(iSt);
453 if (!mLutStorage.count(t_sinoutputName)) {
456 mLutStorage[t_sinoutputName]->setFloatFunction(
457 [ = ](
double aValue) ->
double{
return sin(aValue);},
458 mSignalStorage[t_sininputName],
459 mSignalStorage[
"invPhiAxMin"], mSignalStorage[
"invPhiAxMax"], mSignalStorage[t_sininputName].getToReal(),
463 mLutStorage[t_sinoutputName]->operate(mSignalStorage[t_sininputName], mSignalStorage[t_sinoutputName]);
468 std::map<std::string, Belle2::TRGCDCJLUT* >& mLutStorage)
471 for (
unsigned iSt = 0; iSt < 5; iSt++) {
472 mSignalStorage[
"cosPhiMul_" + to_string(iSt)] <= mSignalStorage[
"cosPhi_" + to_string(iSt)] * mSignalStorage[
"cosPhi_" + to_string(
476 for (
unsigned iSt = 0; iSt < 5; iSt++) {
477 mSignalStorage[
"sinPhiMul_" + to_string(iSt)] <= mSignalStorage[
"sinPhi_" + to_string(iSt)] * mSignalStorage[
"sinPhi_" + to_string(
481 for (
unsigned iSt = 0; iSt < 5; iSt++) {
482 mSignalStorage[
"cossinPhiMul_" + to_string(iSt)] <= mSignalStorage[
"cosPhi_" + to_string(iSt)] * mSignalStorage[
"sinPhi_" +
486 for (
unsigned iSt = 0; iSt < 5; iSt++) {
487 mSignalStorage[
"rcosPhiMul_" + to_string(iSt)] <= mSignalStorage[
"cosPhi_" + to_string(iSt)] * mSignalStorage[
"2Drr_" + to_string(
491 for (
unsigned iSt = 0; iSt < 5; iSt++) {
492 mSignalStorage[
"rsinPhiMul_" + to_string(iSt)] <= mSignalStorage[
"sinPhi_" + to_string(iSt)] * mSignalStorage[
"2Drr_" + to_string(
495 mSignalStorage[
"cossum_p1_0"] <= mSignalStorage[
"cosPhiMul_0"] + mSignalStorage[
"cosPhiMul_1"];
496 mSignalStorage[
"cossum_p1_1"] <= mSignalStorage[
"cosPhiMul_2"] + mSignalStorage[
"cosPhiMul_3"];
497 mSignalStorage[
"cossum_p1_2"] <= mSignalStorage[
"cosPhiMul_4"] ;
498 mSignalStorage[
"cossum_p2_0"] <= mSignalStorage[
"cossum_p1_0"] + mSignalStorage[
"cossum_p1_1"];
499 mSignalStorage[
"cossum_p2_1"] <= mSignalStorage[
"cossum_p1_2"] ;
500 mSignalStorage[
"cossum"] <= mSignalStorage[
"cossum_p2_0"] + mSignalStorage[
"cossum_p2_1"];
501 mSignalStorage[
"sinsum_p1_0"] <= mSignalStorage[
"sinPhiMul_0"] + mSignalStorage[
"sinPhiMul_1"];
502 mSignalStorage[
"sinsum_p1_1"] <= mSignalStorage[
"sinPhiMul_2"] + mSignalStorage[
"sinPhiMul_3"];
503 mSignalStorage[
"sinsum_p1_2"] <= mSignalStorage[
"sinPhiMul_4"] ;
504 mSignalStorage[
"sinsum_p2_0"] <= mSignalStorage[
"sinsum_p1_0"] + mSignalStorage[
"sinsum_p1_1"];
505 mSignalStorage[
"sinsum_p2_1"] <= mSignalStorage[
"sinsum_p1_2"] ;
506 mSignalStorage[
"sinsum"] <= mSignalStorage[
"sinsum_p2_0"] + mSignalStorage[
"sinsum_p2_1"];
507 mSignalStorage[
"cossinsum_p1_0"] <= mSignalStorage[
"cossinPhiMul_0"] + mSignalStorage[
"cossinPhiMul_1"];
508 mSignalStorage[
"cossinsum_p1_1"] <= mSignalStorage[
"cossinPhiMul_2"] + mSignalStorage[
"cossinPhiMul_3"];
509 mSignalStorage[
"cossinsum_p1_2"] <= mSignalStorage[
"cossinPhiMul_4"] ;
510 mSignalStorage[
"cossinsum_p2_0"] <= mSignalStorage[
"cossinsum_p1_0"] + mSignalStorage[
"cossinsum_p1_1"];
511 mSignalStorage[
"cossinsum_p2_1"] <= mSignalStorage[
"cossinsum_p1_2"] ;
512 mSignalStorage[
"cossinsum"] <= mSignalStorage[
"cossinsum_p2_0"] + mSignalStorage[
"cossinsum_p2_1"];
513 mSignalStorage[
"rcossum_p1_0"] <= mSignalStorage[
"rcosPhiMul_0"] + mSignalStorage[
"rcosPhiMul_1"];
514 mSignalStorage[
"rcossum_p1_1"] <= mSignalStorage[
"rcosPhiMul_2"] + mSignalStorage[
"rcosPhiMul_3"];
515 mSignalStorage[
"rcossum_p1_2"] <= mSignalStorage[
"rcosPhiMul_4"] ;
516 mSignalStorage[
"rcossum_p2_0"] <= mSignalStorage[
"rcossum_p1_0"] + mSignalStorage[
"rcossum_p1_1"];
517 mSignalStorage[
"rcossum_p2_1"] <= mSignalStorage[
"rcossum_p1_2"] ;
518 mSignalStorage[
"rcossum"] <= mSignalStorage[
"rcossum_p2_0"] + mSignalStorage[
"rcossum_p2_1"];
519 mSignalStorage[
"rsinsum_p1_0"] <= mSignalStorage[
"rsinPhiMul_0"] + mSignalStorage[
"rsinPhiMul_1"];
520 mSignalStorage[
"rsinsum_p1_1"] <= mSignalStorage[
"rsinPhiMul_2"] + mSignalStorage[
"rsinPhiMul_3"];
521 mSignalStorage[
"rsinsum_p1_2"] <= mSignalStorage[
"rsinPhiMul_4"] ;
522 mSignalStorage[
"rsinsum_p2_0"] <= mSignalStorage[
"rsinsum_p1_0"] + mSignalStorage[
"rsinsum_p1_1"];
523 mSignalStorage[
"rsinsum_p2_1"] <= mSignalStorage[
"rsinsum_p1_2"] ;
524 mSignalStorage[
"rsinsum"] <= mSignalStorage[
"rsinsum_p2_0"] + mSignalStorage[
"rsinsum_p2_1"];
525 mSignalStorage[
"hcx_p1_0"] <= mSignalStorage[
"rcossum"] * mSignalStorage[
"sinsum"];
526 mSignalStorage[
"hcx_p1_1"] <= mSignalStorage[
"rsinsum"] * mSignalStorage[
"cossinsum"];
527 mSignalStorage[
"hcx_p2"] <= mSignalStorage[
"hcx_p1_0"] - mSignalStorage[
"hcx_p1_1"];
528 mSignalStorage[
"hcy_p1_0"] <= mSignalStorage[
"rsinsum"] * mSignalStorage[
"cossum"];
529 mSignalStorage[
"hcy_p1_1"] <= mSignalStorage[
"rcossum"] * mSignalStorage[
"cossinsum"];
530 mSignalStorage[
"hcy_p2"] <= mSignalStorage[
"hcy_p1_0"] - mSignalStorage[
"hcy_p1_1"];
531 mSignalStorage[
"den_p1"] <= mSignalStorage[
"cossum"] * mSignalStorage[
"sinsum"];
532 mSignalStorage[
"den_p2"] <= mSignalStorage[
"cossinsum"] * mSignalStorage[
"cossinsum"];
533 mSignalStorage[
"den"] <= mSignalStorage[
"den_p1"] - mSignalStorage[
"den_p2"];
537 mSignalStorage[
"denMax"] =
Belle2::TRGCDCJSignal(1.332705, mSignalStorage[
"den"].getToReal(), commonData);
543 vector<pair<Belle2::TRGCDCJSignal, vector<pair<Belle2::TRGCDCJSignal*, Belle2::TRGCDCJSignal> > > > t_data;
547 vector<pair<Belle2::TRGCDCJSignal*, Belle2::TRGCDCJSignal> > t_assigns = {
548 make_pair(&mSignalStorage[
"den_c"], mSignalStorage[
"denMax"])
551 t_data.push_back(make_pair(t_compare, t_assigns));
556 make_pair(&mSignalStorage[
"den_c"], mSignalStorage[
"den"].limit(mSignalStorage[
"denMin"], mSignalStorage[
"denMax"]))
559 t_data.push_back(make_pair(t_compare, t_assigns));
564 make_pair(&mSignalStorage[
"den_c"], mSignalStorage[
"denMin"])
567 t_data.push_back(make_pair(t_compare, t_assigns));
572 unsigned long long t_int = mSignalStorage[
"den_c"].getMaxInt();
573 double t_toReal = mSignalStorage[
"den_c"].getToReal();
574 double t_actual = mSignalStorage[
"den_c"].getMaxActual();
575 mSignalStorage[
"iDenMin"] =
Belle2::TRGCDCJSignal(t_int, t_toReal, t_int, t_int, t_actual, t_actual, t_actual, -1, commonData);
576 t_int = mSignalStorage[
"den_c"].getMinInt();
577 t_actual = mSignalStorage[
"den_c"].getMinActual();
578 mSignalStorage[
"iDenMax"] =
Belle2::TRGCDCJSignal(t_int, t_toReal, t_int, t_int, t_actual, t_actual, t_actual, -1, commonData);
581 if (!mLutStorage.count(
"iDen")) {
583 mLutStorage[
"iDen"]->setFloatFunction(
584 [ = ](
double aValue) ->
double{
return 0.5 / aValue;},
585 mSignalStorage[
"den_c"],
586 mSignalStorage[
"iDenMin"], mSignalStorage[
"iDenMax"], pow(1 / mSignalStorage[
"cossinsum"].getToReal(), 2),
591 mLutStorage[
"iDen"]->operate(mSignalStorage[
"den_c"], mSignalStorage[
"iDen"]);
593 mSignalStorage[
"hcx"] <= mSignalStorage[
"hcx_p2"] * mSignalStorage[
"iDen"];
594 mSignalStorage[
"hcy"] <= mSignalStorage[
"hcy_p2"] * mSignalStorage[
"iDen"];
596 mSignalStorage[
"hcx2"] <= mSignalStorage[
"hcx"] * mSignalStorage[
"hcx"];
597 mSignalStorage[
"hcy2"] <= mSignalStorage[
"hcy"] * mSignalStorage[
"hcy"];
598 mSignalStorage[
"sumhc"] <= mSignalStorage[
"hcx2"] + mSignalStorage[
"hcy2"];
600 mSignalStorage[
"sumhcMin"] =
Belle2::TRGCDCJSignal(4489, mSignalStorage[
"sumhc"].getToReal(), commonData);
601 mSignalStorage[
"sumhcMax"] =
Belle2::TRGCDCJSignal(2560000, mSignalStorage[
"sumhc"].getToReal(), commonData);
607 vector<pair<Belle2::TRGCDCJSignal, vector<pair<Belle2::TRGCDCJSignal*, Belle2::TRGCDCJSignal> > > > t_data;
611 vector<pair<Belle2::TRGCDCJSignal*, Belle2::TRGCDCJSignal> > t_assigns = {
612 make_pair(&mSignalStorage[
"sumhc_c"], mSignalStorage[
"sumhcMax"])
615 t_data.push_back(make_pair(t_compare, t_assigns));
620 make_pair(&mSignalStorage[
"sumhc_c"], mSignalStorage[
"sumhc"].limit(mSignalStorage[
"sumhcMin"], mSignalStorage[
"sumhcMax"]))
623 t_data.push_back(make_pair(t_compare, t_assigns));
628 make_pair(&mSignalStorage[
"sumhc_c"], mSignalStorage[
"sumhcMin"])
631 t_data.push_back(make_pair(t_compare, t_assigns));
635 unsigned long long t_int = mSignalStorage[
"sumhc_c"].getMaxInt();
636 double t_toReal = mSignalStorage[
"sumhc_c"].getToReal();
637 double t_actual = mSignalStorage[
"sumhc_c"].getMaxActual();
638 mSignalStorage[
"iRhoMax"] =
Belle2::TRGCDCJSignal(t_int, t_toReal, t_int, t_int, t_actual, t_actual, t_actual, -1, commonData);
639 t_int = mSignalStorage[
"sumhc_c"].getMinInt();
640 t_actual = mSignalStorage[
"sumhc_c"].getMinActual();
641 mSignalStorage[
"iRhoMin"] =
Belle2::TRGCDCJSignal(t_int, t_toReal, t_int, t_int, t_actual, t_actual, t_actual, -1, commonData);
644 if (!mLutStorage.count(
"iRho")) {
646 mLutStorage[
"iRho"]->setFloatFunction(
647 [ = ](
double aValue) ->
double{
return abs(
sqrt(aValue));},
648 mSignalStorage[
"sumhc_c"],
649 mSignalStorage[
"iRhoMin"], mSignalStorage[
"iRhoMax"], mSignalStorage[
"sumhc"].getToReal(),
653 mLutStorage[
"iRho"]->operate(mSignalStorage[
"sumhc_c"], mSignalStorage[
"iRho"]);
655 cout <<
"<<<iRho>>>" << endl;
656 mSignalStorage[
"iRho"].dump();
662 const string& mappingFileMinus)
664 const string sn =
"mappingByFile";
667 const string* fMap[2] = {& mappingFilePlus, & mappingFileMinus};
669 for (
unsigned f = 0; f < 2; f++) {
671 const string& fn = * fMap[f];
672 ifstream infile(fn.c_str(), ios::in);
674 cout <<
" !!! can not open file" << endl
676 <<
" Mapping aborted" << endl;
683 vector<vector<unsigned>> tsf;
687 while (!isdigit(infile.peek())) {
688 getline(infile, ignores);
691 getline(infile, ignores);
695 while (getline(infile, car)) {
698 istringstream in(car);
699 vector<unsigned> slts;
707 if (j != 0 && i != 1)
709 if (j != 1 && i != 2)
716 unsigned axialSuperLayerId = 0;
719 const unsigned nWires = l->
nCells();
721 if (! nWires)
continue;
722 if ((* l)[0]->stereo())
continue;
725 for (
unsigned j = 0; j < nWires; j++) {
730 const unsigned n = x.size();
731 for (
unsigned k = 0; k < n; k++) {
735 for (
unsigned itsf = 0; itsf < tsf[k].size();) {
736 const unsigned sl = tsf[k][itsf];
737 const unsigned id = tsf[k][itsf + 1];
738 if ((i == sl) && (j ==
id)) {
756 std::vector<TRGCDCTrack*>& trackList2DFitted)
767 std::vector<TRGCDCTrack*>& trackList2DFitted)
770 const string sn =
"HoughFinder::doFindingAndFitting (trasan version)";
773 vector<unsigned> peaks[2];
784 std::vector<TRGCDCTrack*>& trackList2DFitted)
787 const string sn =
"HoughFinder::doFindingAndFitting (kaiyu version)";
790 vector<vector<unsigned>> peaks[2];
792 doFitting(trackList2D, trackList2DFitted);
801 std::vector<TRGCDCTrack*>& trackList2D)
804 const string sn =
"HoughFinder::doFinding";
813 for (
unsigned i = 0; i < nLayers; i++) {
815 for (
unsigned j = 0; j < hits.size(); j++) {
816 _plane[0]->
vote(i, hits[j]->cell().localId());
817 _plane[1]->
vote(i, hits[j]->cell().localId());
828 for (
unsigned pm = 0; pm < 2; pm++) {
829 for (
unsigned i = 0; i < peaks[pm].size(); i++) {
833 trackList2D.push_back(t);
835#ifdef TRGCDC_DISPLAY_HOUGH
836 vector<const TCTrack*> cc;
838 const string stg =
"doFinding : track made";
839 const string inf =
" ";
843 D->area().append(cc, Gdk::Color(
"#FF0066009900"));
858 std::vector<TRGCDCTrack*>& trackList2DFitted)
861 const string sn =
"Hough Finder Fitting";
869 bool fEvtTime =
true;
874 for (
unsigned int iTrack = 0; iTrack < trackList2D.size(); iTrack++) {
875 TCTrack& aTrack = *
new TCTrack(* trackList2D[iTrack]);
876 trackList2DFitted.push_back(& aTrack);
879 vector<double> nWires(9);
880 vector<double> rr(9);
882 vector<double> wirePhi2DError(5);
883 wirePhi2DError[0] = 0.0085106;
884 wirePhi2DError[1] = 0.0039841;
885 wirePhi2DError[2] = 0.0025806;
886 wirePhi2DError[3] = 0.0019084;
887 wirePhi2DError[4] = 0.001514;
888 vector<double> driftPhi2DError(5);
889 driftPhi2DError[0] = 0.0085106;
890 driftPhi2DError[1] = 0.0039841;
891 driftPhi2DError[2] = 0.0025806;
892 driftPhi2DError[3] = 0.0019084;
893 driftPhi2DError[4] = 0.001514;
898 const vector<TCLink*>& links = aTrack.links(iSL);
899 const unsigned nSegments = links.size();
903 << nSegments << endl;
907 bool priorityHitTS = 0;
908 for (
unsigned iTS = 0; iTS < nSegments; iTS++) {
909 const TCSegment* _segment =
dynamic_cast<const TCSegment*
>(& links[iTS]->hit()->cell());
910 if (_segment->center().hit() != 0) priorityHitTS = 1;
912 if (nSegments != 1) {
913 if (nSegments == 0) {
917 <<
"=> Not enough TS." << endl;
921 <<
"=> multiple TS are assigned." << endl;
924 if (priorityHitTS == 0) {
928 <<
"=> There are no priority hit TS" << endl;
932 if (trackFull == 0) {
934 CLHEP::HepVector helixParameters(5);
935 helixParameters = aTrack.helix().a();
937 aTrack.setHelix(helix);
941 for (
unsigned iSL = 0; iSL < 9; iSL++) {
946 rr2D = vector<double>({rr[0], rr[2], rr[4], rr[6], rr[8]});
949 vector<double>wirePhi(9);
950 vector<double>driftLength(9);
952 vector<double>lutLR(9);
953 vector<double>mcLR(9);
954 bool fmcLR =
false, fLRLUT =
true;
955 for (
unsigned iSL = 0; iSL < 9; iSL = iSL + 2) {
956 const vector<TCLink*>& links = aTrack.links(iSL);
957 const TCSegment* _segment =
dynamic_cast<const TCSegment*
>(& links[0]->hit()->cell());
958 wirePhi[iSL] = _segment->localId() / nWires[iSL] * 4 * M_PI;
959 lutLR[iSL] = _segment->LUT()->getValue(_segment->lutPattern());
960 mcLR[iSL] = _segment->hit()->mcLR();
961 driftLength[iSL] = _segment->hit()->drift();
963 if (fmcLR) LR[iSL] = mcLR[iSL];
965 else if (fLRLUT) LR[iSL] = lutLR[iSL];
971 vector<double>phi2DError(5);
972 for (
unsigned iAx = 0; iAx < 5; iAx++) {
973 if (LR[2 * iAx] != 2) phi2DError[iAx] = driftPhi2DError[iAx];
974 else phi2DError[iAx] = wirePhi2DError[iAx];
977 vector<double>phi2D(5);
978 for (
unsigned iAx = 0; iAx < 5; iAx++) {
980 driftLength[iAx * 2],
985 cout <<
TRGDebug::tab() <<
"eventTime: " << eventTime << endl;
986 for (
int i = 0 ; i < 5 ; i++) {
991 driftLength[i * 2] << endl;
997 double rho, phi0, pt, chi2;
1001 vector<double>phi2DInvError(5);
1002 for (
unsigned iAx = 0; iAx < 5; iAx++) {
1003 phi2DInvError[iAx] = 1 / phi2DError[iAx];
1008 pt = 0.3 * 1.5 * rho / 100;
1017 CLHEP::HepVector a(5);
1018 a = aTrack.helix().a();
1019 a[1] = (aTrack.charge() < 0) ? fmod(phi0 + M_PI, 2 * M_PI) : phi0;
1020 a[2] = 1 / pt * aTrack.charge();
1022 aTrack.setFitted(1);
1023 aTrack.setHelix(helix);
1024 aTrack.set2DFitChi2(chi2);
1030 double rhoMax = 1600;
1031 int phiBitSize = 12;
1032 int rhoBitSize = 12;
1033 vector<tuple<string, double, int, double, double, int> > t_values = {
1034 make_tuple(
"phi_0", phi2D[0], phiBitSize, -2 * M_PI, 2 * M_PI, 0),
1035 make_tuple(
"phi_1", phi2D[1], phiBitSize, -2 * M_PI, 2 * M_PI, 0),
1036 make_tuple(
"phi_2", phi2D[2], phiBitSize, -2 * M_PI, 2 * M_PI, 0),
1037 make_tuple(
"phi_3", phi2D[3], phiBitSize, -2 * M_PI, 2 * M_PI, 0),
1038 make_tuple(
"phi_4", phi2D[4], phiBitSize, -2 * M_PI, 2 * M_PI, 0),
1039 make_tuple(
"rho", rho, rhoBitSize, rhoMin, rhoMax, 0),
1040 make_tuple(
"2Drr_0", rr[0], 12, 18.8, 105.68, 0),
1041 make_tuple(
"2Drr_1", rr[2], 12, 18.8, 105.68, 0),
1042 make_tuple(
"2Drr_2", rr[4], 12, 18.8, 105.68, 0),
1043 make_tuple(
"2Drr_3", rr[6], 12, 18.8, 105.68, 0),
1044 make_tuple(
"2Drr_4", rr[8], 12, 18.8, 105.68, 0),
1052 bool updateSecondName =
false;
1055 (*it).second.setName((*it).first);
1057 updateSecondName =
true;
1060 if (updateSecondName) {
1075 it->second->makeCOE(it->first +
".coe");
1082#ifdef TRGCDC_DISPLAY_HOUGH
1083 vector<const TCTrack*> cc;
1084 cc.push_back(& aTrack);
1085 const string stg =
"doFitting : Kaiyu method";
1086 const string inf =
" ";
1089 D->information(inf);
1092 D->area().append(cc, Gdk::Color(
"blue"));
1108 const TCHTransformationCircle& tc =
1109 (TCHTransformationCircle&)
_plane[pm]->transformation();
1114 const double r = cc.y();
1115 const double phi = cc.x();
1118 TCCircle c(r, phi, pm ? -1. : 1., *
_plane[pm]);
1119 c.name(
"circle_" + TRGUtil::itostring(
int(peakId) * (pm ? -1 : 1)));
1122 cout <<
TRGDebug::tab() <<
"plane" << pm <<
",serialId=" << peakId
1128 vector<TCLink*> links;
1129 vector<const TCSegment*> segments;
1131 for (
unsigned j = 0; j < nLayers; j++) {
1132 const vector<unsigned>& ptn =
1134 for (
unsigned k = 0; k < ptn.size(); k++) {
1136 segments.push_back(& s);
1138 TCLink* l =
new TCLink(0, s.hit());
1146 cout <<
TRGDebug::tab() <<
"attched segments below" << endl;
1148 for (
unsigned j = 0; j < segments.size(); j++) {
1149 cout << segments[j]->name();
1150 if (j != (segments.size() - 1))
1157 TCTrack* t =
new TCTrack(c);
1158 t->name(
"track_" + TRGUtil::itostring(
int(peakId) * (pm ? -1 : 1)));
The Class for CDC Geometry Parameters.
unsigned nWiresInLayer(int layerId) const
Returns wire numbers in a layer.
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
double senseWireR(int layerId) const
Returns radius of sense wire in each layer.
TRGCDCHelix parameter class.
TRGCDCJSignalData * _commonData
For VHDL code.
TRGCDCHoughPlaneMulti2 * _plane[2]
Hough planes, for + and - charges.
const TRGCDC & _cdc
CDCTRG.
TRGCDCHoughTransformationCircle _circleH
Circle Hough transformtion.
std::map< std::string, TRGCDCJLUT * > m_mLutStorage
Map to hold JLuts.
const unsigned _doit
Doit version.
std::map< std::string, TRGCDCJSignal > m_mSignalStorage
Map to hold JSignals.
TRGCDCPeakFinder _peakFinder
Peak finder.
const unsigned _peakMin
Min. peak height for the peak finder.
const std::vector< unsigned > & patternId(unsigned cellId) const
returns pattern ID which activates specified cell.
A class to use LUTs for TRGCDC.
A class to hold common data for JSignals.
A class to use Signals for TRGCDC 3D tracker.
A class to represent a cell layer.
A class to represent a track segment hit in CDC.
The instance of TRGCDC is a singleton.
A class to represent a point in 2D.
A class to represent a signal timing in the trigger system.
static void rPhiFitter(double *rr, double *phi2, double *invphierror, double &rho, double &myphi0)
A circle fitter with invPhiError without fit chi2 output.
static double calPhi(double wirePhi, double driftLength, double rr, int lr)
Pre 3D fitter functions. rr is in cm scale. driftLength is in cm scale.
double sqrt(double a)
sqrt for double
int doFindingTrasan(std::vector< unsigned > peaks[], std::vector< TRGCDCTrack * > &trackList2D) const
do track finding. (trasan version)
std::vector< const TRGCDCSegmentHit * > segmentHits(void) const
returns a list of TRGCDCSegmentHit.
static void calCosPhi(std::map< std::string, Belle2::TRGCDCJSignal > &mSignalStorage, std::map< std::string, Belle2::TRGCDCJLUT * > &mLutStorage)
Calculate Cos Sin ?
void preparePatterns(unsigned layerId, unsigned nPatterns)
allocate memory for patterns.
static std::string tab(void)
returns tab spaces.
std::vector< const TRGCDCSegmentHit * > axialSegmentHits(unsigned) const
returns a list of TRGCDCSegmentHit in a axial super layer N.
unsigned nSegmentLayers(void) const
returns # of track segment layers.
const TRGCDCWire & center(void) const
returns a center wire.
void printToFile()
Utilities Function to print VHDL code.
static void valuesToMapSignals(std::vector< std::tuple< std::string, double, int, double, double, int > > const &inValues, Belle2::TRGCDCJSignalData *inCommonData, std::map< std::string, Belle2::TRGCDCJSignal > &outMap)
Values => [name, value, bitwidth, min, max, clock] Changes values to signals.
static void calSinPhi(std::map< std::string, Belle2::TRGCDCJSignal > &mSignalStorage, std::map< std::string, Belle2::TRGCDCJLUT * > &mLutStorage)
Calculate Cos Sin ?
static double calPhi(TRGCDCSegmentHit const *segmentHit, double eventTime)
Utility functions.
TRGCDCTrack * makeTrack(const unsigned serialID, const unsigned pm) const
Make a track from serial ID in Hough plane.
int FindAndFit(std::vector< TRGCDCTrack * > &trackList2D, std::vector< TRGCDCTrack * > &trackList2DFitted)
do track finding and fitting (wrapper that can choose between different versions).
void mappingByFile(const std::string &mappingFilePlus, const std::string &mappingFileMinus)
creates mappings by a file.
float charge(void) const
returns charge for this plane.
bool getPrintVhdl() const
Gets the status of m_printVhdl.
static void ifElse(std::vector< std::pair< TRGCDCJSignal, std::vector< std::pair< TRGCDCJSignal *, TRGCDCJSignal > > > > &data, int targetClock)
If else implementation with target clock.
unsigned layerId(void) const
returns layer id.
std::vector< TRGCDCLink * > selectBestHits(const std::vector< TRGCDCLink * > &links) const
selects the best(fastest) hits in each super layer.
void id(unsigned serialId, unsigned &x, unsigned &y) const
returns x and y for serialID.
virtual ~TRGCDCHoughFinder()
Destructor.
int getValue(unsigned) const
get LUT Values
const TRGCDCLayer * segmentLayer(unsigned id) const
returns a pointer to a track segment layer.
bool getPrintedToFile() const
Gets the status of m_printedToFile.
unsigned serialId(unsigned x, unsigned y) const
returns serial ID for position (x, y).
static void enterStage(const std::string &stageName)
Declare that you enter new stage.
unsigned setEntry(unsigned serialId, unsigned layerId, unsigned n)
Sets entry.
void terminate()
termination.
void findPeaks(const TRGCDCHoughPlaneMulti2 &hp, const unsigned threshold, std::vector< std::vector< unsigned > > &peaks) const
do peak finding. Kaiyu's version using p1p2Methode.
unsigned nAxialSuperLayers(void) const
return # of axial super layers.
const CLHEP::HepVector & a(void) const
returns helix parameters.
int doFindingAndFitting(std::vector< TRGCDCTrack * > &trackList2D, std::vector< TRGCDCTrack * > &trackList2DFitted)
do track finding and fitting (Kaiyu version).
unsigned nSuperLayers(void) const
returns # of super layers.
unsigned lutPattern(void) const
hit pattern containing bit for priority position
const TRGCDCSegment & segment(unsigned id) const
returns a track segment.
int doFindingAndFittingTrasan(std::vector< TRGCDCTrack * > &trackList2D, std::vector< TRGCDCTrack * > &trackList2DFitted)
do track finding and fitting (Trasan version).
void merge(void)
Merge layers into one.
void findPeaksTrasan(TRGCDCHoughPlaneMulti2 &hp, const unsigned threshold, const bool centerIsPeak, std::vector< unsigned > &peakSerialIds) const
do peak finding. This is a copy from "trasan".
const TRGCDCLUT * LUT(void) const
returns LUT
unsigned nLayers(void) const
returns # of Hough Boolean layers.
int getEventTime(void) const
returns bad hits(finding invalid hits).
int doFittingTrasan(std::vector< unsigned > peaks[], std::vector< TRGCDCTrack * > &trackList2DFitted) const
do track fitting. (old trasan version)
void clear(void) override
Clears all entries and regions.
int doFitting(std::vector< TRGCDCTrack * > &trackList2D, std::vector< TRGCDCTrack * > &trackList2DFitted)
do track fitting. (kaiyu original)
const TRGCDCSegment & axialSegment(unsigned lyrId, unsigned id) const
returns a track segment in axial layers. (lyrId is axial #)
static void rPhi(std::map< std::string, Belle2::TRGCDCJSignal > &mSignalStorage, std::map< std::string, Belle2::TRGCDCJLUT * > &mLutStorage)
Calculate r * phi ?
static int level(void)
returns the debug level.
std::vector< const TRGCDCWireHit * > hits(void) const
returns a list of TRGCDCWireHit.
float priorityTime(void) const
return priority time in TSHit.
static void leaveStage(const std::string &stageName)
Declare that you leave a stage.
int doFinding(std::vector< std::vector< unsigned > > peaks[], std::vector< TRGCDCTrack * > &trackList2D)
do track finding. (kaiyu version)
void registerPattern(unsigned layerId, unsigned id)
registers a pattern..
TRGCDCHoughFinder(const std::string &name, const TRGCDC &, unsigned peakMin, const std::string &mappingFilePlus, const std::string &mappingFileMinus, unsigned doit)
Contructor.
unsigned nCells(void) const
returns # of cells.
void entryVhdlCode()
Function to print entry VHDL code.
std::string version(void) const
returns version.
const TRGCDCSegment & segment(void) const
returns a pointer to a track segment.
void signalsVhdlCode()
Function to print definition of signal VHDL code.
TRGPoint2D position(unsigned x, unsigned y) const
returns position in Hough plain for a cell (x, y)..
void setVhdlOutputFile(const std::string &)
Sets the filename for VHDL output.
const HepGeom::Point3D< double > ORIGIN
Origin 3D point.
void vote(float rx, float ry, int charge, unsigned layerId, int weight=1)
Voting.
void buffersVhdlCode()
Function to print buffer VHDL code.
void setPrintVhdl(bool)
Sets if to print VHDL output.
unsigned localId(void) const
returns local id in a layer.
static TRGCDCJSignal comp(TRGCDCJSignal const &lhs, const std::string &operate, TRGCDCJSignal const &rhs)
Compare two signals.
Abstract base class for different kinds of events.