/* $Header$ */ #include #include #include #include #include "cctk.h" #include "cctk_Arguments.h" #include "cctk_Parameters.h" #include "util_ErrorCodes.h" #include "util_Table.h" #include "Slab.h" static int convert_index (int const step, int const index, int const * restrict const alldirs, int * restrict const parity) { int srcindex; assert (index>=0 && index<3); assert (alldirs); assert (alldirs[0]>=0 && alldirs[0]<3); assert (alldirs[1]>=0 && alldirs[1]<3); assert (alldirs[1] != alldirs[0]); assert (parity); assert (abs(*parity) == 1); switch (step) { case 0: /* x face: counterclockwise, 90 degrees */ if (index == alldirs[0]) srcindex = alldirs[1], *parity *= -1; else if (index == alldirs[1]) srcindex = alldirs[0], *parity *= +1; else srcindex = index , *parity *= +1; break; case 1: /* y face: clockwise, 90 degrees */ if (index == alldirs[0]) srcindex = alldirs[1], *parity *= +1; else if (index == alldirs[1]) srcindex = alldirs[0], *parity *= -1; else srcindex = index , *parity *= +1; break; case 2: /* xy edge: 180 degrees */ if (index == alldirs[0]) srcindex = alldirs[0], *parity *= -1; else if (index == alldirs[1]) srcindex = alldirs[1], *parity *= -1; else srcindex = index , *parity *= +1; break; default: assert (0); } return srcindex; } int BndRot90VI (cGH const * restrict const cctkGH, int const nvars, CCTK_INT const * restrict const vis) { DECLARE_CCTK_ARGUMENTS; DECLARE_CCTK_PARAMETERS; int * restrict gis; cGroup group; cGroupDynamicData data; char * restrict fullname; void * restrict * restrict varptrs; int * restrict vartypes; int firstvar; char tensortypealias[100]; int numvars; int index; int srcvi; void const * restrict * restrict srcptrs; int * restrict parities; int global_bbox[6]; int global_lbnd[3], global_ubnd[3]; int fake_bbox[6]; CCTK_REAL x0[3], dx[3]; CCTK_REAL origin[3], dorigin[3]; int avoid_origin[3]; int offset[3]; /* offset 0..1 due to avoid_origin */ struct xferinfo * restrict xferinfo; int have_global_bbox, have_local_bbox; int var; int alldirs[2]; int step; int ndirs; /* 1 for face, 2 for edge */ int dir[3]; /* current face or edge */ int otherdir[3]; /* the other direction of the rotation */ int q; /* 0..ndirs-1 */ int d; /* 0..group.dim-1 */ int ierr; /* Check arguments */ assert (cctkGH); assert (nvars==0 || vis); for (var=0; var=0 && vis[var]0); gis = malloc (nvars * sizeof *gis); assert (nvars==0 || gis); varptrs = malloc (nvars * sizeof *varptrs); assert (nvars==0 || varptrs); vartypes = malloc (nvars * sizeof *vartypes); assert (nvars==0 || vartypes); for (var=0; var=0 && gis[var]= 0); } for (d=0; dcctk_bbox[d]; for (d=0; d<3; ++d) local[6+d] = -cctkGH->cctk_lbnd[d]; for (d=0; d<3; ++d) local[9+d] = cctkGH->cctk_ubnd[d]; ierr = CCTK_ReduceLocArrayToArray1D (cctkGH, -1, max_handle, local, global, 12, CCTK_VARIABLE_INT); for (d=0; d<6; ++d) global_bbox[d] = global[ d]; for (d=0; d<3; ++d) global_lbnd[d] = -global[6+d]; for (d=0; d<3; ++d) global_ubnd[d] = global[9+d]; for (d=0; d<3; ++d) { fake_bbox[2*d ] = data.lbnd[d] == global_lbnd[d]; fake_bbox[2*d+1] = data.ubnd[d] == global_ubnd[d]; } } /* directions */ alldirs[0] = 0; alldirs[1] = 1; /* Find grid point that corresponds to the origin */ for (q=0; q<2; ++q) { d = alldirs[q]; /* x0 + dx * origin == 0 */ origin[d] = - x0[d] / dx[d]; dorigin[d] = origin[d] - floor(origin[d]); if (fabs(dorigin[d]) < 1.0e-6 || fabs(dorigin[d] - 1.0) < 1.0e-6) { avoid_origin[d] = 0; } else if (fabs(dorigin[d] - 0.5) < 1.0e-6) { avoid_origin[d] = 1; } else { CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING, "The coordinate origin in the %c-direction falls neither onto a grid point nor into the middle between two grid points.", "xyz"[d]); } offset[d] = avoid_origin[d] ? 0 : 1; } parities = malloc (nvars * sizeof *parities); assert (nvars==0 || parities); srcptrs = malloc (nvars * sizeof *srcptrs); assert (nvars==0 || srcptrs); for (step=0; step<3; ++step) { switch (step) { case 0: /* x face */ ndirs = 1; dir[0] = 0; otherdir[0] = 1; break; case 1: /* y face */ ndirs = 1; dir[0] = 1; otherdir[0] = 0; break; case 2: /* xy edge */ ndirs = 2; dir[0] = 0; dir[1] = 1; otherdir[0] = 1; otherdir[1] = 0; break; default: assert (0); } for (var=0; var0); firstvar = CCTK_FirstVarIndexI(gis[var]); assert (firstvar>=0); index = vis[var] - firstvar; assert (index>=0 && index=0); ierr = Util_TableGetString (table, sizeof tensortypealias, tensortypealias, "tensortypealias"); if (ierr == UTIL_ERROR_TABLE_NO_SUCH_KEY) { /* assume a scalar */ if (numvars != 1) { static int * restrict didwarn = 0; if (! didwarn) { didwarn = calloc (CCTK_NumGroups(), sizeof *didwarn); } if (! didwarn[gis[var]]) { char * groupname = CCTK_GroupName(gis[var]); assert (groupname); CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING, "Group \"%s\" has no tensor type and contains more than one element -- treating these as \"scalar\"", groupname); free (groupname); didwarn[gis[var]] = 1; } } strcpy (tensortypealias, "scalar"); } else if (ierr<0) { char * groupname = CCTK_GroupName(gis[var]); assert (groupname); CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING, "Error in tensor type alias declaration for group \"%s\"", groupname); free (groupname); } if (CCTK_EQUALS (tensortypealias, "scalar")) { /* scalar */ } else if (CCTK_EQUALS (tensortypealias, "4scalar")) { /* 4-scalar */ } else if (CCTK_EQUALS (tensortypealias, "u") || CCTK_EQUALS (tensortypealias, "d")) { /* vector */ assert (numvars == 3); } else if (CCTK_EQUALS (tensortypealias, "4u") || CCTK_EQUALS (tensortypealias, "4d")) { /* 4-vector */ assert (numvars == 4); } else if (CCTK_EQUALS (tensortypealias, "uu_sym") || CCTK_EQUALS (tensortypealias, "dd_sym")) { /* symmetric tensor */ assert (numvars == 6); } else if (CCTK_EQUALS (tensortypealias, "4uu_sym") || CCTK_EQUALS (tensortypealias, "4dd_sym")) { /* symmetric 4-tensor */ assert (numvars == 10); } else { char * groupname = CCTK_GroupName(gis[var]); assert (groupname); CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING, "Illegal tensor type alias for group \"%s\"", groupname); free (groupname); } } parities[var] = +1; if (CCTK_EQUALS (tensortypealias, "scalar")) { /* scalar */ srcvi = vis[var]; /* do nothing */ } else if (CCTK_EQUALS (tensortypealias, "4scalar")) { /* 4-scalar */ srcvi = vis[var]; /* do nothing */ } else if (CCTK_EQUALS (tensortypealias, "u") || CCTK_EQUALS (tensortypealias, "d")) { /* vector */ int srcindex; assert (index>=0 && index<3); srcindex = convert_index (step, index, alldirs, &parities[var]); srcvi = firstvar + srcindex; } else if (CCTK_EQUALS (tensortypealias, "4u") || CCTK_EQUALS (tensortypealias, "4d")) { /* 4-vector */ assert (index>=0 && index<4); if (index==0) { /* temporal component */ srcvi = firstvar; } else { /* spatial components */ int srcindex; int const off = 1; srcindex = off + convert_index (step, index-off, alldirs, &parities[var]); srcvi = firstvar + srcindex; } } else if (CCTK_EQUALS (tensortypealias, "uu_sym") || CCTK_EQUALS (tensortypealias, "dd_sym")) { /* symmetric tensor */ int index1, index2; int srcindex1, srcindex2; int srcindex; assert (index>=0 && index<6); { int const expand1[6] = { 0,0,0,1,1,2 }; int const expand2[6] = { 0,1,2,1,2,2 }; index1 = expand1[index]; index2 = expand2[index]; } srcindex1 = convert_index (step, index1, alldirs, &parities[var]); srcindex2 = convert_index (step, index2, alldirs, &parities[var]); { int const compact[3][3] = { { 0,1,2 }, { 1,3,4 }, { 2,4,5 } }; srcindex = compact[srcindex1][srcindex2]; } srcvi = firstvar + srcindex; } else if (CCTK_EQUALS (tensortypealias, "4uu_sym") || CCTK_EQUALS (tensortypealias, "4dd_sym")) { /* symmetric 4-tensor */ assert (index>=0 && index<10); if (index==0) { /* temporal-temporal component */ srcvi = firstvar; } else if (index<4) { /* temporal-spatial components */ int srcindex; int const off = 1; srcindex = off + convert_index (step, index-off, alldirs, &parities[var]); srcvi = firstvar + srcindex; } else { /* spatial-spatial components */ int index1, index2; int srcindex1, srcindex2; int srcindex; int const off = 4; { int const expand1[6] = { 0,0,0,1,1,2 }; int const expand2[6] = { 0,1,2,1,2,2 }; index1 = expand1[index-off]; index2 = expand2[index-off]; } srcindex1 = convert_index (step, index1, alldirs, &parities[var]); srcindex2 = convert_index (step, index2, alldirs, &parities[var]); { int const compact[3][3] = { { 0,1,2 }, { 1,3,4 }, { 2,4,5 } }; srcindex = off + compact[srcindex1][srcindex2]; } srcvi = firstvar + srcindex; } } else { assert (0); } srcptrs[var] = CCTK_VarDataPtrI (cctkGH, 0, srcvi); assert (srcptrs[var]); } /* for var */ have_global_bbox = 1; for (q=0; q= 0); if (data.gsh[otherdir[q]] < data.nghostzones[dir[q]] + data.nghostzones[otherdir[q]] + offset[otherdir[q]]) { assert (nvars > 0); CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING, "The group \"%s\" has in the %c-direction only %d grid points. " "This is not large enough for a 90 degree rotating symmetry boundary that is %d grid points wide in the %c-direction. " "The group needs to have at least %d grid points in the %c-direction.", CCTK_GroupNameFromVarI(vis[0]), "xyz"[dir[q]], data.gsh[dir[q]], data.nghostzones[dir[q]], "xyz"[otherdir[q]], data.nghostzones[dir[q]] + data.nghostzones[otherdir[q]] + offset[otherdir[q]], "xyz"[dir[q]]); } } /* Allocate slab transfer description */ xferinfo = malloc(group.dim * sizeof *xferinfo); assert (xferinfo); /* Fill in slab transfer description */ for (d=0; d= cctk_lsh[d]) imax[d] = cctk_lsh[d]; } assert (group.dim == 3); assert (group.vartype == CCTK_VARIABLE_REAL); for (k=imin[0]; k=0); indices = malloc (nvars * sizeof *indices); assert (indices); faces = malloc (nvars * sizeof *faces); assert (faces); widths = malloc (nvars * sizeof *widths); assert (widths); tables = malloc (nvars * sizeof *tables); assert (tables); ierr = Boundary_SelectedGVs (cctkGH, nvars, indices, faces, widths, tables, 0); assert (ierr == nvars); for (var=0; var=0 && indices[var]= 0); } ierr = BndRot90VI (cctkGH, nvars, indices); assert (!ierr); free (indices); free (faces); free (widths); free (tables); }