summaryrefslogtreecommitdiff
path: root/src/main/Stagger.c
blob: 67aa8c0fb19f3c1383ba81418d11e1b7291b5fd1 (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
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
 /*@@
   @file      Stagger.c
   @date      Thu Jan 27 15:32:28 2000
   @author    Gerd Lanfermann
   @desc 
   Stuff to deal with staggering.
   @enddesc 
   @version $Header$
 @@*/

#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <ctype.h>

#include "cctk_Flesh.h"
#include "cctk_FortranString.h"
#include "cctk_Groups.h"
#include "cctk_Types.h"
#include "cctk_WarnLevel.h"

#include "cctk_Stagger.h"
#include "cctki_Stagger.h"

static const char *rcsid = "$Header$";

CCTK_FILEVERSION(main_Stagger_c);

/********************************************************************
 *********************     Fortran Wrappers   ***********************
 ********************************************************************/
void CCTK_FCALL CCTK_FNAME (CCTK_GroupStaggerIndexGI)
                           (int *stagcode, int *gindex);
void CCTK_FCALL CCTK_FNAME (CCTK_GroupStaggerIndexGN)
                           (int *scode, ONE_FORTSTRING_ARG);
void CCTK_FCALL CCTK_FNAME (CCTK_StaggerIndex)
                           (int *scode, ONE_FORTSTRING_ARG);
void CCTK_FCALL CCTK_FNAME (CCTK_StaggerDirIndex)
                           (int *dsi, int *dir, int *gsi);
void CCTK_FCALL CCTK_FNAME (CCTK_GroupStaggerDirArray)
                           (int *ierr, int *dindex, int *dim, int *gsc);
void CCTK_FCALL CCTK_FNAME (CCTK_GroupStaggerDirArrayGI)
                           (int *ierr, int *dindex, int *dim, int *gi);
void CCTK_FCALL CCTK_FNAME (CCTK_StaggerDirName)
                           (int *dsc, int *dir, ONE_FORTSTRING_ARG);

/********************************************************************
 *********************     External Routines   **********************
 ********************************************************************/


 /*@@
   @routine    CCTK_GroupStaggerIndexGI
   @date       
   @author     
   @desc       retuns the stagger code for a given group index

   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

int CCTK_GroupStaggerIndexGI(int gindex)
{
  cGroup group;
  int sc;
  CCTK_GroupData(gindex, &group);
  sc = group.stagtype;
  return(sc);
}

void CCTK_FCALL CCTK_FNAME (CCTK_GroupStaggerIndexGI)
                           (int *stagcode, int *gindex) 
{
  *stagcode = CCTK_GroupStaggerIndexGI(*gindex);
}


 /*@@
   @routine    CCTK_GroupStaggerIndexGN
   @date       
   @author     Gerd Lanfermann
   @desc       returns the stagger index for a given group name
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/
 
int CCTK_GroupStaggerIndexGN(const char *gname) 
{
  int gindex;
  gindex = CCTK_GroupIndex(gname);
  return(CCTK_GroupStaggerIndexGI(gindex));
}

void CCTK_FCALL CCTK_FNAME (CCTK_GroupStaggerIndexGN)
                           (int *scode, ONE_FORTSTRING_ARG)
{
  ONE_FORTSTRING_CREATE(gname)
  int gindex;
  gindex = CCTK_GroupIndex(gname);
  *scode = CCTK_GroupStaggerIndexGI(gindex);
  free(gname);
}
  

 /*@@
   @routine    CCTK_StaggerIndex
   @date       
   @author     Gerd Lanfermann
   @desc       returns the stagger index for a given stagger name

   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/
 
int CCTK_StaggerIndex(const char *stype) 
{
  int i,scode,base,dim,m;

  base = 1;
  scode= 0;
  dim  = strlen(stype);

  for (i=0;i<dim;i++) 
  {

    switch (toupper(stype[i]))
    {
      case 'M':m=0; break;
      case 'C':m=1; break;
      case 'P':m=2; break;
      default:
        CCTK_VWarn(1,__LINE__,__FILE__,"Cactus",
                   "CCTK_StaggerIndex: Unknown stagger type %s",stype);
        return(-1);
    }
    scode+= m*base;
    base  = CCTK_NUM_STAGGER * base;
  }
  return(scode);
}

void CCTK_FCALL CCTK_FNAME (CCTK_StaggerIndex)
                           (int *scode, ONE_FORTSTRING_ARG)
{
  ONE_FORTSTRING_CREATE(sname);
  *scode = CCTK_StaggerIndex(sname);
  free(sname);
}


/*@@
   @routine    CCTK_StaggerDirIndex
   @date       
   @author     Gerd Lanfermann
   @desc       returns the stagger index in a direction <dir>
               when given the staggerindex <sc>.

   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/
 
int CCTK_StaggerDirIndex(int dir, int si) 
{
  int val,b,dsi=0;
  static int hash[4],hashed=0;

  if (hashed==0) 
  {
    hash[0]= 1;
    hash[1]= 3;
    hash[2]= 9;
    hash[3]=27;
    hashed = 1;
  }

  for (b=CCTK_NSTAG;b>=0;b--) 
  {
    val = (int)(si / hash[b]);
    si  = si % hash[b];
    if (dir==b) 
    {
      dsi = val;
      break;
    }
  }
  return(dsi);
}



void CCTK_FCALL CCTK_FNAME (CCTK_StaggerDirIndex)
                           (int *dsi, int *dir, int *gsi) 
{
  /* accept fortran indexing [1..]: decrease the directional index
     for the call to the C routine.  */
  *dsi  = CCTK_StaggerDirIndex((*dir)-1, *gsi);
} 



/*@@
   @routine    CCTK_StaggerDirIndexArray
   @date       
   @author     Gerd Lanfermann
   @desc       returns the stagger index for all direction in 
               an array <dindex> of size <dim> when given the staggerindex <sc>.

   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/
 
int CCTK_StaggerDirArray(int *dindex , int dim, int sindex) 
{
  int val,b;
  static int hash[4],hashed=0;

  if (hashed==0) 
  {
    hash[0]= 1;
    hash[1]= 3;
    hash[2]= 9;
    hash[3]=27;
    hashed = 1;
  }

  if (dim>4) 
  {
    CCTK_VWarn(1,__LINE__,__FILE__,"Cactus", 
              "CCTK_StaggerDirArray: Dimension %d exceeds maximum of 4",dim);
    return(-1);
  }

  for (b=CCTK_NSTAG;b>=0;b--) 
  {
    val       = (int)(sindex / hash[b]);
    sindex    = sindex % hash[b];
    if (b<dim) dindex[b] = val;
  }
  return(0);
}


void CCTK_FCALL CCTK_FNAME (CCTK_GroupStaggerDirArray)
                           (int *ierr, int *dindex, int *dim, int *gsc) 
{
  /* accept fortran indexing [1..]: decrease the directional index
     for the call to the C routine.  */
  *ierr = CCTK_StaggerDirArray(dindex, *dim, *gsc);
} 

 /*@@
   @routine    CCTK_GroupStaggerDirArrayGI
   @date       
   @author     Gerd Lanfermann
   @desc       returns the stagger index for all direction in 
               an array <dindex> of size <dim> when given the group
               index <gi>

   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/
 
int CCTK_GroupStaggerDirArrayGI(int *dindex, int dim, int gi) 
{
  int si,ierr;
  si  = CCTK_GroupStaggerIndexGI(gi);
  ierr= CCTK_StaggerDirArray(dindex, dim, si);
  return ierr;
}

void CCTK_FCALL CCTK_FNAME (CCTK_GroupStaggerDirArrayGI)
                           (int *ierr, int *dindex, int *dim, int *gi) 
{
  *ierr = CCTK_GroupStaggerDirArrayGI(dindex, *dim, *gi);
} 


 /*@@
   @routine    CCTK_StaggerDirName
   @date       
   @author     Gerd Lanfermann
   @desc       returns the directional staggering in direction <dir> 
               for a given stagger name <stype>.

   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/
 
int CCTK_StaggerDirName(int dir, const char *stype) 
{
  int scode;
  char hs[7]="MMMMMM";

  sprintf(hs,"%s",stype);

  if (dir> (int) strlen(hs)) 
  {
    CCTK_VWarn(1,__LINE__,__FILE__,"Cactus",
              "CCTK_StaggerDirName: Stagger name too short for direction %d",
              dir);
  }

  switch (toupper(hs[dir]))
  {
    case 'M': scode = 0; break;
    case 'C': scode = 1; break;
    case 'P': scode = 2; break;
    default:
      CCTK_VWarn(1,__LINE__,__FILE__,"Cactus",
                "CCTK_StaggerDirName: Unknown stagger type %s",hs);
      return(-1);
  }
  return(scode);
}

void CCTK_FCALL CCTK_FNAME (CCTK_StaggerDirName)
                           (int *dsc, int *dir, ONE_FORTSTRING_ARG) 
{
  ONE_FORTSTRING_CREATE(sname);

  *dsc = CCTK_StaggerDirName((*dir)-1,sname);

  free(sname);
}






 /*@@
   @routine    CCTKi_ParseStaggerString
   @date       
   @author     Gerd Lanfermann
   @desc       returns the stagger index for a string. Similar routines 
               as CCTK_StaggerIndexName, but does more error checking since 
               it is called during Group setup and if things go wrong, the used
               has a better idea where he specified wrong settings.

   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/
 
int CCTKi_ParseStaggerString(int dim,
                             const char *imp, 
                             const char *gname,
                             const char *stype) 
{
  int i,m;
  int base  = 1;
  int scode = 0;
  char hs[11];

  if (dim>10) 
  {
    CCTK_VWarn(0,__LINE__,__FILE__,"Cactus",
               "CCTKi_ParseStaggerString: Dimension %d exceeds maximum of 10",
               dim);
  }

  /* change possible SHORTCUTS into the official notation*/
  if (CCTK_Equals(stype,"NONE"))
  {
    sprintf(hs,"MMMMMMMMMM");
  }
  else if (CCTK_Equals(stype,"CELL")) 
  {
    sprintf(hs,"CCCCCCCCCC");
  }
  else 
  {  
    if ((int) strlen(stype)!=dim) 
    {  
      CCTK_VWarn(1,__LINE__,__FILE__,"Cactus",
                "CCTKi_ParseStaggerString: Staggering %s for %s unequal to group dimension %d",
                stype,gname,dim);
    }

    sprintf(hs,"%s",stype);
  }
  
  for (i=0;i<dim;i++) 
  {
    switch (toupper(hs[i]))
    {
      case 'M':m=0; break;
      case 'C':m=1; break;
      case 'P':m=2; break;
      default:
        CCTK_VWarn(1,__LINE__,__FILE__,"Cactus",
                   "CCTKi_ParseStaggerString: Unknown stagger type %s for %s::%s",
                   stype,imp,gname);
        return(-1);
    }
    scode+= m*base;
    base  = 3 * base;
  }

  return scode;
}