-
Notifications
You must be signed in to change notification settings - Fork 110
/
Copy pathSBInterpolatedImage.cpp
1475 lines (1309 loc) · 60.9 KB
/
SBInterpolatedImage.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.
*/
//#define DEBUGLOGGING
#include <vector>
#include <map>
#include <algorithm>
#include <cstring> // For memset
#ifdef __SSE2__
#include "xmmintrin.h"
#endif
#include "SBInterpolatedImage.h"
#include "SBInterpolatedImageImpl.h"
#include "Std.h"
namespace galsim {
///////////////////////////////////////////////////////////////////////////////////////////////
// SBInterpolatedImage methods
SBInterpolatedImage::SBInterpolatedImage(
const BaseImage<double>& image,
const Bounds<int>& init_bounds, const Bounds<int>& nonzero_bounds,
const Interpolant& xInterp, const Interpolant& kInterp,
double stepk, double maxk, const GSParams& gsparams) :
SBProfile(new SBInterpolatedImageImpl(
image, init_bounds, nonzero_bounds, xInterp, kInterp, stepk, maxk, gsparams)) {}
SBInterpolatedImage::SBInterpolatedImage(const SBInterpolatedImage& rhs) : SBProfile(rhs) {}
SBInterpolatedImage::~SBInterpolatedImage() {}
const Interpolant& SBInterpolatedImage::getXInterp() const
{
assert(dynamic_cast<const SBInterpolatedImageImpl*>(_pimpl.get()));
return static_cast<const SBInterpolatedImageImpl&>(*_pimpl).getXInterp();
}
const Interpolant& SBInterpolatedImage::getKInterp() const
{
assert(dynamic_cast<const SBInterpolatedImageImpl*>(_pimpl.get()));
return static_cast<const SBInterpolatedImageImpl&>(*_pimpl).getKInterp();
}
void SBInterpolatedImage::calculateStepK(double max_stepk) const
{
assert(dynamic_cast<const SBInterpolatedImageImpl*>(_pimpl.get()));
return static_cast<const SBInterpolatedImageImpl&>(*_pimpl).calculateStepK(max_stepk);
}
void SBInterpolatedImage::calculateMaxK(double max_maxk) const
{
assert(dynamic_cast<const SBInterpolatedImageImpl*>(_pimpl.get()));
return static_cast<const SBInterpolatedImageImpl&>(*_pimpl).calculateMaxK(max_maxk);
}
ConstImageView<double> SBInterpolatedImage::getPaddedImage() const
{
assert(dynamic_cast<const SBInterpolatedImageImpl*>(_pimpl.get()));
return static_cast<const SBInterpolatedImageImpl&>(*_pimpl).getPaddedImage();
}
ConstImageView<double> SBInterpolatedImage::getNonZeroImage() const
{
assert(dynamic_cast<const SBInterpolatedImageImpl*>(_pimpl.get()));
return static_cast<const SBInterpolatedImageImpl&>(*_pimpl).getNonZeroImage();
}
ConstImageView<double> SBInterpolatedImage::getImage() const
{
assert(dynamic_cast<const SBInterpolatedImageImpl*>(_pimpl.get()));
return static_cast<const SBInterpolatedImageImpl&>(*_pimpl).getImage();
}
///////////////////////////////////////////////////////////////////////////////////////////////
// SBInterpolatedImageImpl methods
#define INVALID -1.e300 // dummy value to indicate that flux or centroid not calculated yet.
SBInterpolatedImage::SBInterpolatedImageImpl::SBInterpolatedImageImpl(
const BaseImage<double>& image,
const Bounds<int>& init_bounds, const Bounds<int>& nonzero_bounds,
const Interpolant& xInterp, const Interpolant& kInterp,
double stepk, double maxk, const GSParams& gsparams) :
SBProfileImpl(gsparams),
_image(image.view()), _image_bounds(image.getBounds()),
_init_bounds(init_bounds), _nonzero_bounds(nonzero_bounds),
_xInterp(xInterp), _kInterp(kInterp),
_stepk(stepk), _maxk(maxk),
_flux(INVALID), _xcentroid(INVALID), _ycentroid(INVALID),
_readyToShoot(false)
{
dbg<<"image bounds = "<<image.getBounds()<<std::endl;
dbg<<"init bounds = "<<_init_bounds<<std::endl;
dbg<<"nonzero bounds = "<<_nonzero_bounds<<std::endl;
if (_stepk <= 0.) {
// Calculate stepK:
//
// The amount of flux missed in a circle of radius pi/stepk should be at
// most folding_threshold of the flux.
//
// We add the size of the image and the size of the interpolant in quadrature.
// (Note: Since this isn't a radial profile, R isn't really a radius, but rather
// the size of the square box that is enclosing all the flux.)
double R = std::max((_init_bounds.getXMax()-_init_bounds.getXMin())/2.,
(_init_bounds.getYMax()-_init_bounds.getYMin())/2.);
// Add xInterp range in quadrature just like convolution:
double R2 = _xInterp.xrange();
dbg<<"R(image) = "<<R<<", R(interpolant) = "<<R2<<std::endl;
R = sqrt(R*R + R2*R2);
dbg<<"=> R = "<<R<<std::endl;
_stepk = M_PI / R;
dbg<<"stepk = "<<_stepk<<std::endl;
}
_uscale = 1. / (2.*M_PI);
_maxk1 = _xInterp.urange()/_uscale;
if (_maxk <= 0.) {
// Calculate maxk:
//
// For now, just set this to where the interpolant's FT is <= maxk_threshold.
// Note: since we used kvalue_accuracy for the threshold of _xInterp
// (at least for the default quintic interpolant) rather than maxk_threshold,
// this will probably be larger than we really need.
// We could modify the urange method of Interpolant to take a threshold value
// at that point, rather than just use the constructor's value, but it's
// probably not worth it.
//
// In practice, we will generally call calculateMaxK() after construction to
// refine the value of maxk based on the actual FT of the image.
_maxk = _maxk1;
dbg<<"maxk = "<<_maxk<<std::endl;
}
}
SBInterpolatedImage::SBInterpolatedImageImpl::~SBInterpolatedImageImpl() {}
const Interpolant& SBInterpolatedImage::SBInterpolatedImageImpl::getXInterp() const
{ return _xInterp; }
const Interpolant& SBInterpolatedImage::SBInterpolatedImageImpl::getKInterp() const
{ return _kInterp; }
double SBInterpolatedImage::SBInterpolatedImageImpl::maxSB() const
{ return _image.maxAbsElement(); }
double SBInterpolatedImage::SBInterpolatedImageImpl::xValue(const Position<double>& pos) const
{
double x = pos.x;
double y = pos.y;
int p1, p2, q1, q2; // Range over which we need to sum.
const double SMALL = 10.*std::numeric_limits<double>::epsilon();
// If x or y are integers, only sum 1 element in that direction.
if (std::abs(x-std::floor(x+0.01)) < SMALL*(std::abs(x)+1)) {
p1 = p2 = int(std::floor(x+0.01));
} else {
p1 = int(std::ceil(x-_xInterp.xrange()));
p2 = int(std::floor(x+_xInterp.xrange()));
}
if (std::abs(y-std::floor(y+0.01)) < SMALL*(std::abs(y)+1)) {
q1 = q2 = int(std::floor(y+0.01));
} else {
q1 = int(std::ceil(y-_xInterp.xrange()));
q2 = int(std::floor(y+_xInterp.xrange()));
}
// If either range is not in non-zero part, then value is 0.
if (p2 < _nonzero_bounds.getXMin() ||
p1 > _nonzero_bounds.getXMax() ||
q2 < _nonzero_bounds.getYMin() ||
q1 > _nonzero_bounds.getYMax()) {
return 0.;
}
// Limit to nonzero region
if (p1 < _nonzero_bounds.getXMin()) p1 = _nonzero_bounds.getXMin();
if (p2 > _nonzero_bounds.getXMax()) p2 = _nonzero_bounds.getXMax();
if (q1 < _nonzero_bounds.getYMin()) q1 = _nonzero_bounds.getYMin();
if (q2 > _nonzero_bounds.getYMax()) q2 = _nonzero_bounds.getYMax();
// We'll need these for each row. Save them.
double xwt[p2-p1+1];
for (int p=p1, pp=0; p<=p2; ++p, ++pp) xwt[pp] = _xInterp.xval(p-x);
double sum = 0.;
for (int q=q1; q<=q2; ++q) {
double xsum = 0.;
for (int p=p1, pp=0; p<=p2; ++p, ++pp) {
xsum += xwt[pp] * _image(p,q);
}
sum += xsum * _xInterp.xval(q-y);
}
return sum;
}
int WrapKIndex(int k, int No2, int N)
{
k = (k + No2) % N;
if (k < 0) k += N;
return k - No2;
}
template <bool yn>
struct Maybe // true
{
template <typename T>
static inline void increment(T& p) { ++p; }
template <typename T>
static inline void increment(T& p, int n) { p += n; }
template <typename T>
static inline std::complex<T> conj(const std::complex<T>& x) { return std::conj(x); }
template <typename T, typename T2>
static inline T plus(const T& x, const T2& y) { return x+y; }
};
template <>
struct Maybe<false>
{
template <typename T>
static inline void increment(T& p) { --p; }
template <typename T>
static inline void increment(T& p, int n) { p -= n; }
template <typename T>
static inline std::complex<T> conj(const std::complex<T>& x) { return x; }
template <typename T, typename T2>
static inline T plus(const T& x, const T2& y) { return x-y; }
};
// A helper function for fast calculation of a dot product of real and complex vectors
template <bool c2>
static std::complex<double> ZDot(int n, const double* A, const std::complex<double>* B)
{
if (n) {
#ifdef __SSE2__
std::complex<double> sum(0);
while (n && !IsAligned(A) ) {
sum += *A * *B;
++A;
Maybe<!c2>::increment(B);
--n;
}
int n_2 = (n>>1);
int nb = n-(n_2<<1);
if (n_2) {
union { __m128d xm; double xd[2]; } xsum;
xsum.xm = _mm_set1_pd(0.);
__m128d xsum2 = _mm_set1_pd(0.);
const std::complex<double>* B1 = Maybe<!c2>::plus(B,1);
assert(IsAligned(A));
assert(IsAligned(B));
do {
const __m128d& xA = *(const __m128d*)(A);
const __m128d& xB1 = *(const __m128d*)(B);
const __m128d& xB2 = *(const __m128d*)(B1);
A += 2;
Maybe<!c2>::increment(B,2);
Maybe<!c2>::increment(B1,2);
__m128d xA1 = _mm_shuffle_pd(xA,xA,_MM_SHUFFLE2(0,0));
__m128d xA2 = _mm_shuffle_pd(xA,xA,_MM_SHUFFLE2(1,1));
__m128d x1 = _mm_mul_pd(xA1,xB1);
__m128d x2 = _mm_mul_pd(xA2,xB2);
xsum.xm = _mm_add_pd(xsum.xm,x1);
xsum2 = _mm_add_pd(xsum2,x2);
} while (--n_2);
xsum.xm = _mm_add_pd(xsum.xm,xsum2);
sum += std::complex<double>(xsum.xd[0],xsum.xd[1]);
}
if (nb) {
sum += *A * *B;
++A;
Maybe<!c2>::increment(B);
}
return Maybe<c2>::conj(sum);
#else
std::complex<double> sum = 0.;
do {
sum += *A * *B;
++A;
Maybe<!c2>::increment(B);
} while (--n);
return Maybe<c2>::conj(sum);
#endif
} else {
return 0.;
}
}
// This is the inner loop in all of the KValue calculations, including both the regular
// kValue method and both versions of fillKImage.
std::complex<double> KValueInnerLoop(int n, int p, int q, int No2, int N, double* xwt,
const BaseImage<std::complex<double> >& kimage)
{
std::complex<double> sum = 0.;
#if 0
// This is the more straightforward implementation of this calculation, which we
// preserve here for readability.
for (; n; --n, ++p) {
if (p == No2+1) p -= N;
// _kimage doesn't store p<0 half, so need to use the fact that
// _kimage(p,q) = conj(_kimage(-p,-q)) when p < 0.
if (p < 0)
if (q == -No2)
// N.B. _kimage(p,No2) == _kimage(p,-No2)
sum += *xwt++ * std::conj(kimage(-p,q));
else
sum += *xwt++ * std::conj(kimage(-p,-q));
else
sum += *xwt++ * kimage(p,q);
}
#else
// This version is more efficient, but a little obfuscated.
// It is equivalent in result to the above calculation.
// Note: q=No2 is not stored in the kimage, so when doing negative p values
// we need to use q on the conjugate side, not -q.
// Figure this out now, so we can ignore this subtlety below.
int mq = q == -No2 ? q : -q;
assert(kimage.getStep() == 1);
// First do any negative p values
if (p < 0) {
xdbg<<"Some initial negative p: p = "<<p<<std::endl;
int n1 = std::min(n, -p);
xdbg<<"n1 = "<<n1<<std::endl;
n -= n1;
const std::complex<double>* ptr = &kimage(-p,mq);
sum += ZDot<true>(n1, xwt, ptr);
xwt += n1;
p = 0;
}
// Next do positive p values:
if (n) {
xdbg<<"Positive p: p = "<<p<<std::endl;
const std::complex<double>* ptr = &kimage(p,q);
int n1 = std::min(n, No2+1-p);
xdbg<<"n1 = "<<n1<<std::endl;
n -= n1;
sum += ZDot<false>(n1, xwt, ptr);
xwt += n1;
}
// Finally if we've wrapped around again, do more negative p values:
if (n) {
xdbg<<"More negative p: p = "<<p<<std::endl;
int n1 = std::min(n, No2);
// Note: n1 is always n in practice, but this prevents pointer access going past
// edges of image if there is some unanticipated use case where it isn't.
xassert(n < No2);
xdbg<<"n1 = "<<n1<<std::endl;
p = -No2+1;
const std::complex<double>* ptr = &kimage(-p,mq);
sum += ZDot<true>(n1, xwt, ptr);
}
#endif
return sum;
}
std::complex<double> SBInterpolatedImage::SBInterpolatedImageImpl::kValue(
const Position<double>& kpos) const
{
double kx = kpos.x;
double ky = kpos.y;
dbg<<"evaluating kValue("<<kx<<","<<ky<<")"<<std::endl;
// Don't bother if the desired k value is cut off by the x interpolant:
if (std::abs(kx) > _maxk1 || std::abs(ky) > _maxk1) return std::complex<double>(0.,0.);
checkK();
double xKernelTransform = _xInterp.uval(kx*_uscale) * _xInterp.uval(ky*_uscale);
dbg<<"xKernelTransform = "<<xKernelTransform<<std::endl;
int No2 = _kimage->getBounds().getXMax();
int N = No2 * 2;
double kscale = No2/M_PI; // This is 1/dk
xdbg<<"kimage bounds = "<<_kimage->getBounds()<<", scale = "<<kscale<<std::endl;
kx *= kscale;
ky *= kscale;
int p1, p2, q1, q2; // Range over which we need to sum.
// Note, unlike for xValue, whre we limited the range to the size of the image,
// here, we allow i,j to nominally go off the kimage bounds, in which case we
// wrap around when using it, due to the periodic nature of the fft.
const double SMALL = 10.*std::numeric_limits<double>::epsilon();
if (std::abs(kx-std::floor(kx+0.01)) < SMALL*(std::abs(kx)+1)) {
// If kx or ky are integers, only sum 1 element in that direction.
p1 = p2 = int(std::floor(kx+0.01));
} else {
p1 = int(std::ceil(kx-_kInterp.xrange()));
p2 = int(std::floor(kx+_kInterp.xrange()));
}
if (std::abs(ky-std::floor(ky+0.01)) < SMALL*(std::abs(ky)+1)) {
q1 = q2 = int(std::floor(ky+0.01));
} else {
q1 = int(std::ceil(ky-_kInterp.xrange()));
q2 = int(std::floor(ky+_kInterp.xrange()));
}
dbg<<"p range = "<<p1<<"..."<<p2<<std::endl;
dbg<<"q range = "<<q1<<"..."<<q2<<std::endl;
// We'll need these for each row. Save them.
double xwt[p2-p1+1];
for (int p=p1, pp=0; p<=p2; ++p, ++pp) xwt[pp] = _kInterp.xval(p-kx);
std::complex<double> sum = 0.;
int pwrap1 = WrapKIndex(p1, No2, N);
int qwrap1 = WrapKIndex(q1, No2, N);
dbg<<"pwrap, qwrap = "<<qwrap1<<','<<qwrap1<<std::endl;
dbg<<"kimage bounds = "<<_kimage->getBounds()<<std::endl;
for (int q=q1, qwrap=qwrap1; q<=q2; ++q, ++qwrap) {
if (qwrap == No2) qwrap -= N;
std::complex<double> xsum = KValueInnerLoop(p2-p1+1,pwrap1,qwrap,No2,N,xwt,*_kimage);
sum += xsum * _kInterp.xval(q-ky);
}
dbg<<"sum = "<<sum<<std::endl;
sum *= xKernelTransform;
dbg<<"sum => "<<sum<<std::endl;
return sum;
}
void SBInterpolatedImage::SBInterpolatedImageImpl::checkK() const
{
// Conduct FFT
if (_kimage.get()) return;
int Nx = _image.getXMax()-_image.getXMin()+1;
dbg<<"Nx = "<<Nx<<std::endl;
Bounds<int> b(0,Nx/2,-Nx/2,Nx/2-1);
_kimage.reset(new ImageAlloc<std::complex<double> >(b));
rfft(_image, _kimage->view());
dbg<<"made kimage\n";
dbg<<"kimage bounds = "<<_kimage->getBounds()<<std::endl;
dbg<<"kimage flux = "<<(*_kimage)(0,0).real()<<std::endl;
}
template <typename T>
void SBInterpolatedImage::SBInterpolatedImageImpl::fillXImage(
ImageView<T> im,
double x0, double dx, int izero,
double y0, double dy, int jzero) const
{
dbg<<"SBInterpolatedImage fillXImage\n";
dbg<<"x = "<<x0<<" + i * "<<dx<<", izero = "<<izero<<std::endl;
dbg<<"y = "<<y0<<" + j * "<<dy<<", jzero = "<<jzero<<std::endl;
const int m = im.getNCol();
const int n = im.getNRow();
T* ptr = im.getData();
int skip = im.getNSkip();
assert(im.getStep() == 1);
const double SMALL = 10.*std::numeric_limits<double>::epsilon();
// Notation: We have two images to loop over here, so there are two sets of x,y values.
// To distinguish them, I'll use x,y for the output image,
// and p,q for the image we are interpolating.
// i,j are used for the indices in the output array as usual.
// Find the min/max x and y values where the output can be nonzero.
double minx = _nonzero_bounds.getXMin() - _xInterp.xrange();
double maxx = _nonzero_bounds.getXMax() + _xInterp.xrange();
double miny = _nonzero_bounds.getYMin() - _xInterp.xrange();
double maxy = _nonzero_bounds.getYMax() + _xInterp.xrange();
dbg<<"nonzero bounds = "<<_nonzero_bounds<<std::endl;
dbg<<"min/max x/y = "<<minx<<" "<<maxx<<" "<<miny<<" "<<maxy<<std::endl;
// Figure out what range of i,j these correspond to.
int i1 = int((minx-x0)/dx);
int i2 = int((maxx-x0)/dx);
int j1 = int((miny-y0)/dy);
int j2 = int((maxy-y0)/dy);
if (i2 < i1) std::swap(i1,i2);
if (j2 < j1) std::swap(j1,j2);
++i2; ++j2; // Make them one past the end, rather than last index to use.
if (i1 < 0) i1 = 0;
if (i2 > m) i2 = m;
if (j1 < 0) j1 = 0;
if (j2 > n) j2 = n;
xdbg<<"Output bounds = "<<im.getBounds()<<std::endl;
xdbg<<"Old i,j ranges = "<<0<<" "<<m<<" "<<0<<" "<<n<<std::endl;
xdbg<<"New i,j ranges = "<<i1<<" "<<i2<<" "<<j1<<" "<<j2<<std::endl;
// Fix up x0, y0, ptr, skip to correspond to these i,j ranges.
x0 += i1*dx;
y0 += j1*dy;
if (x0 < minx || x0 > maxx) { x0 += dx; ++i1; } // First points may be able to increase
if (y0 < miny || y0 > maxy) { y0 += dy; ++j1; } // by one more spot.
// If any of these is true, then object is fully off the target image.
if (i1 >= m || i2 < 0 || j1 >= n || j2 < 0 || i1 >= i2 || j1 >= j2) {
im.setZero();
return;
}
ptr += i1 + j1*im.getStride();
int mm = i2-i1; // We'll need this new row length a few times below.
skip += (m - mm);
// Each point in the output image is going to be
//
// I(x,y) = Sum_p,q wt(p-x) wt(q-y) _image(p,q)
//
// We have a lot of points with the same x or the same y, which means that
// xwt and ywt values are repeated a lot. So we want to find a way to reuse
// those _xInterp.xval values.
//
// Finally, it is most efficient if we have the innermost loop be along the
// direction with step=1 in both the input and output images.
//
// The way we try to meet all these goals is to build up the output image one row
// at a time:
//
// I(:,j) = Sum_q wt(q-y_j) rowq(:)
//
// where rowq is an array the same length as the output row.
// Then for each of these y,q values, we can compute rowq as
//
// rowq(i) = Sum_p wt(p-x_i) _image(p,q)
//
// The same set of wt(p-x) values are needed for every rowq calculation, so we
// compute them once and save them. Furthermore, the complete rowq calculation for
// a given q is independent of y, so we save that as well.
double x = x0;
double xwt[_xInterp.ixrange() * mm];
double p1ar[mm];
double p2ar[mm];
int k=0;
for (int i=i1; i<i2; ++i,x+=dx) {
int p1,p2;
// If x is (basically) an integer, only 1 p value.
// Otherwise, have a range based on xInterp.xrange()
if (std::abs(x-std::floor(x+0.01)) < SMALL*(std::abs(x)+1)) {
p1 = p2 = int(std::floor(x+0.01));
} else {
p1 = int(std::ceil(x-_xInterp.xrange()));
p2 = int(std::floor(x+_xInterp.xrange()));
}
// Limit to nonzero region
if (p1 < _nonzero_bounds.getXMin()) p1 = _nonzero_bounds.getXMin();
if (p2 > _nonzero_bounds.getXMax()) p2 = _nonzero_bounds.getXMax();
p1ar[i-i1] = p1;
p2ar[i-i1] = p2;
xdbg<<"i = "<<i<<" x = "<<x<<": p1,p2 = "<<p1<<','<<p2<<std::endl;
assert(p2-p1+1 <= _xInterp.ixrange());
for (int p=p1; p<=p2; ++p) {
xassert(k < _xInterp.ixrange()*mm);
xwt[k++] = _xInterp.xval(p-x);
}
}
// The inner calculation for rowq is the same for multiple y values since it is
// independent of y, so each time we use the same q, the rowq array is the same.
// Therefore we should cache these calculations to reuse when possible.
std::map<int, std::vector<double> > rowq_cache;
im.setZero();
double y = y0;
double temp[mm];
for (int j=j1; j<j2; ++j,y+=dy,ptr+=skip) {
memset(temp, 0, mm * sizeof(double)); // Zero out temp array
xdbg<<"j = "<<j<<", y = "<<y<<std::endl;
// If y is (basically) an integer, only 1 q value.
// Otherwise, have a range based on xInterp.xrange()
// Subtlety: also keep track of the minimum q we want to keep in the cache, which
// may be less than q1 to account for sometimes y being integer, sometimes not.
int q1,q2,qmin;
if (std::abs(y-std::floor(y+0.01)) < SMALL*(std::abs(y)+1)) {
q1 = q2 = int(std::floor(y+0.01));
qmin = int(std::ceil(y-_xInterp.xrange()));
} else {
qmin = q1 = int(std::ceil(y-_xInterp.xrange()));
q2 = int(std::floor(y+_xInterp.xrange()));
}
xdbg<<"q1,q2 = "<<q1<<','<<q2<<std::endl;
// Limit to nonzero region
if (q1 < _nonzero_bounds.getYMin()) q1 = _nonzero_bounds.getYMin();
if (q2 > _nonzero_bounds.getYMax()) q2 = _nonzero_bounds.getYMax();
xdbg<<"q1,q2 => "<<q1<<','<<q2<<std::endl;
// Dump any cached rows we don't need anymore.
while (rowq_cache.size() > 0 && rowq_cache.begin()->first < qmin) {
rowq_cache.erase(rowq_cache.begin());
}
for (int q=q1; q<=q2; ++q) {
// Get rowq from cache. If it isn't there, it will be an empty vector.
std::vector<double>& rowq = rowq_cache[q];
// If this rowq was not in cache, need to make it.
if (rowq.size() == 0) {
rowq.resize(mm);
double x = x0;
int k=0;
std::vector<double>::iterator row_it=rowq.begin();
for (int i=i1; i<i2; ++i,x+=dx,++row_it) {
*row_it = 0.;
int p1 = p1ar[i-i1];
int p2 = p2ar[i-i1];
const double* imptr = &_image(p1,q);
for (int p=p1; p<=p2; ++p) {
*row_it += xwt[k++] * *imptr++;
}
}
}
// Now add that to the output row with the ywt scaling.
double ywt = _xInterp.xval(q-y);
double* tptr = temp;
std::vector<double>::const_iterator row_it = rowq.begin();
for (int i=i1; i<i2; ++i) {
*tptr++ += *row_it++ * ywt;
}
}
// Now finally copy onto the real output image.
// Note: In addition to getting a tiny efficiency gain from having temp on
// the stack while doing the calculation above, this is also important for
// accuracy if the output image is T=float, so we don't gratuitously
// lose precision by adding floats rather than doubles.
double* tptr = temp;
for (int i=i1; i<i2; ++i) *ptr++ = *tptr++;
}
dbg<<"Done SBInterpolatedImage fillXImage\n";
}
template <typename T>
void SBInterpolatedImage::SBInterpolatedImageImpl::fillXImage(
ImageView<T> im,
double x0, double dx, double dxy,
double y0, double dy, double dyx) const
{
dbg<<"SBInterpolatedImage fillXImage\n";
dbg<<"x = "<<x0<<" + i * "<<dx<<" + j * "<<dxy<<std::endl;
dbg<<"y = "<<y0<<" + i * "<<dyx<<" + j * "<<dy<<std::endl;
const int m = im.getNCol();
const int n = im.getNRow();
T* ptr = im.getData();
int skip = im.getNSkip();
assert(im.getStep() == 1);
// In this version every _xInterp.xval call is different, so there's not really any
// way to save any of those calls. The only real optimization that still applies
// is the min/max checks to skip places where the output is zero.
// Find the min/max x and y values where the output can be nonzero.
double minx = _nonzero_bounds.getXMin() - _xInterp.xrange();
double maxx = _nonzero_bounds.getXMax() + _xInterp.xrange();
double miny = _nonzero_bounds.getYMin() - _xInterp.xrange();
double maxy = _nonzero_bounds.getYMax() + _xInterp.xrange();
dbg<<"nonzero bounds = "<<_nonzero_bounds<<std::endl;
dbg<<"min/max x/y = "<<minx<<" "<<maxx<<" "<<miny<<" "<<maxy<<std::endl;
// Figure out what range of i,j these correspond to.
// This is a bit more complicated than the separable case.
// We need to find the i,j corresponding to each corner of the allowed range.
// x = x0 + idx + jdxy
// y = y0 + idyx + jdy
// -> i = ( (x-x0) dy - (y-y0) dxy ) / (dx dy - dxy dyx)
// j = ( (y-y0) dx - (x-x0) dyx ) / (dx dy - dxy dyx)
double denom = dx*dy - dxy*dyx;
int ia = int(((minx-x0)*dy - (miny-y0)*dxy)/denom);
int ja = int(((miny-y0)*dx - (minx-x0)*dyx)/denom);
int ib = int(((minx-x0)*dy - (maxy-y0)*dxy)/denom);
int jb = int(((maxy-y0)*dx - (minx-x0)*dyx)/denom);
int ic = int(((maxx-x0)*dy - (miny-y0)*dxy)/denom);
int jc = int(((miny-y0)*dx - (maxx-x0)*dyx)/denom);
int id = int(((maxx-x0)*dy - (maxy-y0)*dxy)/denom);
int jd = int(((maxy-y0)*dx - (maxx-x0)*dyx)/denom);
dbg<<"Corners at "<<ia<<','<<ja<<" "<<ib<<','<<jb<<" "<<ic<<','<<jc<<" "<<id<<','<<jd<<std::endl;
int i1 = std::min( {ia,ib,ic,id} );
int i2 = std::max( {ia,ib,ic,id} );
int j1 = std::min( {ja,jb,jc,jd} );
int j2 = std::max( {ja,jb,jc,jd} );
dbg<<"i,j ranges = "<<i1<<" "<<i2<<" "<<j1<<" "<<j2<<std::endl;
++i2; ++j2; // Make them one past the end, rather than last index to use.
if (i1 < 0) i1 = 0;
if (i2 > m) i2 = m;
if (j1 < 0) j1 = 0;
if (j2 > n) j2 = n;
dbg<<"Output bounds = "<<im.getBounds()<<std::endl;
dbg<<"Old i,j ranges = "<<0<<" "<<m<<" "<<0<<" "<<n<<std::endl;
dbg<<"New i,j ranges = "<<i1<<" "<<i2<<" "<<j1<<" "<<j2<<std::endl;
if (i1 >= m || i2 < 0 || j1 >= n || j2 < 0 || i1 >= i2 || j1 >= j2) {
im.setZero();
return;
}
// Fix up x0, y0, ptr, skip to correspond to these i,j ranges.
x0 += i1*dx + j1*dxy;
y0 += j1*dy + i1*dyx;
ptr += i1 + j1*im.getStride();
int mm = i2-i1;
skip += (m - mm);
im.setZero();
for (int j=j1; j<j2; ++j,x0+=dxy,y0+=dy,ptr+=skip) {
double x = x0;
double y = y0;
for (int i=i1; i<i2; ++i,x+=dx,y+=dyx,++ptr) {
// Still want this check even with above i1,i2,j1,j2 stuff, since projected
// region is a parallelogram, so some points can still be sipped.
if (y > maxy || y < miny || x > maxx || x < minx) continue;
int p1 = int(std::ceil(x-_xInterp.xrange()));
int p2 = int(std::floor(x+_xInterp.xrange()));
int q1 = int(std::ceil(y-_xInterp.xrange()));
int q2 = int(std::floor(y+_xInterp.xrange()));
if (p1 < _nonzero_bounds.getXMin()) p1 = _nonzero_bounds.getXMin();
if (p2 > _nonzero_bounds.getXMax()) p2 = _nonzero_bounds.getXMax();
if (q1 < _nonzero_bounds.getYMin()) q1 = _nonzero_bounds.getYMin();
if (q2 > _nonzero_bounds.getYMax()) q2 = _nonzero_bounds.getYMax();
double xwt[p2-p1+1];
for (int p=p1, pp=0; p<=p2; ++p, ++pp) {
xwt[pp] = _xInterp.xval(p-x);
}
double sum=0.;
for (int q=q1; q<=q2; ++q) {
double ywt = _xInterp.xval(q-y);
double xsum = 0.;
const double* imptr = &_image(p1,q);
for (int p=p1, pp=0; p<=p2; ++p, ++pp) {
xsum += xwt[pp] * *imptr++;
}
sum += xsum * ywt;
}
xassert(ptr >= im.getData());
xassert(ptr < im.getData() + im.getNElements());
*ptr = sum;
}
}
}
template <typename T>
void SBInterpolatedImage::SBInterpolatedImageImpl::fillKImage(
ImageView<std::complex<T> > im,
double kx0, double dkx, int izero,
double ky0, double dky, int jzero) const
{
dbg<<"SBInterpolatedImage fillKImage\n";
dbg<<"kx = "<<kx0<<" + i * "<<dkx<<", izero = "<<izero<<std::endl;
dbg<<"ky = "<<ky0<<" + j * "<<dky<<", jzero = "<<jzero<<std::endl;
const int m = im.getNCol();
const int n = im.getNRow();
std::complex<T>* ptr = im.getData();
assert(im.getStep() == 1);
int skip = im.getNSkip();
checkK();
const double SMALL = 10.*std::numeric_limits<double>::epsilon();
// Only non-zero where |u| <= maxu
double absdkx = std::abs(dkx);
double absdky = std::abs(dky);
int i1 = std::max( int(-_maxk1/absdkx-kx0/dkx) , 0 );
int i2 = std::min( int(_maxk1/absdkx-kx0/dkx)+1 , m );
int j1 = std::max( int(-_maxk1/absdky-ky0/dky) , 0 );
int j2 = std::min( int(_maxk1/absdky-ky0/dky)+1 , n );
dbg<<"i1,i2,j1,j1 = "<<i1<<','<<i2<<','<<j1<<','<<j2<<std::endl;
if (i1 >= m || i2 < 0 || j1 >= n || j2 < 0 || i1 >= i2 || j1 >= j2) {
im.setZero();
return;
}
kx0 += i1*dkx;
ky0 += j1*dky;
ptr += i1 + j1*im.getStride();
int mm = i2-i1;
skip += (m - mm);
// For the rest of the range, calculate ux, uy values
std::vector<double> ux(i2-i1);
typedef std::vector<double>::iterator It;
It uxit = ux.begin();
double kx = kx0;
for (int i=i1; i<i2; ++i,kx+=dkx) *uxit++ = kx * _uscale;
std::vector<double> uy(j2-j1);
It uyit = uy.begin();
double ky = ky0;
for (int j=j1; j<j2; ++j,ky+=dky) *uyit++ = ky * _uscale;
// Rescale the k values by kscale
int No2 = _kimage->getBounds().getXMax();
int N = No2 * 2;
double kscale = No2/M_PI;
dbg<<"kimage bounds = "<<_kimage->getBounds()<<", scale = "<<kscale<<std::endl;
kx0 *= kscale;
ky0 *= kscale;
dkx *= kscale;
dky *= kscale;
// The caching stuff is the same here as it was for fillXValue. The only difference
// is that we need to wrap around the p,q values and handle the conjugation possibility
// correctly. (cf. comments in kValue method.)
kx = kx0;
double xwt[_kInterp.ixrange() * mm];
double p1ar[mm];
double p2ar[mm];
int k=0;
for (int i=i1; i<i2; ++i,kx+=dkx) {
int p1, p2; // Range over which we need to sum.
if (std::abs(kx-std::floor(kx+0.01)) < SMALL*(std::abs(kx)+1)) {
// If kx is integer, only sum 1 element in that direction.
p1 = p2 = int(std::floor(kx+0.01));
} else {
p1 = int(std::ceil(kx-_kInterp.xrange()));
p2 = int(std::floor(kx+_kInterp.xrange()));
}
p1ar[i-i1] = p1;
p2ar[i-i1] = p2;
xdbg<<"i = "<<i<<" kx = "<<kx<<": p1,p2 = "<<p1<<','<<p2<<std::endl;
assert(p2-p1+1 <= _kInterp.ixrange());
for (int p=p1; p<=p2; ++p) {
xassert(k < _kInterp.ixrange()*mm);
xwt[k++] = _kInterp.xval(p-kx);
}
}
// Again, can cache the rowq vectors.
std::map<int, std::vector<std::complex<double> > > rowq_cache;
// Pre-calculate xInterp factors in place
uxit = ux.begin();
for (int i=i1; i<i2; ++i,++uxit) *uxit = _xInterp.uval(*uxit);
uyit = uy.begin();
for (int j=j1; j<j2; ++j,++uyit) *uyit = _xInterp.uval(*uyit);
im.setZero();
ky = ky0;
uyit = uy.begin();
double temp[2*mm]; // Can't put complex<double> array on stack, so reinterpret_cast below.
for (int j=j1; j<j2; ++j,ky+=dky,ptr+=skip,++uyit) {
memset(temp, 0, 2*mm * sizeof(double));
xdbg<<"j = "<<j<<", ky = "<<ky<<std::endl;
// If y is (basically) an integer, only 1 q value.
int q1,q2,qmin;
if (std::abs(ky-std::floor(ky+0.01)) < SMALL*(std::abs(ky)+1)) {
q1 = q2 = int(std::floor(ky+0.01));
qmin = int(std::ceil(ky-_kInterp.xrange()));
} else {
qmin = q1 = int(std::ceil(ky-_kInterp.xrange()));
q2 = int(std::floor(ky+_kInterp.xrange()));
}
xdbg<<"q1,q2 = "<<q1<<','<<q2<<std::endl;
// Dump any cached rows we don't need anymore.
while (rowq_cache.size() > 0 && rowq_cache.begin()->first < qmin) {
rowq_cache.erase(rowq_cache.begin());
}
int qwrap1 = WrapKIndex(q1, No2, N);
for (int q=q1, qwrap=qwrap1; q<=q2; ++q, ++qwrap) {
if (qwrap == No2) qwrap -= N;
// Get rowq from cache. If it isn't there, it will be an empty vector.
std::vector<std::complex<double> >& rowq = rowq_cache[q];
// If this rowq was not in cache, need to make it.
if (rowq.size() == 0) {
rowq.resize(mm);
kx = kx0;
int k=0;
std::vector<std::complex<double> >::iterator row_it=rowq.begin();
for (int i=i1; i<i2; ++i,kx+=dkx) {
int p1 = p1ar[i-i1];
int p2 = p2ar[i-i1];
int pwrap1 = WrapKIndex(p1, No2, N);
*row_it++ = KValueInnerLoop(p2-p1+1,pwrap1,qwrap,No2,N,&xwt[k],*_kimage);
k += p2-p1+1;
}
}
// Now add that to the output row with the ywt scaling.
double ywt = _kInterp.xval(q-ky);
std::complex<double>* tptr = reinterpret_cast<std::complex<double>*>(temp);
std::vector<std::complex<double> >::const_iterator row_it = rowq.begin();
for (int i=i1; i<i2; ++i) {
xassert(row_it < rowq.end());
xassert((void*)tptr < (void*)(temp + 2*mm));
*tptr++ += *row_it++ * ywt;
}
}
// Now account for the x-interpolant
// And finally copy onto the real output image.
// Note: In addition to getting a tiny efficiency gain from having temp on
// the stack while doing the calculation above, this is also important for
// accuracy if the output image is complex<float>, so we don't gratuitously
// lose precision by adding floats rather than doubles.
uxit = ux.begin();
std::complex<double>* tptr = reinterpret_cast<std::complex<double>*>(temp);
for (int i=i1; i<i2; ++i) {
xassert(ptr < im.getData() + im.getNElements());
xassert(uxit < ux.end());
xassert(uyit < uy.end());
xassert((void*)tptr < (void*)(temp + 2*mm));
*ptr++ = *uxit++ * *uyit * *tptr++;
}
}
dbg<<"Done SBInterpolatedImage fillKImage\n";
}
template <typename T>
void SBInterpolatedImage::SBInterpolatedImageImpl::fillKImage(
ImageView<std::complex<T> > im,
double kx0, double dkx, double dkxy,
double ky0, double dky, double dkyx) const
{
dbg<<"SBInterpolatedImage fillKImage\n";
dbg<<"kx = "<<kx0<<" + i * "<<dkx<<" + j * "<<dkxy<<std::endl;
dbg<<"ky = "<<ky0<<" + i * "<<dkyx<<" + j * "<<dky<<std::endl;
const int m = im.getNCol();
const int n = im.getNRow();
std::complex<T>* ptr = im.getData();
int skip = im.getNSkip();
assert(im.getStep() == 1);
checkK();
double ux0 = kx0 * _uscale;
double uy0 = ky0 * _uscale;
double dux = dkx * _uscale;
double duy = dky * _uscale;
double duxy = dkxy * _uscale;
double duyx = dkyx * _uscale;
int No2 = _kimage->getBounds().getXMax();
int N = No2 * 2;
double kscale = No2/M_PI;
kx0 *= kscale;
dkx *= kscale;
dkxy *= kscale;
ky0 *= kscale;
dky *= kscale;
dkyx *= kscale;
double maxk1 = _maxk1 * kscale;
for (int j=0; j<n; ++j,kx0+=dkxy,ky0+=dky,ux0+=duxy,uy0+=duy,ptr+=skip) {
double kx = kx0;
double ky = ky0;
double ux = ux0;
double uy = uy0;
for (int i=0; i<m; ++i,kx+=dkx,ky+=dkyx,ux+=dux,uy+=duyx) {
if (std::abs(kx) > maxk1 || std::abs(ky) > maxk1) {
*ptr++ = T(0);
} else {
int p1 = int(std::ceil(kx-_kInterp.xrange()));
int p2 = int(std::floor(kx+_kInterp.xrange()));
int q1 = int(std::ceil(ky-_kInterp.xrange()));
int q2 = int(std::floor(ky+_kInterp.xrange()));
double xwt[p2-p1+1];
for (int p=p1, pp=0; p<=p2; ++p, ++pp) xwt[pp] = _kInterp.xval(p-kx);
std::complex<double> sum = 0.;
int pwrap1 = WrapKIndex(p1, No2, N);
int qwrap1 = WrapKIndex(q1, No2, N);
for (int q=q1, qwrap=qwrap1; q<=q2; ++q, ++qwrap) {
if (qwrap == No2) qwrap -= N;
double ywt = _kInterp.xval(q-ky);
sum += ywt * KValueInnerLoop(p2-p1+1,pwrap1,qwrap,No2,N,xwt,*_kimage);
}
*ptr++ = _xInterp.uval(ux) * _xInterp.uval(uy) * sum;
}
}
}
}