GeneralBrokenLines V03-01-01
using EIGEN
Macros | Functions | Variables
JavaNativeAccessWrappers.cpp File Reference
#include "MilleBinary.h"
#include "GblPoint.h"
#include "GblTrajectory.h"
#include "GblUtilities.h"
#include "Eigen/Core"

Go to the source code of this file.

Macros

#define JNA_DO_MONITOR   0
 Do the memory monitoring if either JNA_DEBUG or JNA_MEMORY_MONITOR is defined. More...
 

Functions

std::vector< GblPointptr_array_to_vector (GblPoint *points[], int npoints)
 convert the pointer array of gbl::GblPoint into a vector holding the objects More...
 
MilleBinaryMilleBinaryCtor (const char *fileName, int filenamesize, int doublePrecision, int keepZeros, int aSize)
 Dynamically create new MilleBinary file. More...
 
void MilleBinary_close (MilleBinary *self)
 Closing a gbl::MilleBinary file is the same as destructing it. More...
 
GblPointGblPointCtor (double matrixArray[NROW *NCOL])
 create new gbl::GblPoint from a jacobian More...
 
void GblPoint_delete (GblPoint *self)
 delete gbl::GblPoint More...
 
void GblPoint_printPoint (const GblPoint *self, unsigned int level)
 call gbl::GblPoint::printPoint on self More...
 
int GblPoint_getNumMeasurements (GblPoint *self)
 calculate number of measurements in a gbl::GblPoint More...
 
void GblPoint_addMeasurement2D (GblPoint *self, double *projArray, double *resArray, double *precArray, double minPrecision)
 add a 2D measurement More...
 
void GblPoint_addScatterer (GblPoint *self, double *resArray, double *precArray)
 add a vector-precision scatterer More...
 
void GblPoint_addGlobals (GblPoint *self, int *labels, int nlabels, double *derArray)
 add global derivatives to the point More...
 
void GblPoint_getGlobalLabelsAndDerivatives (GblPoint *self, int *nlabels, int **labels, double **ders)
 get global derivatives and their labels from the point More...
 
GblTrajectoryGblTrajectoryCtorPtrArray (GblPoint *points[], int npoints, int flagCurv, int flagU1dir, int flagU2dir)
 simple construction of a new gbl::GblTrajectory More...
 
GblTrajectoryGblTrajectoryCtorPtrArraySeed (GblPoint *points[], int npoints, int aLabel, double seedArray[], int flagCurv, int flagU1dir, int flagU2dir)
 construct new gbl::GblTrajectory with a seed matrix More...
 
GblTrajectoryGblTrajectoryCtorPtrComposed (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 More...
 
void GblTrajectory_fit (GblTrajectory *self, double *Chi2, int *Ndf, double *lostWeight, char *c_optionList, unsigned int aLabel)
 call gbl::GblTrajectory::fit on self More...
 
void GblTrajectory_delete (GblTrajectory *self)
 delete self More...
 
int GblTrajectory_isValid (GblTrajectory *self)
 call gbl::GblTrajectory::isValid More...
 
int GblTrajectory_getNumPoints (GblTrajectory *self)
 call gbl::GblTrajectory::getNumPoints More...
 
void GblTrajectory_printTrajectory (GblTrajectory *self, int level)
 call gbl::GblTrajectory::printTrajectory More...
 
void GblTrajectory_printData (GblTrajectory *self)
 call gbl::GblTrajectory::printData More...
 
void GblTrajectory_printPoints (GblTrajectory *self, int level)
 call gbl::GblTrajectory::printPoints More...
 
void GblTrajectory_getResults (GblTrajectory *self, int aSignedLabel, double *localPar, int *nLocalPar, double *localCov, int *sizeLocalCov)
 get trajectory results More...
 
int GblTrajectory_getMeasResults (GblTrajectory *self, int aLabel, int *numData, double *aResiduals, double *aMeasErrors, double *aResErrors, double *aDownWeights)
 get measurement results from trajectory More...
 
void GblTrajectory_milleOut (GblTrajectory *self, MilleBinary *millebinary)
 write a trajectory to a gbl::MilleBinary file More...
 
GblDetectorLayerGblDetectorLayerCtor (const char *aName, int aLayer, int aDim, double thickness, double aCenter[], double aResolution[], double aPrecision[], double measTrafo[], double alignTrafo[])
 construct a detector layer More...
 
void GblDetectorLayer_delete (GblDetectorLayer *self)
 delete gbl::GblDetectorLayer More...
 
void GblDetectorLayer_print (GblDetectorLayer *self)
 call gbl::GblDetectorLayer::print More...
 
double GblDetectorLayer_getRadiationLength (GblDetectorLayer *self)
 call gbl::GblDetectorLayer::getRadiationLength More...
 
GblHelixPredictionGblDetectorLayer_intersectWithHelix (GblDetectorLayer *self, GblSimpleHelix *hlx)
 construct a helix prediction from a layer More...
 
GblSimpleHelixGblSimpleHelixCtor (double aRinv, double aPhi0, double aDca, double aDzds, double aZ0)
 construct a new simple helix More...
 
void GblSimpleHelix_delete (GblSimpleHelix *self)
 delete gbl::GblSimpleHelix More...
 
double GblSimpleHelix_getPhi (GblSimpleHelix *self, double aRadius)
 call gbl::GblSimpleHelix::getPhi More...
 
double GblSimpleHelix_getArcLengthR (GblSimpleHelix *self, double aRadius)
 call gbl::GblSimpleHelix::getArcLengthR More...
 
double GblSimpleHelix_getArcLengthXY (GblSimpleHelix *self, double xPos, double yPos)
 call gbl::GblSimpleHelix::getArcLengthXY More...
 
void GblSimpleHelix_moveToXY (GblSimpleHelix *self, double xPos, double yPos, double *newPhi0, double *newDca, double *newZ0)
 call gbl::GblSimpleHelix::moveToXY More...
 
GblHelixPredictionGblSimpleHelix_getPrediction (GblSimpleHelix *self, double refPos[], double uDir[], double vDir[])
 get a helical prediction from a reference coordinate system More...
 
GblHelixPredictionGblHelixPredictionCtor (double sArc, double aPred[], double tDir[], double uDir[], double vDir[], double nDir[], double aPos[])
 create a new helix prediction manually More...
 
void GblHelixPrediction_delete (GblHelixPrediction *self)
 delete a gbl::GblHelixPrediction More...
 
double GblHelixPrediction_getArcLength (GblHelixPrediction *self)
 get the arc length of a helix prediction More...
 
void GblHelixPrediction_getMeasPred (GblHelixPrediction *self, double *prediction)
 get the predicted measurement More...
 
void GblHelixPrediction_getPosition (GblHelixPrediction *self, double *position)
 get the position More...
 
void GblHelixPrediction_getDirection (GblHelixPrediction *self, double direction[])
 get the direction More...
 
double GblHelixPrediction_getCosIncidence (GblHelixPrediction *self)
 get the cosine incidence More...
 
void GblHelixPrediction_getCurvilinearDirs (GblHelixPrediction *self, double curvilinear[])
 get the curvilinear directions More...
 

Variables

const int NROW = 5
 
const int NCOL = 5
 

Detailed Description

Author
Tom Eichlersmith eichl.nosp@m.008@.nosp@m.umn.e.nosp@m.du

Wrappers around construction, getters, setters, and destruction functions of Gbl classes.

These wrappers are put into an extern "C" block so that there names are not mangled and therefore can be accessed by name within a Java Native Access (JNA) class accessing the GBL library.

Note
This file is large however it is not intended to do anything besides interface between JNA and GBL. If any logic in this file affects the calculations done by GBL, that is considered a bug with the implementation.

Information on using these wrappers to access the GBL C++ library from with java using JNA can be found in on the JNA Usage page.

Definition in file JavaNativeAccessWrappers.cpp.

Macro Definition Documentation

◆ JNA_DO_MONITOR

#define JNA_DO_MONITOR   0

Do the memory monitoring if either JNA_DEBUG or JNA_MEMORY_MONITOR is defined.

Definition at line 141 of file JavaNativeAccessWrappers.cpp.

Function Documentation

◆ GblDetectorLayer_delete()

void GblDetectorLayer_delete ( GblDetectorLayer self)

delete gbl::GblDetectorLayer

Definition at line 816 of file JavaNativeAccessWrappers.cpp.

◆ GblDetectorLayer_getRadiationLength()

double GblDetectorLayer_getRadiationLength ( GblDetectorLayer self)

◆ GblDetectorLayer_intersectWithHelix()

GblHelixPrediction * GblDetectorLayer_intersectWithHelix ( GblDetectorLayer self,
GblSimpleHelix hlx 
)

construct a helix prediction from a layer

Parameters
[in]selfgbl::GblDetectorLayer to operate on
[in]hlxgbl::GblSimpleHelix to do helical calculations
Returns
newly constructed gbl::GblHelixPrediction

Definition at line 864 of file JavaNativeAccessWrappers.cpp.

References gbl::GblDetectorLayer::getCenter(), gbl::GblDetectorLayer::getMeasSystemDirs(), and gbl::GblSimpleHelix::getPrediction().

◆ GblDetectorLayer_print()

void GblDetectorLayer_print ( GblDetectorLayer self)

call gbl::GblDetectorLayer::print

Parameters
[in]selfgbl::GblDetectorLayer to operate on

Definition at line 830 of file JavaNativeAccessWrappers.cpp.

References gbl::GblDetectorLayer::print().

◆ GblDetectorLayerCtor()

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

See also
gbl::GblDetectorLayer::GblDetectorLayer
Parameters
[in]aNameC-style string name
[in]aLayerinteger ID for layer
[in]aDimdimension of layer
[in]thicknessthickness of alyer
[in]aCenterdouble array of length 3 defining center of layer
[in]aResolutiondouble array of length 2
[in]aPrecisiondouble array of length 2
[in]measTrafodouble array defining 3x3 matrix row wise
[in]alignTrafodouble array defining 3x3 matrix row wise
Returns
newly constructed detector layer

Definition at line 788 of file JavaNativeAccessWrappers.cpp.

◆ GblHelixPrediction_delete()

void GblHelixPrediction_delete ( GblHelixPrediction self)

delete a gbl::GblHelixPrediction

Parameters
[in]selfgbl::GblHelixPrediction to delete

Definition at line 1010 of file JavaNativeAccessWrappers.cpp.

◆ GblHelixPrediction_getArcLength()

double GblHelixPrediction_getArcLength ( GblHelixPrediction self)

get the arc length of a helix prediction

See also
gbl::GblHelixPrediction::getArcLength
Parameters
[in]selfgbl::GblHelixPrediction to operate on
Returns
arc length

Definition at line 1023 of file JavaNativeAccessWrappers.cpp.

References gbl::GblHelixPrediction::getArcLength().

◆ GblHelixPrediction_getCosIncidence()

double GblHelixPrediction_getCosIncidence ( GblHelixPrediction self)

get the cosine incidence

See also
gbl::GblHelixPrediction::getCosIncidence
Parameters
[in]selfgbl::GblHelixPrediction to operate on
Returns
value of cosine incidence

Definition at line 1074 of file JavaNativeAccessWrappers.cpp.

References gbl::GblHelixPrediction::getCosIncidence().

◆ GblHelixPrediction_getCurvilinearDirs()

void GblHelixPrediction_getCurvilinearDirs ( GblHelixPrediction self,
double  curvilinear[] 
)

get the curvilinear directions

See also
gbl::GblHelixPrediction::getCurvilinearDirs
Parameters
[in]selfgbl::GblHelixPrediction to operate on
[out]curvilinearlength-6 double array that will hold the two direction vectors

Definition at line 1084 of file JavaNativeAccessWrappers.cpp.

References gbl::GblHelixPrediction::getCurvilinearDirs().

◆ GblHelixPrediction_getDirection()

void GblHelixPrediction_getDirection ( GblHelixPrediction self,
double  direction[] 
)

get the direction

See also
gbl::GblHelixPrediction::getDirection
Parameters
[in]selfgbl::GblHelixPrediction to operate on
[out]directionlength-3 double array that will hold direction

Definition at line 1060 of file JavaNativeAccessWrappers.cpp.

References gbl::GblHelixPrediction::getDirection().

◆ GblHelixPrediction_getMeasPred()

void GblHelixPrediction_getMeasPred ( GblHelixPrediction self,
double *  prediction 
)

get the predicted measurement

See also
gbl::GblHelixPrediction::getMeadPred
Parameters
[in]selfgbl::GblHelixPrediction to operate on
[out]predictionlength-2 double array that will hold predicted measurement

Definition at line 1033 of file JavaNativeAccessWrappers.cpp.

References gbl::GblHelixPrediction::getMeasPred().

◆ GblHelixPrediction_getPosition()

void GblHelixPrediction_getPosition ( GblHelixPrediction self,
double *  position 
)

get the position

See also
gbl::GblHelixPrediction::getPosition
Parameters
[in]selfgbl::GblHelixPrediction to operate on
[out]positionlength-3 double array that will hold position

Definition at line 1046 of file JavaNativeAccessWrappers.cpp.

References gbl::GblHelixPrediction::getPosition().

◆ GblHelixPredictionCtor()

GblHelixPrediction * GblHelixPredictionCtor ( double  sArc,
double  aPred[],
double  tDir[],
double  uDir[],
double  vDir[],
double  nDir[],
double  aPos[] 
)

create a new helix prediction manually

Parameters
[in]sArcarc length
[in]aPredlength-2 double array predicted measurement
[in]tDirlength-3 double array defining t direction
[in]uDirlength-3 double array defining u direction
[in]vDirlength-3 double array defining v direction
[in]nDirlength-3 double array defining n direction
[in]aPoslength-3 double array defining position

Definition at line 989 of file JavaNativeAccessWrappers.cpp.

◆ GblPoint_addGlobals()

void GblPoint_addGlobals ( GblPoint self,
int *  labels,
int  nlabels,
double *  derArray 
)

add global derivatives to the point

See also
gbl::GblPoint::addGlobals
Parameters
[in]selfgbl::GblPoint to operate on
[in]labelsint array of labels nlabels long
[in]nlabelsnumber of labels
[in]derArraydouble array of derivatives nlabels long

Definition at line 399 of file JavaNativeAccessWrappers.cpp.

References gbl::GblPoint::addGlobals().

◆ GblPoint_addMeasurement2D()

void GblPoint_addMeasurement2D ( GblPoint self,
double *  projArray,
double *  resArray,
double *  precArray,
double  minPrecision 
)

add a 2D measurement

Note
We are only supporting a 2D position residual and a 2x2 projection matrix!
See also
gbl::GblPoint::addMeasurement
Parameters
[in]selfgbl::GblPoint to operate on
[in]projArraylength-4 array holding the entries in the 2x2 proj matrix
[in]resArraylength-2 array holding residuals
[in]precArraylength-2 array holding the precisions
[in]minPrecisionMinimal precision to accept measurement

Definition at line 352 of file JavaNativeAccessWrappers.cpp.

References gbl::GblPoint::addMeasurement().

◆ GblPoint_addScatterer()

void GblPoint_addScatterer ( GblPoint self,
double *  resArray,
double *  precArray 
)

add a vector-precision scatterer

Note
only supporting vector-precision at this time
See also
gbl::GblPoint::addScatterer
Parameters
[in]selfgbl::GblPoint to operate on
[in]resArraylength-2 residuals array
[in]precArraylength-2 precision array

Definition at line 378 of file JavaNativeAccessWrappers.cpp.

References gbl::GblPoint::addScatterer().

◆ GblPoint_delete()

void GblPoint_delete ( GblPoint self)

delete gbl::GblPoint

provided so that users of JNA can clean up the memory that is not handled automatically by java itself

Parameters
[in]selfgbl::GblPoint to delete

Definition at line 300 of file JavaNativeAccessWrappers.cpp.

◆ GblPoint_getGlobalLabelsAndDerivatives()

void GblPoint_getGlobalLabelsAndDerivatives ( GblPoint self,
int *  nlabels,
int **  labels,
double **  ders 
)

get global derivatives and their labels from the point

See also
gbl::GblPoint::getGlobalLabelsAndDerivatives
Parameters
[in]selfGblPoint to operate on
[out]nlabelspointer to int where number of labels will be stored
[out]labelspointer to array where labels will be stored
[out]derspointer to array where derivatives will be stored

Definition at line 420 of file JavaNativeAccessWrappers.cpp.

References gbl::GblPoint::getGlobalLabelsAndDerivatives(), gbl::GblPoint::getMeasBegin(), and gbl::GblPoint::getMeasEnd().

◆ GblPoint_getNumMeasurements()

int GblPoint_getNumMeasurements ( GblPoint self)

calculate number of measurements in a gbl::GblPoint

Since java struggles to handle the C++ iterators, we calculate the difference between the iterators here so that the java side can access how many measurements a GblPoint has.

Parameters
[in]selfgbl::GblPoint to operate on
Returns
number of measurements

Definition at line 331 of file JavaNativeAccessWrappers.cpp.

References gbl::GblPoint::getMeasBegin(), and gbl::GblPoint::getMeasEnd().

◆ GblPoint_printPoint()

void GblPoint_printPoint ( const GblPoint self,
unsigned int  level 
)

call gbl::GblPoint::printPoint on self

Definition at line 314 of file JavaNativeAccessWrappers.cpp.

References gbl::GblPoint::printPoint().

◆ GblPointCtor()

GblPoint * GblPointCtor ( double  matrixArray[NROW *NCOL])

create new gbl::GblPoint from a jacobian

Parameters
[in]matrixArrayC-style double array listing the jacobian row-by-row.
Returns
dynamically created GblPoint

Definition at line 280 of file JavaNativeAccessWrappers.cpp.

◆ GblSimpleHelix_delete()

void GblSimpleHelix_delete ( GblSimpleHelix self)

delete gbl::GblSimpleHelix

Definition at line 897 of file JavaNativeAccessWrappers.cpp.

◆ GblSimpleHelix_getArcLengthR()

double GblSimpleHelix_getArcLengthR ( GblSimpleHelix self,
double  aRadius 
)

call gbl::GblSimpleHelix::getArcLengthR

Parameters
[in]selfgbl::GblSimpleHelix to operate on
[in]aRadiusradius to calculate arc length of
Returns
arc length calculation

Definition at line 920 of file JavaNativeAccessWrappers.cpp.

References gbl::GblSimpleHelix::getArcLengthR().

◆ GblSimpleHelix_getArcLengthXY()

double GblSimpleHelix_getArcLengthXY ( GblSimpleHelix self,
double  xPos,
double  yPos 
)

call gbl::GblSimpleHelix::getArcLengthXY

Parameters
[in]selfgbl::GblSimpleHelix to operate on
[in]xPosx-position
[in]yPosy-position
Returns
arc length calculation

Definition at line 931 of file JavaNativeAccessWrappers.cpp.

References gbl::GblSimpleHelix::getArcLengthXY().

◆ GblSimpleHelix_getPhi()

double GblSimpleHelix_getPhi ( GblSimpleHelix self,
double  aRadius 
)

call gbl::GblSimpleHelix::getPhi

Parameters
[in]selfgbl::GblSimpleHelix to operate on
[in]aRadiusradius to calculate phi at
Returns
phi calculation

Definition at line 910 of file JavaNativeAccessWrappers.cpp.

References gbl::GblSimpleHelix::getPhi().

◆ GblSimpleHelix_getPrediction()

GblHelixPrediction * GblSimpleHelix_getPrediction ( GblSimpleHelix self,
double  refPos[],
double  uDir[],
double  vDir[] 
)

get a helical prediction from a reference coordinate system

See also
gbl::GblSimpleHelix::getPrediction
Parameters
[in]selfgbl::GblSimpleHelix to operate on
[in]refPosdouble array of length three containing reference position
[in]uDir3-length double array defining u direction
[in]vDir3-length double array defining v direction
Returns
new gbl::GblHelixPrediction from this coordinate system

Definition at line 959 of file JavaNativeAccessWrappers.cpp.

References gbl::GblSimpleHelix::getPrediction().

◆ GblSimpleHelix_moveToXY()

void GblSimpleHelix_moveToXY ( GblSimpleHelix self,
double  xPos,
double  yPos,
double *  newPhi0,
double *  newDca,
double *  newZ0 
)

call gbl::GblSimpleHelix::moveToXY

Parameters
[in]selfgbl::GblSimpleHelix to operate on
[in]xPosx-position
[in]yPosy-position
[out]newPhi0double address to store resulting phi
[out]newDcadouble address to store resulting Dca
[out]newZ0double address to store resulting Z0

Definition at line 944 of file JavaNativeAccessWrappers.cpp.

References gbl::GblSimpleHelix::moveToXY().

◆ GblSimpleHelixCtor()

GblSimpleHelix * GblSimpleHelixCtor ( double  aRinv,
double  aPhi0,
double  aDca,
double  aDzds,
double  aZ0 
)

construct a new simple helix

The input parameters are all the same as the C++ constructor gbl::GblSimpleHelix::GblSimpleHelix.

Returns
newly constructed gbl::GblSimpleHelix

Definition at line 887 of file JavaNativeAccessWrappers.cpp.

◆ GblTrajectory_delete()

void GblTrajectory_delete ( GblTrajectory self)

delete self

Cleanup for JNA, handles cleaning up all GblPoints within it!

Definition at line 611 of file JavaNativeAccessWrappers.cpp.

References gbl::GblTrajectory::getNumPoints().

◆ GblTrajectory_fit()

void GblTrajectory_fit ( GblTrajectory self,
double *  Chi2,
int *  Ndf,
double *  lostWeight,
char *  c_optionList,
unsigned int  aLabel 
)

call gbl::GblTrajectory::fit on self

Parameters
[in]selfgbl::GblTrajectory to operate on
[out]Chi2pointer to double where Chi2 result will be stored
[out]Ndfpointer to double where Ndf result will be stored
[out]lostWeightpointer to double where lostWeight result will be stored
[in]c_optionListC-style string listing options
[in]aLabelinteger label for fit

Definition at line 595 of file JavaNativeAccessWrappers.cpp.

References gbl::GblTrajectory::fit().

◆ GblTrajectory_getMeasResults()

int GblTrajectory_getMeasResults ( GblTrajectory self,
int  aLabel,
int *  numData,
double *  aResiduals,
double *  aMeasErrors,
double *  aResErrors,
double *  aDownWeights 
)

get measurement results from trajectory

See also
gbl::GblTrajectory::getMeasResults
Parameters
[in]selfgbl::GblTrajectory to operate on
[in]aLabelinteger label of results
[out]numDatalength of resulting arrays (expect 2)
[out]aResidualsdouble array of residuals
[out]aMeasErrorsdouble array of measurement errors
[out]aResErrorsdouble array of residual errors
[out]aDownWeightsdouble array of down weights

Definition at line 727 of file JavaNativeAccessWrappers.cpp.

References gbl::GblTrajectory::getMeasResults().

◆ GblTrajectory_getNumPoints()

int GblTrajectory_getNumPoints ( GblTrajectory self)

call gbl::GblTrajectory::getNumPoints

Parameters
[in]selfgbl::GblTrajectory to operate on
Returns
number of points

Definition at line 641 of file JavaNativeAccessWrappers.cpp.

References gbl::GblTrajectory::getNumPoints().

◆ GblTrajectory_getResults()

void GblTrajectory_getResults ( GblTrajectory self,
int  aSignedLabel,
double *  localPar,
int *  nLocalPar,
double *  localCov,
int *  sizeLocalCov 
)

get trajectory results

Note
only supports 5-vector and 5x5 covariance matrix!
Parameters
[in]selfgbl::GblTrajectory to operate on
[in]aSignedLabelinteger label for results
[out]localParlength 5 double array where local results will be stored
[out]nLocalParinteger length of double array (always set to 5)
[out]localCovlength 25 double array where local covariance results will be stored
[out]sizeLocalCovinteger dimension of covariance matrix (always set to 5)

Definition at line 693 of file JavaNativeAccessWrappers.cpp.

References gbl::GblTrajectory::getResults().

◆ GblTrajectory_isValid()

int GblTrajectory_isValid ( GblTrajectory self)

call gbl::GblTrajectory::isValid

Parameters
[in]selfgbl::GblTrajectory to operate on
Returns
0 if true, 1 if false

Definition at line 629 of file JavaNativeAccessWrappers.cpp.

References gbl::GblTrajectory::isValid().

◆ GblTrajectory_milleOut()

void GblTrajectory_milleOut ( GblTrajectory self,
MilleBinary millebinary 
)

write a trajectory to a gbl::MilleBinary file

See also
gbl::GblTrajectory::milleOut
Parameters
[in]selfgbl::GblTrajectory to write out
[in]millebinarygbl::MilleBinary to write to

Definition at line 766 of file JavaNativeAccessWrappers.cpp.

References gbl::GblTrajectory::milleOut().

◆ GblTrajectory_printData()

void GblTrajectory_printData ( GblTrajectory self)

call gbl::GblTrajectory::printData

Parameters
[in]selfgbl::GblTrajectory to operate on

Definition at line 664 of file JavaNativeAccessWrappers.cpp.

References gbl::GblTrajectory::printData().

◆ GblTrajectory_printPoints()

void GblTrajectory_printPoints ( GblTrajectory self,
int  level 
)

call gbl::GblTrajectory::printPoints

Parameters
[in]selfgbl::GblTrajectory to operate on
[in]levelinteger level for printing

Definition at line 676 of file JavaNativeAccessWrappers.cpp.

References gbl::GblTrajectory::printPoints().

◆ GblTrajectory_printTrajectory()

void GblTrajectory_printTrajectory ( GblTrajectory self,
int  level 
)

call gbl::GblTrajectory::printTrajectory

Parameters
[in]selfgbl::GblTrajectory to operate on
[in]levelinteger level for printing

Definition at line 653 of file JavaNativeAccessWrappers.cpp.

References gbl::GblTrajectory::printTrajectory().

◆ GblTrajectoryCtorPtrArray()

GblTrajectory * GblTrajectoryCtorPtrArray ( GblPoint points[],
int  npoints,
int  flagCurv,
int  flagU1dir,
int  flagU2dir 
)

simple construction of a new gbl::GblTrajectory

See also
gbl::GblTrajectory::GblTrajectory

Just constructing a trajectory from a set of points. No seeding matrix.

Parameters
[in]pointsarray of pointers to gbl::GblPoint to put into trajectory
[in]npointsnumber of points in array
[in]flagCurvuse q/p - non-zero for true, zero for false
[in]flagU1diruse in u1 direction - non-zero for true, zero for false
[in]flagU2diruse in u2 direction - non-zero for true, zero for false
Returns
dynamically allocated trajectory wrapping input points

Definition at line 486 of file JavaNativeAccessWrappers.cpp.

References ptr_array_to_vector().

◆ GblTrajectoryCtorPtrArraySeed()

GblTrajectory * GblTrajectoryCtorPtrArraySeed ( GblPoint points[],
int  npoints,
int  aLabel,
double  seedArray[],
int  flagCurv,
int  flagU1dir,
int  flagU2dir 
)

construct new gbl::GblTrajectory with a seed matrix

See also
gbl::GblTrajectory::GblTrajectory
Parameters
[in]pointsarray of pointers to gbl::GblPoint to put into trajectory
[in]npointsnumber of points in array
[in]aLabelinteger label for seed
[in]seedArraydouble-array of length 25 listing seed matrix elements row-wise
[in]flagCurvuse q/p - non-zero for true, zero for false
[in]flagU1diruse in u1 direction - non-zero for true, zero for false
[in]flagU2diruse in u2 direction - non-zero for true, zero for false
Returns
dynamically allocated trajectory wrapping input points

Definition at line 514 of file JavaNativeAccessWrappers.cpp.

References ptr_array_to_vector().

◆ GblTrajectoryCtorPtrComposed()

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

See also
gbl::GblTrajectory::GblTrajectory
Parameters
[in]points_1array of gbl::GblPoint for one track
[in]npoints_1number of points
[in]trafo_1double array listing elements of 2x3 track matrix
[in]points_2array of gbl::GblPoint for one track
[in]npoints_2number of points
[in]trafo_2double array listing elements of 2x3 track matrix
Returns
dynamically created trajectory composed of two tracks

Definition at line 545 of file JavaNativeAccessWrappers.cpp.

References ptr_array_to_vector().

◆ MilleBinary_close()

void MilleBinary_close ( MilleBinary self)

Closing a gbl::MilleBinary file is the same as destructing it.

Since the destructor of the gbl::MilleBinary class is what handles performing the final write operations, we simply delete the object pointed to by the passed pointer.

Note
This means the object on the JNA side will be invalid and will cause a program crash if it is accessed after using this function on it!
Parameters
[in]selfgbl::MilleBinary to delete

Definition at line 264 of file JavaNativeAccessWrappers.cpp.

◆ MilleBinaryCtor()

MilleBinary * MilleBinaryCtor ( const char *  fileName,
int  filenamesize,
int  doublePrecision,
int  keepZeros,
int  aSize 
)

Dynamically create new MilleBinary file.

We return the raw pointer to this object so that JNA can effectively access it.

Unfortunately the translation of booleans to/from java is not very stable, so it is safer to simply pass integers and check if they are 0 (for false) or nonzero (for true).

Parameters
[in]fileNamename of file in a C-style string
[in]filenamesizelength of filename for accessing C-style string
[in]doublePrecisionnon-zero if should use doublePrecision, zero if should use float precision
[in]keepZerosnon-zero to keep zeros, zero to not keep them
[in]aSizesize of buffer to keep in memory

Definition at line 235 of file JavaNativeAccessWrappers.cpp.

◆ ptr_array_to_vector()

std::vector< GblPoint > ptr_array_to_vector ( GblPoint points[],
int  npoints 
)

convert the pointer array of gbl::GblPoint into a vector holding the objects

This is a helper function and should not be bound to a function by JNA, so we are keeping it outside the extern "C" block so its name is mangled.

Note
We copy the data pointed to into the vector, so the points input into this function still need to be deleted.
Parameters
[in]pointsarray of pointers to gblGblPoint to put into vector
[in]npointsnumber of points (size of array)
Returns
vector of GblPoints with same content as array

Definition at line 194 of file JavaNativeAccessWrappers.cpp.

Referenced by GblTrajectoryCtorPtrArray(), GblTrajectoryCtorPtrArraySeed(), and GblTrajectoryCtorPtrComposed().

Variable Documentation

◆ NCOL

const int NCOL = 5

Definition at line 129 of file JavaNativeAccessWrappers.cpp.

◆ NROW

const int NROW = 5

Definition at line 128 of file JavaNativeAccessWrappers.cpp.