aboutsummaryrefslogtreecommitdiff
path: root/src/gr/expansion.cc
blob: 5a0378354748911b53f359f17090cd657cacde7f (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
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
// expansion.cc -- evaluate expansion Theta of trial horizon surface
// $Header$
//
// <<<prototypes for functions local to this file>>>
// expansion - top-level driver
//
/// setup_xyz_posns - setup global xyz posns of grid points
/// interpolate_geometry - interpolate g_ij and K_ij from Cactus 3-D grid
/// convert_conformal_to_physical - convert conformal gij to physical
///
/// h_is_finite - does function h contain any NaNs/infinities?
/// geometry_is_finite - do geometry vars contain NaN/infty?
///
/// compute_Theta - compute Theta(h) given earlier setup
///

//
// debugging flags
//
#undef GEOMETRY_INTERP_DEBUG	// define this for verbose debugging
				// of geometry interpolator calls/results
#undef GEOMETRY_INTERP_DEBUG2	// define this for even more verbose debugging
				// of geometry interpolator calls/results
#undef GEOMETRY_INTERP_SYNC_SLEEP	// use sleep() to try to ensure we're
					// synchronized across all processors
					// before calling geometry interpolator
#undef GEOMETRY_INTERP_SYNC_BARRIER	// use CCTK_Barrier() to ensure we're
					// synchronized across all processors
					// before calling geometry interpolator
#undef COMPUTE_THETA_DEBUG	// define this for verbose debugging
				// in compute_Theta()

#ifdef GEOMETRY_INTERP_SYNC_SLEEP
#include <unistd.h>		// sleep()
#endif

#include <stdio.h>
#include <assert.h>
#include <math.h>

#include "util_Table.h"
#include "cctk.h"

#include "config.h"
#include "stdc.h"
#include "../jtutil/util.hh"
#include "../jtutil/array.hh"
#include "../jtutil/cpm_map.hh"
#include "../jtutil/linear_map.hh"

#include "../patch/coords.hh"
#include "../patch/grid.hh"
#include "../patch/fd_grid.hh"
#include "../patch/patch.hh"
#include "../patch/patch_edge.hh"
#include "../patch/patch_interp.hh"
#include "../patch/ghost_zone.hh"
#include "../patch/patch_system.hh"

#include "../elliptic/Jacobian.hh"

#include "gfns.hh"
#include "gr.hh"

// all the code in this file is inside this namespace
namespace AHFinderDirect
	  {
using jtutil::error_exit;
using jtutil::pow2;
using jtutil::pow4;

//******************************************************************************
//******************************************************************************
//******************************************************************************

//
// ***** prototypes for functions local to this file *****
//

namespace {
void setup_xyz_posns(patch_system& ps, bool print_msg_flag);
enum expansion_status
  interpolate_geometry(patch_system* ps_ptr,
		       const struct cactus_grid_info& cgi,
		       const struct geometry_info& gi,
		       const struct error_info& error_info, bool initial_flag,
		       bool print_msg_flag);
void convert_conformal_to_physical(patch_system& ps,
				   bool print_msg_flag);
void Schwarzschild_EF_geometry(patch_system& ps,
			       const struct cactus_grid_info& cgi,
			       const struct geometry_info& gi,
			       bool print_msg_flag);

bool h_is_finite(patch_system& ps,
		 const struct error_info& error_info, bool initial_flag,
		 bool print_msg_flag);
bool geometry_is_finite(patch_system& ps,
			const struct error_info& error_info, bool initial_flag,
			bool print_msg_flag);

bool compute_Theta(patch_system& ps, fp add_to_expansion,
		   bool Jacobian_flag, jtutil::norm<fp>* Theta_norms_ptr,
		   const struct error_info& error_info, bool initial_flag,
		   bool print_msg_flag);
	  }

//******************************************************************************
//******************************************************************************
//******************************************************************************

//
// If ps_ptr != NULL, this function computes the LHS function Theta(h),
// and optionally also its Jacobian coefficients (from which the Jacobian
// matrix may be computed later).
//
// If ps_ptr == NULL, this function does a dummy computation, described
// below.
//
// Inputs (angular gridfns, on ghosted grid):
// ... defined on ghosted grid
// ... only values on nominal grid are actually used as input
//	h				# shape of trial surface
//
// Inputs (Cactus 3-D gridfns):
//	gxx,gxy,gxz,gyy,gyz,gzz		# 3-metric $g_{ij}$
//	kxx,kxy,kxz,kyy,kyz,kzz		# extrinsic curvature $K_{ij}$
//
// Outputs (temporaries computed at each grid point)
//	## computed by hand-written code
//	global_[xyz]			# xyz positions of grid points
//	X_ud_*, X_udd_*			# xyz derivative coefficients
//	## computed by Maple-generated code
//	g_uu_{11,12,13,22,23,33}	# $g^{ij}$
//	K				# $K$
//	K_dd_{11,12,13,22,23,33}	# $K^{ij}$
//	partial_d_ln_sqrt_g_d		# $\partial_i \ln \sqrt{g}$
//	partial_d_g_uu_{1,2,3}{11,12,13,22,23,33}	# $\partial_k g^{ij}$
//
// Outputs (angular gridfns, all on the nominal grid):
//	## interpolated from 3-D Cactus grid
//	g_dd_{11,12,13,22,23,33}			# $g_{ij}$
//	K_dd_{11,12,13,22,23,33}			# $K_{ij}$
//	partial_d_g_dd_{1,2,3}{11,12,13,22,23,33}	# $\partial_k g_{ij}$
//	## computed at the nominal grid points
//	Theta						# $\Theta = \Theta(h)$
//
// Arguments:
// ps_ptr --> The patch system, or == NULL to do (only) a dummy computation
//	      in which only the parameter-table setup and a dummy geometry
//	      interpolator call are done, the latter with the number of
//	      interpolation points is set to 0 and all the output array
//	      pointers set to NULL.
// add_to_expansion = A real number which is added to the computed expansion
//		      at each grid point.
// initial_flag = true if this is the first evaluation of  expansion()
//		       for this horizon,
//		  false otherwise;
//		  this is used (only) to select which elements of  error_info
//		  are relevant
// Jacobian_flag = true to compute the Jacobian coefficients,
//		   false to skip this.
// print_msg_flag = true to print status messages,
//		    false to skip this.
// Theta_norms_ptr = (out) If this pointer is non-NULL, the norm object it
//			   points to is updated with all the Theta values
//			   in the grid.  This norm object can then be used
//			   to compute various (gridwise) norms of Theta.
//
// Results:
// This function returns a status code indicating whether the computation
// succeeded or failed, and if the latter, what caused the failure.
//
enum expansion_status
  expansion(patch_system* ps_ptr, fp add_to_expansion,
	    const struct cactus_grid_info& cgi,
	    const struct geometry_info& gi,
	    const struct error_info& error_info, bool initial_flag,
	    bool Jacobian_flag /* = false */,
	    bool print_msg_flag /* = false */,
	    jtutil::norm<fp>* Theta_norms_ptr /* = NULL */)
{
const bool active_flag = (ps_ptr != NULL);
if (print_msg_flag)
   then CCTK_VInfo(CCTK_THORNSTRING,
		   "   %sexpansion",
		   active_flag ? "" : "dummy ");

if (active_flag)
   then {
	//
	// normal computation
	//

	// fill in values of all ghosted gridfns in ghost zones
	ps_ptr->synchronize();

	if (gi.check_that_h_is_finite && !h_is_finite(*ps_ptr,
						      error_info, initial_flag,
						      print_msg_flag))
	   then return expansion_failure__surface_nonfinite;
							// *** ERROR RETURN ***

	// set up xyz positions of grid points
	setup_xyz_posns(*ps_ptr, print_msg_flag);
	}

// compute the "geometry" g_ij, K_ij, and partial_k g_ij
if (gi.hardwire_Schwarzschild_EF_geometry)
   then {
	if (active_flag)
	   then Schwarzschild_EF_geometry(*ps_ptr,
					  gi,
					  print_msg_flag);
	}
   else {
	// this is the only function we call unconditionally; it looks at
	// ps_ptr (non-NULL vs NULL) to choose a normal vs dummy computation
	const enum expansion_status status
		= interpolate_geometry(ps_ptr,
				       cgi, gi,
				       error_info, initial_flag,
				       print_msg_flag);
	if (status != expansion_success)
	   then return status;				// *** ERROR RETURN ***
	if (active_flag && cgi.use_Cactus_conformal_metric)
	   then convert_conformal_to_physical(*ps_ptr,
					      print_msg_flag);
	}

if (active_flag)
   then {
	if (gi.check_that_geometry_is_finite
	    && !geometry_is_finite(*ps_ptr,
				   error_info, initial_flag,
				   print_msg_flag))
	   then return expansion_failure__geometry_nonfinite;
							// *** ERROR RETURN ***

	// compute remaining gridfns --> $\Theta$
	// and optionally also the Jacobian coefficients
	// by algebraic ops and angular finite differencing
	if (!compute_Theta(*ps_ptr, add_to_expansion,
			   Jacobian_flag, Theta_norms_ptr,
			   error_info, initial_flag,
			   print_msg_flag))
	   then return expansion_failure__gij_not_positive_definite;
							// *** ERROR RETURN ***
	}

return expansion_success;				// *** NORMAL RETURN ***
}

//******************************************************************************
//******************************************************************************
//******************************************************************************

//
// This function sets up the xyz-position gridfns:
// * global_[xyz] (used by  interplate_geometry() )
// * global_{xx,xy,xz,yy,yz,zz} (used by higher-level code to compute
//                               quadrupole moments of horizons)
//
// Bugs:
// * We initialize the global_{xx,xy,xz,yy,yz,zz} gridfns every time
//   this function is called, i.e. at each horizon-finding Newton
//   iteration, even though they're only needed at the end of the
//   horizon-finding process.  In practice the extra cost is small,
//   though, so it's probably not worth fixing this...
//
namespace {
void setup_xyz_posns(patch_system& ps, bool print_msg_flag)
{
if (print_msg_flag)
   then CCTK_VInfo(CCTK_THORNSTRING,
		   "      xyz positions and derivative coefficients");

	for (int pn = 0 ; pn < ps.N_patches() ; ++pn)
	{
	patch& p = ps.ith_patch(pn);

		for (int irho = p.min_irho() ; irho <= p.max_irho() ; ++irho)
		{
		for (int isigma = p.min_isigma() ;
		     isigma <= p.max_isigma() ;
		     ++isigma)
		{
		const fp r = p.ghosted_gridfn(gfns::gfn__h, irho,isigma);
		const fp rho = p.rho_of_irho(irho);
		const fp sigma = p.sigma_of_isigma(isigma);

		fp local_x, local_y, local_z;
		p.xyz_of_r_rho_sigma(r,rho,sigma, local_x,local_y,local_z);

		const fp global_x = ps.origin_x() + local_x;
		const fp global_y = ps.origin_y() + local_y;
		const fp global_z = ps.origin_z() + local_z;

		p.gridfn(gfns::gfn__global_x, irho,isigma) = global_x;
		p.gridfn(gfns::gfn__global_y, irho,isigma) = global_y;
		p.gridfn(gfns::gfn__global_z, irho,isigma) = global_z;

		const fp global_xx = global_x * global_x;
		const fp global_xy = global_x * global_y;
		const fp global_xz = global_x * global_z;
		const fp global_yy = global_y * global_y;
		const fp global_yz = global_y * global_z;
		const fp global_zz = global_z * global_z;

		p.gridfn(gfns::gfn__global_xx, irho,isigma) = global_xx;
		p.gridfn(gfns::gfn__global_xy, irho,isigma) = global_xy;
		p.gridfn(gfns::gfn__global_xz, irho,isigma) = global_xz;
		p.gridfn(gfns::gfn__global_yy, irho,isigma) = global_yy;
		p.gridfn(gfns::gfn__global_yz, irho,isigma) = global_yz;
		p.gridfn(gfns::gfn__global_zz, irho,isigma) = global_zz;
		}
		}
	}
}
	  }

//******************************************************************************

//
// If ps_ptr != NULL, this function interpolates the Cactus gridfns
//	gxx...gzz
//	kxx...kzz
//	psi			# optional
// to determine the nominal-grid angular gridfns
//	g_dd_ij
//	partial_d_g_dd_kij
//	K_dd_ij
//	psi			# optional
//	partial_d_psi_k		# optional
// at the nominal-grid trial horizon surface positions given by the
// global_(x,y,z) angular gridfns in the patch system *ps_ptr.  The psi
// interpolation is only done if the cgi.use_Cactus_conformal_metric flag
// is set.  Note that this function ignores the physical-vs-conformal
// semantics of the gridfns; it just interpolates and takes derivatives
// of the stored gridfn values.
//
// If ps_ptr == NULL, this function does (only) the parameter-table
// setup and a a dummy interpolator call, as described in the comments
// to  expansion()  above.
//
// The interpolation is done via  CCTK_InterpGridArrays() .  This has the
// option to return both an overall interpolation status, and a "local"
// status which gives the results of interpolating only the points requested
// on *this* processor; if the local status is available we use it, otherwise
// we fall back to the overall status.
//
// Inputs (angular gridfns, all on the nominal grid):
//	global_[xyz]			# xyz positions of grid points
//
// Inputs (Cactus 3-D gridfns):
//	gxx,gxy,gxz,gyy,gyz,gzz		# 3-metric $g_{ij}$
//					# (may be either physical or conformal)
//	kxx,kxy,kxz,kyy,kyz,kzz		# extrinsic curvature $K_{ij}$
//	psi				# optional conformal factor $\psi$
//
// Outputs (angular gridfns, all on the nominal grid):
//	g_dd_{11,12,13,22,23,33}		# $\stored{g}_{ij}$
//	K_dd_{11,12,13,22,23,33}		# $K_{ij}$
//	partial_d_g_dd_[123]{11,12,13,22,23,33}	# $\partial_k \stored{g}_{ij}$
//	psi					# (optional) $\psi$
//	partial_d_psi_[123]			# (optional) $\partial_k \psi$
//
// This function may also modify the interpolator parameter table.
//
// Results:
// This function returns a status code indicating whether the computation
// succeeded or failed, and if the latter, what caused the failure.  Possible
// failure codes are
// * expansion_failure__surface_outside_grid
// * expansion_failure__surface_in_excised_region	// not implemented yet
//
namespace {
enum expansion_status
  interpolate_geometry(patch_system* ps_ptr,
			  const struct cactus_grid_info& cgi,
			  const struct geometry_info& gi,
			  const struct error_info& error_info, bool initial_flag,
			  bool print_msg_flag)
{
const bool active_flag = (ps_ptr != NULL);
const bool psi_flag = cgi.use_Cactus_conformal_metric;

//
// Implementation Notes:
//
// To handle the optional interpolation of psi, we set up all the data
// type and pointer arrays to include psi, but with the psi entries at
// the end.  We then choose the array sizes passed to the interpolator
// to either include or exclude the psi entries as appropriate.
//
// We remember whether or not psi was interpolated on the previous call,
// and only modify the interpolator parameter table if this changes (or
// if this is our first call).
//

if (print_msg_flag)
   then CCTK_VInfo(CCTK_THORNSTRING,
		   "      interpolating %s from Cactus grid",
		   (psi_flag ? "{g_ij, K_ij, psi}" : "{g_ij, K_ij}"));

int status;

#define CAST_PTR_OR_NULL(type_,ptr_)	\
	(ps_ptr == NULL) ? NULL : static_cast<type_>(ptr_)


//
// ***** interpolation points *****
//
const int N_interp_points = (ps_ptr == NULL) ? 0 : ps_ptr->N_grid_points();
const int interp_coords_type_code = CCTK_VARIABLE_REAL;
const void* const interp_coords[N_GRID_DIMS]
  = {
    CAST_PTR_OR_NULL(const void*, ps_ptr->gridfn_data(gfns::gfn__global_x)),
    CAST_PTR_OR_NULL(const void*, ps_ptr->gridfn_data(gfns::gfn__global_y)),
    CAST_PTR_OR_NULL(const void*, ps_ptr->gridfn_data(gfns::gfn__global_z)),
    };


//
// ***** input arrays *****
//

const CCTK_INT input_array_variable_indices[]
	= {
	  cgi.g_dd_11_varindex, cgi.g_dd_12_varindex, cgi.g_dd_13_varindex,
				cgi.g_dd_22_varindex, cgi.g_dd_23_varindex,
						      cgi.g_dd_33_varindex,
	  cgi.K_dd_11_varindex, cgi.K_dd_12_varindex, cgi.K_dd_13_varindex,
				cgi.K_dd_22_varindex, cgi.K_dd_23_varindex,
						      cgi.K_dd_33_varindex,
	  cgi.psi_varindex,
	  };
const CCTK_INT input_array_type_codes[]
	= {
	  // g_ij
	  CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
			      CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
						  CCTK_VARIABLE_REAL,
	  // K_ij
	  CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
			      CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
						  CCTK_VARIABLE_REAL,
	  // psi
	  CCTK_VARIABLE_REAL,
	  };
const int N_input_arrays_for_psi = 1;
const int N_input_arrays_dim =   sizeof(input_array_variable_indices)
			       / sizeof(input_array_variable_indices[0]);
const int N_input_arrays_use
	= psi_flag ? N_input_arrays_dim
		   : N_input_arrays_dim - N_input_arrays_for_psi;


//
// ***** output arrays *****
//

const CCTK_INT output_array_type_codes[]
	= {
 // $g_{ij}$         $\partial_x g_{ij}$ $\partial_y g_{ij}$ $\partial_z g_{ij}$
 CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
 CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
 CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
 CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
 CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
 CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
 // $K_{ij}$
 CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
		     CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
					 CCTK_VARIABLE_REAL,
 // $\psi$           $\partial_x \psi$   $\partial_y \psi$   $\partial_z \psi$
 CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL,
	  };

const CCTK_INT operand_indices[]
	= {
	  0, 0, 0, 0,		// g_dd_11, partial_[xyz] g_dd_11
	  1, 1, 1, 1,		// g_dd_12, partial_[xyz] g_dd_12
	  2, 2, 2, 2,		// g_dd_13, partial_[xyz] g_dd_13
	  3, 3, 3, 3,		// g_dd_22, partial_[xyz] g_dd_22
	  4, 4, 4, 4,		// g_dd_23, partial_[xyz] g_dd_23
	  5, 5, 5, 5,		// g_dd_33, partial_[xyz] g_dd_33
	  6, 7, 8, 9, 10, 11,	// K_dd_{11,12,13,22,23,33}
	  12, 12, 12, 12,	// psi, partial_[xyz] psi
	  };
#define DERIV(x)	x
const CCTK_INT operation_codes[]
  = {
    DERIV(0), DERIV(1), DERIV(2), DERIV(3), // g_dd_11, partial_[xyz] g_dd_11
    DERIV(0), DERIV(1), DERIV(2), DERIV(3), // g_dd_12, partial_[xyz] g_dd_12
    DERIV(0), DERIV(1), DERIV(2), DERIV(3), // g_dd_13, partial_[xyz] g_dd_13
    DERIV(0), DERIV(1), DERIV(2), DERIV(3), // g_dd_22, partial_[xyz] g_dd_22
    DERIV(0), DERIV(1), DERIV(2), DERIV(3), // g_dd_23, partial_[xyz] g_dd_23
    DERIV(0), DERIV(1), DERIV(2), DERIV(3), // g_dd_33, partial_[xyz] g_dd_33
    DERIV(0), DERIV(0), DERIV(0), DERIV(0), DERIV(0), DERIV(0),
					    // K_dd_{11,12,13,22,23,33}
    DERIV(0), DERIV(1), DERIV(2), DERIV(3), // psi, partial_[xyz] psi
    };

void* const output_arrays[]
  = {
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__g_dd_11)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_111)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_211)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_311)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__g_dd_12)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_112)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_212)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_312)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__g_dd_13)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_113)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_213)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_313)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__g_dd_22)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_122)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_222)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_322)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__g_dd_23)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_123)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_223)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_323)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__g_dd_33)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_133)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_233)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_g_dd_333)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__K_dd_11)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__K_dd_12)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__K_dd_13)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__K_dd_22)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__K_dd_23)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__K_dd_33)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__psi)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_psi_1)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_psi_2)),
    CAST_PTR_OR_NULL(void*, ps_ptr->gridfn_data(gfns::gfn__partial_d_psi_3)),
    };

const int N_output_arrays_for_psi = 4;
const int N_output_arrays_dim
	= sizeof(output_arrays) / sizeof(output_arrays[0]);
const int N_output_arrays_use
	= psi_flag ? N_output_arrays_dim
		   : N_output_arrays_dim - N_output_arrays_for_psi;


//
// ***** parameter table *****
//

// this flag is true if and only if the parameter table already has the
//	suppress_warnings
//	operand_indices
//	operand_codes
// entries for  psi_flag == par_table_psi_flag .
static bool par_table_setup = false;

// if  par_table_setup,
//    this flag is the value of  psi_flag  for the parameter table entries
// otherwise this flag is ignored
static bool par_table_psi_flag = false;

if (par_table_setup && (psi_flag == par_table_psi_flag))
   then {
	// parameter table is already set to just what we need
	// ==> no-op here
	}
   else {
	// store derivative info in interpolator parameter table
	if (print_msg_flag)
	   then CCTK_VInfo(CCTK_THORNSTRING,
			   "         setting up interpolator derivative info");

	// set a dummy value for the key "suppress_warnings"
	// to tell CCTK_InterpGridArrays() not to print warning messages
	// for points outside the grid
	status = Util_TableSetInt(gi.param_table_handle,
				  0,
				  "suppress_warnings");
	if (status < 0)
	   then error_exit(ERROR_EXIT,
"***** interpolate_geometry():\n"
"        unable to set \"suppress_warnings\" key in interpolator parameter table!\n"
"        Util_TableSetInt() status=%d\n"
			   ,
			   status);				/*NOTREACHED*/

	status = Util_TableSetIntArray(gi.param_table_handle,
				       N_output_arrays_use, operand_indices,
				       "operand_indices");
	if (status < 0)
	   then error_exit(ERROR_EXIT,
"***** interpolate_geometry():\n"
"        unable to set operand_indices in interpolator parameter table!\n"
"        Util_TableSetIntArray() status=%d\n"
			   ,
			   status);				/*NOTREACHED*/

	status = Util_TableSetIntArray(gi.param_table_handle,
				       N_output_arrays_use, operation_codes,
				       "operation_codes");
	if (status < 0)
	   then error_exit(ERROR_EXIT,
"***** interpolate_geometry():\n"
"        unable to set operation_codes in interpolator parameter table!\n"
"        Util_TableSetIntArray() status=%d\n"
			   ,
			   status);				/*NOTREACHED*/

	par_table_setup = true;
	par_table_psi_flag = psi_flag;
	}


//
// ***** the actual interpolation *****
//
if (print_msg_flag)
   then CCTK_VInfo(CCTK_THORNSTRING,
		   "         calling geometry interpolator (%s%d points)",
		   (active_flag ? "" : "dummy: "), N_interp_points);

#ifdef GEOMETRY_INTERP_DEBUG
printf("AHFinderDirect:: proc %d: initializing interpolator outputs to 999.999\n",
       int(CCTK_MyProc(cgi.GH)));
	  {
	for (int pt = 0 ; pt < N_interp_points ; ++pt)
	{
		for (int out = 0 ; out < N_output_arrays_use ; ++out)
		{
		CCTK_REAL* const out_ptr
			= static_cast<CCTK_REAL*>(output_arrays[out]);
		out_ptr[pt] = 999.999;
		}
	}
	  }
#endif

#ifdef GEOMETRY_INTERP_DEBUG
printf("AHFinderDirect:: proc %d: calling CCTK_InterpGridArrays(N_interp_points=%d)\n",
       int(CCTK_MyProc(cgi.GH)), N_interp_points);
fflush(stdout);
#endif

#ifdef GEOMETRY_INTERP_SYNC_SLEEP
sleep(1);
#endif
#ifdef GEOMETRY_INTERP_SYNC_BARRIER
CCTK_Barrier(cgi.GH);
#endif

status = CCTK_InterpGridArrays(cgi.GH, N_GRID_DIMS,
			       gi.operator_handle, gi.param_table_handle,
			       cgi.coord_system_handle,
			       N_interp_points,
				  interp_coords_type_code,
				  interp_coords,
			       N_input_arrays_use,
				  input_array_variable_indices,
			       N_output_arrays_use,
				  output_array_type_codes,
				  output_arrays);

#ifdef GEOMETRY_INTERP_DEBUG
printf("AHFinderDirect:: proc %d: CCTK_InterpGridArrays() returned status=%d\n",
       int(CCTK_MyProc(cgi.GH)), status);
fflush(stdout);
#endif

#ifdef GEOMETRY_INTERP_DEBUG2
	  {
	for (int pt = 0 ; pt < N_interp_points ; pt = 2*pt + (pt == 0))
	{
	printf("AHFinderDirect:: proc %d: CCTK_InterpGridArrays() results for pt=%d:\n",
	       int(CCTK_MyProc(cgi.GH)), pt);
		for (int out = 0 ; out < N_output_arrays_use ; ++out)
		{
		const CCTK_REAL* const out_ptr
			= static_cast<const CCTK_REAL*>(output_arrays[out]);
		printf("   out=%d   result=%g\n", out, double(out_ptr[pt]));
		}
	}
	  }
#endif	/* GEOMETRY_INTERP_DEBUG2 */


//
// ***** get the local interpolator status (if available)
//
// CCTK_InterpGridArrays() returns a status that reflects the interpolation
// on *all* processors, but what we really want here is just the status
// for the interpolation of the points that *we* (this processor) asked for
// ==> get the local status if it's available
//
CCTK_INT local_interpolator_status;
const int table_status = Util_TableGetInt(gi.param_table_handle,
					  &local_interpolator_status,
					  "local_interpolator_status");
if	(table_status == 1)
   then {
	// we got the local interpolator status successfully
	// ==> use it for all further checks
	#ifdef GEOMETRY_INTERP_DEBUG
	printf("AHFinderDirect:: proc %d: replacing status=%d with local_interpolator_status=%d\n",
	       int(CCTK_MyProc(cgi.GH)), status, int(local_interpolator_status));
	fflush(stdout);
	#endif
	status = local_interpolator_status;
	}
else if (table_status == UTIL_ERROR_TABLE_NO_SUCH_KEY)
   then {
	// evidently the interpolators don't provide the
	// local interpolator status
	// ==> stick with the status have
	// ==> no-op here
	}
else	error_exit(ERROR_EXIT,
"***** interpolate_geometry():\n"
"        error return %d trying to get local interpolator status\n"
"        from parameter table!  (CCTK_InterpGridArrays() status=%d)\n"
		   ,
		   table_status, status);		/*NOTREACHED*/


//
// ***** handle any interpolation errors *****
//
if (status == CCTK_ERROR_INTERP_POINT_OUTSIDE)
   then {
	if (print_msg_flag)
	   then {
		// see if we can get further info
		const int warn_level
		   = initial_flag
		     ? error_info.warn_level__point_outside__initial
		     : error_info.warn_level__point_outside__subsequent;


		CCTK_VWarn(warn_level, __LINE__, __FILE__, CCTK_THORNSTRING,
"interpolate_geometry():\n"
"        one or more points on the trial horizon surface point"
"        is/are outside the grid (or too close to the grid boundary)"
"        (in a single-processor run, this may also mean that\n"
"         driver::ghost_size is too small for this geometry interpolator)\n");
		}

	return expansion_failure__surface_outside_grid;	// *** ERROR RETURN ***
	}

else if (status == CCTK_ERROR_INTERP_GHOST_SIZE_TOO_SMALL)
   then error_exit(ERROR_EXIT,
"***** interpolate_geometry(): driver::ghost_size is too small\n"
"                              for this geometry interpolator!\n");
								/*NOTREACHED*/

else if (status < 0)
   then error_exit(ERROR_EXIT,
"***** interpolate_geometry(): error return %d from interpolator!\n",
		   status);					/*NOTREACHED*/

return expansion_success;				// *** NORMAL RETURN ***
}
	  }

//******************************************************************************

//
// This function converts the g_dd_ij and partial_d_g_dd_kij gridfns
// from the Cactus conformal semantics to the physical semantics.  As
// documented in the CactusEinstein/ConformalState thorn guide, the
// Cactus conformal semantics are
//	physical_gij = psi^4 * stored_gij
// From this it's trivial to derive the transformation for the derivatives,
//	partial_k physical_gij = 4 psi^3 (partial_k psi) stored_gij
//				 + psi^4 partial_k stored_gij
//
// Inputs (angular gridfns, all on the nominal grid):
//   psi					# $\psi$
//   partial_d_psi_{1,2,3}			# $\partial_k \psi$
//   g_dd_{11,12,13,22,23,33}			# $\stored{g}_{ij}$
//   partial_d_g_dd_{1,2,3}{11,12,13,22,23,33}	# $\partial_k \stored{g}_{ij}$
//
//
// Outputs (angular gridfns, all on the nominal grid):
//   g_dd_{11,12,13,22,23,33}			# $g_{ij}$
//   partial_d_g_dd_{1,2,3}{11,12,13,22,23,33}	# $\partial_k g_{ij}$
//
namespace {
void convert_conformal_to_physical(patch_system& ps, bool print_msg_flag)
{
if (print_msg_flag)
   then CCTK_VInfo(CCTK_THORNSTRING,
		   "      converting Cactus conformal gij --> physical gij");

    for (int pn = 0 ; pn < ps.N_patches() ; ++pn)
    {
    patch& p = ps.ith_patch(pn);

	for (int irho = p.min_irho() ; irho <= p.max_irho() ; ++irho)
	{
	for (int isigma = p.min_isigma() ; isigma <= p.max_isigma() ; ++isigma)
	{
	const fp psi = p.gridfn(gfns::gfn__psi, irho,isigma);
	const fp psi3 = jtutil::pow3(psi);
	const fp psi4 = jtutil::pow4(psi);

	const fp partial_d_psi_1
		= p.gridfn(gfns::gfn__partial_d_psi_1, irho,isigma);
	const fp partial_d_psi_2
		= p.gridfn(gfns::gfn__partial_d_psi_2, irho,isigma);
	const fp partial_d_psi_3
		= p.gridfn(gfns::gfn__partial_d_psi_3, irho,isigma);

	const fp stored_g_dd_11 = p.gridfn(gfns::gfn__g_dd_11, irho, isigma);
	const fp stored_g_dd_12 = p.gridfn(gfns::gfn__g_dd_12, irho, isigma);
	const fp stored_g_dd_13 = p.gridfn(gfns::gfn__g_dd_13, irho, isigma);
	const fp stored_g_dd_22 = p.gridfn(gfns::gfn__g_dd_22, irho, isigma);
	const fp stored_g_dd_23 = p.gridfn(gfns::gfn__g_dd_23, irho, isigma);
	const fp stored_g_dd_33 = p.gridfn(gfns::gfn__g_dd_33, irho, isigma);

	p.gridfn(gfns::gfn__g_dd_11, irho, isigma) *= psi4;
	p.gridfn(gfns::gfn__g_dd_12, irho, isigma) *= psi4;
	p.gridfn(gfns::gfn__g_dd_13, irho, isigma) *= psi4;
	p.gridfn(gfns::gfn__g_dd_22, irho, isigma) *= psi4;
	p.gridfn(gfns::gfn__g_dd_23, irho, isigma) *= psi4;
	p.gridfn(gfns::gfn__g_dd_33, irho, isigma) *= psi4;

	p.gridfn(gfns::gfn__partial_d_g_dd_111, irho,isigma)
		= 4.0*psi3*partial_d_psi_1*stored_g_dd_11
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_111, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_112, irho,isigma)
		= 4.0*psi3*partial_d_psi_1*stored_g_dd_12
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_112, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_113, irho,isigma)
		= 4.0*psi3*partial_d_psi_1*stored_g_dd_13
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_113, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_122, irho,isigma)
		= 4.0*psi3*partial_d_psi_1*stored_g_dd_22
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_122, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_123, irho,isigma)
		= 4.0*psi3*partial_d_psi_1*stored_g_dd_23
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_123, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_133, irho,isigma)
		= 4.0*psi3*partial_d_psi_1*stored_g_dd_33
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_133, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_211, irho,isigma)
		= 4.0*psi3*partial_d_psi_2*stored_g_dd_11
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_211, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_212, irho,isigma)
		= 4.0*psi3*partial_d_psi_2*stored_g_dd_12
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_212, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_213, irho,isigma)
		= 4.0*psi3*partial_d_psi_2*stored_g_dd_13
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_213, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_222, irho,isigma)
		= 4.0*psi3*partial_d_psi_2*stored_g_dd_22
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_222, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_223, irho,isigma)
		= 4.0*psi3*partial_d_psi_2*stored_g_dd_23
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_223, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_233, irho,isigma)
		= 4.0*psi3*partial_d_psi_2*stored_g_dd_33
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_233, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_311, irho,isigma)
		= 4.0*psi3*partial_d_psi_3*stored_g_dd_11
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_311, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_312, irho,isigma)
		= 4.0*psi3*partial_d_psi_3*stored_g_dd_12
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_312, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_313, irho,isigma)
		= 4.0*psi3*partial_d_psi_3*stored_g_dd_13
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_313, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_322, irho,isigma)
		= 4.0*psi3*partial_d_psi_3*stored_g_dd_22
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_322, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_323, irho,isigma)
		= 4.0*psi3*partial_d_psi_3*stored_g_dd_23
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_323, irho,isigma);
	p.gridfn(gfns::gfn__partial_d_g_dd_333, irho,isigma)
		= 4.0*psi3*partial_d_psi_3*stored_g_dd_33
		  + psi4*p.gridfn(gfns::gfn__partial_d_g_dd_333, irho,isigma);
	}
	}
    }
}
	  }

//******************************************************************************
//******************************************************************************
//******************************************************************************

//
// This function whether the horizon shape function  h  (nominal grid only)
// contains all finite floating-point numbers, or whether it contain one
// or more NaNs or infinities.
//
// Results:
#ifdef HAVE_FINITE
// This function returns true if all the h values are finite, false
// otherwise (i.e. if it contains any NaNs or infinities).
#else
// This function is a no-op, and just prints a warning message and
// returns true.
#endif
//
namespace {
bool h_is_finite(patch_system& ps,
		 const struct error_info& error_info, bool initial_flag,
		 bool print_msg_flag)
{
if (print_msg_flag)
   then CCTK_VInfo(CCTK_THORNSTRING, "      checking that h is finite");

#ifdef HAVE_FINITE
	for (int pn = 0 ; pn < ps.N_patches() ; ++pn)
	{
	patch& p = ps.ith_patch(pn);

		for (int irho = p.min_irho() ; irho <= p.max_irho() ; ++irho)
		{
		for (int isigma = p.min_isigma() ;
		     isigma <= p.max_isigma() ;
		     ++isigma)
		{
		const fp h = p.ghosted_gridfn(gfns::gfn__h, irho,isigma);
		if (!finite(h))
		   then {
			const fp rho = p.rho_of_irho(irho);
			const fp sigma = p.sigma_of_isigma(isigma);
			const fp drho   = jtutil::degrees_of_radians(rho);
			const fp dsigma = jtutil::degrees_of_radians(sigma);
			CCTK_VWarn(error_info.warn_level__nonfinite_geometry,
				   __LINE__, __FILE__, CCTK_THORNSTRING,
"\n"
"   h=%g isn't finite!\n"
"   %s patch (rho,sigma)=(%g,%g) (drho,dsigma)=(%g,%g)\n"
				   ,
				   double(h),
				   p.name(), double(rho), double(sigma),
					     double(drho), double(dsigma));
			return false;			// *** found a NaN ***
			}
		}
		}
	}
return true;					// *** all values finite ***
#else
CCTK_VWarn(error_info.warn_level__skipping_finite_check,
	   __LINE__, __FILE__, CCTK_THORNSTRING,
           "      no finite() fn ==>  skipping is-h-finite check");
return true;					// *** no check possible ***
#endif
}
	  }

//******************************************************************************

//
// This function whether the geometry variables
//	g_dd_ij
//	partial_d_g_dd_kij
//	K_dd_ij
// are all finite floating-point numbers, or whether they contain one
// or more NaNs or infinities.
//
// Results:
#ifdef HAVE_FINITE
// This function returns true if all the geometry variables are finite,
// false otherwise (i.e. if they contain any NaNs or infinities).
#else
// This function is a no-op, and just prints a warning message and
// returns true.
#endif
//
namespace {
bool geometry_is_finite(patch_system& ps,
			const struct error_info& error_info, bool initial_flag,
			bool print_msg_flag)
{
if (print_msg_flag)
   then CCTK_VInfo(CCTK_THORNSTRING, "      checking that geometry is finite");

#ifdef HAVE_FINITE
    for (int pn = 0 ; pn < ps.N_patches() ; ++pn)
    {
    patch& p = ps.ith_patch(pn);

	for (int irho = p.min_irho() ; irho <= p.max_irho() ; ++irho)
	{
	for (int isigma = p.min_isigma() ;
	     isigma <= p.max_isigma() ;
	     ++isigma)
	{
	const fp g_dd_11 = p.gridfn(gfns::gfn__g_dd_11, irho,isigma);
	const fp g_dd_12 = p.gridfn(gfns::gfn__g_dd_12, irho,isigma);
	const fp g_dd_13 = p.gridfn(gfns::gfn__g_dd_13, irho,isigma);
	const fp g_dd_22 = p.gridfn(gfns::gfn__g_dd_22, irho,isigma);
	const fp g_dd_23 = p.gridfn(gfns::gfn__g_dd_23, irho,isigma);
	const fp g_dd_33 = p.gridfn(gfns::gfn__g_dd_33, irho,isigma);

	const fp K_dd_11 = p.gridfn(gfns::gfn__K_dd_11, irho,isigma);
	const fp K_dd_12 = p.gridfn(gfns::gfn__K_dd_12, irho,isigma);
	const fp K_dd_13 = p.gridfn(gfns::gfn__K_dd_13, irho,isigma);
	const fp K_dd_22 = p.gridfn(gfns::gfn__K_dd_22, irho,isigma);
	const fp K_dd_23 = p.gridfn(gfns::gfn__K_dd_23, irho,isigma);
	const fp K_dd_33 = p.gridfn(gfns::gfn__K_dd_33, irho,isigma);

	const fp partial_d_g_dd_111
		= p.gridfn(gfns::gfn__partial_d_g_dd_111, irho,isigma);
	const fp partial_d_g_dd_112
		= p.gridfn(gfns::gfn__partial_d_g_dd_112, irho,isigma);
	const fp partial_d_g_dd_113
		= p.gridfn(gfns::gfn__partial_d_g_dd_113, irho,isigma);
	const fp partial_d_g_dd_122
		= p.gridfn(gfns::gfn__partial_d_g_dd_122, irho,isigma);
	const fp partial_d_g_dd_123
		= p.gridfn(gfns::gfn__partial_d_g_dd_123, irho,isigma);
	const fp partial_d_g_dd_133
		= p.gridfn(gfns::gfn__partial_d_g_dd_133, irho,isigma);
	const fp partial_d_g_dd_211
		= p.gridfn(gfns::gfn__partial_d_g_dd_211, irho,isigma);
	const fp partial_d_g_dd_212
		= p.gridfn(gfns::gfn__partial_d_g_dd_212, irho,isigma);
	const fp partial_d_g_dd_213
		= p.gridfn(gfns::gfn__partial_d_g_dd_213, irho,isigma);
	const fp partial_d_g_dd_222
		= p.gridfn(gfns::gfn__partial_d_g_dd_222, irho,isigma);
	const fp partial_d_g_dd_223
		= p.gridfn(gfns::gfn__partial_d_g_dd_223, irho,isigma);
	const fp partial_d_g_dd_233
		= p.gridfn(gfns::gfn__partial_d_g_dd_233, irho,isigma);
	const fp partial_d_g_dd_311
		= p.gridfn(gfns::gfn__partial_d_g_dd_311, irho,isigma);
	const fp partial_d_g_dd_312
		= p.gridfn(gfns::gfn__partial_d_g_dd_312, irho,isigma);
	const fp partial_d_g_dd_313
		= p.gridfn(gfns::gfn__partial_d_g_dd_313, irho,isigma);
	const fp partial_d_g_dd_322
		= p.gridfn(gfns::gfn__partial_d_g_dd_322, irho,isigma);
	const fp partial_d_g_dd_323
		= p.gridfn(gfns::gfn__partial_d_g_dd_323, irho,isigma);
	const fp partial_d_g_dd_333
		= p.gridfn(gfns::gfn__partial_d_g_dd_333, irho,isigma);

	if (    !finite(g_dd_11) || !finite(g_dd_12) || !finite(g_dd_13)
				 || !finite(g_dd_22) || !finite(g_dd_23)
						     || !finite(g_dd_33)
	     || !finite(K_dd_11) || !finite(K_dd_12) || !finite(K_dd_13)
				 || !finite(K_dd_22) || !finite(K_dd_23)
						     || !finite(K_dd_33)
	     || !finite(partial_d_g_dd_111)
	     || !finite(partial_d_g_dd_112)
	     || !finite(partial_d_g_dd_113)
	     || !finite(partial_d_g_dd_122)
	     || !finite(partial_d_g_dd_123)
	     || !finite(partial_d_g_dd_133)
	     || !finite(partial_d_g_dd_211)
	     || !finite(partial_d_g_dd_212)
	     || !finite(partial_d_g_dd_213)
	     || !finite(partial_d_g_dd_222)
	     || !finite(partial_d_g_dd_223)
	     || !finite(partial_d_g_dd_233)
	     || !finite(partial_d_g_dd_311)
	     || !finite(partial_d_g_dd_312)
	     || !finite(partial_d_g_dd_313)
	     || !finite(partial_d_g_dd_322)
	     || !finite(partial_d_g_dd_323)
	     || !finite(partial_d_g_dd_333)    )
	   then {
		const fp h = p.ghosted_gridfn(gfns::gfn__h, irho,isigma);
		const fp rho = p.rho_of_irho(irho);
		const fp sigma = p.sigma_of_isigma(isigma);
		const fp drho   = jtutil::degrees_of_radians(rho);
		const fp dsigma = jtutil::degrees_of_radians(sigma);
		fp local_x, local_y, local_z;
		p.xyz_of_r_rho_sigma(h,rho,sigma, local_x,local_y,local_z);
		const fp global_x = ps.origin_x() + local_x;
		const fp global_y = ps.origin_y() + local_y;
		const fp global_z = ps.origin_z() + local_z;
		CCTK_VWarn(error_info.warn_level__nonfinite_geometry,
			   __LINE__, __FILE__, CCTK_THORNSTRING,
"\n"
"   geometry isn't finite at %s patch\n"
"   h=%g (rho,sigma)=(%g,%g) (drho,dsigma)=(%g,%g)\n"
"   local_(x,y,z)=(%g,%g,%g)\n"
"   global_(x,y,z)=(%g,%g,%g)\n"
"   g_dd_11=%g   _12=%g   _13=%g\n"
"       _22=%g   _23=%g   _33=%g\n"
"   K_dd_11=%g   _12=%g   _13=%g\n"
"       _22=%g   _23=%g   _33=%g\n"
"   partial_d_g_dd_111=%g   _112=%g   _113=%g\n"
"                 _122=%g   _123=%g   _133=%g\n"
"   partial_d_g_dd_211=%g   _212=%g   _213=%g\n"
"                 _222=%g   _223=%g   _233=%g\n"
"   partial_d_g_dd_311=%g   _312=%g   _313=%g\n"
"                 _322=%g   _323=%g   _333=%g\n"
			   ,
			   p.name(),
			   double(h), double(rho), double(sigma),
				      double(drho), double(dsigma),
			   double(local_x), double(local_y), double(local_z),
			   double(global_x), double(global_y), double(global_z),
			   double(g_dd_11), double(g_dd_12), double(g_dd_13),
			   double(g_dd_22), double(g_dd_23), double(g_dd_33),
			   double(K_dd_11), double(K_dd_12), double(K_dd_13),
			   double(K_dd_22), double(K_dd_23), double(K_dd_33),
			   double(partial_d_g_dd_111),
			   double(partial_d_g_dd_112),
			   double(partial_d_g_dd_113),
			   double(partial_d_g_dd_122),
			   double(partial_d_g_dd_123),
			   double(partial_d_g_dd_133),
			   double(partial_d_g_dd_211),
			   double(partial_d_g_dd_212),
			   double(partial_d_g_dd_213),
			   double(partial_d_g_dd_222),
			   double(partial_d_g_dd_223),
			   double(partial_d_g_dd_233),
			   double(partial_d_g_dd_311),
			   double(partial_d_g_dd_312),
			   double(partial_d_g_dd_313),
			   double(partial_d_g_dd_322),
			   double(partial_d_g_dd_323),
			   double(partial_d_g_dd_333));
		return false;			// *** found a NaN ***
		}
	}
	}
    }
return true;						// *** no NaNs found ***
#else
CCTK_VWarn(error_info.warn_level__skipping_finite_check,
	   __LINE__, __FILE__, CCTK_THORNSTRING,
           "      no finite() ==>  skipping is-geometry-finite check");
return true;					// *** no check possible ***
#endif
}
	  }

//******************************************************************************
//******************************************************************************
//******************************************************************************

//
// This function computes the expansion Theta(h), and optionally also
// its Jacobian coefficients, (from which the Jacobian matrix may be
// computed later).  This function uses a mixture of algebraic operations
// and (rho,sigma) finite differencing.  The computation is done entirely
// on the nominal angular grid.
//
// N.b. This function #includes "cg.hh", which defines "dangerous" macros
//      which will stay in effect for the rest of this compilation unit!
//
// Arguments:
// Jacobian_flag = true to compute the Jacobian coefficients,
//		   false to skip this.
//
// Results:
// This function returns true for a successful computation, or false
// if the computation failed because Theta_D <= 0 (this means the interpolated
// g_ij isn't positive definite).
//
namespace {
bool compute_Theta(patch_system& ps, fp add_to_expansion,
		   bool Jacobian_flag, jtutil::norm<fp>* Theta_norms_ptr,
		   const struct error_info& error_info, bool initial_flag,
		   bool print_msg_flag)
{
if (print_msg_flag)
   then CCTK_VInfo(CCTK_THORNSTRING, "      computing Theta(h)");

	for (int pn = 0 ; pn < ps.N_patches() ; ++pn)
	{
	patch& p = ps.ith_patch(pn);

		for (int irho = p.min_irho() ; irho <= p.max_irho() ; ++irho)
		{
		for (int isigma = p.min_isigma() ;
		     isigma <= p.max_isigma() ;
		     ++isigma)
		{
		//
		// compute the X_ud and X_udd derivative coefficients
		// ... n.b. this uses the *local* (x,y,z) coordinates
		//
		const fp r = p.ghosted_gridfn(gfns::gfn__h, irho,isigma);
		const fp rho = p.rho_of_irho(irho);
		const fp sigma = p.sigma_of_isigma(isigma);
		fp xx, yy, zz;
		p.xyz_of_r_rho_sigma(r,rho,sigma, xx, yy, zz);
		#ifdef COMPUTE_THETA_DEBUG
		if (    ((irho == 0) && (isigma == 0))
		     || ((irho == 5) && (isigma == 5))    )
		   then {
			printf("AHFinderDirect:: computing at %s patch (%d,%d) r=%g ==> local xyz=(%g,%g,%g)\n",
			       p.name(), irho,isigma, double(r),
			       double(xx), double(yy), double(zz));
			printf("AHFinderDirect:: got g_dd_11=%g g_dd_12=%g g_dd_13=%g\n",
			       p.gridfn(gfns::gfn__g_dd_11, irho,isigma),
			       p.gridfn(gfns::gfn__g_dd_12, irho,isigma),
			       p.gridfn(gfns::gfn__g_dd_13, irho,isigma));
			printf("                                g_dd_22=%g g_dd_23=%g\n",
			       p.gridfn(gfns::gfn__g_dd_22, irho,isigma),
			       p.gridfn(gfns::gfn__g_dd_23, irho,isigma));
			printf("                                           g_dd_33=%g\n",
			       p.gridfn(gfns::gfn__g_dd_33, irho,isigma));
			fflush(stdout);
			}
		#endif

		// 1st derivative coefficients X_ud
		const fp X_ud_11 = p.partial_rho_wrt_x(xx, yy, zz);
		const fp X_ud_12 = p.partial_rho_wrt_y(xx, yy, zz);
		const fp X_ud_13 = p.partial_rho_wrt_z(xx, yy, zz);
		const fp X_ud_21 = p.partial_sigma_wrt_x(xx, yy, zz);
		const fp X_ud_22 = p.partial_sigma_wrt_y(xx, yy, zz);
		const fp X_ud_23 = p.partial_sigma_wrt_z(xx, yy, zz);

		// 2nd derivative coefficient gridfns X_udd
		const fp X_udd_111 = p.partial2_rho_wrt_xx(xx, yy, zz);
		const fp X_udd_112 = p.partial2_rho_wrt_xy(xx, yy, zz);
		const fp X_udd_113 = p.partial2_rho_wrt_xz(xx, yy, zz);
		const fp X_udd_122 = p.partial2_rho_wrt_yy(xx, yy, zz);
		const fp X_udd_123 = p.partial2_rho_wrt_yz(xx, yy, zz);
		const fp X_udd_133 = p.partial2_rho_wrt_zz(xx, yy, zz);
		const fp X_udd_211 = p.partial2_sigma_wrt_xx(xx, yy, zz);
		const fp X_udd_212 = p.partial2_sigma_wrt_xy(xx, yy, zz);
		const fp X_udd_213 = p.partial2_sigma_wrt_xz(xx, yy, zz);
		const fp X_udd_222 = p.partial2_sigma_wrt_yy(xx, yy, zz);
		const fp X_udd_223 = p.partial2_sigma_wrt_yz(xx, yy, zz);
		const fp X_udd_233 = p.partial2_sigma_wrt_zz(xx, yy, zz);

		//
		// "call" the Maple-generated code
		// ... each cg/*.c file has a separate set of temp variables,
		//     and so must be inside its own set of { } braces
		//

		// gridfn #defines
		#include "cg.hh"

		  {
		// g_uu
		#include "../gr.cg/inverse_metric.c"
		  }

		  {
		// K, K_uu
		#include "../gr.cg/extrinsic_curvature_trace_raise.c"
		  }

		  {
		// partial_d_g_uu
		#include "../gr.cg/inverse_metric_gradient.c"
		  }

		  {
		// partial_d_ln_sqrt_g
		#include "../gr.cg/metric_det_gradient.c"
		  }

		  {
		// Theta_A, Theta_B, Theta_C, Theta_D
		#include "../gr.cg/expansion.c"
		  }

		if (Theta_D <= 0)
		   then {
const int warn_level
  = initial_flag ? error_info.warn_level__gij_not_positive_definite__initial
		 : error_info.warn_level__gij_not_positive_definite__subsequent;
CCTK_VWarn(warn_level, __LINE__, __FILE__, CCTK_THORNSTRING,
"\n"
"   compute_Theta(): Theta_D = $g^{ij} s_i s_j$ = %g <= 0\n"
"                    at %s patch rho=%g sigma=%g!\n"
"                    (i.e. the interpolated g_ij isn't positive definite)",
	   double(Theta_D),
	   p.name(), double(rho), double(sigma));
			return false;			// *** ERROR RETURN ***
			}

		// compute H via equation (14) of my 1996 horizon finding paper
		const fp sqrt_Theta_D = sqrt(Theta_D);
		Theta = + Theta_A/(Theta_D*sqrt_Theta_D)
			+ Theta_B/sqrt_Theta_D
			+ Theta_C/Theta_D
			- K
			+ add_to_expansion;

		// update running norms of Theta(h) function
		if (Theta_norms_ptr != NULL)
		   then Theta_norms_ptr->data(Theta);

		if (Jacobian_flag)
		   then {
			// partial_Theta_wrt_partial_d_h,
			// partial_Theta_wrt_partial_dd_h
			#include "../gr.cg/expansion_Jacobian.c"
			}
		}
		}
	}

return true;						// *** NORMAL RETURN ***
}
	  }

//******************************************************************************
//******************************************************************************
//******************************************************************************

	  }	// namespace AHFinderDirect