-
Notifications
You must be signed in to change notification settings - Fork 110
/
Copy pathSilicon.cpp
1594 lines (1411 loc) · 68.9 KB
/
Silicon.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/* -*- c++ -*-
* Copyright (c) 2012-2023 by the GalSim developers team on GitHub
* https://github.com/GalSim-developers
*
* This file is part of GalSim: The modular galaxy image simulation toolkit.
* https://github.com/GalSim-developers/GalSim
*
* GalSim is free software: redistribution and use in source and binary forms,
* with or without modification, are permitted provided that the following
* conditions are met:
*
* 1. Redistributions of source code must retain the above copyright notice, this
* list of conditions, and the disclaimer given in the accompanying LICENSE
* file.
* 2. Redistributions in binary form must reproduce the above copyright notice,
* this list of conditions, and the disclaimer given in the documentation
* and/or other materials provided with the distribution.
*/
/*
* ------------------------------------------------------------------------------
* Author: Craig Lage, UC Davis
* Date: Feb 17, 2017
* Routines for integrating the CCD simulations into GalSim
*/
#include <cmath>
#include <string>
#include <fstream>
#include <sstream>
#include <vector>
#include <algorithm>
#include <climits>
// Uncomment this for debugging output
//#define DEBUGLOGGING
#ifdef DEBUGLOGGING
#undef _OPENMP
#endif
#ifdef _OPENMP
#include <omp.h>
#endif
#include "Std.h"
#include "Silicon.h"
#include "Image.h"
#include "PhotonArray.h"
#ifdef DEBUGLOGGING
// This function was used to help track down a memory leak. Saving it for future resuse.
#include<mach/mach.h>
double RSS() {
struct task_basic_info t_info;
mach_msg_type_number_t t_info_count = TASK_BASIC_INFO_COUNT;
if (KERN_SUCCESS != task_info(mach_task_self(),
TASK_BASIC_INFO, (task_info_t)&t_info,
&t_info_count))
{
return -1;
}
// resident size is in t_info.resident_size;
// virtual size is in t_info.virtual_size;
return double(t_info.resident_size) / (1 << 30); // 2^30 == GB
}
#endif
namespace galsim {
// std::min and std::max are not supported in GPU offloaded code, so define our own
// integer versions here. (C standard library provides floating point versions which do
// work on GPU).
int imin(int a, int b) { return a < b ? a : b; }
int imax(int a, int b) { return a > b ? a : b; }
// Helper function used in a few places below.
void buildEmptyPoly(Polygon& poly, int numVertices)
{
dbg<<"buildEmptyPoly\n";
double dtheta = M_PI / (2.0 * (numVertices + 1.0));
double theta0 = - M_PI / 4.0;
poly.reserve(numVertices*4 + 8);
// First the corners
dbg<<"corners:\n";
for (int xpix=0; xpix<2; xpix++) {
for (int ypix=0; ypix<2; ypix++) {
poly.add(Position<double>(xpix, ypix));
// Two copies of the corner to be consistent with new code
// that has two corner points.
poly.add(Position<double>(xpix, ypix));
}
}
// Next the edges
dbg<<"x edges:\n";
for (int xpix=0; xpix<2; xpix++) {
for (int n=0; n<numVertices; n++) {
double theta = theta0 + (n + 1.0) * dtheta;
poly.add(Position<double>(xpix, (std::tan(theta) + 1.0) / 2.0));
}
}
dbg<<"y edges:\n";
for (int ypix=0; ypix<2; ypix++) {
for (int n=0; n<numVertices; n++) {
double theta = theta0 + (n + 1.0) * dtheta;
poly.add(Position<double>((std::tan(theta) + 1.0) / 2.0, ypix));
}
}
poly.sort();
}
Silicon::Silicon(int numVertices, double numElec, int nx, int ny, int qDist,
double diffStep, double pixelSize,
double sensorThickness, double* vertex_data,
const Table& tr_radial_table, Position<double> treeRingCenter,
const Table& abs_length_table, bool transpose) :
_numVertices(numVertices), _nx(nx), _ny(ny), _qDist(qDist),
_diffStep(diffStep), _pixelSize(pixelSize),
_sensorThickness(sensorThickness),
_tr_radial_table(tr_radial_table), _treeRingCenter(treeRingCenter),
_abs_length_table(abs_length_table), _transpose(transpose),
_targetData(nullptr)
{
dbg<<"Silicon constructor\n";
// This constructor reads in the distorted pixel shapes from the Poisson solver
// and builds arrays of points for calculating the distorted pixel shapes
// as a function of charge in the surrounding pixels.
// First build the distorted points. We have linear boundary arrays,
// an undistorted polygon, and a polygon for test.
int nv1 = 4 * _numVertices + 4; // Number of vertices in each pixel in input file
_nv = 4 * _numVertices + 8; // Number of vertices in each pixel
dbg<<"_numVertices = "<<_numVertices<<", _nv = "<<_nv<<std::endl;
dbg<<"nx,ny = "<<nx<<", "<<ny<<" ntot = "<<nx*ny<<std::endl;
dbg<<"total memory = "<<nx*ny*_nv*sizeof(Position<float>)/(1024.*1024.)<<" MBytes"<<std::endl;
buildEmptyPoly(_emptypoly, _numVertices);
// Next, we read in the pixel distortions from the Poisson_CCD simulations
if (_transpose) std::swap(_nx,_ny);
_horizontalDistortions.resize(horizontalRowStride(_nx) * (_ny + 1));
_verticalDistortions.resize(verticalColumnStride(_ny) * (_nx + 1));
_horizontalDistortions.shrink_to_fit();
_verticalDistortions.shrink_to_fit();
for (int index=0; index < nv1*_nx*_ny; index++) {
int n1 = index % nv1;
int j = (index / nv1) % _ny;
int i = index / (nv1 * _ny);
xdbg<<"index = "<<index<<std::endl;
xdbg<<"i,j = "<<i<<','<<j<<std::endl;
xassert(index == (i * _ny + j) * nv1 + n1);
double x0 = vertex_data[5*index+0];
double y0 = vertex_data[5*index+1];
double x1 = vertex_data[5*index+3];
double y1 = vertex_data[5*index+4];
if (_transpose) {
xdbg<<"Original i,j,n = "<<i<<','<<j<<','<<n1<<std::endl;
std::swap(i,j);
std::swap(x0,y0);
std::swap(x1,y1);
n1 = (_numVertices - n1 + nv1) % nv1;
xdbg<<" => "<<i<<','<<j<<','<<n1<<std::endl;
}
// Figure out the new index around the pixel polygon when there are two corner points.
int n = n1;
if (n >= cornerIndexBottomLeft()) ++n;
if (n > cornerIndexBottomRight()) ++n;
if (n >= cornerIndexTopRight()) ++n;
if (n > cornerIndexTopLeft()) ++n;
xdbg<<"n1 = "<<n1<<", n = "<<n<<std::endl;
// The following captures the pixel displacement. These are translated into
// coordinates compatible with (x,y). These are per electron.
double x = _emptypoly[n].x;
x = ((x1 - x0) / _pixelSize + 0.5 - x) / numElec;
double y = _emptypoly[n].y;
y = ((y1 - y0) / _pixelSize + 0.5 - y) / numElec;
// populate the linear distortions arrays
// make sure to always use values from closest to center pixel
if ((((n < cornerIndexBottomLeft()) || (n > cornerIndexTopLeft())) && (i <= (_nx / 2))) || // LHS
(((n > cornerIndexBottomRight()) && (n < cornerIndexTopRight())) && (i >= (_nx / 2))) || // RHS
(((n >= cornerIndexBottomLeft()) && (n <= cornerIndexBottomRight())) && (j <= (_ny / 2))) || // bottom
(((n >= cornerIndexTopRight()) && (n <= cornerIndexTopLeft())) && (j >= (_ny / 2)))) { // top
bool horiz = false;
int bidx = getBoundaryIndex(i, j, n, &horiz);
if (horiz) {
_horizontalDistortions[bidx].x = x;
_horizontalDistortions[bidx].y = y;
}
else {
_verticalDistortions[bidx].x = x;
_verticalDistortions[bidx].y = y;
}
}
// If this is a corner pixel, we need to increment n and add it again.
// For the old method (using _distortions), this just duplicates the point.
// But in the new method, this will do different things in the horizontal and
// vertical directions to make sure we get both possible corner locations
// in the two directions if appropriate.
bool corner = ((n == cornerIndexBottomLeft() - 1) ||
(n == cornerIndexBottomRight()) ||
(n == cornerIndexTopRight() - 1) ||
(n == cornerIndexTopLeft()));
if (corner) {
// Increment n to the next location around the polygon.
++n;
xdbg<<"corner. n => "<<n<<std::endl;
// Do all the same stuff as above again. (Could consider pulling this little
// section out into a function that we call twice.)
if ((((n < cornerIndexBottomLeft()) || (n > cornerIndexTopLeft())) && (i <= (_nx / 2))) || // LHS
(((n > cornerIndexBottomRight()) && (n < cornerIndexTopRight())) && (i >= (_nx / 2))) || // RHS
(((n >= cornerIndexBottomLeft()) && (n <= cornerIndexBottomRight())) && (j <= (_ny / 2))) || // bottom
(((n >= cornerIndexTopRight()) && (n <= cornerIndexTopLeft())) && (j >= (_ny / 2)))) { // top
bool horiz = false;
int bidx = getBoundaryIndex(i, j, n, &horiz);
if (horiz) {
_horizontalDistortions[bidx].x = x;
_horizontalDistortions[bidx].y = y;
}
else {
_verticalDistortions[bidx].x = x;
_verticalDistortions[bidx].y = y;
}
}
}
}
// Process _abs_length_table and _emptypoly ready for GPU
// this will only be fully accurate for cases where the table uses linear
// interpolation, and the data points are evenly spaced. Currently this is
// always the case for _abs_length_table.
_abs_length_arg_min = _abs_length_table.argMin();
_abs_length_arg_max = _abs_length_table.argMax();
_abs_length_size = _abs_length_table.size();
_abs_length_table_GPU.resize(_abs_length_size);
_abs_length_increment = (_abs_length_arg_max - _abs_length_arg_min) /
(double)(_abs_length_size - 1);
for (int i = 0; i < _abs_length_size; i++) {
_abs_length_table_GPU[i] =
_abs_length_table.lookup(_abs_length_arg_min + (((double)i) * _abs_length_increment));
}
_emptypolyGPU.resize(_emptypoly.size());
for (int i=0; i<int(_emptypoly.size()); i++) {
_emptypolyGPU[i].x = _emptypoly[i].x;
_emptypolyGPU[i].y = _emptypoly[i].y;
}
}
Silicon::~Silicon()
{
if (_targetData != nullptr) {
finalize();
}
}
void Silicon::updatePixelBounds(int nx, int ny, size_t k,
Bounds<double>* pixelInnerBoundsData,
Bounds<double>* pixelOuterBoundsData,
Position<float>* horizontalBoundaryPointsData,
Position<float>* verticalBoundaryPointsData)
{
// update the bounding rectangles for pixel k
// get pixel co-ordinates
int x = k / ny;
int y = k % ny;
// compute outer bounds first
pixelOuterBoundsData[k] = Bounds<double>();
// iterate over pixel boundary
int n, idx;
// LHS lower half
for (n = 0; n < cornerIndexBottomLeft(); n++) {
idx = verticalPixelIndex(x, y, ny) + n + cornerIndexBottomLeft();
pixelOuterBoundsData[k] += verticalBoundaryPointsData[idx];
}
// bottom row including corners
for (; n <= cornerIndexBottomRight(); n++) {
idx = horizontalPixelIndex(x, y, nx) + (n - cornerIndexBottomLeft());
pixelOuterBoundsData[k] += horizontalBoundaryPointsData[idx];
}
// RHS
for (; n < cornerIndexTopRight(); n++) {
idx = verticalPixelIndex(x + 1, y, ny) + (cornerIndexTopRight() - n - 1);
Position<double> p = verticalBoundaryPointsData[idx];
p.x += 1.0;
pixelOuterBoundsData[k] += p;
}
// top row including corners
for (; n <= cornerIndexTopLeft(); n++) {
idx = horizontalPixelIndex(x, y + 1, nx) + (cornerIndexTopLeft() - n);
Position<double> p = horizontalBoundaryPointsData[idx];
p.y += 1.0;
pixelOuterBoundsData[k] += p;
}
// LHS upper half
for (; n < _nv; n++) {
idx = verticalPixelIndex(x, y, ny) + (n - cornerIndexTopLeft() - 1);
pixelOuterBoundsData[k] += verticalBoundaryPointsData[idx];
}
Position<double> center = pixelOuterBoundsData[k].center();
// compute inner bounds
// initialize inner from outer
double ibxmin = pixelOuterBoundsData[k].getXMin();
double ibxmax = pixelOuterBoundsData[k].getXMax();
double ibymin = pixelOuterBoundsData[k].getYMin();
double ibymax = pixelOuterBoundsData[k].getYMax();
// iterate over pixel boundary
// LHS lower half
for (n = 0; n < cornerIndexBottomLeft(); n++) {
idx = verticalPixelIndex(x, y, ny) + n + cornerIndexBottomLeft();
double px = verticalBoundaryPointsData[idx].x;
double py = verticalBoundaryPointsData[idx].y;
if (px-center.x >= std::abs(py-center.y) && px < ibxmax) ibxmax = px;
if (px-center.x <= -std::abs(py-center.y) && px > ibxmin) ibxmin = px;
if (py-center.y >= std::abs(px-center.x) && py < ibymax) ibymax = py;
if (py-center.y <= -std::abs(px-center.x) && py > ibymin) ibymin = py;
}
// bottom row including corners
for (; n <= cornerIndexBottomRight(); n++) {
idx = horizontalPixelIndex(x, y, nx) + (n - cornerIndexBottomLeft());
double px = horizontalBoundaryPointsData[idx].x;
double py = horizontalBoundaryPointsData[idx].y;
if (px-center.x >= std::abs(py-center.y) && px < ibxmax) ibxmax = px;
if (px-center.x <= -std::abs(py-center.y) && px > ibxmin) ibxmin = px;
if (py-center.y >= std::abs(px-center.x) && py < ibymax) ibymax = py;
if (py-center.y <= -std::abs(px-center.x) && py > ibymin) ibymin = py;
}
// RHS
for (; n < cornerIndexTopRight(); n++) {
idx = verticalPixelIndex(x + 1, y, ny) + (cornerIndexTopRight() - n - 1);
double px = verticalBoundaryPointsData[idx].x + 1.0;
double py = verticalBoundaryPointsData[idx].y;
if (px-center.x >= std::abs(py-center.y) && px < ibxmax) ibxmax = px;
if (px-center.x <= -std::abs(py-center.y) && px > ibxmin) ibxmin = px;
if (py-center.y >= std::abs(px-center.x) && py < ibymax) ibymax = py;
if (py-center.y <= -std::abs(px-center.x) && py > ibymin) ibymin = py;
}
// top row including corners
for (; n <= cornerIndexTopLeft(); n++) {
idx = horizontalPixelIndex(x, y + 1, nx) + (cornerIndexTopLeft() - n);
double px = horizontalBoundaryPointsData[idx].x;
double py = horizontalBoundaryPointsData[idx].y + 1.0;
if (px-center.x >= std::abs(py-center.y) && px < ibxmax) ibxmax = px;
if (px-center.x <= -std::abs(py-center.y) && px > ibxmin) ibxmin = px;
if (py-center.y >= std::abs(px-center.x) && py < ibymax) ibymax = py;
if (py-center.y <= -std::abs(px-center.x) && py > ibymin) ibymin = py;
}
// LHS upper half
for (; n < _nv; n++) {
idx = verticalPixelIndex(x, y, ny) + (n - cornerIndexTopLeft() - 1);
double px = verticalBoundaryPointsData[idx].x;
double py = verticalBoundaryPointsData[idx].y;
if (px-center.x >= std::abs(py-center.y) && px < ibxmax) ibxmax = px;
if (px-center.x <= -std::abs(py-center.y) && px > ibxmin) ibxmin = px;
if (py-center.y >= std::abs(px-center.x) && py < ibymax) ibymax = py;
if (py-center.y <= -std::abs(px-center.x) && py > ibymin) ibymin = py;
}
// store results in actual bound structure
pixelInnerBoundsData[k].setXMin(ibxmin);
pixelInnerBoundsData[k].setXMax(ibxmax);
pixelInnerBoundsData[k].setYMin(ibymin);
pixelInnerBoundsData[k].setYMax(ibymax);
}
template <typename T>
void Silicon::updatePixelDistortions(ImageView<T> target)
{
dbg<<"updatePixelDistortions\n";
// This updates the pixel distortions in the linear boundary arrays
// based on the amount of additional charge in each pixel
// This distortion assumes the electron is created at the
// top of the silicon. It mus be scaled based on the conversion depth
// This is handled in insidePixel.
int nxCenter = (_nx - 1) / 2;
int nyCenter = (_ny - 1) / 2;
// Now add in the displacements
const int nx = target.getNCol();
const int ny = target.getNRow();
const int step = target.getStep();
const int stride = target.getStride();
T* targetData = target.getData();
const int npix = nx * ny;
Position<float>* horizontalBoundaryPointsData = _horizontalBoundaryPoints.data();
Position<float>* verticalBoundaryPointsData = _verticalBoundaryPoints.data();
Position<float>* horizontalDistortionsData = _horizontalDistortions.data();
Position<float>* verticalDistortionsData = _verticalDistortions.data();
bool* changedData = _changed.get();
// Loop through the boundary arrays and update any points affected by nearby pixels
// Horizontal array first
// map image data and changed array throughout all GPU loops
#ifdef _OPENMP
#ifndef GALSIM_USE_GPU
#pragma omp parallel for
#else
#pragma omp target teams distribute parallel for
#endif
#endif
for (int p=0; p < npix; p++) {
// Calculate which pixel we are currently below
int x = p % nx;
int y = p / nx;
// Loop over rectangle of pixels that could affect this row of points
int polyi1 = imax(x - _qDist, 0);
int polyi2 = imin(x + _qDist, nx - 1);
// NB. We are working between rows y and y-1, so need polyj1 = y-1 - _qDist.
int polyj1 = imax(y - (_qDist + 1), 0);
int polyj2 = imin(y + _qDist, ny - 1);
bool change = false;
for (int j=polyj1; j <= polyj2; j++) {
for (int i=polyi1; i <= polyi2; i++) {
// Check whether this pixel has charge on it
double charge = targetData[(j * stride) + (i * step)];
if (charge != 0.0) {
change = true;
// Work out corresponding index into distortions array
int dist_index = (((y - j + nyCenter) * _nx) + (x - i + nxCenter)) * horizontalPixelStride();
int index = p * horizontalPixelStride();
// Loop over boundary points and update them
for (int n=0; n < horizontalPixelStride(); ++n, ++index, ++dist_index) {
horizontalBoundaryPointsData[index].x =
double(horizontalBoundaryPointsData[index].x) +
horizontalDistortionsData[dist_index].x * charge;
horizontalBoundaryPointsData[index].y =
double(horizontalBoundaryPointsData[index].y) +
horizontalDistortionsData[dist_index].y * charge;
}
}
}
}
// update changed array
if (change) {
if (y < ny) changedData[(x * ny) + y] = true; // pixel above
if (y > 0) changedData[(x * ny) + (y - 1)] = true; // pixel below
}
}
// Now vertical array
#ifdef _OPENMP
#ifndef GALSIM_USE_GPU
#pragma omp parallel for
#else
#pragma omp target teams distribute parallel for
#endif
#endif
for (int p=0; p < (nx * ny); p++) {
// Calculate which pixel we are currently on
int x = p / ny;
int y = (ny - 1) - (p % ny); // remember vertical points run top-to-bottom
// Loop over rectangle of pixels that could affect this column of points
int polyi1 = imax(x - (_qDist + 1), 0);
int polyi2 = imin(x + _qDist, nx - 1);
int polyj1 = imax(y - _qDist, 0);
int polyj2 = imin(y + _qDist, ny - 1);
bool change = false;
for (int j=polyj1; j <= polyj2; j++) {
for (int i=polyi1; i <= polyi2; i++) {
// Check whether this pixel has charge on it
double charge = targetData[(j * stride) + (i * step)];
if (charge != 0.0) {
change = true;
// Work out corresponding index into distortions array
int dist_index = (((x - i + nxCenter) * _ny) + ((_ny - 1) - (y - j + nyCenter))) * verticalPixelStride() + (verticalPixelStride() - 1);
int index = p * verticalPixelStride() + (verticalPixelStride() - 1);
// Loop over boundary points and update them
for (int n=0; n < verticalPixelStride(); ++n, --index, --dist_index) {
verticalBoundaryPointsData[index].x =
double(verticalBoundaryPointsData[index].x) +
verticalDistortionsData[dist_index].x * charge;
verticalBoundaryPointsData[index].y =
double(verticalBoundaryPointsData[index].y) +
verticalDistortionsData[dist_index].y * charge;
}
}
}
}
// update changed array
if (change) {
if (x < nx) changedData[(x * ny) + y] = true;
if (x > 0) changedData[((x - 1) * ny) + y] = true;
}
}
Bounds<double>* pixelInnerBoundsData = _pixelInnerBounds.data();
Bounds<double>* pixelOuterBoundsData = _pixelOuterBounds.data();
#ifdef _OPENMP
#ifndef GALSIM_USE_GPU
#pragma omp parallel for
#else
#pragma omp target teams distribute parallel for
#endif
#endif
for (int k=0; k<npix; ++k) {
if (changedData[k]) {
updatePixelBounds(nx, ny, k, pixelInnerBoundsData,
pixelOuterBoundsData,
horizontalBoundaryPointsData,
verticalBoundaryPointsData);
changedData[k] = false;
}
}
}
// This version of calculateTreeRingDistortion only distorts a polygon.
// Used in the no-flux pixel area calculation.
void Silicon::calculateTreeRingDistortion(int i, int j, Position<int> orig_center,
Polygon& poly) const
{
for (int n=0; n<_nv; n++) {
xdbg<<"i,j,n = "<<i<<','<<j<<','<<n<<": x,y = "<<
poly[n].x <<" "<< poly[n].y<<std::endl;
double tx = (double)i + poly[n].x - _treeRingCenter.x + (double)orig_center.x;
double ty = (double)j + poly[n].y - _treeRingCenter.y + (double)orig_center.y;
xdbg<<"tx,ty = "<<tx<<','<<ty<<std::endl;
double r = sqrt(tx * tx + ty * ty);
if (r > 0 && r < _tr_radial_table.argMax()) {
double shift = _tr_radial_table.lookup(r);
xdbg<<"r = "<<r<<", shift = "<<shift<<std::endl;
// Shifts are along the radial vector in direction of the doping gradient
double dx = shift * tx / r;
double dy = shift * ty / r;
xdbg<<"dx,dy = "<<dx<<','<<dy<<std::endl;
poly[n].x = double(poly[n].x) + dx;
poly[n].y = double(poly[n].y) + dy;
xdbg<<" x,y => "<<poly[n].x <<" "<< poly[n].y;
}
}
}
// This version updates the linear boundary
void Silicon::calculateTreeRingDistortion(int i, int j, Position<int> orig_center,
int nx, int ny, int i1, int j1)
{
iteratePixelBoundary(
i-i1, j-j1, nx, ny, [&](int n, Position<float>& pt, bool rhs, bool top) {
Position<double> p = pt;
// only do bottom and left points unless we're on top/right edge
if ((rhs) && ((i - i1) < (nx - 1))) return;
if ((top) && ((j - j1) < (ny - 1))) return;
if (rhs) p.x += 1.0;
if (top) p.y += 1.0;
//xdbg<<"x,y = "<<p.x<<','<<p.y<<std::endl;
double tx = (double)i + p.x - _treeRingCenter.x + (double)orig_center.x;
double ty = (double)j + p.y - _treeRingCenter.y + (double)orig_center.y;
//xdbg<<"tx,ty = "<<tx<<','<<ty<<std::endl;
double r = sqrt(tx * tx + ty * ty);
if (r > 0 && r < _tr_radial_table.argMax()) {
double shift = _tr_radial_table.lookup(r);
//xdbg<<"r = "<<r<<", shift = "<<shift<<std::endl;
// Shifts are along the radial vector in direction of the doping gradient
double dx = shift * tx / r;
double dy = shift * ty / r;
//xdbg<<"dx,dy = "<<dx<<','<<dy<<std::endl;
pt.x += dx;
pt.y += dy;
}
}
);
}
template <typename T>
void Silicon::addTreeRingDistortions(ImageView<T> target, Position<int> orig_center)
{
if (_tr_radial_table.size() == 2) {
//dbg<<"Trivial radial table\n";
// The no tree rings case is indicated with a table of size 2, which
// wouldn't make any sense as a user input.
return;
}
dbg<<"addTreeRings\n";
// This updates the pixel distortions in the _imagepolys
// pixel list based on a model of tree rings.
// The coordinates _treeRingCenter are the coordinates
// of the tree ring center, shifted to compensate for the
// fact that target has its origin shifted to (0,0).
Bounds<int> b = target.getBounds();
const int i1 = b.getXMin();
const int i2 = b.getXMax();
const int j1 = b.getYMin();
const int j2 = b.getYMax();
const int nx = target.getNCol();
const int ny = target.getNRow();
// Now we cycle through the pixels in the target image and add
// the (small) distortions due to tree rings
std::vector<bool> changed(nx * ny, false);
for (int i=i1; i<=i2; ++i) {
for (int j=j1; j<=j2; ++j) {
int index = (i - i1) * ny + (j - j1);
calculateTreeRingDistortion(i, j, orig_center, nx, ny, i1, j1);
changed[index] = true;
}
}
for (size_t k=0; k<changed.size(); ++k) {
if (changed[k]) {
updatePixelBounds(nx, ny, k, _pixelInnerBounds.data(),
_pixelOuterBounds.data(),
_horizontalBoundaryPoints.data(),
_verticalBoundaryPoints.data());
}
}
}
// Scales a linear pixel boundary into a polygon object.
void Silicon::scaleBoundsToPoly(int i, int j, int nx, int ny,
const Polygon& emptypoly, Polygon& result,
double factor) const
{
result = emptypoly;
iteratePixelBoundary(
i, j, nx, ny,
[&](int n, const Position<float>& pt, bool rhs, bool top) {
Position<double> p = pt;
if (rhs) p.x += 1.0;
if (top) p.y += 1.0;
result[n].x += (p.x - emptypoly[n].x) * factor;
result[n].y += (p.y - emptypoly[n].y) * factor;
}
);
result.updateBounds();
}
bool Silicon::insidePixel(int ix, int iy, double x, double y, double zconv,
Bounds<int>& targetBounds, bool* off_edge,
int emptypolySize,
Bounds<double>* pixelInnerBoundsData,
Bounds<double>* pixelOuterBoundsData,
Position<float>* horizontalBoundaryPointsData,
Position<float>* verticalBoundaryPointsData,
Position<double>* emptypolyData) const
{
// This scales the pixel distortion based on the zconv, which is the depth
// at which the electron is created, and then tests to see if the delivered
// point is inside the pixel.
// (ix,iy) is the pixel being tested, and (x,y) is the coordinate of the
// photon within the pixel, with (0,0) in the lower left
// If test pixel is off the image, return false. (Avoids seg faults!)
if (!targetBounds.includes(ix, iy)) {
if (off_edge) *off_edge = true;
return false;
}
const int i1 = targetBounds.getXMin();
const int i2 = targetBounds.getXMax();
const int j1 = targetBounds.getYMin();
const int j2 = targetBounds.getYMax();
const int nx = i2-i1+1;
const int ny = j2-j1+1;
int index = (ix - i1) * ny + (iy - j1);
// First do some easy checks if the point isn't terribly close to the boundary.
bool inside;
if (pixelInnerBoundsData[index].includes(x, y)) {
inside = true;
} else if (!pixelOuterBoundsData[index].includes(x, y)) {
inside = false;
} else {
// OK, it must be near the boundary, so now be careful.
// The term zfactor decreases the pixel shifts as we get closer to the bottom
// It is an empirical fit to the Poisson solver simulations, and only matters
// when we get quite close to the bottom. This could be more accurate by making
// the Vertices files have an additional look-up variable (z), but this doesn't
// seem necessary at this point
const double zfit = 12.0;
const double zfactor = std::tanh(zconv / zfit);
#if 0
// Old version that used a temporary polygon per thread (_testpoly)
#ifdef _OPENMP
int t = omp_get_thread_num();
#else
int t = 0;
#endif
// Scale the testpoly vertices by zfactor
scaleBoundsToPoly(ix - i1, iy - j1, nx, ny, _emptypoly, _testpoly[t],
zfactor);
// Now test to see if the point is inside
Position<double> p(x, y);
inside = _testpoly[t].contains(p);
#else
// New version that doesn't use a temporary polygon object
// This is required for GPU as due to the high number of threads,
// having a temporary polygon per thread is not practical
// compute first point of polygon (index 0)
double x1=0, y1=0, xfirst, yfirst, xinters = 0.0;
inside = false;
for (int n = 0; n <= _nv; n++) {
double xp = 0.0, yp = 0.0;
double epx = 0.0, epy = 0.0;
if (n < _nv) {
epx = emptypolyData[n].x;
epy = emptypolyData[n].y;
}
xp = epx;
yp = epy;
int idx;
// compute this point
if (n < cornerIndexBottomLeft()) {
idx = verticalPixelIndex(ix - i1, iy - j1, ny) + n + cornerIndexBottomLeft();
xp += (verticalBoundaryPointsData[idx].x - epx) * zfactor;
yp += (verticalBoundaryPointsData[idx].y - epy) * zfactor;
}
else if (n <= cornerIndexBottomRight()) {
// bottom row including corners
idx = horizontalPixelIndex(ix - i1, iy - j1, nx) + (n - cornerIndexBottomLeft());
double px = horizontalBoundaryPointsData[idx].x;
xp += (px - epx) * zfactor;
yp += (horizontalBoundaryPointsData[idx].y - epy) * zfactor;
}
// RHS
else if (n < cornerIndexTopRight()) {
idx = verticalPixelIndex(ix - i1 + 1, iy - j1, ny) + (cornerIndexTopRight() - n - 1);
xp += ((verticalBoundaryPointsData[idx].x + 1.0) - epx) * zfactor;
yp += (verticalBoundaryPointsData[idx].y - epy) * zfactor;
}
// top row including corners
else if (n <= cornerIndexTopLeft()) {
idx = horizontalPixelIndex(ix - i1, iy - j1 + 1, nx) + (cornerIndexTopLeft() - n);
double px = horizontalBoundaryPointsData[idx].x;
xp += (px - epx) * zfactor;
yp += ((horizontalBoundaryPointsData[idx].y + 1.0) - epy) * zfactor;
}
else if (n < _nv) {
// LHS upper half
idx = verticalPixelIndex(ix - i1, iy - j1, ny) + (n - cornerIndexTopLeft() - 1);
xp += (verticalBoundaryPointsData[idx].x - epx) * zfactor;
yp += (verticalBoundaryPointsData[idx].y - epy) * zfactor;
}
if (n == 0) {
// save first point for later
x1 = xp;
y1 = yp;
xfirst = xp;
yfirst = yp;
}
else {
// shoelace algorithm
double x2 = xp;
double y2 = yp;
if (n == _nv) {
x2 = xfirst;
y2 = yfirst;
}
double ymin = y1 < y2 ? y1 : y2;
double ymax = y1 > y2 ? y1 : y2;
double xmax = x1 > x2 ? x1 : x2;
if (y > ymin) {
if (y <= ymax) {
if (x <= xmax) {
if (y1 != y2) {
xinters = (y - y1) * (x2 - x1) / (y2 - y1) + x1;
}
if ((x1 == x2) || (x <= xinters)) {
inside = !inside;
}
}
}
}
x1 = x2;
y1 = y2;
}
}
#endif
}
// If the nominal pixel is on the edge of the image and the photon misses in the
// direction of falling off the image, (possibly) report that in off_edge.
if (!inside && off_edge) {
*off_edge = false;
if ((ix == i1) && (x < pixelInnerBoundsData[index].getXMin())) *off_edge = true;
if ((ix == i2) && (x > pixelInnerBoundsData[index].getXMax())) *off_edge = true;
if ((iy == j1) && (y < pixelInnerBoundsData[index].getYMin())) *off_edge = true;
if ((iy == j2) && (y > pixelInnerBoundsData[index].getYMax())) *off_edge = true;
}
return inside;
}
// Helper function to calculate how far down into the silicon the photon converts into
// an electron.
double Silicon::calculateConversionDepth(bool photonsHasAllocatedWavelengths,
const double* photonsWavelength,
const double* abs_length_table_data,
bool photonsHasAllocatedAngles,
const double* photonsDXDZ,
const double* photonsDYDZ, int i,
double randomNumber) const
{
// Determine the distance the photon travels into the silicon
double si_length;
if (photonsHasAllocatedWavelengths) {
double lambda = photonsWavelength[i]; // in nm
// Lookup the absorption length in the imported table
// perform abs_length_table lookup with linear interpolation
int tableIdx = int((lambda - _abs_length_arg_min) / _abs_length_increment);
double alpha = (lambda - ((float(tableIdx) * _abs_length_increment) +
_abs_length_arg_min)) / _abs_length_increment;
int tableIdx1 = tableIdx + 1;
if (tableIdx <= 0) tableIdx = tableIdx1 = 0;
if (tableIdx >= _abs_length_size-1) tableIdx = tableIdx1 = _abs_length_size - 1;
double abs_length = (abs_length_table_data[tableIdx] * (1.0 - alpha)) +
(abs_length_table_data[tableIdx1] * alpha); // in microns
si_length = -abs_length * log(1.0 - randomNumber); // in microns
#ifdef DEBUGLOGGING
if (i % 1000 == 0) {
xdbg<<"lambda = "<<lambda<<std::endl;
xdbg<<"si_length = "<<si_length<<std::endl;
}
#endif
} else {
// If no wavelength info, assume conversion takes place near the top.
si_length = 1.0;
}
// Next we partition the si_length into x,y,z. Assuming dz is positive downward
if (photonsHasAllocatedAngles) {
double dxdz = photonsDXDZ[i];
double dydz = photonsDYDZ[i];
double dz = si_length / std::sqrt(1.0 + dxdz*dxdz + dydz*dydz); // in microns
#ifdef DEBUGLOGGING
if (i % 1000 == 0) {
xdbg<<"dxdz = "<<dxdz<<std::endl;
xdbg<<"dydz = "<<dydz<<std::endl;
xdbg<<"dz = "<<dz<<std::endl;
}
#endif
return std::fmin(_sensorThickness - 1.0, dz); // max 1 micron from bottom
} else {
return si_length;
}
}
bool searchNeighbors(const Silicon& silicon, int& ix, int& iy, double x, double y, double zconv,
Bounds<int>& targetBounds, int& step, int emptypolysize,
Bounds<double>* pixelInnerBoundsData,
Bounds<double>* pixelOuterBoundsData,
Position<float>* horizontalBoundaryPointsData,
Position<float>* verticalBoundaryPointsData,
Position<double>* emptypolyData)
{
const int xoff[9] = {0,1,1,0,-1,-1,-1,0,1}; // Displacements to neighboring pixels
const int yoff[9] = {0,0,1,1,1,0,-1,-1,-1}; // Displacements to neighboring pixels
// The following code finds which pixel we are in given
// pixel distortion due to the brighter-fatter effect
// The following are set up to start the search in the undistorted pixel, then
// search in the nearest neighbor first if it's not in the undistorted pixel.
if ((x > y) && (x > 1.0 - y)) step = 1;
else if ((x < y) && (x < 1.0 - y)) step = 7;
else if ((x < y) && (x > 1.0 - y)) step = 3;
else step = 5;
int n=step;
for (int m=1; m<9; m++) {
int ix_off = ix + xoff[n];
int iy_off = iy + yoff[n];
double x_off = x - xoff[n];
double y_off = y - yoff[n];
if (silicon.insidePixel(ix_off, iy_off, x_off, y_off, zconv,
targetBounds, nullptr, emptypolysize,
pixelInnerBoundsData, pixelOuterBoundsData,
horizontalBoundaryPointsData,
verticalBoundaryPointsData,
emptypolyData)) {
ix = ix_off;
iy = iy_off;
return true;
}
n = ((n-1) + step) % 8 + 1;
// This is intended to start with the nearest neighbor, then cycle through others.
}
return false;
}
// Calculates the area of a pixel based on the linear boundaries.
double Silicon::pixelArea(int i, int j, int nx, int ny) const
{
double area = 0.0;
bool horizontal1, horizontal2;
// compute sum of triangle areas using cross-product rule (shoelace formula)
for (int n = 0; n < _nv; n++) {
int pi1 = getBoundaryIndex(i, j, n, &horizontal1, nx, ny);
Position<double> p1 = horizontal1 ? _horizontalBoundaryPoints[pi1] :
_verticalBoundaryPoints[pi1];
if ((n > cornerIndexBottomRight()) && (n < cornerIndexTopRight())) p1.x += 1.0;
if ((n >= cornerIndexTopRight()) && (n <= cornerIndexTopLeft())) p1.y += 1.0;
int n2 = (n + 1) % _nv;
int pi2 = getBoundaryIndex(i, j, n2, &horizontal2, nx, ny);
Position<double> p2 = horizontal2 ? _horizontalBoundaryPoints[pi2] :
_verticalBoundaryPoints[pi2];
if ((n2 > cornerIndexBottomRight()) && (n2 < cornerIndexTopRight())) p2.x += 1.0;
if ((n2 >= cornerIndexTopRight()) && (n2 <= cornerIndexTopLeft())) p2.y += 1.0;
area += p1.x * p2.y;
area -= p2.x * p1.y;
}
return std::abs(area) / 2.0;
}
template <typename T>
void Silicon::fillWithPixelAreas(ImageView<T> target, Position<int> orig_center,
bool use_flux)
{
Bounds<int> b = target.getBounds();
if (!b.isDefined())
throw std::runtime_error("Attempting to PhotonArray::addTo an Image with"
" undefined Bounds");
const int i1 = b.getXMin();
const int i2 = b.getXMax();
const int j1 = b.getYMin();
const int j2 = b.getYMax();
const int nx = target.getNCol();
const int ny = target.getNRow();
const int npix = nx * ny;
if (use_flux) {
dbg<<"Start full pixel area calculation\n";
dbg<<"nx,ny = "<<nx<<','<<ny<<std::endl;
dbg<<"total memory = "<<npix*_nv*sizeof(Position<float>)/(1024.*1024.)<<" MBytes"<<std::endl;
// This will add distortions according to the current flux in the image, on the
// GPU where appropriate.
initialize(target, orig_center);
// Copy the distorted pixel boundaries from GPU back to CPU if necessary.
#ifdef _OPENMP
#ifdef GALSIM_USE_GPU
Position<float>* horizontalBoundaryPointsData = _horizontalBoundaryPoints.data();
Position<float>* verticalBoundaryPointsData = _verticalBoundaryPoints.data();
int hbpSize = _horizontalBoundaryPoints.size();
int vbpSize = _verticalBoundaryPoints.size();
#pragma omp target update from(horizontalBoundaryPointsData[0:hbpSize], verticalBoundaryPointsData[0:vbpSize])
#endif
#endif
// Fill target with the area in each pixel.
const int skip = target.getNSkip();
const int step = target.getStep();
T* ptr = target.getData();
for (int j=j1; j<=j2; ++j, ptr+=skip) {
for (int i=i1; i<=i2; ++i, ptr+=step) {
double newArea = pixelArea(i - i1, j - j1, nx, ny);
*ptr = newArea;
}
}
} else {
// If we don't care about respecting the flux in the image (which we usually don't
// since this is generally used for making sky images with the tree ring patterns),