Belle II Software development
quadrilateral_endcap_t Struct Referenceabstract

quadrilateral struct for end cap More...

Inheritance diagram for quadrilateral_endcap_t:
quadrilateral_t shape_t

Public Member Functions

bool istrap () const override
 is trapped?
 
map< int, G4ThreeVector > make_verticies (double wrapthick) const override
 create map of vertices
 
G4VSolid * get_tesselatedsolid (const string &prefix, double wrapthick, G4Translate3D &shift UNUSED) const override
 get tesselated solid
 
virtual G4VSolid * get_tesselatedsolid (const std::string &prefix, double wrapthick, G4Translate3D &shift) const =0
 get tesselated solid
 
G4VSolid * get_extrudedsolid (const string &prefix, double wrapthick, G4Translate3D &shift UNUSED) const override
 get extruded solid
 
virtual G4VSolid * get_extrudedsolid (const std::string &prefix, double wrapthick, G4Translate3D &shift) const =0
 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
 
virtual G4VSolid * get_bellecrystal (const std::string &prefix, double wrapthick, G4Translate3D &shift) const =0
 get Belle crystal
 
G4VSolid * get_solid (const std::string &prefix, double wrapthick, G4Translate3D &shift) const
 get solid
 

Public Attributes

union {
   struct {
      double   A
 
      double   B
 
      double   C
 
      double   D
 
      double   a
 
      double   b
 
      double   c
 
      double   d
 
      double   H_aA
 
      double   H_dD
 
      double   dg13
 
      double   dg24
 
      double   dg57
 
      double   dg68
 
      double   a1
 
      double   a2
 
      double   a3
 
      double   a4
 
      double   Volume
 
      double   Weight
 
   } 
 
   double   t [20]
 
}; 
 
int nshape
 shapes
 

Detailed Description

quadrilateral struct for end cap

Definition at line 349 of file shapes.cc.

Constructor & Destructor Documentation

◆ quadrilateral_endcap_t()

Definition at line 356 of file shapes.cc.

356{}

◆ ~quadrilateral_endcap_t()

virtual ~quadrilateral_endcap_t ( )
inlinevirtual

Definition at line 357 of file shapes.cc.

357{}

Member Function Documentation

◆ get_bellecrystal()

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

get Belle crystal

Definition at line 257 of file shapes.cc.

258 {
259 map<int, G4ThreeVector> v = make_verticies(wrapthick);
260
261 string name(prefix);
262 name += to_string(nshape);
263
264 G4ThreeVector pt[8];
265 pt[0] = v[4];
266 pt[1] = v[1];
267 pt[2] = v[3];
268 pt[3] = v[2];
269 pt[4] = v[8];
270 pt[5] = v[5];
271 pt[6] = v[7];
272 pt[7] = v[6];
273
274 for (int i = 0; i < 8; i++) pt[i] = v[i + 1];
275 G4VSolid* shape = new BelleCrystal(name.c_str(), 4, pt);
276 // cout<<name<<" "<<shape->GetCubicVolume()*1e-3<<" "<<Volume<<endl;
277 // for(int i=0;i<100000;i++){
278 // G4ThreeVector a = shape->GetPointOnSurface();
279 // cout<<a.x()<<" "<<a.y()<<" "<<a.z()<<endl;
280 // }
281 // exit(0);
282 return shape;
283 }
virtual map< int, G4ThreeVector > make_verticies(double wrapthick) const =0
create map of vertices
int nshape
shapes
Definition: shapes.h:30

◆ get_extrudedsolid()

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

get extruded solid

Definition at line 146 of file shapes.cc.

147 {
148 map<int, G4ThreeVector> v = make_verticies(wrapthick);
149
150 string name(prefix);
151 name += to_string(nshape);
152
153 std::vector<G4TwoVector> p1, p2;
154 p1.push_back(G4TwoVector(v[1].x(), v[1].y()));
155 p1.push_back(G4TwoVector(v[2].x(), v[2].y()));
156 p1.push_back(G4TwoVector(v[3].x(), v[3].y()));
157 p1.push_back(G4TwoVector(v[4].x(), v[4].y()));
158
159 p2.push_back(G4TwoVector(v[1 + 4].x(), v[1 + 4].y()));
160 p2.push_back(G4TwoVector(v[2 + 4].x(), v[2 + 4].y()));
161 p2.push_back(G4TwoVector(v[3 + 4].x(), v[3 + 4].y()));
162 p2.push_back(G4TwoVector(v[4 + 4].x(), v[4 + 4].y()));
163
164 double sum = 0, smin = 1e9, smax = -1e9;
165 for (int i = 0; i < 4; i++) {
166 for (int j = i + 1; j < 4; j++) {
167 double s2 = (p2[j] - p2[i]).mag2() / (p1[j] - p1[i]).mag2();
168 double s = sqrt(s2);
169 sum += s;
170 if (s > smax) smax = s;
171 if (s < smin) smin = s;
172 }
173 }
174 double ave = sum / 6;
175
176 double scale = ave;
177 G4TwoVector off1(0, 0), off2(scale * p1[0].x() - p2[0].x(), scale * p1[0].y() - p2[0].y());
178
179 return new G4ExtrudedSolid(name, p1, abs(v[1].z()), off1, 1, -off2, scale);
180 }
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
inlineoverrideinherited

get tesselated solid

Definition at line 121 of file shapes.cc.

122 {
123 map<int, G4ThreeVector> v = make_verticies(wrapthick);
124
125 string name(prefix);
126 name += to_string(nshape);
127 G4TessellatedSolid* s = new G4TessellatedSolid(name.c_str());
128
129 // Now add the facets to the solid
130 // top plane
131 s->AddFacet(new G4QuadrangularFacet(v[1], v[4], v[3], v[2], ABSOLUTE));
132 // bottom plane
133 s->AddFacet(new G4QuadrangularFacet(v[5], v[6], v[7], v[8], ABSOLUTE));
134 // lateral sides
135 s->AddFacet(new G4QuadrangularFacet(v[1], v[2], v[6], v[5], ABSOLUTE));
136 s->AddFacet(new G4QuadrangularFacet(v[2], v[3], v[7], v[6], ABSOLUTE));
137 s->AddFacet(new G4QuadrangularFacet(v[3], v[4], v[8], v[7], ABSOLUTE));
138 s->AddFacet(new G4QuadrangularFacet(v[4], v[1], v[5], v[8], ABSOLUTE));
139
140 // Finally declare the solid is complete
141 s->SetSolidClosed(true);
142 return s;
143 }

◆ get_trapezoid()

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

get trapezoid

Implements shape_t.

Definition at line 183 of file shapes.cc.

184 {
185 map<int, G4ThreeVector> v = make_verticies(wrapthick);
186
187 // make sides @ +-Y parallel
188 G4ThreeVector b = v[4] - v[1];
189 G4ThreeVector d12 = v[2] - v[1], d43 = v[3] - v[4];
190 double b2 = b.mag2();
191 double h1 = b.cross(d12).mag2();
192 double h4 = b.cross(d43).mag2();
193 // if(abs(h1/h4-1)<0.02) cout<<"Shape = "<<nshape<<" "<<h1<<" "<<h4<<" "<<h1/h4-1<<endl;
194 if (h1 < h4) {
195 G4ThreeVector d24 = v[4] - v[2];
196 double d43b = d43 * b, s = (b2 * (d24 * d43) - (d24 * b) * d43b) / (d43b * d43b - b2 * d43.mag2());
197 v[3] = v[4] + s * d43;
198 v[7] = v[8] + s * (v[7] - v[8]);
199 } else {
200 G4ThreeVector d31 = v[1] - v[3];
201 double d12b = d12 * b, s = (b2 * (d31 * d12) - (d31 * b) * d12b) / (d12b * d12b - b2 * d12.mag2());
202 v[2] = v[1] + s * d12;
203 v[6] = v[5] + s * (v[6] - v[5]);
204 }
205
206 string name(prefix);
207 name += to_string(nshape);
208
209 G4ThreeVector pt[8];
210 pt[0] = v[4];
211 pt[1] = v[1];
212 pt[2] = v[3];
213 pt[3] = v[2];
214 pt[4] = v[8];
215 pt[5] = v[5];
216 pt[6] = v[7];
217 pt[7] = v[6];
218
219 auto alignz = [&](int i, int j) {pt[j].setZ(pt[i].z());};
220 auto aligny = [&](int i, int j) {pt[j].setY(pt[i].y());};
221
222 alignz(0, 1);
223 alignz(0, 2);
224 alignz(0, 3);
225
226 alignz(4, 1 + 4);
227 alignz(4, 2 + 4);
228 alignz(4, 3 + 4);
229
230 aligny(0, 1);
231 aligny(2, 3);
232
233 aligny(4, 5);
234 aligny(6, 7);
235
236 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();
237 dx /= 8;
238 double dy = pt[0].y() + pt[2].y() + pt[4].y() + pt[6].y();
239 dy /= 4;
240 for (int j = 0; j < 8; j++) {
241 pt[j].setX(pt[j].x() - dx);
242 pt[j].setY(pt[j].y() - dy);
243 }
244
245 // int oprec = cout.precision(17);
246 // cout<<dx<<" "<<dy<<endl;
247 // cout.precision(oprec);
248
249 shift = G4Translate3D(dx, dy, 0);
250 G4VSolid* shape = new G4Trap(name.c_str(), pt);
251 // cout<<name<<" "<<shape->GetCubicVolume()*1e-3<<" "<<Volume<<endl;
252 // G4VSolid *shape = new BelleCrystal(name.c_str(), 4, pt);
253 return shape;
254 }

◆ istrap()

bool istrap ( ) const
inlineoverridevirtual

is trapped?

Implements shape_t.

Definition at line 360 of file shapes.cc.

361 {
362 double h1 = sind(a1) * D, h4 = sind(a4) * C;
363 return abs(h1 - h4) < 0.01 * h1;
364 }

◆ make_verticies()

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

create map of vertices

Implements quadrilateral_t.

Definition at line 367 of file shapes.cc.

368 {
369 double minh = std::min(sind(a1) * D, sind(a4) * C);
370 double h2 = minh / 2;
371 double db = (tand(a1 - 90) - tand(a4 - 90)) * h2 / 2;
372
373 map<int, G4ThreeVector> v;
374 v[5] = G4ThreeVector(-A / 2 + db, -h2, -150);
375 v[6] = v[5] + moveto(D, a1);
376 v[8] = G4ThreeVector(A / 2 + db, -h2, -150);
377 v[7] = v[8] + moveto(C, 180 - a4);
378
379 // adjust position of v[2],v[3],v[7],v[6] to have all 4 points in a plane, other combination already in a plane by the construction
380 G4ThreeVector vB = v[7] - v[6], va(a, 0, 0), vd = moveto(d, a1), vc = moveto(c, 180 - a4);
381 double delta = vB.cross(va + vc - vd).z() / vB.cross(vc + vd).z(); // delta should be very small ~10^-6 or less
382
383 vd *= 1 + delta;
384 vc *= 1 - delta;
385
386 v[1] = v[5] + G4ThreeVector((H_aA * cosd(a1) + H_dD) / sind(a1), H_aA, 300);
387 v[2] = v[1] + vd;
388 v[4] = v[1] + va;
389 v[3] = v[4] + vc;
390
391 for (int j = 1; j <= 8; j++) v[j] = G4ThreeVector(v[j].x(), -v[j].y(), -v[j].z());
392
393 // G4ThreeVector c0 = centerofgravity(v, 1, 4);
394 // G4ThreeVector c1 = centerofgravity(v, 5, 4);
395 // G4ThreeVector cz = 0.5*(c0+c1);
396 // cout<<c0<<" "<<c1<<" "<<cz<<endl;
397
398 if (wrapthick != 0) {
399 map<int, G4ThreeVector> nv;
400 nv[1] = newvertex(wrapthick, v[1], v[5], v[2], v[4]);
401 nv[2] = newvertex(wrapthick, v[2], v[6], v[3], v[1]);
402 nv[3] = newvertex(wrapthick, v[3], v[7], v[4], v[2]);
403 nv[4] = newvertex(wrapthick, v[4], v[8], v[1], v[3]);
404 nv[5] = newvertex(wrapthick, v[5], v[1], v[8], v[6]);
405 nv[6] = newvertex(wrapthick, v[6], v[2], v[5], v[7]);
406 nv[7] = newvertex(wrapthick, v[7], v[3], v[6], v[8]);
407 nv[8] = newvertex(wrapthick, v[8], v[4], v[7], v[5]);
408 std::swap(nv, v);
409 }
410 // if(nshape==2){ for(int j=1;j<=8;j++) cout<<v[j]<<" "; cout<<endl;}
411 return v;
412 }

Member Data Documentation

◆ a

double a

Definition at line 352 of file shapes.cc.

◆ A

double A

Definition at line 352 of file shapes.cc.

◆ a1

double a1

Definition at line 352 of file shapes.cc.

◆ a2

double a2

Definition at line 352 of file shapes.cc.

◆ a3

double a3

Definition at line 352 of file shapes.cc.

◆ a4

double a4

Definition at line 352 of file shapes.cc.

◆ B

double B

Definition at line 352 of file shapes.cc.

◆ b

double b

Definition at line 352 of file shapes.cc.

◆ C

double C

Definition at line 352 of file shapes.cc.

◆ c

double c

Definition at line 352 of file shapes.cc.

◆ D

double D

Definition at line 352 of file shapes.cc.

◆ d

double d

Definition at line 352 of file shapes.cc.

◆ dg13

double dg13

Definition at line 352 of file shapes.cc.

◆ dg24

double dg24

Definition at line 352 of file shapes.cc.

◆ dg57

double dg57

Definition at line 352 of file shapes.cc.

◆ dg68

double dg68

Definition at line 352 of file shapes.cc.

◆ H_aA

double H_aA

Definition at line 352 of file shapes.cc.

◆ H_dD

double H_dD

Definition at line 352 of file shapes.cc.

◆ nshape

int nshape
inherited

shapes

Definition at line 30 of file shapes.h.

◆ t

double t[20]

Definition at line 354 of file shapes.cc.

◆ Volume

double Volume

Definition at line 352 of file shapes.cc.

◆ Weight

double Weight

Definition at line 352 of file shapes.cc.


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