aboutsummaryrefslogtreecommitdiff
path: root/src/GeneralizedPolynomial-Uniform/InterpLocalUniform.h
blob: 4ff2e5287958b7b807f530a5f560128cf5e952de (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
/* InterpLocalUniform.h -- private stuff for this interpolator */
/* $Header$ */

#ifdef LOCALINTERP_STANDALONE_BUILD
  typedef int    CCTK_INT;
  typedef double CCTK_REAL;
#endif

typedef CCTK_REAL fp;

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

/* misc stuff that Jonathan Thornburg likes to use */
typedef int bool;
#define false	0
#define true	1

#define then	/* empty */

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

/* number of integers in the range [x,y] inclusive */
#define HOW_MANY_IN_RANGE(x,y)  ((y) - (x) + 1)

/* is the integer x even/odd? */
#define IS_EVEN(x)	(((x) & 0x1) == 0)
#define IS_ODD (x)	(((x) & 0x1) != 0)

/* round floating-point value to nearest integer */
/* ... result is expressed as floating point! */
/* ... needs <math.h> for floor() */
#define JT_ROUND(x)	floor((x) + 0.5)

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

/*
 * Cactus has some odd anomolies between real and complex datatypes
 * as to whether they're #define or typedef, and how their presence
 * should be detected at preprocessor-time.  These macros isolate this
 * and provide a clean way for the rest of the code to test for which
 * types are defined:
 */

#ifdef CCTK_REAL4
  #define HAVE_CCTK_REAL4
  #define HAVE_CCTK_COMPLEX8
#endif

#ifdef CCTK_REAL8
  #define HAVE_CCTK_REAL8
  #define HAVE_CCTK_COMPLEX16
#endif

#ifdef CCTK_REAL16
  #define HAVE_CCTK_REAL16
  #define HAVE_CCTK_COMPLEX32
#endif

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

/*
 * compile-time upper bounds for sizing arrays etc
 */

/*
 * We must have 1 <= N_dims <= MAX_N_DIMS.  Alas, there are lots of places
 * in "template.c" where code explicitly enumerates all possible N_DIMS
 * values at compile-time, so changing (eg raising) this limit requires
 * checking all preprocessor uses of N_DIMS as well as MAX_N_DIMS. :( :(
 */
#define MAX_N_DIMS	3

/* a "boundary" is the combination of a dimension and a min/max "side" */
#define MAX_N_BOUNDARIES	(2*MAX_N_DIMS)

/*
 * if molecule_family_string is a C string specifying a molecule family
 * (i.e. value associated with the "molecule_family" key in the parameter
 * table), we must have
 *	strlen(molecule_family_string) <= MAX_MOLECULE_FAMILY_STRLEN
 * n.b. exceeding this won't cause a buffer overflow, it will "just"
 *      cause the string to be truncated (and probably not recognized
 *      by the interpolator)
 */
#define MAX_MOLECULE_FAMILY_STRLEN	20

/*
 * N_MOLECULE_FAMILIES is the number of distinct molecule families
 */
enum	molecule_family
	{
	molecule_family_cube = 0,
	N_MOLECULE_FAMILIES	/* this must be the last entry in the enum */
	};

/*
 * We must have 1 <= order <= MAX_ORDER.
 */
#define MAX_ORDER	6

/*
 * We must have 0 <= smoothing <= MAX_SMOOTHING
 */
#define MAX_SMOOTHING	0

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

/*
 * other compile-time settings
 */

/* defaults for boundary_{off_centering,extrapolation}_tolerance[] */
#ifdef NOT_YET
  #define LAGRANGE_BNDRY_OFF_CNTR_TOL_DEF	999.0
  #define LAGRANGE_BNDRY_EXTRAP_TOL_DEF		1.0e-10
  #define HERMITE_BNDRY_OFF_CNTR_TOL_DEF	1.0e-10
  #define HERMITE_BNDRY_EXTRAP_TOL_DEF		0.0
#else
  #define LAGRANGE_BNDRY_OFF_CNTR_TOL_DEF	999.0
  #define LAGRANGE_BNDRY_EXTRAP_TOL_DEF		1.0e-10
  #define HERMITE_BNDRY_OFF_CNTR_TOL_DEF	999.0
  #define HERMITE_BNDRY_EXTRAP_TOL_DEF		1.0e-10
#endif

/* CCTK_VWarn() severity level for error/warning messages */
#define BUG_MSG_SEVERITY_LEVEL			0
#define ERROR_MSG_SEVERITY_LEVEL		0
#define WARNING_MSG_SEVERITY_LEVEL		1

/* CCTK_Abort() exit code for internal error (interpolator bug) aborts */
#define BUG_ABORT_CODE				42

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

/*
 * data structures used in multiple files
 */

/* number of real "parts" in a complex number */
#define COMPLEX_N_PARTS 2

struct	error_flags
	{
	int error_pt;
	int error_ibndry;
	int error_axis;
	int error_direction;
	};

struct	molecule_structure_flags
	{
	bool MSS_is_fn_of_interp_coords;
	bool MSS_is_fn_of_which_operation;
	bool MSS_is_fn_of_input_array_values;
	bool Jacobian_is_fn_of_input_array_values;
	};

struct	molecule_min_max_m_info
	{
	CCTK_INT molecule_min_m[MAX_N_DIMS];
	CCTK_INT molecule_max_m[MAX_N_DIMS];
	};

struct	Jacobian_info
	{
	CCTK_REAL** Jacobian_pointer;
	CCTK_INT*   Jacobian_offset;
	CCTK_INT Jacobian_interp_point_stride;
	CCTK_INT Jacobian_m_strides[MAX_N_DIMS];
	CCTK_INT Jacobian_part_stride;
	};

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

/*
 * error codes for LocalInterp_molecule_posn()
 */

/* x < minimum allowable x in grid */
#define MOLECULE_POSN_ERROR_X_LT_MIN	(-1)

/* x > maximum allowable x in grid */
#define MOLECULE_POSN_ERROR_X_GT_MAX	(-2)

/* grid is smaller than molecule */
#define MOLECULE_POSN_ERROR_GRID_TINY	(-3)

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

/*
 * prototypes for functions visible from multiple files
 */

/* functions in InterpLocalUniform.c */
int LocalInterp_ILU_Lagrange_TP(int N_dims,
				int param_table_handle,
				/***** coordinate system *****/
				const CCTK_REAL coord_origin[],
				const CCTK_REAL coord_delta[],
				/***** interpolation points *****/
				int N_interp_points,
				int interp_coords_type_code,
				const void *const interp_coords[],
				/***** input arrays *****/
				int N_input_arrays,
				const CCTK_INT input_array_dims[],
				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[]);
int LocalInterp_ILU_Lagrange_MD(int N_dims,
				int param_table_handle,
				/***** coordinate system *****/
				const CCTK_REAL coord_origin[],
				const CCTK_REAL coord_delta[],
				/***** interpolation points *****/
				int N_interp_points,
				int interp_coords_type_code,
				const void *const interp_coords[],
				/***** input arrays *****/
				int N_input_arrays,
				const CCTK_INT input_array_dims[],
				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[]);
int LocalInterp_ILU_Hermite(int N_dims,
			    int param_table_handle,
			    /***** coordinate system *****/
			    const CCTK_REAL coord_origin[],
			    const CCTK_REAL coord_delta[],
			    /***** interpolation points *****/
			    int N_interp_points,
			    int interp_coords_type_code,
			    const void *const interp_coords[],
			    /***** input arrays *****/
			    int N_input_arrays,
			    const CCTK_INT input_array_dims[],
			    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[]);

/* functions in molecule_posn.c */
int LocalInterp_molecule_posn(fp grid_origin, fp grid_delta,
			      int grid_i_min, int grid_i_max,
			      int molecule_size,
			      fp boundary_off_centering_tolerance_min,
			      fp boundary_off_centering_tolerance_max,
			      fp boundary_extrapolation_tolerance_min,
			      fp boundary_extrapolation_tolerance_max,
			      fp x,
			      int* i_center, fp* x_rel);

/* functions in util.c */
int LocalInterp_decode_N_parts(int type_code);