gridder_cxx.h 41.7 KB
Newer Older
Martin Reinecke's avatar
Martin Reinecke committed
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
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
#ifndef GRIDDER_CXX_H
#define GRIDDER_CXX_H

/*
 *  This file is part of nifty_gridder.
 *
 *  nifty_gridder is free software; you can redistribute it and/or modify
 *  it under the terms of the GNU General Public License as published by
 *  the Free Software Foundation; either version 2 of the License, or
 *  (at your option) any later version.
 *
 *  nifty_gridder is distributed in the hope that it will be useful,
 *  but WITHOUT ANY WARRANTY; without even the implied warranty of
 *  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 *  GNU General Public License for more details.
 *
 *  You should have received a copy of the GNU General Public License
 *  along with nifty_gridder; if not, write to the Free Software
 *  Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
 */

/* Copyright (C) 2019 Max-Planck-Society
   Author: Martin Reinecke */

#include <iostream>
#include <algorithm>
#include <cstdlib>
#include <cmath>
#include <vector>
#include <array>

#include "mr_util/error_handling.h"
#include "mr_util/fft.h"
#include "mr_util/threading.h"
#include "mr_util/useful_macros.h"
#include "mr_util/mav.h"
#include "mr_util/gl_integrator.h"

#if defined(__GNUC__)
#define ALIGNED(align) __attribute__ ((aligned(align)))
#else
#define ALIGNED(align)
#endif

namespace gridder {

namespace detail {

using namespace std;
using namespace mr;

template<size_t ndim> void checkShape
  (const array<size_t, ndim> &shp1, const array<size_t, ndim> &shp2)
  {
  for (size_t i=0; i<ndim; ++i)
    MR_assert(shp1[i]==shp2[i], "shape mismatch");
  }

template<typename T> inline T fmod1 (T v)
  { return v-floor(v); }

template<typename T, size_t ndim> class tmpStorage
  {
  private:
    vector<T> d;
    mav<T,ndim> mav_;

    static size_t prod(const array<size_t,ndim> &shp)
      {
      size_t res=1;
      for (auto v: shp) res*=v;
      return res;
      }

  public:
    tmpStorage(const array<size_t,ndim> &shp)
      : d(prod(shp)), mav_(d.data(), shp) {}
    mav<T,ndim> &getMav() { return mav_; }
    const_mav<T,ndim> getCmav() { return cmav(mav_); }
    void fill(const T & val)
      { std::fill(d.begin(), d.end(), val); }
  };

//
// Start of real gridder functionality
//

template<typename T> void complex2hartley
  (const const_mav<complex<T>, 2> &grid, const mav<T,2> &grid2, size_t nthreads)
  {
  checkShape(grid.shape(), grid2.shape());
  size_t nu=grid.shape(0), nv=grid.shape(1);

  execStatic(nu, nthreads, 0, [&](Scheduler &sched)
    {
    while (auto rng=sched.getNext()) for(auto u=rng.lo; u<rng.hi; ++u)
      {
      size_t xu = (u==0) ? 0 : nu-u;
      for (size_t v=0; v<nv; ++v)
        {
        size_t xv = (v==0) ? 0 : nv-v;
        grid2(u,v) += T(0.5)*(grid( u, v).real()+grid( u, v).imag()+
                              grid(xu,xv).real()-grid(xu,xv).imag());
        }
      }
    });
  }

template<typename T> void hartley2complex
  (const const_mav<T,2> &grid, const mav<complex<T>,2> &grid2, size_t nthreads)
  {
  checkShape(grid.shape(), grid2.shape());
  size_t nu=grid.shape(0), nv=grid.shape(1);

  execStatic(nu, nthreads, 0, [&](Scheduler &sched)
    {
    while (auto rng=sched.getNext()) for(auto u=rng.lo; u<rng.hi; ++u)
      {
      size_t xu = (u==0) ? 0 : nu-u;
      for (size_t v=0; v<nv; ++v)
        {
        size_t xv = (v==0) ? 0 : nv-v;
        T v1 = T(0.5)*grid( u, v);
        T v2 = T(0.5)*grid(xu,xv);
        grid2(u,v) = std::complex<T>(v1+v2, v1-v2);
        }
      }
    });
  }

template<typename T> void hartley2_2D(const const_mav<T,2> &in,
  const mav<T,2> &out, size_t nthreads)
  {
  checkShape(in.shape(), out.shape());
  size_t nu=in.shape(0), nv=in.shape(1);
  ptrdiff_t sz=ptrdiff_t(sizeof(T));
  stride_t stri{sz*in.stride(0), sz*in.stride(1)};
  stride_t stro{sz*out.stride(0), sz*out.stride(1)};
  auto d_i = in.data();
  auto ptmp = out.data();
  r2r_separable_hartley({nu, nv}, stri, stro, {0,1}, d_i, ptmp, T(1),
    nthreads);
  execStatic((nu+1)/2-1, nthreads, 0, [&](Scheduler &sched)
    {
    while (auto rng=sched.getNext()) for(auto i=rng.lo+1; i<rng.hi+1; ++i)
      for(size_t j=1; j<(nv+1)/2; ++j)
         {
         T a = ptmp[i*nv+j];
         T b = ptmp[(nu-i)*nv+j];
         T c = ptmp[i*nv+nv-j];
         T d = ptmp[(nu-i)*nv+nv-j];
         ptmp[i*nv+j] = T(0.5)*(a+b+c-d);
         ptmp[(nu-i)*nv+j] = T(0.5)*(a+b+d-c);
         ptmp[i*nv+nv-j] = T(0.5)*(a+c+d-b);
         ptmp[(nu-i)*nv+nv-j] = T(0.5)*(b+c+d-a);
         }
     });
  }

class ES_Kernel
  {
  private:
    static constexpr double pi = 3.141592653589793238462643383279502884197;
    double beta;
    int p;
    vector<double> x, wgt, psi;
    size_t supp;

  public:
    ES_Kernel(size_t supp_, double ofactor, size_t nthreads)
      : beta(get_beta(supp_,ofactor)*supp_), p(int(1.5*supp_+2)), supp(supp_)
      {
      GL_Integrator integ(2*p,nthreads);
      x = integ.coordsSymmetric();
      wgt = integ.weightsSymmetric();
      psi=x;
      for (auto &v:psi)
        v=operator()(v);
      }
    ES_Kernel(size_t supp_, size_t nthreads)
      : ES_Kernel(supp_, 2., nthreads){}

    double operator()(double v) const { return exp(beta*(sqrt(1.-v*v)-1.)); }
    /* Compute correction factors for the ES gridding kernel
       This implementation follows eqs. (3.8) to (3.10) of Barnett et al. 2018 */
    double corfac(double v) const
      {
      double tmp=0;
      for (int i=0; i<p; ++i)
        tmp += wgt[i]*psi[i]*cos(pi*supp*v*x[i]);
      return 2./(supp*tmp);
      }
    static double get_beta(size_t supp, double ofactor=2)
      {
      MR_assert((supp>=2) && (supp<=15), "unsupported support size");
      if (ofactor>=2)
        {
        static const vector<double> opt_beta {-1, 0.14, 1.70, 2.08, 2.205, 2.26,
          2.29, 2.307, 2.316, 2.3265, 2.3324, 2.282, 2.294, 2.304, 2.3138, 2.317};
        MR_assert(supp<opt_beta.size(), "bad support size");
        return opt_beta[supp];
        }
      if (ofactor>=1.2)
        {
        // empirical, but pretty accurate approximation
        static const array<double,16> betacorr{0,0,-0.51,-0.21,-0.1,-0.05,-0.025,-0.0125,0,0,0,0,0,0,0,0};
        auto x0 = 1./(2*ofactor);
        auto bcstrength=1.+(x0-0.25)*2.5;
        return 2.32+bcstrength*betacorr[supp]+(0.25-x0)*3.1;
        }
      MR_fail("oversampling factor is too small");
      }

    static size_t get_supp(double epsilon, double ofactor=2)
      {
      double epssq = epsilon*epsilon;
      if (ofactor>=2)
        {
        static const vector<double> maxmaperr { 1e8, 0.19, 2.98e-3, 5.98e-5,
          1.11e-6, 2.01e-8, 3.55e-10, 5.31e-12, 8.81e-14, 1.34e-15, 2.17e-17,
          2.12e-19, 2.88e-21, 3.92e-23, 8.21e-25, 7.13e-27 };

        for (size_t i=2; i<maxmaperr.size(); ++i)
          if (epssq>maxmaperr[i]) return i;
        MR_fail("requested epsilon too small - minimum is 1e-13");
        }
      if (ofactor>=1.2)
        {
        for (size_t w=2; w<16; ++w)
          {
          auto estimate = 12*exp(-2.*w*ofactor); // empirical, not very good approximation
          if (epssq>estimate) return w;
          }
        MR_fail("requested epsilon too small");
        }
      MR_fail("oversampling factor is too small");
      }
  };

/* Compute correction factors for the ES gridding kernel
   This implementation follows eqs. (3.8) to (3.10) of Barnett et al. 2018 */
vector<double> correction_factors(size_t n, double ofactor, size_t nval, size_t supp,
  size_t nthreads)
  {
  ES_Kernel kernel(supp, ofactor, nthreads);
  vector<double> res(nval);
  double xn = 1./n;
  execStatic(nval, nthreads, 0, [&](Scheduler &sched)
    {
    while (auto rng=sched.getNext()) for(auto i=rng.lo; i<rng.hi; ++i)
      res[i] = kernel.corfac(i*xn);
    });
  return res;
  }

using idx_t = uint32_t;

struct RowChan
  {
  idx_t row, chan;
  };

struct UVW
  {
  double u, v, w;
  UVW() {}
  UVW(double u_, double v_, double w_) : u(u_), v(v_), w(w_) {}
  UVW operator* (double fct) const
    { return UVW(u*fct, v*fct, w*fct); }
  void Flip() { u=-u; v=-v; w=-w; }
  bool FixW()
    {
    bool flip = w<0;
    if (flip) Flip();
    return flip;
    }
  };

class Baselines
  {
  protected:
    vector<UVW> coord;
    vector<double> f_over_c;
    idx_t nrows, nchan;
    idx_t shift, mask;

  public:
    template<typename T> Baselines(const const_mav<T,2> &coord_,
      const const_mav<T,1> &freq, bool negate_v=false)
      {
      constexpr double speedOfLight = 299792458.;
      MR_assert(coord_.shape(1)==3, "dimension mismatch");
      auto hugeval = size_t(~(idx_t(0)));
      MR_assert(coord_.shape(0)<hugeval, "too many entries in MS");
      MR_assert(coord_.shape(1)<hugeval, "too many entries in MS");
      MR_assert(coord_.size()<hugeval, "too many entries in MS");
      nrows = coord_.shape(0);
      nchan = freq.shape(0);
      shift=0;
      while((idx_t(1)<<shift)<nchan) ++shift;
      mask=(idx_t(1)<<shift)-1;
      MR_assert(nrows*(mask+1)<hugeval, "too many entries in MS");
      f_over_c.resize(nchan);
      for (size_t i=0; i<nchan; ++i)
        {
        MR_assert(freq[i]>0, "negative channel frequency encountered");
        f_over_c[i] = freq(i)/speedOfLight;
        }
      coord.resize(nrows);
      if (negate_v)
        for (size_t i=0; i<coord.size(); ++i)
          coord[i] = UVW(coord_(i,0), -coord_(i,1), coord_(i,2));
      else
        for (size_t i=0; i<coord.size(); ++i)
          coord[i] = UVW(coord_(i,0), coord_(i,1), coord_(i,2));
      }

    RowChan getRowChan(idx_t index) const
      { return RowChan{index>>shift, index&mask}; }

    UVW effectiveCoord(const RowChan &rc) const
      { return coord[rc.row]*f_over_c[rc.chan]; }
    UVW effectiveCoord(idx_t index) const
      { return effectiveCoord(getRowChan(index)); }
    size_t Nrows() const { return nrows; }
    size_t Nchannels() const { return nchan; }
    idx_t getIdx(idx_t irow, idx_t ichan) const
      { return ichan+(irow<<shift); }
  };

class GridderConfig
  {
  protected:
    size_t nx_dirty, ny_dirty, nu, nv;
    double ofactor, eps, psx, psy;
    size_t supp, nsafe;
    double beta;
    size_t nthreads;
    double ushift, vshift;
    int maxiu0, maxiv0;

    template<typename T> complex<T> wscreen(T x, T y, T w, bool adjoint) const
      {
      constexpr T pi = T(3.141592653589793238462643383279502884197);
      T tmp = 1-x-y;
      if (tmp<=0) return 1; // no phase factor beyond the horizon
      T nm1 = (-x-y)/(sqrt(tmp)+1); // more accurate form of sqrt(1-x-y)-1
      T phase = 2*pi*w*nm1;
      if (adjoint) phase *= -1;
      return complex<T>(cos(phase), sin(phase));
      }

  public:
    GridderConfig(size_t nxdirty, size_t nydirty, size_t nu_, size_t nv_,
      double epsilon, double pixsize_x, double pixsize_y, size_t nthreads_)
      : nx_dirty(nxdirty), ny_dirty(nydirty), nu(nu_), nv(nv_),
        ofactor(min(double(nu)/nxdirty, double(nv)/nydirty)),
        eps(epsilon),
        psx(pixsize_x), psy(pixsize_y),
        supp(ES_Kernel::get_supp(epsilon, ofactor)), nsafe((supp+1)/2),
        beta(ES_Kernel::get_beta(supp, ofactor)*supp),
        nthreads(nthreads_),
        ushift(supp*(-0.5)+1+nu), vshift(supp*(-0.5)+1+nv),
        maxiu0((nu+nsafe)-supp), maxiv0((nv+nsafe)-supp)
      {
      MR_assert(nu>=2*nsafe, "nu too small");
      MR_assert(nv>=2*nsafe, "nv too small");
      MR_assert((nx_dirty&1)==0, "nx_dirty must be even");
      MR_assert((ny_dirty&1)==0, "ny_dirty must be even");
      MR_assert((nu&1)==0, "nu must be even");
      MR_assert((nv&1)==0, "nv must be even");
      MR_assert(epsilon>0, "epsilon must be positive");
      MR_assert(pixsize_x>0, "pixsize_x must be positive");
      MR_assert(pixsize_y>0, "pixsize_y must be positive");
      MR_assert(ofactor>=1.2, "oversampling factor smaller than 1.2");
      }
    GridderConfig(size_t nxdirty, size_t nydirty,
      double epsilon, double pixsize_x, double pixsize_y, size_t nthreads_)
      : GridderConfig(nxdirty, nydirty, max<size_t>(30,2*nxdirty),
                      max<size_t>(30,2*nydirty), epsilon, pixsize_x,
                      pixsize_y, nthreads_) {}
    size_t Nxdirty() const { return nx_dirty; }
    size_t Nydirty() const { return ny_dirty; }
    double Epsilon() const { return eps; }
    double Pixsize_x() const { return psx; }
    double Pixsize_y() const { return psy; }
    size_t Nu() const { return nu; }
    size_t Nv() const { return nv; }
    size_t Supp() const { return supp; }
    size_t Nsafe() const { return nsafe; }
    double Beta() const { return beta; }
    size_t Nthreads() const { return nthreads; }
    double Ofactor() const{ return ofactor; }

    template<typename T> void grid2dirty_post(const mav<T,2> &tmav,
      const mav<T,2> &dirty) const
      {
      checkShape(dirty.shape(), {nx_dirty,ny_dirty});
      auto cfu = correction_factors(nu, ofactor, nx_dirty/2+1, supp, nthreads);
      auto cfv = correction_factors(nv, ofactor, ny_dirty/2+1, supp, nthreads);
      execStatic(nx_dirty, nthreads, 0, [&](Scheduler &sched)
        {
        while (auto rng=sched.getNext()) for(auto i=rng.lo; i<rng.hi; ++i)
          {
          int icfu = abs(int(nx_dirty/2)-int(i));
          for (size_t j=0; j<ny_dirty; ++j)
            {
            int icfv = abs(int(ny_dirty/2)-int(j));
            size_t i2 = nu-nx_dirty/2+i;
            if (i2>=nu) i2-=nu;
            size_t j2 = nv-ny_dirty/2+j;
            if (j2>=nv) j2-=nv;
            // FIXME: for some reason g++ warns about double-to-float conversion
            // here, even though there is an explicit cast...
            dirty(i,j) = tmav(i2,j2)*T(cfu[icfu]*cfv[icfv]);
            }
          }
        });
      }
    template<typename T> void grid2dirty_post2(
      const mav<complex<T>,2> &tmav, const mav<T,2> &dirty, T w) const
      {
      checkShape(dirty.shape(), {nx_dirty,ny_dirty});
      double x0 = -0.5*nx_dirty*psx,
             y0 = -0.5*ny_dirty*psy;
      execStatic(nx_dirty/2+1, nthreads, 0, [&](Scheduler &sched)
        {
        while (auto rng=sched.getNext()) for(auto i=rng.lo; i<rng.hi; ++i)
          {
          T fx = T(x0+i*psx);
          fx *= fx;
          for (size_t j=0; j<=ny_dirty/2; ++j)
            {
            T fy = T(y0+j*psy);
            auto ws = wscreen(fx, fy*fy, w, true);
            size_t ix = nu-nx_dirty/2+i;
            if (ix>=nu) ix-=nu;
            size_t jx = nv-ny_dirty/2+j;
            if (jx>=nv) jx-=nv;
            dirty(i,j) += (tmav(ix,jx)*ws).real(); // lower left
            size_t i2 = nx_dirty-i, j2 = ny_dirty-j;
            size_t ix2 = nu-nx_dirty/2+i2;
            if (ix2>=nu) ix2-=nu;
            size_t jx2 = nv-ny_dirty/2+j2;
            if (jx2>=nv) jx2-=nv;
            if ((i>0)&&(i<i2))
              {
              dirty(i2,j) += (tmav(ix2,jx)*ws).real(); // lower right
              if ((j>0)&&(j<j2))
                dirty(i2,j2) += (tmav(ix2,jx2)*ws).real(); // upper right
              }
            if ((j>0)&&(j<j2))
              dirty(i,j2) += (tmav(ix,jx2)*ws).real(); // upper left
            }
          }
        });
      }

    template<typename T> void grid2dirty(const const_mav<T,2> &grid,
      const mav<T,2> &dirty) const
      {
      checkShape(grid.shape(), {nu,nv});
      tmpStorage<T,2> tmpdat({nu,nv});
      auto tmav = tmpdat.getMav();
      hartley2_2D<T>(grid, tmav, nthreads);
      grid2dirty_post(tmav, dirty);
      }

    template<typename T> void grid2dirty_c_overwrite_wscreen_add
      (const mav<complex<T>,2> &grid, const mav<T,2> &dirty, T w) const
      {
      checkShape(grid.shape(), {nu,nv});
      constexpr auto sc = ptrdiff_t(sizeof(complex<T>));
      c2c({nu,nv},{grid.stride(0)*sc,grid.stride(1)*sc},
        {grid.stride(0)*sc, grid.stride(1)*sc}, {0,1}, BACKWARD,
        grid.data(), grid.data(), T(1), nthreads);
      grid2dirty_post2(grid, dirty, w);
      }

    template<typename T> void dirty2grid_pre(const const_mav<T,2> &dirty,
      const mav<T,2> &grid) const
      {
      checkShape(dirty.shape(), {nx_dirty, ny_dirty});
      checkShape(grid.shape(), {nu, nv});
      auto cfu = correction_factors(nu, ofactor, nx_dirty/2+1, supp, nthreads);
      auto cfv = correction_factors(nv, ofactor, ny_dirty/2+1, supp, nthreads);
      grid.fill(0);
      execStatic(nx_dirty, nthreads, 0, [&](Scheduler &sched)
        {
        while (auto rng=sched.getNext()) for(auto i=rng.lo; i<rng.hi; ++i)
          {
          int icfu = abs(int(nx_dirty/2)-int(i));
          for (size_t j=0; j<ny_dirty; ++j)
            {
            int icfv = abs(int(ny_dirty/2)-int(j));
            size_t i2 = nu-nx_dirty/2+i;
            if (i2>=nu) i2-=nu;
            size_t j2 = nv-ny_dirty/2+j;
            if (j2>=nv) j2-=nv;
            grid(i2,j2) = dirty(i,j)*T(cfu[icfu]*cfv[icfv]);
            }
          }
        });
      }
    template<typename T> void dirty2grid_pre2(const const_mav<T,2> &dirty,
      const mav<complex<T>,2> &grid, T w) const
      {
      checkShape(dirty.shape(), {nx_dirty, ny_dirty});
      checkShape(grid.shape(), {nu, nv});
      grid.fill(0);

      double x0 = -0.5*nx_dirty*psx,
             y0 = -0.5*ny_dirty*psy;
      execStatic(nx_dirty/2+1, nthreads, 0, [&](Scheduler &sched)
        {
        while (auto rng=sched.getNext()) for(auto i=rng.lo; i<rng.hi; ++i)
          {
          T fx = T(x0+i*psx);
          fx *= fx;
          for (size_t j=0; j<=ny_dirty/2; ++j)
            {
            T fy = T(y0+j*psy);
            auto ws = wscreen(fx, fy*fy, w, false);
            size_t ix = nu-nx_dirty/2+i;
            if (ix>=nu) ix-=nu;
            size_t jx = nv-ny_dirty/2+j;
            if (jx>=nv) jx-=nv;
            grid(ix,jx) = dirty(i,j)*ws; // lower left
            size_t i2 = nx_dirty-i, j2 = ny_dirty-j;
            size_t ix2 = nu-nx_dirty/2+i2;
            if (ix2>=nu) ix2-=nu;
            size_t jx2 = nv-ny_dirty/2+j2;
            if (jx2>=nv) jx2-=nv;
            if ((i>0)&&(i<i2))
              {
              grid(ix2,jx) = dirty(i2,j)*ws; // lower right
              if ((j>0)&&(j<j2))
                grid(ix2,jx2) = dirty(i2,j2)*ws; // upper right
              }
            if ((j>0)&&(j<j2))
              grid(ix,jx2) = dirty(i,j2)*ws; // upper left
            }
          }
        });
      }

    template<typename T> void dirty2grid(const const_mav<T,2> &dirty,
      const mav<T,2> &grid) const
      {
      dirty2grid_pre(dirty, grid);
      hartley2_2D<T>(cmav(grid), grid, nthreads);
      }

    template<typename T> void dirty2grid_c_wscreen(const const_mav<T,2> &dirty,
      const mav<complex<T>,2> &grid, T w) const
      {
      dirty2grid_pre2(dirty, grid, w);
      constexpr auto sc = ptrdiff_t(sizeof(complex<T>));
      stride_t strides{grid.stride(0)*sc,grid.stride(1)*sc};
      c2c({nu,nv}, strides, strides, {0,1}, FORWARD,
        grid.data(), grid.data(), T(1), nthreads);
      }

    void getpix(double u_in, double v_in, double &u, double &v, int &iu0, int &iv0) const
      {
      u=fmod1(u_in*psx)*nu;
      iu0 = min(int(u+ushift)-int(nu), maxiu0);
      v=fmod1(v_in*psy)*nv;
      iv0 = min(int(v+vshift)-int(nv), maxiv0);
      }

    template<typename T> void apply_wscreen(const const_mav<complex<T>,2> &dirty,
      mav<complex<T>,2> &dirty2, double w, bool adjoint) const
      {
      checkShape(dirty.shape(), {nx_dirty, ny_dirty});
      checkShape(dirty2.shape(), {nx_dirty, ny_dirty});

      double x0 = -0.5*nx_dirty*psx,
             y0 = -0.5*ny_dirty*psy;
      execStatic(nx_dirty/2+1, nthreads, 0, [&](Scheduler &sched)
        {
        while (auto rng=sched.getNext()) for(auto i=rng.lo; i<rng.hi; ++i)
          {
          T fx = T(x0+i*psx);
          fx *= fx;
          for (size_t j=0; j<=ny_dirty/2; ++j)
            {
            T fy = T(y0+j*psy);
            auto ws = wscreen(fx, fy*fy, T(w), adjoint);
            dirty2(i,j) = dirty(i,j)*ws; // lower left
            size_t i2 = nx_dirty-i, j2 = ny_dirty-j;
            if ((i>0)&&(i<i2))
              {
              dirty2(i2,j) = dirty(i2,j)*ws; // lower right
              if ((j>0)&&(j<j2))
                dirty2(i2,j2) = dirty(i2,j2)*ws; // upper right
              }
            if ((j>0)&&(j<j2))
              dirty2(i,j2) = dirty(i,j2)*ws; // upper left
            }
          }
        });
      }
  };

constexpr int logsquare=4;

template<typename T, typename T2=complex<T>> class Helper
  {
  private:
    const GridderConfig &gconf;
    int nu, nv, nsafe, supp;
    T beta;
    const T2 *grid_r;
    T2 *grid_w;
    int su, sv;
    int iu0, iv0; // start index of the current visibility
    int bu0, bv0; // start index of the current buffer

    vector<T2> rbuf, wbuf;
    bool do_w_gridding;
    double w0, xdw;
    size_t nexp;
    size_t nvecs;
    vector<std::mutex> &locks;

    void dump() const
      {
      if (bu0<-nsafe) return; // nothing written into buffer yet

      int idxu = (bu0+nu)%nu;
      int idxv0 = (bv0+nv)%nv;
      for (int iu=0; iu<su; ++iu)
        {
        int idxv = idxv0;
        {
        std::lock_guard<std::mutex> lock(locks[idxu]);
        for (int iv=0; iv<sv; ++iv)
          {
          grid_w[idxu*nv + idxv] += wbuf[iu*sv + iv];
          if (++idxv>=nv) idxv=0;
          }
        }
        if (++idxu>=nu) idxu=0;
        }
      }

    void load()
      {
      int idxu = (bu0+nu)%nu;
      int idxv0 = (bv0+nv)%nv;
      for (int iu=0; iu<su; ++iu)
        {
        int idxv = idxv0;
        for (int iv=0; iv<sv; ++iv)
          {
          rbuf[iu*sv + iv] = grid_r[idxu*nv + idxv];
          if (++idxv>=nv) idxv=0;
          }
        if (++idxu>=nu) idxu=0;
        }
      }

  public:
    const T2 *p0r;
    T2 *p0w;
    T kernel[64] ALIGNED(64);
    static constexpr size_t vlen=native_simd<T>::size();

    Helper(const GridderConfig &gconf_, const T2 *grid_r_, T2 *grid_w_,
      vector<std::mutex> &locks_, double w0_=-1, double dw_=-1)
      : gconf(gconf_), nu(gconf.Nu()), nv(gconf.Nv()), nsafe(gconf.Nsafe()),
        supp(gconf.Supp()), beta(T(gconf.Beta())), grid_r(grid_r_),
        grid_w(grid_w_), su(2*nsafe+(1<<logsquare)), sv(2*nsafe+(1<<logsquare)),
        bu0(-1000000), bv0(-1000000),
        rbuf(su*sv*(grid_r!=nullptr),T(0)),
        wbuf(su*sv*(grid_w!=nullptr),T(0)),
        do_w_gridding(dw_>0),
        w0(w0_),
        xdw(T(1)/dw_),
        nexp(2*supp + do_w_gridding),
        nvecs(vlen*((nexp+vlen-1)/vlen)),
        locks(locks_)
      {}
    ~Helper() { if (grid_w) dump(); }

    int lineJump() const { return sv; }
    T Wfac() const { return kernel[2*supp]; }
    void prep(const UVW &in)
      {
      double u, v;
      gconf.getpix(in.u, in.v, u, v, iu0, iv0);
      double xsupp=2./supp;
      double x0 = xsupp*(iu0-u);
      double y0 = xsupp*(iv0-v);
      for (int i=0; i<supp; ++i)
        {
        kernel[i  ] = T(x0+i*xsupp);
        kernel[i+supp] = T(y0+i*xsupp);
        }
      if (do_w_gridding)
        kernel[2*supp] = min(T(1), T(xdw*xsupp*abs(w0-in.w)));
      for (size_t i=nexp; i<nvecs; ++i)
        kernel[i]=0;
      for (size_t i=0; i<nvecs; ++i)
        kernel[i] = exp(beta*(sqrt(T(1)-kernel[i]*kernel[i])-T(1)));
      if ((iu0<bu0) || (iv0<bv0) || (iu0+supp>bu0+su) || (iv0+supp>bv0+sv))
        {
        if (grid_w) { dump(); fill(wbuf.begin(), wbuf.end(), T(0)); }
        bu0=((((iu0+nsafe)>>logsquare)<<logsquare))-nsafe;
        bv0=((((iv0+nsafe)>>logsquare)<<logsquare))-nsafe;
        if (grid_r) load();
        }
      p0r = grid_r ? rbuf.data() + sv*(iu0-bu0) + iv0-bv0 : nullptr;
      p0w = grid_w ? wbuf.data() + sv*(iu0-bu0) + iv0-bv0 : nullptr;
      }
  };

template<class T, class Serv> class SubServ
  {
  private:
    const Serv &srv;
    const_mav<idx_t,1> subidx;

  public:
    SubServ(const Serv &orig, const const_mav<idx_t,1> &subidx_)
      : srv(orig), subidx(subidx_){}
    size_t Nvis() const { return subidx.size(); }
    const Baselines &getBaselines() const { return srv.getBaselines(); }
    UVW getCoord(size_t i) const
      { return srv.getCoord(subidx[i]); }
    complex<T> getVis(size_t i) const
      { return srv.getVis(subidx[i]); }
    idx_t getIdx(size_t i) const { return srv.getIdx(subidx[i]); }
    void setVis (size_t i, const complex<T> &v) const
      { srv.setVis(subidx[i], v); }
    void addVis (size_t i, const complex<T> &v) const
      { srv.addVis(subidx[i], v); }
  };

template<class T, class T2> class MsServ
  {
  private:
    const Baselines &baselines;
    const_mav<idx_t,1> idx;
    T2 ms;
    const_mav<T,2> wgt;
    size_t nvis;
    bool have_wgt;

  public:
    using Tsub = SubServ<T, MsServ>;

    MsServ(const Baselines &baselines_,
    const const_mav<idx_t,1> &idx_, T2 ms_, const const_mav<T,2> &wgt_)
      : baselines(baselines_), idx(idx_), ms(ms_), wgt(wgt_),
        nvis(idx.shape(0)), have_wgt(wgt.size()!=0)
      {
      checkShape(ms.shape(), {baselines.Nrows(), baselines.Nchannels()});
      if (have_wgt) checkShape(wgt.shape(), ms.shape());
      }
    Tsub getSubserv(const const_mav<idx_t,1> &subidx) const
      { return Tsub(*this, subidx); }
    size_t Nvis() const { return nvis; }
    const Baselines &getBaselines() const { return baselines; }
    UVW getCoord(size_t i) const
      { return baselines.effectiveCoord(idx(i)); }
    complex<T> getVis(size_t i) const
      {
      auto rc = baselines.getRowChan(idx(i));
      return have_wgt ? ms(rc.row, rc.chan)*wgt(rc.row, rc.chan)
                      : ms(rc.row, rc.chan);
      }
    idx_t getIdx(size_t i) const { return idx[i]; }
    void setVis (size_t i, const complex<T> &v) const
      {
      auto rc = baselines.getRowChan(idx(i));
      ms(rc.row, rc.chan) = have_wgt ? v*wgt(rc.row, rc.chan) : v;
      }
    void addVis (size_t i, const complex<T> &v) const
      {
      auto rc = baselines.getRowChan(idx(i));
      ms(rc.row, rc.chan) += have_wgt ? v*wgt(rc.row, rc.chan) : v;
      }
  };
template<class T, class T2> MsServ<T, T2> makeMsServ
  (const Baselines &baselines,
   const const_mav<idx_t,1> &idx, const T2 &ms, const const_mav<T,2> &wgt)
  { return MsServ<T, T2>(baselines, idx, ms, wgt); }

template<typename T, typename Serv> void x2grid_c
  (const GridderConfig &gconf, const Serv &srv, mav<complex<T>,2> &grid,
  double w0=-1, double dw=-1)
  {
  checkShape(grid.shape(), {gconf.Nu(), gconf.Nv()});
  MR_assert(grid.contiguous(), "grid is not contiguous");
  size_t supp = gconf.Supp();
  size_t nthreads = gconf.Nthreads();
  bool do_w_gridding = dw>0;
  vector<std::mutex> locks(gconf.Nu());

  size_t np = srv.Nvis();
  execGuided(np, nthreads, 100, 0.2, [&](Scheduler &sched)
    {
    Helper<T> hlp(gconf, nullptr, grid.data(), locks, w0, dw);
    int jump = hlp.lineJump();
    const T * MRUTIL_RESTRICT ku = hlp.kernel;
    const T * MRUTIL_RESTRICT kv = hlp.kernel+supp;

    while (auto rng=sched.getNext()) for(auto ipart=rng.lo; ipart<rng.hi; ++ipart)
      {
      UVW coord = srv.getCoord(ipart);
      auto flip = coord.FixW();
      hlp.prep(coord);
      auto * MRUTIL_RESTRICT ptr = hlp.p0w;
      auto v(srv.getVis(ipart));
      if (do_w_gridding) v*=hlp.Wfac();
      if (flip) v=conj(v);
      for (size_t cu=0; cu<supp; ++cu)
        {
        complex<T> tmp(v*ku[cu]);
        size_t cv=0;
        for (; cv+3<supp; cv+=4)
          {
          ptr[cv  ] += tmp*kv[cv  ];
          ptr[cv+1] += tmp*kv[cv+1];
          ptr[cv+2] += tmp*kv[cv+2];
          ptr[cv+3] += tmp*kv[cv+3];
          }
        for (; cv<supp; ++cv)
          ptr[cv] += tmp*kv[cv];
        ptr+=jump;
        }
      }
    });
  }

template<typename T, typename Serv> void grid2x_c
  (const GridderConfig &gconf, const const_mav<complex<T>,2> &grid,
  const Serv &srv, double w0=-1, double dw=-1)
  {
  checkShape(grid.shape(), {gconf.Nu(), gconf.Nv()});
  MR_assert(grid.contiguous(), "grid is not contiguous");
  size_t supp = gconf.Supp();
  size_t nthreads = gconf.Nthreads();
  bool do_w_gridding = dw>0;
  vector<std::mutex> locks(gconf.Nu());

  // Loop over sampling points
  size_t np = srv.Nvis();
  execGuided(np, nthreads, 1000, 0.5, [&](Scheduler &sched)
    {
    Helper<T> hlp(gconf, grid.data(), nullptr, locks, w0, dw);
    int jump = hlp.lineJump();
    const T * MRUTIL_RESTRICT ku = hlp.kernel;
    const T * MRUTIL_RESTRICT kv = hlp.kernel+supp;

    while (auto rng=sched.getNext()) for(auto ipart=rng.lo; ipart<rng.hi; ++ipart)
      {
      UVW coord = srv.getCoord(ipart);
      auto flip = coord.FixW();
      hlp.prep(coord);
      complex<T> r = 0;
      const auto * MRUTIL_RESTRICT ptr = hlp.p0r;
      for (size_t cu=0; cu<supp; ++cu)
        {
        complex<T> tmp(0);
        size_t cv=0;
        for (; cv+3<supp; cv+=4)
          tmp += ptr[cv  ]*kv[cv  ]
               + ptr[cv+1]*kv[cv+1]
               + ptr[cv+2]*kv[cv+2]
               + ptr[cv+3]*kv[cv+3];
        for (; cv<supp; ++cv)
          tmp += ptr[cv] * kv[cv];
        r += tmp*ku[cu];
        ptr += jump;
        }
      if (flip) r=conj(r);
      if (do_w_gridding) r*=hlp.Wfac();
      srv.addVis(ipart, r);
      }
    });
  }

template<typename T> void apply_global_corrections(const GridderConfig &gconf,
  const mav<T,2> &dirty, const ES_Kernel &kernel, double dw, bool divide_by_n)
  {
  auto nx_dirty=gconf.Nxdirty();
  auto ny_dirty=gconf.Nydirty();
  size_t nthreads = gconf.Nthreads();
  auto psx=gconf.Pixsize_x();
  auto psy=gconf.Pixsize_y();
  double x0 = -0.5*nx_dirty*psx,
         y0 = -0.5*ny_dirty*psy;
  auto cfu = correction_factors(gconf.Nu(), gconf.Ofactor(),
                                nx_dirty/2+1, gconf.Supp(), nthreads);
  auto cfv = correction_factors(gconf.Nv(), gconf.Ofactor(),
                                ny_dirty/2+1, gconf.Supp(), nthreads);
  execStatic(nx_dirty/2+1, nthreads, 0, [&](Scheduler &sched)
    {
    while (auto rng=sched.getNext()) for(auto i=rng.lo; i<rng.hi; ++i)
      {
      auto fx = T(x0+i*psx);
      fx *= fx;
      for (size_t j=0; j<=ny_dirty/2; ++j)
        {
        auto fy = T(y0+j*psy);
        fy*=fy;
        T fct = 0;
        auto tmp = 1-fx-fy;
        if (tmp>=0)
          {
          auto nm1 = (-fx-fy)/(sqrt(tmp)+1); // accurate form of sqrt(1-x-y)-1
          fct = T(kernel.corfac(nm1*dw));
          if (divide_by_n)
            fct /= nm1+1;
          }
        else // beyond the horizon, don't really know what to do here
          {
          if (divide_by_n)
            fct=0;
          else
            {
            auto nm1 = sqrt(-tmp)-1;
            fct = T(kernel.corfac(nm1*dw));
            }
          }
        fct *= T(cfu[nx_dirty/2-i]*cfv[ny_dirty/2-j]);
        size_t i2 = nx_dirty-i, j2 = ny_dirty-j;
        dirty(i,j)*=fct;
        if ((i>0)&&(i<i2))
          {
          dirty(i2,j)*=fct;
          if ((j>0)&&(j<j2))
            dirty(i2,j2)*=fct;
          }
        if ((j>0)&&(j<j2))
          dirty(i,j2)*=fct;
        }
      }
    });
  }

template<typename Serv> class WgridHelper
  {
  private:
    const Serv &srv;
    double wmin, dw;
    size_t nplanes, supp;
    vector<vector<idx_t>> minplane;
    size_t verbosity;

    int curplane;
    vector<idx_t> subidx;

    static void wminmax(const Serv &srv, double &wmin, double &wmax)
      {
      size_t nvis = srv.Nvis();

      wmin= 1e38;
      wmax=-1e38;
      // FIXME maybe this can be done more intelligently
      for (size_t ipart=0; ipart<nvis; ++ipart)
        {
        auto wval = abs(srv.getCoord(ipart).w);
        wmin = min(wmin,wval);
        wmax = max(wmax,wval);
        }
      }

    template<typename T> static void update_idx(vector<T> &v, const vector<T> &add,
      const vector<T> &del)
      {
      MR_assert(v.size()>=del.size(), "must not happen");
      vector<T> res;
      res.reserve((v.size()+add.size())-del.size());
      auto iin=v.begin(), ein=v.end();
      auto iadd=add.begin(), eadd=add.end();
      auto irem=del.begin(), erem=del.end();

      while(iin!=ein)
        {
        if ((irem!=erem) && (*iin==*irem))
          {  ++irem; ++iin; } // skip removed entry
        else if ((iadd!=eadd) && (*iadd<*iin))
           res.push_back(*(iadd++)); // add new entry
        else
          res.push_back(*(iin++));
        }
      MR_assert(irem==erem, "must not happen");
      while(iadd!=eadd)
        res.push_back(*(iadd++));
      MR_assert(res.size()==(v.size()+add.size())-del.size(), "must not happen");
      v.swap(res);
      }

  public:
    WgridHelper(const GridderConfig &gconf, const Serv &srv_, size_t verbosity_)
      : srv(srv_), verbosity(verbosity_), curplane(-1)
      {
      size_t nvis = srv.Nvis();
      size_t nthreads = gconf.Nthreads();
      double wmax;

      wminmax(srv, wmin, wmax);
      if (verbosity>0) cout << "Using " << nthreads << " thread"
                            << ((nthreads!=1) ? "s" : "") << endl;
      if (verbosity>0) cout << "W range: " << wmin << " to " << wmax << endl;

      double x0 = -0.5*gconf.Nxdirty()*gconf.Pixsize_x(),
             y0 = -0.5*gconf.Nydirty()*gconf.Pixsize_y();
      double nmin = sqrt(max(1.-x0*x0-y0*y0,0.))-1.;
      if (x0*x0+y0*y0>1.)
        nmin = -sqrt(abs(1.-x0*x0-y0*y0))-1.;
      dw = 0.25/abs(nmin);
      nplanes = size_t((wmax-wmin)/dw+2);
      dw = (1.+1e-13)*(wmax-wmin)/(nplanes-1);

      supp = gconf.Supp();
      wmin -= (0.5*supp-1)*dw;
      wmax += (0.5*supp-1)*dw;
      nplanes += supp-2;
      if (verbosity>0) cout << "Kernel support: " << supp << endl;
      if (verbosity>0) cout << "nplanes: " << nplanes << endl;

      minplane.resize(nplanes);
#if 0
      // extra short, but potentially inefficient version:
      for (size_t ipart=0; ipart<nvis; ++ipart)
        {
        int plane0 = max(0,int(1+(abs(srv.getCoord(ipart).w)-(0.5*supp*dw)-wmin)/dw));
        minplane[plane0].push_back(idx_t(ipart));
        }
#else
      // more efficient: precalculate final vector sizes and avoid reallocations
      vector<size_t> cnt(nplanes,0);
      for(size_t ipart=0; ipart<nvis; ++ipart)
        {
        int plane0 = max(0,int(1+(abs(srv.getCoord(ipart).w)-(0.5*supp*dw)-wmin)/dw));
        ++cnt[plane0];
        }

      // fill minplane
      for (size_t j=0; j<nplanes; ++j)
        minplane[j].resize(cnt[j]);
      vector<size_t> ofs(nplanes, 0);
      for (size_t ipart=0; ipart<nvis; ++ipart)
        {
        int plane0 = max(0,int(1+(abs(srv.getCoord(ipart).w)-(0.5*supp*dw)-wmin)/dw));
        minplane[plane0][ofs[plane0]++]=idx_t(ipart);
        }
#endif
      }

    typename Serv::Tsub getSubserv() const
      {
      auto subidx2 = const_mav<idx_t, 1>(subidx.data(), {subidx.size()});
      return srv.getSubserv(subidx2);
      }
    double W() const { return wmin+curplane*dw; }
    size_t Nvis() const { return subidx.size(); }
    double DW() const { return dw; }
    bool advance()
      {
      if (++curplane>=int(nplanes)) return false;
      update_idx(subidx, minplane[curplane], curplane>=int(supp) ? minplane[curplane-supp] : vector<idx_t>());
      if (verbosity>1)
        cout << "Working on plane " << curplane << " containing " << subidx.size()
             << " visibilities" << endl;
      return true;
      }
  };

template<typename T, typename Serv> void x2dirty(
  const GridderConfig &gconf, const Serv &srv, const mav<T,2> &dirty,
  bool do_wstacking, size_t verbosity)
  {
  if (do_wstacking)
    {
    size_t nthreads = gconf.Nthreads();
    if (verbosity>0) cout << "Gridding using improved w-stacking" << endl;
    WgridHelper<Serv> hlp(gconf, srv, verbosity);
    double dw = hlp.DW();
    dirty.fill(0);
    tmpStorage<complex<T>,2> grid_({gconf.Nu(),gconf.Nv()});
    auto grid=grid_.getMav();
    while(hlp.advance())  // iterate over w planes
      {
      if (hlp.Nvis()==0) continue;
      grid.fill(0);
      x2grid_c(gconf, hlp.getSubserv(), grid, hlp.W(), dw);
      gconf.grid2dirty_c_overwrite_wscreen_add(grid, dirty, T(hlp.W()));
      }
    // correct for w gridding etc.
    apply_global_corrections(gconf, dirty, ES_Kernel(gconf.Supp(), gconf.Ofactor(), nthreads), dw, true);
    }
  else
    {
    if (verbosity>0)
      cout << "Gridding without w-stacking: " << srv.Nvis()
           << " visibilities" << endl;
    if (verbosity>0) cout << "Using " << gconf.Nthreads() << " threads" << endl;

    tmpStorage<complex<T>,2> grid_({gconf.Nu(), gconf.Nv()});
    auto grid=grid_.getMav();
    grid_.fill(0.);
    x2grid_c(gconf, srv, grid);
    tmpStorage<T,2> rgrid_(grid.shape());
    auto rgrid=rgrid_.getMav();
    complex2hartley(cmav(grid), rgrid, gconf.Nthreads());
    gconf.grid2dirty(cmav(rgrid), dirty);
    }
  }

template<typename T, typename Serv> void dirty2x(
  const GridderConfig &gconf,  const const_mav<T,2> &dirty,
  const Serv &srv, bool do_wstacking, size_t verbosity)
  {
  if (do_wstacking)
    {
    size_t nx_dirty=gconf.Nxdirty(), ny_dirty=gconf.Nydirty();
    size_t nthreads = gconf.Nthreads();
    if (verbosity>0) cout << "Degridding using improved w-stacking" << endl;
    WgridHelper<Serv> hlp(gconf, srv, verbosity);
    double dw = hlp.DW();
    tmpStorage<T,2> tdirty_({nx_dirty,ny_dirty});
    auto tdirty=tdirty_.getMav();
    for (size_t i=0; i<nx_dirty; ++i)
      for (size_t j=0; j<ny_dirty; ++j)
        tdirty(i,j) = dirty(i,j);
    // correct for w gridding etc.
    apply_global_corrections(gconf, tdirty, ES_Kernel(gconf.Supp(), gconf.Ofactor(), nthreads), dw, true);
    tmpStorage<complex<T>,2> grid_({gconf.Nu(),gconf.Nv()});
    auto grid=grid_.getMav();
    while(hlp.advance())  // iterate over w planes
      {
      if (hlp.Nvis()==0) continue;
      gconf.dirty2grid_c_wscreen(cmav(tdirty), grid, T(hlp.W()));
      grid2x_c(gconf, cmav(grid), hlp.getSubserv(), hlp.W(), dw);
      }
    }
  else
    {
    if (verbosity>0)
      cout << "Degridding without w-stacking: " << srv.Nvis()
           << " visibilities" << endl;
    if (verbosity>0) cout << "Using " << gconf.Nthreads() << " threads" << endl;

    tmpStorage<T,2> grid_({gconf.Nu(), gconf.Nv()});
    auto grid=grid_.getMav();
    gconf.dirty2grid(dirty, grid);
    tmpStorage<complex<T>,2> grid2_(grid.shape());
    auto grid2=grid2_.getMav();
    hartley2complex(cmav(grid), grid2, gconf.Nthreads());
    grid2x_c(gconf, cmav(grid2), srv);
    }
  }

void calc_share(size_t nshares, size_t myshare, size_t nwork, size_t &lo,
  size_t &hi)
  {
  size_t nbase = nwork/nshares;
  size_t additional = nwork%nshares;
  lo = myshare*nbase + ((myshare<additional) ? myshare : additional);
  hi = lo+nbase+(myshare<additional);
  }


template<typename T> vector<idx_t> getWgtIndices(const Baselines &baselines,
  const GridderConfig &gconf, const const_mav<T,2> &wgt,
  const const_mav<complex<T>,2> &ms)
  {
  size_t nrow=baselines.Nrows(),
         nchan=baselines.Nchannels(),
         nsafe=gconf.Nsafe();
  bool have_wgt=wgt.size()!=0;
  if (have_wgt) checkShape(wgt.shape(),{nrow,nchan});
  bool have_ms=ms.size()!=0;
  if (have_ms) checkShape(ms.shape(), {nrow,nchan});
  constexpr int side=1<<logsquare;
  size_t nbu = (gconf.Nu()+1+side-1) >> logsquare,
         nbv = (gconf.Nv()+1+side-1) >> logsquare;
  vector<idx_t> acc(nbu*nbv+1,0);
  vector<idx_t> tmp(nrow*nchan);

  for (idx_t irow=0, idx=0; irow<nrow; ++irow)
    for (idx_t ichan=0; ichan<nchan; ++ichan, ++idx)
      if (((!have_ms ) || (norm(ms (irow,ichan))!=0)) &&
          ((!have_wgt) || (wgt(irow,ichan)!=0)))
        {
        auto uvw = baselines.effectiveCoord(RowChan{irow,idx_t(ichan)});
        if (uvw.w<0) uvw.Flip();
        double u, v;
        int iu0, iv0;
        gconf.getpix(uvw.u, uvw.v, u, v, iu0, iv0);
        iu0 = (iu0+nsafe)>>logsquare;
        iv0 = (iv0+nsafe)>>logsquare;
        ++acc[nbv*iu0 + iv0 + 1];
        tmp[idx] = nbv*iu0 + iv0;
        }
      else
        tmp[idx] = ~idx_t(0);

  for (size_t i=1; i<acc.size(); ++i)
    acc[i] += acc[i-1];

  vector<idx_t> res(acc.back());
  for (size_t irow=0, idx=0; irow<nrow; ++irow)
    for (size_t ichan=0; ichan<nchan; ++ichan, ++idx)
      if (tmp[idx]!=(~idx_t(0)))
        res[acc[tmp[idx]]++] = baselines.getIdx(irow, ichan);
  return res;
  }

template<typename T> void ms2dirty_general(const const_mav<double,2> &uvw,
  const const_mav<double,1> &freq, const const_mav<complex<T>,2> &ms,
  const const_mav<T,2> &wgt, double pixsize_x, double pixsize_y, size_t nu, size_t nv, double epsilon,
  bool do_wstacking, size_t nthreads, const mav<T,2> &dirty, size_t verbosity,
  bool negate_v=false)
  {
  Baselines baselines(uvw, freq, negate_v);
  GridderConfig gconf(dirty.shape(0), dirty.shape(1), nu, nv, epsilon, pixsize_x, pixsize_y, nthreads);
  auto idx = getWgtIndices(baselines, gconf, wgt, ms);
  auto idx2 = const_mav<idx_t,1>(idx.data(),{idx.size()});
  x2dirty(gconf, makeMsServ(baselines,idx2,ms,wgt), dirty, do_wstacking, verbosity);
  }
template<typename T> void ms2dirty(const const_mav<double,2> &uvw,
  const const_mav<double,1> &freq, const const_mav<complex<T>,2> &ms,
  const const_mav<T,2> &wgt, double pixsize_x, double pixsize_y, double epsilon,
  bool do_wstacking, size_t nthreads, const mav<T,2> &dirty, size_t verbosity)
  {
  ms2dirty_general(uvw, freq, ms, wgt, pixsize_x, pixsize_y,
    2*dirty.shape(0), 2*dirty.shape(1), epsilon, do_wstacking, nthreads,
    dirty, verbosity);
  }

template<typename T> void dirty2ms_general(const const_mav<double,2> &uvw,
  const const_mav<double,1> &freq, const const_mav<T,2> &dirty,
  const const_mav<T,2> &wgt, double pixsize_x, double pixsize_y, size_t nu, size_t nv,double epsilon,
  bool do_wstacking, size_t nthreads, const mav<complex<T>,2> &ms,
  size_t verbosity, bool negate_v=false)
  {
  Baselines baselines(uvw, freq, negate_v);
  GridderConfig gconf(dirty.shape(0), dirty.shape(1), nu, nv, epsilon, pixsize_x, pixsize_y, nthreads);
  const_mav<complex<T>,2> null_ms(nullptr, {0,0});
  auto idx = getWgtIndices(baselines, gconf, wgt, null_ms);
  auto idx2 = const_mav<idx_t,1>(idx.data(),{idx.size()});
  ms.fill(0);
  dirty2x(gconf, dirty, makeMsServ(baselines,idx2,ms,wgt), do_wstacking, verbosity);
  }
template<typename T> void dirty2ms(const const_mav<double,2> &uvw,
  const const_mav<double,1> &freq, const const_mav<T,2> &dirty,
  const const_mav<T,2> &wgt, double pixsize_x, double pixsize_y, double epsilon,
  bool do_wstacking, size_t nthreads, const mav<complex<T>,2> &ms, size_t verbosity)
  {
  dirty2ms_general(uvw, freq, dirty, wgt, pixsize_x, pixsize_y,
    2*dirty.shape(0), 2*dirty.shape(1), epsilon, do_wstacking, nthreads, ms,
    verbosity);
  }

} // namespace detail

// public names
using detail::ms2dirty_general;
using detail::dirty2ms_general;
using detail::ms2dirty;
using detail::dirty2ms;
} // namespace gridder

#endif