aboutsummaryrefslogtreecommitdiff
path: root/src/periodic.c
blob: 7339d9ed531df490ffbb6ef47cf0b45e2d676fc6 (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
/* $Header$ */

#include <assert.h>
#include <stdlib.h>
#include <stdio.h>
#include "cctk.h"
#include "cctk_Parameters.h"
#include "Slab.h"
#include "periodic.h"

static const char * restrict const rcsid = "$Header$";
CCTK_FILEVERSION(TAT_Periodic_periodic_c);



int
BndPeriodicVI (cGH const * restrict const cctkGH,
               int const * restrict const stencil,
               int const vi)
{
  DECLARE_CCTK_PARAMETERS;
  
  cGroup group;
  cGroupDynamicData data;
  void * restrict varptr;
  struct xferinfo * restrict xferinfo;
  int do_periodic[3];
  int gi;
  int dir;
  int d, f;
  int ierr;
  
  /* Check arguments */
  assert (cctkGH);
  assert (stencil);
  assert (vi>=0 && vi<CCTK_NumVars());
  
  /* Get and check group info */
  gi = CCTK_GroupIndexFromVarI (vi);
  assert (gi>=0 && gi<CCTK_NumGroups());
  
  ierr = CCTK_GroupData (gi, &group);
  assert (!ierr);
  assert (group.grouptype == CCTK_GF);
  assert (group.disttype == CCTK_DISTRIB_DEFAULT);
  assert (group.stagtype == 0);
  
  ierr = CCTK_GroupDynamicData (cctkGH, gi, &data);
  assert (!ierr);
  
  varptr = CCTK_VarDataPtrI (cctkGH, 0, vi);
  assert (varptr);
  
  /* Condition periodicity information */
  do_periodic[0] = periodic || periodic_x;
  do_periodic[1] = periodic || periodic_y;
  do_periodic[2] = periodic || periodic_z;
  
  printf ("BndPeriodicVI var=%s\n", CCTK_FullName(vi));
  for (d=0; d<3; ++d) {
    printf ("  d=%d per=%d ste=%d\n", d, do_periodic[d], stencil[d]);
  }
  
  /* Allocate slab transfer description */
  xferinfo = malloc(group.dim * sizeof *xferinfo);
  assert (xferinfo);
  
  for (dir=0; dir<group.dim; ++dir) {
    if (dir<3 && do_periodic[dir]) {
      
      /* Loop over faces */
      for (f=0; f<2; ++f) {
        
        /* Fill in slab transfer description */
        for (d=0; d<group.dim; ++d) {
          xferinfo[d].src.gsh         = data.gsh        [d];
          xferinfo[d].src.lbnd        = data.lbnd       [d];
          xferinfo[d].src.lsh         = data.lsh        [d];
          xferinfo[d].src.lbbox       = data.bbox       [2*d  ];
          xferinfo[d].src.ubbox       = data.bbox       [2*d+1];
          xferinfo[d].src.nghostzones = data.nghostzones[d];
          xferinfo[d].src.off         = 0;
          xferinfo[d].src.str         = 1;
          xferinfo[d].src.len         = data.lsh        [d];
          
          xferinfo[d].dst.gsh         = data.gsh        [d];
          xferinfo[d].dst.lbnd        = data.lbnd       [d];
          xferinfo[d].dst.lsh         = data.lsh        [d];
          xferinfo[d].dst.lbbox       = data.bbox       [2*d  ];
          xferinfo[d].dst.ubbox       = data.bbox       [2*d+1];
          xferinfo[d].dst.nghostzones = data.nghostzones[d];
          xferinfo[d].dst.off         = 0;
          xferinfo[d].dst.str         = 1;
          xferinfo[d].dst.len         = data.lsh        [d];
          
          xferinfo[d].xpose = d;
          xferinfo[d].flip  = 0;
        }
        
        if (f==0) {
          /* Fill in lower face */
          xferinfo[dir].src.off = data.gsh[dir] - 2 * stencil[dir];
          xferinfo[dir].src.len = stencil[dir];
          xferinfo[dir].dst.len = stencil[dir];
        } else {
          /* Fill in upper face */
          xferinfo[dir].src.off = stencil[dir];
          xferinfo[dir].src.len = stencil[dir];
          xferinfo[dir].dst.off = data.gsh[dir] - stencil[dir];
          xferinfo[dir].dst.len = stencil[dir];
        }
        
        ierr = Slab_Transfer (cctkGH, group.dim, xferinfo, -1,
                              group.vartype, varptr, group.vartype, varptr);
        assert (!ierr);
        
      } /* for f */
      
    } /* if dir */
  } /* for dir */
  
  free (xferinfo);
  
  return 0;
}



int
BndPeriodicVN (cGH const * restrict const cctkGH,
               int const * restrict const stencil,
               char const * restrict const vn)
{
  int const vi = CCTK_VarIndex (vn);
  assert (vi>=0 && vi<CCTK_NumVars());
  return BndPeriodicVI (cctkGH, stencil, vi);
}



int
BndPeriodicGI (cGH const * restrict const cctkGH,
               int const * restrict const stencil,
               int const gi)
{
  int vi1, vi2;
  int vi;
  int ierr;
  assert (gi>=0 && gi<CCTK_NumGroups());
  vi1 = CCTK_FirstVarIndexI(gi);
  assert (vi1>=0 && vi1<CCTK_NumVars());
  vi2 = vi1 + CCTK_NumVarsInGroupI(gi);
  assert (vi2>=vi1 && vi2<CCTK_NumVars());
  assert (vi>=0 && vi<CCTK_NumVars());
  for (vi=vi1; vi<vi2; ++vi) {
    ierr = BndPeriodicVI (cctkGH, stencil, vi);
    if (ierr) return ierr;
  }
  return 0;
}



int
BndPeriodicGN (cGH const * restrict const cctkGH,
               int const * restrict const stencil,
               char const * restrict const gn)
{
  int const gi = CCTK_GroupIndex (gn);
  assert (gi>=0 && gi<CCTK_NumGroups());
  return BndPeriodicGI (cctkGH, stencil, gi);
}



void
Periodic_ApplyBC (cGH const * restrict const cctkGH)
{
  int nvars;
  CCTK_INT * restrict indices;
  CCTK_INT * restrict faces;
  CCTK_INT * restrict widths;
  CCTK_INT * restrict tables;
  int vi;
  int dim;
  int * restrict stencil;
  int i;
  int d;
  int ierr;
  
  assert (cctkGH);
  
  nvars = Boundary_SelectedGVs (cctkGH, 0, 0, 0, 0, 0, 0);
  assert (nvars>=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 (i=0; i<nvars; ++i) {
    vi = indices[i];
    assert (vi>=0 && vi<CCTK_NumVars());
    
    assert (widths[i] >= 0);
    
    dim = CCTK_GroupDimFromVarI (vi);
    assert (dim>=0);
    
    stencil = malloc (dim * sizeof *stencil);
    assert (stencil);
#if 0
    for (d=0; d<dim; ++d) {
      stencil[d] = widths[i];
    }
#endif
    ierr = CCTK_GroupnghostzonesVI (cctkGH, dim, stencil, vi);
    assert (!ierr);
    
    ierr = BndPeriodicVI (cctkGH, stencil, vi);
    assert (!ierr);
    
    free (stencil);
  }
  
  free (indices);
  free (faces);
  free (widths);
  free (tables);
}