Garfield++ 3.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ComponentElmer.cc
Go to the documentation of this file.
1// Copied and modified ComponentAnsys123.cc
2
3#include <math.h>
4#include <stdlib.h>
5#include <fstream>
6#include <iostream>
7
9
10namespace {
11
12void PrintErrorReadingFile(const std::string& hdr, const std::string& file,
13 const int line) {
14 std::cerr << hdr << "\n Error reading file " << file << " (line " << line
15 << ").\n";
16}
17
18void PrintErrorOpeningFile(const std::string& hdr, const std::string& filetype,
19 const std::string& filename) {
20 std::cerr << hdr << "\n Could not open " << filetype << " file "
21 << filename << " for reading.\n";
22 std::cerr << " The file perhaps does not exist.\n";
23}
24}
25
26namespace Garfield {
27
29 m_className = "ComponentElmer";
30}
31
32ComponentElmer::ComponentElmer(const std::string& header,
33 const std::string& elist,
34 const std::string& nlist,
35 const std::string& mplist,
36 const std::string& volt, const std::string& unit)
38 m_className = "ComponentElmer";
39 Initialise(header, elist, nlist, mplist, volt, unit);
40}
41
42bool ComponentElmer::Initialise(const std::string& header,
43 const std::string& elist,
44 const std::string& nlist,
45 const std::string& mplist,
46 const std::string& volt,
47 const std::string& unit) {
48 const std::string hdr = m_className + "::Initialise:";
49 m_debug = false;
50 m_ready = false;
51 m_warning = false;
52 m_nWarnings = 0;
53
54 // Keep track of the success.
55 bool ok = true;
56
57 // Buffer for reading
58 const int size = 100;
59 char line[size];
60
61 // Open the header.
62 std::ifstream fheader;
63 fheader.open(header.c_str(), std::ios::in);
64 if (fheader.fail()) {
65 PrintErrorOpeningFile(hdr, "header", header);
66 }
67
68 // Temporary variables for use in file reading
69 char* token = NULL;
70 bool readerror = false;
71 bool readstop = false;
72 int il = 0;
73
74 // Read the header to get the number of nodes and elements.
75 fheader.getline(line, size, '\n');
76 token = strtok(line, " ");
77 nNodes = ReadInteger(token, 0, readerror);
78 token = strtok(NULL, " ");
79 nElements = ReadInteger(token, 0, readerror);
80 std::cout << hdr << "\n Read " << nNodes << " nodes and " << nElements
81 << " elements from file " << header << ".\n";
82 if (readerror) {
83 PrintErrorReadingFile(hdr, header, il);
84 fheader.close();
85 return false;
86 }
87
88 // Close the header file.
89 fheader.close();
90
91 // Open the nodes list.
92 std::ifstream fnodes;
93 fnodes.open(nlist.c_str(), std::ios::in);
94 if (fnodes.fail()) {
95 PrintErrorOpeningFile(hdr, "nodes", nlist);
96 }
97
98 // Check the value of the unit.
99 double funit;
100 if (strcmp(unit.c_str(), "mum") == 0 || strcmp(unit.c_str(), "micron") == 0 ||
101 strcmp(unit.c_str(), "micrometer") == 0) {
102 funit = 0.0001;
103 } else if (strcmp(unit.c_str(), "mm") == 0 ||
104 strcmp(unit.c_str(), "millimeter") == 0) {
105 funit = 0.1;
106 } else if (strcmp(unit.c_str(), "cm") == 0 ||
107 strcmp(unit.c_str(), "centimeter") == 0) {
108 funit = 1.0;
109 } else if (strcmp(unit.c_str(), "m") == 0 ||
110 strcmp(unit.c_str(), "meter") == 0) {
111 funit = 100.0;
112 } else {
113 std::cerr << hdr << " Unknown length unit " << unit << ".\n";
114 ok = false;
115 funit = 1.0;
116 }
117 if (m_debug) std::cout << hdr << " Unit scaling factor = " << funit << ".\n";
118
119 // Read the nodes from the file.
120 for (il = 0; il < nNodes; il++) {
121 // Get a line from the nodes file.
122 fnodes.getline(line, size, '\n');
123
124 // Ignore the first two characters.
125 token = strtok(line, " ");
126 token = strtok(NULL, " ");
127
128 // Get the node coordinates.
129 token = strtok(NULL, " ");
130 double xnode = ReadDouble(token, -1, readerror);
131 token = strtok(NULL, " ");
132 double ynode = ReadDouble(token, -1, readerror);
133 token = strtok(NULL, " ");
134 double znode = ReadDouble(token, -1, readerror);
135 if (readerror) {
136 PrintErrorReadingFile(hdr, nlist, il);
137 fnodes.close();
138 return false;
139 }
140
141 // Set up and create a new node.
142 Node newNode;
143 newNode.w.clear();
144 newNode.x = xnode * funit;
145 newNode.y = ynode * funit;
146 newNode.z = znode * funit;
147 nodes.push_back(std::move(newNode));
148 }
149
150 // Close the nodes file.
151 fnodes.close();
152
153 // Open the potential file.
154 std::ifstream fvolt;
155 fvolt.open(volt.c_str(), std::ios::in);
156 if (fvolt.fail()) {
157 PrintErrorOpeningFile(hdr, "potentials", volt);
158 }
159
160 // Reset the line counter.
161 il = 1;
162
163 // Read past the header.
164 while (!readstop && fvolt.getline(line, size, '\n')) {
165 token = strtok(line, " ");
166 if (strcmp(token, "Perm:") == 0) readstop = true;
167 il++;
168 }
169
170 // Should have stopped: if not, print error message.
171 if (!readstop) {
172 std::cerr << hdr << "\n Error reading past header of potentials file "
173 << volt << ".\n";
174 fvolt.close();
175 return false;
176 }
177
178 // Read past the permutation map (number of lines = nNodes).
179 for (int tl = 0; tl < nNodes; tl++) {
180 fvolt.getline(line, size, '\n');
181 il++;
182 }
183
184 // Read the potentials.
185 for (int tl = 0; tl < nNodes; tl++) {
186 fvolt.getline(line, size, '\n');
187 token = strtok(line, " ");
188 double v = ReadDouble(token, -1, readerror);
189 if (readerror) {
190 PrintErrorReadingFile(hdr, volt, il);
191 fvolt.close();
192 return false;
193 }
194 // Place the voltage in its appropriate node.
195 nodes[tl].v = v;
196 }
197
198 // Close the potentials file.
199 fvolt.close();
200
201 // Open the materials file.
202 std::ifstream fmplist;
203 fmplist.open(mplist.c_str(), std::ios::in);
204 if (fmplist.fail()) {
205 PrintErrorOpeningFile(hdr, "materials", mplist);
206 }
207
208 // Read the dielectric constants from the materials file.
209 fmplist.getline(line, size, '\n');
210 token = strtok(line, " ");
211 if (readerror) {
212 std::cerr << hdr << "\n Error reading number of materials from "
213 << mplist << ".\n";
214 fmplist.close();
215 ok = false;
216 return false;
217 }
218 m_nMaterials = ReadInteger(token, 0, readerror);
219 materials.resize(m_nMaterials);
220 for (unsigned int i = 0; i < m_nMaterials; ++i) {
221 materials[i].ohm = -1;
222 materials[i].eps = -1;
223 materials[i].medium = nullptr;
224 }
225 for (il = 2; il < ((int)m_nMaterials + 2); il++) {
226 fmplist.getline(line, size, '\n');
227 token = strtok(line, " ");
228 ReadInteger(token, -1, readerror);
229 token = strtok(NULL, " ");
230 double dc = ReadDouble(token, -1.0, readerror);
231 if (readerror) {
232 PrintErrorReadingFile(hdr, mplist, il);
233 fmplist.close();
234 ok = false;
235 return false;
236 }
237 materials[il - 2].eps = dc;
238 std::cout << hdr << "\n Set material " << il - 2 << " of "
239 << m_nMaterials << " to eps " << dc << ".\n";
240 }
241
242 // Close the materials file.
243 fmplist.close();
244
245 // Find the lowest epsilon, check for eps = 0, set default drift media.
246 double epsmin = -1.;
247 unsigned int iepsmin = 0;
248 for (unsigned int imat = 0; imat < m_nMaterials; ++imat) {
249 if (materials[imat].eps < 0) continue;
250 if (materials[imat].eps == 0) {
251 std::cerr << hdr << "\n Material " << imat
252 << " has been assigned a permittivity equal to zero in\n "
253 << mplist << ".\n";
254 ok = false;
255 } else if (epsmin < 0. || epsmin > materials[imat].eps) {
256 epsmin = materials[imat].eps;
257 iepsmin = imat;
258 }
259 }
260
261 if (epsmin < 0.) {
262 std::cerr << hdr << "\n No material with positive permittivity found \n"
263 << " in material list " << mplist << ".\n";
264 ok = false;
265 } else {
266 for (unsigned int imat = 0; imat < m_nMaterials; ++imat) {
267 materials[imat].driftmedium = imat == iepsmin ? true : false;
268 }
269 }
270
271 // Open the elements file.
272 std::ifstream felems;
273 felems.open(elist.c_str(), std::ios::in);
274 if (felems.fail()) {
275 PrintErrorOpeningFile(hdr, "elements", elist);
276 }
277
278 // Read the elements and their material indices.
279 elements.clear();
280 int highestnode = 0;
281 Element newElement;
282 for (il = 0; il < nElements; il++) {
283 // Get a line
284 felems.getline(line, size, '\n');
285
286 // Split into tokens.
287 token = strtok(line, " ");
288 // Read the 2nd-order element
289 // Note: Ordering of Elmer elements can be described in the
290 // ElmerSolver manual (appendix D. at the time of this comment)
291 // If the order read below is compared to the shape functions used
292 // eg. in ElectricField, the order is wrong, but note at the
293 // end of this function the order of elements 5,6,7 will change to
294 // 7,5,6 when actually recorded in newElement.emap to correct for this
295 token = strtok(NULL, " ");
296 int imat = ReadInteger(token, -1, readerror) - 1;
297 token = strtok(NULL, " ");
298 token = strtok(NULL, " ");
299 int in0 = ReadInteger(token, -1, readerror);
300 token = strtok(NULL, " ");
301 int in1 = ReadInteger(token, -1, readerror);
302 token = strtok(NULL, " ");
303 int in2 = ReadInteger(token, -1, readerror);
304 token = strtok(NULL, " ");
305 int in3 = ReadInteger(token, -1, readerror);
306 token = strtok(NULL, " ");
307 int in4 = ReadInteger(token, -1, readerror);
308 token = strtok(NULL, " ");
309 int in5 = ReadInteger(token, -1, readerror);
310 token = strtok(NULL, " ");
311 int in6 = ReadInteger(token, -1, readerror);
312 token = strtok(NULL, " ");
313 int in7 = ReadInteger(token, -1, readerror);
314 token = strtok(NULL, " ");
315 int in8 = ReadInteger(token, -1, readerror);
316 token = strtok(NULL, " ");
317 int in9 = ReadInteger(token, -1, readerror);
318
319 if (m_debug && il < 10) {
320 std::cout << " Read nodes " << in0 << ", " << in1 << ", " << in2
321 << ", " << in3 << ", ... from element " << il + 1 << " of "
322 << nElements << " with mat " << imat << ".\n";
323 }
324
325 // Check synchronisation.
326 if (readerror) {
327 PrintErrorReadingFile(hdr, elist, il);
328 felems.close();
329 ok = false;
330 return false;
331 }
332
333 // Check the material number and ensure that epsilon is non-negative.
334 if (imat < 0 || imat > (int)m_nMaterials) {
335 std::cerr << hdr << "\n Out-of-range material number on file " << elist
336 << " (line " << il << ").\n";
337 std::cerr << " Element: " << il << ", material: " << imat << "\n";
338 std::cerr << " nodes: (" << in0 << ", " << in1 << ", " << in2 << ", "
339 << in3 << ", " << in4 << ", " << in5 << ", " << in6 << ", "
340 << in7 << ", " << in8 << ", " << in9 << ")\n";
341 ok = false;
342 }
343 if (materials[imat].eps < 0) {
344 std::cerr << hdr << "\n Element " << il << " in element list " << elist
345 << "\n uses material " << imat
346 << " which has not been assigned a positive permittivity in "
347 << mplist << ".\n";
348 ok = false;
349 }
350
351 // Check the node numbers.
352 if (in0 < 1 || in1 < 1 || in2 < 1 || in3 < 1 || in4 < 1 || in5 < 1 ||
353 in6 < 1 || in7 < 1 || in8 < 1 || in9 < 1) {
354 std::cerr << hdr << "\n Found a node number < 1 on file " << elist
355 << " (line " << il << ").\n Element: " << il
356 << ", material: " << imat << "\n nodes: (" << in0 << ", "
357 << in1 << ", " << in2 << ", " << in3 << ", " << in4 << ", "
358 << in5 << ", " << in6 << ", " << in7 << ", " << in8 << ", "
359 << in9 << ")\n";
360 ok = false;
361 }
362 if (in0 > highestnode) highestnode = in0;
363 if (in1 > highestnode) highestnode = in1;
364 if (in2 > highestnode) highestnode = in2;
365 if (in3 > highestnode) highestnode = in3;
366 if (in4 > highestnode) highestnode = in4;
367 if (in5 > highestnode) highestnode = in5;
368 if (in6 > highestnode) highestnode = in6;
369 if (in7 > highestnode) highestnode = in7;
370 if (in8 > highestnode) highestnode = in8;
371 if (in9 > highestnode) highestnode = in9;
372
373 // These elements must not be degenerate.
374 if (in0 == in1 || in0 == in2 || in0 == in3 || in0 == in4 || in0 == in5 ||
375 in0 == in6 || in0 == in7 || in0 == in8 || in0 == in9 || in1 == in2 ||
376 in1 == in3 || in1 == in4 || in1 == in5 || in1 == in6 || in1 == in7 ||
377 in1 == in8 || in1 == in9 || in2 == in3 || in2 == in4 || in2 == in5 ||
378 in2 == in6 || in2 == in7 || in2 == in8 || in2 == in9 || in3 == in4 ||
379 in3 == in5 || in3 == in6 || in3 == in7 || in3 == in8 || in3 == in9 ||
380 in4 == in5 || in4 == in6 || in4 == in7 || in4 == in8 || in4 == in9 ||
381 in5 == in6 || in5 == in7 || in5 == in8 || in5 == in9 || in6 == in7 ||
382 in6 == in8 || in6 == in9 || in7 == in8 || in7 == in9 || in8 == in9) {
383 std::cerr << hdr << "\n Element " << il << " of file " << elist
384 << " is degenerate,\n"
385 << " no such elements are allowed in this type of map.\n";
386 ok = false;
387 }
388
389 newElement.degenerate = false;
390
391 // Store the material reference.
392 newElement.matmap = imat;
393
394 // Node references
395 newElement.emap[0] = in0 - 1;
396 newElement.emap[1] = in1 - 1;
397 newElement.emap[2] = in2 - 1;
398 newElement.emap[3] = in3 - 1;
399 newElement.emap[4] = in4 - 1;
400 newElement.emap[7] = in5 - 1;
401 newElement.emap[5] = in6 - 1;
402 newElement.emap[6] = in7 - 1;
403 newElement.emap[8] = in8 - 1;
404 newElement.emap[9] = in9 - 1;
405 elements.push_back(newElement);
406 }
407
408 // Close the elements file.
409 felems.close();
410
411 // Set the ready flag.
412 if (ok) {
413 m_ready = true;
414 } else {
415 std::cerr << hdr << "\n Field map could not be "
416 << "read and cannot be interpolated.\n";
417 return false;
418 }
419
420 std::cout << hdr << " Finished.\n";
421
422 // Remove weighting fields (if any).
423 wfields.clear();
424 wfieldsOk.clear();
426
427 // Establish the ranges.
428 SetRange();
430 return true;
431}
432
433bool ComponentElmer::SetWeightingField(std::string wvolt, std::string label) {
434 const std::string hdr = m_className + "::SetWeightingField:";
435 if (!m_ready) {
436 PrintNotReady("SetWeightingField");
437 std::cerr << " Weighting field cannot be added.\n";
438 return false;
439 }
440
441 // Keep track of the success.
442 bool ok = true;
443
444 // Open the voltage list.
445 std::ifstream fwvolt;
446 fwvolt.open(wvolt.c_str(), std::ios::in);
447 if (fwvolt.fail()) {
448 PrintErrorOpeningFile(hdr, "potential", wvolt);
449 return false;
450 }
451
452 // Check if a weighting field with the same label already exists.
453 int iw = nWeightingFields;
454 for (int i = nWeightingFields; i--;) {
455 if (wfields[i] == label) {
456 iw = i;
457 break;
458 }
459 }
460 if (iw == nWeightingFields) {
464 for (int j = nNodes; j--;) {
465 nodes[j].w.resize(nWeightingFields);
466 }
467 } else {
468 std::cout << hdr << "\n Replacing existing weighting field " << label
469 << ".\n";
470 }
471 wfields[iw] = label;
472 wfieldsOk[iw] = false;
473
474 // Temporary variables for use in file reading
475 const int size = 100;
476 char line[size];
477 char* token = NULL;
478 bool readerror = false;
479 bool readstop = false;
480 int il = 1;
481
482 // Read past the header.
483 while (!readstop && fwvolt.getline(line, size, '\n')) {
484 token = strtok(line, " ");
485 if (strcmp(token, "Perm:") == 0) readstop = true;
486 il++;
487 }
488
489 // Should have stopped: if not, print error message.
490 if (!readstop) {
491 std::cerr << hdr << "\n Error reading past header of potentials file "
492 << wvolt << ".\n";
493 fwvolt.close();
494 ok = false;
495 return false;
496 }
497
498 // Read past the permutation map (number of lines = nNodes).
499 for (int tl = 0; tl < nNodes; tl++) {
500 fwvolt.getline(line, size, '\n');
501 il++;
502 }
503
504 // Read the potentials.
505 for (int tl = 0; tl < nNodes; tl++) {
506 double v;
507 fwvolt.getline(line, size, '\n');
508 token = strtok(line, " ");
509 v = ReadDouble(token, -1, readerror);
510 if (readerror) {
511 PrintErrorReadingFile(hdr, wvolt, il);
512 fwvolt.close();
513 ok = false;
514 return false;
515 }
516 // Place the weighting potential at its appropriate node and index.
517 nodes[tl].w[iw] = v;
518 }
519
520 // Close the potentials file.
521 fwvolt.close();
522 std::cout << hdr << "\n Read potentials from file " << wvolt << ".\n";
523
524 // Set the ready flag.
525 wfieldsOk[iw] = ok;
526 if (!ok) {
527 std::cerr << hdr << "\n Field map could not "
528 << "be read and cannot be interpolated.\n";
529 return false;
530 }
531 return true;
532}
533
534void ComponentElmer::ElectricField(const double x, const double y,
535 const double z, double& ex, double& ey,
536 double& ez, Medium*& m, int& status) {
537 double v = 0.;
538 ElectricField(x, y, z, ex, ey, ez, v, m, status);
539}
540
541void ComponentElmer::ElectricField(const double xin, const double yin,
542 const double zin, double& ex, double& ey,
543 double& ez, double& volt, Medium*& m,
544 int& status) {
545 // Copy the coordinates
546 double x = xin, y = yin, z = zin;
547
548 // Map the coordinates onto field map coordinates
549 bool xmirr, ymirr, zmirr;
550 double rcoordinate, rotation;
551 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
552
553 // Initial values
554 ex = ey = ez = volt = 0.;
555 status = 0;
556 m = nullptr;
557
558 // Do not proceed if not properly initialised.
559 if (!m_ready) {
560 status = -10;
561 PrintNotReady("ElectricField");
562 return;
563 }
564
565 if (m_warning) PrintWarning("ElectricField");
566
567 // Find the element that contains this point
568 double t1, t2, t3, t4, jac[4][4], det;
569 const int imap = FindElement13(x, y, z, t1, t2, t3, t4, jac, det);
570 if (imap < 0) {
571 if (m_debug) {
572 std::cout << m_className << "::ElectricField:\n Point (" << x << ", "
573 << y << ", " << z << ") is not in the mesh.\n";
574 }
575 status = -6;
576 return;
577 }
578
579 const Element& element = elements[imap];
580 if (m_debug) {
581 PrintElement("ElectricField", x, y, z, t1, t2, t3, t4, element, 10);
582 }
583 const Node& n0 = nodes[element.emap[0]];
584 const Node& n1 = nodes[element.emap[1]];
585 const Node& n2 = nodes[element.emap[2]];
586 const Node& n3 = nodes[element.emap[3]];
587 const Node& n4 = nodes[element.emap[4]];
588 const Node& n5 = nodes[element.emap[5]];
589 const Node& n6 = nodes[element.emap[6]];
590 const Node& n7 = nodes[element.emap[7]];
591 const Node& n8 = nodes[element.emap[8]];
592 const Node& n9 = nodes[element.emap[9]];
593 // Shorthands.
594 const double fourt1 = 4 * t1;
595 const double fourt2 = 4 * t2;
596 const double fourt3 = 4 * t3;
597 const double fourt4 = 4 * t4;
598 const double invdet = 1. / det;
599 // Tetrahedral field
600 volt = n0.v * t1 * (2 * t1 - 1) + n1.v * t2 * (2 * t2 - 1) +
601 n2.v * t3 * (2 * t3 - 1) + n3.v * t4 * (2 * t4 - 1) +
602 4 * n4.v * t1 * t2 + 4 * n5.v * t1 * t3 + 4 * n6.v * t1 * t4 +
603 4 * n7.v * t2 * t3 + 4 * n8.v * t2 * t4 + 4 * n9.v * t3 * t4;
604 ex = -(n0.v * (fourt1 - 1) * jac[0][1] + n1.v * (fourt2 - 1) * jac[1][1] +
605 n2.v * (fourt3 - 1) * jac[2][1] + n3.v * (fourt4 - 1) * jac[3][1] +
606 n4.v * (fourt2 * jac[0][1] + fourt1 * jac[1][1]) +
607 n5.v * (fourt3 * jac[0][1] + fourt1 * jac[2][1]) +
608 n6.v * (fourt4 * jac[0][1] + fourt1 * jac[3][1]) +
609 n7.v * (fourt3 * jac[1][1] + fourt2 * jac[2][1]) +
610 n8.v * (fourt4 * jac[1][1] + fourt2 * jac[3][1]) +
611 n9.v * (fourt4 * jac[2][1] + fourt3 * jac[3][1])) *
612 invdet;
613 ey = -(n0.v * (fourt1 - 1) * jac[0][2] + n1.v * (fourt2 - 1) * jac[1][2] +
614 n2.v * (fourt3 - 1) * jac[2][2] + n3.v * (fourt4 - 1) * jac[3][2] +
615 n4.v * (fourt2 * jac[0][2] + fourt1 * jac[1][2]) +
616 n5.v * (fourt3 * jac[0][2] + fourt1 * jac[2][2]) +
617 n6.v * (fourt4 * jac[0][2] + fourt1 * jac[3][2]) +
618 n7.v * (fourt3 * jac[1][2] + fourt2 * jac[2][2]) +
619 n8.v * (fourt4 * jac[1][2] + fourt2 * jac[3][2]) +
620 n9.v * (fourt4 * jac[2][2] + fourt3 * jac[3][2])) *
621 invdet;
622 ez = -(n0.v * (fourt1 - 1) * jac[0][3] + n1.v * (fourt2 - 1) * jac[1][3] +
623 n2.v * (fourt3 - 1) * jac[2][3] + n3.v * (fourt4 - 1) * jac[3][3] +
624 n4.v * (fourt2 * jac[0][3] + fourt1 * jac[1][3]) +
625 n5.v * (fourt3 * jac[0][3] + fourt1 * jac[2][3]) +
626 n6.v * (fourt4 * jac[0][3] + fourt1 * jac[3][3]) +
627 n7.v * (fourt3 * jac[1][3] + fourt2 * jac[2][3]) +
628 n8.v * (fourt4 * jac[1][3] + fourt2 * jac[3][3]) +
629 n9.v * (fourt4 * jac[2][3] + fourt3 * jac[3][3])) *
630 invdet;
631
632 // Transform field to global coordinates
633 UnmapFields(ex, ey, ez, x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
634
635 // Drift medium?
636 const Material& mat = materials[element.matmap];
637 if (m_debug) {
638 std::cout << m_className << "::ElectricField:\n Material "
639 << element.matmap << ", drift flag " << mat.driftmedium << ".\n";
640 }
641 m = mat.medium;
642 status = -5;
643 if (mat.driftmedium && m && m->IsDriftable()) status = 0;
644}
645
646void ComponentElmer::WeightingField(const double xin, const double yin,
647 const double zin, double& wx, double& wy,
648 double& wz, const std::string& label) {
649 // Initial values
650 wx = wy = wz = 0;
651
652 // Do not proceed if not properly initialised.
653 if (!m_ready) return;
654
655 // Look for the label.
656 int iw = 0;
657 bool found = false;
658 for (int i = nWeightingFields; i--;) {
659 if (wfields[i] == label) {
660 iw = i;
661 found = true;
662 break;
663 }
664 }
665
666 // Do not proceed if the requested weighting field does not exist.
667 if (!found) return;
668 // Check if the weighting field is properly initialised.
669 if (!wfieldsOk[iw]) return;
670
671 // Copy the coordinates.
672 double x = xin, y = yin, z = zin;
673
674 // Map the coordinates onto field map coordinates
675 bool xmirr, ymirr, zmirr;
676 double rcoordinate, rotation;
677 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
678
679 if (m_warning) PrintWarning("WeightingField");
680
681 // Find the element that contains this point.
682 double t1, t2, t3, t4, jac[4][4], det;
683 const int imap = FindElement13(x, y, z, t1, t2, t3, t4, jac, det);
684 // Check if the point is in the mesh.
685 if (imap < 0) return;
686
687 const Element& element = elements[imap];
688 if (m_debug) {
689 PrintElement("WeightingField", x, y, z, t1, t2, t3, t4, element, 10, iw);
690 }
691 const Node& n0 = nodes[element.emap[0]];
692 const Node& n1 = nodes[element.emap[1]];
693 const Node& n2 = nodes[element.emap[2]];
694 const Node& n3 = nodes[element.emap[3]];
695 const Node& n4 = nodes[element.emap[4]];
696 const Node& n5 = nodes[element.emap[5]];
697 const Node& n6 = nodes[element.emap[6]];
698 const Node& n7 = nodes[element.emap[7]];
699 const Node& n8 = nodes[element.emap[8]];
700 const Node& n9 = nodes[element.emap[9]];
701 // Shorthands.
702 const double fourt1 = 4 * t1;
703 const double fourt2 = 4 * t2;
704 const double fourt3 = 4 * t3;
705 const double fourt4 = 4 * t4;
706 const double invdet = 1. / det;
707 // Tetrahedral field
708 wx = -(n0.w[iw] * (fourt1 - 1) * jac[0][1] +
709 n1.w[iw] * (fourt2 - 1) * jac[1][1] +
710 n2.w[iw] * (fourt3 - 1) * jac[2][1] +
711 n3.w[iw] * (fourt4 - 1) * jac[3][1] +
712 n4.w[iw] * (fourt2 * jac[0][1] + fourt1 * jac[1][1]) +
713 n5.w[iw] * (fourt3 * jac[0][1] + fourt1 * jac[2][1]) +
714 n6.w[iw] * (fourt4 * jac[0][1] + fourt1 * jac[3][1]) +
715 n7.w[iw] * (fourt3 * jac[1][1] + fourt2 * jac[2][1]) +
716 n8.w[iw] * (fourt4 * jac[1][1] + fourt2 * jac[3][1]) +
717 n9.w[iw] * (fourt4 * jac[2][1] + fourt3 * jac[3][1])) *
718 invdet;
719 wy = -(n0.w[iw] * (fourt1 - 1) * jac[0][2] +
720 n1.w[iw] * (fourt2 - 1) * jac[1][2] +
721 n2.w[iw] * (fourt3 - 1) * jac[2][2] +
722 n3.w[iw] * (fourt4 - 1) * jac[3][2] +
723 n4.w[iw] * (fourt2 * jac[0][2] + fourt1 * jac[1][2]) +
724 n5.w[iw] * (fourt3 * jac[0][2] + fourt1 * jac[2][2]) +
725 n6.w[iw] * (fourt4 * jac[0][2] + fourt1 * jac[3][2]) +
726 n7.w[iw] * (fourt3 * jac[1][2] + fourt2 * jac[2][2]) +
727 n8.w[iw] * (fourt4 * jac[1][2] + fourt2 * jac[3][2]) +
728 n9.w[iw] * (fourt4 * jac[2][2] + fourt3 * jac[3][2])) *
729 invdet;
730 wz = -(n0.w[iw] * (fourt1 - 1) * jac[0][3] +
731 n1.w[iw] * (fourt2 - 1) * jac[1][3] +
732 n2.w[iw] * (fourt3 - 1) * jac[2][3] +
733 n3.w[iw] * (fourt4 - 1) * jac[3][3] +
734 n4.w[iw] * (fourt2 * jac[0][3] + fourt1 * jac[1][3]) +
735 n5.w[iw] * (fourt3 * jac[0][3] + fourt1 * jac[2][3]) +
736 n6.w[iw] * (fourt4 * jac[0][3] + fourt1 * jac[3][3]) +
737 n7.w[iw] * (fourt3 * jac[1][3] + fourt2 * jac[2][3]) +
738 n8.w[iw] * (fourt4 * jac[1][3] + fourt2 * jac[3][3]) +
739 n9.w[iw] * (fourt4 * jac[2][3] + fourt3 * jac[3][3])) *
740 invdet;
741
742 // Transform field to global coordinates
743 UnmapFields(wx, wy, wz, x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
744}
745
746double ComponentElmer::WeightingPotential(const double xin, const double yin,
747 const double zin,
748 const std::string& label) {
749 // Do not proceed if not properly initialised.
750 if (!m_ready) return 0.;
751
752 // Look for the label.
753 int iw = 0;
754 bool found = false;
755 for (int i = nWeightingFields; i--;) {
756 if (wfields[i] == label) {
757 iw = i;
758 found = true;
759 break;
760 }
761 }
762
763 // Do not proceed if the requested weighting field does not exist.
764 if (!found) return 0.;
765 // Check if the weighting field is properly initialised.
766 if (!wfieldsOk[iw]) return 0.;
767
768 // Copy the coordinates.
769 double x = xin, y = yin, z = zin;
770
771 // Map the coordinates onto field map coordinates.
772 bool xmirr, ymirr, zmirr;
773 double rcoordinate, rotation;
774 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
775
776 if (m_warning) PrintWarning("WeightingPotential");
777
778 // Find the element that contains this point.
779 double t1, t2, t3, t4, jac[4][4], det;
780 const int imap = FindElement13(x, y, z, t1, t2, t3, t4, jac, det);
781 if (imap < 0) return 0.;
782
783 const Element& element = elements[imap];
784 if (m_debug) {
785 PrintElement("WeightingPotential", x, y, z, t1, t2, t3, t4, element, 10,
786 iw);
787 }
788 const Node& n0 = nodes[element.emap[0]];
789 const Node& n1 = nodes[element.emap[1]];
790 const Node& n2 = nodes[element.emap[2]];
791 const Node& n3 = nodes[element.emap[3]];
792 const Node& n4 = nodes[element.emap[4]];
793 const Node& n5 = nodes[element.emap[5]];
794 const Node& n6 = nodes[element.emap[6]];
795 const Node& n7 = nodes[element.emap[7]];
796 const Node& n8 = nodes[element.emap[8]];
797 const Node& n9 = nodes[element.emap[9]];
798 // Tetrahedral field
799 return n0.w[iw] * t1 * (2 * t1 - 1) + n1.w[iw] * t2 * (2 * t2 - 1) +
800 n2.w[iw] * t3 * (2 * t3 - 1) + n3.w[iw] * t4 * (2 * t4 - 1) +
801 4 * n4.w[iw] * t1 * t2 + 4 * n5.w[iw] * t1 * t3 +
802 4 * n6.w[iw] * t1 * t4 + 4 * n7.w[iw] * t2 * t3 +
803 4 * n8.w[iw] * t2 * t4 + 4 * n9.w[iw] * t3 * t4;
804}
805
806Medium* ComponentElmer::GetMedium(const double xin, const double yin,
807 const double zin) {
808 // Copy the coordinates
809 double x = xin, y = yin, z = zin;
810
811 // Map the coordinates onto field map coordinates
812 bool xmirr, ymirr, zmirr;
813 double rcoordinate, rotation;
814 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
815
816 // Do not proceed if not properly initialised.
817 if (!m_ready) {
818 PrintNotReady("GetMedium");
819 return nullptr;
820 }
821 if (m_warning) PrintWarning("GetMedium");
822
823 // Find the element that contains this point
824 double t1, t2, t3, t4, jac[4][4], det;
825 const int imap = FindElement13(x, y, z, t1, t2, t3, t4, jac, det);
826 if (imap < 0) {
827 if (m_debug) {
828 std::cout << m_className << "::GetMedium:\n Point (" << x << ", " << y
829 << ", " << z << ") is not in the mesh.\n";
830 }
831 return nullptr;
832 }
833 const Element& element = elements[imap];
834 if (element.matmap >= m_nMaterials) {
835 if (m_debug) {
836 std::cerr << m_className << "::GetMedium:\n Point (" << x << ", " << y
837 << ", " << z << ") has out of range material number " << imap
838 << ".\n";
839 }
840 return nullptr;
841 }
842
843 if (m_debug) PrintElement("GetMedium", x, y, z, t1, t2, t3, t4, element, 10);
844
845 return materials[element.matmap].medium;
846}
847
848double ComponentElmer::GetElementVolume(const unsigned int i) {
849 if (i >= elements.size()) return 0.;
850 const Element& element = elements[i];
851 const Node& n0 = nodes[element.emap[0]];
852 const Node& n1 = nodes[element.emap[1]];
853 const Node& n2 = nodes[element.emap[2]];
854 const Node& n3 = nodes[element.emap[3]];
855
856 // Uses formula V = |a (dot) b x c|/6
857 // with a => "3", b => "1", c => "2" and origin "0"
858 const double vol =
859 fabs((n3.x - n0.x) *
860 ((n1.y - n0.y) * (n2.z - n0.z) - (n2.y - n0.y) * (n1.z - n0.z)) +
861 (n3.y - n0.y) *
862 ((n1.z - n0.z) * (n2.x - n0.x) - (n2.z - n0.z) * (n1.x - n0.x)) +
863 (n3.z - n0.z) * ((n1.x - n0.x) * (n2.y - n0.y) -
864 (n3.x - n0.x) * (n1.y - n0.y))) /
865 6.;
866 return vol;
867}
868
869void ComponentElmer::GetAspectRatio(const unsigned int i, double& dmin,
870 double& dmax) {
871 if (i >= elements.size()) {
872 dmin = dmax = 0.;
873 return;
874 }
875
876 const Element& element = elements[i];
877 const int np = 4;
878 // Loop over all pairs of vertices.
879 for (int j = 0; j < np - 1; ++j) {
880 const Node& nj = nodes[element.emap[j]];
881 for (int k = j + 1; k < np; ++k) {
882 const Node& nk = nodes[element.emap[k]];
883 // Compute distance.
884 const double dx = nj.x - nk.x;
885 const double dy = nj.y - nk.y;
886 const double dz = nj.z - nk.z;
887 const double dist = sqrt(dx * dx + dy * dy + dz * dz);
888 if (k == 1) {
889 dmin = dmax = dist;
890 } else {
891 if (dist < dmin) dmin = dist;
892 if (dist > dmax) dmax = dist;
893 }
894 }
895 }
896}
897
898} // namespace Garfield
bool m_ready
Ready for use?
std::string m_className
Class name.
bool m_debug
Switch on/off debugging messages.
void WeightingField(const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label) override
bool Initialise(const std::string &header="mesh.header", const std::string &elist="mesh.elements", const std::string &nlist="mesh.nodes", const std::string &mplist="dielectrics.dat", const std::string &volt="out.result", const std::string &unit="cm")
void UpdatePeriodicity() override
Verify periodicities.
bool SetWeightingField(std::string prnsol, std::string label)
Import a list of voltages to be used as weighting field.
ComponentElmer()
Default constructor.
void ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status) override
double GetElementVolume(const unsigned int i) override
void GetAspectRatio(const unsigned int i, double &dmin, double &dmax) override
double WeightingPotential(const double x, const double y, const double z, const std::string &label) override
Medium * GetMedium(const double x, const double y, const double z) override
Get the medium at a given location (x, y, z).
Base class for components based on finite-element field maps.
double ReadDouble(char *token, double def, bool &error)
void MapCoordinates(double &xpos, double &ypos, double &zpos, bool &xmirrored, bool &ymirrored, bool &zmirrored, double &rcoordinate, double &rotation) const
Move (xpos, ypos, zpos) to field map coordinates.
void PrintWarning(const std::string &header)
int ReadInteger(char *token, int def, bool &error)
std::vector< Material > materials
int FindElement13(const double x, const double y, const double z, double &t1, double &t2, double &t3, double &t4, double jac[4][4], double &det)
Find the element for a point in curved quadratic tetrahedra.
void PrintElement(const std::string &header, const double x, const double y, const double z, const double t1, const double t2, const double t3, const double t4, const Element &element, const unsigned int n, const int iw=-1) const
virtual void SetRange()
Calculate x, y, z, V and angular ranges.
std::vector< Element > elements
std::vector< std::string > wfields
void UnmapFields(double &ex, double &ey, double &ez, double &xpos, double &ypos, double &zpos, bool &xmirrored, bool &ymirrored, bool &zmirrored, double &rcoordinate, double &rotation) const
Move (ex, ey, ez) to global coordinates.
void PrintNotReady(const std::string &header) const
Abstract base class for media.
Definition: Medium.hh:13
bool IsDriftable() const
Is charge carrier transport enabled in this medium?
Definition: Medium.hh:74