Belle II Software development
pent_t Struct Referenceabstract

pentagon shape More...

Inheritance diagram for pent_t:
shape_t

Public Member Functions

void adjust ()
 adjust sizes to have flat sides
 
bool istrap () const override
 is trapped?
 
map< int, G4ThreeVector > make_verticies (double wrapthick) const
 create map of vertices
 
G4VSolid * get_tesselatedsolid (const string &prefix, double wrapthick, G4Translate3D &shift UNUSED) const override
 get tesselated solid
 
G4VSolid * get_extrudedsolid (const string &prefix, double wrapthick, G4Translate3D &shift UNUSED) const override
 get extruded solid
 
G4VSolid * get_trapezoid (const string &prefix, double wrapthick, G4Translate3D &shift) const override
 get trapezoid
 
G4VSolid * get_bellecrystal (const string &prefix, double wrapthick, G4Translate3D &shift UNUSED) const override
 get Belle crystal
 
G4VSolid * get_solid (const std::string &prefix, double wrapthick, G4Translate3D &shift) const
 get solid
 
virtual G4VSolid * get_tesselatedsolid (const std::string &prefix, double wrapthick, G4Translate3D &shift) const =0
 get tesselated solid
 
virtual G4VSolid * get_extrudedsolid (const std::string &prefix, double wrapthick, G4Translate3D &shift) const =0
 get extruded solid
 
virtual G4VSolid * get_bellecrystal (const std::string &prefix, double wrapthick, G4Translate3D &shift) const =0
 get Belle crystal
 

Public Attributes

union {
   struct {
      double   A
 
      double   C
 
      double   D
 
      double   a
 
      double   c
 
      double   d
 
      double   B
 
      double   b
 
      double   H_aA
 
      double   H_dD
 
      double   dg13
 
      double   dg24
 
      double   dg57
 
      double   dg68
 
      double   a1
 
      double   a4
 
      double   a2
 
      double   a3
 
      double   a9
 
      double   Volume
 
      double   Weight
 
   } 
 
   double   t [21]
 
}; 
 
bool _adjusted
 are sizes adjusted?
 
int nshape
 shapes
 

Detailed Description

pentagon shape

Definition at line 416 of file shapes.cc.

Constructor & Destructor Documentation

◆ pent_t()

pent_t ( )
inline

Definition at line 424 of file shapes.cc.

424: _adjusted(false) {}
bool _adjusted
are sizes adjusted?
Definition: shapes.cc:423

◆ ~pent_t()

virtual ~pent_t ( )
inlinevirtual

Definition at line 425 of file shapes.cc.

425{}

Member Function Documentation

◆ adjust()

void adjust ( )
inline

adjust sizes to have flat sides

Definition at line 428 of file shapes.cc.

429 {
430 if (!_adjusted) {
431 // adjust sizes to have flat sides
432 H_dD -= 0.00005551197484235;
433 B += 0.0011245236213532729;
434 b += -0.00044853029662963;
435 _adjusted = true;
436 }
437 }

◆ get_bellecrystal()

G4VSolid * get_bellecrystal ( const string &  prefix,
double  wrapthick,
G4Translate3D &shift  UNUSED 
) const
inlineoverride

get Belle crystal

Definition at line 608 of file shapes.cc.

609 {
610 map<int, G4ThreeVector> v = make_verticies(wrapthick);
611
612 string name(prefix);
613 name += to_string(nshape);
614
615 G4ThreeVector pt[10];
616 pt[0] = v[1];
617 pt[1] = v[2];
618 pt[2] = v[9];
619 pt[3] = v[3];
620 pt[4] = v[4];
621 pt[5] = v[5];
622 pt[6] = v[6];
623 pt[7] = v[10];
624 pt[8] = v[7];
625 pt[9] = v[8];
626
627 G4VSolid* shape = new BelleCrystal(name.c_str(), 5, pt);
628 // cout<<name<<" "<<shape->GetCubicVolume()*1e-3<<" "<<Volume<<endl;
629 // for(int i=0;i<100000;i++){
630 // G4ThreeVector a = shape->GetPointOnSurface();
631 // cout<<a.x()<<" "<<a.y()<<" "<<a.z()<<endl;
632 // }
633 // exit(0);
634 // if(prefix.find("crystal")!=string::npos){
635 // pentspeed(shape);
636 // exit(0);
637 // }
638 // if(prefix.find("wrap")!=string::npos){
639 // pentspeed2(shape);
640 // exit(0);
641 // }
642 return shape;
643 }
map< int, G4ThreeVector > make_verticies(double wrapthick) const
create map of vertices
Definition: shapes.cc:443
int nshape
shapes
Definition: shapes.h:30

◆ get_extrudedsolid()

G4VSolid * get_extrudedsolid ( const string &  prefix,
double  wrapthick,
G4Translate3D &shift  UNUSED 
) const
inlineoverride

get extruded solid

Definition at line 512 of file shapes.cc.

513 {
514 map<int, G4ThreeVector> v = make_verticies(wrapthick);
515
516 string name(prefix);
517 name += to_string(nshape);
518
519 std::vector<G4TwoVector> p1, p2;
520 p1.push_back(G4TwoVector(v[1].x(), v[1].y()));
521 p1.push_back(G4TwoVector(v[2].x(), v[2].y()));
522 p1.push_back(G4TwoVector(v[9].x(), v[9].y()));
523 p1.push_back(G4TwoVector(v[3].x(), v[3].y()));
524 p1.push_back(G4TwoVector(v[4].x(), v[4].y()));
525
526 p2.push_back(G4TwoVector(v[1 + 4].x(), v[1 + 4].y()));
527 p2.push_back(G4TwoVector(v[2 + 4].x(), v[2 + 4].y()));
528 p2.push_back(G4TwoVector(v[ 10].x(), v[ 10].y()));
529 p2.push_back(G4TwoVector(v[3 + 4].x(), v[3 + 4].y()));
530 p2.push_back(G4TwoVector(v[4 + 4].x(), v[4 + 4].y()));
531
532 double sum = 0, smin = 1e9, smax = -1e9;
533 int count = 0;
534 for (int i = 0; i < 5; i++) {
535 for (int j = i + 1; j < 5; j++) {
536 double s2 = (p2[j] - p2[i]).mag2() / (p1[j] - p1[i]).mag2();
537 double s = sqrt(s2);
538 sum += s;
539 if (s > smax) smax = s;
540 if (s < smin) smin = s;
541 count++;
542 }
543 }
544 double ave = sum / count;
545
546 double scale = ave;
547 G4TwoVector off1(0, 0), off2(scale * p1[0].x() - p2[0].x(), scale * p1[0].y() - p2[0].y());
548
549 return new G4ExtrudedSolid(name, p1, abs(v[1].z()), off1, 1, -off2, scale);
550 }
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28

◆ get_solid()

G4VSolid * get_solid ( const std::string &  prefix,
double  wrapthick,
G4Translate3D &  shift 
) const
inherited

get solid

Definition at line 73 of file shapes.cc.

74 {
75 return get_bellecrystal(prefix, wrapthick, shift);
76 }
virtual G4VSolid * get_bellecrystal(const std::string &prefix, double wrapthick, G4Translate3D &shift) const =0
get Belle crystal

◆ get_tesselatedsolid()

G4VSolid * get_tesselatedsolid ( const string &  prefix,
double  wrapthick,
G4Translate3D &shift  UNUSED 
) const
inlineoverride

get tesselated solid

Definition at line 480 of file shapes.cc.

481 {
482 if (nshape != 36) return nullptr; // only one crystal has pentagon shape
483
484 map<int, G4ThreeVector> v = make_verticies(wrapthick);
485
486 string name(prefix);
487 name += to_string(nshape);
488 G4TessellatedSolid* s = new G4TessellatedSolid(name.c_str());
489
490 // Now add the facets to the solid
491 // top plane
492 s->AddFacet(new G4QuadrangularFacet(v[1], v[4], v[3], v[2], ABSOLUTE));
493 s->AddFacet(new G4TriangularFacet(v[2], v[3], v[9], ABSOLUTE));
494
495 //bottom plane
496 s->AddFacet(new G4QuadrangularFacet(v[5], v[6], v[7], v[8], ABSOLUTE));
497 s->AddFacet(new G4TriangularFacet(v[6], v[10], v[7], ABSOLUTE));
498
499 //sides
500 s->AddFacet(new G4QuadrangularFacet(v[1], v[2], v[6], v[5], ABSOLUTE));
501 s->AddFacet(new G4QuadrangularFacet(v[2], v[9], v[10], v[6], ABSOLUTE));
502 s->AddFacet(new G4QuadrangularFacet(v[9], v[3], v[7], v[10], ABSOLUTE));
503 s->AddFacet(new G4QuadrangularFacet(v[3], v[4], v[8], v[7], ABSOLUTE));
504 s->AddFacet(new G4QuadrangularFacet(v[4], v[1], v[5], v[8], ABSOLUTE));
505
506 // Finally declare the solid is complete
507 s->SetSolidClosed(true);
508 return s;
509 }

◆ get_trapezoid()

G4VSolid * get_trapezoid ( const string &  prefix,
double  wrapthick,
G4Translate3D &  shift 
) const
inlineoverridevirtual

get trapezoid

Implements shape_t.

Definition at line 553 of file shapes.cc.

554 {
555 if (nshape != 36) return nullptr; // only one crystal has pentagon shape
556
557 map<int, G4ThreeVector> v = make_verticies(wrapthick);
558
559 string name(prefix);
560 name += to_string(nshape);
561
562 G4ThreeVector pt[8];
563 pt[0] = v[4];
564 pt[1] = v[1];
565 pt[2] = v[3];
566 pt[3] = v[2];
567 pt[4] = v[8];
568 pt[5] = v[5];
569 pt[6] = v[7];
570 pt[7] = v[6];
571
572 auto alignz = [&](int i, int j) {pt[j].setZ(pt[i].z());};
573 auto aligny = [&](int i, int j) {pt[j].setY(pt[i].y());};
574
575 alignz(0, 1);
576 alignz(0, 2);
577 alignz(0, 3);
578
579 alignz(4, 1 + 4);
580 alignz(4, 2 + 4);
581 alignz(4, 3 + 4);
582
583 aligny(0, 1);
584 aligny(2, 3);
585
586 aligny(4, 5);
587 aligny(6, 7);
588
589
590 double dx = pt[0].x() + pt[1].x() + pt[4].x() + pt[5].x() + pt[2].x() + pt[3].x() + pt[6].x() + pt[7].x();
591 dx /= 8;
592 double dy = pt[0].y() + pt[2].y() + pt[4].y() + pt[6].y();
593 dy /= 4;
594 for (int j = 0; j < 8; j++) {
595 pt[j].setX(pt[j].x() - dx);
596 pt[j].setY(pt[j].y() - dy);
597 }
598 shift = G4Translate3D(dx, dy, 0);
599
600 // cout<<name<<" "<<dx<<" "<<dy<<endl;
601
602 G4VSolid* shape = new G4Trap(name.c_str(), pt);
603 // G4VSolid *shape = new BelleCrystal(name.c_str(), 4, pt);
604 return shape;
605 }

◆ istrap()

bool istrap ( ) const
inlineoverridevirtual

is trapped?

Implements shape_t.

Definition at line 440 of file shapes.cc.

440{ return false;}

◆ make_verticies()

map< int, G4ThreeVector > make_verticies ( double  wrapthick) const
inline

create map of vertices

Definition at line 443 of file shapes.cc.

444 {
445 assert(_adjusted);
446
447 double h2 = cosd(a1 - 90) * D / 2;
448 map<int, G4ThreeVector> v;
449 v[5] = G4ThreeVector(-A / 2, -h2, -150);
450 v[6] = v[5] + moveto(D, a1);
451 v[10] = v[6] + moveto(B, a1 + a2 - 180);
452 v[8] = G4ThreeVector(A / 2, -h2, -150);
453 v[7] = v[8] + moveto(D, 180 - a1);
454
455 v[1] = v[5] + G4ThreeVector((H_aA * cosd(a1) + H_dD) / sind(a1), H_aA, 300);
456 v[2] = v[1] + moveto(d, a1);
457 v[9] = v[2] + moveto(b, a1 + a2 - 180);
458 v[4] = v[1] + moveto(a, 0);
459 v[3] = v[4] + moveto(d, 180 - a1);
460
461 for (int j = 1; j <= 10; j++) v[j] = G4ThreeVector(v[j].x(), -v[j].y(), -v[j].z());
462 if (wrapthick != 0) {
463 map<int, G4ThreeVector> nv;
464 nv[ 1] = newvertex(wrapthick, v[1], v[5], v[2], v[4]);
465 nv[ 2] = newvertex(wrapthick, v[2], v[6], v[9], v[1]);
466 nv[ 9] = newvertex(wrapthick, v[9], v[10], v[3], v[2]);
467 nv[ 3] = newvertex(wrapthick, v[3], v[7], v[4], v[9]);
468 nv[ 4] = newvertex(wrapthick, v[4], v[8], v[1], v[3]);
469 nv[ 5] = newvertex(wrapthick, v[5], v[1], v[8], v[6]);
470 nv[ 6] = newvertex(wrapthick, v[6], v[2], v[5], v[10]);
471 nv[10] = newvertex(wrapthick, v[10], v[9], v[6], v[7]);
472 nv[ 7] = newvertex(wrapthick, v[7], v[3], v[10], v[8]);
473 nv[ 8] = newvertex(wrapthick, v[8], v[4], v[7], v[5]);
474 std::swap(nv, v);
475 }
476 return v;
477 }

Member Data Documentation

◆ _adjusted

bool _adjusted

are sizes adjusted?

Definition at line 423 of file shapes.cc.

◆ a

double a

Definition at line 419 of file shapes.cc.

◆ A

double A

Definition at line 419 of file shapes.cc.

◆ a1

double a1

Definition at line 419 of file shapes.cc.

◆ a2

double a2

Definition at line 419 of file shapes.cc.

◆ a3

double a3

Definition at line 419 of file shapes.cc.

◆ a4

double a4

Definition at line 419 of file shapes.cc.

◆ a9

double a9

Definition at line 419 of file shapes.cc.

◆ b

double b

Definition at line 419 of file shapes.cc.

◆ B

double B

Definition at line 419 of file shapes.cc.

◆ c

double c

Definition at line 419 of file shapes.cc.

◆ C

double C

Definition at line 419 of file shapes.cc.

◆ d

double d

Definition at line 419 of file shapes.cc.

◆ D

double D

Definition at line 419 of file shapes.cc.

◆ dg13

double dg13

Definition at line 419 of file shapes.cc.

◆ dg24

double dg24

Definition at line 419 of file shapes.cc.

◆ dg57

double dg57

Definition at line 419 of file shapes.cc.

◆ dg68

double dg68

Definition at line 419 of file shapes.cc.

◆ H_aA

double H_aA

Definition at line 419 of file shapes.cc.

◆ H_dD

double H_dD

Definition at line 419 of file shapes.cc.

◆ nshape

int nshape
inherited

shapes

Definition at line 30 of file shapes.h.

◆ t

double t[21]

Definition at line 421 of file shapes.cc.

◆ Volume

double Volume

Definition at line 419 of file shapes.cc.

◆ Weight

double Weight

Definition at line 419 of file shapes.cc.


The documentation for this struct was generated from the following file: