97   std::lock_guard<std::mutex> guard(
fMutex);
 
   98   std::vector<ThreadData_t *>::iterator i = 
fThreadData.begin();
 
 
  114   std::lock_guard<std::mutex> guard(
fMutex);
 
 
  220         Fatal(
"ComputeBBox", 
"Wrong section order");
 
  227      Fatal(
"ComputeBBox", 
"Shape %s at index %d: Not allowed first two or last two sections at same Z", 
GetName(),
 
  280   fOrigin[2] = 0.5 * (zmax + zmin);
 
  283   fDZ = 0.5 * (zmax - zmin);
 
 
  320   if ((
fZ[
ipl + 1] - point[2]) < (point[2] - 
fZ[
ipl]))
 
  377   if (
norm[0] * dir[0] + 
norm[1] * dir[1] + 
norm[2] * dir[2] < 0) {
 
 
  390   if (point[2] < 
fZ[0])
 
  392   if (point[2] > 
fZ[
fNz - 1])
 
 
  471      ((
TGeoPgon *)
this)->CreateThreadData(1);
 
  481      if ((point[0] * dir[1] - point[1] * dir[0]) > 0) {
 
 
  591      sphi[0] = 
TMath::Sqrt((point[0] * point[0] + point[1] * point[1]) / (1. - dir[2] * dir[2]));
 
 
  725      for (i = 0; i < 3; i++)
 
  726         pt[i] = point[i] + step * dir[i];
 
 
  783         for (i = 0; i < 3; i++)
 
  814      for (i = 0; i < 3; i++)
 
  815         pt[i] = point[i] + step * dir[i];
 
 
  939            for (i = 0; i < 3; i++)
 
 
 1030         for (i = 0; i < 3; i++)
 
 
 1082      for (
Int_t i = 0; i < 3; i++)
 
 
 1230      for (i = 0; i < 3; i++)
 
 1258      ((
TGeoPgon *)
this)->CreateThreadData(1);
 
 1283   if (
fDphi < 360.0) {
 
 1354                     if (
pt[2] * dir[2] < 0)
 
 1399      if (
iph[0] >= 0 && 
sph[0] > 1.E-8)
 
 
 1440   Double_t zmax = start + ndiv * step;
 
 1445      Error(
"Divide", 
"makes no sense dividing a pgon on radius");
 
 1449         Error(
"Divide", 
"ndiv should divide number of pgon edges");
 
 1463      for (
id = 0; 
id < ndiv; 
id++) {
 
 1464         voldiv->AddNodeOffset(vol, 
id, start + 
id * step + step / 2, opt.
Data());
 
 1471         if (start < 
fZ[
ipl])
 
 1474            if ((start + ndiv * step) > 
fZ[
ipl + 1])
 
 1483         Error(
"Divide", 
"cannot divide pcon on Z if divided region is not between 2 consecutive planes");
 
 1491      for (
id = 0; 
id < ndiv; 
id++) {
 
 1503         voldiv->AddNodeOffset(vol, 
id, start + 
id * step + step / 2, opt.
Data());
 
 1507   default: 
Error(
"Divide", 
"Wrong axis type for division"); 
return nullptr;
 
 
 1517   param[0] = 
fRmin[0]; 
 
 1518   param[1] = 
fRmax[0]; 
 
 1520      if (
fRmin[i] < param[0])
 
 1521         param[0] = 
fRmin[i];
 
 1522      if (
fRmax[i] > param[1])
 
 1523         param[1] = 
fRmax[i];
 
 1527   param[0] *= param[0];
 
 1528   param[1] *= param[1];
 
 1535   param[3] = param[2] + 
fDphi;                     
 
 
 1597   for (i = 0; i < 
nz * 2; i++) {
 
 1599      for (
j = 1; 
j < 
n; 
j++) {
 
 1612   for (i = 0; i < 2; i++) {
 
 1614      for (
j = 0; 
j < 
n; 
j++) {
 
 1622   for (i = 0; i < (
nz - 1); i++) {
 
 1625      for (
j = 0; 
j < 
n; 
j++) {
 
 1632      for (
j = 0; 
j < 
n; 
j++) {
 
 1642      for (i = 1; i < (
nz - 1); i++) {
 
 1643         for (
j = 0; 
j < 2; 
j++) {
 
 1646            buff.fSegs[
indx++] = (2 * i + 1) * 
n + 
j * (
n - 1);
 
 1657   for (
j = 0; 
j < 
n - 1; 
j++) {
 
 1674   for (
j = 0; 
j < 
n - 1; 
j++) {
 
 1692   for (k = 0; k < (
nz - 1); k++) {
 
 1694      for (
j = 0; 
j < 
n - 1; 
j++) {
 
 1697         buff.fPols[
indx++] = 
nz * 2 * 
m + (2 * k + i * 1 + 2) * 
n + 
j + 1;
 
 1698         buff.fPols[
indx++] = (2 * k + i * 1 + 2) * 
m + 
j;
 
 1699         buff.fPols[
indx++] = 
nz * 2 * 
m + (2 * k + i * 1 + 2) * 
n + 
j;
 
 1700         buff.fPols[
indx++] = (2 * k + i * 1) * 
m + 
j;
 
 1705         buff.fPols[
indx++] = 
nz * 2 * 
m + (2 * k + i * 1 + 2) * 
n;
 
 1706         buff.fPols[
indx++] = (2 * k + i * 1 + 2) * 
m + 
j;
 
 1707         buff.fPols[
indx++] = 
nz * 2 * 
m + (2 * k + i * 1 + 2) * 
n + 
j;
 
 1708         buff.fPols[
indx++] = (2 * k + i * 1) * 
m + 
j;
 
 1711      for (
j = 0; 
j < 
n - 1; 
j++) {
 
 1714         buff.fPols[
indx++] = (2 * k + i * 1) * 
m + 
j;
 
 1715         buff.fPols[
indx++] = 
nz * 2 * 
m + (2 * k + i * 1 + 2) * 
n + 
j;
 
 1716         buff.fPols[
indx++] = (2 * k + i * 1 + 2) * 
m + 
j;
 
 1717         buff.fPols[
indx++] = 
nz * 2 * 
m + (2 * k + i * 1 + 2) * 
n + 
j + 1;
 
 1722         buff.fPols[
indx++] = (2 * k + i * 1) * 
m + 
j;
 
 1723         buff.fPols[
indx++] = 
nz * 2 * 
m + (2 * k + i * 1 + 2) * 
n + 
j;
 
 1724         buff.fPols[
indx++] = (2 * k + i * 1 + 2) * 
m + 
j;
 
 1725         buff.fPols[
indx++] = 
nz * 2 * 
m + (2 * k + i * 1 + 2) * 
n;
 
 1733      for (k = 0; k < (
nz - 1); k++) {
 
 
 1762   if ((
nz < 2) || (
nbPnts <= 0) || (
n < 2))
 
 1770   for (i = 0; i < 
nz; i++) {
 
 1772      for (
j = 1; 
j < 
n; 
j++) {
 
 1781   for (
j = 0; 
j < 
n; 
j++) {
 
 1789   for (
j = 0; 
j < 
n; 
j++) {
 
 1796   for (i = 0; i < (
nz - 1); i++) {
 
 1799      for (
j = 0; 
j < 
n; 
j++) {
 
 1811   for (
j = 0; 
j < 
n - 1; 
j++) {
 
 1822   for (
j = 0; 
j < 
n - 1; 
j++) {
 
 1831   for (
Int_t k = 0; k < (
nz - 1); k++) {
 
 1834      for (
j = 0; 
j < 
n - 1; 
j++) {
 
 
 1856      Fatal(
"Rpg", 
"Plane index parameter ipl=%i out of range\n", 
ipl);
 
 
 1894   a = ((point[0] * dir[2] - point[2] * dir[0]) * 
cphi + (point[1] * dir[2] - point[2] * dir[1]) * 
sphi) * 
invdirz;
 
 
 1928      r = 
TMath::Sqrt(point[0] * point[0] + point[1] * point[1]);
 
 1939      for (i = 0; i < 3; i++)
 
 1965      for (i = 0; i < 3; i++)
 
 
 2069   out << 
"   // Shape: " << 
GetName() << 
" type: " << 
ClassName() << std::endl;
 
 2070   out << 
"   phi1    = " << 
fPhi1 << 
";" << std::endl;
 
 2071   out << 
"   dphi    = " << 
fDphi << 
";" << std::endl;
 
 2072   out << 
"   nedges = " << 
fNedges << 
";" << std::endl;
 
 2073   out << 
"   nz      = " << 
fNz << 
";" << std::endl;
 
 2074   out << 
"   auto " << 
GetPointerName() << 
" = new TGeoPgon(\"" << 
GetName() << 
"\", phi1, dphi, nedges, nz);" 
 2077      out << 
"      z     = " << 
fZ[i] << 
";" << std::endl;
 
 2078      out << 
"      rmin  = " << 
fRmin[i] << 
";" << std::endl;
 
 2079      out << 
"      rmax  = " << 
fRmax[i] << 
";" << std::endl;
 
 2080      out << 
"   " << 
GetPointerName() << 
"->DefineSection(" << i << 
", z, rmin, rmax);" << std::endl;
 
 
 2095      Error(
"SetDimensions", 
"Pgon %s: Number of Z sections must be > 2", 
GetName());
 
 2111      DefineSection(i, param[4 + 3 * i], param[5 + 3 * i], param[6 + 3 * i]);
 
 
 2129      for (i = 0; i < 
GetNz(); i++) {
 
 2131            for (
j = 0; 
j < 
n; 
j++) {
 
 2137         for (
j = 0; 
j < 
n; 
j++) {
 
 
 2172      for (i = 0; i < 
fNz; i++) {
 
 2174            for (
j = 0; 
j < 
n; 
j++) {
 
 2180         for (
j = 0; 
j < 
n; 
j++) {
 
 
 2221      npols = 2 * (
n - 1) + (
nz - 1) * (
n - 1);
 
 
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
void Fatal(const char *location, const char *msgfmt,...)
Use this function in case of a fatal error. It will abort the program.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize id
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
R__EXTERN TGeoManager * gGeoManager
Generic 3D primitive description class.
Bool_t SectionsValid(UInt_t mask) const
void SetSectionsValid(UInt_t mask)
Bool_t SetRawSizes(UInt_t reqPnts, UInt_t reqPntsCapacity, UInt_t reqSegs, UInt_t reqSegsCapacity, UInt_t reqPols, UInt_t reqPolsCapacity)
Set kRaw tessellation section of buffer with supplied sizes.
void FillBuffer3D(TBuffer3D &buffer, Int_t reqSections, Bool_t localFrame) const override
Fills the supplied buffer, with sections in desired frame See TBuffer3D.h for explanation of sections...
Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const override
Compute distance from outside point to surface of the box.
TGeoVolumeMulti * MakeVolumeMulti(const char *name, TGeoMedium *medium)
Make a TGeoVolumeMulti handling a list of volumes.
TObjArray * GetListOfShapes() const
static Int_t ThreadId()
Translates the current thread id to an ordinal number.
Node containing an offset.
a cylindrical phi divison pattern
base finder class for patterns. A pattern is specifying a division type
A polycone is represented by a sequence of tubes/cones, glued together at defined Z planes.
virtual void DefineSection(Int_t snum, Double_t z, Double_t rmin, Double_t rmax)
Defines z position of a section plane, rmin and rmax at this z.
void InspectShape() const override
print shape parameters
Bool_t HasInsideSurface() const
Returns true when pgon has internal surface It will be only disabled when all Rmin values are 0.
Polygons are defined in the same way as polycones, the difference being just that the segments betwee...
void Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const override
Compute safe distance from each of the points in the input array.
TBuffer3D * MakeBuffer3D() const override
Creates a TBuffer3D describing this shape.
void SetPoints(Double_t *points) const override
create polygone mesh points
Bool_t Contains(const Double_t *point) const override
test if point is inside this shape check total z range
Bool_t SliceCrossingInZ(const Double_t *point, const Double_t *dir, Int_t nphi, Int_t *iphi, Double_t *sphi, Double_t &snext, Double_t stepmax) const
Performs ray propagation between Z segments.
~TGeoPgon() override
destructor
Bool_t SliceCrossing(const Double_t *point, const Double_t *dir, Int_t nphi, Int_t *iphi, Double_t *sphi, Double_t &snext, Double_t stepmax) const
Check boundary crossing inside phi slices.
void CreateThreadData(Int_t nthreads) override
Create thread data for n threads max.
void Sizeof3D() const override
fill size of this 3-D object
Int_t GetNmeshVertices() const override
Return number of vertices of the mesh representation.
void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm) const override
Compute normal to closest surface from POINT.
void GetBoundingCylinder(Double_t *param) const override
Fill vector param[4] with the bounding cylinder parameters.
void DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const override
Compute distance from array of input points having directions specified by dirs. Store output in dist...
Bool_t SliceCrossingZ(const Double_t *point, const Double_t *dir, Int_t nphi, Int_t *iphi, Double_t *sphi, Double_t &snext, Double_t stepmax) const
Performs ray propagation between Z segments.
void InspectShape() const override
Inspect the PGON parameters.
void LocatePhi(const Double_t *point, Int_t &ipsec) const
Locates index IPSEC of the phi sector containing POINT.
std::mutex fMutex
Size for the navigation data array.
TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step) override
Divide this polygone shape belonging to volume "voldiv" into ndiv volumes called divname,...
std::vector< ThreadData_t * > fThreadData
void GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const override
Returns numbers of vertices, segments and polygons composing the shape mesh.
ThreadData_t & GetThreadData() const
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save a primitive as a C++ statement(s) on output stream "out".
Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const override
computes the closest distance from given point to this shape, according to option.
void ClearThreadData() const override
Double_t SafetyToSegment(const Double_t *point, Int_t ipl, Int_t iphi, Bool_t in, Double_t safphi, Double_t safmin=TGeoShape::Big()) const
Compute safety from POINT to segment between planes ipl, ipl+1 within safmin.
Double_t Capacity() const override
Computes capacity of the shape in [length^3].
Double_t Rpg(Double_t z, Int_t ipl, Bool_t inner, Double_t &a, Double_t &b) const
Computes projected pgon radius (inner or outer) corresponding to a given Z value.
void SetDimensions(Double_t *param) override
Set PGON dimensions starting from an array.
void SetSegsAndPolsNoInside(TBuffer3D &buff) const
Fill TBuffer3D structure for segments and polygons, when no inner surface exists.
Bool_t SliceCrossingIn(const Double_t *point, const Double_t *dir, Int_t ipl, Int_t nphi, Int_t *iphi, Double_t *sphi, Double_t &snext, Double_t stepmax) const
Check boundary crossing inside phi slices.
void Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const override
Check the inside status for each of the points in the array.
void DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const override
Compute distance from array of input points having directions specified by dirs. Store output in dist...
const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const override
Fills a static 3D buffer and returns a reference.
Int_t fThreadSize
Navigation data per thread.
void ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize) override
Compute the normal for an array o points so that norm.dot.dir is positive Input: Arrays of point coor...
void ComputeBBox() override
compute bounding box for a polygone Check if the sections are in increasing Z order
Double_t Rproj(Double_t z, const Double_t *point, const Double_t *dir, Double_t cphi, Double_t sphi, Double_t &a, Double_t &b) const
Computes projected distance at a given Z for a given ray inside a given sector and fills coefficients...
Double_t DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const override
compute distance from inside point to surface of the polygone first find out in which Z section the p...
Bool_t IsCrossingSlice(const Double_t *point, const Double_t *dir, Int_t iphi, Double_t sstart, Int_t &ipl, Double_t &snext, Double_t stepmax) const
Check crossing of a given pgon slice, from a starting point inside the slice.
Int_t GetPhiCrossList(const Double_t *point, const Double_t *dir, Int_t istart, Double_t *sphi, Int_t *iphi, Double_t stepmax=TGeoShape::Big()) const
Mutex for thread data.
void SetSegsAndPols(TBuffer3D &buff) const override
Fill TBuffer3D structure for segments and polygons.
Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const override
Compute distance from outside point to surface of the polygone.
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
compute closest distance from point px,py to each corner
Base abstract class for all shapes.
Int_t GetBasicColor() const
Get the basic color (0-7).
void TransformPoints(Double_t *points, UInt_t NbPoints) const
Tranform a set of points (LocalToMaster)
void SetShapeBit(UInt_t f, Bool_t set)
Equivalent of TObject::SetBit.
static Double_t SafetyPhi(const Double_t *point, Bool_t in, Double_t phi1, Double_t phi2)
Static method to compute safety w.r.t a phi corner defined by cosines/sines of the angles phi1,...
static Bool_t IsSameWithinTolerance(Double_t a, Double_t b)
Check if two numbers differ with less than a tolerance.
const char * GetPointerName() const
Provide a pointer name containing uid.
Int_t ShapeDistancetoPrimitive(Int_t numpoints, Int_t px, Int_t py) const
Returns distance to shape primitive mesh.
static void NormalPhi(const Double_t *point, const Double_t *dir, Double_t *norm, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
Static method to compute normal to phi planes.
static Bool_t IsCrossingSemiplane(const Double_t *point, const Double_t *dir, Double_t cphi, Double_t sphi, Double_t &snext, Double_t &rxy)
Compute distance from POINT to semiplane defined by PHI angle along DIR.
const char * GetName() const override
Get the shape name.
static Double_t Tolerance()
static Bool_t IsCloseToPhi(Double_t epsil, const Double_t *point, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
True if point is closer than epsil to one of the phi planes defined by c1,s1 or c2,...
static Double_t DistFromOutsideS(const Double_t *point, const Double_t *dir, Double_t rmin, Double_t rmax, Double_t dz)
Static method to compute distance from outside point to a tube with given parameters Boundary safe al...
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Int_t IndexOf(const TObject *obj) const override
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
virtual const char * ClassName() const
Returns name of class to which the object belongs.
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
const char * Data() const
Long64_t LocMin(Long64_t n, const T *a)
Returns index of array with the minimum element.
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
T1 Sign(T1 a, T2 b)
Returns a value with the magnitude of a and the sign of b.
Double_t ATan2(Double_t y, Double_t x)
Returns the principal value of the arc tangent of y/x, expressed in radians.
Long64_t LocMax(Long64_t n, const T *a)
Returns index of array with the maximum element.
constexpr Double_t DegToRad()
Conversion from degree to radian: .
Double_t Sqrt(Double_t x)
Returns the square root of x.
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Double_t Tan(Double_t)
Returns the tangent of an angle of x radians.
Long64_t BinarySearch(Long64_t n, const T *array, T value)
Binary search in an array of n values to locate value.
constexpr Double_t RadToDeg()
Conversion from radian to degree: .
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
ThreadData_t()
[fNedges+4] temporary double buffer array
~ThreadData_t()
Destructor.