Garfield++ 5.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 if (m_nebem) {
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 const bool empty = !RangeSet(gPad);
102 const double bm = canvas->GetBottomMargin();
103 const double lm = canvas->GetLeftMargin();
104 const double rm = canvas->GetRightMargin();
105 const double tm = canvas->GetTopMargin();
106 if (!empty) {
107 TPad* pad = new TPad("cell", "", 0, 0, 1, 1);
108 pad->SetFillStyle(0);
109 pad->SetFrameFillStyle(0);
110 pad->Draw();
111 pad->cd();
112 }
113 gPad->Range(x0 - (x1 - x0) * (lm / (1. - rm - lm)),
114 y0 - (y1 - y0) * (bm / (1. - tm - lm)),
115 x1 + (x1 - x0) * (rm / (1. - rm - lm)),
116 y1 + (y1 - y0) * (tm / (1. - tm - lm)));
117 }
118
119 if (m_nebem) return PlotNeBem(twod);
120
121 // Get the cell type.
122 const std::string cellType = m_component->GetCellType();
123
124 // Get the x/y periodicities.
125 double sx = 0., sy = 0.;
126 const bool perX = m_component->GetPeriodicityX(sx);
127 const bool perY = m_component->GetPeriodicityY(sy);
128 // Determine the number of periods present in the cell.
129 const int nMinX = perX ? int(x0 / sx) - 1 : 0;
130 const int nMaxX = perX ? int(x1 / sx) + 1 : 0;
131 const int nMinY = perY ? int(y0 / sy) - 1 : 0;
132 const int nMaxY = perY ? int(y1 / sy) + 1 : 0;
133 // Get the phi periodicity.
134 double sphi = 0.;
135 const bool perPhi = m_component->GetPeriodicityPhi(sphi);
136 const int nPhi = perPhi ? int(360. / sphi) : 0;
137 sphi *= DegreeToRad;
138 if (!twod) SetupGeo(dx, dy, dz);
139 const bool polar = m_component->IsPolar();
140
141 // Get the number of wires.
142 const unsigned int nWires = m_component->GetNumberOfWires();
143 std::vector<std::string> wireTypes;
144 // Loop over the wires.
145 for (unsigned int i = 0; i < nWires; ++i) {
146 double xw = 0., yw = 0., dw = 0., vw = 0., lw = 0., qw = 0.;
147 std::string lbl;
148 int nTrap;
149 m_component->GetWire(i, xw, yw, dw, vw, lbl, lw, qw, nTrap);
150 auto it = std::find(wireTypes.begin(), wireTypes.end(), lbl);
151 const int type = std::distance(wireTypes.begin(), it);
152 if (it == wireTypes.end()) wireTypes.push_back(lbl);
153 if (polar) {
154 const double r = xw;
155 for (int j = 0; j <= nPhi; ++j) {
156 const double phi = yw * Garfield::DegreeToRad + j * sphi;
157 const double x = r * cos(phi);
158 const double y = r * sin(phi);
159 if (twod) {
160 PlotWire(x, y, dw, type);
161 } else {
162 PlotWire(x, y, dw, type, std::min(0.5 * lw, dz));
163 }
164 }
165 continue;
166 }
167 for (int nx = nMinX; nx <= nMaxX; ++nx) {
168 const double x = xw + nx * sx;
169 if (x + 0.5 * dw <= x0 || x - 0.5 * dw >= x1) continue;
170 for (int ny = nMinY; ny <= nMaxY; ++ny) {
171 const double y = yw + ny * sy;
172 if (y + 0.5 * dw <= y0 || y - 0.5 * dw >= y1) continue;
173 if (twod) {
174 PlotWire(x, y, dw, type);
175 } else {
176 PlotWire(x, y, dw, type, std::min(0.5 * lw, dz));
177 }
178 }
179 }
180 }
181
182 // Draw the x planes.
183 const unsigned int nPlanesX = m_component->GetNumberOfPlanesX();
184 for (unsigned int i = 0; i < nPlanesX; ++i) {
185 double xp = 0., vp = 0.;
186 std::string lbl;
187 m_component->GetPlaneX(i, xp, vp, lbl);
188 for (int nx = nMinX; nx <= nMaxX; ++nx) {
189 const double x = xp + nx * sx;
190 if (x < x0 || x > x1) continue;
191 if (twod) {
192 PlotPlane(x, y0, x, y1);
193 } else {
194 const double width = std::min(0.01 * dx, 0.01 * dy);
195 PlotPlane(width, dy, dz, x, 0);
196 }
197 }
198 }
199
200 // Draw the y planes.
201 const unsigned int nPlanesY = m_component->GetNumberOfPlanesY();
202 for (unsigned int i = 0; i < nPlanesY; ++i) {
203 double yp = 0., vp = 0.;
204 std::string lbl;
205 m_component->GetPlaneY(i, yp, vp, lbl);
206 for (int ny = nMinY; ny <= nMaxY; ++ny) {
207 const double y = yp + ny * sy;
208 if (y < y0 || y > y1) continue;
209 if (twod) {
210 PlotPlane(x0, y, x1, y);
211 } else {
212 const double width = std::min(0.01 * dx, 0.01 * dy);
213 PlotPlane(dx, width, dz, 0, y);
214 }
215 }
216 }
217
218 // Draw the r and phi planes.
219 const unsigned int nPlanesR = m_component->GetNumberOfPlanesR();
220 const unsigned int nPlanesPhi = m_component->GetNumberOfPlanesPhi();
221 if (nPlanesR > 0 || nPlanesPhi > 0) {
222 double vp = 0.;
223 std::string lbl;
224 std::array<double, 2> phi = {0., 360.};
225 double dphi = 360.;
226 if (nPlanesPhi == 2 && !m_component->GetPeriodicityPhi(dphi)) {
227 m_component->GetPlanePhi(0, phi[0], vp, lbl);
228 m_component->GetPlanePhi(1, phi[1], vp, lbl);
229 }
230 double w = 0.01 * std::min(dx, dy);
231 std::array<double, 2> r = {0., 0.};
232 if (nPlanesR > 0) {
233 m_component->GetPlaneR(0, r[0], vp, lbl);
234 w = 0.02 * r[0];
235 }
236 if (nPlanesR == 2) {
237 m_component->GetPlaneR(1, r[1], vp, lbl);
238 w = 0.01 * (r[1] - r[0]);
239 }
240 for (unsigned int i = 0; i < nPlanesR; ++i) {
241 if (twod) {
242 TEllipse circle;
243 circle.SetDrawOption("same");
244 circle.SetFillStyle(0);
245 circle.SetNoEdges();
246 circle.DrawEllipse(0, 0, r[i], r[i], phi[0], phi[1], 0);
247 continue;
248 }
249 const auto med = m_geo->GetMedium("Metal");
250 const auto tubs = m_geo->MakeTubs("Plane", med, r[i] - w, r[i] + w,
251 dz, phi[0], phi[1]);
252 tubs->SetLineColor(kGreen - 5);
253 tubs->SetTransparency(75);
254 m_geo->GetTopVolume()->AddNode(tubs, 1);
255 }
256 if (nPlanesR == 1) {
257 std::swap(r[0], r[1]);
258 } else if (nPlanesR == 0) {
259 r[1] = std::max(dx, dy);
260 }
261 for (unsigned int i = 0; i < nPlanesPhi; ++i) {
262 const double cp = cos(phi[i] * DegreeToRad);
263 const double sp = sin(phi[i] * DegreeToRad);
264 if (twod) {
265 PlotPlane(r[0] * cp, r[0] * sp, r[1] * cp, r[1] * sp);
266 continue;
267 }
268 const auto med = m_geo->GetMedium("Metal");
269 const double dr = 0.5 * (r[1] - r[0]);
270 TGeoVolume* plane = m_geo->MakeBox("Plane", med, dr, w, dz);
271 plane->SetLineColor(kGreen - 5);
272 plane->SetTransparency(75);
273 TGeoRotation* rot = new TGeoRotation("", phi[i], 0, 0);
274 const double rm = 0.5 * (r[0] + r[1]);
275 TGeoCombiTrans* tr = new TGeoCombiTrans(cp * rm, sp * rm, 0, rot);
276 m_geo->GetTopVolume()->AddNode(plane, 1, tr);
277 }
278 }
279 double rt = 0., vt = 0.;
280 int nt = 0;
281 std::string lbl;
282 if (m_component->GetTube(rt, vt, nt, lbl)) {
283 if (twod) {
284 PlotTube(0., 0., rt, nt);
285 } else {
286 PlotTube(0., 0., 0.98 * rt, 1.02 * rt, nt, dz);
287 }
288 }
289
290 if (!twod) {
291 m_geo->CloseGeometry();
292 m_geo->GetTopNode()->Draw("ogl");
293 } else {
294 gPad->Update();
295 }
296
297 return true;
298}
299
300bool ViewCell::PlotNeBem(const bool twod) {
301
302 if (!twod) {
303 std::cerr << m_className << "::PlotNeBem: 3D plot not implemented yet.\n";
304 return false;
305 }
306
307 // Draw the regions.
308 const unsigned int nRegions = m_nebem->GetNumberOfRegions();
309 for (unsigned int i = nRegions; i-- > 0;) {
310 std::vector<double> xv;
311 std::vector<double> yv;
312 Medium* medium = nullptr;
313 unsigned int bctype = 1;
314 double v = 0.;
315 if (!m_nebem->GetRegion(i, xv, yv, medium, bctype, v)) continue;
316 const unsigned int n = xv.size();
317 if (n < 3) continue;
318 TLine line;
319 line.SetDrawOption("same");
320 if (bctype == 4) {
321 line.SetLineStyle(2);
322 } else {
323 line.SetLineStyle(1);
324 }
325 for (unsigned int j = 0; j < n; ++j) {
326 const unsigned int k = j < n - 1 ? j + 1 : 0;
327 line.DrawLine(xv[j], yv[j], xv[k], yv[k]);
328 }
329 }
330
331 // Draw the wires.
332 const unsigned int nWires = m_nebem->GetNumberOfWires();
333 for (unsigned int i = 0; i < nWires; ++i) {
334 double x = 0., y = 0., d = 0., v = 0., q = 0.;
335 if (!m_nebem->GetWire(i, x, y, d, v, q)) continue;
336 PlotWire(x, y, d, 0);
337 }
338
339 // Draw the straight-line segments.
340 const unsigned int nSegments = m_nebem->GetNumberOfSegments();
341 for (unsigned int i = 0; i < nSegments; ++i) {
342 double x0 = 0., y0 = 0., x1 = 0., y1 = 0., v = 0.;
343 if (!m_nebem->GetSegment(i, x0, y0, x1, y1, v)) continue;
344 PlotPlane(x0, y0, x1, y1);
345 }
346
347 gPad->Update();
348 return true;
349}
350
351void ViewCell::SetupGeo(const double dx, const double dy, const double dz) {
352
353 if (!m_geo) {
354 gGeoManager = nullptr;
355 m_geo.reset(new TGeoManager("ViewCellGeoManager", "Cell layout"));
356 TGeoMaterial* matVacuum = new TGeoMaterial("Vacuum", 0., 0., 0.);
357 TGeoMaterial* matMetal = new TGeoMaterial("Metal", 63.546, 29., 8.92);
358 TGeoMedium* medVacuum = new TGeoMedium("Vacuum", 0, matVacuum);
359 TGeoMedium* medMetal = new TGeoMedium("Metal", 1, matMetal);
360 m_geo->AddMaterial(matVacuum);
361 m_geo->AddMaterial(medMetal->GetMaterial());
362 TGeoVolume* world =
363 m_geo->MakeBox("World", medVacuum, 1.05 * dx, 1.05 * dy, 1.05 * dz);
364 m_geo->SetTopVolume(world);
365 } else {
366 TGeoVolume* top = m_geo->GetTopVolume();
367 TGeoBBox* box = dynamic_cast<TGeoBBox*>(top);
368 double halfLenghts[3] = {1.05 * dx, 1.05 * dy, 1.05 * dz};
369 if (box) box->SetDimensions(halfLenghts);
370 }
371}
372
373void ViewCell::PlotWire(const double x, const double y, const double d,
374 const int type) {
375 if (m_useWireMarker) {
376 int markerStyle = 1;
377 if (type == 0) {
378 markerStyle = kFullCircle;
379 } else if (type == 1) {
380 markerStyle = kOpenCircle;
381 } else if (type == 2) {
382 markerStyle = kFullSquare;
383 } else if (type == 3) {
384 markerStyle = kOpenSquare;
385 } else {
386 markerStyle = 26 + type;
387 }
388 TMarker marker;
389 marker.SetMarkerStyle(markerStyle);
390 marker.SetDrawOption("Psame");
391 marker.DrawMarker(x, y);
392 return;
393 }
394
395 TEllipse circle;
396 circle.SetDrawOption("same");
397 circle.SetFillColor(kWhite);
398 circle.DrawEllipse(x, y, 0.5 * d, 0.5 * d, 0, 360, 0);
399}
400
401void ViewCell::PlotWire(const double x, const double y, const double d,
402 const int type, const double lz) {
403
404 const auto medium = m_geo->GetMedium("Metal");
405 TGeoVolume* wire = m_geo->MakeTube("Wire", medium, 0., 0.5 * d, lz);
406 switch (type) {
407 case 0:
408 wire->SetLineColor(kGray + 2);
409 break;
410 case 1:
411 wire->SetLineColor(kRed + 2);
412 break;
413 case 2:
414 wire->SetLineColor(kPink + 3);
415 break;
416 case 3:
417 wire->SetLineColor(kCyan + 3);
418 break;
419 default:
420 wire->SetLineColor(kBlue + type);
421 break;
422 }
423 m_geo->GetTopVolume()->AddNode(wire, 1, new TGeoTranslation(x, y, 0.));
424}
425
426void ViewCell::PlotTube(const double x0, const double y0, const double r,
427 const int n) {
428 if (n <= 0) {
429 TEllipse circle;
430 circle.SetDrawOption("same");
431 circle.SetFillStyle(0);
432 circle.DrawEllipse(x0, y0, r, r, 0, 360, 0);
433 return;
434 }
435
436 std::vector<double> x(n + 1);
437 std::vector<double> y(n + 1);
438 for (int i = 0; i <= n; ++i) {
439 const double phi = i * TwoPi / double(n);
440 x[i] = x0 + r * cos(phi);
441 y[i] = y0 + r * sin(phi);
442 }
443 TPolyLine pline;
444 pline.SetDrawOption("same");
445 pline.DrawPolyLine(n + 1, x.data(), y.data());
446}
447
448void ViewCell::PlotTube(const double x0, const double y0,
449 const double r1, const double r2, const int n,
450 const double lz) {
451
452 TGeoVolume* tube = nullptr;
453 if (n <= 0) {
454 // Round tube.
455 tube = m_geo->MakeTube("Tube", m_geo->GetMedium("Metal"), r1, r2, lz);
456 } else {
457 tube = m_geo->MakePgon("Tube", m_geo->GetMedium("Metal"), 0., 360., n, 2);
458 TGeoPgon* pgon = dynamic_cast<TGeoPgon*>(tube->GetShape());
459 if (pgon) {
460 pgon->DefineSection(0, -lz, r1, r2);
461 pgon->DefineSection(1, +lz, r1, r2);
462 }
463 }
464 tube->SetLineColor(kGreen + 2);
465 tube->SetTransparency(75);
466 m_geo->GetTopVolume()->AddNode(tube, 1, new TGeoTranslation(x0, y0, 0));
467}
468
469void ViewCell::PlotPlane(const double x0, const double y0,
470 const double x1, const double y1) {
471
472 TLine line;
473 line.SetDrawOption("same");
474 line.DrawLine(x0, y0, x1, y1);
475}
476
477void ViewCell::PlotPlane(const double dx, const double dy, const double dz,
478 const double x0, const double y0) {
479
480 const auto medium = m_geo->GetMedium("Metal");
481 TGeoVolume* plane = m_geo->MakeBox("Plane", medium, dx, dy, dz);
482 plane->SetLineColor(kGreen - 5);
483 plane->SetTransparency(75);
484 m_geo->GetTopVolume()->AddNode(plane, 1, new TGeoTranslation(x0, y0, 0));
485}
486
487}
bool GetVoltageRange(double &pmin, double &pmax) override
Calculate the voltage range [V].
Two-dimensional implementation of the nearly exact Boundary Element Method.
std::string m_className
Definition ViewBase.hh:82
TPad * GetCanvas()
Retrieve the canvas.
Definition ViewBase.cc:74
static bool RangeSet(TVirtualPad *)
Definition ViewBase.cc:83
ViewBase()=delete
Default constructor.
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