Garfield++ 3.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ComponentAnsys121.cc
Go to the documentation of this file.
1#include <math.h>
2#include <stdlib.h>
3#include <fstream>
4#include <iostream>
5
7
8namespace Garfield {
9
11 m_className = "ComponentAnsys121";
12 m_is3d = false;
13 // Default bounding box
14 m_minBoundingBox[2] = -50;
15 m_maxBoundingBox[2] = 50;
16}
17
18bool ComponentAnsys121::Initialise(std::string elist, std::string nlist,
19 std::string mplist, std::string prnsol,
20 std::string unit) {
21 m_ready = false;
22 m_warning = false;
23 m_nWarnings = 0;
24 // Keep track of the success.
25 bool ok = true;
26
27 // Buffer for reading
28 const int size = 100;
29 char line[size];
30
31 // Open the material list.
32 std::ifstream fmplist;
33 fmplist.open(mplist.c_str(), std::ios::in);
34 if (fmplist.fail()) {
35 std::cerr << m_className << "::Initialise:\n"
36 << " Could not open material file " << mplist
37 << " for reading. The file perhaps does not exist.\n";
38 return false;
39 }
40
41 // Read the material list.
42 m_nMaterials = 0;
43 int il = 0;
44 unsigned int icurrmat = 0;
45 bool readerror = false;
46 while (fmplist.getline(line, size, '\n')) {
47 il++;
48 // Skip page feed
49 if (strcmp(line, "1") == 0) {
50 fmplist.getline(line, size, '\n');
51 il++;
52 fmplist.getline(line, size, '\n');
53 il++;
54 fmplist.getline(line, size, '\n');
55 il++;
56 fmplist.getline(line, size, '\n');
57 il++;
58 fmplist.getline(line, size, '\n');
59 il++;
60 continue;
61 }
62 // Split the line in tokens
63 char* token = NULL;
64 token = strtok(line, " ");
65 // Skip blank lines and headers
66 if (!token || strcmp(token, " ") == 0 || strcmp(token, "\n") == 0 ||
67 strcmp(token, "TEMPERATURE") == 0 || strcmp(token, "PROPERTY=") == 0 ||
68 int(token[0]) == 10 || int(token[0]) == 13)
69 continue;
70 // Read number of materials,
71 // ensure it does not exceed the maximum and initialise the list
72 if (strcmp(token, "LIST") == 0) {
73 token = strtok(NULL, " ");
74 token = strtok(NULL, " ");
75 token = strtok(NULL, " ");
76 token = strtok(NULL, " ");
77 m_nMaterials = ReadInteger(token, -1, readerror);
78 if (readerror) {
79 std::cerr << m_className << "::Initialise:\n";
80 std::cerr << " Error reading file " << mplist << " (line " << il
81 << ").\n";
82 fmplist.close();
83 ok = false;
84 return false;
85 }
86 materials.resize(m_nMaterials);
87 for (unsigned int i = 0; i < m_nMaterials; ++i) {
88 materials[i].ohm = -1;
89 materials[i].eps = -1;
90 materials[i].medium = nullptr;
91 }
92 if (m_debug) {
93 std::cout << m_className << "::Initialise:\n";
94 std::cout << " Number of materials: " << m_nMaterials << "\n";
95 }
96 } else if (strcmp(token, "MATERIAL") == 0) {
97 // Version 12 format: read material number
98 token = strtok(NULL, " ");
99 token = strtok(NULL, " ");
100 const int imat = ReadInteger(token, -1, readerror);
101 if (readerror || imat < 0) {
102 std::cerr << m_className << "::Initialise:\n";
103 std::cerr << " Error reading file " << mplist << " (line " << il
104 << ").\n";
105 fmplist.close();
106 ok = false;
107 return false;
108 }
109 icurrmat = imat;
110 } else if (strcmp(token, "TEMP") == 0) {
111 // Version 12 format: read property tag and value
112 token = strtok(NULL, " ");
113 int itype = 0;
114 if (strncmp(token, "PERX", 4) == 0) {
115 itype = 1;
116 } else if (strncmp(token, "RSVX", 4) == 0) {
117 itype = 2;
118 } else {
119 std::cerr << m_className << "::Initialise:\n";
120 std::cerr << " Found unknown material property flag " << token
121 << "\n";
122 std::cerr << " on material properties file " << mplist << " (line "
123 << il << ").\n";
124 ok = false;
125 }
126 fmplist.getline(line, size, '\n');
127 il++;
128 token = NULL;
129 token = strtok(line, " ");
130 if (icurrmat < 1 || icurrmat > m_nMaterials) {
131 std::cerr << m_className << "::Initialise:\n";
132 std::cerr << " Found out-of-range current material index "
133 << icurrmat << "\n";
134 std::cerr << " in material properties file " << mplist << ".\n";
135 ok = false;
136 readerror = false;
137 } else if (itype == 1) {
138 materials[icurrmat - 1].eps = ReadDouble(token, -1, readerror);
139 } else if (itype == 2) {
140 materials[icurrmat - 1].ohm = ReadDouble(token, -1, readerror);
141 }
142 if (readerror) {
143 std::cerr << m_className << "::Initialise:\n";
144 std::cerr << " Error reading file " << mplist << " (line " << il
145 << ").\n";
146 fmplist.close();
147 ok = false;
148 return false;
149 }
150 } else if (strcmp(token, "PROPERTY") == 0) {
151 // Version 11 format
152 token = strtok(NULL, " ");
153 token = strtok(NULL, " ");
154 int itype = 0;
155 if (strcmp(token, "PERX") == 0) {
156 itype = 1;
157 } else if (strcmp(token, "RSVX") == 0) {
158 itype = 2;
159 } else {
160 std::cerr << m_className << "::Initialise:\n";
161 std::cerr << m_className << "::Initialise:\n";
162 std::cerr << " Found unknown material property flag " << token
163 << "\n";
164 std::cerr << " on material properties file " << mplist << " (line "
165 << il << ").\n";
166 ok = false;
167 }
168 token = strtok(NULL, " ");
169 token = strtok(NULL, " ");
170 int imat = ReadInteger(token, -1, readerror);
171 if (readerror) {
172 std::cerr << m_className << "::Initialise:\n";
173 std::cerr << " Error reading file " << mplist << " (line " << il
174 << ").\n";
175 fmplist.close();
176 ok = false;
177 return false;
178 } else if (imat < 1 || imat > (int)m_nMaterials) {
179 std::cerr << m_className << "::Initialise:\n";
180 std::cerr << " Found out-of-range current material index " << imat
181 << "\n";
182 std::cerr << " in material properties file " << mplist << ".\n";
183 ok = false;
184 } else {
185 fmplist.getline(line, size, '\n');
186 il++;
187 fmplist.getline(line, size, '\n');
188 il++;
189 token = NULL;
190 token = strtok(line, " ");
191 token = strtok(NULL, " ");
192 if (itype == 1) {
193 materials[imat - 1].eps = ReadDouble(token, -1, readerror);
194 } else if (itype == 2) {
195 materials[imat - 1].ohm = ReadDouble(token, -1, readerror);
196 }
197 if (readerror) {
198 std::cerr << m_className << "::Initialise:\n";
199 std::cerr << " Error reading file " << mplist << " (line " << il
200 << ").\n";
201 fmplist.close();
202 ok = false;
203 return false;
204 }
205 }
206 }
207 }
208
209 // Close the file
210 fmplist.close();
211
212 // Find the lowest epsilon, check for eps = 0, set default drift media
213 double epsmin = -1;
214 unsigned int iepsmin = 0;
215 for (unsigned int imat = 0; imat < m_nMaterials; ++imat) {
216 if (materials[imat].eps < 0) continue;
217 if (materials[imat].eps == 0) {
218 std::cerr << m_className << "::Initialise:\n";
219 std::cerr << " Material " << imat
220 << " has been assigned a permittivity\n";
221 std::cerr << " equal to zero in " << mplist << ".\n";
222 ok = false;
223 } else if (epsmin < 0. || epsmin > materials[imat].eps) {
224 epsmin = materials[imat].eps;
225 iepsmin = imat;
226 }
227 }
228
229 if (epsmin < 0.) {
230 std::cerr << m_className << "::Initialise:\n";
231 std::cerr << " No material with positive permittivity found \n";
232 std::cerr << " in material list " << mplist << ".\n";
233 ok = false;
234 } else {
235 for (unsigned int imat = 0; imat < m_nMaterials; ++imat) {
236 if (imat == iepsmin) {
237 materials[imat].driftmedium = true;
238 } else {
239 materials[imat].driftmedium = false;
240 }
241 }
242 }
243
244 // Tell how many lines read
245 std::cout << m_className << "::Initialise:\n";
246 std::cout << " Read properties of " << m_nMaterials
247 << " materials from file " << mplist << ".\n";
248 if (m_debug) PrintMaterials();
249
250 // Open the element list
251 std::ifstream felist;
252 felist.open(elist.c_str(), std::ios::in);
253 if (felist.fail()) {
254 std::cerr << m_className << "::Initialise:\n";
255 std::cerr << " Could not open element file " << elist
256 << " for reading.\n";
257 std::cerr << " The file perhaps does not exist.\n";
258 return false;
259 }
260 // Read the element list
261 elements.clear();
262 nElements = 0;
263 Element newElement;
264 int ndegenerate = 0;
265 int nbackground = 0;
266 il = 0;
267 int highestnode = 0;
268 while (felist.getline(line, size, '\n')) {
269 il++;
270 // Split the line in tokens
271 char* token = NULL;
272 // Split into tokens
273 token = strtok(line, " ");
274 // Skip blank lines and headers
275 if (!token || strcmp(token, " ") == 0 || strcmp(token, "\n") == 0 ||
276 int(token[0]) == 10 || int(token[0]) == 13 ||
277 strcmp(token, "LIST") == 0 || strcmp(token, "ELEM") == 0)
278 continue;
279 // Read the element
280 int ielem = ReadInteger(token, -1, readerror);
281 token = strtok(NULL, " ");
282 int imat = ReadInteger(token, -1, readerror);
283 token = strtok(NULL, " ");
284 token = strtok(NULL, " ");
285 token = strtok(NULL, " ");
286 token = strtok(NULL, " ");
287 token = strtok(NULL, " ");
288 if (!token) std::cerr << "No token available\n";
289 int in0 = ReadInteger(token, -1, readerror);
290 token = strtok(NULL, " ");
291 int in1 = ReadInteger(token, -1, readerror);
292 token = strtok(NULL, " ");
293 int in2 = ReadInteger(token, -1, readerror);
294 token = strtok(NULL, " ");
295 int in3 = ReadInteger(token, -1, readerror);
296 token = strtok(NULL, " ");
297 int in4 = ReadInteger(token, -1, readerror);
298 token = strtok(NULL, " ");
299 int in5 = ReadInteger(token, -1, readerror);
300 token = strtok(NULL, " ");
301 int in6 = ReadInteger(token, -1, readerror);
302 token = strtok(NULL, " ");
303 int in7 = ReadInteger(token, -1, readerror);
304
305 // Check synchronisation
306 if (readerror) {
307 std::cerr << m_className << "::Initialise:\n";
308 std::cerr << " Error reading file " << elist << " (line " << il
309 << ").\n";
310 felist.close();
311 ok = false;
312 return false;
313 } else if (ielem - 1 != nElements + nbackground) {
314 std::cerr << m_className << "::Initialise:\n";
315 std::cerr << " Synchronisation lost on file " << elist << " (line "
316 << il << ").\n";
317 std::cerr << " Element: " << ielem << " (expected " << nElements
318 << "), material: " << imat << ", nodes: (" << in0 << ", " << in1
319 << ", " << in2 << ", " << in3 << ", " << in4 << ", " << in5
320 << ", " << in6 << ", " << in7 << ")\n";
321 ok = false;
322 }
323 // Check the material number and ensure that epsilon is non-negative
324 if (imat < 1 || imat > (int)m_nMaterials) {
325 std::cerr << m_className << "::Initialise:\n";
326 std::cerr << " Out-of-range material number on file " << elist
327 << " (line " << il << ").\n";
328 std::cerr << " Element: " << ielem << ", material: " << imat
329 << ", nodes: (" << in0 << ", " << in1 << ", " << in2 << ", "
330 << in3 << ", " << in4 << ", " << in5 << ", " << in6 << ", "
331 << in7 << ")\n";
332 ok = false;
333 }
334 if (materials[imat - 1].eps < 0) {
335 std::cerr << m_className << "::Initialise:\n";
336 std::cerr << " Element " << ielem << " in element list " << elist
337 << " uses material " << imat << " which\n",
338 std::cerr << " has not been assigned a positive permittivity\n";
339 std::cerr << " in material list " << mplist << ".\n";
340 ok = false;
341 }
342 // Check the node numbers
343 if (in0 < 1 || in1 < 1 || in2 < 1 || in3 < 1 || in4 < 1 || in5 < 1 ||
344 in6 < 1 || in7 < 1) {
345 std::cerr << m_className << "::Initialise:\n";
346 std::cerr << " Found a node number < 1 on file " << elist << " (line "
347 << il << ").\n";
348 std::cerr << " Element: " << ielem << ", material: " << imat << "\n";
349 ok = false;
350 }
351 if (in0 > highestnode) highestnode = in0;
352 if (in1 > highestnode) highestnode = in1;
353 if (in2 > highestnode) highestnode = in2;
354 if (in3 > highestnode) highestnode = in3;
355 if (in4 > highestnode) highestnode = in4;
356 if (in5 > highestnode) highestnode = in5;
357 if (in6 > highestnode) highestnode = in6;
358 if (in7 > highestnode) highestnode = in7;
359 // Skip quadrilaterals which are background.
360 if (m_deleteBackground && materials[imat - 1].ohm == 0) {
361 nbackground++;
362 continue;
363 }
364 // Store the element, degeneracy
365 if (in2 == in3 && in3 == in6) {
366 ndegenerate++;
367 newElement.degenerate = true;
368 } else {
369 newElement.degenerate = false;
370 }
371 // Store the material reference
372 newElement.matmap = imat - 1;
373 // Node references
374 if (newElement.degenerate) {
375 newElement.emap[0] = in0 - 1;
376 newElement.emap[1] = in1 - 1;
377 newElement.emap[2] = in2 - 1;
378 newElement.emap[3] = in4 - 1;
379 newElement.emap[4] = in7 - 1;
380 newElement.emap[5] = in5 - 1;
381 newElement.emap[6] = in3 - 1;
382 newElement.emap[7] = in6 - 1;
383 } else {
384 newElement.emap[0] = in0 - 1;
385 newElement.emap[1] = in1 - 1;
386 newElement.emap[2] = in2 - 1;
387 newElement.emap[3] = in3 - 1;
388 newElement.emap[4] = in4 - 1;
389 newElement.emap[5] = in5 - 1;
390 newElement.emap[6] = in6 - 1;
391 newElement.emap[7] = in7 - 1;
392 }
393 elements.push_back(newElement);
394 ++nElements;
395 }
396 // Close the file
397 felist.close();
398 // Tell how many lines read
399 std::cout << m_className << "::Initialise:\n";
400 std::cout << " Read " << nElements << " elements from file " << elist
401 << ",\n";
402 std::cout << " highest node number: " << highestnode << ",\n";
403 std::cout << " degenerate elements: " << ndegenerate << ",\n";
404 std::cout << " background elements skipped: " << nbackground << ".\n";
405 // Check the value of the unit
406 double funit;
407 if (strcmp(unit.c_str(), "mum") == 0 || strcmp(unit.c_str(), "micron") == 0 ||
408 strcmp(unit.c_str(), "micrometer") == 0) {
409 funit = 0.0001;
410 } else if (strcmp(unit.c_str(), "mm") == 0 ||
411 strcmp(unit.c_str(), "millimeter") == 0) {
412 funit = 0.1;
413 } else if (strcmp(unit.c_str(), "cm") == 0 ||
414 strcmp(unit.c_str(), "centimeter") == 0) {
415 funit = 1.0;
416 } else if (strcmp(unit.c_str(), "m") == 0 ||
417 strcmp(unit.c_str(), "meter") == 0) {
418 funit = 100.0;
419 } else {
420 std::cerr << m_className << "::Initialise:\n";
421 std::cerr << " Unknown length unit " << unit << ".\n";
422 ok = false;
423 funit = 1.0;
424 }
425 if (m_debug) {
426 std::cout << m_className << "::Initialise:\n";
427 std::cout << " Unit scaling factor = " << funit << ".\n";
428 }
429
430 // Open the node list
431 std::ifstream fnlist;
432 fnlist.open(nlist.c_str(), std::ios::in);
433 if (fnlist.fail()) {
434 std::cerr << m_className << "::Initialise:\n";
435 std::cerr << " Could not open nodes file " << nlist << " for reading.\n";
436 std::cerr << " The file perhaps does not exist.\n";
437 return false;
438 }
439 // Read the node list
440 nodes.clear();
441 nNodes = 0;
442 Node newNode;
443 newNode.w.clear();
444 il = 0;
445 while (fnlist.getline(line, size, '\n')) {
446 il++;
447 // Split the line in tokens
448 char* token = NULL;
449 // Split into tokens
450 token = strtok(line, " ");
451 // Skip blank lines and headers
452 if (!token || strcmp(token, " ") == 0 || strcmp(token, "\n") == 0 ||
453 int(token[0]) == 10 || int(token[0]) == 13 ||
454 strcmp(token, "LIST") == 0 || strcmp(token, "NODE") == 0)
455 continue;
456 // Read the element
457 int inode = ReadInteger(token, -1, readerror);
458 token = strtok(NULL, " ");
459 double xnode = ReadDouble(token, -1, readerror);
460 token = strtok(NULL, " ");
461 double ynode = ReadDouble(token, -1, readerror);
462 token = strtok(NULL, " ");
463 double znode = ReadDouble(token, -1, readerror);
464 // Check syntax
465 if (readerror) {
466 std::cerr << m_className << "::Initialise:\n";
467 std::cerr << " Error reading file " << nlist << " (line " << il
468 << ").\n";
469 fnlist.close();
470 ok = false;
471 return false;
472 }
473 // Check synchronisation
474 if (inode - 1 != nNodes) {
475 std::cerr << m_className << "::Initialise:\n";
476 std::cerr << " Synchronisation lost on file " << nlist << " (line "
477 << il << ").\n";
478 std::cerr << " Node: " << inode << " (expected " << nNodes
479 << "), (x,y,z) = (" << xnode << ", " << ynode << ", " << znode
480 << ")\n";
481 ok = false;
482 }
483
484 // Store the point coordinates
485 newNode.x = xnode * funit;
486 newNode.y = ynode * funit;
487 newNode.z = znode * funit;
488 nodes.push_back(newNode);
489 ++nNodes;
490 }
491 // Close the file
492 fnlist.close();
493 // Tell how many lines read
494 std::cout << m_className << "::Initialise:\n";
495 std::cout << " Read " << nNodes << " nodes from file " << nlist << ".\n";
496 // Check number of nodes
497 if (nNodes != highestnode) {
498 std::cerr << m_className << "::Initialise:\n";
499 std::cerr << " Number of nodes read (" << nNodes << ") on " << nlist
500 << "\n";
501 std::cerr << " does not match element list (" << highestnode << ").\n";
502 ok = false;
503 }
504
505 // Open the voltage list
506 std::ifstream fprnsol;
507 fprnsol.open(prnsol.c_str(), std::ios::in);
508 if (fprnsol.fail()) {
509 std::cerr << m_className << "::Initialise:\n";
510 std::cerr << " Could not open potential file " << prnsol
511 << " for reading.\n";
512 std::cerr << " The file perhaps does not exist.\n";
513 return false;
514 }
515 // Read the voltage list
516 il = 0;
517 int nread = 0;
518 readerror = false;
519 while (fprnsol.getline(line, size, '\n')) {
520 il++;
521 // Split the line in tokens
522 char* token = NULL;
523 token = strtok(line, " ");
524 // Skip blank lines and headers
525 if (!token || strcmp(token, " ") == 0 || strcmp(token, "\n") == 0 ||
526 int(token[0]) == 10 || int(token[0]) == 13 ||
527 strcmp(token, "PRINT") == 0 || strcmp(token, "*****") == 0 ||
528 strcmp(token, "LOAD") == 0 || strcmp(token, "TIME=") == 0 ||
529 strcmp(token, "MAXIMUM") == 0 || strcmp(token, "VALUE") == 0 ||
530 strcmp(token, "NODE") == 0)
531 continue;
532 // Read the element
533 int inode = ReadInteger(token, -1, readerror);
534 token = strtok(NULL, " ");
535 double volt = ReadDouble(token, -1, readerror);
536 // Check syntax
537 if (readerror) {
538 std::cerr << m_className << "::Initialise:\n";
539 std::cerr << " Error reading file " << prnsol << " (line " << il
540 << ").\n";
541 fprnsol.close();
542 ok = false;
543 return false;
544 }
545 // Check node number and store if OK
546 if (inode < 1 || inode > nNodes) {
547 std::cerr << m_className << "::Initialise:\n";
548 std::cerr << " Node number " << inode << " out of range\n";
549 std::cerr << " on potential file " << prnsol << " (line " << il
550 << ").\n";
551 ok = false;
552 } else {
553 nodes[inode - 1].v = volt;
554 nread++;
555 }
556 }
557 // Close the file
558 fprnsol.close();
559 // Tell how many lines read
560 std::cout << m_className << "::Initialise:\n";
561 std::cout << " Read " << nread << " potentials from file " << prnsol
562 << ".\n";
563 // Check number of nodes
564 if (nread != nNodes) {
565 std::cerr << m_className << "::Initialise:\n";
566 std::cerr << " Number of nodes read (" << nread << ") on potential file "
567 << prnsol << " does not\n",
568 std::cerr << " match the node list (" << nNodes << ").\n";
569 ok = false;
570 }
571 // Set the ready flag
572 if (ok) {
573 m_ready = true;
574 } else {
575 std::cerr << m_className << "::Initialise:\n";
576 std::cerr
577 << " Field map could not be read and cannot be interpolated.\n";
578 return false;
579 }
580
581 // Remove weighting fields (if any).
582 wfields.clear();
583 wfieldsOk.clear();
585
586 // Establish the ranges
587 SetRange();
589 return true;
590}
591
593 std::string label) {
594 if (!m_ready) {
595 PrintNotReady("SetWeightingField");
596 std::cerr << " Weighting field cannot be added.\n";
597 return false;
598 }
599
600 // Open the voltage list.
601 std::ifstream fprnsol;
602 fprnsol.open(prnsol.c_str(), std::ios::in);
603 if (fprnsol.fail()) {
604 std::cerr << m_className << "::SetWeightingField:\n";
605 std::cerr << " Could not open potential file " << prnsol
606 << " for reading.\n";
607 std::cerr << " The file perhaps does not exist.\n";
608 return false;
609 }
610
611 // Check if a weighting field with the same label already exists.
612 int iw = nWeightingFields;
613 for (int i = nWeightingFields; i--;) {
614 if (wfields[i] == label) {
615 iw = i;
616 break;
617 }
618 }
619 if (iw == nWeightingFields) {
623 for (int j = nNodes; j--;) {
624 nodes[j].w.resize(nWeightingFields);
625 }
626 } else {
627 std::cout << m_className << "::SetWeightingField:\n";
628 std::cout << " Replacing existing weighting field " << label << ".\n";
629 }
630 wfields[iw] = label;
631 wfieldsOk[iw] = false;
632
633 // Buffer for reading
634 const int size = 100;
635 char line[size];
636
637 bool ok = true;
638 // Read the voltage list.
639 int il = 0;
640 int nread = 0;
641 bool readerror = false;
642 while (fprnsol.getline(line, size, '\n')) {
643 il++;
644 // Split the line in tokens.
645 char* token = NULL;
646 token = strtok(line, " ");
647 // Skip blank lines and headers.
648 if (!token || strcmp(token, " ") == 0 || strcmp(token, "\n") == 0 ||
649 int(token[0]) == 10 || int(token[0]) == 13 ||
650 strcmp(token, "PRINT") == 0 || strcmp(token, "*****") == 0 ||
651 strcmp(token, "LOAD") == 0 || strcmp(token, "TIME=") == 0 ||
652 strcmp(token, "MAXIMUM") == 0 || strcmp(token, "VALUE") == 0 ||
653 strcmp(token, "NODE") == 0)
654 continue;
655 // Read the element.
656 int inode = ReadInteger(token, -1, readerror);
657 token = strtok(NULL, " ");
658 double volt = ReadDouble(token, -1, readerror);
659 // Check the syntax.
660 if (readerror) {
661 std::cerr << m_className << "::SetWeightingField:\n";
662 std::cerr << " Error reading file " << prnsol << " (line " << il
663 << ").\n";
664 fprnsol.close();
665 return false;
666 }
667 // Check node number and store if OK.
668 if (inode < 1 || inode > nNodes) {
669 std::cerr << m_className << "::SetWeightingField:\n";
670 std::cerr << " Node number " << inode << " out of range\n";
671 std::cerr << " on potential file " << prnsol << " (line " << il
672 << ").\n";
673 ok = false;
674 } else {
675 nodes[inode - 1].w[iw] = volt;
676 nread++;
677 }
678 }
679 // Close the file.
680 fprnsol.close();
681 // Tell how many lines read.
682 std::cout << m_className << "::SetWeightingField:\n";
683 std::cout << " Read " << nread << " potentials from file " << prnsol
684 << ".\n";
685 // Check the number of nodes.
686 if (nread != nNodes) {
687 std::cerr << m_className << "::SetWeightingField:\n";
688 std::cerr << " Number of nodes read (" << nread << ") "
689 << " on potential file " << prnsol << "\n";
690 std::cerr << " does not match the node list (" << nNodes << ").\n";
691 ok = false;
692 }
693
694 // Set the ready flag.
695 wfieldsOk[iw] = ok;
696 if (!ok) {
697 std::cerr << m_className << "::SetWeightingField:\n";
698 std::cerr << " Field map could not be read "
699 << "and cannot be interpolated.\n";
700 return false;
701 }
702 return true;
703}
704
705void ComponentAnsys121::ElectricField(const double x, const double y,
706 const double z, double& ex, double& ey,
707 double& ez, Medium*& m, int& status) {
708 double v;
709 ElectricField(x, y, z, ex, ey, ez, v, m, status);
710}
711
712void ComponentAnsys121::ElectricField(const double xin, const double yin,
713 const double zin, double& ex, double& ey,
714 double& ez, double& volt, Medium*& m,
715 int& status) {
716 // Copy the coordinates
717 double x = xin, y = yin, z = 0.;
718
719 // Map the coordinates onto field map coordinates
720 bool xmirr, ymirr, zmirr;
721 double rcoordinate, rotation;
722 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
723
724 // Initial values
725 ex = ey = ez = volt = 0;
726 m = nullptr;
727 status = 0;
728
729 // Do not proceed if not properly initialised.
730 if (!m_ready) {
731 status = -10;
732 PrintNotReady("ElectricField");
733 return;
734 }
735
736 if (m_warning) PrintWarning("ElectricField");
737
738 if (zin < m_minBoundingBox[2] || zin > m_maxBoundingBox[2]) {
739 status = -5;
740 return;
741 }
742
743 // Find the element that contains this point
744 double t1, t2, t3, t4, jac[4][4], det;
745 const int imap = FindElement5(x, y, z, t1, t2, t3, t4, jac, det);
746 if (imap < 0) {
747 if (m_debug) {
748 std::cout << m_className << "::ElectricField:\n";
749 std::cout << " Point (" << x << ", " << y << ") not in the mesh.\n";
750 }
751 status = -6;
752 return;
753 }
754
755 const Element& element = elements[imap];
756 if (m_debug) {
757 PrintElement("ElectricField", x, y, z, t1, t2, t3, t4, element, 8);
758 }
759
760 const Node& n0 = nodes[element.emap[0]];
761 const Node& n1 = nodes[element.emap[1]];
762 const Node& n2 = nodes[element.emap[2]];
763 const Node& n3 = nodes[element.emap[3]];
764 const Node& n4 = nodes[element.emap[4]];
765 const Node& n5 = nodes[element.emap[5]];
766 // Calculate quadrilateral field, which can degenerate to a triangular field
767 const double invdet = 1. / det;
768 if (element.degenerate) {
769 volt = n0.v * t1 * (2 * t1 - 1) + n1.v * t2 * (2 * t2 - 1) +
770 n2.v * t3 * (2 * t3 - 1) + 4 * n3.v * t1 * t2 + 4 * n4.v * t1 * t3 +
771 4 * n5.v * t2 * t3;
772 ex = -(n0.v * (4 * t1 - 1) * jac[0][1] + n1.v * (4 * t2 - 1) * jac[1][1] +
773 n2.v * (4 * t3 - 1) * jac[2][1] +
774 n3.v * (4 * t2 * jac[0][1] + 4 * t1 * jac[1][1]) +
775 n4.v * (4 * t3 * jac[0][1] + 4 * t1 * jac[2][1]) +
776 n5.v * (4 * t3 * jac[1][1] + 4 * t2 * jac[2][1])) *
777 invdet;
778 ey = -(n0.v * (4 * t1 - 1) * jac[0][2] + n1.v * (4 * t2 - 1) * jac[1][2] +
779 n2.v * (4 * t3 - 1) * jac[2][2] +
780 n3.v * (4 * t2 * jac[0][2] + 4 * t1 * jac[1][2]) +
781 n4.v * (4 * t3 * jac[0][2] + 4 * t1 * jac[2][2]) +
782 n5.v * (4 * t3 * jac[1][2] + 4 * t2 * jac[2][2])) *
783 invdet;
784 } else {
785 const Node& n6 = nodes[element.emap[6]];
786 const Node& n7 = nodes[element.emap[7]];
787 volt = -n0.v * (1 - t1) * (1 - t2) * (1 + t1 + t2) * 0.25 -
788 n1.v * (1 + t1) * (1 - t2) * (1 - t1 + t2) * 0.25 -
789 n2.v * (1 + t1) * (1 + t2) * (1 - t1 - t2) * 0.25 -
790 n3.v * (1 - t1) * (1 + t2) * (1 + t1 - t2) * 0.25 +
791 n4.v * (1 - t1) * (1 + t1) * (1 - t2) * 0.5 +
792 n5.v * (1 + t1) * (1 + t2) * (1 - t2) * 0.5 +
793 n6.v * (1 - t1) * (1 + t1) * (1 + t2) * 0.5 +
794 n7.v * (1 - t1) * (1 + t2) * (1 - t2) * 0.5;
795 ex = -(n0.v * ((1 - t2) * (2 * t1 + t2) * jac[0][0] +
796 (1 - t1) * (t1 + 2 * t2) * jac[1][0]) *
797 0.25 +
798 n1.v * ((1 - t2) * (2 * t1 - t2) * jac[0][0] -
799 (1 + t1) * (t1 - 2 * t2) * jac[1][0]) *
800 0.25 +
801 n2.v * ((1 + t2) * (2 * t1 + t2) * jac[0][0] +
802 (1 + t1) * (t1 + 2 * t2) * jac[1][0]) *
803 0.25 +
804 n3.v * ((1 + t2) * (2 * t1 - t2) * jac[0][0] -
805 (1 - t1) * (t1 - 2 * t2) * jac[1][0]) *
806 0.25 +
807 n4.v * (t1 * (t2 - 1) * jac[0][0] +
808 (t1 - 1) * (t1 + 1) * jac[1][0] * 0.5) +
809 n5.v * ((1 - t2) * (1 + t2) * jac[0][0] * 0.5 -
810 (1 + t1) * t2 * jac[1][0]) +
811 n6.v * (-t1 * (1 + t2) * jac[0][0] +
812 (1 - t1) * (1 + t1) * jac[1][0] * 0.5) +
813 n7.v * ((t2 - 1) * (t2 + 1) * jac[0][0] * 0.5 +
814 (t1 - 1) * t2 * jac[1][0])) *
815 invdet;
816 ey = -(n0.v * ((1 - t2) * (2 * t1 + t2) * jac[0][1] +
817 (1 - t1) * (t1 + 2 * t2) * jac[1][1]) *
818 0.25 +
819 n1.v * ((1 - t2) * (2 * t1 - t2) * jac[0][1] -
820 (1 + t1) * (t1 - 2 * t2) * jac[1][1]) *
821 0.25 +
822 n2.v * ((1 + t2) * (2 * t1 + t2) * jac[0][1] +
823 (1 + t1) * (t1 + 2 * t2) * jac[1][1]) *
824 0.25 +
825 n3.v * ((1 + t2) * (2 * t1 - t2) * jac[0][1] -
826 (1 - t1) * (t1 - 2 * t2) * jac[1][1]) *
827 0.25 +
828 n4.v * (t1 * (t2 - 1) * jac[0][1] +
829 (t1 - 1) * (t1 + 1) * jac[1][1] * 0.5) +
830 n5.v * ((1 - t2) * (1 + t2) * jac[0][1] * 0.5 -
831 (1 + t1) * t2 * jac[1][1]) +
832 n6.v * (-t1 * (1 + t2) * jac[0][1] +
833 (1 - t1) * (1 + t1) * jac[1][1] * 0.5) +
834 n7.v * ((t2 - 1) * (t2 + 1) * jac[0][1] * 0.5 +
835 (t1 - 1) * t2 * jac[1][1])) *
836 invdet;
837 }
838
839 // Transform field to global coordinates
840 UnmapFields(ex, ey, ez, x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
841
842 // Drift medium?
843 if (m_debug) {
844 std::cout << m_className << "::ElectricField:\n";
845 std::cout << " Material " << element.matmap << ", drift flag "
846 << materials[element.matmap].driftmedium << ".\n";
847 }
848 m = materials[element.matmap].medium;
849 status = -5;
850 if (materials[element.matmap].driftmedium) {
851 if (m && m->IsDriftable()) status = 0;
852 }
853}
854
855void ComponentAnsys121::WeightingField(const double xin, const double yin,
856 const double zin, double& wx, double& wy,
857 double& wz, const std::string& label) {
858 // Initial values
859 wx = wy = wz = 0;
860
861 // Do not proceed if not properly initialised.
862 if (!m_ready) return;
863
864 // Look for the label.
865 int iw = 0;
866 bool found = false;
867 for (int i = nWeightingFields; i--;) {
868 if (wfields[i] == label) {
869 iw = i;
870 found = true;
871 break;
872 }
873 }
874
875 // Do not proceed if the requested weighting field does not exist.
876 if (!found) return;
877 // Check if the weighting field is properly initialised.
878 if (!wfieldsOk[iw]) return;
879
880 // Copy the coordinates.
881 double x = xin, y = yin, z = zin;
882
883 // Map the coordinates onto field map coordinates
884 bool xmirr, ymirr, zmirr;
885 double rcoordinate, rotation;
886 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
887
888 if (m_warning) PrintWarning("WeightingField");
889
890 // Find the element that contains this point.
891 double t1, t2, t3, t4, jac[4][4], det;
892 const int imap = FindElement5(x, y, z, t1, t2, t3, t4, jac, det);
893 // Check if the point is in the mesh.
894 if (imap < 0) return;
895
896 const Element& element = elements[imap];
897 if (m_debug) {
898 PrintElement("WeightingField", x, y, z, t1, t2, t3, t4, element, 8, iw);
899 }
900 const Node& n0 = nodes[element.emap[0]];
901 const Node& n1 = nodes[element.emap[1]];
902 const Node& n2 = nodes[element.emap[2]];
903 const Node& n3 = nodes[element.emap[3]];
904 const Node& n4 = nodes[element.emap[4]];
905 const Node& n5 = nodes[element.emap[5]];
906 // Calculate quadrilateral field, which can degenerate to a triangular field
907 const double invdet = 1. / det;
908 if (elements[imap].degenerate) {
909 wx = -(n0.w[iw] * (4 * t1 - 1) * jac[0][1] +
910 n1.w[iw] * (4 * t2 - 1) * jac[1][1] +
911 n2.w[iw] * (4 * t3 - 1) * jac[2][1] +
912 n3.w[iw] * (4 * t2 * jac[0][1] + 4 * t1 * jac[1][1]) +
913 n4.w[iw] * (4 * t3 * jac[0][1] + 4 * t1 * jac[2][1]) +
914 n5.w[iw] * (4 * t3 * jac[1][1] + 4 * t2 * jac[2][1])) *
915 invdet;
916 wy = -(n0.w[iw] * (4 * t1 - 1) * jac[0][2] +
917 n1.w[iw] * (4 * t2 - 1) * jac[1][2] +
918 n2.w[iw] * (4 * t3 - 1) * jac[2][2] +
919 n3.w[iw] * (4 * t2 * jac[0][2] + 4 * t1 * jac[1][2]) +
920 n4.w[iw] * (4 * t3 * jac[0][2] + 4 * t1 * jac[2][2]) +
921 n5.w[iw] * (4 * t3 * jac[1][2] + 4 * t2 * jac[2][2])) *
922 invdet;
923 } else {
924 const Node& n6 = nodes[element.emap[6]];
925 const Node& n7 = nodes[element.emap[7]];
926 wx = -(n0.w[iw] * ((1 - t2) * (2 * t1 + t2) * jac[0][0] +
927 (1 - t1) * (t1 + 2 * t2) * jac[1][0]) *
928 0.25 +
929 n1.w[iw] * ((1 - t2) * (2 * t1 - t2) * jac[0][0] -
930 (1 + t1) * (t1 - 2 * t2) * jac[1][0]) *
931 0.25 +
932 n2.w[iw] * ((1 + t2) * (2 * t1 + t2) * jac[0][0] +
933 (1 + t1) * (t1 + 2 * t2) * jac[1][0]) *
934 0.25 +
935 n3.w[iw] * ((1 + t2) * (2 * t1 - t2) * jac[0][0] -
936 (1 - t1) * (t1 - 2 * t2) * jac[1][0]) *
937 0.25 +
938 n4.w[iw] * (t1 * (t2 - 1) * jac[0][0] +
939 (t1 - 1) * (t1 + 1) * jac[1][0] * 0.5) +
940 n5.w[iw] * ((1 - t2) * (1 + t2) * jac[0][0] * 0.5 -
941 (1 + t1) * t2 * jac[1][0]) +
942 n6.w[iw] * (-t1 * (1 + t2) * jac[0][0] +
943 (1 - t1) * (1 + t1) * jac[1][0] * 0.5) +
944 n7.w[iw] * ((t2 - 1) * (1 + t2) * jac[0][0] * 0.5 +
945 (t1 - 1) * t2 * jac[1][0])) *
946 invdet;
947 wy = -(n0.w[iw] * ((1 - t2) * (2 * t1 + t2) * jac[0][1] +
948 (1 - t1) * (t1 + 2 * t2) * jac[1][1]) *
949 0.25 +
950 n1.w[iw] * ((1 - t2) * (2 * t1 - t2) * jac[0][1] -
951 (1 + t1) * (t1 - 2 * t2) * jac[1][1]) *
952 0.25 +
953 n2.w[iw] * ((1 + t2) * (2 * t1 + t2) * jac[0][1] +
954 (1 + t1) * (t1 + 2 * t2) * jac[1][1]) *
955 0.25 +
956 n3.w[iw] * ((1 + t2) * (2 * t1 - t2) * jac[0][1] -
957 (1 - t1) * (t1 - 2 * t2) * jac[1][1]) *
958 0.25 +
959 n4.w[iw] * (t1 * (t2 - 1) * jac[0][1] +
960 (t1 - 1) * (t1 + 1) * jac[1][1] * 0.5) +
961 n5.w[iw] * ((1 - t2) * (1 + t2) * jac[0][1] * 0.5 -
962 (1 + t1) * t2 * jac[1][1]) +
963 n6.w[iw] * (-t1 * (1 + t2) * jac[0][1] +
964 (1 - t1) * (1 + t1) * jac[1][1] * 0.5) +
965 n7.w[iw] * ((t2 - 1) * (t2 + 1) * jac[0][1] * 0.5 +
966 (t1 - 1) * t2 * jac[1][1])) *
967 invdet;
968 }
969
970 // Transform field to global coordinates
971 UnmapFields(wx, wy, wz, x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
972}
973
974double ComponentAnsys121::WeightingPotential(const double xin, const double yin,
975 const double zin,
976 const std::string& label) {
977 // Do not proceed if not properly initialised.
978 if (!m_ready) return 0.;
979
980 // Look for the label.
981 int iw = 0;
982 bool found = false;
983 for (int i = nWeightingFields; i--;) {
984 if (wfields[i] == label) {
985 iw = i;
986 found = true;
987 break;
988 }
989 }
990
991 // Do not proceed if the requested weighting field does not exist.
992 if (!found) return 0.;
993 // Check if the weighting field is properly initialised.
994 if (!wfieldsOk[iw]) return 0.;
995
996 // Copy the coordinates.
997 double x = xin, y = yin, z = zin;
998
999 // Map the coordinates onto field map coordinates.
1000 bool xmirr, ymirr, zmirr;
1001 double rcoordinate, rotation;
1002 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
1003
1004 if (m_warning) PrintWarning("WeightingPotential");
1005
1006 // Find the element that contains this point.
1007 double t1, t2, t3, t4, jac[4][4], det;
1008 const int imap = FindElement5(x, y, z, t1, t2, t3, t4, jac, det);
1009 // Check if the point is in the mesh
1010 if (imap < 0) return 0.;
1011
1012 const Element& element = elements[imap];
1013 if (m_debug) {
1014 PrintElement("WeightingPotential", x, y, z, t1, t2, t3, t4, element, 8, iw);
1015 }
1016 // Calculate quadrilateral field, which can degenerate to a triangular field
1017 const Node& n0 = nodes[element.emap[0]];
1018 const Node& n1 = nodes[element.emap[1]];
1019 const Node& n2 = nodes[element.emap[2]];
1020 const Node& n3 = nodes[element.emap[3]];
1021 const Node& n4 = nodes[element.emap[4]];
1022 const Node& n5 = nodes[element.emap[5]];
1023 if (element.degenerate) {
1024 return n0.w[iw] * t1 * (2 * t1 - 1) + n1.w[iw] * t2 * (2 * t2 - 1) +
1025 n2.w[iw] * t3 * (2 * t3 - 1) + 4 * n3.w[iw] * t1 * t2 +
1026 4 * n4.w[iw] * t1 * t3 + 4 * n5.w[iw] * t2 * t3;
1027 }
1028
1029 const Node& n6 = nodes[element.emap[6]];
1030 const Node& n7 = nodes[element.emap[7]];
1031 return -n0.w[iw] * (1 - t1) * (1 - t2) * (1 + t1 + t2) * 0.25 -
1032 n1.w[iw] * (1 + t1) * (1 - t2) * (1 - t1 + t2) * 0.25 -
1033 n2.w[iw] * (1 + t1) * (1 + t2) * (1 - t1 - t2) * 0.25 -
1034 n3.w[iw] * (1 - t1) * (1 + t2) * (1 + t1 - t2) * 0.25 +
1035 n4.w[iw] * (1 - t1) * (1 + t1) * (1 - t2) * 0.5 +
1036 n5.w[iw] * (1 + t1) * (1 + t2) * (1 - t2) * 0.5 +
1037 n6.w[iw] * (1 - t1) * (1 + t1) * (1 + t2) * 0.5 +
1038 n7.w[iw] * (1 - t1) * (1 + t2) * (1 - t2) * 0.5;
1039}
1040
1041Medium* ComponentAnsys121::GetMedium(const double xin, const double yin,
1042 const double zin) {
1043 // Copy the coordinates.
1044 double x = xin, y = yin, z = 0.;
1045
1046 // Map the coordinates onto field map coordinates.
1047 bool xmirr, ymirr, zmirr;
1048 double rcoordinate, rotation;
1049 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
1050
1051 if (zin < m_minBoundingBox[2] || z > m_maxBoundingBox[2]) {
1052 return nullptr;
1053 }
1054
1055 // Do not proceed if not properly initialised.
1056 if (!m_ready) {
1057 PrintNotReady("GetMedium");
1058 return nullptr;
1059 }
1060 if (m_warning) PrintWarning("GetMedium");
1061
1062 // Find the element that contains this point.
1063 double t1, t2, t3, t4, jac[4][4], det;
1064 const int imap = FindElement5(x, y, z, t1, t2, t3, t4, jac, det);
1065 if (imap < 0) {
1066 if (m_debug) {
1067 std::cerr << m_className << "::GetMedium:\n";
1068 std::cerr << " Point (" << x << ", " << y << ") not in the mesh.\n";
1069 }
1070 return nullptr;
1071 }
1072 const Element& element = elements[imap];
1073 if (element.matmap >= m_nMaterials) {
1074 if (m_debug) {
1075 std::cerr << m_className << "::GetMedium:\n";
1076 std::cerr << " Point (" << x << ", " << y << ")"
1077 << " has out of range material number " << imap << ".\n";
1078 }
1079 return nullptr;
1080 }
1081
1082 if (m_debug) {
1083 PrintElement("GetMedium", x, y, z, t1, t2, t3, t4, element, 8);
1084 }
1085
1086 // Assign a medium.
1087 return materials[element.matmap].medium;
1088}
1089
1090void ComponentAnsys121::SetRangeZ(const double zmin, const double zmax) {
1091 if (fabs(zmax - zmin) <= 0.) {
1092 std::cerr << m_className << "::SetRangeZ: Zero range is not permitted.\n";
1093 return;
1094 }
1095 m_minBoundingBox[2] = m_mapmin[2] = std::min(zmin, zmax);
1096 m_maxBoundingBox[2] = m_mapmax[2] = std::max(zmin, zmax);
1097}
1098
1102}
1103
1104double ComponentAnsys121::GetElementVolume(const unsigned int i) {
1105 if (i >= elements.size()) return 0.;
1106 const Element& element = elements[i];
1107 const Node& n0 = nodes[element.emap[0]];
1108 const Node& n1 = nodes[element.emap[1]];
1109 const Node& n2 = nodes[element.emap[2]];
1110 const Node& n3 = nodes[element.emap[3]];
1111 const double surf =
1112 0.5 *
1113 (fabs((n1.x - n0.x) * (n2.y - n0.y) - (n2.x - n0.x) * (n1.y - n0.y)) +
1114 fabs((n3.x - n0.x) * (n2.y - n0.y) - (n2.x - n0.x) * (n3.y - n0.y)));
1115 return surf;
1116}
1117
1118void ComponentAnsys121::GetAspectRatio(const unsigned int i, double& dmin,
1119 double& dmax) {
1120 if (i >= elements.size()) {
1121 dmin = dmax = 0.;
1122 return;
1123 }
1124
1125 const Element& element = elements[i];
1126 const int np = 8;
1127 // Loop over all pairs of vertices.
1128 for (int j = 0; j < np - 1; ++j) {
1129 const Node& nj = nodes[element.emap[j]];
1130 for (int k = j + 1; k < np; ++k) {
1131 const Node& nk = nodes[element.emap[k]];
1132 // Compute distance.
1133 const double dx = nj.x - nk.x;
1134 const double dy = nj.y - nk.y;
1135 const double dist = sqrt(dx * dx + dy * dy);
1136 if (k == 1) {
1137 dmin = dmax = dist;
1138 } else {
1139 if (dist < dmin) dmin = dist;
1140 if (dist > dmax) dmax = dist;
1141 }
1142 }
1143 }
1144}
1145}
bool SetWeightingField(std::string prnsol, std::string label)
void GetAspectRatio(const unsigned int i, double &dmin, double &dmax) override
bool Initialise(std::string elist="ELIST.lis", std::string nlist="NLIST.lis", std::string mplist="MPLIST.lis", std::string prnsol="PRNSOL.lis", std::string unit="cm")
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 UpdatePeriodicity() override
Verify periodicities.
void SetRangeZ(const double zmin, const double zmax)
double WeightingPotential(const double x, const double y, const double z, const std::string &label) override
void WeightingField(const double x, const double y, const double z, double &wx, double &wy, double &wz, 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).
bool m_ready
Ready for use?
std::string m_className
Class name.
bool m_debug
Switch on/off debugging messages.
Base class for components based on finite-element field maps.
void PrintMaterials()
List all currently defined materials.
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)
std::array< double, 3 > m_maxBoundingBox
std::array< double, 3 > m_mapmin
std::array< double, 3 > m_minBoundingBox
int ReadInteger(char *token, int def, bool &error)
int FindElement5(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 quadrilaterals.
std::vector< Material > materials
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::array< double, 3 > m_mapmax
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