Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ViewCell.cc
Go to the documentation of this file.
1#include <cmath>
2#include <iostream>
3
4#include <TEllipse.h>
5#include <TGeoBBox.h>
6#include <TGeoPgon.h>
7#include <TLine.h>
8#include <TMarker.h>
9#include <TPolyLine.h>
10
13#include "Garfield/Plotting.hh"
14#include "Garfield/ViewCell.hh"
15
16namespace Garfield {
17
18ViewCell::ViewCell() : ViewBase("ViewCell") {}
19
21 if (!cmp) {
22 std::cerr << m_className << "::SetComponent: Null pointer.\n";
23 return;
24 }
25 m_component = cmp;
26}
27
29 if (!cmp) {
30 std::cerr << m_className << "::SetComponent: Null pointer.\n";
31 return;
32 }
33 m_nebem = cmp;
34}
35
37 if (!Plot(true)) {
38 std::cerr << m_className << "::Plot2d: Error creating plot.\n";
39 }
40}
41
43 if (!Plot(false)) {
44 std::cerr << m_className << "::Plot3d: Error creating plot.\n";
45 }
46}
47
48bool ViewCell::Plot(const bool twod) {
49 if (!m_component && !m_nebem) {
50 std::cerr << m_className << "::Plot: Component is not defined.\n";
51 return false;
52 }
53
54 double pmin = 0., pmax = 0.;
55 if (m_component) {
56 if (!m_component->GetVoltageRange(pmin, pmax)) {
57 std::cerr << m_className << "::Plot: Component is not ready.\n";
58 return false;
59 }
60 } else {
61 if (!m_nebem->GetVoltageRange(pmin, pmax)) {
62 std::cerr << m_className << "::Plot: Component is not ready.\n";
63 return false;
64 }
65 }
66
67 // Get the bounding box.
68 double x0 = m_xMinBox, y0 = m_yMinBox, z0 = m_zMinBox;
69 double x1 = m_xMaxBox, y1 = m_yMaxBox, z1 = m_zMaxBox;
70 if (twod && m_userPlotLimits) {
71 x0 = m_xMinPlot;
72 y0 = m_yMinPlot;
73 x1 = m_xMaxPlot;
74 y1 = m_yMaxPlot;
75 } else if (!m_userBox) {
76 if (m_component) {
77 if (!m_component->GetBoundingBox(x0, y0, z0, x1, y1, z1)) {
78 std::cerr << m_className << "::Plot:\n"
79 << " Bounding box cannot be determined.\n"
80 << " Call SetArea first.\n";
81 return false;
82 }
83 } else {
84 if (!m_nebem->GetBoundingBox(x0, y0, z0, x1, y1, z1)) {
85 std::cerr << m_className << "::Plot:\n"
86 << " Bounding box cannot be determined.\n"
87 << " Call SetArea first.\n";
88 return false;
89 }
90 }
91 }
92 const double dx = std::max(fabs(x0), fabs(x1));
93 const double dy = std::max(fabs(y0), fabs(y1));
94 const double dz = std::max(fabs(z0), fabs(z1));
95
96 auto canvas = GetCanvas();
97 canvas->cd();
98 canvas->SetTitle("Cell layout");
99
100 if (twod) {
101 bool empty = false;
102 if (!gPad ||
103 (gPad->GetListOfPrimitives()->GetSize() == 0 && gPad->GetX1() == 0 &&
104 gPad->GetX2() == 1 && gPad->GetY1() == 0 && gPad->GetY2() == 1)) {
105 empty = true;
106 }
107 const double bm = canvas->GetBottomMargin();
108 const double lm = canvas->GetLeftMargin();
109 const double rm = canvas->GetRightMargin();
110 const double tm = canvas->GetTopMargin();
111 if (!empty) {
112 TPad* pad = new TPad("cell", "", 0, 0, 1, 1);
113 pad->SetFillStyle(0);
114 pad->SetFrameFillStyle(0);
115 pad->Draw();
116 pad->cd();
117 }
118 gPad->Range(x0 - (x1 - x0) * (lm / (1. - rm - lm)),
119 y0 - (y1 - y0) * (bm / (1. - tm - lm)),
120 x1 + (x1 - x0) * (rm / (1. - rm - lm)),
121 y1 + (y1 - y0) * (tm / (1. - tm - lm)));
122 }
123
124 if (m_nebem) return PlotNeBem(twod);
125
126 // Get the cell type.
127 const std::string cellType = m_component->GetCellType();
128
129 // Get the x/y periodicities.
130 double sx = 0., sy = 0.;
131 const bool perX = m_component->GetPeriodicityX(sx);
132 const bool perY = m_component->GetPeriodicityY(sy);
133 // Determine the number of periods present in the cell.
134 const int nMinX = perX ? int(x0 / sx) - 1 : 0;
135 const int nMaxX = perX ? int(x1 / sx) + 1 : 0;
136 const int nMinY = perY ? int(y0 / sy) - 1 : 0;
137 const int nMaxY = perY ? int(y1 / sy) + 1 : 0;
138 // Get the phi periodicity.
139 double sphi = 0.;
140 const bool perPhi = m_component->GetPeriodicityPhi(sphi);
141 const int nPhi = perPhi ? int(360. / sphi) : 0;
142 sphi *= DegreeToRad;
143 if (!twod) SetupGeo(dx, dy, dz);
144 const bool polar = m_component->IsPolar();
145
146 // Get the number of wires.
147 const unsigned int nWires = m_component->GetNumberOfWires();
148 std::vector<std::string> wireTypes;
149 // Loop over the wires.
150 for (unsigned int i = 0; i < nWires; ++i) {
151 double xw = 0., yw = 0., dw = 0., vw = 0., lw = 0., qw = 0.;
152 std::string lbl;
153 int nTrap;
154 m_component->GetWire(i, xw, yw, dw, vw, lbl, lw, qw, nTrap);
155 auto it = std::find(wireTypes.begin(), wireTypes.end(), lbl);
156 const int type = std::distance(wireTypes.begin(), it);
157 if (it == wireTypes.end()) wireTypes.push_back(lbl);
158 if (polar) {
159 const double r = xw;
160 for (int j = 0; j <= nPhi; ++j) {
161 const double phi = yw * Garfield::DegreeToRad + j * sphi;
162 const double x = r * cos(phi);
163 const double y = r * sin(phi);
164 if (twod) {
165 PlotWire(x, y, dw, type);
166 } else {
167 PlotWire(x, y, dw, type, std::min(0.5 * lw, dz));
168 }
169 }
170 continue;
171 }
172 for (int nx = nMinX; nx <= nMaxX; ++nx) {
173 const double x = xw + nx * sx;
174 if (x + 0.5 * dw <= x0 || x - 0.5 * dw >= x1) continue;
175 for (int ny = nMinY; ny <= nMaxY; ++ny) {
176 const double y = yw + ny * sy;
177 if (y + 0.5 * dw <= y0 || y - 0.5 * dw >= y1) continue;
178 if (twod) {
179 PlotWire(x, y, dw, type);
180 } else {
181 PlotWire(x, y, dw, type, std::min(0.5 * lw, dz));
182 }
183 }
184 }
185 }
186
187 // Draw the x planes.
188 const unsigned int nPlanesX = m_component->GetNumberOfPlanesX();
189 for (unsigned int i = 0; i < nPlanesX; ++i) {
190 double xp = 0., vp = 0.;
191 std::string lbl;
192 m_component->GetPlaneX(i, xp, vp, lbl);
193 for (int nx = nMinX; nx <= nMaxX; ++nx) {
194 const double x = xp + nx * sx;
195 if (x < x0 || x > x1) continue;
196 if (twod) {
197 PlotPlane(x, y0, x, y1);
198 } else {
199 const double width = std::min(0.01 * dx, 0.01 * dy);
200 PlotPlane(width, dy, dz, x, 0);
201 }
202 }
203 }
204
205 // Draw the y planes.
206 const unsigned int nPlanesY = m_component->GetNumberOfPlanesY();
207 for (unsigned int i = 0; i < nPlanesY; ++i) {
208 double yp = 0., vp = 0.;
209 std::string lbl;
210 m_component->GetPlaneY(i, yp, vp, lbl);
211 for (int ny = nMinY; ny <= nMaxY; ++ny) {
212 const double y = yp + ny * sy;
213 if (y < y0 || y > y1) continue;
214 if (twod) {
215 PlotPlane(x0, y, x1, y);
216 } else {
217 const double width = std::min(0.01 * dx, 0.01 * dy);
218 PlotPlane(dx, width, dz, 0, y);
219 }
220 }
221 }
222
223 // Draw the r and phi planes.
224 const unsigned int nPlanesR = m_component->GetNumberOfPlanesR();
225 const unsigned int nPlanesPhi = m_component->GetNumberOfPlanesPhi();
226 if (nPlanesR > 0 || nPlanesPhi > 0) {
227 double vp = 0.;
228 std::string lbl;
229 std::array<double, 2> phi = {0., 360.};
230 double dphi = 360.;
231 if (nPlanesPhi == 2 && !m_component->GetPeriodicityPhi(dphi)) {
232 m_component->GetPlanePhi(0, phi[0], vp, lbl);
233 m_component->GetPlanePhi(1, phi[1], vp, lbl);
234 }
235 double w = 0.01 * std::min(dx, dy);
236 std::array<double, 2> r = {0., 0.};
237 if (nPlanesR > 0) {
238 m_component->GetPlaneR(0, r[0], vp, lbl);
239 w = 0.02 * r[0];
240 }
241 if (nPlanesR == 2) {
242 m_component->GetPlaneR(1, r[1], vp, lbl);
243 w = 0.01 * (r[1] - r[0]);
244 }
245 for (unsigned int i = 0; i < nPlanesR; ++i) {
246 if (twod) {
247 TEllipse circle;
248 circle.SetDrawOption("same");
249 circle.SetFillStyle(0);
250 circle.SetNoEdges();
251 circle.DrawEllipse(0, 0, r[i], r[i], phi[0], phi[1], 0);
252 continue;
253 }
254 const auto med = m_geo->GetMedium("Metal");
255 const auto tubs = m_geo->MakeTubs("Plane", med, r[i] - w, r[i] + w,
256 dz, phi[0], phi[1]);
257 tubs->SetLineColor(kGreen - 5);
258 tubs->SetTransparency(75);
259 m_geo->GetTopVolume()->AddNode(tubs, 1);
260 }
261 if (nPlanesR == 1) {
262 std::swap(r[0], r[1]);
263 } else if (nPlanesR == 0) {
264 r[1] = std::max(dx, dy);
265 }
266 for (unsigned int i = 0; i < nPlanesPhi; ++i) {
267 const double cp = cos(phi[i] * DegreeToRad);
268 const double sp = sin(phi[i] * DegreeToRad);
269 if (twod) {
270 PlotPlane(r[0] * cp, r[0] * sp, r[1] * cp, r[1] * sp);
271 continue;
272 }
273 const auto med = m_geo->GetMedium("Metal");
274 const double dr = 0.5 * (r[1] - r[0]);
275 TGeoVolume* plane = m_geo->MakeBox("Plane", med, dr, w, dz);
276 plane->SetLineColor(kGreen - 5);
277 plane->SetTransparency(75);
278 TGeoRotation* rot = new TGeoRotation("", phi[i], 0, 0);
279 const double rm = 0.5 * (r[0] + r[1]);
280 TGeoCombiTrans* tr = new TGeoCombiTrans(cp * rm, sp * rm, 0, rot);
281 m_geo->GetTopVolume()->AddNode(plane, 1, tr);
282 }
283 }
284 double rt = 0., vt = 0.;
285 int nt = 0;
286 std::string lbl;
287 if (m_component->GetTube(rt, vt, nt, lbl)) {
288 if (twod) {
289 PlotTube(0., 0., rt, nt);
290 } else {
291 PlotTube(0., 0., 0.98 * rt, 1.02 * rt, nt, dz);
292 }
293 }
294
295 if (!twod) {
296 m_geo->CloseGeometry();
297 m_geo->GetTopNode()->Draw("ogl");
298 } else {
299 gPad->Update();
300 }
301
302 return true;
303}
304
305bool ViewCell::PlotNeBem(const bool twod) {
306
307 if (!twod) {
308 std::cerr << m_className << "::PlotNeBem: 3D plot not implemented yet.\n";
309 return false;
310 }
311
312 // Draw the regions.
313 const unsigned int nRegions = m_nebem->GetNumberOfRegions();
314 for (unsigned int i = nRegions; i-- > 0;) {
315 std::vector<double> xv;
316 std::vector<double> yv;
317 Medium* medium = nullptr;
318 unsigned int bctype = 1;
319 double v = 0.;
320 if (!m_nebem->GetRegion(i, xv, yv, medium, bctype, v)) continue;
321 const unsigned int n = xv.size();
322 if (n < 3) continue;
323 TLine line;
324 line.SetDrawOption("same");
325 if (bctype == 4) {
326 line.SetLineStyle(2);
327 } else {
328 line.SetLineStyle(1);
329 }
330 for (unsigned int j = 0; j < n; ++j) {
331 const unsigned int k = j < n - 1 ? j + 1 : 0;
332 line.DrawLine(xv[j], yv[j], xv[k], yv[k]);
333 }
334 }
335
336 // Draw the wires.
337 const unsigned int nWires = m_nebem->GetNumberOfWires();
338 for (unsigned int i = 0; i < nWires; ++i) {
339 double x = 0., y = 0., d = 0., v = 0., q = 0.;
340 if (!m_nebem->GetWire(i, x, y, d, v, q)) continue;
341 PlotWire(x, y, d, 0);
342 }
343
344 // Draw the straight-line segments.
345 const unsigned int nSegments = m_nebem->GetNumberOfSegments();
346 for (unsigned int i = 0; i < nSegments; ++i) {
347 double x0 = 0., y0 = 0., x1 = 0., y1 = 0., v = 0.;
348 if (!m_nebem->GetSegment(i, x0, y0, x1, y1, v)) continue;
349 PlotPlane(x0, y0, x1, y1);
350 }
351
352 gPad->Update();
353 return true;
354}
355
356void ViewCell::SetupGeo(const double dx, const double dy, const double dz) {
357
358 if (!m_geo) {
359 gGeoManager = nullptr;
360 m_geo.reset(new TGeoManager("ViewCellGeoManager", "Cell layout"));
361 TGeoMaterial* matVacuum = new TGeoMaterial("Vacuum", 0., 0., 0.);
362 TGeoMaterial* matMetal = new TGeoMaterial("Metal", 63.546, 29., 8.92);
363 TGeoMedium* medVacuum = new TGeoMedium("Vacuum", 0, matVacuum);
364 TGeoMedium* medMetal = new TGeoMedium("Metal", 1, matMetal);
365 m_geo->AddMaterial(matVacuum);
366 m_geo->AddMaterial(medMetal->GetMaterial());
367 TGeoVolume* world =
368 m_geo->MakeBox("World", medVacuum, 1.05 * dx, 1.05 * dy, 1.05 * dz);
369 m_geo->SetTopVolume(world);
370 } else {
371 TGeoVolume* top = m_geo->GetTopVolume();
372 TGeoBBox* box = dynamic_cast<TGeoBBox*>(top);
373 double halfLenghts[3] = {1.05 * dx, 1.05 * dy, 1.05 * dz};
374 if (box) box->SetDimensions(halfLenghts);
375 }
376}
377
378void ViewCell::PlotWire(const double x, const double y, const double d,
379 const int type) {
380 if (m_useWireMarker) {
381 int markerStyle = 1;
382 if (type == 0) {
383 markerStyle = kFullCircle;
384 } else if (type == 1) {
385 markerStyle = kOpenCircle;
386 } else if (type == 2) {
387 markerStyle = kFullSquare;
388 } else if (type == 3) {
389 markerStyle = kOpenSquare;
390 } else {
391 markerStyle = 26 + type;
392 }
393 TMarker marker;
394 marker.SetMarkerStyle(markerStyle);
395 marker.SetDrawOption("Psame");
396 marker.DrawMarker(x, y);
397 return;
398 }
399
400 TEllipse circle;
401 circle.SetDrawOption("same");
402 circle.SetFillColor(kWhite);
403 circle.DrawEllipse(x, y, 0.5 * d, 0.5 * d, 0, 360, 0);
404}
405
406void ViewCell::PlotWire(const double x, const double y, const double d,
407 const int type, const double lz) {
408
409 const auto medium = m_geo->GetMedium("Metal");
410 TGeoVolume* wire = m_geo->MakeTube("Wire", medium, 0., 0.5 * d, lz);
411 switch (type) {
412 case 0:
413 wire->SetLineColor(kGray + 2);
414 break;
415 case 1:
416 wire->SetLineColor(kRed + 2);
417 break;
418 case 2:
419 wire->SetLineColor(kPink + 3);
420 break;
421 case 3:
422 wire->SetLineColor(kCyan + 3);
423 break;
424 default:
425 wire->SetLineColor(kBlue + type);
426 break;
427 }
428 m_geo->GetTopVolume()->AddNode(wire, 1, new TGeoTranslation(x, y, 0.));
429}
430
431void ViewCell::PlotTube(const double x0, const double y0, const double r,
432 const int n) {
433 if (n <= 0) {
434 TEllipse circle;
435 circle.SetDrawOption("same");
436 circle.SetFillStyle(0);
437 circle.DrawEllipse(x0, y0, r, r, 0, 360, 0);
438 return;
439 }
440
441 std::vector<double> x(n + 1);
442 std::vector<double> y(n + 1);
443 for (int i = 0; i <= n; ++i) {
444 const double phi = i * TwoPi / double(n);
445 x[i] = x0 + r * cos(phi);
446 y[i] = y0 + r * sin(phi);
447 }
448 TPolyLine pline;
449 pline.SetDrawOption("same");
450 pline.DrawPolyLine(n + 1, x.data(), y.data());
451}
452
453void ViewCell::PlotTube(const double x0, const double y0,
454 const double r1, const double r2, const int n,
455 const double lz) {
456
457 TGeoVolume* tube = nullptr;
458 if (n <= 0) {
459 // Round tube.
460 tube = m_geo->MakeTube("Tube", m_geo->GetMedium("Metal"), r1, r2, lz);
461 } else {
462 tube = m_geo->MakePgon("Tube", m_geo->GetMedium("Metal"), 0., 360., n, 2);
463 TGeoPgon* pgon = dynamic_cast<TGeoPgon*>(tube->GetShape());
464 if (pgon) {
465 pgon->DefineSection(0, -lz, r1, r2);
466 pgon->DefineSection(1, +lz, r1, r2);
467 }
468 }
469 tube->SetLineColor(kGreen + 2);
470 tube->SetTransparency(75);
471 m_geo->GetTopVolume()->AddNode(tube, 1, new TGeoTranslation(x0, y0, 0));
472}
473
474void ViewCell::PlotPlane(const double x0, const double y0,
475 const double x1, const double y1) {
476
477 TLine line;
478 line.SetDrawOption("same");
479 line.DrawLine(x0, y0, x1, y1);
480}
481
482void ViewCell::PlotPlane(const double dx, const double dy, const double dz,
483 const double x0, const double y0) {
484
485 const auto medium = m_geo->GetMedium("Metal");
486 TGeoVolume* plane = m_geo->MakeBox("Plane", medium, dx, dy, dz);
487 plane->SetLineColor(kGreen - 5);
488 plane->SetTransparency(75);
489 m_geo->GetTopVolume()->AddNode(plane, 1, new TGeoTranslation(x0, y0, 0));
490}
491
492}
bool GetVoltageRange(double &pmin, double &pmax) override
Calculate the voltage range [V].
bool GetPeriodicityY(double &s)
Get the periodic length in the y-direction.
bool IsPolar() const
Are polar coordinates being used?
unsigned int GetNumberOfWires() const
Get the number of wires.
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.
unsigned int GetNumberOfPlanesY() const
Get the number of equipotential planes at constant y.
unsigned int GetNumberOfPlanesX() const
Get the number of equipotential planes at constant x.
bool GetPeriodicityX(double &s)
Get the periodic length in the x-direction.
bool GetTube(double &r, double &voltage, int &nEdges, std::string &label) const
Retrieve the tube parameters.
bool GetPlaneX(const unsigned int i, double &x, double &voltage, std::string &label) const
Retrieve the parameters of a plane at constant x.
bool GetPeriodicityPhi(double &s)
Get the periodicity [degree] in phi.
bool GetPlanePhi(const unsigned int i, double &phi, double &voltage, std::string &label) const
Retrieve the parameters of a plane at constant phi.
bool GetBoundingBox(double &x0, double &y0, double &z0, double &x1, double &y1, double &z1) override
Get the bounding box coordinates.
unsigned int GetNumberOfPlanesPhi() const
Get the number of equipotential planes at constant phi.
bool GetPlaneY(const unsigned int i, double &y, double &voltage, std::string &label) const
Retrieve the parameters of a plane at constant y.
unsigned int GetNumberOfPlanesR() const
Get the number of equipotential planes at constant radius.
Two-dimensional implementation of the nearly exact Boundary Element Method.
unsigned int GetNumberOfSegments() const
Return the number of conducting straight-line segments.
bool GetVoltageRange(double &vmin, double &vmax) override
Calculate the voltage range [V].
bool GetSegment(const unsigned int i, double &x0, double &y0, double &x1, double &x2, double &v) const
Return the coordinates and voltage of a given straight-line segment.
bool GetBoundingBox(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax) override
Get the bounding box coordinates.
unsigned int GetNumberOfRegions() const
Return the number of regions.
bool GetRegion(const unsigned int i, std::vector< double > &xv, std::vector< double > &yv, Medium *&medium, unsigned int &bctype, double &v)
Return the properties of a given region.
bool GetWire(const unsigned int i, double &x, double &y, double &d, double &v, double &q) const
Return the coordinates, diameter, potential and charge of a given wire.
unsigned int GetNumberOfWires() const
Return the number of wires.
Base class for visualization classes.
Definition: ViewBase.hh:18
std::string m_className
Definition: ViewBase.hh:78
TPad * GetCanvas()
Retrieve the canvas.
Definition: ViewBase.cc:74
void Plot3d()
Make a three-dimensional drawing of the cell layout (using TGeo).
Definition: ViewCell.cc:42
void SetComponent(ComponentAnalyticField *comp)
Set the component for which to draw the cell geometry.
Definition: ViewCell.cc:20
ViewCell()
Constructor.
Definition: ViewCell.cc:18
void Plot2d()
Make a two-dimensional drawing of the cell layout.
Definition: ViewCell.cc:36
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:432
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:384