1#ifndef G_COMPONENT_ANALYTIC_FIELD_H
2#define G_COMPONENT_ANALYTIC_FIELD_H
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 =
"");
49 const std::string& label =
"");
59 const double smax,
const std::string& label,
60 const double gap = -1.);
63 const double smax,
const std::string& label,
64 const double gap = -1.);
68 const double smax,
const std::string& label,
69 const double gap = -1.);
72 const double smin,
const double smax,
73 const std::string& label,
const double gap = -1.);
85 const double zmin,
const double zmax,
86 const std::string& label,
const double gap = -1.,
87 const double rot = 0.);
90 const double zmin,
const double zmax,
91 const std::string& label,
const double gap = -1.,
92 const double rot = 0.);
95 const double phimax,
const double zmin,
96 const double zmax,
const std::string& label,
97 const double gap = -1.);
100 const double rmax,
const double zmin,
101 const double zmax,
const std::string& label,
102 const double gap = -1.);
130 void AddCharge(
const double x,
const double y,
const double z,
165 void AddReadout(
const std::string& label,
const bool silent =
false);
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);
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;
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;
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);
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);
247 const std::string& field_function,
const double target,
248 const std::vector<unsigned int>& wires,
249 const bool print =
true);
258 const double eps = 1.e-4,
259 const unsigned int nMaxIter = 10);
283 void SetGravity(
const double dx,
const double dy,
const double dz);
300 std::vector<double>& yMap,
301 std::vector<std::vector<double> >& fxMap,
302 std::vector<std::vector<double> >& fyMap);
314 std::vector<double>& csag, std::vector<double>& xsag,
315 std::vector<double>& ysag,
double& stretch,
316 const bool print =
true);
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;
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;
331 using Component::ElectricField;
335 double& wx,
double& wy,
double& wz,
336 const std::string& label)
override {
339 Wfield(x, y, z, wx, wy, wz, label);
342 const std::string& label)
override {
344 return Wpot(x, y, z, label);
348 double& y1,
double& z1)
override;
350 double& y1,
double& z1)
override;
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;
358 const double z0,
double& xw,
double& yx,
359 double& rw)
override;
362 const double x1,
const double y1,
const double z1,
363 double& xc,
double& yc,
double& zc)
override;
457 std::vector<std::complex<double> >
m_zw;
513 std::vector<std::vector<std::vector<double> > >
m_qwire;
515 std::vector<std::vector<std::vector<double> > >
m_qplane;
579 bool Update(
const std::vector<double>& vw,
const std::array<double, 5>& vp);
594 std::vector<std::vector<std::complex<double> > >& mat);
596 std::vector<std::vector<std::complex<double> > >& mat);
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);
608 bool Charge(std::vector<std::vector<double> >& mat);
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;
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;
648 bool Wfield(
const double x,
const double y,
const double z,
double& ex,
649 double& ey,
double& ez,
const std::string& label)
const;
651 const int mx,
const int my,
652 const std::vector<double>& qw)
const;
654 const int my,
const std::vector<double>& qw)
const;
656 const int mx,
const std::vector<double>& qw)
const;
658 const std::vector<double>& qw)
const;
660 const std::vector<double>& qw)
const;
662 const std::vector<double>& qw)
const;
664 const std::vector<double>& qw)
const;
666 const std::vector<double>& qw)
const;
668 const int mx,
const int my,
669 const std::vector<double>& qp)
const;
671 const int my,
const std::vector<double>& qp)
const;
673 const int mx,
const std::vector<double>& qp)
const;
675 const std::vector<double>& qp)
const;
677 const std::vector<double>& qp)
const;
679 const std::vector<double>& qp)
const;
681 const std::vector<double>& qp)
const;
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;
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;
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;
701 const std::vector<double>& qw)
const;
703 const std::vector<double>& qw)
const;
705 const std::vector<double>& qw)
const;
707 const std::vector<double>& qw)
const;
709 const std::vector<double>& qw)
const;
711 const std::vector<double>& qw)
const;
713 const std::vector<double>& qw)
const;
715 const int my,
const std::vector<double>& qp)
const;
717 const std::vector<double>& qp)
const;
719 const std::vector<double>& qp)
const;
721 const std::vector<double>& qp)
const;
723 const std::vector<double>& qp)
const;
725 const std::vector<double>& qp)
const;
727 const std::vector<double>& qp)
const;
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;
741 double& ey,
const std::vector<bool>& cnalso)
const;
743 double& ey,
const std::vector<bool>& cnalso)
const;
745 double& ey,
const std::vector<bool>& cnalso)
const;
747 double& ey,
const std::vector<bool>& cnalso)
const;
749 double& ey,
const std::vector<bool>& cnalso)
const;
751 double& ey,
const std::vector<bool>& cnalso)
const;
753 double& ey,
const std::vector<bool>& cnalso)
const;
755 double& ey,
const std::vector<bool>& cnalso)
const;
757 double& ey,
const std::vector<bool>& cnalso)
const;
759 double& ey,
const std::vector<bool>& cnalso)
const;
761 double& ey,
const std::vector<bool>& cnalso)
const;
763 double& ey,
const std::vector<bool>& cnalso)
const;
766 double& ey,
double& volt,
const bool opt)
const;
768 double& ey,
double& volt,
const bool opt)
const;
770 double& ey,
double& volt,
const bool opt)
const;
772 double& ey,
double& volt,
const bool opt)
const;
774 double& ey,
double& volt,
const bool opt)
const;
777 double Ph2(
const double xpos,
const double ypos)
const;
778 double Ph2Lim(
const double radius)
const {
781 void E2Sum(
const double xpos,
const double ypos,
double& ex,
785 void ConformalMap(
const std::complex<double>& z, std::complex<double>& ww,
786 std::complex<double>& wd)
const;
788 static bool InTube(
const double x0,
const double y0,
const double a,
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;
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;
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;
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;
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);
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.
std::complex< double > m_zmult
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.
unsigned int m_nTermBessel
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
std::vector< double > m_amp2
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
bool m_useGravitationalForce
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
std::array< Plane, 5 > m_planes
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).
bool m_useElectrostaticForce
double WpotStripXy(const double x, const double y, const double z, const int ip, const Strip &strip) const
double StepSizeHint() override
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.
ComponentAnalyticField()
Constructor.
bool Wfield(const double x, const double y, const double z, double &ex, double &ey, double &ez, const std::string &label) const
std::vector< double > m_cosph2
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
std::array< double, 3 > m_down
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.
std::vector< double > m_b2sin
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
ScanningRange m_scanRange
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.
std::vector< std::string > m_readout
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
std::vector< Charge3d > m_ch3d
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.
~ComponentAnalyticField()
Destructor.
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
void UpdatePeriodicity() override
double Wpot(const double x, const double y, const double z, const std::string &label) const
std::vector< double > m_sinph2
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 smax
Coordinates in x/y.
double zmax
Coordinates in z.
double gap
Distance to the opposite electrode.
std::vector< Strip > strips1
x/y strips.
double ewycor
Background weighting fields.
std::vector< Pixel > pixels
Pixels.
std::vector< Strip > strips2
z strips.
double gap
Distance to the opposite electrode.
double tension
Stretching weight.
int nTrap
Trap radius. Particle is "trapped" if within nTrap * radius of wire.