ATLAS Offline Software
Loading...
Searching...
No Matches
LArWheelSliceSolidTests.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include <iostream>
6#include <stdexcept>
7#include "boost/io/ios_state.hpp"
8#include <map>
9
10#include "TRandom3.h"
11#include "TF1.h"
12#include "TNtupleD.h"
13#include "TFile.h"
14// For root version ifdef
15#include "TROOT.h"
16
18#include "LArWheelSliceSolid.h"
19#include "CLHEP/Units/SystemOfUnits.h"
20#include "CLHEP/Units/PhysicalConstants.h"
21//#define LOCAL_DEBUG 1
22#include<stdlib.h>
23
26
27static double IntPrecision = 0.0001;
28
29EInside LArWheelSliceSolid::Inside_accordion(const G4ThreeVector &p) const
30{
31 G4ThreeVector p1 = p;
32 int p1_fan = 0;
33 G4double d = m_FanHalfThickness - fabs(GetCalculator()->DistanceToTheNearestFan(p1, p1_fan));
34 if(d > s_Tolerance) return kInside;
35 if(d > -s_Tolerance) return kSurface;
36 return kOutside;
37}
38
40#ifdef LOCAL_DEBUG
41 const char *m) const
42#else
43 const char *) const
44#endif
45 {
46 p[0] = 0.; p[1] = m_Rmin;
47 if(m_Zmin < s_Tolerance) p[2] = 0.;
48 else if(m_Calculator->GetWheelThickness() - m_Zmax < s_Tolerance) p[2] = m_Zmax;
49 else p[2] = (m_Zmin + m_Zmax) * 0.5;
50 int p_fan = 0;
51 GetCalculator()->DistanceToTheNearestFan(p, p_fan);
52
53#ifdef LOCAL_DEBUG
54 if(m) std::cout << m << std::endl;
55#endif
56 }
57
58static void get_r(const G4VSolid *p, G4double z,
59 G4double &rmin, G4double &rmax
60 )
61{
62 G4ThreeVector from(10.*CLHEP::m, 0., z);
63 rmax = from[0] - p->DistanceToIn(from, G4ThreeVector(-1., 0., 0.));
64 from[0] = 0.;
65 rmin = p->DistanceToIn(from, G4ThreeVector(1., 0., 0.));
66}
67
68static TRandom *rnd = 0;
70{
71 if(rnd == 0) rnd = new TRandom3(0);
72 G4double r = rnd->Uniform();
73 G4ThreeVector p(0., 0., 0.);
74
75 G4double level1 = .980;
76 G4double level2 = .993;
77 const char *v = getenv("LARWHEELSLICESOLID_TEST_MODE_LEVEL1");
78 if(v) level1 = atof(v);
79 const char *v1 = getenv("LARWHEELSLICESOLID_TEST_MODE_LEVEL2");
80 if(v1) level2 = atof(v1);
81
82#if LOCAL_DEBUG > 1
83 std::cout << "LWS::GPOS " << r << std::endl;
84#endif
85
86 if(r <= level1){
88 } else if(r <= level2){
90 } else if(r <= 1.){
92 } else {
93 G4Exception(
94 "LArWheelSliceSolid", "Rnd generator error",
95 FatalException, "GetPointOnSurface: Wrong data from rnd generator"
96 );
97 }
98 return p;
99}
100
102{
103 p[0] = 0.; p[1] = 0.; p[2] = rnd->Uniform(m_Zmin, m_Zmax);
104
105 G4double rmin, rmax;
106 get_r(m_BoundingShape, p[2], rmin, rmax);
107
108 p[1] = rnd->Uniform(rmin, rmax);
109 p.setPhi(rnd->Uniform(0., CLHEP::twopi));
110 G4double dphi = p.phi();
111 int p_fan = 0;
112 G4double d = GetCalculator()->DistanceToTheNearestFan(p, p_fan);
113 dphi -= p.phi();
114
115 G4int side = 0;
116 if(d < 0.) side = -1;
117 if(d >= 0.) side = 1;
118
119 G4double a = GetCalculator()->AmplitudeOfSurface(p, side, p_fan);
120 p[0] = a;
121
122 p.rotateZ(dphi);
123
124 if(m_BoundingShape->Inside(p) == kOutside){
125 G4ThreeVector D = p; D[2] = 0.;
126 G4double d1 = m_BoundingShape->DistanceToIn(p, D);
127 if(d1 > 10.*CLHEP::m){
128 D *= -1;
129 d1 = m_BoundingShape->DistanceToIn(p, D);
130 }
131 if(d1 > 10.*CLHEP::m){
132 set_failover_point(p, "acc fail0");
133 return;
134 }
135 d1 *= 2.;
136
137 G4ThreeVector B = p + D * d1;
138 G4double dphi = B.phi();
139 int B_fan = 0;
141 dphi -= B.phi();
142
143 B[0] = GetCalculator()->AmplitudeOfSurface(B, side, B_fan);
144 B.rotateZ(dphi);
145 EInside Bi = m_BoundingShape->Inside(B);
146 if(Bi == kSurface){
147 p = B;
148 return;
149 }
150 if(Bi == kOutside){
151 set_failover_point(p, "acc fail1");
152 return;
153 }
154 G4ThreeVector D1 = (p - B).unit();
155 G4ThreeVector X = B + D1 * m_BoundingShape->DistanceToOut(B, D1);
156 if(Inside(X) == kSurface){
157 p = X;
158 } else { // failed
159 set_failover_point(p, "acc fail2");
160 return;
161 }
162 }
163}
164
166{
167 const G4double z = rnd->Uniform(m_Zmin, m_Zmax);
168 G4double rmin, rmax;
169 get_r(m_BoundingShape, z, rmin, rmax);
170 const bool inner = rnd->Uniform() > 0.5? true: false;
171
172 p[0] = 0.; p[1] = inner? rmin: rmax; p[2] = z;
173 p.setPhi(rnd->Uniform(0., CLHEP::twopi));
174 G4double dphi = p.phi();
175 int p_fan = 0;
177 dphi -= p.phi();
178
179 const G4double r = p[1];
180
181 G4ThreeVector A1(0., r, z);
182 A1[0] = GetCalculator()->AmplitudeOfSurface(A1, -1, p_fan);
183 A1.rotateZ(dphi);
184 EInside A1i = m_BoundingShape->Inside(A1);
185 // EInside A1a = Inside_accordion(A1);
186 //std::cout << "A1: " << A1i << " " << A1a << std::endl;
187 if(A1i == kSurface){
188 //std::cout << "got A1" << std::endl;
189 p = A1;
190 return;
191 }
192
193 G4ThreeVector A2(0., r, z);
194 A2[0] = GetCalculator()->AmplitudeOfSurface(A2, 1, p_fan);
195 A2.rotateZ(dphi);
196 EInside A2i = m_BoundingShape->Inside(A2);
197 // EInside A2a = Inside_accordion(A2);
198 //std::cout << "A2: " << A2i << " " << A2a << std::endl;
199 if(A2i == kSurface){
200 //std::cout << "got A2" << std::endl;
201 p = A2;
202 return;
203 }
204
205 if(A1i != A2i){
206 if(A2i == kOutside){
207 std::swap(A1, A2);
208 std::swap(A1i, A2i);
209 }
210 // here A1 is outside BP, A2 is inside BP
211 G4ThreeVector d = (A2 - A1).unit();
212 p = A1 + d * m_BoundingShape->DistanceToIn(A1, d);
213 //std::cout << "got A1<->A2" << std::endl;
214 return;
215 }
216 // here A1i == A2i
217
218 G4double step;
219 if(A1i == kInside){
220 G4double d1 = m_BoundingShape->DistanceToOut(A1);
221 G4double d2 = m_BoundingShape->DistanceToOut(A2);
222 step = d1 > d2? d1 : d2;
223 if(inner) step *= -2;
224 else step *= 2;
225 } else {
226 G4double d1 = m_BoundingShape->DistanceToIn(A1);
227 G4double d2 = m_BoundingShape->DistanceToIn(A2);
228 step = d1 > d2? d1 : d2;
229 if(inner) step *= 2;
230 else step *= -2;
231 }
232
233 G4ThreeVector B1(0., r + step, z);
234 B1[0] = GetCalculator()->AmplitudeOfSurface(B1, -1, p_fan);
235 B1.rotateZ(dphi);
236 EInside B1i = m_BoundingShape->Inside(B1);
237 // EInside B1a = Inside_accordion(B1);
238 //std::cout << "B1: " << B1i << " " << B1a << std::endl;
239 if(B1i == kSurface){
240 //std::cout << "got B1" << std::endl;
241 p = B1;
242 return;
243 }
244 G4ThreeVector B2(0., r + step, z);
245 B2[0] = GetCalculator()->AmplitudeOfSurface(B2, 1, p_fan);
246 B2.rotateZ(dphi);
247 EInside B2i = m_BoundingShape->Inside(B2);
248 // EInside B2a = Inside_accordion(B2);
249 //std::cout << "B2: " << B2i << " " << B2a << std::endl;
250 if(B2i == kSurface){
251 //std::cout << "got B2" << std::endl;
252 p = B2;
253 return;
254 }
255
256 if(B1i == A1i || B2i == A1i){ // failed
257 set_failover_point(p, "pol fail1");
258 return;
259 }
260 if(A1i == kInside){
261 std::swap(A1, B1);
262 std::swap(A2, B2);
263 std::swap(A1i, B1i);
264 std::swap(A2i, B2i);
265 }
266 // here A* outside, B* inside, all on accordion surface
267
268 G4ThreeVector d1 = (A1 - B1).unit();
269 G4ThreeVector X1 = B1 + d1 * m_BoundingShape->DistanceToOut(B1, d1);
270 G4ThreeVector d2 = (A2 - B2).unit();
271 G4ThreeVector X2 = B2 + d2 * m_BoundingShape->DistanceToOut(B2, d2);
272
273 G4ThreeVector X = X1;
274 G4double phi1 = X1.phi(), phi2 = X2.phi();
275 // X1 corresponds to side = -, X2 to side = +
276 if(phi1 > 0. && phi2 < 0.) phi2 += CLHEP::twopi;
277 G4double phiX = rnd->Uniform(phi1, phi2);
278 if(phiX > CLHEP::pi) phiX -= CLHEP::twopi;
279 X.setPhi(phiX);
280
281 if(Inside(X) == kSurface){
282 p = X;
283 } else { // failed
284 set_failover_point(p, "pol fail2");
285 }
286}
287
289{
290 p[0] = 0.;
291 p[1] = 0.;
292 p[2] = rnd->Uniform() > 0.5? m_Zmin: m_Zmax;
293
294 G4double rmin, rmax;
295 get_r(m_BoundingShape, p[2], rmin, rmax);
296
297 p[1] = rnd->Uniform(rmin, rmax);
298 p.setPhi(rnd->Uniform(0., CLHEP::twopi));
299 G4double dphi = p.phi();
300 int p_fan = 0;
302 dphi -= p.phi();
303
304 p[0] = rnd->Uniform(
305 GetCalculator()->AmplitudeOfSurface(p, -1, p_fan),
306 GetCalculator()->AmplitudeOfSurface(p, 1, p_fan)
307 );
308
309 p.rotateZ(dphi);
310
311 if(m_BoundingShape->Inside(p) != kSurface){
312 set_failover_point(p, "flat fail");
313 }
314}
315
317{
318 // sagging ignored, effect should be negligible
319 double result =
320 m_f_vol->Integral(m_Rmin, m_Rmax, IntPrecision)
321
322#ifndef LOCAL_DEBUG
324#endif
325 ;
326 return result;
327}
328
330{
331 return m_f_area_on_pc->Integral(m_Zmin, m_Zmax);
332}
333
335{
336 G4double result = 0.;
337 G4double rmin, rmax;
338 get_r(m_BoundingShape, m_Zmin, rmin, rmax);
339 result += rmax - rmin;
340 get_r(m_BoundingShape, m_Zmax, rmin, rmax);
341 result += rmax - rmin;
343 return result;
344}
345
347{
348 m_f_length->SetParameter(0, r);
349
350 double zmin = m_BoundingShape->DistanceToIn(
351 G4ThreeVector(0., r, -m_Zmin), G4ThreeVector(0., 0., 1.)
352 );
353 zmin -= m_Zmin * 2.;
354 double zmax = m_BoundingShape->DistanceToIn(
355 G4ThreeVector(0., r, m_Zmax), G4ThreeVector(0., 0., -1.)
356 );
357 zmax = m_Zmax - zmax;
358 double result = m_f_length->Integral(zmin, zmax);
359 return result;
360}
361
363{
364 return m_f_side_area->Integral(m_Rmin, m_Rmax, IntPrecision);
365}
366
368{
369 double result = 0.;
370
371 double a1 = get_area_on_polycone();
372 result += a1;
373#ifdef LOCAL_DEBUG
374 std::cout << "get_area_on_polycone: " << a1/mm2 << std::endl;
375#endif
376
377 double a2 = get_area_on_face();
378 result += a2;
379#ifdef LOCAL_DEBUG
380 std::cout << "get_area_on_face: " << a2/mm2 << std::endl;
381#endif
382
383 double a3 = get_area_on_side();
384 result += a3;
385#ifdef LOCAL_DEBUG
386 std::cout << "get_area_on_side: " << a3/mm2 << std::endl;
387#endif
388
389 // sagging ignored, effect should be negligible
390 return result
391#ifndef LOCAL_DEBUG
393#endif
394 ;
395}
396
398{
399 boost::io::ios_all_saver ias(std::cout);
400 const char *on = getenv("LARWHEELSLICESOLID_TEST");
401 if(on == 0) return;
402 std::string test_mode = on;
403
404 std::cout << "============| LArWheelSliceSolid test() routine |=============="
405 << std::endl;
406 std::cout << "Solid of type " << TypeStr() << std::endl;
407 std::cout.precision(6);
408 std::cout << std::fixed;
409 const char *prec = getenv("LARWHEELSLICESOLID_TEST_INTPRECISION");
410 if(prec) IntPrecision = atof(prec);
411 std::cout << "Int. precision " << IntPrecision << std::endl;
412 std::cout << "test mode " << test_mode << std::endl;
413
414 std::cout << "m_Rmin = " << m_Rmin << " m_Rmax = " << m_Rmax << std::endl
415 << "m_Zmin = " << m_Zmin << " m_Zmax = " << m_Zmax << std::endl;
416
417 // TFile *F = new TFile("LArWheelSliceSolid_test.root", "RECREATE");
418 TFile *F = 0;
419 TNtupleD *T = 0;
420 if(test_mode.find("root") != std::string::npos){
421 F = new TFile("LArWheelSliceSolid_test.root", "UPDATE");
422 T = new TNtupleD(GetName(), GetName(), "x:y:z");
423 }
424 int N = 1000000;
425 const int Nmax(1000000000);
426 char *NN = getenv("LARWHEELSLICESOLID_TEST_NPOINTS");
427
428 if(NN) {
429 char *endptr;
430 N = strtol(NN, &endptr, 0);
431 if (endptr[0] != '\0') {
432 throw std::invalid_argument("Could not convert string to int: " + std::string(NN));
433 }
434 }
435 if (Nmax<N) {
436 std::cout << "Number of points from LARWHEELSLICESOLID_TEST_NPOINTS environment variable ("<<N<<") is too large. Using " << Nmax << " instead." << std::endl;
437 N=Nmax;
438 }
439 if (N<0) {
440 std::cout << "Number of points from LARWHEELSLICESOLID_TEST_NPOINTS environment variable ("<<N<<") is negative!!. Using 0 instead." << std::endl;
441 N=0;
442 }
443 if(test_mode.find("points") == std::string::npos){
444 N = 0;
445 } else {
446 std::cout << N << " points" << std::endl;
447 }
448 for(int i = 0; i < N; ++ i){
449 G4ThreeVector p = GetPointOnSurface();
450#ifdef LOCAL_DEBUG
451 EInside ii = Inside(p);
452 if(ii != kSurface){
453 std::cout << i << " "
454 << (ii == kInside? "inside": "outside")
455 << std::endl;
456 }
457#endif
458 if(T) T->Fill(p[0], p[1], p[2]);
459 }
460 if(F){
461 T->Write();
462 F->Write();
463 F->Close();
464 delete F;
465 }
466
467 if(test_mode.find("volume") != std::string::npos){
468 double cv = GetCubicVolume();
469 std::cout << "GetCubicVolume: " << cv/CLHEP::mm3 << " mm^3" << std::endl;
470 }
471
472 if(test_mode.find("area") != std::string::npos){
473 double sa = GetSurfaceArea();
474 std::cout << "GetSurfaceArea: " << sa/CLHEP::mm2 << " mm^2" << std::endl;
475 }
476
477 std::cout << "======= end of ArWheelSolid test() routine ============="
478 << std::endl;
479
480 if(test_mode.find("once") != std::string::npos) exit(0);
481 ias.restore();
482}
483
485 if(m_f_area) {
486 delete m_f_area;
487 m_f_area = 0;
488 }
489 if(m_f_vol) {
490 delete m_f_vol;
491 m_f_vol = 0;
492 }
493
494 if(m_f_area_on_pc) {
495 delete m_f_area_on_pc;
496 m_f_area_on_pc = 0;
497 }
498
499 if(m_f_length) {
500 delete m_f_length;
501 m_f_length = 0;
502 }
503 if(m_f_side_area) {
504 delete m_f_side_area;
505 m_f_side_area = 0;
506 }
507}
508
509G4double LArWheelSliceSolid::get_area_at_r(G4double r) const
510{
511 m_f_area->SetParameter(0, r);
512
513 double zmin = m_BoundingShape->DistanceToIn(
514 G4ThreeVector(0., r, m_Zmin), G4ThreeVector(0., 0., 1.)
515 );
516 double zmax = m_BoundingShape->DistanceToIn(
517 G4ThreeVector(0., r, m_Zmax), G4ThreeVector(0., 0., -1.)
518 );
519 zmax = m_Zmax - zmax;
520
521 double result = m_f_area->Integral(zmin, zmax);
522
523 return result;
524}
525
526static std::map<double, LArWheelSliceSolid *> solid;
527
528double LArWheelSliceSolid_fcn_area(double *x, double *p)
529{
530 const double &z = x[0];
531 const double &r = p[0];
532 const double &index = p[1];
533
534 G4ThreeVector a(0., r, z);
535 double b = solid[index]->GetCalculator()->AmplitudeOfSurface(a, -1, 121) // sagging ignored, effect should be negligible, use arbitrary fan number
536 - solid[index]->GetCalculator()->AmplitudeOfSurface(a, 1, 121);
537 return b;
538}
539
540double LArWheelSliceSolid_fcn_vol(double *x, double *p)
541{
542 const double &r = x[0];
543 const double &index = p[0];
544
545 return solid[index]->get_area_at_r(r);
546}
547
548double LArWheelSliceSolid_fcn_area_on_pc(double *x, double *p)
549{
550 const double &z = x[0];
551 const double &index = p[0];
552
553 G4double rmin, rmax;
554 get_r(solid[index]->m_BoundingShape, z, rmin, rmax);
555
556 double result = 0.;
557 G4ThreeVector a(0., rmin, z);
558 result += solid[index]->GetCalculator()->AmplitudeOfSurface(a, -1, 232) // sagging ignored, effect should be negligible, use arbitrary fan number
559 - solid[index]->GetCalculator()->AmplitudeOfSurface(a, 1, 232);
560 a[1] = rmax;
561 result += solid[index]->GetCalculator()->AmplitudeOfSurface(a, -1, 343)
562 - solid[index]->GetCalculator()->AmplitudeOfSurface(a, 1, 343);
563
564 return result;
565}
566
567
568double LArWheelSliceSolid_get_dl(double *x, double *par, G4int side)
569{
570 const double &z = x[0];
571 const double &r = par[0];
572 const double &index = par[1];
573
574 const double h = 0.001;
575
576 //check what happens if z+h > m_Zmax etc
577 G4ThreeVector p(0., r, z + h);
578 G4double D1 = solid[index]->GetCalculator()->AmplitudeOfSurface(p, side, 5665); // sagging ignored, effect should be negligible, use arbitrary fan number
579 p[2] = z - h;
580 D1 -= solid[index]->GetCalculator()->AmplitudeOfSurface(p, side, 5665);
581 D1 /= 2 * h;
582
583 p[2] = z + h / 2;
584 G4double D2 = solid[index]->GetCalculator()->AmplitudeOfSurface(p, side, 5665);
585 p[2] = z - h / 2;
586 D2 -= solid[index]->GetCalculator()->AmplitudeOfSurface(p, side, 5665);
587 D2 /= h;
588
589 G4double D = (D2 * 4 - D1) / 3.;
590 G4double dl = sqrt(1 + D * D);
591
592 return dl;
593}
594
595static double fcn_length(double *x, double *p)
596{
598}
599
600double LArWheelSliceSolid_fcn_side_area(double *x, double *p)
601{
602 const double &r = x[0];
603 const double &index = p[0];
604
605 return solid[index]->get_length_at_r(r);
606}
607
609{
610 m_test_index = double(solid.size());
611 solid[m_test_index] = this;
612
613#ifdef DEBUG_LARWHEELSLICESOLID
614 if(Verbose > 0)
615 std::cout << "LArWheelSliceSolid::init_tests: put " << this
616 << " with index " << m_test_index << std::endl;
617#endif
618
619 m_f_area = new TF1(
620 (GetName() + "_f_area").c_str(), &LArWheelSliceSolid_fcn_area,
621 m_Zmin, m_Zmax, 2
622 );
623 m_f_area->FixParameter(1, m_test_index);
624
625 m_f_vol = new TF1(
626 (GetName() + "_f_vol").c_str(), &LArWheelSliceSolid_fcn_vol,
627 m_Rmin, m_Rmax, 1
628 );
629 m_f_vol->FixParameter(0, m_test_index);
630
631 m_f_area_on_pc = new TF1(
632 (GetName() + "_f_area_pc").c_str(), &LArWheelSliceSolid_fcn_area_on_pc,
633 m_Zmin, m_Zmax, 1
634 );
635 m_f_area_on_pc->FixParameter(0, m_test_index);
636
637 m_f_length = new TF1(
638 (GetName() + "_f_length").c_str(), &fcn_length,
639 m_Zmin, m_Zmax, 2
640 );
641 m_f_length->FixParameter(1, m_test_index);
642
643 m_f_side_area = new TF1(
644 (GetName() + "_f_side_area").c_str(), &LArWheelSliceSolid_fcn_side_area,
645 m_Rmin, m_Rmax, 1
646 );
647 m_f_side_area->FixParameter(0, m_test_index);
648}
649
650#ifdef DEBUG_LARWHEELSLICESOLID
651G4bool LArWheelSliceSolid::test_dti(
652 const G4ThreeVector &inputP, const G4ThreeVector &inputV,
653 const G4double distance
654) const
655{
656 if(distance > 10.*CLHEP::m){
657 LWSDBG(1, std::cout << "DTI test skipped, distance > 10 m"
658 << std::endl);
659 return false;
660 }
661 unsigned long counter = 0;
662 counter ++;
663 G4ThreeVector p;
664 if(m_BoundingShape->Inside(inputP) == kOutside){
665 p = inputP + inputV * m_BoundingShape->DistanceToIn(inputP, inputV);
666 } else p = inputP;
667 const G4double phi0 = p.phi();
668 int p_fan = 0;
669 const G4double d = GetCalculator()->DistanceToTheNearestFan(p, p_fan);
670 if(fabs(d) < m_FHTplusT){
671 std::cout << "DTI test already inside" << MSG_VECTOR(p) << std::endl;
672 return false;
673 }
674 G4ThreeVector v( inputV );
675 v.rotateZ(p.phi() - phi0);
676 const G4double dd = s_IterationPrecision;
677 LWSDBG(1, std::cout << "Start DTI test, expect "
678 << long(distance / dd) << " iterations"
679 << std::endl);
680
681 G4int V = Verbose;
682 Verbose = 0;
683
684 const G4double d1 = distance - s_IterationPrecision;
685 bool first = true;
686 for(G4double t = s_IterationPrecision; t < d1; t += dd){
687 G4ThreeVector p1 = p + v * t;
688 if(fabs(GetCalculator()->DistanceToTheNeutralFibre(p1, p_fan)) < m_FHTplusT){
689 std::cout << "DTI test at " << MSG_VECTOR(inputP) << " -> "
690 << MSG_VECTOR(inputV) << ", translated to "
691 << MSG_VECTOR(p) << " - > " << MSG_VECTOR(v)
692 << " in range of "
693 << distance << ": found nearer intersection at local point"
694 << MSG_VECTOR(p1) << ", distance " << t
695 << ", call " << counter
696 << std::endl;
697 Verbose = V;
698
699 if(first){
700 first = false;
701 FILE *F = fopen("dti_error.dat", "w");
702 if(F){
703 fprintf(F, "%10e %10e %10e\n", p.x(), p.y(), p.z());
704 fprintf(F, "%10e %10e %10e\n", v.x(), v.y(), v.z());
705 for(G4double e = s_IterationPrecision; e < d1; e += dd){
706 p1 = p + v * e;
707 G4double f = fabs(GetCalculator()->DistanceToTheNeutralFibre(p1, p_fan)) - m_FanHalfThickness;
708 fprintf(F, "%10e %10e\n", e, f);
709 }
710 fclose(F);
711 }
712 }
713
714 return true;
715 }
716 }
717 Verbose = V;
718 LWSDBG(1, std::cout << "DTI test at " << MSG_VECTOR(p) << " -> "
719 << MSG_VECTOR(v) << " in range of "
720 << distance << ": allright" << std::endl);
721 return false;
722}
723
724G4bool LArWheelSliceSolid::test_dto(
725 const G4ThreeVector &inputP, const G4ThreeVector &inputV,
726 const G4double distance
727) const
728{
729 if(distance > 10.*CLHEP::m){
730 LWSDBG(1, std::cout << "DTO test skipped, distance > 10 m"
731 << std::endl);
732 return false;
733 }
734 unsigned long counter = 0;
735 counter ++;
736 G4ThreeVector p( inputP );
737 const G4double phi0 = p.phi();
738 int p_fan = 0;
739 const G4double d = GetCalculator()->DistanceToTheNearestFan(p, p_fan);
740 if(fabs(d) > m_FHTplusT){
741 std::cout << "DTO test already outside" << MSG_VECTOR(p) << std::endl;
742 return false;
743 }
744 G4ThreeVector v( inputV );
745 v.rotateZ(p.phi() - phi0);
746 const G4double dd = s_IterationPrecision;
747 LWSDBG(1, std::cout << "Start DTO test, expect "
748 << long(distance / dd) << " iterations"
749 << std::endl);
750
751 G4int V = Verbose;
752 Verbose = 0;
753
754 const G4double d1 = distance - s_IterationPrecision;
755 static bool first = true;
756 for(G4double t = s_IterationPrecision; t < d1; t += dd){
757 G4ThreeVector p1 = p + v * t;
758 if(fabs(GetCalculator()->DistanceToTheNeutralFibre(p1, p_fan)) > m_FHTplusT){
759 std::cout << "DTO test at " << MSG_VECTOR(inputP) << " -> "
760 << MSG_VECTOR(inputV) << ", translated to "
761 << MSG_VECTOR(p) << " - > " << MSG_VECTOR(v)
762 << " in range of "
763 << distance << ": found nearer intersection at local point"
764 << MSG_VECTOR(p1) << ", distance " << t
765 << ", call " << counter
766 << std::endl;
767 Verbose = V;
768
769 if(first){
770 first = false;
771 FILE *F = fopen("dto_error.dat", "w");
772 if(F){
773 fprintf(F, "%10e %10e %10e\n", p.x(), p.y(), p.z());
774 fprintf(F, "%10e %10e %10e\n", v.x(), v.y(), v.z());
775 for(G4double e = s_IterationPrecision; e < d1; e += dd){
776 p1 = p + v * e;
777 G4double f = fabs(GetCalculator()->DistanceToTheNeutralFibre(p1, p_fan)) - m_FanHalfThickness;
778 fprintf(F, "%10e %10e\n", e, f);
779 }
780 fclose(F);
781 }
782 }
783
784 return true;
785 }
786 }
787 Verbose = V;
788 LWSDBG(1, std::cout << "DTO test at " << MSG_VECTOR(p) << " -> "
789 << MSG_VECTOR(v) << " in range of "
790 << distance << ": allright" << std::endl);
791 return false;
792}
793#endif
const PlainObject unit() const
This is a plugin that makes Eigen look like CLHEP & defines some convenience methods.
static Double_t a
double LArWheelSliceSolid_fcn_vol(double *x, double *p)
double LArWheelSliceSolid_fcn_area_on_pc(double *x, double *p)
double LArWheelSliceSolid_fcn_side_area(double *x, double *p)
double LArWheelSliceSolid_fcn_area(double *x, double *p)
static std::map< double, LArWheelSliceSolid * > solid
static double fcn_length(double *x, double *p)
static void get_r(const G4VSolid *p, G4double z, G4double &rmin, G4double &rmax)
static TRandom * rnd
static double IntPrecision
double LArWheelSliceSolid_get_dl(double *x, double *par, G4int side)
#define LWSDBG(a, b)
#define F(x, y, z)
Definition MD5.cxx:112
#define x
#define z
Define macros for attributes used to control the static checker.
#define ATLAS_NO_CHECK_FILE_THREAD_SAFETY
Header file for AthHistogramAlgorithm.
double AmplitudeOfSurface(const CLHEP::Hep3Vector &P, int side, int fan_number) const
double GetFanHalfThickness(LArG4::LArWheelCalculator_t) const
double DistanceToTheNearestFan(CLHEP::Hep3Vector &p, int &out_fan_number) const
Determines the nearest to the input point fan.
friend double LArWheelSliceSolid_fcn_vol(double *, double *)
static const G4double s_IterationPrecision
friend double LArWheelSliceSolid_fcn_area_on_pc(double *, double *)
const LArWheelCalculator * m_Calculator
friend double LArWheelSliceSolid_fcn_side_area(double *, double *)
G4double get_area_on_side(void) const
G4double get_area_at_r(G4double r) const
void get_point_on_polycone_surface(G4ThreeVector &) const
const LArWheelCalculator * GetCalculator(void) const
friend double LArWheelSliceSolid_fcn_area(double *, double *)
static const G4double s_Tolerance
void set_failover_point(G4ThreeVector &p, const char *m=0) const
EInside Inside_accordion(const G4ThreeVector &) const
G4double get_length_at_r(G4double r) const
void get_point_on_flat_surface(G4ThreeVector &) const
EInside Inside(const G4ThreeVector &) const
G4double get_area_on_polycone(void) const
G4double get_area_on_face(void) const
G4String TypeStr(void) const
G4ThreeVector GetPointOnSurface(void) const
void get_point_on_accordion_surface(G4ThreeVector &) const
int r
Definition globals.cxx:22
float distance(const Amg::Vector3D &p1, const Amg::Vector3D &p2)
calculates the distance between two point in 3D space
bool first
Definition DeMoScan.py:534
@ Verbose
Definition ZDCMsg.h:18
Definition index.py:1
void swap(ElementLinkVector< DOBJ > &lhs, ElementLinkVector< DOBJ > &rhs)