Geant4  10.02.p03
G4GenericPolycone.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4GenericPolycone.cc 72131 2013-07-11 07:12:24Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class source file
32 //
33 //
34 // G4GenericPolycone.cc
35 //
36 // Implementation of a CSG polycone
37 //
38 // --------------------------------------------------------------------
39 
40 #include "G4GenericPolycone.hh"
41 
42 #if !defined(G4GEOM_USE_UGENERICPOLYCONE)
43 
44 #include "G4PolyconeSide.hh"
45 #include "G4PolyPhiFace.hh"
46 
47 #include "Randomize.hh"
48 
49 #include "G4Polyhedron.hh"
50 #include "G4EnclosingCylinder.hh"
51 #include "G4ReduciblePolygon.hh"
52 #include "G4VPVParameterisation.hh"
53 
54 using namespace CLHEP;
55 
56 //
57 // Constructor (generic parameters)
58 //
60  G4double phiStart,
61  G4double phiTotal,
62  G4int numRZ,
63  const G4double r[],
64  const G4double z[] )
65  : G4VCSGfaceted( name )
66 {
67 
68  G4ReduciblePolygon *rz = new G4ReduciblePolygon( r, z, numRZ );
69 
70  Create( phiStart, phiTotal, rz );
71 
72  // Set original_parameters struct for consistency
73  //
74  //SetOriginalParameters(rz);
75 
76  delete rz;
77 }
78 
79 //
80 // Create
81 //
82 // Generic create routine, called by each constructor after
83 // conversion of arguments
84 //
86  G4double phiTotal,
87  G4ReduciblePolygon *rz )
88 {
89  //
90  // Perform checks of rz values
91  //
92  if (rz->Amin() < 0.0)
93  {
94  std::ostringstream message;
95  message << "Illegal input parameters - " << GetName() << G4endl
96  << " All R values must be >= 0 !";
97  G4Exception("G4GenericPolycone::Create()", "GeomSolids0002",
98  FatalErrorInArgument, message);
99  }
100 
101  G4double rzArea = rz->Area();
102  if (rzArea < -kCarTolerance)
103  {
104  rz->ReverseOrder();
105  }
106  else if (rzArea < kCarTolerance)
107  {
108  std::ostringstream message;
109  message << "Illegal input parameters - " << GetName() << G4endl
110  << " R/Z cross section is zero or near zero: " << rzArea;
111  G4Exception("G4GenericPolycone::Create()", "GeomSolids0002",
112  FatalErrorInArgument, message);
113  }
114 
116  || (!rz->RemoveRedundantVertices( kCarTolerance )) )
117  {
118  std::ostringstream message;
119  message << "Illegal input parameters - " << GetName() << G4endl
120  << " Too few unique R/Z values !";
121  G4Exception("G4GenericPolycone::Create()", "GeomSolids0002",
122  FatalErrorInArgument, message);
123  }
124 
125  if (rz->CrossesItself(1/kInfinity))
126  {
127  std::ostringstream message;
128  message << "Illegal input parameters - " << GetName() << G4endl
129  << " R/Z segments cross !";
130  G4Exception("G4GenericPolycone::Create()", "GeomSolids0002",
131  FatalErrorInArgument, message);
132  }
133 
134  numCorner = rz->NumVertices();
135 
136  //
137  // Phi opening? Account for some possible roundoff, and interpret
138  // nonsense value as representing no phi opening
139  //
140  if (phiTotal <= 0 || phiTotal > twopi-1E-10)
141  {
142  phiIsOpen = false;
143  startPhi = 0;
144  endPhi = twopi;
145  }
146  else
147  {
148  phiIsOpen = true;
149 
150  //
151  // Convert phi into our convention
152  //
153  startPhi = phiStart;
154  while( startPhi < 0 ) // Loop checking, 13.08.2015, G.Cosmo
155  startPhi += twopi;
156 
157  endPhi = phiStart+phiTotal;
158  while( endPhi < startPhi ) // Loop checking, 13.08.2015, G.Cosmo
159  endPhi += twopi;
160  }
161 
162  //
163  // Allocate corner array.
164  //
166 
167  //
168  // Copy corners
169  //
170  G4ReduciblePolygonIterator iterRZ(rz);
171 
172  G4PolyconeSideRZ *next = corners;
173  iterRZ.Begin();
174  do // Loop checking, 13.08.2015, G.Cosmo
175  {
176  next->r = iterRZ.GetA();
177  next->z = iterRZ.GetB();
178  } while( ++next, iterRZ.Next() );
179 
180  //
181  // Allocate face pointer array
182  //
184  faces = new G4VCSGface*[numFace];
185 
186  //
187  // Construct conical faces
188  //
189  // But! Don't construct a face if both points are at zero radius!
190  //
191  G4PolyconeSideRZ *corner = corners,
192  *prev = corners + numCorner-1,
193  *nextNext;
194  G4VCSGface **face = faces;
195  do // Loop checking, 13.08.2015, G.Cosmo
196  {
197  next = corner+1;
198  if (next >= corners+numCorner) next = corners;
199  nextNext = next+1;
200  if (nextNext >= corners+numCorner) nextNext = corners;
201 
202  if (corner->r < 1/kInfinity && next->r < 1/kInfinity) continue;
203 
204  //
205  // We must decide here if we can dare declare one of our faces
206  // as having a "valid" normal (i.e. allBehind = true). This
207  // is never possible if the face faces "inward" in r.
208  //
209  G4bool allBehind;
210  if (corner->z > next->z)
211  {
212  allBehind = false;
213  }
214  else
215  {
216  //
217  // Otherwise, it is only true if the line passing
218  // through the two points of the segment do not
219  // split the r/z cross section
220  //
221  allBehind = !rz->BisectedBy( corner->r, corner->z,
222  next->r, next->z, kCarTolerance );
223  }
224 
225  *face++ = new G4PolyconeSide( prev, corner, next, nextNext,
226  startPhi, endPhi-startPhi, phiIsOpen, allBehind );
227  } while( prev=corner, corner=next, corner > corners );
228 
229  if (phiIsOpen)
230  {
231  //
232  // Construct phi open edges
233  //
234  *face++ = new G4PolyPhiFace( rz, startPhi, 0, endPhi );
235  *face++ = new G4PolyPhiFace( rz, endPhi, 0, startPhi );
236  }
237 
238  //
239  // We might have dropped a face or two: recalculate numFace
240  //
241  numFace = face-faces;
242 
243  //
244  // Make enclosingCylinder
245  //
247  new G4EnclosingCylinder( rz, phiIsOpen, phiStart, phiTotal );
248 }
249 
250 
251 //
252 // Fake default constructor - sets only member data and allocates memory
253 // for usage restricted to object persistency.
254 //
256  : G4VCSGfaceted(a), startPhi(0.), endPhi(0.), phiIsOpen(false),
258 {
259 }
260 
261 
262 //
263 // Destructor
264 //
266 {
267  delete [] corners;
268  delete enclosingCylinder;
269 }
270 
271 
272 //
273 // Copy constructor
274 //
276  : G4VCSGfaceted( source )
277 {
278  CopyStuff( source );
279 }
280 
281 
282 //
283 // Assignment operator
284 //
286 {
287  if (this == &source) return *this;
288 
289  G4VCSGfaceted::operator=( source );
290 
291  delete [] corners;
292  // if (original_parameters) delete original_parameters;
293 
294  delete enclosingCylinder;
295 
296  CopyStuff( source );
297 
298  return *this;
299 }
300 
301 
302 //
303 // CopyStuff
304 //
306 {
307  //
308  // Simple stuff
309  //
310  startPhi = source.startPhi;
311  endPhi = source.endPhi;
312  phiIsOpen = source.phiIsOpen;
313  numCorner = source.numCorner;
314 
315  //
316  // The corner array
317  //
319 
320  G4PolyconeSideRZ *corn = corners,
321  *sourceCorn = source.corners;
322  do // Loop checking, 13.08.2015, G.Cosmo
323  {
324  *corn = *sourceCorn;
325  } while( ++sourceCorn, ++corn < corners+numCorner );
326 
327  //
328  // Enclosing cylinder
329  //
331 
332  fRebuildPolyhedron = false;
333  fpPolyhedron = 0;
334 }
335 
336 
337 //
338 // Reset
339 //
341 {
342 
343  std::ostringstream message;
344  message << "Solid " << GetName() << " built using generic construct."
345  << G4endl << "Not applicable to the generic construct !";
346  G4Exception("G4GenericPolycone::Reset()", "GeomSolids1001",
347  JustWarning, message, "Parameters NOT resetted.");
348  return 1;
349 
350 }
351 
352 
353 //
354 // Inside
355 //
356 // This is an override of G4VCSGfaceted::Inside, created in order
357 // to speed things up by first checking with G4EnclosingCylinder.
358 //
360 {
361  //
362  // Quick test
363  //
364  if (enclosingCylinder->MustBeOutside(p)) return kOutside;
365 
366  //
367  // Long answer
368  //
369  return G4VCSGfaceted::Inside(p);
370 }
371 
372 
373 //
374 // DistanceToIn
375 //
376 // This is an override of G4VCSGfaceted::Inside, created in order
377 // to speed things up by first checking with G4EnclosingCylinder.
378 //
380  const G4ThreeVector &v ) const
381 {
382  //
383  // Quick test
384  //
385  if (enclosingCylinder->ShouldMiss(p,v))
386  return kInfinity;
387 
388  //
389  // Long answer
390  //
391  return G4VCSGfaceted::DistanceToIn( p, v );
392 }
393 
394 
395 //
396 // DistanceToIn
397 //
399 {
400  return G4VCSGfaceted::DistanceToIn(p);
401 }
402 
403 
404 //
405 // ComputeDimensions
406 //
407 /*void G4GenericPolycone::ComputeDimensions( G4VPVParameterisation* p,
408  const G4int n,
409  const G4VPhysicalVolume* pRep )
410 {
411  p->ComputeDimensions(*this,n,pRep);
412 }
413 */
414 //
415 // GetEntityType
416 //
418 {
419  return G4String("G4GenericPolycone");
420 }
421 
422 
423 //
424 // Make a clone of the object
425 //
427 {
428  return new G4GenericPolycone(*this);
429 }
430 
431 //
432 // Stream object contents to an output stream
433 //
434 std::ostream& G4GenericPolycone::StreamInfo( std::ostream& os ) const
435 {
436  G4int oldprc = os.precision(16);
437  os << "-----------------------------------------------------------\n"
438  << " *** Dump for solid - " << GetName() << " ***\n"
439  << " ===================================================\n"
440  << " Solid type: G4GenericPolycone\n"
441  << " Parameters: \n"
442  << " starting phi angle : " << startPhi/degree << " degrees \n"
443  << " ending phi angle : " << endPhi/degree << " degrees \n";
444  G4int i=0;
445 
446  os << " number of RZ points: " << numCorner << "\n"
447  << " RZ values (corners): \n";
448  for (i=0; i<numCorner; i++)
449  {
450  os << " "
451  << corners[i].r << ", " << corners[i].z << "\n";
452  }
453  os << "-----------------------------------------------------------\n";
454  os.precision(oldprc);
455 
456  return os;
457 }
458 
459 
460 
461 //
462 // GetPointOnSurface
463 //
465 {
466  return GetPointOnSurfaceGeneric();
467 
468 }
469 
470 //
471 // CreatePolyhedron
472 //
474 {
475  // The following code prepares for:
476  // HepPolyhedron::createPolyhedron(int Nnodes, int Nfaces,
477  // const double xyz[][3],
478  // const int faces_vec[][4])
479  // Here is an extract from the header file HepPolyhedron.h:
497  const G4int numSide =
499  * (endPhi - startPhi) / twopi) + 1;
500  G4int nNodes;
501  G4int nFaces;
502  typedef G4double double3[3];
503  double3* xyz;
504  typedef G4int int4[4];
505  int4* faces_vec;
506  if (phiIsOpen)
507  {
508  // Triangulate open ends. Simple ear-chopping algorithm...
509  // I'm not sure how robust this algorithm is (J.Allison).
510  //
511  std::vector<G4bool> chopped(numCorner, false);
512  std::vector<G4int*> triQuads;
513  G4int remaining = numCorner;
514  G4int iStarter = 0;
515  while (remaining >= 3) // Loop checking, 13.08.2015, G.Cosmo
516  {
517  // Find unchopped corners...
518  //
519  G4int A = -1, B = -1, C = -1;
520  G4int iStepper = iStarter;
521  do // Loop checking, 13.08.2015, G.Cosmo
522  {
523  if (A < 0) { A = iStepper; }
524  else if (B < 0) { B = iStepper; }
525  else if (C < 0) { C = iStepper; }
526  do // Loop checking, 13.08.2015, G.Cosmo
527  {
528  if (++iStepper >= numCorner) { iStepper = 0; }
529  }
530  while (chopped[iStepper]);
531  }
532  while (C < 0 && iStepper != iStarter);
533 
534  // Check triangle at B is pointing outward (an "ear").
535  // Sign of z cross product determines...
536  //
537  G4double BAr = corners[A].r - corners[B].r;
538  G4double BAz = corners[A].z - corners[B].z;
539  G4double BCr = corners[C].r - corners[B].r;
540  G4double BCz = corners[C].z - corners[B].z;
541  if (BAr * BCz - BAz * BCr < kCarTolerance)
542  {
543  G4int* tq = new G4int[3];
544  tq[0] = A + 1;
545  tq[1] = B + 1;
546  tq[2] = C + 1;
547  triQuads.push_back(tq);
548  chopped[B] = true;
549  --remaining;
550  }
551  else
552  {
553  do // Loop checking, 13.08.2015, G.Cosmo
554  {
555  if (++iStarter >= numCorner) { iStarter = 0; }
556  }
557  while (chopped[iStarter]);
558  }
559  }
560  // Transfer to faces...
561  //
562  nNodes = (numSide + 1) * numCorner;
563  nFaces = numSide * numCorner + 2 * triQuads.size();
564  faces_vec = new int4[nFaces];
565  G4int iface = 0;
566  G4int addition = numCorner * numSide;
567  G4int d = numCorner - 1;
568  for (G4int iEnd = 0; iEnd < 2; ++iEnd)
569  {
570  for (size_t i = 0; i < triQuads.size(); ++i)
571  {
572  // Negative for soft/auxiliary/normally invisible edges...
573  //
574  G4int a, b, c;
575  if (iEnd == 0)
576  {
577  a = triQuads[i][0];
578  b = triQuads[i][1];
579  c = triQuads[i][2];
580  }
581  else
582  {
583  a = triQuads[i][0] + addition;
584  b = triQuads[i][2] + addition;
585  c = triQuads[i][1] + addition;
586  }
587  G4int ab = std::abs(b - a);
588  G4int bc = std::abs(c - b);
589  G4int ca = std::abs(a - c);
590  faces_vec[iface][0] = (ab == 1 || ab == d)? a: -a;
591  faces_vec[iface][1] = (bc == 1 || bc == d)? b: -b;
592  faces_vec[iface][2] = (ca == 1 || ca == d)? c: -c;
593  faces_vec[iface][3] = 0;
594  ++iface;
595  }
596  }
597 
598  // Continue with sides...
599 
600  xyz = new double3[nNodes];
601  const G4double dPhi = (endPhi - startPhi) / numSide;
602  G4double phi = startPhi;
603  G4int ixyz = 0;
604  for (G4int iSide = 0; iSide < numSide; ++iSide)
605  {
606  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
607  {
608  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
609  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
610  xyz[ixyz][2] = corners[iCorner].z;
611  if (iSide == 0) // startPhi
612  {
613  if (iCorner < numCorner - 1)
614  {
615  faces_vec[iface][0] = ixyz + 1;
616  faces_vec[iface][1] = -(ixyz + numCorner + 1);
617  faces_vec[iface][2] = ixyz + numCorner + 2;
618  faces_vec[iface][3] = ixyz + 2;
619  }
620  else
621  {
622  faces_vec[iface][0] = ixyz + 1;
623  faces_vec[iface][1] = -(ixyz + numCorner + 1);
624  faces_vec[iface][2] = ixyz + 2;
625  faces_vec[iface][3] = ixyz - numCorner + 2;
626  }
627  }
628  else if (iSide == numSide - 1) // endPhi
629  {
630  if (iCorner < numCorner - 1)
631  {
632  faces_vec[iface][0] = ixyz + 1;
633  faces_vec[iface][1] = ixyz + numCorner + 1;
634  faces_vec[iface][2] = ixyz + numCorner + 2;
635  faces_vec[iface][3] = -(ixyz + 2);
636  }
637  else
638  {
639  faces_vec[iface][0] = ixyz + 1;
640  faces_vec[iface][1] = ixyz + numCorner + 1;
641  faces_vec[iface][2] = ixyz + 2;
642  faces_vec[iface][3] = -(ixyz - numCorner + 2);
643  }
644  }
645  else
646  {
647  if (iCorner < numCorner - 1)
648  {
649  faces_vec[iface][0] = ixyz + 1;
650  faces_vec[iface][1] = -(ixyz + numCorner + 1);
651  faces_vec[iface][2] = ixyz + numCorner + 2;
652  faces_vec[iface][3] = -(ixyz + 2);
653  }
654  else
655  {
656  faces_vec[iface][0] = ixyz + 1;
657  faces_vec[iface][1] = -(ixyz + numCorner + 1);
658  faces_vec[iface][2] = ixyz + 2;
659  faces_vec[iface][3] = -(ixyz - numCorner + 2);
660  }
661  }
662  ++iface;
663  ++ixyz;
664  }
665  phi += dPhi;
666  }
667 
668  // Last corners...
669 
670  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
671  {
672  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
673  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
674  xyz[ixyz][2] = corners[iCorner].z;
675  ++ixyz;
676  }
677  }
678  else // !phiIsOpen - i.e., a complete 360 degrees.
679  {
680  nNodes = numSide * numCorner;
681  nFaces = numSide * numCorner;;
682  xyz = new double3[nNodes];
683  faces_vec = new int4[nFaces];
684  const G4double dPhi = (endPhi - startPhi) / numSide;
685  G4double phi = startPhi;
686  G4int ixyz = 0, iface = 0;
687  for (G4int iSide = 0; iSide < numSide; ++iSide)
688  {
689  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
690  {
691  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
692  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
693  xyz[ixyz][2] = corners[iCorner].z;
694 
695  if (iSide < numSide - 1)
696  {
697  if (iCorner < numCorner - 1)
698  {
699  faces_vec[iface][0] = ixyz + 1;
700  faces_vec[iface][1] = -(ixyz + numCorner + 1);
701  faces_vec[iface][2] = ixyz + numCorner + 2;
702  faces_vec[iface][3] = -(ixyz + 2);
703  }
704  else
705  {
706  faces_vec[iface][0] = ixyz + 1;
707  faces_vec[iface][1] = -(ixyz + numCorner + 1);
708  faces_vec[iface][2] = ixyz + 2;
709  faces_vec[iface][3] = -(ixyz - numCorner + 2);
710  }
711  }
712  else // Last side joins ends...
713  {
714  if (iCorner < numCorner - 1)
715  {
716  faces_vec[iface][0] = ixyz + 1;
717  faces_vec[iface][1] = -(ixyz + numCorner - nFaces + 1);
718  faces_vec[iface][2] = ixyz + numCorner - nFaces + 2;
719  faces_vec[iface][3] = -(ixyz + 2);
720  }
721  else
722  {
723  faces_vec[iface][0] = ixyz + 1;
724  faces_vec[iface][1] = -(ixyz - nFaces + numCorner + 1);
725  faces_vec[iface][2] = ixyz - nFaces + 2;
726  faces_vec[iface][3] = -(ixyz - numCorner + 2);
727  }
728  }
729  ++ixyz;
730  ++iface;
731  }
732  phi += dPhi;
733  }
734  }
735  G4Polyhedron* polyhedron = new G4Polyhedron;
736  G4int problem = polyhedron->createPolyhedron(nNodes, nFaces, xyz, faces_vec);
737  delete [] faces_vec;
738  delete [] xyz;
739  if (problem)
740  {
741  std::ostringstream message;
742  message << "Problem creating G4Polyhedron for: " << GetName();
743  G4Exception("G4GenericPolycone::CreatePolyhedron()", "GeomSolids1002",
744  JustWarning, message);
745  delete polyhedron;
746  return 0;
747  }
748  else
749  {
750  return polyhedron;
751  }
752 }
753 
754 #endif
G4bool CrossesItself(G4double tolerance)
G4VCSGfaceted & operator=(const G4VCSGfaceted &source)
G4bool fRebuildPolyhedron
G4Polyhedron * fpPolyhedron
static const G4double kInfinity
Definition: geomdefs.hh:42
G4PolyconeSideRZ * corners
G4VSolid * Clone() const
Float_t d
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])
G4bool ShouldMiss(const G4ThreeVector &p, const G4ThreeVector &v) const
G4String name
Definition: TRTMaterials.hh:40
std::ostream & StreamInfo(std::ostream &os) const
double C(double temp)
int G4int
Definition: G4Types.hh:78
G4GenericPolycone(const G4String &name, G4double phiStart, G4double phiTotal, G4int numRZ, const G4double r[], const G4double z[])
G4String GetName() const
G4bool RemoveDuplicateVertices(G4double tolerance)
virtual EInside Inside(const G4ThreeVector &p) const
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
G4EnclosingCylinder * enclosingCylinder
void Create(G4double phiStart, G4double phiTotal, G4ReduciblePolygon *rz)
G4bool RemoveRedundantVertices(G4double tolerance)
G4Polyhedron * CreatePolyhedron() const
double A(double temperature)
G4bool MustBeOutside(const G4ThreeVector &p) const
void CopyStuff(const G4GenericPolycone &source)
bool G4bool
Definition: G4Types.hh:79
G4GeometryType GetEntityType() const
static const double twopi
Definition: G4SIunits.hh:75
G4VCSGface ** faces
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4GenericPolycone & operator=(const G4GenericPolycone &source)
G4ThreeVector GetPointOnSurfaceGeneric() const
static const G4double ab
EInside Inside(const G4ThreeVector &p) const
static G4int GetNumberOfRotationSteps()
EInside
Definition: geomdefs.hh:58
Definition: DoubConv.h:18
G4double Amin() const
G4int NumVertices() const
static const double degree
Definition: G4SIunits.hh:143
G4bool BisectedBy(G4double a1, G4double b1, G4double a2, G4double b2, G4double tolerance)
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
#define G4endl
Definition: G4ios.hh:61
G4double kCarTolerance
Definition: G4VSolid.hh:304
double G4double
Definition: G4Types.hh:76
G4ThreeVector GetPointOnSurface() const