aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetInterp2/src/fasterp.cc
blob: 6c5880f44973ad92d2e09ba8c197108d973d3172 (plain)
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
1272
1273
1274
1275
1276
#include <cctk.h>
#include <cctk_Parameters.h>

#include <algorithm>
#include <cassert>
#include <cmath>
#include <cstddef>
#include <cstdlib>
#include <iostream>
#include <sstream>

#ifdef CCTK_MPI
#  include <mpi.h>
#else
#  include "nompi.h"
#endif

#include <carpet.hh>
#include <vect.hh>

#include "fasterp.hh"



namespace CarpetInterp2 {
  
  
  
  int const ipoison = -1234567890;
  CCTK_REAL const poison = -1.0e+12;
  
  int get_poison (int const &) CCTK_ATTRIBUTE_CONST;
  int get_poison (int const &) { return ipoison; }
  CCTK_REAL get_poison (CCTK_REAL const &) CCTK_ATTRIBUTE_CONST;
  CCTK_REAL get_poison (CCTK_REAL const &) { return poison; }
  
  template <typename T>
  void fill (vector<T> & v, T const & val)
  {
    // fill (v.begin(), v.end(), val);
#pragma omp parallel for
    for (int n=0; n<int(v.size()); ++n) {
      v.AT(n) = val;
    }
  }
  
  template <typename T>
  void fill_with_poison (vector<T> & v)
  {
#ifndef NDEBUG
    T dummy;
    fill (v, get_poison(dummy));
    // fill (v.begin(), v.end(), get_poison(dummy));
#endif
  }
  
  template <typename T>
  void check_for_poison (vector<T> const & v)
  {
#pragma omp parallel for
    for (int n=0; n<int(v.size()); ++n) {
      T dummy;
      assert (v.AT(n) != get_poison(dummy));
    }
  }
  
  
  
  // Create an MPI datatype for location information
  MPI_Datatype
  fasterp_iloc_t::mpi_datatype ()
  {
    static bool initialised = false;
    static MPI_Datatype newtype;
    if (not initialised) {
      static fasterp_iloc_t s;
#define ENTRY(type, name)                                               \
      {                                                                 \
        sizeof s.name / sizeof(type), /* count elements */              \
          (char*)&s.name - (char*)&s, /* offsetof doesn't work (why?) */ \
          dist::mpi_datatype<type>(), /* find MPI datatype */           \
          STRINGIFY(name), /* field name */                             \
          STRINGIFY(type), /* type name */                              \
          }
      dist::mpi_struct_descr_t const descr[] = {
        ENTRY(int, mrc),
#ifdef CARPETINTERP2_CHECK
        ENTRY(int, pn),
        ENTRY(int, ipos),
        ENTRY(int, ind),
#endif
        ENTRY(int, ind3d),
        ENTRY(CCTK_REAL, offset),
        {1, sizeof(s), MPI_UB, "MPI_UB", "MPI_UB"}
      };
#undef ENTRY
      newtype =
        dist::create_mpi_datatype (sizeof descr / sizeof descr[0], descr,
                                   "fasterp_iloc_t::mpi_datatype", sizeof s);
      initialised = true;
    }
    return newtype;
  }
  
  void
  fasterp_iloc_t::
  output (ostream& os)
    const
  {
    os << "fasterp_iloc_t{"
       << "mrc=" << mrc << ","
#ifdef CARPETINTERP2_CHECK
       << "pn=" << pn << ","
       << "ipos=" << ipos << ","
       << "ind=" << ind << ","
#endif
       << "ind3d=" << ind3d << ","
       << "offset=" << offset << "}";
  }
  
  
  
  //////////////////////////////////////////////////////////////////////////////
  //////////////////////////////////////////////////////////////////////////////
  //////////////////////////////////////////////////////////////////////////////
  
  
  
  int mrc_t::maps       = -1;
  int mrc_t::reflevels  = -1;
  int mrc_t::components = -1;
  
  void
  mrc_t::
  determine_mrc_info ()
  {
    maps = Carpet::maps;
    reflevels = Carpet::reflevels;
    components = 0;
    for (int m=0; m<maps; ++m) {
      for (int rl=0; rl<reflevels; ++rl) {
        int const ncomps = Carpet::vhh.AT(m)->components (rl);
        components = max (components, ncomps);
      }
    }
  }
  
  mrc_t::
  mrc_t (int ind)
  {
    assert (ind >= 0);
    int const ind0 = ind;
    c = ind % components;
    ind /= components;
    rl = ind % reflevels;
    ind /= reflevels;
    m = ind % maps;
    ind /= maps;
    assert (ind == 0);
    assert (get_ind() == ind0);
  }
  
  void
  mrc_t::
  output (ostream& os)
    const
  {
    os << "mrc{m=" << m << ",rl=" << rl << ",c=" << c << "}";
  }
  
  
  
  void
  pn_t::
  output (ostream& os)
    const
  {
    os << "pn{p=" << p << ",n=" << n << "}";
  }
  
  
  
  //////////////////////////////////////////////////////////////////////////////
  //////////////////////////////////////////////////////////////////////////////
  //////////////////////////////////////////////////////////////////////////////
  
  
  
  // Calculate the coefficients for one interpolation stencil
  // TODO: Could templatify this function on the order to improve
  // efficiency
  int
  fasterp_src_loc_t::
  calc_stencil (fasterp_iloc_t const & iloc,
                ivect const & lsh,
                int const order)
  {
    assert (order <= max_order);
    CCTK_REAL const eps = 1.0e-12;
    
#ifdef CARPETINTERP2_CHECK
    mrc  = iloc.mrc;
    pn   = iloc.pn;
    ipos = iloc.ipos;
    
    // Save lsh
    saved_lsh = lsh;
#endif
    
#ifndef NDEBUG
    // Poison all coefficients
    for (int d=0; d<dim; ++d) {
      for (int n=0; n<=max_order; ++n) {
        coeffs[d][n] = poison;
      }
    }
#endif
    
    if (order == 0) {
      // Special case
      for (int d=0; d<dim; ++d) {
        exact[d] = true;
      }
      
#ifdef CARPETINTERP2_CHECK
      ind   = iloc.ind;
#endif
      ind3d = iloc.ind3d;
      return 0;
    }
    
    // Choose stencil anchor (first stencil point)
    ivect iorigin;
    if (order % 2 == 0) {
      iorigin = - ivect(order/2);
    } else {
      // Potentially shift the stencil anchor for odd interpolation
      // orders (i.e., for even numbers of stencil points)
      ivect const ioffset (iloc.offset < CCTK_REAL(0.0));
      iorigin = - ivect((order-1)/2) - ioffset;
    }
    rvect const offset = iloc.offset - rvect(iorigin);
    // Ensure that the stencil is centred
    assert (all (offset >= CCTK_REAL(0.5)*(order-1) and
                 offset < CCTK_REAL(0.5)*(order+1)));
    
    for (int d=0; d<dim; ++d) {
      // C_n = PRODUCT_m,m!=n [(x - x_m) / (x_n - x_m)]
      CCTK_REAL const x = offset[d];
      // round is not available with PGI compilers
      // CCTK_REAL const rx = round(x);
      CCTK_REAL const rx = floor(x+0.5);
      if (abs(x - rx) < eps * (1.0 + abs(x))) {
        // The interpolation point coincides with a grid point; no
        // interpolation is necessary (this is a special case)
        iorigin[d] += int(rx);
        exact[d] = true;
      } else {
        for (int n=0; n<=order; ++n) {
          CCTK_REAL const xn = n;
          CCTK_REAL c = 1.0;
          for (int m=0; m<=order; ++m) {
            if (m != n) {
              CCTK_REAL const xm = m;
              c *= (x - xm) / (xn - xm);
            }
          }
          coeffs[d][n] = c;
        }
        exact[d] = false;
        // The sum of the coefficients must be one
        CCTK_REAL s = 0.0;
        for (int n=0; n<=order; ++n) {
          s += coeffs[d][n];
        }
        assert (abs(s - 1.0) <= eps);
      }
    }
    
    // Set 3D location of stencil anchor
#ifdef CARPETINTERP2_CHECK
    ind = iloc.ind + iorigin;
    if (not (all (ind>=0 and ind+either(exact,0,order)<lsh))) {
      stringstream buf;
      buf << "*this=" << *this << " iloc=" << iloc << " "
          << "lsh=" << lsh << " order=" << order;
      CCTK_VWarn (CCTK_WARN_ALERT, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Array access out of bounds: %s. Most likely, the interpolation point is too close to an outer or to a symmetry/inter-patch boundary. More information follows...",
                  buf.str().c_str());
      return -1;
    }
#endif
    ind3d = iloc.ind3d + index(lsh, iorigin);
#ifdef CARPETINTERP2_CHECK
    assert (index(lsh,ind) == ind3d);
#endif
    return 0;
  }
  
  
  
  // Interpolate a set of variables at the same location, with a given
  // set of interpolation coefficients.  The template mechanism allows
  // the order in each direction to be different; in particular, the
  // order can be 0 if the interpolation location coincides with a
  // source grid point.  For interpolation to order O, this function
  // should be instantiated eight times, with On=O and On=0, for
  // n=[0,1,2].  We hope that the compiler optimises the for loops
  // away if On=0.
  template <int O0, int O1, int O2>
  void
  fasterp_src_loc_t::
  interpolate (ivect const & lsh,
               vector<CCTK_REAL const *> const & varptrs,
               CCTK_REAL * restrict const vals)
    const
  {
#ifdef CARPETINTERP2_CHECK
    assert (all (lsh == saved_lsh));
#endif
    assert (O0 == 0 or not exact[0]);
    assert (O1 == 0 or not exact[1]);
    assert (O2 == 0 or not exact[2]);
    
    size_t const di = 1;
    size_t const dj = di * lsh[0];
    size_t const dk = dj * lsh[1];
    
#ifdef CARPETINTERP2_CHECK
    assert (all (ind>=0 and ind+ivect(O0,O1,O2)<lsh));
    assert (ind3d == index(lsh,ind));
    assert (int(di*ind[0] + dj*ind[1] + dk*ind[2]) == index(lsh,ind));
#endif
    
    for (size_t v=0; v<varptrs.size(); ++v) {
      vals[v] = 0.0;
    }
    
    for (size_t k=0; k<=O2; ++k) {
      assert (O2 == 0 or coeffs[2][k] != poison);
      CCTK_REAL const coeff_k = (O2==0 ? 1.0 : coeffs[2][k]);
      for (size_t j=0; j<=O1; ++j) {
        assert (O1 == 0 or coeffs[1][j] != poison);
        CCTK_REAL const coeff_jk = coeff_k * (O1==0 ? 1.0 : coeffs[1][j]);
        for (size_t i=0; i<=O0; ++i) {
          assert (O0 == 0 or coeffs[0][i] != poison);
          CCTK_REAL const coeff_ijk = coeff_jk * (O0==0 ? 1.0 : coeffs[0][i]);
          
          for (size_t v=0; v<varptrs.size(); ++v) {
            vals[v] += coeff_ijk * varptrs.AT(v)[ind3d + i*di + j*dj + k*dk];
          } // for v
          
        }
      }
    }
#ifdef CARPETINTERP2_CHECK
#  if 0
    for (size_t v=0; v<varptrs.size(); ++v) {
      vals[v] = ind[0] + 1000 * (ind[1] + 1000 * ind[2]);
    } // for v
#  endif
#  if 0
    for (size_t v=0; v<varptrs.size(); ++v) {
      vals[v] = pn.p * 1000000 + pn.n;
    } // for v
#  endif
#endif
  }
  
  
  
  // Interpolate a set of variables at the same location, with a given
  // set of interpolation coefficients.  This calls the specialised
  // interpolation function, depending on whether interpolation is
  // exact in each of the directions.
  template <int O>
  void
  fasterp_src_loc_t::
  interpolate (ivect const & lsh,
               vector<CCTK_REAL const *> const & varptrs,
               CCTK_REAL * restrict const vals)
    const
  {
    int const Z = 0;
    if (exact[2]) {
      if (exact[1]) {
        if (exact[0]) {
          interpolate<Z,Z,Z> (lsh, varptrs, vals);
        } else {
          interpolate<O,Z,Z> (lsh, varptrs, vals);
        }
      } else {
        if (exact[0]) {
          interpolate<Z,O,Z> (lsh, varptrs, vals);
        } else {
          interpolate<O,O,Z> (lsh, varptrs, vals);
        }
      }
    } else {
      if (exact[1]) {
        if (exact[0]) {
          interpolate<Z,Z,O> (lsh, varptrs, vals);
        } else {
          interpolate<O,Z,O> (lsh, varptrs, vals);
        }
      } else {
        if (exact[0]) {
          interpolate<Z,O,O> (lsh, varptrs, vals);
        } else {
          interpolate<O,O,O> (lsh, varptrs, vals);
        }
      }
    }
  }
  
  
  
  // Interpolate a set of variables at the same location, with a given
  // set of interpolation coefficients.  This calls a specialised
  // interpolation function, depending on whether interpolation is
  // exact in each of the directions.
  void
  fasterp_src_loc_t::
  interpolate (ivect const & lsh,
               int const order,
               vector<CCTK_REAL const *> const & varptrs,
               CCTK_REAL * restrict const vals)
    const
  {
    switch (order) {
    case  0: interpolate< 0> (lsh, varptrs, vals); break;
    case  1: interpolate< 1> (lsh, varptrs, vals); break;
    case  2: interpolate< 2> (lsh, varptrs, vals); break;
    case  3: interpolate< 3> (lsh, varptrs, vals); break;
    case  4: interpolate< 4> (lsh, varptrs, vals); break;
    case  5: interpolate< 5> (lsh, varptrs, vals); break;
    case  6: interpolate< 6> (lsh, varptrs, vals); break;
    case  7: interpolate< 7> (lsh, varptrs, vals); break;
    case  8: interpolate< 8> (lsh, varptrs, vals); break;
    case  9: interpolate< 9> (lsh, varptrs, vals); break;
    case 10: interpolate<10> (lsh, varptrs, vals); break;
    case 11: interpolate<11> (lsh, varptrs, vals); break;
    default:
      // Add higher orders here as desired
      CCTK_WARN (CCTK_WARN_ABORT,
                 "Interpolation orders larger than 11 are not yet implemented");
      assert (0);
    }
  }
  
  
  
  void
  fasterp_src_loc_t::
  output (ostream& os)
    const
  {
    os << "fasterp_src_loc_t{";
    os << "coeffs=[";
    for (int d=0; d<dim; ++d) {
      if (d>0) os << ",";
      os << "[";
      for (int n=0; n<=max_order; ++n) {
        if (n>0) os << ",";
        os << coeffs[d][n];
      }
      os << "]";
    }
    os << "],";
    os << "exact=" << exact << ",";
#ifdef CARPETINTERP2_CHECK
    os << "pn=" << pn << ",";
    os << "mrc=" << mrc << ",";
    os << "ipos=" << ipos << ",";
    os << "ind=" << ind << ",";
#endif
    os << "ind3d=" << ind3d;
#ifdef CARPETINTERP2_CHECK
    os << "," << "saved_lsh=" << saved_lsh;
#endif
    os << "}";
  }
  
  
  
  //////////////////////////////////////////////////////////////////////////////
  //////////////////////////////////////////////////////////////////////////////
  //////////////////////////////////////////////////////////////////////////////
  
  
  
  // Set up an interpolation starting from global coordinates
  fasterp_setup_t::
  fasterp_setup_t (cGH const * restrict const cctkGH,
                   fasterp_glocs_t const & locations,
                   int const order_)
    : order (order_),
      reflevel (Carpet::reflevel),
      regridding_epoch (reflevel == -1 ?
                        Carpet::regridding_epoch :
                        Carpet::level_regridding_epochs.AT(reflevel))
  {
    // Some global properties
    int const npoints = locations.size();
    
    
    
    // Calculate patch numbers and local coordinates
    fasterp_llocs_t local_locations (npoints);
    
    if (CCTK_IsFunctionAliased ("MultiPatch_GlobalToLocal")) {
      // This is a multi-patch simulation: convert global to local
      // coordinates
      
      CCTK_REAL const * coords[dim];
      CCTK_REAL       * local_coords[dim];
      for (int d=0; d<dim; ++d) {
        coords[d]       = &locations.coords[d].front();
        local_coords[d] = &local_locations.coords[d].front();
      }
      
      int const ierr = MultiPatch_GlobalToLocal
        (cctkGH, dim, npoints,
         coords,
         &local_locations.maps.front(), local_coords, NULL, NULL);
      assert (not ierr);
      
    } else {
      // This is a single-patch simulation
      
      // TODO: Optimise this: don't copy coordinates, and don't store
      // map numbers if there is only one map.
      for (int n=0; n<npoints; ++n) {
        int const m = 0;
        local_locations.maps.AT(n) = m;
        for (int d=0; d<dim; ++d) {
          local_locations.coords[d].AT(n) = locations.coords[d].AT(n);
        }
      }
      
    } // if not multi-patch
    
    setup (cctkGH, local_locations);
  }
  
  
  
  // Set up an interpolation starting from local coordinates
  fasterp_setup_t::
  fasterp_setup_t (cGH const * restrict const cctkGH,
                   fasterp_llocs_t const & locations,
                   int const order_)
    : order (order_)
  {
    setup (cctkGH, locations);
  }
  
  
  
  // Helper for setting up an interpolation
  void
  fasterp_setup_t::
  setup (cGH const * restrict const cctkGH,
         fasterp_llocs_t const & locations)
  {
    DECLARE_CCTK_PARAMETERS;
    
    if (verbose) CCTK_VInfo (CCTK_THORNSTRING,
                             "Setting up interpolation for %d grid points",
                             int(locations.size()));
    
    if (order < 0) {
      CCTK_WARN (CCTK_WARN_ABORT,
                 "Interpolation order must be non-negative");
    }
    if (order > max_order) {
      CCTK_VWarn (CCTK_WARN_ABORT, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Interpolation order cannot be larger than max_order=%d; "
                  "order=%d was requested.  "
                  "(You can increase the compile time constant max_order "
                  "in thorn CarpetInterp2.)",
                  max_order, order);
    }
    
    // Some global properties
    int const npoints = locations.size();
    int const nprocs = CCTK_nProcs (cctkGH);
    // int const myproc = CCTK_MyProc (cctkGH);
    
    mrc_t::determine_mrc_info();
    int const maxmrc = mrc_t::get_max_ind();
    
    MPI_Comm & comm_world = * (MPI_Comm *) GetMPICommWorld (cctkGH);
    
    
    
    // Obtain the coordinate ranges for all patches
    vector<rvect> lower (Carpet::maps);
    vector<rvect> upper (Carpet::maps);
    vector<rvect> delta (Carpet::maps); // spacing on finest possible grid
    vector<rvect> idelta (Carpet::maps); // inverse spacing
    for (int m=0; m<Carpet::maps; ++m) {
      jvect gsh;
      int const ierr =
        GetCoordRange (cctkGH, m, Carpet::mglevel, dim,
                       & gsh[0],
                       & lower.AT(m)[0], & upper.AT(m)[0], & delta.AT(m)[0]);
      assert (not ierr);
      //delta.AT(m) /= Carpet::maxspacereflevelfact;
      gh const * const hh = Carpet::vhh.AT(m);
      ibbox const & baseext = hh->baseextent(Carpet::mglevel, 0);
      delta.AT(m) /= baseext.stride();
      idelta.AT(m) = CCTK_REAL(1.0) / delta.AT(m);
      if (veryverbose) {
        cout << "GetCoordRange[" << m << "]: lower=" << lower.AT(m) << " upper=" << upper.AT(m) << " delta=" << delta.AT(m) << endl;
      }
    }
    
    // Calculate refinement levels, components, and integer grid point
    // indices
    if (verbose) CCTK_INFO ("Mapping points onto components");
    vector<fasterp_iloc_t> ilocs (npoints);
    vector<int> proc (npoints);
    fill_with_poison (proc);
    vector<int> nlocs (nprocs, 0);
    int min_rl, max_rl;
    if (Carpet::is_level_mode()) {
      min_rl = Carpet::reflevel;
      max_rl = Carpet::reflevel + 1;
    } else if (Carpet::is_global_mode()) {
      min_rl = 0;
      max_rl = Carpet::reflevels;
    }
#pragma omp parallel for
    for (int n=0; n<npoints; ++n) {
      int const m = locations.maps.AT(n);
      rvect const pos (locations.coords[0].AT(n),
                       locations.coords[1].AT(n),
                       locations.coords[2].AT(n));
      
      gh const * const hh = Carpet::vhh.AT(m);
      // ibbox const & baseext = hh->baseextent(Carpet::mglevel, 0);
      // rvect const rpos =
      //   (pos - lower.AT(m)) / (upper.AT(m) - lower.AT(m)) *
      //   rvect (baseext.upper() - baseext.lower());
      rvect const rpos = (pos - lower.AT(m)) * idelta.AT(m);
      
      // Find refinement level and component
      int rl, c;
      ivect ipos;
      hh->locate_position (rpos, Carpet::mglevel, min_rl, max_rl,
                           rl, c, ipos);
      if (not (rl>=0 and c>=0)) {
#pragma omp critical
        {
          ostringstream msg;
          msg << "Interpolation point " << n << " on map " << m << " "
              << "at " << pos << " is outside of the grid hierarchy";
          msg << "\n"
              << "rl=" << rl << " c=" << c << "\n"
              << "rpos=" << rpos << "\n"
              << "ipos=" << ipos << "\n"
              << "lower=" << lower << "\n"
              << "upper=" << upper << "\n"
              << "delta=" << delta << "\n"
              << "idelta=" << idelta << "\n"
              << "hh=" << *hh << "\n";
          CCTK_WARN (CCTK_WARN_ABORT, msg.str().c_str());
        }
      }
      assert (rl>=0 and c>=0);
      
      ibbox const & ext =
        Carpet::vdd.AT(m)
        ->light_boxes.AT(Carpet::mglevel).AT(rl).AT(c).exterior;
      rvect dpos = rpos - rvect(ipos);
      
      // Convert from Carpet indexing to grid point indexing
      assert (all (ipos % ext.stride() == ivect(0)));
      ipos /= ext.stride();
      dpos /= rvect(ext.stride());
      if (not (all (abs(dpos) <= rvect(0.5)))) {
        cout << "fasterp.cc:659\n"
             << "   dpos=" << dpos << "\n"
             << "   ext=" << ext << "\n";
      }
      assert (all (abs(dpos) <= rvect(0.5)));
      
      ivect const ind = ipos - ext.lower() / ext.stride();
      ivect const lsh = ext.shape() / ext.stride();
      int const ind3d = index(lsh, ind);
#if 0
      ENTER_SINGLEMAP_MODE (cctkGH, m, CCTK_GF) {
        ENTER_LOCAL_MODE (cctkGH, c, CCTK_GF) {
          CCTK_REAL const * restrict const xptr = (CCTK_REAL const *) CCTK_VarDataPtr (cctkGH, 0, "grid::x");
          CCTK_REAL const * restrict const yptr = (CCTK_REAL const *) CCTK_VarDataPtr (cctkGH, 0, "grid::y");
          CCTK_REAL const * restrict const zptr = (CCTK_REAL const *) CCTK_VarDataPtr (cctkGH, 0, "grid::z");
          assert (xptr);
          assert (yptr);
          assert (zptr);
          cout << "CI2 map=" << m << " pos=" << pos << " ind=" << ind << " x=" << xptr[ind3d] << " y=" << yptr[ind3d] << " z=" << zptr[ind3d] << endl;
        } LEAVE_LOCAL_MODE;
      } LEAVE_SINGLEMAP_MODE;
#endif
      
      // TODO: assert that there are enough ghost zones
      
      // TODO: store for every face/direction/component/map/reflevel
      // how wide the boundaries are in every direction. Then check
      // against this when the stencils are set up.
      
      // Store result
      fasterp_iloc_t & iloc = ilocs.AT(n);
      iloc.mrc    = mrc_t(m, rl, c);
#ifdef CARPETINTERP2_CHECK
      iloc.pn.p   = dist::rank();
      iloc.pn.n   = n;
      iloc.ipos   = ipos * ext.stride();
      iloc.ind    = ind;
#endif
      iloc.ind3d  = ind3d;
      iloc.offset = dpos;
      
      // Find source processor
      int const p = Carpet::vhh.AT(m)->processor(rl,c);
      proc.AT(n) = p;
      int& nlocs_p = nlocs.AT(p);
#pragma omp atomic
      ++ nlocs_p;
      
      // Output
      if (veryverbose) {
#pragma omp critical
        {
          cout << "Point #" << n << " at " << pos << ": iloc " << iloc << endl;
        }
      }
    }
    
    
    
    // Find mapping from processors to "processor indices": It may be
    // too expensive to store data for all processors, so we store
    // only data about those processors with which we actually
    // communicate.
    if (verbose) CCTK_INFO ("Determine set of communicating processors");
    {
      int n_nz_nlocs = 0;
      for (int p=0; p<nprocs; ++p) {
        if (nlocs.AT(p) > 0) {
          ++ n_nz_nlocs;
        }
      }
      recv_descr.procs.resize (n_nz_nlocs);
      recv_descr.procinds.resize (nprocs, -1);
      int pp = 0;
      int offset = 0;
      for (int p=0; p<nprocs; ++p) {
        if (nlocs.AT(p) > 0) {
          recv_descr.procs.AT(pp).p       = p;
          recv_descr.procs.AT(pp).offset  = offset;
          recv_descr.procs.AT(pp).npoints = nlocs.AT(p);
          recv_descr.procinds.AT(p) = pp;
          ++ pp;
          offset += nlocs.AT(p);
        }
      }
      assert (pp == n_nz_nlocs);
      assert (offset == npoints);
      recv_descr.npoints = npoints;
    }
    
    // Create a mapping "index" from location index, as specified by
    // the user, to the index order in which the data are received
    // from the other processors.
    if (verbose) {
      CCTK_INFO ("Determine inter-processor gather index permutation");
    }
    {
      vector<int> index (recv_descr.procs.size());
      fill_with_poison (index);
      for (int pp=0; pp<int(recv_descr.procs.size()); ++pp) {
        index.AT(pp) = recv_descr.procs.AT(pp).offset;
      }
      recv_descr.index.resize (npoints);
      // TODO: parallelise with OpenMP
      for (int n=0; n<npoints; ++n) {
        int const p = proc.AT(n);
        int const pp = recv_descr.procinds.AT(p);
        assert (pp >= 0);
        recv_descr.index.AT(n) = index.AT(pp);
        ++index.AT(pp);
      }
      for (int pp=0; pp<int(recv_descr.procs.size()); ++pp) {
        int const recv_npoints = index.AT(pp) - recv_descr.procs.AT(pp).offset;
        assert (recv_npoints == recv_descr.procs.AT(pp).npoints);
      }
#ifndef NDEBUG
      vector<bool> received (npoints);
      fill (received, false);
      // NOTE: Can't use OMP parallel here -- vector<bool> has the
      // wrong memory layout for this
      for (int n=0; n<npoints; ++n) {
        assert (not received.AT(n));
        received.AT(n) = true;
      }
      for (int n=0; n<npoints; ++n) {
        assert (received.AT(n));
      }
#endif
    }
    
    
    
    // Count the number of points which have to be sent to other
    // processors, and exchange this information with MPI
    if (verbose) {
      CCTK_INFO ("Count and exchange number of communicated grid points");
    }
    vector<int> recv_npoints (nprocs), recv_offsets (nprocs);
    fill (recv_npoints, 0);
    fill_with_poison (recv_offsets);
    {
      int offset = 0;
      for (int pp=0; pp<int(recv_descr.procs.size()); ++pp) {
        int const p = recv_descr.procs.AT(pp).p;
        recv_npoints.AT(p) = recv_descr.procs.AT(pp).npoints;
        recv_offsets.AT(p) = offset;
        offset += recv_npoints.AT(p);
      }
      assert (offset == npoints);
    }
    vector<int> send_npoints (nprocs), send_offsets (nprocs);
    fill_with_poison (send_npoints);
    fill_with_poison (send_offsets);
    MPI_Alltoall (&recv_npoints.front(), 1, MPI_INT,
                  &send_npoints.front(), 1, MPI_INT,
                  comm_world);
    int npoints_send = 0;
    for (int p=0; p<nprocs; ++p) {
      send_offsets.AT(p) = npoints_send;
      npoints_send += send_npoints.AT(p);
    }
    
    // Scatter the location information into a send-array, and
    // exchange it with MPI
    vector<fasterp_iloc_t> scattered_ilocs(npoints);
#pragma omp parallel for
    for (int n=0; n<npoints; ++n) {
      scattered_ilocs.AT(recv_descr.index.AT(n)) = ilocs.AT(n);
    }
    vector<fasterp_iloc_t> gathered_ilocs(npoints_send);
    MPI_Alltoallv
      (&scattered_ilocs.front(), &recv_npoints.front(), &recv_offsets.front(),
       fasterp_iloc_t::mpi_datatype(),
       &gathered_ilocs.front(),  &send_npoints.front(), &send_offsets.front(),
       fasterp_iloc_t::mpi_datatype(),
       comm_world);
    
#ifdef CARPETINTERP2_CHECK
    // Ensure that the ilocs were sent to the right processor
    for (int p=0; p<nprocs; ++p) {
#pragma omp parallel for
      for (int n=0; n<send_npoints.at(p); ++n) {
        fasterp_iloc_t const & iloc = gathered_ilocs.at(send_offsets.at(p) + n);
        mrc_t const & mrc = iloc.mrc;
        int const m  = mrc.m;
        int const rl = mrc.rl;
        int const c  = mrc.c;
        assert (Carpet::vhh.AT(m)->is_local(rl,c));
      }
    }
#endif
    
    
    
    // Fill in send descriptors
    send_descr.npoints = npoints_send;
    
    {
      int n_nz_send_npoints = 0;
      for (int p=0; p<nprocs; ++p) {
        if (send_npoints.AT(p) > 0) ++n_nz_send_npoints;
      }
      send_descr.procs.resize (n_nz_send_npoints);
      int pp = 0;
      for (int p=0; p<nprocs; ++p) {
        if (send_npoints.AT(p) > 0) {
          send_proc_t & send_proc = send_descr.procs.AT(pp);
          send_proc.p       = p;
          send_proc.offset  = send_offsets.AT(p);
          send_proc.npoints = send_npoints.AT(p);
          ++pp;
        }
      }
      assert (pp == n_nz_send_npoints);
    }
    
    
    
    // Calculate stencil coefficients
    if (verbose) CCTK_INFO ("Calculate stencil coefficients");
    
    for (size_t pp=0; pp<send_descr.procs.size(); ++pp) {
      send_proc_t & send_proc = send_descr.procs.AT(pp);
      
      vector<int> mrc2comp (maxmrc, -1);
      vector<int> comp2mrc (maxmrc);
      fill_with_poison (comp2mrc);
      int ncomps = 0;
      
      vector<int> npoints_comp (maxmrc);
      fill (npoints_comp, 0);
      
      // TODO: parallelise with OpenMP
      for (int n=0; n<send_proc.npoints; ++n) {
        fasterp_iloc_t const & iloc = gathered_ilocs.AT(send_proc.offset + n);
        int const mrc = iloc.mrc.get_ind();
        if (mrc2comp.AT(mrc) == -1) {
          mrc2comp.AT(mrc) = ncomps;
          comp2mrc.AT(ncomps) = mrc;
          ++ ncomps;
        }
        ++ npoints_comp.AT(mrc);
      }
      assert (ncomps <= maxmrc);
      send_proc.comps.resize(ncomps);
      
      int offset = 0;
      for (int comp=0; comp<ncomps; ++comp) {
        send_comp_t & send_comp = send_proc.comps.AT(comp);
        int const mrc = comp2mrc.AT(comp);
        send_comp.mrc = mrc;
        send_comp.locs.reserve (npoints_comp.AT(mrc));
        
        mrc_t const themrc (mrc);
        int const m  = themrc.m;
        int const rl = themrc.rl;
        int const c  = themrc.c;
        assert (Carpet::vhh.AT(m)->is_local(rl,c));
        ibbox const & ext =
          Carpet::vdd.AT(m)->light_boxes.AT(Carpet::mglevel).AT(rl).AT(c).exterior;
        send_comp.lsh = ext.shape() / ext.stride();
        
        send_comp.offset = offset;
        send_comp.npoints = npoints_comp.AT(mrc);
        offset += send_comp.npoints;
      }
      assert (offset == send_proc.npoints);
      
      send_proc.index.resize (send_proc.npoints);
      fill_with_poison (send_proc.index);
      // TODO: This is not parallel!  Loop over comps instead?
      // #pragma omp parallel for
      for (int n=0; n<send_proc.npoints; ++n) {
        fasterp_iloc_t const & iloc = gathered_ilocs.AT(send_proc.offset + n);
        int const mrc = iloc.mrc.get_ind();
        int const comp = mrc2comp.AT(mrc);
        send_comp_t & send_comp = send_proc.comps.AT(comp);
        
        send_proc.index.AT(n) = send_comp.offset + send_comp.locs.size();
        
        fasterp_src_loc_t sloc;
        int const ierr = sloc.calc_stencil (iloc, send_comp.lsh, order);
        if (ierr) {
          CCTK_VWarn (CCTK_WARN_ABORT, __LINE__, __FILE__, CCTK_THORNSTRING,
                      "Could not determine valid interpolation stencil for point %d on map %d, refinement level %d, component %d",
                      n, iloc.mrc.m, iloc.mrc.rl, iloc.mrc.c);
        }
        send_comp.locs.push_back (sloc);
      }
      
      for (int comp=0; comp<ncomps; ++comp) {
        send_comp_t & send_comp = send_proc.comps.AT(comp);
        assert (int(send_comp.locs.size()) == send_comp.npoints);
      }
      
#ifndef NDEBUG
      {
        vector<bool> used(send_proc.npoints, false);
        for (int n=0; n<send_proc.npoints; ++n) {
          assert (not used.AT(send_proc.index.AT(n)));
          used.AT(send_proc.index.AT(n)) = true;
        }
        for (int n=0; n<send_proc.npoints; ++n) {
          assert (used.AT(send_proc.index.AT(n)));
        }
      }
#endif
      
#ifdef CARPETINTERP2_CHECK
      for (int comp=0; comp<ncomps; ++comp) {
        send_comp_t const & send_comp = send_proc.comps.at(comp);
        assert (int(send_comp.locs.size()) == send_comp.npoints);
        for (int n=0; n<send_comp.npoints; ++n) {
          fasterp_src_loc_t const & sloc = send_comp.locs.at(n);
          assert (sloc.mrc == send_comp.mrc);
          assert (all (sloc.saved_lsh == send_comp.lsh));
        }
      }
#endif
      
    } // for pp
    
#ifdef CARPETINTERP2_CHECK
    {
      if (verbose) CCTK_INFO ("Compare send and receive counts");
      
      vector<int> recv_count (dist::size());
      fill (recv_count, 0);
      for (size_t pp=0; pp<recv_descr.procs.size(); ++pp) {
        recv_proc_t const & recv_proc = recv_descr.procs.AT(pp);
        assert (recv_count.AT(recv_proc.p) == 0);
        assert (recv_proc.npoints > 0);
        recv_count.AT(recv_proc.p) = recv_proc.npoints;
      }
      
      vector<int> send_count (dist::size());
      fill (send_count, 0);
      for (size_t pp=0; pp<send_descr.procs.size(); ++pp) {
        send_proc_t const & send_proc = send_descr.procs.AT(pp);
        assert (send_count.AT(send_proc.p) == 0);
        assert (send_proc.npoints > 0);
        send_count.AT(send_proc.p) = send_proc.npoints;
      }
      
      {
        vector<int> tmp_count (dist::size());
        MPI_Alltoall (&send_count.front(), 1, MPI_INT,
                      &tmp_count .front(), 1, MPI_INT,
                      dist::comm());
        swap (send_count, tmp_count);
      }
      bool error = false;
      for (int p=0; p<dist::size(); ++p) {
        if (not (send_count.AT(p) == recv_count.AT(p))) {
          ostringstream buf;
          buf << "Error: processor " << p << " "
              << "sends me " << send_count.AT(p) << " points, "
              << "while I receive " << recv_count.AT(p) << " from it";
          CCTK_WARN (CCTK_WARN_ALERT, buf.str().c_str());
          error = true;
        }
      }
      if (error) {
        CCTK_WARN (CCTK_WARN_ABORT, "Internal error in interpolation setup");
      }
    }
#endif
    
    if (verbose) CCTK_INFO ("Done.");
  }
  
  
  
  // Free the setup for one interpolation
  fasterp_setup_t::
  ~fasterp_setup_t ()
  {
    // do nothing -- C++ calls destructors automatically
  }
  
  
  
  // Interpolate
  void 
  fasterp_setup_t::
  interpolate (cGH const * restrict const cctkGH,
               vector<int> const & varinds,
               vector<CCTK_REAL *> & values)
    const
  {
    DECLARE_CCTK_PARAMETERS;
    
    // Check regridding epoch
    if (outofdate())
    {
      if (reflevel == -1) {
        CCTK_VWarn (CCTK_WARN_ABORT, __LINE__, __FILE__, CCTK_THORNSTRING,
                    "The Carpet grid structure was changed since this fasterp_setup was created");
      } else {
        CCTK_VWarn (CCTK_WARN_ABORT, __LINE__, __FILE__, CCTK_THORNSTRING,
                    "The Carpet grid structure on level %d was changed since this fasterp_setup was created",
                    reflevel);
      }
    }
    
    // Desired time level
    int const tl = 0;           // current time
    
    size_t const nvars = varinds.size();
    if (verbose) CCTK_VInfo (CCTK_THORNSTRING,
                             "Interpolating %d variables", int(nvars));
    assert (values.size() == nvars);
    
    if (nvars == 0) return;
    
    for (size_t v=0; v<values.size(); ++v) {
      int const vi = varinds.AT(v);
      assert (vi >= 0 and vi <= CCTK_NumVars());
      int const gi = CCTK_GroupIndexFromVarI (vi);
      assert (gi >= 0);
      cGroup group;
      check (not CCTK_GroupData (gi, &group));
      assert (group.grouptype == CCTK_GF);
      assert (group.vartype == CCTK_VARIABLE_REAL);
      assert (group.dim == dim);
#if 0
      // Cannot be called in level mode
      cGroupDynamicData dyndata;
      {
        int const ierr = CCTK_GroupDynamicData (cctkGH, gi, &dyndata);
        assert (not ierr);
      }
      assert (dyndata.activetimelevels > tl);
#endif
      assert (recv_descr.npoints == 0 or values.AT(v) != NULL);
    }
    
    MPI_Comm & comm_world = * (MPI_Comm *) GetMPICommWorld (cctkGH);
    int const mpi_tag = 0;
    
    // Post Irecvs
    if (verbose) CCTK_INFO ("Posting MPI_Irecvs");
    vector<CCTK_REAL> recv_points (recv_descr.npoints * nvars);
    fill_with_poison (recv_points);
    vector<MPI_Request> recv_reqs (recv_descr.procs.size());
#ifdef CARPETINTERP2_CHECK
    vector<pn_t> recv_pn (recv_descr.npoints);
    vector<MPI_Request> recv_reqs_pn (recv_descr.procs.size());
#endif
    for (size_t pp=0; pp<recv_descr.procs.size(); ++pp) {
      recv_proc_t const & recv_proc = recv_descr.procs.AT(pp);
      
      MPI_Irecv (& recv_points.AT(recv_proc.offset * nvars),
                 recv_proc.npoints * nvars,
                 dist::mpi_datatype<CCTK_REAL>(), recv_proc.p, mpi_tag,
                 comm_world, & recv_reqs.AT(pp));
#ifdef CARPETINTERP2_CHECK
      MPI_Irecv (& recv_pn.AT(recv_proc.offset),
                 recv_proc.npoints * 2,
                 dist::mpi_datatype<int>(), recv_proc.p, mpi_tag,
                 comm_world, & recv_reqs_pn.AT(pp));
#endif
    }
    
    // Interpolate data and post Isends
    if (verbose) CCTK_INFO ("Interpolating and posting MPI_Isends");
    // TODO: Use one array per processor?
    vector<CCTK_REAL> send_points (send_descr.npoints * nvars);
    fill_with_poison (send_points);
    vector<MPI_Request> send_reqs (send_descr.procs.size());
#ifdef CARPETINTERP2_CHECK
    vector<pn_t> send_pn (send_descr.npoints);
    vector<MPI_Request> send_reqs_pn (send_descr.procs.size());
#endif
    for (size_t pp=0; pp<send_descr.procs.size(); ++pp) {
      send_proc_t const & send_proc = send_descr.procs.AT(pp);
      
      vector<CCTK_REAL> computed_points (send_proc.npoints * nvars);
      fill_with_poison (computed_points);
#ifdef CARPETINTERP2_CHECK
      vector<pn_t> computed_pn (send_descr.npoints);
#endif
      for (size_t comp=0; comp<send_proc.comps.size(); ++comp) {
        send_comp_t const & send_comp = send_proc.comps.AT(comp);
        
        int const m  = send_comp.mrc.m;
        int const rl = send_comp.mrc.rl;
        int const c  = send_comp.mrc.c;
        
        int const lc = Carpet::vhh.AT(m)->get_local_component(rl,c);
        
        // Get pointers to 3D data
        vector<CCTK_REAL const *> varptrs (nvars);
        for (size_t v=0; v<nvars; ++v) {
          int const gi = CCTK_GroupIndexFromVarI (varinds.AT(v));
          assert (gi >= 0);
          int const vi = varinds.AT(v) - CCTK_FirstVarIndexI (gi);
          assert (vi >= 0);
          ggf const *const ff = Carpet::arrdata.AT(gi).AT(m).data.AT(vi);
          void const *const ptr =
            ff->data_pointer(tl, rl, lc, Carpet::mglevel)->storage();
          varptrs.AT(v) = (CCTK_REAL const *)ptr;
          assert (varptrs.AT(v));
        }
        
        // TODO: This loops seems unbalanced.  Maybe the different
        // interpolations have different costs.
#pragma omp parallel for schedule (dynamic, 1000)
        for (int n=0; n<int(send_comp.locs.size()); ++n) {
          size_t const ind = (send_comp.offset + n) * nvars;
          send_comp.locs.AT(n).interpolate
            (send_comp.lsh, order, varptrs, &computed_points.AT(ind));
#ifdef CARPETINTERP2_CHECK
          computed_pn.AT(send_comp.offset + n) = send_comp.locs.AT(n).pn;
#endif
        }
        
      } // for comp
      
      // Gather send points
#pragma omp parallel for
      for (int n=0; n<send_proc.npoints; ++n) {
        size_t const nn = send_proc.index.AT(n);
        for (size_t v=0; v<nvars; ++v) {
          send_points.AT((send_proc.offset + n) * nvars + v) =
            computed_points.AT(nn * nvars + v);
        }
#ifdef CARPETINTERP2_CHECK
        send_pn.AT(send_proc.offset + n) = computed_pn.AT(nn);
#endif
      }
      // TODO: Use a scatter index list instead of a gather index
      // list, and combine the scattering with the calculation
      
#ifdef CARPETINTERP2_CHECK
#pragma omp parallel for
      for (int n=0; n<send_proc.npoints; ++n) {
        assert (send_pn.AT(send_proc.offset + n).p == send_proc.p);
        if (n>0) {
          assert (send_pn.AT(send_proc.offset + n  ).n >
                  send_pn.AT(send_proc.offset + n-1).n);
        }
      }
#endif
      
      MPI_Isend (& send_points.AT(send_proc.offset * nvars),
                 send_proc.npoints * nvars,
                 dist::mpi_datatype<CCTK_REAL>(), send_proc.p, mpi_tag,
                 comm_world, & send_reqs.AT(pp));
#ifdef CARPETINTERP2_CHECK
      MPI_Isend (& send_pn.AT(send_proc.offset),
                 send_proc.npoints * 2,
                 dist::mpi_datatype<int>(), send_proc.p, mpi_tag,
                 comm_world, & send_reqs_pn.AT(pp));
#endif
    } // for pp
    
    // Wait for Irecvs to complete
    if (verbose) CCTK_INFO ("Waiting for MPI_Irevcs to complete");
    MPI_Waitall (recv_reqs.size(), & recv_reqs.front(), MPI_STATUSES_IGNORE);
#ifdef CARPETINTERP2_CHECK
    MPI_Waitall (recv_reqs.size(), & recv_reqs_pn.front(), MPI_STATUSES_IGNORE);
#endif
    
    // Gather data
    if (verbose) CCTK_INFO ("Gathering data");
#pragma omp parallel for
    for (int n=0; n<recv_descr.npoints; ++n) {
      size_t const nn = recv_descr.index.AT(n);
      for (size_t v=0; v<nvars; ++v) {
        values.AT(v)[n] = recv_points.AT(nn * nvars + v);
#ifndef NDEBUG
        recv_points.AT(nn * nvars + v) = poison;
#endif
      }
#ifdef CARPETINTERP2_CHECK
      assert (recv_pn.AT(nn).p == dist::rank());
      assert (recv_pn.AT(nn).n == n);
#endif
    }
    
    // Wait for Isends to complete
    if (verbose) CCTK_INFO ("Waiting for MPI_Isends to complete");
    MPI_Waitall (send_reqs.size(), & send_reqs.front(), MPI_STATUSES_IGNORE);
#ifdef CARPETINTERP2_CHECK
    MPI_Waitall (send_reqs.size(), & send_reqs_pn.front(), MPI_STATUSES_IGNORE);
#endif
    
    if (verbose) CCTK_INFO ("Done.");
  }
  
  
  
} // namespace CarpetInterp2