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
1277
1278
1279
1280
1281
1282
1283
1284
|
/*@@
@file template.c
@date 23 Oct 2001
@author Jonathan Thornburg <jthorn@aei.mpg.de>
@desc
This file is a template for generalized interpolation for
1-, 2-, or 3-d molecules. It's used by defining various
C preprocessor macros, then #including this file. Each
such inclusion defines a single interpolation function,
which does generalized interpolation for a single
combination of (N_dims, molecule_family, order, smoothing).
The following header files must be #included before
#including this file:
<math.h>
<limits.h>
<stdlib.h>
<string.h>
<stdio.h>
"util_ErrorCodes.h"
"cctk.h"
"InterpLocalUniform.h"
The following preprocessor macros must be defined before
#including this file (note the macros which are file names
must all include the proper quotes for a #include in their
expansion, e.g.
#define DATA_VAR_DCL_FILE_NAME \
"coeffs/2d.cube.size3/data-var.dcl.c"
@enddesc
@var FUNCTION_NAME
@vdesc The name of the interpolation function, e.g.
#define FUNCTION_NAME LocalInterp_ILA_2d_cube_ord2_sm0
@endvar
@var N_DIMS
@vdesc The number of dimensions in which to interpolate, e.g.
#define N_DIMS 2
The present implementation restricts this to 1, 2,
or 3, but this could easily be changed if needed.
Note that MAX_N_DIMS (defined in "InterpLocalUniform.c")
is a compile-time upper bound for N_DIMS, useful for
sizing arrays etc.
@endvar
@var MOLECULE_MIN_M
@vdesc The minimum m coordinate in the molecule, e.g.
#define MOLECULE_MIN_M -1
The present implementation takes this to be the same in each
dimension, but this could easily be changed if needed.
@endvar
@var MOLECULE_MAX_M
@vdesc The maximum m coordinate in the molecule, e.g.
#define MOLECULE_MAX_M 1
The present implementation takes this to be the same in each
dimension, but this could easily be changed if needed.
@endvar
@var MOLECULE_SIZE
@vdesc The diameter of (number of points in) the molecules to be used,
e.g.
#define MOLECULE_SIZE 3
The present implementation takes this to be the same in each
dimension, but this could easily be changed if needed.
@endvar
@var HAVE_OP_{I,DX,DY,DXX,DXY,DYY,...}
Each of these symbols should be defined or not defined
according as if the corresponding derivative operator is
to be supported or not supported by this function.
@endvar
@var DATA_VAR_DCL_FILE_NAME
@vdesc The name of a file (presumably machine-generated) containing
a sequence of one or more C declarations for a molecule-sized
set of "data variables", eg (for 2D size-3 molecules)
fp data_m1_m1;
fp data_0_m1;
fp data_p1_m1;
fp data_m1_0;
fp data_0_0;
fp data_p1_0;
fp data_m1_p1;
fp data_0_p1;
fp data_p1_p1;
@endvar
@var COEFF_{I,DX,DY,DXX,DXY,DYY,...}_DCL_FILE_NAME
@vdesc Each of these macros should be the name of a file
(presumably machine-generated) containing a sequence
of one or more C declarations for a molecule-sized set
of coefficients for the corresponding derivative operator,
eg (for dx with 2D size-3 molecules)
fp coeff_dx_m1_m1;
fp coeff_dx_0_m1;
fp coeff_dx_p1_m1;
fp coeff_dx_m1_0;
fp coeff_dx_0_0;
fp coeff_dx_p1_0;
fp coeff_dx_m1_p1;
fp coeff_dx_0_p1;
fp coeff_dx_p1_p1;
@endvar
@var COEFF_{I,DX,DY,DXX,DXY,DYY,...}_VAR_STORE_FILE_NAME
@vdesc The name of a file (presumably machine-generated) containing
a sequence of C assignment statements to assign the interpolation
coefficients to the corresponding COEFF(...) expressions,
eg (for 2D size-3 molecules)
COEFF(-1,-1) = coeff_dx_m1_m1;
COEFF(0,-1) = coeff_dx_0_m1;
COEFF(1,-1) = coeff_dx_p1_m1;
COEFF(-1,0) = coeff_dx_m1_0;
COEFF(0,0) = coeff_dx_0_0;
COEFF(1,0) = coeff_dx_p1_0;
COEFF(-1,1) = coeff_dx_m1_p1;
COEFF(0,1) = coeff_dx_0_p1;
COEFF(1,1) = coeff_dx_p1_p1;
@endvar
@var DATA_VAR_ASSIGN_FILE_NAME
@vdesc The name of a file (presumably machine-generated) containing
a sequence of C assignment statements to assign DATA(...) to
the corresponding data variables for each molecule component,
eg (for 2D size-3 molecules)
data_m1_m1 = DATA(-1,-1);
data_0_m1 = DATA(0,-1);
data_p1_m1 = DATA(1,-1);
data_m1_0 = DATA(-1,0);
data_0_0 = DATA(0,0);
data_p1_0 = DATA(1,0);
data_m1_p1 = DATA(-1,+1);
data_0_p1 = DATA(0,+1);
data_p1_p1 = DATA(1,+1);
@endvar
@var INTERP_{I,DX,DY,DXX,DXY,DYY,...}_COMPUTE_FILE_NAME
@vdesc Each of these macros should be the name of a file
(presumably machine-generated) containing a C assignment
statement (or a sequence of such statements) to compute
the variable result as the molecule-sized linear combination
of the data variables (cf DATA_VAR_{DCL,ASSIGN}_FILE_NAME)
and the interpolation coefficients
(cf COEFF_{I,DX,DY,DXX,DXY,DYY,...}_COMPUTE_FILE_NAME),
eg (for 2D size-3 molecules, dx operator)
result =
coeff_dx_m1_m1*data_m1_m1
+ coeff_dx_0_m1*data_0_m1
+ coeff_dx_p1_m1*data_p1_m1
+ coeff_dx_p2_m1*data_p2_m1
+ coeff_dx_m1_0*data_m1_0
+ coeff_dx_0_0*data_0_0
+ coeff_dx_p1_0*data_p1_0
+ coeff_dx_p2_0*data_p2_0
+ coeff_dx_m1_p1*data_m1_p1
+ coeff_dx_0_p1*data_0_p1
+ coeff_dx_p1_p1*data_p1_p1
+ coeff_dx_p2_p1*data_p2_p1
+ coeff_dx_m1_p2*data_m1_p2
+ coeff_dx_0_p2*data_0_p2
+ coeff_dx_p1_p2*data_p1_p2
+ coeff_dx_p2_p2*data_p2_p2;
Note that this may *NOT* include any variable declarations;
it should be valid to appear in the middle of a sequence of
C statements.
@endvar
@var COEFF_{I,DX,DY,DXX,DXY,DYY,...}_COMPUTE_FILE_NAME
@vdesc Each of these macros should be the name of a file
(presumably machine-generated) containing a sequence of C
assignment statements to compute all the (molecule-sized
set of) interpolation coefficients for the specified
derivative operator, as polynomials in the variables
(x,y,z), eg (for 2D size-3 molecules, dx operator)
fp t35,
t42,
t41,
t40,
t39,
t36;
t35 = RATIONAL(1.0,3.0)*x;
t42 = t35+RATIONAL(-1.0,4.0)*y;
t41 = t35+RATIONAL(1.0,4.0)*y;
t40 = RATIONAL(-1.0,6.0);
t39 = RATIONAL(1.0,6.0);
t36 = RATIONAL(-2.0,3.0)*x;
coeff_dx_m1_m1 = t40+t41;
coeff_dx_0_m1 = t36;
coeff_dx_p1_m1 = t39+t42;
coeff_dx_m1_0 = t40+t35;
coeff_dx_0_0 = t36;
coeff_dx_p1_0 = t35+t39;
coeff_dx_m1_p1 = t40+t42;
coeff_dx_0_p1 = t36;
coeff_dx_p1_p1 = t39+t41;
As illustrated, the code may use the macro RATIONAL
(defined later in this file) to represent rational-number
coefficients, and it may also declare temporary variables
as needed.
Note that these are the only coefficients which depend
on the actual interpolation scheme used; all the others
just depend on the interpolation dimension and molecule
family and size.
@endvar
@version $Id$
@@*/
/******************************************************************************/
/*@@
@routine FUNCTION_NAME
@date 23 Oct 2001
@author Jonathan Thornburg <jthorn@aei.mpg.de>
@desc
This function does generalized interpolation of one or more
2d arrays to arbitrary points. For details, see the header
comments for InterpLocalUniform() (in "InterpLocalUniform.c"
in this same directory).
This function's arguments are mostly all a subset of those of
InterpLocalUniform() ; the difference is that this function
takes all its arguments explicitly, whereas InputLocalArrays()
takes some of them indirectly via a key/value parameter table.
Here we document only the "new" explicit arguments, and these
only briefly; see the InterpLocalUniform() documentation
and/or the thorn guide for this thorn for further details.
@var error_flags
@vdesc If we encounter an out-of-range point and this pointer
is non-NULL, we store a description of the out-of-range
point in the pointed-to structure.
@vtype struct error_flags *error_flags;
@vio out
@var molecule_structure_flags
@vdesc If this pointer is non-NULL, we store flags describing
the interpolation molecule's structure in the pointed-to
structure.
@vtype struct molecule_structure_flags *molecule_structure_flags;
@vio out
@var molecule_min_max_m_info
@vdesc If this pointer is non-NULL, we store the interpolation
molecule's min/max m in the pointed-to structure.
@vtype struct molecule_min_max_m_info *molecule_min_max_m_info;
@vio out
@var molecule_positions
@vdesc If this pointer is non-NULL, we store the interpolation
molecule's positions in the (caller-supplied) arrays
pointed to by this pointer.
@vtype CCTK_INT* const molecule_positions[];
@vio out
@returntype int
@returndesc This function's return result is the same as that of
InterpLocalUniform():
0 success
UTIL_ERROR_BAD_INPUT one of the input arguments is invalid
UTIL_ERROR_NO_MEMORY unable to malloc() temporary memory
CCTK_ERROR_INTERP_POINT_X_RANGE
interpolation point is out of range
@endreturndesc
@@*/
int FUNCTION_NAME(/***** coordinate system *****/
const CCTK_REAL coord_system_origin[],
const CCTK_REAL grid_spacing[],
/***** interpolation points *****/
int N_interp_points,
int interp_coords_type_code,
const void* const interp_coords[],
const CCTK_REAL out_of_range_tolerance[],
/***** input arrays *****/
int N_input_arrays,
const CCTK_INT input_array_offsets[],
const CCTK_INT input_array_strides[],
const CCTK_INT input_array_min_subscripts[],
const CCTK_INT input_array_max_subscripts[],
const CCTK_INT input_array_type_codes[],
const void* const input_arrays[],
/***** output arrays *****/
int N_output_arrays,
const CCTK_INT output_array_type_codes[],
void* const output_arrays[],
/***** operation info */
const CCTK_INT operand_indices[],
const CCTK_INT operation_codes[],
/***** other return results *****/
struct error_flags* error_flags,
struct molecule_structure_flags* molecule_structure_flags,
struct molecule_min_max_m_info* molecule_min_max_m_info,
CCTK_INT* const molecule_positions[],
struct Jacobian_info* Jacobian_info)
{
/*
* Implementation notes:
*
* The basic outline of this function is as follows:
*
* store molecule structure flags
* if (querying molecule min/max m)
* then store molecule min/max m
* compute "which derivatives are wanted" flags
* precompute 1/dx factors
* for each interpolation point
* {
* declare all the coefficients
* declare all the data-values variables
* ***fetch*** interpolation point coordinates
* compute position of interpolation molecules
* with respect to the grid
* if (querying molecule positions)
* then store this molecule position
* compute coefficients for all derivatives which are wanted
* for each output array
* {
* ***decode*** the input/output array datatypes
* to determine whether they're real or complex
* (they must both be the same in this regard),
* and define
* int N_parts = data is complex ? 2 : 1;
* int part;
* for (part = 0 ; part < N_parts ; ++part)
* {
* if (this output array is computed
* using a different input array
* than the previous one || part != 0)
* then ***fetch*** the input array values
* into local "data" variables
* {
* fp result;
* switch (operation_code)
* {
* case 0:
* result = compute the interpolant
* itself as a linear combination
* of data variables & op=0 coeffs
* break;
* case 1:
* result = compute the interpolant
* itself as a linear combination
* of data variables & op=1 coeffs
* result *= inverse_dx;
* break;
* case ...
* }
* if (querying Jacobian)
* then {
* fp factor;
* switch (operation_code)
* {
* case 0:
* store the op=0 Jacobian values
* break;
* case 1:
* factor = inverse_dx;
* store the op=1 Jacobian values
* (multiplying by factor )
* break;
* case ...
* }
* }
* ***store*** result in output array
* }
* }
* }
* }
*
* Here "***decode***", "***fetch***", and "***store***" are all
* actually switches on the various array datatypes. For complex
* datatypes "***fetch***" and "***store***" pointer-alias the
* N-dimensional input array to a 2-dimensional array where the 1st axis
* is the 1-D subscript corresponding to the N input axes, and the 2nd
* axes has 2 elements subscripted by [part] for the (real,imaginary)
* components of the complex values.
*
* At present we do all floating-point computations in type "fp"
* (typically a typedef for CCTK_REAL), so arrays of higher precision
* than this will incur extra rounding errors. In practice these should
* be negligible compared to the "truncation" interpolation errors.
*/
/*
* Naming conventions:
* in, out = 0-origin indices each selecting an input/output array
* pt = 0-origin index selecting an interpolation point
*/
/* number of real "parts" in a complex number */
#define COMPLEX_N_PARTS 2
/* layout of axes in N_dims-element arrays */
#define X_AXIS 0
#define Y_AXIS 1
#define Z_AXIS 2
#if (N_DIMS > 3)
#error "N_DIMS may not be > 3!"
#endif
/* basic sanity check on molecule size */
#define MOLECULE_M_COUNT (MOLECULE_MAX_M - MOLECULE_MIN_M + 1)
#if (MOLECULE_SIZE != MOLECULE_M_COUNT)
#error "MOLECULE_SIZE inconsistent with MOLECULE_{MIN,MAX}_M!"
#endif
/* input array size, strides, and subscripting computation */
#if (N_DIMS >= 1)
const int stride_i = input_array_strides[X_AXIS];
#endif
#if (N_DIMS >= 2)
const int stride_j = input_array_strides[Y_AXIS];
#endif
#if (N_DIMS >= 3)
const int stride_k = input_array_strides[Z_AXIS];
#endif
#if (N_DIMS == 1)
#define SUB1(i) (i*stride_i)
#elif (N_DIMS == 2)
#define SUB2(i,j) (i*stride_i + j*stride_j)
#elif (N_DIMS == 3)
#define SUB3(i,j,k) (i*stride_i + j*stride_j + k*stride_k)
#else
#error "N_DIMS must be 1, 2, or 3!"
#endif
/* macros used by machine-generated interpolation coefficient expressions */
/*
* FIXME: right now this is used as (eg) RATIONAL(1.0,2.0);
* it might be cleaner if it were RATIONAL(1,2) with the
* preprocessor ## operator used to glue on the .0
* (I _think_ this is portable, but is it really??)
*/
#define RATIONAL(num,den) (num/den)
/*
***** execution begins here *****
*/
/*
* Jacobian structure info
*/
if (molecule_structure_flags != NULL)
then {
molecule_structure_flags->MSS_is_fn_of_interp_coords = 0;
molecule_structure_flags->MSS_is_fn_of_which_operation = 0;
molecule_structure_flags->MSS_is_fn_of_input_array_values = 0;
molecule_structure_flags->Jacobian_is_fn_of_input_array_values = 0;
}
if (molecule_min_max_m_info != NULL)
then {
int axis;
for (axis = 0 ; axis < N_DIMS ; ++axis)
{
molecule_min_max_m_info->molecule_min_m[axis] = MOLECULE_MIN_M;
molecule_min_max_m_info->molecule_max_m[axis] = MOLECULE_MAX_M;
}
}
/*
* compute flags specifying which derivatives are wanted
*/
{
#ifdef HAVE_OP_I
bool want_I = false;
#endif
#ifdef HAVE_OP_DX
bool want_dx = false;
#endif
#ifdef HAVE_OP_DY
bool want_dy = false;
#endif
#ifdef HAVE_OP_DZ
bool want_dz = false;
#endif
#ifdef HAVE_OP_DXX
bool want_dxx = false;
#endif
#ifdef HAVE_OP_DXY
bool want_dxy = false;
#endif
#ifdef HAVE_OP_DXZ
bool want_dxz = false;
#endif
#ifdef HAVE_OP_DYY
bool want_dyy = false;
#endif
#ifdef HAVE_OP_DYZ
bool want_dyz = false;
#endif
#ifdef HAVE_OP_DZZ
bool want_dzz = false;
#endif
{
int out;
for (out = 0 ; out < N_output_arrays ; ++out)
{
switch (operation_codes[out])
{
#ifdef HAVE_OP_I
case 0: want_I = true; break;
#endif
#ifdef HAVE_OP_DX
case 1: want_dx = true; break;
#endif
#ifdef HAVE_OP_DY
case 2: want_dy = true; break;
#endif
#ifdef HAVE_OP_DZ
case 3: want_dz = true; break;
#endif
#ifdef HAVE_OP_DXX
case 11: want_dxx = true; break;
#endif
#ifdef HAVE_OP_DXY
case 12:
case 21: want_dxy = true; break;
#endif
#ifdef HAVE_OP_DXZ
case 13:
case 31: want_dxz = true; break;
#endif
#ifdef HAVE_OP_DYY
case 22: want_dyy = true; break;
#endif
#ifdef HAVE_OP_DYZ
case 23:
case 32: want_dyz = true; break;
#endif
#ifdef HAVE_OP_DZZ
case 33: want_dzz = true; break;
#endif
default:
CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
"Generalized interpolation operation_code %d not supported!",
operation_codes[out]); /*NOTREACHED*/
return UTIL_ERROR_BAD_INPUT; /*** ERROR RETURN ***/
}
}
}
/*
* save origin/delta variables, precompute 1/delta factors
* ... in theory we could compute only those factors we're going to use,
* but it's not worth the trouble, so we just compute them all
*/
{
#if N_DIMS >= 1
const fp origin_x = coord_system_origin[X_AXIS];
const fp delta_x = grid_spacing[X_AXIS];
#if defined(HAVE_OP_DX) || defined(HAVE_OP_DXY) || defined(HAVE_OP_DXZ)
const fp inverse_delta_x = 1.0 / delta_x;
#endif
#ifdef HAVE_OP_DXX
const fp inverse_delta_x2 = 1.0 / (delta_x*delta_x);
#endif
#endif
#if N_DIMS >= 2
const fp origin_y = coord_system_origin[Y_AXIS];
const fp delta_y = grid_spacing[Y_AXIS];
#if defined(HAVE_OP_DY) || defined(HAVE_OP_DXY) || defined(HAVE_OP_DYZ)
const fp inverse_delta_y = 1.0 / delta_y;
#endif
#ifdef HAVE_OP_DYY
const fp inverse_delta_y2 = 1.0 / (delta_y*delta_y);
#endif
#endif
#if N_DIMS >= 3
const fp origin_z = coord_system_origin[Z_AXIS];
const fp delta_z = grid_spacing[Z_AXIS];
#if defined(HAVE_OP_DZ) || defined(HAVE_OP_DXZ) || defined(HAVE_OP_DYZ)
const fp inverse_delta_z = 1.0 / delta_z;
#endif
#ifdef HAVE_OP_DZZ
const fp inverse_delta_z2 = 1.0 / (delta_z*delta_z);
#endif
#endif
/*
* interpolate at each point
*/
int pt;
for (pt = 0 ; pt < N_interp_points ; ++pt)
{
/* declare all the data-values variables */
#include DATA_VAR_DCL_FILE_NAME
/* declare all the interpolation coefficients */
#ifdef HAVE_OP_I
#include COEFF_I_DCL_FILE_NAME
#endif
#ifdef HAVE_OP_DX
#include COEFF_DX_DCL_FILE_NAME
#endif
#ifdef HAVE_OP_DY
#include COEFF_DY_DCL_FILE_NAME
#endif
#ifdef HAVE_OP_DZ
#include COEFF_DZ_DCL_FILE_NAME
#endif
#ifdef HAVE_OP_DXX
#include COEFF_DXX_DCL_FILE_NAME
#endif
#ifdef HAVE_OP_DXY
#include COEFF_DXY_DCL_FILE_NAME
#endif
#ifdef HAVE_OP_DXZ
#include COEFF_DXZ_DCL_FILE_NAME
#endif
#ifdef HAVE_OP_DYY
#include COEFF_DYY_DCL_FILE_NAME
#endif
#ifdef HAVE_OP_DYZ
#include COEFF_DYZ_DCL_FILE_NAME
#endif
#ifdef HAVE_OP_DZZ
#include COEFF_DZZ_DCL_FILE_NAME
#endif
/*
* ***fetch*** interpolation point coordinates
*/
fp interp_coords_fp[N_DIMS];
int axis;
for (axis = 0 ; axis < N_DIMS ; ++axis)
{
/* pointer to array of interp coords for this axis */
const void *const interp_coords_ptr = interp_coords[axis];
switch (interp_coords_type_code)
{
case CCTK_VARIABLE_REAL:
{
const CCTK_REAL *const interp_coords_ptr_real
= (const CCTK_REAL *) interp_coords_ptr;
interp_coords_fp[axis] = interp_coords_ptr_real[pt];
break;
}
#ifdef HAVE_CCTK_REAL4
case CCTK_VARIABLE_REAL4:
{
const CCTK_REAL4 *const interp_coords_ptr_real4
= (const CCTK_REAL4 *) interp_coords_ptr;
interp_coords_fp[axis] = interp_coords_ptr_real4[pt];
break;
}
#endif
#ifdef HAVE_CCTK_REAL8
case CCTK_VARIABLE_REAL8:
{
const CCTK_REAL8 *const interp_coords_ptr_real8
= (const CCTK_REAL8 *) interp_coords_ptr;
interp_coords_fp[axis] = interp_coords_ptr_real8[pt];
break;
}
#endif
#ifdef HAVE_CCTK_REAL16
case CCTK_VARIABLE_REAL16:
{
/* FIXME: maybe we should warn (once per cactus run) */
/* that we may be doing arithmetic in lower */
/* precision than the interp coords? */
const CCTK_REAL16 *const interp_coords_ptr_real16
= (const CCTK_REAL16 *) interp_coords_ptr;
interp_coords_fp[axis]
= interp_coords_ptr_real16[pt];
break;
}
#endif
default:
CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
"interp-coords datatype %d not supported!",
interp_coords_type_code);
/*NOTREACHED*/
return UTIL_ERROR_BAD_INPUT; /*** ERROR RETURN ***/
/* end of switch (interp_coords_type_code) */
}
/* end of for (axis = ...) loop */
}
/*
* compute position of interpolation molecules with respect to
* the grid, i.e. compute (x,y,z) coordinates of interpolation
* point relative to molecule center, in units of the grid spacing
*
* n.b. we need the final answers in variables with the magic
* names (x,y,z) (machine-generated code uses these names),
* but we use temp variables as intermediates for (likely)
* better performance: the temp variables have their addresses
* taken and so may not be register-allocated, whereas the
* final (x,y,z) are "clean" and thus more likely to be
* register-allocated
*/
{
#if (N_DIMS >= 1)
fp x_temp;
const int center_i
= LocalInterp_molecule_posn(origin_x, delta_x,
input_array_min_subscripts[X_AXIS],
input_array_max_subscripts[X_AXIS],
out_of_range_tolerance[X_AXIS],
MOLECULE_SIZE,
interp_coords_fp[X_AXIS],
&x_temp,
(int *) NULL, (int *) NULL);
const fp x = x_temp;
#endif
#if (N_DIMS >= 2)
fp y_temp;
const int center_j
= LocalInterp_molecule_posn(origin_y, delta_y,
input_array_min_subscripts[Y_AXIS],
input_array_max_subscripts[Y_AXIS],
out_of_range_tolerance[Y_AXIS],
MOLECULE_SIZE,
interp_coords_fp[Y_AXIS],
&y_temp,
(int *) NULL, (int *) NULL);
const fp y = y_temp;
#endif
#if (N_DIMS >= 3)
fp z_temp;
const int center_k
= LocalInterp_molecule_posn(origin_z, delta_z,
input_array_min_subscripts[Z_AXIS],
input_array_max_subscripts[Z_AXIS],
out_of_range_tolerance[Z_AXIS],
MOLECULE_SIZE,
interp_coords_fp[Z_AXIS],
&z_temp,
(int *) NULL, (int *) NULL);
const fp z = z_temp;
#endif
/* is the interpolation point out-of-range? */
#if (N_DIMS >= 1)
if ((center_i == INT_MIN) || (center_i == INT_MAX))
then {
if (error_flags != NULL)
then {
error_flags->error_pt = pt;
error_flags->error_axis = X_AXIS;
error_flags->error_end = (center_i > 0) ? +1 : -1;
}
return CCTK_ERROR_INTERP_POINT_X_RANGE; /*** ERROR RETURN ***/
}
#endif
#if (N_DIMS >= 2)
if ((center_j == INT_MIN) || (center_j == INT_MAX))
then {
if (error_flags != NULL)
then {
error_flags->error_pt = pt;
error_flags->error_axis = Y_AXIS;
error_flags->error_end = (center_j > 0) ? +1 : -1;
}
return CCTK_ERROR_INTERP_POINT_X_RANGE; /*** ERROR RETURN ***/
}
#endif
#if (N_DIMS >= 3)
if ((center_k == INT_MIN) || (center_k == INT_MAX))
then {
if (error_flags != NULL)
then {
error_flags->error_pt = pt;
error_flags->error_axis = Z_AXIS;
error_flags->error_end = (center_k > 0) ? +1 : -1;
}
return CCTK_ERROR_INTERP_POINT_X_RANGE; /*** ERROR RETURN ***/
}
#endif
/* 1D subscript of molecule center in input arrays */
/* (this is used in DATA data-fetching macros below) */
{
#if (N_DIMS == 1)
const int center_sub = SUB1(center_i);
#elif (N_DIMS == 2)
const int center_sub = SUB2(center_i, center_j);
#elif (N_DIMS == 3)
const int center_sub = SUB3(center_i, center_j, center_k);
#else
#error "N_DIMS must be 1, 2, or 3!"
#endif
/* molecule position queries */
if (molecule_positions != NULL)
then {
#if (N_DIMS >= 1)
molecule_positions[X_AXIS][pt] = center_i;
#endif
#if (N_DIMS >= 2)
molecule_positions[Y_AXIS][pt] = center_j;
#endif
#if (N_DIMS >= 3)
molecule_positions[Z_AXIS][pt] = center_k;
#endif
}
/*
* compute the coefficients for whichever derivatives are wanted
* using machine-generated coefficient expressions
* (these are polynomials in the variables (x,y,z))
*/
#ifdef HAVE_OP_I
if (want_I)
then {
#include COEFF_I_COMPUTE_FILE_NAME
}
#endif
#ifdef HAVE_OP_DX
if (want_dx)
then {
#include COEFF_DX_COMPUTE_FILE_NAME
}
#endif
#ifdef HAVE_OP_DY
if (want_dy)
then {
#include COEFF_DY_COMPUTE_FILE_NAME
}
#endif
#ifdef HAVE_OP_DZ
if (want_dz)
then {
#include COEFF_DZ_COMPUTE_FILE_NAME
}
#endif
#ifdef HAVE_OP_DXX
if (want_dxx)
then {
#include COEFF_DXX_COMPUTE_FILE_NAME
}
#endif
#ifdef HAVE_OP_DXY
if (want_dxy)
then {
#include COEFF_DXY_COMPUTE_FILE_NAME
}
#endif
#ifdef HAVE_OP_DXZ
if (want_dxz)
then {
#include COEFF_DXZ_COMPUTE_FILE_NAME
}
#endif
#ifdef HAVE_OP_DYY
if (want_dyy)
then {
#include COEFF_DYY_COMPUTE_FILE_NAME
}
#endif
#ifdef HAVE_OP_DYZ
if (want_dyz)
then {
#include COEFF_DYZ_COMPUTE_FILE_NAME
}
#endif
#ifdef HAVE_OP_DZZ
if (want_dzz)
then {
#include COEFF_DZZ_COMPUTE_FILE_NAME
}
#endif
/*
* compute each output array at this point
*/
{
int out;
const void *input_array_ptr = NULL;
for (out = 0 ; out < N_output_arrays ; ++out)
{
const int in = operand_indices[out];
const int input_offset = input_array_offsets[in];
const int sub_temp = input_offset + center_sub;
/*
* ***decode*** the input/output array datatypes
* to determine whether they're real or complex,
* and verify that they're both the same in this regard
* ==> define
* const int N_parts = data is complex ? 2 : 1;
*/
const int N_input_parts
= LocalInterp_decode_N_parts(input_array_type_codes[in]);
const int N_output_parts
= LocalInterp_decode_N_parts(output_array_type_codes[out]);
if (N_input_parts != N_output_parts)
then {
CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
"\n"
" can't do real input --> complex output\n"
" or complex input --> real output interpolation!\n"
" (0-origin) input #%d datatype = %d, output #%d datatype = %d\n"
,
in, (int) input_array_type_codes[in],
out, (int) output_array_type_codes[out]); /*NOTREACHED*/
return UTIL_ERROR_BAD_INPUT; /*** ERROR RETURN ***/
}
{
const int N_parts = N_input_parts;
int part;
for (part = 0 ; part < N_parts ; ++part)
{
if ( (input_arrays[in] != input_array_ptr)
|| (part != 0) )
then {
/*
* ***fetch*** the input array values
* into local variables
*/
input_array_ptr = input_arrays[in];
switch (input_array_type_codes[in])
{
case CCTK_VARIABLE_REAL:
{
const CCTK_REAL *const input_array_ptr_real
= (const CCTK_REAL *) input_array_ptr;
#undef DATA
#if (N_DIMS == 1)
#define DATA(i) input_array_ptr_real[sub_temp + SUB1(i)]
#elif (N_DIMS == 2)
#define DATA(i,j) input_array_ptr_real[sub_temp + SUB2(i,j)]
#elif (N_DIMS == 3)
#define DATA(i,j,k) input_array_ptr_real[sub_temp + SUB3(i,j,k)]
#else
#error "N_DIMS must be 1, 2, or 3!"
#endif
#include DATA_VAR_ASSIGN_FILE_NAME
break;
}
#ifdef HAVE_CCTK_REAL4
case CCTK_VARIABLE_REAL4:
{
const CCTK_REAL4 *const input_array_ptr_real4
= (const CCTK_REAL4 *) input_array_ptr;
#undef DATA
#if (N_DIMS == 1)
#define DATA(i) input_array_ptr_real4[sub_temp + SUB1(i)]
#elif (N_DIMS == 2)
#define DATA(i,j) input_array_ptr_real4[sub_temp + SUB2(i,j)]
#elif (N_DIMS == 3)
#define DATA(i,j,k) input_array_ptr_real4[sub_temp + SUB3(i,j,k)]
#else
#error "N_DIMS must be 1, 2, or 3!"
#endif
#include DATA_VAR_ASSIGN_FILE_NAME
break;
}
#endif /* HAVE_CCTK_REAL4 */
#ifdef HAVE_CCTK_REAL8
case CCTK_VARIABLE_REAL8:
{
const CCTK_REAL8 *const input_array_ptr_real8
= (const CCTK_REAL8 *) input_array_ptr;
#undef DATA
#if (N_DIMS == 1)
#define DATA(i) input_array_ptr_real8[sub_temp + SUB1(i)]
#elif (N_DIMS == 2)
#define DATA(i,j) input_array_ptr_real8[sub_temp + SUB2(i,j)]
#elif (N_DIMS == 3)
#define DATA(i,j,k) input_array_ptr_real8[sub_temp + SUB3(i,j,k)]
#else
#error "N_DIMS must be 1, 2, or 3!"
#endif
#include DATA_VAR_ASSIGN_FILE_NAME
break;
}
#endif /* HAVE_CCTK_REAL8 */
#ifdef HAVE_CCTK_REAL16
case CCTK_VARIABLE_REAL16:
{
/* FIXME: maybe we should warn (once per cactus run) that we may be */
/* doing arithmetic in lower precision than the data type? */
const CCTK_REAL16 *const input_array_ptr_real16
= (const CCTK_REAL16 *) input_array_ptr;
#undef DATA
#if (N_DIMS == 1)
#define DATA(i) input_array_ptr_real16[sub_temp + SUB1(i)]
#elif (N_DIMS == 2)
#define DATA(i,j) input_array_ptr_real16[sub_temp + SUB2(i,j)]
#elif (N_DIMS == 3)
#define DATA(i,j,k) input_array_ptr_real16[sub_temp + SUB3(i,j,k)]
#else
#error "N_DIMS must be 1, 2, or 3!"
#endif
#include DATA_VAR_ASSIGN_FILE_NAME
break;
}
#endif /* HAVE_CCTK_REAL16 */
case CCTK_VARIABLE_COMPLEX:
{
const CCTK_REAL (*const input_array_ptr_complex)[COMPLEX_N_PARTS]
= (const CCTK_REAL (*)[COMPLEX_N_PARTS]) input_array_ptr;
#undef DATA
#if (N_DIMS == 1)
#define DATA(i) input_array_ptr_complex[sub_temp + SUB1(i)][part]
#elif (N_DIMS == 2)
#define DATA(i,j) input_array_ptr_complex[sub_temp + SUB2(i,j)][part]
#elif (N_DIMS == 3)
#define DATA(i,j,k) input_array_ptr_complex[sub_temp + SUB3(i,j,k)][part]
#else
#error "N_DIMS must be 1, 2, or 3!"
#endif
#include DATA_VAR_ASSIGN_FILE_NAME
break;
}
#ifdef HAVE_CCTK_COMPLEX8
case CCTK_VARIABLE_COMPLEX8:
{
const CCTK_REAL4 (*const input_array_ptr_complex8)[COMPLEX_N_PARTS]
= (const CCTK_REAL4 (*)[COMPLEX_N_PARTS]) input_array_ptr;
#undef DATA
#if (N_DIMS == 1)
#define DATA(i) input_array_ptr_complex8[sub_temp + SUB1(i)][part]
#elif (N_DIMS == 2)
#define DATA(i,j) input_array_ptr_complex8[sub_temp + SUB2(i,j)][part]
#elif (N_DIMS == 3)
#define DATA(i,j,k) input_array_ptr_complex8[sub_temp + SUB3(i,j,k)][part]
#else
#error "N_DIMS must be 1, 2, or 3!"
#endif
#include DATA_VAR_ASSIGN_FILE_NAME
break;
}
#endif /* HAVE_CCTK_COMPLEX8 */
#ifdef HAVE_CCTK_COMPLEX16
case CCTK_VARIABLE_COMPLEX16:
{
const CCTK_REAL8 (*const input_array_ptr_complex16)[COMPLEX_N_PARTS]
= (const CCTK_REAL8 (*)[COMPLEX_N_PARTS]) input_array_ptr;
#undef DATA
#if (N_DIMS == 1)
#define DATA(i) input_array_ptr_complex16[sub_temp + SUB1(i)][part]
#elif (N_DIMS == 2)
#define DATA(i,j) input_array_ptr_complex16[sub_temp + SUB2(i,j)][part]
#elif (N_DIMS == 3)
#define DATA(i,j,k) input_array_ptr_complex16[sub_temp + SUB3(i,j,k)][part]
#else
#error "N_DIMS must be 1, 2, or 3!"
#endif
#include DATA_VAR_ASSIGN_FILE_NAME
break;
}
#endif /* HAVE_CCTK_COMPLEX16 */
#ifdef HAVE_CCTK_COMPLEX32
case CCTK_VARIABLE_COMPLEX32:
{
const CCTK_REAL16 (*const input_array_ptr_complex32)[COMPLEX_N_PARTS]
= (const CCTK_REAL16 (*)[COMPLEX_N_PARTS]) input_array_ptr;
#undef DATA
#if (N_DIMS == 1)
#define DATA(i) input_array_ptr_complex32[sub_temp + SUB1(i)][part]
#elif (N_DIMS == 2)
#define DATA(i,j) input_array_ptr_complex32[sub_temp + SUB2(i,j)][part]
#elif (N_DIMS == 3)
#define DATA(i,j,k) input_array_ptr_complex32[sub_temp + SUB3(i,j,k)][part]
#else
#error "N_DIMS must be 1, 2, or 3!"
#endif
#include DATA_VAR_ASSIGN_FILE_NAME
break;
}
#endif /* HAVE_CCTK_COMPLEX32 */
default:
CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
"input datatype %d not supported!",
input_array_type_codes[in]); /*NOTREACHED*/
return UTIL_ERROR_BAD_INPUT; /*** ERROR RETURN ***/
}
}
/*
* compute the interpolant itself
* as a linear combination of the data variables
*/
{
fp result;
switch (operation_codes[out])
{
#ifdef HAVE_OP_I
case 0:
#include INTERP_I_COMPUTE_FILE_NAME
break;
#endif
#ifdef HAVE_OP_DX
case 1:
#include INTERP_DX_COMPUTE_FILE_NAME
result *= inverse_delta_x;
break;
#endif
#ifdef HAVE_OP_DY
case 2:
#include INTERP_DY_COMPUTE_FILE_NAME
result *= inverse_delta_y;
break;
#endif
#ifdef HAVE_OP_DZ
case 3:
#include INTERP_DZ_COMPUTE_FILE_NAME
result *= inverse_delta_z;
break;
#endif
#ifdef HAVE_OP_DXX
case 11:
#include INTERP_DXX_COMPUTE_FILE_NAME
result *= inverse_delta_x2;
break;
#endif
#ifdef HAVE_OP_DXY
case 12:
case 21:
#include INTERP_DXY_COMPUTE_FILE_NAME
result *= inverse_delta_x * inverse_delta_y;
break;
#endif
#ifdef HAVE_OP_DXZ
case 13:
case 31:
#include INTERP_DXZ_COMPUTE_FILE_NAME
result *= inverse_delta_x * inverse_delta_z;
break;
#endif
#ifdef HAVE_OP_DYY
case 22:
#include INTERP_DYY_COMPUTE_FILE_NAME
result *= inverse_delta_y2;
break;
#endif
#ifdef HAVE_OP_DYZ
case 23:
case 32:
#include INTERP_DYZ_COMPUTE_FILE_NAME
result *= inverse_delta_y * inverse_delta_z;
break;
#endif
#ifdef HAVE_OP_DZZ
case 33:
#include INTERP_DZZ_COMPUTE_FILE_NAME
result *= inverse_delta_z2;
break;
#endif
default:
CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
"Generalized interpolation operation_code %d not supported!",
operation_codes[out]); /*NOTREACHED*/
return UTIL_ERROR_BAD_INPUT;
/*** ERROR RETURN ***/
/* end of switch (operation_codes[out]) */
}
/*
* ***store*** the result in the output array
*/
switch (output_array_type_codes[out])
{
case CCTK_VARIABLE_REAL:
{
CCTK_REAL *const output_array_ptr_real
= (CCTK_REAL *) output_arrays[out];
output_array_ptr_real[pt] = result;
break;
}
#ifdef HAVE_CCTK_REAL4
case CCTK_VARIABLE_REAL4:
{
CCTK_REAL4 *const output_array_ptr_real4
= (CCTK_REAL4 *) output_arrays[out];
output_array_ptr_real4[pt] = result;
break;
}
#endif
#ifdef HAVE_CCTK_REAL8
case CCTK_VARIABLE_REAL8:
{
CCTK_REAL8 *const output_array_ptr_real8
= (CCTK_REAL8 *) output_arrays[out];
output_array_ptr_real8[pt] = result;
break;
}
#endif
#ifdef HAVE_CCTK_REAL16
case CCTK_VARIABLE_REAL16:
{
CCTK_REAL16 *const output_array_ptr_real16
= (CCTK_REAL16 *) output_arrays[out];
output_array_ptr_real16[pt] = result;
break;
}
#endif
case CCTK_VARIABLE_COMPLEX:
{
CCTK_REAL (*const output_array_ptr_complex)[COMPLEX_N_PARTS]
= (CCTK_REAL (*)[COMPLEX_N_PARTS]) output_arrays[out];
output_array_ptr_complex[pt][part] = result;
break;
}
#ifdef HAVE_CCTK_COMPLEX8
case CCTK_VARIABLE_COMPLEX8:
{
CCTK_REAL4 (*const output_array_ptr_complex8)[COMPLEX_N_PARTS]
= (CCTK_REAL4 (*)[COMPLEX_N_PARTS]) output_arrays[out];
output_array_ptr_complex8[pt][part] = result;
break;
}
#endif /* HAVE_CCTK_COMPLEX8 */
#ifdef HAVE_CCTK_COMPLEX16
case CCTK_VARIABLE_COMPLEX16:
{
CCTK_REAL8 (*const output_array_ptr_complex16)[COMPLEX_N_PARTS]
= (CCTK_REAL8 (*)[COMPLEX_N_PARTS]) output_arrays[out];
output_array_ptr_complex16[pt][part] = result;
break;
}
#endif /* HAVE_CCTK_COMPLEX16 */
#ifdef HAVE_CCTK_COMPLEX32
case CCTK_VARIABLE_COMPLEX32:
{
CCTK_REAL16 (*const output_array_ptr_complex32)[COMPLEX_N_PARTS]
= (CCTK_REAL16 (*)[COMPLEX_N_PARTS]) output_arrays[out];
output_array_ptr_complex32[pt][part] = result;
break;
}
#endif /* HAVE_CCTK_COMPLEX32 */
default:
CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
"output datatype %d not supported",
output_array_type_codes[out]); /*NOTREACHED*/
return UTIL_ERROR_BAD_INPUT; /*** ERROR RETURN ***/
/* end of switch (output type code) */
}
}
/* end of for (part = ...) loop */
}
}
/* end of for (out = ...) loop */
}
}
}
}
/* end of for (pt = ...) loop */
}
return 0; /*** NORMAL RETURN ***/
}
}
}
|