Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
Loading...
Searching...
No Matches
ComponentAnalyticField.hh
Go to the documentation of this file.
1#ifndef G_COMPONENT_ANALYTIC_FIELD_H
2#define G_COMPONENT_ANALYTIC_FIELD_H
3
4#include <array>
5#include <complex>
6#include <mutex>
7#include <string>
8#include <vector>
9
10#include "Garfield/Component.hh"
11
12class TPad;
13
14namespace Garfield {
15
16class Medium;
17
20
21class ComponentAnalyticField : public Component {
22 public:
27
29 void SetMedium(Medium* medium) { m_medium = medium; }
31 void AddWire(const double x, const double y, const double diameter,
32 const double voltage, const std::string& label = "",
33 const double length = 100., const double tension = 50.,
34 const double rho = 19.3, const int ntrap = 5);
36 void AddTube(const double radius, const double voltage, const int nEdges,
37 const std::string& label = "");
39 void AddPlaneX(const double x, const double voltage,
40 const std::string& label = "");
42 void AddPlaneY(const double y, const double voltage,
43 const std::string& label = "");
45 void AddPlaneR(const double r, const double voltage,
46 const std::string& label = "");
48 void AddPlanePhi(const double phi, const double voltage,
49 const std::string& label = "");
58 void AddStripOnPlaneX(const char direction, const double x, const double smin,
59 const double smax, const std::string& label,
60 const double gap = -1.);
62 void AddStripOnPlaneY(const char direction, const double y, const double smin,
63 const double smax, const std::string& label,
64 const double gap = -1.);
67 void AddStripOnPlaneR(const char direction, const double r, const double smin,
68 const double smax, const std::string& label,
69 const double gap = -1.);
71 void AddStripOnPlanePhi(const char direction, const double phi,
72 const double smin, const double smax,
73 const std::string& label, const double gap = -1.);
84 void AddPixelOnPlaneX(const double x, const double ymin, const double ymax,
85 const double zmin, const double zmax,
86 const std::string& label, const double gap = -1.,
87 const double rot = 0.);
89 void AddPixelOnPlaneY(const double y, const double xmin, const double xmax,
90 const double zmin, const double zmax,
91 const std::string& label, const double gap = -1.,
92 const double rot = 0.);
94 void AddPixelOnPlaneR(const double r, const double phimin,
95 const double phimax, const double zmin,
96 const double zmax, const std::string& label,
97 const double gap = -1.);
99 void AddPixelOnPlanePhi(const double phi, const double rmin,
100 const double rmax, const double zmin,
101 const double zmax, const std::string& label,
102 const double gap = -1.);
103
105 void SetPeriodicityX(const double s);
107 void SetPeriodicityY(const double s);
109 void SetPeriodicityPhi(const double phi);
111 bool GetPeriodicityX(double& s);
113 bool GetPeriodicityY(double& s);
115 bool GetPeriodicityPhi(double& s);
116
122 bool IsPolar() const { return m_polar; }
123
125 void PrintCell();
127 void PlotCell(TPad* pad);
128
130 void AddCharge(const double x, const double y, const double z,
131 const double q);
135 void PrintCharges() const;
136
156 std::string GetCellType() {
157 if (!m_cellset) {
158 if (CellCheck()) CellType();
159 }
160 return GetCellType(m_cellType);
161 }
162
165 void AddReadout(const std::string& label, const bool silent = false);
166
167 void SetNumberOfCellCopies(const unsigned int nfourier);
168
179 bool MultipoleMoments(const unsigned int iw, const unsigned int order = 4,
180 const bool print = false, const bool plot = false,
181 const double rmult = 1., const double eps = 1.e-4,
182 const unsigned int nMaxIter = 20);
184 void EnableDipoleTerms(const bool on = true);
186 void EnableChargeCheck(const bool on = true) { m_chargeCheck = on; }
187
189 unsigned int GetNumberOfWires() const { return m_nWires; }
191 bool GetWire(const unsigned int i, double& x, double& y, double& diameter,
192 double& voltage, std::string& label, double& length,
193 double& charge, int& ntrap) const;
194
196 unsigned int GetNumberOfPlanesX() const;
198 unsigned int GetNumberOfPlanesY() const;
200 unsigned int GetNumberOfPlanesR() const;
202 unsigned int GetNumberOfPlanesPhi() const;
204 bool GetPlaneX(const unsigned int i, double& x, double& voltage,
205 std::string& label) const;
207 bool GetPlaneY(const unsigned int i, double& y, double& voltage,
208 std::string& label) const;
210 bool GetPlaneR(const unsigned int i, double& r, double& voltage,
211 std::string& label) const;
213 bool GetPlanePhi(const unsigned int i, double& phi, double& voltage,
214 std::string& label) const;
216 bool GetTube(double& r, double& voltage, int& nEdges,
217 std::string& label) const;
218
232 bool OptimiseOnTrack(const std::vector<std::string>& groups,
233 const std::string& field_function, const double target,
234 const double x0, const double y0, const double x1,
235 const double y1, const unsigned int nP = 20,
236 const bool print = true);
239 bool OptimiseOnGrid(const std::vector<std::string>& groups,
240 const std::string& field_function, const double target,
241 const double x0, const double y0, const double x1,
242 const double y1, const unsigned int nX = 10,
243 const unsigned int nY = 10, const bool print = true);
246 bool OptimiseOnWires(const std::vector<std::string>& groups,
247 const std::string& field_function, const double target,
248 const std::vector<unsigned int>& wires,
249 const bool print = true);
250
257 void SetOptimisationParameters(const double dist = 1.,
258 const double eps = 1.e-4,
259 const unsigned int nMaxIter = 10);
260
263 bool ElectricFieldAtWire(const unsigned int iw, double& ex, double& ey);
264
267 void SetScanningGrid(const unsigned int nX, const unsigned int nY);
270 void SetScanningArea(const double xmin, const double xmax, const double ymin,
271 const double ymax);
277 void SetScanningAreaFirstOrder(const double scale = 2.);
280 void EnableExtrapolation(const bool on = true) { m_extrapolateForces = on; }
281
283 void SetGravity(const double dx, const double dy, const double dz);
285 void GetGravity(double& dx, double& dy, double& dz) const;
287 void EnableGravity(const bool on = true) { m_useGravitationalForce = on; }
289 void EnableElectrostaticForce(const bool on = true) {
291 }
292
299 bool ForcesOnWire(const unsigned int iw, std::vector<double>& xMap,
300 std::vector<double>& yMap,
301 std::vector<std::vector<double> >& fxMap,
302 std::vector<std::vector<double> >& fyMap);
313 bool WireDisplacement(const unsigned int iw, const bool detailed,
314 std::vector<double>& csag, std::vector<double>& xsag,
315 std::vector<double>& ysag, double& stretch,
316 const bool print = true);
319 void SetNumberOfShots(const unsigned int n) { m_nShots = n; }
321 void SetNumberOfSteps(const unsigned int n);
322
323 Medium* GetMedium(const double x, const double y, const double z) override;
324 void ElectricField(const double x, const double y, const double z, double& ex,
325 double& ey, double& ez, Medium*& m, int& status) override;
326
327 void ElectricField(const double x, const double y, const double z, double& ex,
328 double& ey, double& ez, double& v, Medium*& m,
329 int& status) override;
330
331 using Component::ElectricField;
332 bool GetVoltageRange(double& pmin, double& pmax) override;
333
334 void WeightingField(const double x, const double y, const double z,
335 double& wx, double& wy, double& wz,
336 const std::string& label) override {
337 wx = wy = wz = 0.;
338 if (!m_sigset) PrepareSignals();
339 Wfield(x, y, z, wx, wy, wz, label);
340 }
341 double WeightingPotential(const double x, const double y, const double z,
342 const std::string& label) override {
343 if (!m_sigset) PrepareSignals();
344 return Wpot(x, y, z, label);
345 }
346
347 bool GetBoundingBox(double& x0, double& y0, double& z0, double& x1,
348 double& y1, double& z1) override;
349 bool GetElementaryCell(double& x0, double& y0, double& z0, double& x1,
350 double& y1, double& z1) override;
351
352 bool CrossedWire(const double x0, const double y0, const double z0,
353 const double x1, const double y1, const double z1,
354 double& xc, double& yc, double& zc, const bool centre,
355 double& rc) override;
356
357 bool InTrapRadius(const double q0, const double x0, const double y0,
358 const double z0, double& xw, double& yx,
359 double& rw) override;
360
361 bool CrossedPlane(const double x0, const double y0, const double z0,
362 const double x1, const double y1, const double z1,
363 double& xc, double& yc, double& zc) override;
364
381
382 double StepSizeHint() override;
383
384 private:
385 std::mutex m_mutex;
386
387 Medium* m_medium = nullptr;
388
389 bool m_chargeCheck = false;
390
391 bool m_cellset = false;
392 bool m_sigset = false;
393
394 bool m_polar = false;
395
396 // Cell type.
398
399 // Bounding box
400 double m_xmin, m_xmax;
401 double m_ymin, m_ymax;
402 double m_zmin, m_zmax;
403
404 // Voltage range
405 double m_vmin, m_vmax;
406
407 // Periodicities
408 bool m_perx = false;
409 bool m_pery = false;
410 double m_sx, m_sy;
411
412 // Signals
413 int m_nFourier = 1;
415 bool m_fperx = false;
416 bool m_fpery = false;
417 int m_mxmin = 0;
418 int m_mxmax = 0;
419 int m_mymin = 0;
420 int m_mymax = 0;
421 int m_mfexp = 0;
422
423 std::vector<std::string> m_readout;
424
425 // Wires
426 unsigned int m_nWires;
427 struct Wire {
428 double x, y;
429 double r;
430 double v;
431 double e;
432 std::string type;
433 double u;
434 int ind;
436 int nTrap;
437 double tension;
438 double density;
439 };
440 std::vector<Wire> m_w;
441
442 // Option for computation of dipole terms
443 bool m_dipole = false;
444 // Dipole angle and amplitude
445 std::vector<double> m_cosph2;
446 std::vector<double> m_sinph2;
447 std::vector<double> m_amp2;
448
449 // Parameters for B2 type cells
450 std::vector<double> m_b2sin;
451 // Parameters for C type cells
453 std::complex<double> m_zmult;
454 double m_p1, m_p2, m_c1;
455 // Parameters for D3 type cells
456 // Wire positions in conformal mapping
457 std::vector<std::complex<double> > m_zw;
458 double m_kappa;
459
460 // Reference potential
461 double m_v0;
463
464 // Planes
465 // Existence
466 bool m_ynplan[4];
468 // Coordinates
469 double m_coplan[4];
471 // Voltages
472 double m_vtplan[4];
473
474 struct Strip {
475 std::string type;
476 int ind;
477 double smin, smax;
478 double gap;
479 };
480
481 struct Pixel {
482 std::string type;
483 int ind = 0;
484 double smin = 0., smax = 0.;
485 double zmin = 0., zmax = 0.;
486 double gap = -1.;
487 double cphi = 1.;
488 double sphi = 0.;
489 };
490
491 struct Plane {
492 std::string type;
493 int ind;
494 double ewxcor, ewycor;
495 std::vector<Strip> strips1;
496 std::vector<Strip> strips2;
497 std::vector<Pixel> pixels;
498 };
499 std::array<Plane, 5> m_planes;
500
501 // Tube
502 bool m_tube = false;
503 int m_mtube = 0;
504 int m_ntube = 1;
505 double m_cotube = 1.;
506 double m_cotube2 = 1.;
507 double m_vttube = 0.;
508
509 // Smallest dimension.
510 double m_dmin = -1.;
511
512 // Wire weighting charges.
513 std::vector<std::vector<std::vector<double> > > m_qwire;
514 // Plane weighting charges.
515 std::vector<std::vector<std::vector<double> > > m_qplane;
516
517 // Point charges
518 struct Charge3d {
519 double x, y, z;
520 double e;
521 };
522 std::vector<Charge3d> m_ch3d;
523 unsigned int m_nTermBessel = 10;
524 unsigned int m_nTermPoly = 100;
525
528 // Gravity
529 std::array<double, 3> m_down{{0, 0, 1}};
530 // Number of shots used for solving the wire sag differential equations
531 unsigned int m_nShots = 2;
532 // Number of integration steps within each shot.
533 unsigned int m_nSteps = 20;
534 // Options for setting the range of wire shifts
535 // for which the forces are computed.
536 enum class ScanningRange { Largest = 0, FirstOrder, User };
538 // Limits of the user-specified scanning range.
539 double m_xScanMin = 0.;
540 double m_xScanMax = 0.;
541 double m_yScanMin = 0.;
542 double m_yScanMax = 0.;
543 // Scaling factor for first-order estimate of the scanning range.
544 double m_scaleRange = 2.;
545 // Number of grid lines at which the forces are stored.
546 unsigned int m_nScanX = 11;
547 unsigned int m_nScanY = 11;
548 // Extrapolate beyond the scanning range or not.
550
551 // Maximum deviation between target and field function at which
552 // to allow the iteration to stop.
553 double m_optDist = 1.;
554 // Relative change in Euclidean distance between target and
555 // field function at which to allow the iteration to stop.
556 double m_optEps = 1.e-4;
557 // Maximum number of iterations in the optimisation fit.
558 unsigned int m_optNitmax = 10;
559
560 void UpdatePeriodicity() override;
561 void Reset() override {
562 CellInit();
563 m_medium = nullptr;
564 }
565
566 void CellInit();
567 bool Prepare();
568 bool CellCheck();
569 bool WireCheck() const;
570 bool CellType();
571 std::string GetCellType(const Cell) const;
576
577 // Calculation of charges
578 bool Setup();
579 bool Update(const std::vector<double>& vw, const std::array<double, 5>& vp);
580 bool SetupA00();
581 bool SetupB1X();
582 bool SetupB1Y();
583 bool SetupB2X();
584 bool SetupB2Y();
585 bool SetupC10();
586 bool SetupC2X();
587 bool SetupC2Y();
588 bool SetupC30();
589 bool SetupD10();
590 bool SetupD20();
591 bool SetupD30();
592
593 bool IprA00(const int mx, const int my,
594 std::vector<std::vector<std::complex<double> > >& mat);
595 bool IprB2X(const int my,
596 std::vector<std::vector<std::complex<double> > >& mat);
597 bool IprB2Y(const int mx,
598 std::vector<std::vector<std::complex<double> > >& mat);
599 bool IprC2X(std::vector<std::vector<std::complex<double> > >& mat);
600 bool IprC2Y(std::vector<std::vector<std::complex<double> > >& mat);
601 bool IprC30(std::vector<std::vector<std::complex<double> > >& mat);
602 bool IprD10(std::vector<std::vector<std::complex<double> > >& mat);
603 bool IprD30(std::vector<std::vector<std::complex<double> > >& mat);
604
606
607 // Inversion of the capacitance matrix
608 bool Charge(std::vector<std::vector<double> >& mat);
609
610 // Evaluation of the electric field
611 int Field(const double xin, const double yin, const double zin, double& ex,
612 double& ey, double& ez, double& volt, const bool opt);
613 void FieldA00(const double xpos, const double ypos, double& ex, double& ey,
614 double& volt, const bool opt) const;
615 void FieldB1X(const double xpos, const double ypos, double& ex, double& ey,
616 double& volt, const bool opt) const;
617 void FieldB1Y(const double xpos, const double ypos, double& ex, double& ey,
618 double& volt, const bool opt) const;
619 void FieldB2X(const double xpos, const double ypos, double& ex, double& ey,
620 double& volt, const bool opt) const;
621 void FieldB2Y(const double xpos, const double ypos, double& ex, double& ey,
622 double& volt, const bool opt) const;
623 void FieldC10(const double xpos, const double ypos, double& ex, double& ey,
624 double& volt, const bool opt) const;
625 void FieldC2X(const double xpos, const double ypos, double& ex, double& ey,
626 double& volt, const bool opt) const;
627 void FieldC2Y(const double xpos, const double ypos, double& ex, double& ey,
628 double& volt, const bool opt) const;
629 void FieldC30(const double xpos, const double ypos, double& ex, double& ey,
630 double& volt, const bool opt) const;
631 void FieldD10(const double xpos, const double ypos, double& ex, double& ey,
632 double& volt, const bool opt) const;
633 void FieldD20(const double xpos, const double ypos, double& ex, double& ey,
634 double& volt, const bool opt) const;
635 void FieldD30(const double xpos, const double ypos, double& ex, double& ey,
636 double& volt, const bool opt) const;
637
638 // Field due to point charges
639 void Field3dA00(const double x, const double y, const double z, double& ex,
640 double& ey, double& ez, double& volt) const;
641 void Field3dB2X(const double x, const double y, const double z, double& ex,
642 double& ey, double& ez, double& volt) const;
643 void Field3dB2Y(const double x, const double y, const double z, double& ex,
644 double& ey, double& ez, double& volt) const;
645 void Field3dD10(const double x, const double y, const double z, double& ex,
646 double& ey, double& ez, double& volt) const;
647 // Evaluation of the weighting field
648 bool Wfield(const double x, const double y, const double z, double& ex,
649 double& ey, double& ez, const std::string& label) const;
650 void WfieldWireA00(const double x, const double y, double& ex, double& ey,
651 const int mx, const int my,
652 const std::vector<double>& qw) const;
653 void WfieldWireB2X(const double x, const double y, double& ex, double& ey,
654 const int my, const std::vector<double>& qw) const;
655 void WfieldWireB2Y(const double x, const double y, double& ex, double& ey,
656 const int mx, const std::vector<double>& qw) const;
657 void WfieldWireC2X(const double x, const double y, double& ex, double& ey,
658 const std::vector<double>& qw) const;
659 void WfieldWireC2Y(const double x, const double y, double& ex, double& ey,
660 const std::vector<double>& qw) const;
661 void WfieldWireC30(const double x, const double y, double& ex, double& ey,
662 const std::vector<double>& qw) const;
663 void WfieldWireD10(const double x, const double y, double& ex, double& ey,
664 const std::vector<double>& qw) const;
665 void WfieldWireD30(const double x, const double y, double& ex, double& ey,
666 const std::vector<double>& qw) const;
667 void WfieldPlaneA00(const double x, const double y, double& ex, double& ey,
668 const int mx, const int my,
669 const std::vector<double>& qp) const;
670 void WfieldPlaneB2X(const double x, const double y, double& ex, double& ey,
671 const int my, const std::vector<double>& qp) const;
672 void WfieldPlaneB2Y(const double x, const double ypos, double& ex, double& ey,
673 const int mx, const std::vector<double>& qp) const;
674 void WfieldPlaneC2X(const double x, const double y, double& ex, double& ey,
675 const std::vector<double>& qp) const;
676 void WfieldPlaneC2Y(const double x, const double y, double& ex, double& ey,
677 const std::vector<double>& qp) const;
678 void WfieldPlaneC30(const double x, const double y, double& ex, double& ey,
679 const std::vector<double>& qp) const;
680 void WfieldPlaneD10(const double x, const double y, double& ex, double& ey,
681 const std::vector<double>& qp) const;
682 void WfieldPlaneD30(const double x, const double y, double& ex, double& ey,
683 const std::vector<double>& qp) const;
684 void WfieldStripZ(const double x, const double y, double& ex, double& ey,
685 const int ip, const Strip& strip) const;
686 void WfieldStripXy(const double x, const double y, const double z, double& ex,
687 double& ey, double& ez, const int ip,
688 const Strip& strip) const;
689 void WfieldStrip(const double x, const double y, const double g,
690 const double w, double& fx, double& fy) const;
691 void WfieldPixel(const double x, const double y, const double z, double& ex,
692 double& ey, double& ez, const int ip,
693 const Pixel& pixel) const;
694
695 // Evaluation of the weighting potential.
696 double Wpot(const double x, const double y, const double z,
697 const std::string& label) const;
698 double WpotWireA00(const double x, const double y, const int mx, const int my,
699 const std::vector<double>& qw) const;
700 double WpotWireB2X(const double x, const double y, const int my,
701 const std::vector<double>& qw) const;
702 double WpotWireB2Y(const double x, const double y, const int mx,
703 const std::vector<double>& qw) const;
704 double WpotWireC2X(const double x, const double y,
705 const std::vector<double>& qw) const;
706 double WpotWireC2Y(const double x, const double y,
707 const std::vector<double>& qw) const;
708 double WpotWireC30(const double x, const double y,
709 const std::vector<double>& qw) const;
710 double WpotWireD10(const double x, const double y,
711 const std::vector<double>& qw) const;
712 double WpotWireD30(const double x, const double y,
713 const std::vector<double>& qw) const;
714 double WpotPlaneA00(const double x, const double y, const int mx,
715 const int my, const std::vector<double>& qp) const;
716 double WpotPlaneB2X(const double x, const double y, const int my,
717 const std::vector<double>& qp) const;
718 double WpotPlaneB2Y(const double x, const double y, const int mx,
719 const std::vector<double>& qp) const;
720 double WpotPlaneC2X(const double x, const double y,
721 const std::vector<double>& qp) const;
722 double WpotPlaneC2Y(const double x, const double y,
723 const std::vector<double>& qp) const;
724 double WpotPlaneC30(const double x, const double y,
725 const std::vector<double>& qp) const;
726 double WpotPlaneD10(const double x, const double y,
727 const std::vector<double>& qp) const;
728 double WpotPlaneD30(const double x, const double y,
729 const std::vector<double>& qp) const;
730 double WpotStripZ(const double x, const double y, const int ip,
731 const Strip& strip) const;
732 double WpotStripXy(const double x, const double y, const double z,
733 const int ip, const Strip& strip) const;
734 double WpotPixel(const double x, const double y, const double z, const int ip,
735 const Pixel& pixel) const;
736
737 // Functions for calculating the electric field at a given wire position,
738 // as if the wire itself were not there but with the presence
739 // of its mirror images.
740 void FieldAtWireA00(const double xpos, const double ypos, double& ex,
741 double& ey, const std::vector<bool>& cnalso) const;
742 void FieldAtWireB1X(const double xpos, const double ypos, double& ex,
743 double& ey, const std::vector<bool>& cnalso) const;
744 void FieldAtWireB1Y(const double xpos, const double ypos, double& ex,
745 double& ey, const std::vector<bool>& cnalso) const;
746 void FieldAtWireB2X(const double xpos, const double ypos, double& ex,
747 double& ey, const std::vector<bool>& cnalso) const;
748 void FieldAtWireB2Y(const double xpos, const double ypos, double& ex,
749 double& ey, const std::vector<bool>& cnalso) const;
750 void FieldAtWireC10(const double xpos, const double ypos, double& ex,
751 double& ey, const std::vector<bool>& cnalso) const;
752 void FieldAtWireC2X(const double xpos, const double ypos, double& ex,
753 double& ey, const std::vector<bool>& cnalso) const;
754 void FieldAtWireC2Y(const double xpos, const double ypos, double& ex,
755 double& ey, const std::vector<bool>& cnalso) const;
756 void FieldAtWireC30(const double xpos, const double ypos, double& ex,
757 double& ey, const std::vector<bool>& cnalso) const;
758 void FieldAtWireD10(const double xpos, const double ypos, double& ex,
759 double& ey, const std::vector<bool>& cnalso) const;
760 void FieldAtWireD20(const double xpos, const double ypos, double& ex,
761 double& ey, const std::vector<bool>& cnalso) const;
762 void FieldAtWireD30(const double xpos, const double ypos, double& ex,
763 double& ey, const std::vector<bool>& cnalso) const;
764
765 void DipoleFieldA00(const double xpos, const double ypos, double& ex,
766 double& ey, double& volt, const bool opt) const;
767 void DipoleFieldB1X(const double xpos, const double ypos, double& ex,
768 double& ey, double& volt, const bool opt) const;
769 void DipoleFieldB1Y(const double xpos, const double ypos, double& ex,
770 double& ey, double& volt, const bool opt) const;
771 void DipoleFieldB2X(const double xpos, const double ypos, double& ex,
772 double& ey, double& volt, const bool opt) const;
773 void DipoleFieldB2Y(const double xpos, const double ypos, double& ex,
774 double& ey, double& volt, const bool opt) const;
775
776 // Auxiliary functions for C type cells
777 double Ph2(const double xpos, const double ypos) const;
778 double Ph2Lim(const double radius) const {
779 return -log(abs(m_zmult) * radius * (1. - 3. * m_p1 + 5. * m_p2));
780 }
781 void E2Sum(const double xpos, const double ypos, double& ex,
782 double& ey) const;
783
784 // Mapping function for D30 type cells
785 void ConformalMap(const std::complex<double>& z, std::complex<double>& ww,
786 std::complex<double>& wd) const;
787
788 static bool InTube(const double x0, const double y0, const double a,
789 const int n);
790
791 bool SagDetailed(const Wire& wire, const std::vector<double>& xMap,
792 const std::vector<double>& yMap,
793 const std::vector<std::vector<double> >& fxMap,
794 const std::vector<std::vector<double> >& fyMap,
795 std::vector<double>& csag, std::vector<double>& xsag,
796 std::vector<double>& ysag) const;
797 bool GetForceRatio(const Wire& wire, const double coor,
798 const std::array<double, 2>& bend,
799 const std::array<double, 2>& dbend,
800 std::array<double, 2>& f, const std::vector<double>& xMap,
801 const std::vector<double>& yMap,
802 const std::vector<std::vector<double> >& fxMap,
803 const std::vector<std::vector<double> >& fyMap) const;
804 bool FindZeroes(const Wire& wire, const double h, std::vector<double>& x,
805 const std::vector<double>& xMap,
806 const std::vector<double>& yMap,
807 const std::vector<std::vector<double> >& fxMap,
808 const std::vector<std::vector<double> >& fyMap) const;
809 bool StepRKN(const Wire& wire, const double h, double& x,
810 std::array<double, 2>& y, std::array<double, 2>& yp,
811 const std::vector<double>& xMap, const std::vector<double>& yMap,
812 const std::vector<std::vector<double> >& fxMap,
813 const std::vector<std::vector<double> >& fyMap) const;
814 bool Trace(const Wire& wire, const double h, const std::vector<double>& xx,
815 std::vector<double>& f, const std::vector<double>& xMap,
816 const std::vector<double>& yMap,
817 const std::vector<std::vector<double> >& fxMap,
818 const std::vector<std::vector<double> >& fyMap) const;
819 size_t SignalLayer(const int mx, const int my) const;
820
822 const std::vector<std::string>& groups, std::vector<double>& vw0,
823 std::array<double, 5>& vp0, std::vector<double>& aFit,
824 std::vector<std::vector<unsigned int> >& wiresInGroup,
825 std::vector<std::vector<unsigned int> >& planesInGroup);
826};
827} // namespace Garfield
828
829#endif
double WpotPlaneD30(const double x, const double y, const std::vector< double > &qp) const
bool CrossedPlane(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, double &xc, double &yc, double &zc) override
void AddPlaneX(const double x, const double voltage, const std::string &label="")
Add a plane at constant x.
void PrintCell()
Print all available information on the cell.
void SetGravity(const double dx, const double dy, const double dz)
Set the gravity orientation.
bool GetVoltageRange(double &pmin, double &pmax) override
void SetPolarCoordinates()
Use polar coordinates.
bool FindZeroes(const Wire &wire, const double h, std::vector< double > &x, const std::vector< double > &xMap, const std::vector< double > &yMap, const std::vector< std::vector< double > > &fxMap, const std::vector< std::vector< double > > &fyMap) const
void InitialiseFitParameters(const std::vector< std::string > &groups, std::vector< double > &vw0, std::array< double, 5 > &vp0, std::vector< double > &aFit, std::vector< std::vector< unsigned int > > &wiresInGroup, std::vector< std::vector< unsigned int > > &planesInGroup)
void FieldB2X(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
void AddReadout(const std::string &label, const bool silent=false)
Request calculation of weighting field and potential for a given group of wires or planes.
void SetPeriodicityX(const double s)
Set the periodic length [cm] in the x-direction.
void FieldD20(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
bool IprB2X(const int my, std::vector< std::vector< std::complex< double > > > &mat)
void WfieldWireC2X(const double x, const double y, double &ex, double &ey, const std::vector< double > &qw) const
bool IprC30(std::vector< std::vector< std::complex< double > > > &mat)
void EnableDipoleTerms(const bool on=true)
Request dipole terms be included for each of the wires (default: off).
bool GetPeriodicityY(double &s)
Get the periodic length in the y-direction.
void SetCartesianCoordinates()
Use Cartesian coordinates (default).
bool CrossedWire(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, double &xc, double &yc, double &zc, const bool centre, double &rc) override
void DipoleFieldB1X(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
void FieldD10(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
void FieldAtWireD30(const double xpos, const double ypos, double &ex, double &ey, const std::vector< bool > &cnalso) const
void WfieldPlaneB2Y(const double x, const double ypos, double &ex, double &ey, const int mx, const std::vector< double > &qp) const
void AddPixelOnPlanePhi(const double phi, const double rmin, const double rmax, const double zmin, const double zmax, const std::string &label, const double gap=-1.)
Add a pixel on an existing plane at constant phi.
bool IprB2Y(const int mx, std::vector< std::vector< std::complex< double > > > &mat)
void WfieldStrip(const double x, const double y, const double g, const double w, double &fx, double &fy) const
std::vector< std::vector< std::vector< double > > > m_qplane
static bool InTube(const double x0, const double y0, const double a, const int n)
void DipoleFieldB2Y(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
void EnableChargeCheck(const bool on=true)
Check the quality of the capacitance matrix inversion (default: off).
void AddPixelOnPlaneX(const double x, const double ymin, const double ymax, const double zmin, const double zmax, const std::string &label, const double gap=-1., const double rot=0.)
Add a pixel on an existing plane at constant x.
bool IsPolar() const
Are polar coordinates being used?
std::string GetCellType(const Cell) const
void SetMedium(Medium *medium)
Set the medium inside the cell.
bool GetForceRatio(const Wire &wire, const double coor, const std::array< double, 2 > &bend, const std::array< double, 2 > &dbend, std::array< double, 2 > &f, const std::vector< double > &xMap, const std::vector< double > &yMap, const std::vector< std::vector< double > > &fxMap, const std::vector< std::vector< double > > &fyMap) const
void WfieldWireD30(const double x, const double y, double &ex, double &ey, const std::vector< double > &qw) const
void WfieldPlaneC30(const double x, const double y, double &ex, double &ey, const std::vector< double > &qp) const
void AddCharge(const double x, const double y, const double z, const double q)
Add a point charge.
void SetNumberOfCellCopies(const unsigned int nfourier)
bool Trace(const Wire &wire, const double h, const std::vector< double > &xx, std::vector< double > &f, const std::vector< double > &xMap, const std::vector< double > &yMap, const std::vector< std::vector< double > > &fxMap, const std::vector< std::vector< double > > &fyMap) const
double WpotWireB2X(const double x, const double y, const int my, const std::vector< double > &qw) const
double WpotPlaneD10(const double x, const double y, const std::vector< double > &qp) const
unsigned int GetNumberOfWires() const
Get the number of wires.
void FieldAtWireA00(const double xpos, const double ypos, double &ex, double &ey, const std::vector< bool > &cnalso) const
void SetScanningGrid(const unsigned int nX, const unsigned int nY)
Set the number of grid lines at which the electrostatic force as function of the wire displacement is...
void SetNumberOfSteps(const unsigned int n)
Set the number of integration steps within each shot (must be >= 1).
double WpotStripXy(const double x, const double y, const double z, const int ip, const Strip &strip) const
void WfieldWireB2X(const double x, const double y, double &ex, double &ey, const int my, const std::vector< double > &qw) const
std::string GetCellType()
Return the cell type.
void SetNumberOfShots(const unsigned int n)
Set the number of shots used for numerically solving the wire sag differential equation.
void FieldAtWireB1X(const double xpos, const double ypos, double &ex, double &ey, const std::vector< bool > &cnalso) const
void ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, double &v, Medium *&m, int &status) override
double Ph2Lim(const double radius) const
void AddPlaneY(const double y, const double voltage, const std::string &label="")
Add a plane at constant y.
size_t SignalLayer(const int mx, const int my) const
void FieldAtWireB2Y(const double xpos, const double ypos, double &ex, double &ey, const std::vector< bool > &cnalso) const
void E2Sum(const double xpos, const double ypos, double &ex, double &ey) const
void ConformalMap(const std::complex< double > &z, std::complex< double > &ww, std::complex< double > &wd) const
void AddPixelOnPlaneR(const double r, const double phimin, const double phimax, const double zmin, const double zmax, const std::string &label, const double gap=-1.)
Add a pixel on an existing plane at constant radius.
double WpotPlaneA00(const double x, const double y, const int mx, const int my, const std::vector< double > &qp) const
void WfieldWireB2Y(const double x, const double y, double &ex, double &ey, const int mx, const std::vector< double > &qw) const
bool GetWire(const unsigned int i, double &x, double &y, double &diameter, double &voltage, std::string &label, double &length, double &charge, int &ntrap) const
Retrieve the parameters of a wire.
bool GetPlaneR(const unsigned int i, double &r, double &voltage, std::string &label) const
Retrieve the parameters of a plane at constant radius.
bool WireDisplacement(const unsigned int iw, const bool detailed, std::vector< double > &csag, std::vector< double > &xsag, std::vector< double > &ysag, double &stretch, const bool print=true)
Compute the sag profile of a wire.
void WfieldPixel(const double x, const double y, const double z, double &ex, double &ey, double &ez, const int ip, const Pixel &pixel) const
void EnableGravity(const bool on=true)
Include gravity in the sag computation or not.
void FieldAtWireB1Y(const double xpos, const double ypos, double &ex, double &ey, const std::vector< bool > &cnalso) const
void GetGravity(double &dx, double &dy, double &dz) const
Get the gravity orientation.
bool Wfield(const double x, const double y, const double z, double &ex, double &ey, double &ez, const std::string &label) const
void AddStripOnPlaneX(const char direction, const double x, const double smin, const double smax, const std::string &label, const double gap=-1.)
Add a strip in the y or z direction on an existing plane at constant x.
void WfieldPlaneA00(const double x, const double y, double &ex, double &ey, const int mx, const int my, const std::vector< double > &qp) const
void FieldD30(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
void ClearCharges()
Remove all point charges.
void FieldB1Y(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
unsigned int GetNumberOfPlanesY() const
Get the number of equipotential planes at constant y.
void PlotCell(TPad *pad)
Make a plot of the cell layout.
double WpotPixel(const double x, const double y, const double z, const int ip, const Pixel &pixel) const
void Field3dB2X(const double x, const double y, const double z, double &ex, double &ey, double &ez, double &volt) const
double WpotPlaneC2Y(const double x, const double y, const std::vector< double > &qp) const
double WpotWireC2Y(const double x, const double y, const std::vector< double > &qw) const
bool IprA00(const int mx, const int my, std::vector< std::vector< std::complex< double > > > &mat)
double Ph2(const double xpos, const double ypos) const
bool IprC2X(std::vector< std::vector< std::complex< double > > > &mat)
std::vector< std::complex< double > > m_zw
bool ElectricFieldAtWire(const unsigned int iw, double &ex, double &ey)
Calculate the electric field at a given wire position, as if the wire itself were not there,...
bool SagDetailed(const Wire &wire, const std::vector< double > &xMap, const std::vector< double > &yMap, const std::vector< std::vector< double > > &fxMap, const std::vector< std::vector< double > > &fyMap, std::vector< double > &csag, std::vector< double > &xsag, std::vector< double > &ysag) const
bool MultipoleMoments(const unsigned int iw, const unsigned int order=4, const bool print=false, const bool plot=false, const double rmult=1., const double eps=1.e-4, const unsigned int nMaxIter=20)
Calculate multipole moments for a given wire.
bool Update(const std::vector< double > &vw, const std::array< double, 5 > &vp)
bool InTrapRadius(const double q0, const double x0, const double y0, const double z0, double &xw, double &yx, double &rw) override
void AddWire(const double x, const double y, const double diameter, const double voltage, const std::string &label="", const double length=100., const double tension=50., const double rho=19.3, const int ntrap=5)
Add a wire at (x, y) .
void EnableExtrapolation(const bool on=true)
Switch on/off extrapolation of electrostatic forces beyond the scanning area (default: off).
void WfieldWireC2Y(const double x, const double y, double &ex, double &ey, const std::vector< double > &qw) const
unsigned int GetNumberOfPlanesX() const
Get the number of equipotential planes at constant x.
void SetScanningAreaLargest()
Set the scanning area to the largest area around each wire which is free from other cell elements.
int Field(const double xin, const double yin, const double zin, double &ex, double &ey, double &ez, double &volt, const bool opt)
bool GetPeriodicityX(double &s)
Get the periodic length in the x-direction.
void WfieldPlaneC2Y(const double x, const double y, double &ex, double &ey, const std::vector< double > &qp) const
bool GetTube(double &r, double &voltage, int &nEdges, std::string &label) const
Retrieve the tube parameters.
void FieldB1X(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
bool GetPlaneX(const unsigned int i, double &x, double &voltage, std::string &label) const
Retrieve the parameters of a plane at constant x.
void AddTube(const double radius, const double voltage, const int nEdges, const std::string &label="")
Add a tube.
void DipoleFieldB2X(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
void Field3dB2Y(const double x, const double y, const double z, double &ex, double &ey, double &ez, double &volt) const
void WfieldPlaneD30(const double x, const double y, double &ex, double &ey, const std::vector< double > &qp) const
void FieldAtWireC10(const double xpos, const double ypos, double &ex, double &ey, const std::vector< bool > &cnalso) const
void AddPlanePhi(const double phi, const double voltage, const std::string &label="")
Add a plane at constant phi.
double WpotPlaneB2Y(const double x, const double y, const int mx, const std::vector< double > &qp) const
void FieldAtWireB2X(const double xpos, const double ypos, double &ex, double &ey, const std::vector< bool > &cnalso) const
bool IprD10(std::vector< std::vector< std::complex< double > > > &mat)
void WfieldStripZ(const double x, const double y, double &ex, double &ey, const int ip, const Strip &strip) const
bool OptimiseOnWires(const std::vector< std::string > &groups, const std::string &field_function, const double target, const std::vector< unsigned int > &wires, const bool print=true)
Vary the potential of selected electrodes to match (a function of) the potential or field on the surf...
void Field3dA00(const double x, const double y, const double z, double &ex, double &ey, double &ez, double &volt) const
bool GetElementaryCell(double &x0, double &y0, double &z0, double &x1, double &y1, double &z1) override
void SetPeriodicityPhi(const double phi)
Set the periodicity [degree] in phi.
void FieldAtWireD20(const double xpos, const double ypos, double &ex, double &ey, const std::vector< bool > &cnalso) const
void FieldC2Y(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
bool GetPeriodicityPhi(double &s)
Get the periodicity [degree] in phi.
void AddPlaneR(const double r, const double voltage, const std::string &label="")
Add a plane at constant radius.
std::vector< std::vector< std::vector< double > > > m_qwire
bool GetPlanePhi(const unsigned int i, double &phi, double &voltage, std::string &label) const
Retrieve the parameters of a plane at constant phi.
void DipoleFieldA00(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
void WfieldPlaneB2X(const double x, const double y, double &ex, double &ey, const int my, const std::vector< double > &qp) const
void AddStripOnPlaneY(const char direction, const double y, const double smin, const double smax, const std::string &label, const double gap=-1.)
Add a strip in the x or z direction on an existing plane at constant y.
void DipoleFieldB1Y(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
double WpotWireB2Y(const double x, const double y, const int mx, const std::vector< double > &qw) const
void FieldAtWireC30(const double xpos, const double ypos, double &ex, double &ey, const std::vector< bool > &cnalso) const
void WfieldWireA00(const double x, const double y, double &ex, double &ey, const int mx, const int my, const std::vector< double > &qw) const
void FieldC10(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
double WpotPlaneC2X(const double x, const double y, const std::vector< double > &qp) const
bool GetBoundingBox(double &x0, double &y0, double &z0, double &x1, double &y1, double &z1) override
void WfieldWireD10(const double x, const double y, double &ex, double &ey, const std::vector< double > &qw) const
bool IprC2Y(std::vector< std::vector< std::complex< double > > > &mat)
void FieldC30(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
void WfieldStripXy(const double x, const double y, const double z, double &ex, double &ey, double &ez, const int ip, const Strip &strip) const
bool ForcesOnWire(const unsigned int iw, std::vector< double > &xMap, std::vector< double > &yMap, std::vector< std::vector< double > > &fxMap, std::vector< std::vector< double > > &fyMap)
Calculate a table of the forces acting on a wire.
void WfieldWireC30(const double x, const double y, double &ex, double &ey, const std::vector< double > &qw) const
void FieldC2X(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
void WeightingField(const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label) override
bool StepRKN(const Wire &wire, const double h, double &x, std::array< double, 2 > &y, std::array< double, 2 > &yp, const std::vector< double > &xMap, const std::vector< double > &yMap, const std::vector< std::vector< double > > &fxMap, const std::vector< std::vector< double > > &fyMap) const
void ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status) override
void FieldA00(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
void SetOptimisationParameters(const double dist=1., const double eps=1.e-4, const unsigned int nMaxIter=10)
Set the conditions at which to allow the iteration to stop.
void Field3dD10(const double x, const double y, const double z, double &ex, double &ey, double &ez, double &volt) const
void WfieldPlaneD10(const double x, const double y, double &ex, double &ey, const std::vector< double > &qp) const
void SetScanningAreaFirstOrder(const double scale=2.)
Set the scanning area based on the zeroth-order estimates of the wire shift, enlarged by a scaling fa...
double WpotStripZ(const double x, const double y, const int ip, const Strip &strip) const
double Wpot(const double x, const double y, const double z, const std::string &label) const
void FieldAtWireC2Y(const double xpos, const double ypos, double &ex, double &ey, const std::vector< bool > &cnalso) const
void AddPixelOnPlaneY(const double y, const double xmin, const double xmax, const double zmin, const double zmax, const std::string &label, const double gap=-1., const double rot=0.)
Add a pixel on an existing plane at constant y.
void EnableElectrostaticForce(const bool on=true)
Include the electrostatic force in the sag computation or not.
double WeightingPotential(const double x, const double y, const double z, const std::string &label) override
void SetPeriodicityY(const double s)
Set the periodic length [cm] in the y-direction.
double WpotWireD30(const double x, const double y, const std::vector< double > &qw) const
Medium * GetMedium(const double x, const double y, const double z) override
unsigned int GetNumberOfPlanesPhi() const
Get the number of equipotential planes at constant phi.
void FieldAtWireC2X(const double xpos, const double ypos, double &ex, double &ey, const std::vector< bool > &cnalso) const
bool Charge(std::vector< std::vector< double > > &mat)
double WpotWireC2X(const double x, const double y, const std::vector< double > &qw) const
double WpotPlaneC30(const double x, const double y, const std::vector< double > &qp) const
double WpotPlaneB2X(const double x, const double y, const int my, const std::vector< double > &qp) const
void AddStripOnPlaneR(const char direction, const double r, const double smin, const double smax, const std::string &label, const double gap=-1.)
Add a strip in the phi or z direction on an existing plane at constant radius.
void PrintCharges() const
Print a list of the point charges.
double WpotWireC30(const double x, const double y, const std::vector< double > &qw) const
bool OptimiseOnGrid(const std::vector< std::string > &groups, const std::string &field_function, const double target, const double x0, const double y0, const double x1, const double y1, const unsigned int nX=10, const unsigned int nY=10, const bool print=true)
Vary the potential of selected electrodes to match (a function of) the potential or field on an x-y (...
void SetScanningArea(const double xmin, const double xmax, const double ymin, const double ymax)
Specify explicitly the boundaries of the the scanning area (i.
bool OptimiseOnTrack(const std::vector< std::string > &groups, const std::string &field_function, const double target, const double x0, const double y0, const double x1, const double y1, const unsigned int nP=20, const bool print=true)
Vary the potential of selected electrodes to match (a function of) the potential or field along a giv...
double WpotWireD10(const double x, const double y, const std::vector< double > &qw) const
bool GetPlaneY(const unsigned int i, double &y, double &voltage, std::string &label) const
Retrieve the parameters of a plane at constant y.
double WpotWireA00(const double x, const double y, const int mx, const int my, const std::vector< double > &qw) const
void FieldB2Y(const double xpos, const double ypos, double &ex, double &ey, double &volt, const bool opt) const
bool IprD30(std::vector< std::vector< std::complex< double > > > &mat)
void FieldAtWireD10(const double xpos, const double ypos, double &ex, double &ey, const std::vector< bool > &cnalso) const
unsigned int GetNumberOfPlanesR() const
Get the number of equipotential planes at constant radius.
void AddStripOnPlanePhi(const char direction, const double phi, const double smin, const double smax, const std::string &label, const double gap=-1.)
Add a strip in the r or z direction on an existing plane at constant phi.
void WfieldPlaneC2X(const double x, const double y, double &ex, double &ey, const std::vector< double > &qp) const
double gap
Distance to the opposite electrode.
double ewycor
Background weighting fields.
double gap
Distance to the opposite electrode.
int nTrap
Trap radius. Particle is "trapped" if within nTrap * radius of wire.