GeneralBrokenLines V03-01-03
using EIGEN
JavaNativeAccessWrappers.cpp
Go to the documentation of this file.
1
121#include "MilleBinary.h"
122#include "GblPoint.h"
123#include "GblTrajectory.h"
124#include "GblUtilities.h"
125
126#include "Eigen/Core"
127
128const int NROW = 5;
129const int NCOL = 5;
130
131using namespace gbl;
132using namespace Eigen;
133
138#if (defined (JNA_DEBUG) || defined (JNA_MEMORY_MONITOR))
139#define JNA_DO_MONITOR 1
140#else
141#define JNA_DO_MONITOR 0
142#endif
143
144#if JNA_DO_MONITOR
145#include <iostream>
146
147long int num_gbl_point = 0;
148long int num_gbl_traj = 0;
149long int num_mille_bin = 0;
150long int num_gbl_det_layer = 0;
151long int num_gbl_simple_helix = 0;
152long int num_gbl_helix_prediction = 0;
153
167__attribute__((destructor))
168void print_status() {
169 std::cout
170 << " GBL Structures Left\n"
171 << "GBL Points: " << num_gbl_point << "\n"
172 << "GBL Trajectories: " << num_gbl_traj << "\n"
173 << "Mille Binaries: " << num_mille_bin << "\n"
174 << "GBL Det Layers: " << num_gbl_det_layer << "\n"
175 << "GBL Simple Helix: " << num_gbl_simple_helix << "\n"
176 << "GBL Helix Pred: " << num_gbl_helix_prediction << "\n"
177 << std::flush;
178}
179#endif
180
194std::vector<GblPoint> ptr_array_to_vector(GblPoint* points[], int npoints) {
195 std::vector<GblPoint> points_vec;
196 // since we already know the size of the vector, we 'reserve' the size
197 // so that the vector doesn't need to waste time copying/moving the GblPoints
198 // around as it grows in size
199 // we do *not* use 'resize' since that would involve default-constructing
200 // all the GblPoints
201 points_vec.reserve(npoints);
202
203 for (int i{0}; i < npoints; ++i) {
204 // get the pointer
205 GblPoint* gblpoint = points[i];
206 // COPY the data into the vector,
207 // this copy-constructs a /new/ gbl point
208 points_vec.emplace_back(*(gblpoint));
209#if JNA_DO_MONITOR
210 ++num_gbl_point;
211#endif
212#ifdef JNA_DEBUG
213 std::cout << "COPY GblPoint " << gblpoint << " -> " << &(points_vec.back()) << std::endl;
214#endif
215 }
216
217 return points_vec;
218}
219
220extern "C" {
235MilleBinary* MilleBinaryCtor(const char* fileName, int filenamesize, int doublePrecision, int keepZeros, int aSize) {
236#ifdef JNA_DEBUG
237 std::cout << "MilleBinaryCtor(" << fileName << ", " << filenamesize << ", "
238 << doublePrecision << ", " << keepZeros << ", " << aSize << ")" << std::endl;
239#endif
240 std::string binName(fileName,filenamesize);
241 MilleBinary* mb = new MilleBinary(binName,doublePrecision!=0,keepZeros!=0,aSize);
242#ifdef JNA_DEBUG
243 std::cout << "MilleBinary created at " << mb << std::endl;
244#endif
245#if JNA_DO_MONITOR
246 ++num_mille_bin;
247#endif
248 return mb;
249}
250
265#ifdef JNA_DEBUG
266 std::cout << "MilleBinary_close(" << self << ")" << std::endl;
267#endif
268#if JNA_DO_MONITOR
269 --num_mille_bin;
270#endif
271 if (self) delete self;
272}
273
280GblPoint* GblPointCtor(double matrixArray[NROW*NCOL]) {
281 Map<Matrix5d> jacobian(matrixArray,5,5);
282 GblPoint* self = new GblPoint(jacobian);
283#if JNA_DO_MONITOR
284 ++num_gbl_point;
285#endif
286#ifdef JNA_DEBUG
287 std::cout << "GblPointCtor at " << self << " " << num_gbl_point << std::endl;
288#endif
289 return self;
290}
291
301#if JNA_DO_MONITOR
302 --num_gbl_point;
303#endif
304#ifdef JNA_DEBUG
305 std::cout << "GblPoint_delete(" << self << ") " << num_gbl_point << std::endl;
306#endif
307 if (self) delete self;
308}
309
310
314void GblPoint_printPoint(const GblPoint* self, unsigned int level) {
315#ifdef JNA_DEBUG
316 std::cout << "GblPoint_printPoint(" << self << ", " << level << ")" << std::endl;
317#endif
318 self->printPoint(level);
319}
320
332#ifdef JNA_DEBUG
333 std::cout << "GblPoint_getNumMeasurements(" << self << ")" << std::endl;
334#endif
335 return (self->getMeasEnd() - self->getMeasBegin());
336}
337
353 double *projArray,
354 double *resArray,
355 double *precArray,
356 double minPrecision) {
357#ifdef JNA_DEBUG
358 std::cout << "GblPoint_addMeasurement2D("
359 << self << ", " << projArray << ", " << resArray << ", " << precArray << ", " << minPrecision
360 << ")" << std::endl;
361#endif
362
363 Map<Matrix2d> aProjection(projArray,2,2);
364 Map<Vector2d> aResiduals(resArray, 2);
365 Map<Vector2d> aPrecision(precArray,2);
366
367 self->addMeasurement(aProjection, aResiduals, aPrecision, minPrecision);
368}
369
378void GblPoint_addScatterer(GblPoint* self, double *resArray, double *precArray) {
379#ifdef JNA_DEBUG
380 std::cout << "GblPoint_addScatterer(" << self << ", "
381 << resArray << ", " << precArray << ")" << std::endl;
382#endif
383 // chose to do the Vector2d addScatterer since
384 // PF's original comment "only support vector precision"
385 Eigen::Vector2d aResiduals(resArray);
386 Eigen::Vector2d aPrecision(precArray);
387
388 self->addScatterer(aResiduals,aPrecision);
389}
390
399void GblPoint_addGlobals(GblPoint* self, int *labels, int nlabels, double* derArray) {
400#ifdef JNA_DEBUG
401 std::cout << "GblPoint_addGlobals(" << self
402 << ", " << labels << ", " << nlabels << ", " << derArray << ")" << std::endl;
403#endif
404 std::vector<int> aLabels;
405 for (int i=0; i<nlabels; i++) {
406 aLabels.push_back(labels[i]);
407 }
408 Map<Eigen::MatrixXd> derivatives(derArray,1,nlabels);
409 self->addGlobals(aLabels, derivatives);
410}
411
420void GblPoint_getGlobalLabelsAndDerivatives(GblPoint* self, int* nlabels, int** labels, double** ders) {
421#ifdef JNA_DEBUG
422 std::cout << "GblPoint_getGlobalLabelsAndDerivatives("
423 << self << ", " << labels << ", " << ders
424 << ")" << std::endl;
425#endif
426 std::vector<int> glabels;
427 std::vector<double> gders;
428
429#ifdef JNA_DEBUG
430 std::cout << " GblPoint has " << std::flush
431 << self->getMeasEnd() - self->getMeasBegin()
432 << " measurements." << std::endl;
433#endif
434
435 //Should I add the number of derivatives? - Row/Col? CHECK CHECK CHECK
437 0 /* aMeas */, 0 /* aRow */,
438 glabels, gders);
439
440#ifdef JNA_DEBUG
441 std::cout <<" Aquired " << glabels.size() << " labels and " << gders.size() << " derivatives." <<std::endl;
442#endif
443
444 *nlabels = glabels.size();
445 *labels = new int[*nlabels];
446 *ders = new double[*nlabels];
447
448#ifdef JNA_DEBUG
449 std::cout <<" Allocated return arrays." << std::endl;
450#endif
451
452
453 for (std::size_t il{0}; il < glabels.size(); ++il) {
454 (*labels)[il] = glabels.at(il);
455 //std::cout<<glabels.at(il)<<std::endl;
456 }
457
458 //std::cout<<"GblPointWrapper::gders"<<std::endl;
459
460 for (std::size_t id{0}; id < glabels.size(); ++id) {
461 (*ders)[id] = gders.at(id);
462 //std::cout<<gders.at(il)<<std::endl;
463 }
464
465 // set array using Eigen
466 //Map<MatrixXd>(ders,1,gders.size()) = gders;
467#ifdef JNA_DEBUG
468 std::cout << " Done with assignment and leaving." << std::endl;
469#endif
470}
471
487 int flagCurv, int flagU1dir, int flagU2dir) {
488#ifdef JNA_DEBUG
489 std::cout << "GblTracjectoryCtorPtrArray("
490 << points << ", " << npoints << ", "
491 << flagCurv << ", " << flagU1dir << ", " << flagU2dir
492 << ")" << std::endl;
493#endif
494#if JNA_DO_MONITOR
495 ++num_gbl_traj;
496#endif
497
498 return new GblTrajectory(ptr_array_to_vector(points, npoints),
499 flagCurv!=0, flagU1dir!=0, flagU2dir!=0);
500}
501
515 int aLabel, double seedArray[],
516 int flagCurv, int flagU1dir, int flagU2dir) {
517#ifdef JNA_DEBUG
518 std::cout << "GblTrajectoryCtorPtrArraySeed("
519 << points << ", " << npoints << ", "
520 << aLabel << ", " << seedArray << ", "
521 << flagCurv << ", " << flagU1dir << ", " << flagU2dir
522 << ")" << std::endl;
523#endif
524#if JNA_DO_MONITOR
525 ++num_gbl_traj;
526#endif
527
528 Map<Matrix5d> seed(seedArray,5,5);
529
530 return new GblTrajectory(ptr_array_to_vector(points, npoints), aLabel, seed,
531 flagCurv!=0, flagU1dir!=0, flagU2dir!=0);
532}
533
545GblTrajectory* GblTrajectoryCtorPtrComposed(GblPoint* points_1[], int npoints_1, double trafo_1[],
546 GblPoint* points_2[], int npoints_2, double trafo_2[]) {
547#ifdef JNA_DEBUG
548 std::cout << "GblTrajectoryCtorPtrComposed("
549 << points_1 << ", " << npoints_1 << ", " << trafo_1 << ", "
550 << points_2 << ", " << npoints_2 << ", " << trafo_2 << ")"
551 << std::endl;
552#endif
553#if JNA_DO_MONITOR
554 ++num_gbl_traj;
555#endif
556
557
558 MatrixXd inner_1(2,3);
559 inner_1(0,0)=trafo_1[0];
560 inner_1(0,1)=trafo_1[1];
561 inner_1(0,2)=trafo_1[2];
562
563 inner_1(1,0)=trafo_1[3];
564 inner_1(1,1)=trafo_1[4];
565 inner_1(1,2)=trafo_1[5];
566
567 std::pair<std::vector<GblPoint>, MatrixXd> track_trafo_1 = std::make_pair(ptr_array_to_vector(points_1, npoints_1), inner_1);
568
569 //second track
570
571 MatrixXd inner_2(2,3);
572 inner_2(0,0)=trafo_2[0];
573 inner_2(0,1)=trafo_2[1];
574 inner_2(0,2)=trafo_2[2];
575
576 inner_2(1,0)=trafo_2[3];
577 inner_2(1,1)=trafo_2[4];
578 inner_2(1,2)=trafo_2[5];
579
580 std::pair<std::vector<GblPoint>, MatrixXd> track_trafo_2 = std::make_pair(ptr_array_to_vector(points_2, npoints_2), inner_2);
581
582 return new GblTrajectory({track_trafo_1, track_trafo_2});
583
584}
585
595void GblTrajectory_fit(GblTrajectory* self, double* Chi2, int* Ndf, double* lostWeight, char* c_optionList, unsigned int aLabel) {
596#ifdef JNA_DEBUG
597 std::cout << "GblTrajectory_fit("
598 << self << ", " << Chi2 << ", " << Ndf << ", " << lostWeight << ", " << c_optionList << ", " << aLabel << ")"
599 << std::endl;
600#endif
601
602 std::string optionList(c_optionList);
603 self->fit(*Chi2, *Ndf, *lostWeight, optionList,aLabel);
604}
605
612#ifdef JNA_DEBUG
613 std::cout << "GblTrajectory_delete(" << self << ")" << std::endl;
614#endif
615#if JNA_DO_MONITOR
616 if (self) {
617 --num_gbl_traj;
618 num_gbl_point -= self->getNumPoints();
619 }
620#endif
621 if (self) delete self;
622}
623
630#ifdef JNA_DEBUG
631 std::cout << "GblTrajectory_isValid(" << self << ")" << std::endl;
632#endif
633 return self->isValid() ? 0 : 1;
634}
635
642#ifdef JNA_DEBUG
643 std::cout << "GblTrajectory_getNumPoints(" << self << ")" << std::endl;
644#endif
645 return (int) self->getNumPoints();
646}
647
654#ifdef JNA_DEBUG
655 std::cout << "GblTrajectory_printTrajectory(" << self << ", " << level << ")" << std::endl;
656#endif
657 return self->printTrajectory(level);
658}
659
665#ifdef JNA_DEBUG
666 std::cout << "GblTrajectory_printData(" << self << ")" << std::endl;
667#endif
668 return self->printData();
669}
670
677#ifdef JNA_DEBUG
678 std::cout << "GblTrajectory_printPoints(" << self << ", " << level << ")" << std::endl;
679#endif
680 return self->printPoints(level);
681}
682
693void GblTrajectory_getResults(GblTrajectory* self, int aSignedLabel, double* localPar, int* nLocalPar,
694 double * localCov, int* sizeLocalCov) {
695#ifdef JNA_DEBUG
696 std::cout << "GblTrajectory_getResults(" << self
697 << ", " << aSignedLabel << ", " << localPar << ", " << nLocalPar
698 << ", " << localCov << ", " << sizeLocalCov << ")" << std::endl;
699#endif
700
701 Eigen::VectorXd e_localPar(5);
702 Eigen::MatrixXd e_localCov(5,5);
703
704 self->getResults(aSignedLabel, e_localPar, e_localCov);
705
706 //std::cout<<"gblTrajectoryWrapper::getResults"<<std::endl;
707 //std::cout<<e_localPar<<std::endl;
708 //std::cout<<e_localCov<<std::endl;
709
710 Map<Vector5d>(localPar, 5) = e_localPar;
711 Map<Matrix5d>(localCov, 5, 5) = e_localCov;
712 *nLocalPar = 5;
713 *sizeLocalCov = 5;
714}
715
727int GblTrajectory_getMeasResults(GblTrajectory* self, int aLabel, int* numData,
728 double* aResiduals, double* aMeasErrors, double* aResErrors,
729 double* aDownWeights) {
730#ifdef JNA_DEBUG
731 std::cout << "GblTrajectory_getMeasResults(" << self
732 << ", " << aLabel << ", " << numData << ", " << aResiduals
733 << ", " << aMeasErrors << ", " << aResErrors << ", " << aDownWeights << ")" << std::endl;
734#endif
735
736 Eigen::VectorXd e_aResiduals(2);
737 Eigen::VectorXd e_aMeasErrors(2);
738 Eigen::VectorXd e_aResErrors(2);
739 Eigen::VectorXd e_aDownWeights(2);
740 unsigned int num_data = 0;
741
742 unsigned int out = self->getMeasResults(aLabel, num_data, e_aResiduals, e_aMeasErrors,
743 e_aResErrors, e_aDownWeights);
744#ifdef JNA_DEBUG
745 std::cout << " getMeasResults returned the status code " << out << std::endl;
746#endif
747
748 *numData = num_data;
749
750 for (unsigned int i = 0; i < num_data; i++) {
751 aResiduals[i] = e_aResiduals(i);
752 aMeasErrors[i] = e_aMeasErrors(i);
753 aResErrors[i] = e_aResErrors(i);
754 aDownWeights[i] = e_aDownWeights(i);
755 }
756
757 return out;
758}
759
767#ifdef JNA_DEBUG
768 std::cout << "GblTrajectory_milleOut(" << self
769 << ", " << millebinary << ")" << std::endl;
770#endif
771 self->milleOut(*millebinary);
772}
773
788GblDetectorLayer* GblDetectorLayerCtor(const char* aName, int aLayer, int aDim, double thickness,
789 double aCenter[], double aResolution[], double aPrecision[],
790 double measTrafo[], double alignTrafo[]) {
791#ifdef JNA_DEBUG
792 std::cout << "GblDetectorLayerCtor(" << aName << ", " << aLayer
793 << ", " << aDim << ", " << thickness << ", " << aCenter << ", " << aResolution
794 << ", " << aPrecision << ", " << measTrafo << ", " << alignTrafo << ")" << std::endl;
795#endif
796#if JNA_DO_MONITOR
797 ++num_gbl_det_layer;
798#endif
799
800 // uses eigen's Map structure to decompose an array into our type of Vector/Matrix
801 // need to then copy that object into a local variable since Map's only allow const
802 // reference access
803 Vector3d e_aCenter(Map<Vector3d>(aCenter,3));
804 Vector2d e_aResolution(Map<Vector2d>(aResolution,2));
805 Vector2d e_aPrecision(Map<Vector2d>(aPrecision,2));
806 Matrix3d e_measTrafo(Map<Matrix3d>(measTrafo,3,3));
807 Matrix3d e_alignTrafo(Map<Matrix3d>(alignTrafo,3,3));
808
809 GblDetectorLayer* layer = new GblDetectorLayer(aName, aLayer, aDim, thickness, e_aCenter, e_aResolution, e_aPrecision, e_measTrafo, e_alignTrafo);
810 return layer;
811}
812
817#ifdef JNA_DEBUG
818 std::cout << "GblDetectorLayer_delete(" << self << ")" << std::endl;
819#endif
820#if JNA_DO_MONITOR
821 --num_gbl_det_layer;
822#endif
823 if (self) delete self;
824}
825
831#ifdef JNA_DEBUG
832 std::cout << "GblDetectorLayer_print(" << self << ")" << std::endl;
833#endif
834 self->print();
835}
836
842#ifdef JNA_DEBUG
843 std::cout << "GblDetectorLayer_getRadiationLength(" << self << ")" << std::endl;
844#endif
845 return self->getRadiationLength();
846}
847
848//TODO Implement wrappers for
849//getResolution
850//getPrecision
851//getCenter
852//getMeasSystemDirs
853//getAlignSystemDirs
854
855
856//Helix prediction on layer
857
865#ifdef JNA_DEBUG
866 std::cout << "GblDetectorLayer_intersectWithHelix(" << self << ", " << hlx << ")" << std::endl;
867#endif
868#if JNA_DO_MONITOR
869 ++num_gbl_helix_prediction;
870#endif
871
872 Vector3d center = self->getCenter();
873 Vector3d udir = (self->getMeasSystemDirs()).row(0);
874 Vector3d vdir = (self->getMeasSystemDirs()).row(1);
875
876 return new GblHelixPrediction(hlx->getPrediction(center, udir, vdir));
877}
878
887GblSimpleHelix* GblSimpleHelixCtor(double aRinv, double aPhi0, double aDca, double aDzds, double aZ0) {
888#if JNA_DO_MONITOR
889 ++num_gbl_simple_helix;
890#endif
891 return new GblSimpleHelix(aRinv, aPhi0, aDca, aDzds, aZ0);
892}
893
898#if JNA_DO_MONITOR
899 --num_gbl_simple_helix;
900#endif
901 if (self) delete self;
902}
903
910double GblSimpleHelix_getPhi(GblSimpleHelix* self, double aRadius) {
911 return self->getPhi(aRadius);
912}
913
920double GblSimpleHelix_getArcLengthR(GblSimpleHelix* self, double aRadius) {
921 return self->getArcLengthR(aRadius);
922}
923
931double GblSimpleHelix_getArcLengthXY(GblSimpleHelix* self, double xPos, double yPos) {
932 return self->getArcLengthXY(xPos, yPos);
933}
934
944void GblSimpleHelix_moveToXY(GblSimpleHelix* self, double xPos, double yPos,
945 double* newPhi0, double* newDca, double* newZ0) {
946 self->moveToXY(xPos, yPos,
947 *newPhi0, *newDca, *newZ0);
948}
949
959GblHelixPrediction* GblSimpleHelix_getPrediction(GblSimpleHelix* self, double refPos[], double uDir[], double vDir[]) {
960
961 Map<Vector3d> e_refPos(refPos,3);
962 Map<Vector3d> e_uDir(uDir,3);
963 Map<Vector3d> e_vDir(vDir,3);
964 GblHelixPrediction prediction = self->getPrediction(e_refPos,e_uDir,e_vDir);
965
966 /*std::cout<<"Cross Check GBL predicted position!"<<std::endl;
967 std::cout<<prediction->getPosition()<<std::endl;
968
969 std::cout<<"Cross Check GBL meas predicted position!"<<std::endl;
970 std::cout<<prediction->getMeasPred()<<std::endl;*/
971
972#if JNA_DO_MONITOR
973 ++num_gbl_helix_prediction;
974#endif
975 // JNA only deals with pointers so we need to dynamically create a new copy
976 return new GblHelixPrediction(prediction);
977}
978
989GblHelixPrediction* GblHelixPredictionCtor(double sArc, double aPred[], double tDir[], double uDir[], double vDir[],
990 double nDir[], double aPos[]) {
991 Map<Vector2d> e_aPred(aPred,2);
992 Map<Vector3d> e_tDir(tDir,3);
993 Map<Vector3d> e_uDir(uDir,3);
994 Map<Vector3d> e_vDir(vDir,3);
995 Map<Vector3d> e_nDir(nDir,3);
996 Map<Vector3d> e_aPos(aPos,3);
997
998
999#if JNA_DO_MONITOR
1000 ++num_gbl_helix_prediction;
1001#endif
1002 return new GblHelixPrediction(sArc, e_aPred, e_tDir, e_uDir, e_vDir,
1003 e_nDir, e_aPos);
1004}
1005
1011#if JNA_DO_MONITOR
1012 --num_gbl_helix_prediction;
1013#endif
1014 if (self) delete self;
1015}
1016
1024 return self->getArcLength();
1025}
1026
1034 Vector2d e_pred = self->getMeasPred();
1035
1036 prediction[0] = e_pred(0);
1037 prediction[1] = e_pred(1);
1038}
1039
1047 Vector3d e_pos = self->getPosition();
1048
1049 position[0] = e_pos(0);
1050 position[1] = e_pos(1);
1051 position[2] = e_pos(2);
1052}
1053
1061 Vector3d e_dir = self->getDirection();
1062
1063 direction[0] = e_dir(0);
1064 direction[1] = e_dir(1);
1065 direction[2] = e_dir(2);
1066}
1067
1075 return self->getCosIncidence();
1076}
1077
1085 Matrix<double,2,3> curDirs = self->getCurvilinearDirs();
1086
1087 curvilinear[0] = curDirs(0,0);
1088 curvilinear[1] = curDirs(0,1);
1089 curvilinear[2] = curDirs(0,2);
1090 curvilinear[3] = curDirs(1,0);
1091 curvilinear[4] = curDirs(1,1);
1092 curvilinear[5] = curDirs(1,2);
1093}
1094}
1095
GblPoint definition.
GblTrajectory definition.
Definitions for GBL utilities.
void GblTrajectory_printTrajectory(GblTrajectory *self, int level)
call gbl::GblTrajectory::printTrajectory
void GblHelixPrediction_getPosition(GblHelixPrediction *self, double *position)
get the position
void GblHelixPrediction_getCurvilinearDirs(GblHelixPrediction *self, double curvilinear[])
get the curvilinear directions
GblHelixPrediction * GblDetectorLayer_intersectWithHelix(GblDetectorLayer *self, GblSimpleHelix *hlx)
construct a helix prediction from a layer
void GblHelixPrediction_delete(GblHelixPrediction *self)
delete a gbl::GblHelixPrediction
double GblSimpleHelix_getPhi(GblSimpleHelix *self, double aRadius)
call gbl::GblSimpleHelix::getPhi
GblTrajectory * GblTrajectoryCtorPtrComposed(GblPoint *points_1[], int npoints_1, double trafo_1[], GblPoint *points_2[], int npoints_2, double trafo_2[])
compose trajectory for a 2-body decay
const int NCOL
void GblTrajectory_getResults(GblTrajectory *self, int aSignedLabel, double *localPar, int *nLocalPar, double *localCov, int *sizeLocalCov)
get trajectory results
void MilleBinary_close(MilleBinary *self)
Closing a gbl::MilleBinary file is the same as destructing it.
double GblSimpleHelix_getArcLengthXY(GblSimpleHelix *self, double xPos, double yPos)
call gbl::GblSimpleHelix::getArcLengthXY
void GblTrajectory_delete(GblTrajectory *self)
delete self
int GblPoint_getNumMeasurements(GblPoint *self)
calculate number of measurements in a gbl::GblPoint
void GblHelixPrediction_getDirection(GblHelixPrediction *self, double direction[])
get the direction
MilleBinary * MilleBinaryCtor(const char *fileName, int filenamesize, int doublePrecision, int keepZeros, int aSize)
Dynamically create new MilleBinary file.
int GblTrajectory_getNumPoints(GblTrajectory *self)
call gbl::GblTrajectory::getNumPoints
GblTrajectory * GblTrajectoryCtorPtrArraySeed(GblPoint *points[], int npoints, int aLabel, double seedArray[], int flagCurv, int flagU1dir, int flagU2dir)
construct new gbl::GblTrajectory with a seed matrix
std::vector< GblPoint > ptr_array_to_vector(GblPoint *points[], int npoints)
convert the pointer array of gbl::GblPoint into a vector holding the objects
GblTrajectory * GblTrajectoryCtorPtrArray(GblPoint *points[], int npoints, int flagCurv, int flagU1dir, int flagU2dir)
simple construction of a new gbl::GblTrajectory
double GblHelixPrediction_getArcLength(GblHelixPrediction *self)
get the arc length of a helix prediction
int GblTrajectory_isValid(GblTrajectory *self)
call gbl::GblTrajectory::isValid
void GblSimpleHelix_moveToXY(GblSimpleHelix *self, double xPos, double yPos, double *newPhi0, double *newDca, double *newZ0)
call gbl::GblSimpleHelix::moveToXY
GblHelixPrediction * GblHelixPredictionCtor(double sArc, double aPred[], double tDir[], double uDir[], double vDir[], double nDir[], double aPos[])
create a new helix prediction manually
void GblDetectorLayer_print(GblDetectorLayer *self)
call gbl::GblDetectorLayer::print
GblHelixPrediction * GblSimpleHelix_getPrediction(GblSimpleHelix *self, double refPos[], double uDir[], double vDir[])
get a helical prediction from a reference coordinate system
void GblPoint_addGlobals(GblPoint *self, int *labels, int nlabels, double *derArray)
add global derivatives to the point
void GblPoint_delete(GblPoint *self)
delete gbl::GblPoint
void GblTrajectory_milleOut(GblTrajectory *self, MilleBinary *millebinary)
write a trajectory to a gbl::MilleBinary file
void GblPoint_addScatterer(GblPoint *self, double *resArray, double *precArray)
add a vector-precision scatterer
double GblSimpleHelix_getArcLengthR(GblSimpleHelix *self, double aRadius)
call gbl::GblSimpleHelix::getArcLengthR
GblSimpleHelix * GblSimpleHelixCtor(double aRinv, double aPhi0, double aDca, double aDzds, double aZ0)
construct a new simple helix
GblPoint * GblPointCtor(double matrixArray[NROW *NCOL])
create new gbl::GblPoint from a jacobian
void GblHelixPrediction_getMeasPred(GblHelixPrediction *self, double *prediction)
get the predicted measurement
double GblHelixPrediction_getCosIncidence(GblHelixPrediction *self)
get the cosine incidence
GblDetectorLayer * GblDetectorLayerCtor(const char *aName, int aLayer, int aDim, double thickness, double aCenter[], double aResolution[], double aPrecision[], double measTrafo[], double alignTrafo[])
construct a detector layer
void GblTrajectory_fit(GblTrajectory *self, double *Chi2, int *Ndf, double *lostWeight, char *c_optionList, unsigned int aLabel)
call gbl::GblTrajectory::fit on self
void GblSimpleHelix_delete(GblSimpleHelix *self)
delete gbl::GblSimpleHelix
void GblDetectorLayer_delete(GblDetectorLayer *self)
delete gbl::GblDetectorLayer
void GblPoint_getGlobalLabelsAndDerivatives(GblPoint *self, int *nlabels, int **labels, double **ders)
get global derivatives and their labels from the point
void GblTrajectory_printPoints(GblTrajectory *self, int level)
call gbl::GblTrajectory::printPoints
void GblPoint_printPoint(const GblPoint *self, unsigned int level)
call gbl::GblPoint::printPoint on self
void GblPoint_addMeasurement2D(GblPoint *self, double *projArray, double *resArray, double *precArray, double minPrecision)
add a 2D measurement
double GblDetectorLayer_getRadiationLength(GblDetectorLayer *self)
call gbl::GblDetectorLayer::getRadiationLength
void GblTrajectory_printData(GblTrajectory *self)
call gbl::GblTrajectory::printData
int GblTrajectory_getMeasResults(GblTrajectory *self, int aLabel, int *numData, double *aResiduals, double *aMeasErrors, double *aResErrors, double *aDownWeights)
get measurement results from trajectory
const int NROW
MilleBinary definition.
Detector layer.
Definition: GblUtilities.h:107
Eigen::Vector3d getCenter() const
Get center.
double getRadiationLength() const
Get radiation length.
void print() const
Print GblDetectorLayer.
Eigen::Matrix3d getMeasSystemDirs() const
Get directions of measurement system.
Prediction on helix.
Definition: GblUtilities.h:49
Eigen::Matrix< double, 2, 3 > getCurvilinearDirs() const
Get curvilinear directions (U,V)
double getArcLength() const
Get arc-length.
const Eigen::Vector2d & getMeasPred() const
Get (measurement) prediction.
double getCosIncidence() const
Get cosine of incidence.
const Eigen::Vector3d & getPosition() const
Get position.
const Eigen::Vector3d & getDirection() const
Get position.
Point on trajectory.
Definition: GblPoint.h:62
void addMeasurement(const Eigen::MatrixBase< Projection > &aProjection, const Eigen::MatrixBase< Residuals > &aResiduals, const Eigen::MatrixBase< Precision > &aPrecision, double minPrecision=0.)
Add a measurement to a point.
Definition: GblPoint.h:238
void getGlobalLabelsAndDerivatives(unsigned int aMeas, unsigned int aRow, std::vector< int > &aLabels, std::vector< double > &aDerivatives) const
Retrieve global derivatives from a measurement at a point for a single row.
Definition: GblPoint.cpp:625
void printPoint(unsigned int level=0) const
Print GblPoint.
Definition: GblPoint.cpp:554
std::vector< GblMeasurement >::iterator getMeasBegin()
Get GblMeasurement iterator for begin.
Definition: GblPoint.cpp:609
void addGlobals(const std::vector< int > &aLabels, const Eigen::MatrixBase< Derivative > &aDerivatives)
Add global derivatives to a point.
Definition: GblPoint.h:293
std::vector< GblMeasurement >::iterator getMeasEnd()
Get GblMeasurement iterator for end.
Definition: GblPoint.cpp:614
void addScatterer(const Eigen::Vector2d &aResiduals, const Eigen::Matrix2d &aPrecision)
Add a thin scatterer to a point.
Definition: GblPoint.cpp:159
Simple helix.
Definition: GblUtilities.h:78
double getArcLengthR(double aRadius) const
Get (2D) arc length for given radius (to ref. point)
void moveToXY(double xPos, double yPos, double &newPhi0, double &newDca, double &newZ0) const
Move to new reference point (X,Y)
GblHelixPrediction getPrediction(const Eigen::Vector3d &refPos, const Eigen::Vector3d &uDir, const Eigen::Vector3d &vDir) const
Get prediction.
double getPhi(double aRadius) const
Get phi (of point on circle) for given radius (to ref. point)
double getArcLengthXY(double xPos, double yPos) const
Get (2D) arc length for given point.
GBL trajectory.
Definition: GblTrajectory.h:50
void printPoints(unsigned int level=0) const
Print GblPoints on trajectory.
unsigned int getMeasResults(unsigned int aLabel, unsigned int &numData, Eigen::VectorXd &aResiduals, Eigen::VectorXd &aMeasErrors, Eigen::VectorXd &aResErrors, Eigen::VectorXd &aDownWeights)
Get residuals from fit at point for measurement "long list".
void printTrajectory(unsigned int level=0) const
Print GblTrajectory.
void milleOut(MilleBinary &aMille)
Write valid trajectory to Millepede-II binary file.
bool isValid() const
Retrieve validity of trajectory.
unsigned int getResults(int aSignedLabel, Eigen::VectorXd &localPar, Eigen::MatrixXd &localCov) const
Get fit results at point.
unsigned int fit(double &Chi2, int &Ndf, double &lostWeight, const std::string &optionList="", unsigned int aLabel=0)
Perform fit of (valid) trajectory.
void printData() const
Print GblData blocks for trajectory.
unsigned int getNumPoints() const
Retrieve number of point from trajectory.
Millepede-II (binary) record.
Definition: MilleBinary.h:81
Namespace for the general broken lines package.