aboutsummaryrefslogtreecommitdiff
path: root/src/gr/driver.cc
blob: f63bcd476da7440135993908bd5e674d0392065b (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
// driver.cc -- top level driver for finding apparent horizons
// $Id$
//
// AHFinderDirect_driver - top-level driver
//

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

#include "util_Table.h"
#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"

#include "jt/stdc.h"
#include "jt/util.hh"
#include "jt/array.hh"
#include "jt/cpm_map.hh"
#include "jt/linear_map.hh"
using jtutil::error_exit;

#include "../config.hh"

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

#include "gfn.hh"
#include "AHFinderDirect.hh"

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

//
// This function is the Cactus interface for the test driver.
//
extern "C"
  void AHFinderDirect_driver(CCTK_ARGUMENTS)
{
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS

CCTK_VInfo(CCTK_THORNSTRING, "initializing AHFinderDirect data structures");

//
// set up the interpatch interpolator
//
CCTK_VInfo(CCTK_THORNSTRING, "   setting up interpatch interpolator");
const int interp_handle = CCTK_InterpHandle(interpatch_interpolator_name);
if (interp_handle < 0)
   then CCTK_VWarn(-1, __LINE__, __FILE__, CCTK_THORNSTRING,
		   "couldn't find interpolator \"%s\"!",
		   interpatch_interpolator_name);		/*NOTREACHED*/
const int interp_par_table_handle
	= Util_TableCreateFromString(interpatch_interpolator_pars);
if (interp_par_table_handle < 0)
   then CCTK_VWarn(-1, __LINE__, __FILE__, CCTK_THORNSTRING,
		   "bad interpatch-interpolator parameter(s) \"%s\"!",
		   interpatch_interpolator_pars);		/*NOTREACHED*/

//
// set up the Cactus grid info
//
CCTK_VInfo(CCTK_THORNSTRING, "   setting up Cactus grid info");
struct cactus_grid_info cgi;
cgi.GH = cctkGH;
cgi.coord_origin[0] = cctk_origin_space[0];
cgi.coord_origin[1] = cctk_origin_space[1];
cgi.coord_origin[2] = cctk_origin_space[2];
cgi.coord_delta[0] = cctk_delta_space[0];
cgi.coord_delta[1] = cctk_delta_space[1];
cgi.coord_delta[2] = cctk_delta_space[2];
cgi.gridfn_dims[0] = cctk_lsh[0];
cgi.gridfn_dims[1] = cctk_lsh[1];
cgi.gridfn_dims[2] = cctk_lsh[2];
cgi.g_dd_11_data = static_cast<const fp *>(
			CCTK_VarDataPtr(cctkGH, 0, "einstein::gxx")
					  );
cgi.g_dd_12_data = static_cast<const fp *>(
			CCTK_VarDataPtr(cctkGH, 0, "einstein::gxy")
					  );
cgi.g_dd_13_data = static_cast<const fp *>(
			CCTK_VarDataPtr(cctkGH, 0, "einstein::gxz")
					  );
cgi.g_dd_22_data = static_cast<const fp *>(
			CCTK_VarDataPtr(cctkGH, 0, "einstein::gyy")
					  );
cgi.g_dd_23_data = static_cast<const fp *>(
			CCTK_VarDataPtr(cctkGH, 0, "einstein::gyz")
					  );
cgi.g_dd_33_data = static_cast<const fp *>(
			CCTK_VarDataPtr(cctkGH, 0, "einstein::gzz")
					  );
cgi.K_dd_11_data = static_cast<const fp *>(
			CCTK_VarDataPtr(cctkGH, 0, "einstein::kxx")
					  );
cgi.K_dd_12_data = static_cast<const fp *>(
			CCTK_VarDataPtr(cctkGH, 0, "einstein::kxy")
					  );
cgi.K_dd_13_data = static_cast<const fp *>(
			CCTK_VarDataPtr(cctkGH, 0, "einstein::kxz")
					  );
cgi.K_dd_22_data = static_cast<const fp *>(
			CCTK_VarDataPtr(cctkGH, 0, "einstein::kyy")
					  );
cgi.K_dd_23_data = static_cast<const fp *>(
			CCTK_VarDataPtr(cctkGH, 0, "einstein::kyz")
					  );
cgi.K_dd_33_data = static_cast<const fp *>(
			CCTK_VarDataPtr(cctkGH, 0, "einstein::kzz")
					  );

//
// set up the geometry interpolator
//
struct geometry_interpolator_info gii;
CCTK_VInfo(CCTK_THORNSTRING, "   setting up geometry interpolator");
gii.operator_handle = CCTK_InterpHandle(geometry_interpolator_name);
if (interp_handle < 0)
   then CCTK_VWarn(-1, __LINE__, __FILE__, CCTK_THORNSTRING,
		   "couldn't find interpolator \"%s\"!",
		   geometry_interpolator_name);		/*NOTREACHED*/
gii.param_table_handle
	= Util_TableCreateFromString(geometry_interpolator_pars);
if (interp_par_table_handle < 0)
   then CCTK_VWarn(-1, __LINE__, __FILE__, CCTK_THORNSTRING,
		   "bad geometry-interpolator parameter(s) \"%s\"!",
		   geometry_interpolator_pars);		/*NOTREACHED*/

//
// create the patch system and initialize the xyz derivative coefficients
//
CCTK_VInfo(CCTK_THORNSTRING, "   creating patch system");
patch_system ps(origin_x, origin_y, origin_z,
		patch_system::type_of_name(patch_system_type),
		N_ghost_points, N_overlap_points, delta_drho_dsigma,
		nominal_gfns::min_gfn, nominal_gfns::max_gfn,
		ghosted_gfns::min_gfn, ghosted_gfns::max_gfn,
		interp_handle, interp_par_table_handle);

//
// find the apparent horizon
//
if      (STRING_EQUAL(method, "horizon"))
   then {
	ps.read_ghosted_gridfn(ghosted_gfns::gfn__h, "h.dat", false);
	const int timer_handle = CCTK_TimerCreateI();
	CCTK_TimerStartI(timer_handle);
	horizon_function(ps, cgi, gii);
	CCTK_TimerStopI(timer_handle);
	CCTK_TimerPrintDataI(timer_handle, -1);
	ps.print_gridfn(nominal_gfns::gfn__H, "H.dat");
	}
else	CCTK_VWarn(-1, __LINE__, __FILE__, CCTK_THORNSTRING,
		   "unknown method=\"%s\"!",
		   method);					/*NOTREACHED*/
}