diff options
author | schnetter <schnetter@6ca9aeac-0e4f-0410-a746-fe4df63e9d0c> | 2005-08-11 16:30:21 +0000 |
---|---|---|
committer | schnetter <schnetter@6ca9aeac-0e4f-0410-a746-fe4df63e9d0c> | 2005-08-11 16:30:21 +0000 |
commit | 3b17f6feb768ef00a03bfb1ac12c39eb46b5f238 (patch) | |
tree | 5215b074f346937899341006d12842c59a468632 | |
parent | 835882819951638883a2f47036a17bd2bbb5c092 (diff) |
Support taking space or time derivatives while interpolating. So far
implemented only for interpolating into 1D arrays.
git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/InterpToArray/trunk@5 6ca9aeac-0e4f-0410-a746-fe4df63e9d0c
-rw-r--r-- | param.ccl | 10 | ||||
-rw-r--r-- | src/interp.c | 35 |
2 files changed, 44 insertions, 1 deletions
@@ -59,6 +59,16 @@ STRING array1d_vars[100] "Names of the grid functions that should be interpolate "^[A-Za-z][A-Za-z0-9_]*[:][:][A-Za-z][A-Za-z0-9_]*$" :: "grid function name" } "" +INT array1d_spacederivs[100] "Space derivative orders for each grid function" +{ + 0:* :: "" +} 0 + +INT array1d_timederivs[100] "Time derivative order for each grid function" +{ + 0:* :: "" +} 0 + REAL array1d_x0 "Origin" STEERABLE=always { *:* :: "" diff --git a/src/interp.c b/src/interp.c index badaead..91a84b8 100644 --- a/src/interp.c +++ b/src/interp.c @@ -28,6 +28,8 @@ InterpToArray (CCTK_ARGUMENTS) CCTK_INT * restrict inputs; CCTK_INT * restrict output_types; CCTK_POINTER * restrict outputs; + CCTK_INT * restrict operation_codes; + CCTK_INT * restrict time_deriv_order; int nvars; int npoints; @@ -97,7 +99,7 @@ InterpToArray (CCTK_ARGUMENTS) for (n=0; n<nvars; ++n) { outputs[n] = &scalars[npoints * n]; } - + ierr = CCTK_InterpGridArrays (cctkGH, 3, interpolator, options_table, coord_handle, npoints, CCTK_VARIABLE_REAL, coords, @@ -165,19 +167,50 @@ InterpToArray (CCTK_ARGUMENTS) for (n=0; n<nvars; ++n) { outputs[n] = &arrays1d[npoints * n]; } + + operation_codes = malloc (nvars * sizeof * operation_codes); + assert (operation_codes); + + for (n=0; n<nvars; ++n) { + operation_codes[n] = array1d_spacederivs[n]; + } + + ierr = Util_TableSetIntArray + (options_table, nvars, operation_codes, "operation_codes"); + assert (! ierr); + + time_deriv_order = malloc (nvars * sizeof * time_deriv_order); + assert (time_deriv_order); + + for (n=0; n<nvars; ++n) { + time_deriv_order[n] = array1d_timederivs[n]; + } + ierr = Util_TableSetIntArray + (options_table, nvars, time_deriv_order, "time_deriv_order"); + assert (! ierr); + ierr = CCTK_InterpGridArrays (cctkGH, 3, interpolator, options_table, coord_handle, npoints, CCTK_VARIABLE_REAL, coords, nvars, inputs, nvars, output_types, outputs); + assert (! ierr); + + ierr = Util_TableDeleteKey (options_table, "operation_codes"); + assert (! ierr); + ierr = Util_TableDeleteKey (options_table, "time_deriv_order"); + assert (! ierr); + free (coordsx); free (coordsy); free (coordsz); free (inputs); free (output_types); free (outputs); + free (operation_codes); + free (time_deriv_order); } } |