aboutsummaryrefslogtreecommitdiff
path: root/src/PughUtils.c
blob: 621cc48d0d6b757e26b913ec156dc2f093ea786e (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
 /*@@
   @file      PUGH_utils.c
   @date      Sunday 12th September 1999
   @author    Gabrielle Allen
   @desc
   @version   $Id$
   @enddesc
 @@*/

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

#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"

#include "pugh.h"

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

CCTK_FILEVERSION(CactusPUGH_PUGH_PughUtils_c);


/********************************************************************
 ********************    Static Variables   *************************
 ********************************************************************/

/********************************************************************
 ********************    Internal Routines   ************************
 ********************************************************************/


/********************************************************************
 ********************    External Routines   ************************
 ********************************************************************/
int PUGH_QueryGroupStorage (const cGH *GH, int group, const char *groupname);
void PUGH_Report(cGH *GH);
void PUGH_PrintStorageReport (cGH *GH);
void PUGHi_PrintStorageReport (void);
void PUGH_PrintFinalStorageReport (cGH *GH);
void PUGH_PrintStorage(const cGH *GH);
void CCTK_FCALL CCTK_FNAME (PUGH_PrintStorage)
                           (const cGH **GH);

 /*@@
   @routine    PUGH_Topology
   @date       Sunday 28 October 2001
   @author     Gabrielle Allen
   @desc
               Return the processor topology for grid functions
               for a given dimension
   @enddesc
@@*/
const int *PUGH_Topology (const cGH *GH, int dim)
{
  pGH *pughGH;


  if (dim <= 0)
  {
    CCTK_WARN (0, "Illegal dimension");
  }

  pughGH = PUGH_pGH (GH);

  return pughGH->Connectivity[dim-1]->nprocs;
}

 /*@@
   @routine    PUGH_Report
   @date       Sunday 12 September 1999
   @author     Gabrielle Allen
   @desc
               Report on PUGH set up
   @enddesc
@@*/
void PUGH_Report (cGH *GH)
{
  pGH *pughGH;
  int i,gi,dim;
  int *havedims;
  char *mess;
#ifdef CCTK_MPI
  int in;
#endif
  DECLARE_CCTK_PARAMETERS

  pughGH = PUGH_pGH(GH);

  havedims = (int *)calloc(GH->cctk_dim,sizeof(int));

  mess = (char *)malloc(200*sizeof(char));

#ifdef CCTK_MPI
  CCTK_VInfo(CCTK_THORNSTRING, "MPI Evolution on %d processors", pughGH->nprocs);
#else
  CCTK_VInfo(CCTK_THORNSTRING, "Single processor evolution");
#endif

  /* Report on grid decomposition for any grid functions */
  for (gi=0;gi<CCTK_NumGroups();gi++)
  {
    if (CCTK_GroupTypeI(gi) == CCTK_GF)
    {
      dim = CCTK_GroupDimI(gi);
      if (dim > 0)
      {
        havedims[dim-1]=1;
      }
    }
  }

  for (dim=0;dim<GH->cctk_dim;dim++)
  {
    if (havedims[dim])
    {
      CCTK_VInfo(CCTK_THORNSTRING, "%d-dimensional grid functions",dim+1);
      sprintf(mess,"  Size:");
      for (i=0;i<dim+1;i++)
      {
        sprintf(mess,"%s %d",mess,pughGH->GFExtras[dim]->nsize[i]);
      }
      CCTK_INFO(mess);
#ifdef CCTK_MPI
      sprintf(mess,"  Processor topology:");
      for (i=0;i<dim;i++)
      {
        sprintf(mess,"%s %d x",mess,pughGH->Connectivity[dim]->nprocs[i]);
      }
      sprintf(mess,"%s %d",mess,pughGH->Connectivity[dim]->nprocs[dim]);
      CCTK_INFO(mess);

      if (CCTK_Equals(partition, "automatic"))
      {
        sprintf(mess,"  Local load: %d   [",pughGH->GFExtras[dim]->npoints);
        for (i=0;i<dim;i++)
        {
          sprintf(mess,"%s%d x ",mess,pughGH->GFExtras[dim]->lnsize[i]);
        }
        sprintf(mess,"%s%d]",mess,pughGH->GFExtras[dim]->lnsize[i]);
        CCTK_INFO(mess);

        if (CCTK_Equals(info,"load"))
        {
          for (in=0; in<pughGH->nprocs; in++)
          {
            sprintf(mess,"  Local load on processor %d:  %d  [",
                    in,pughGH->GFExtras[dim]->rnpoints[in]);
            for (i=0;i<dim;i++)
            {
              sprintf(mess,"%s%d x ",mess,
                      pughGH->GFExtras[dim]->rnsize[in][i]);
            }
            sprintf(mess,"%s%d]",mess,pughGH->GFExtras[dim]->rnsize[in][i]);
            CCTK_INFO(mess);
          }
        }

      }
      else
      { /* manual partition */
        for (in=0; in<pughGH->nprocs; in++)
        {
          sprintf(mess,"  Local load on proc %d: %d   [",
                  in,pughGH->GFExtras[dim]->rnpoints[in]);
          for (i=0;i<dim;i++)
          {
            sprintf(mess,"%s%d x ",mess,pughGH->GFExtras[dim]->rnsize[in][i]);
          }
          sprintf(mess,"%s%d]",mess,pughGH->GFExtras[dim]->rnsize[in][i]);
          CCTK_INFO(mess);
        }
      }
      CCTK_VInfo(CCTK_THORNSTRING, "  Maximum load skew: %f",
               pughGH->GFExtras[dim]->maxskew);

#endif

    }
  }

  free(havedims);
  free(mess);

}

 /*@@
   @routine    PUGH_pGH
   @date       Wed Feb  2 13:27:41 2000
   @author     Tom Goodale
   @desc
               This takes a cGH and returns the active pGH.
               Needed once the pugh GH Extension changes from being
               a pGH to a structure containing pGHs for multi-patch stuff.
   @enddesc
   @calls      CCTK_GHExtension

   @var        GH
   @vdesc      pointer to grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar

   @returntype pGH *
   @returndesc
               the PUGH GH extension pointer,<BR>
               or NULL if no PUGH GH extension was registered yet
   @endreturndesc
@@*/
pGH *PUGH_pGH (const cGH *GH)
{
  return ((pGH *) CCTK_GHExtension (GH, "PUGH"));
}


 /*@@
   @routine    PUGH_MPIDataType
   @date       Fri 26 Jan 2001
   @author     Thomas Radke
   @desc
               This routine returns the MPI data type to use
               for communicating a given CCTK datatype.
   @enddesc

   @var        pughGH
   @vdesc      Pointer to PUGH extensions
   @vtype      const pGH *
   @vio        in
   @endvar
   @var        cctk_type
   @vdesc      CCTK data type identifier
   @vtype      int
   @vio        in
   @endvar

   @returntype MPI_Datatype
   @returndesc
               the appropriate datatype for success
               otherwise MPI_DATATYPE_NULL
   @endreturndesc
@@*/
#ifdef CCTK_MPI
MPI_Datatype PUGH_MPIDataType (const pGH *pughGH, int cctk_type)
{
  MPI_Datatype retval;


  switch (cctk_type)
  {
    case CCTK_VARIABLE_BYTE:      retval = PUGH_MPI_BYTE; break;
    case CCTK_VARIABLE_INT:       retval = PUGH_MPI_INT; break;
    case CCTK_VARIABLE_REAL:      retval = PUGH_MPI_REAL; break;
    case CCTK_VARIABLE_COMPLEX:   retval = pughGH->PUGH_mpi_complex; break;
#ifdef CCTK_INT1
    case CCTK_VARIABLE_INT1:      retval = PUGH_MPI_INT1; break;
#endif
#ifdef CCTK_INT2
    case CCTK_VARIABLE_INT2:      retval = PUGH_MPI_INT2; break;
#endif
#ifdef CCTK_INT4
    case CCTK_VARIABLE_INT4:      retval = PUGH_MPI_INT4; break;
#endif
#ifdef CCTK_INT8
    case CCTK_VARIABLE_INT8:      retval = PUGH_MPI_INT8; break;
#endif
#ifdef CCTK_REAL4
    case CCTK_VARIABLE_REAL4:     retval = PUGH_MPI_REAL4; break;
    case CCTK_VARIABLE_COMPLEX8:  retval = pughGH->PUGH_mpi_complex8; break;
#endif
#ifdef CCTK_REAL8
    case CCTK_VARIABLE_REAL8:     retval = PUGH_MPI_REAL8; break;
    case CCTK_VARIABLE_COMPLEX16: retval = pughGH->PUGH_mpi_complex16; break;
#endif
#ifdef CCTK_REAL16
    case CCTK_VARIABLE_REAL16:    retval = PUGH_MPI_REAL16; break;
    case CCTK_VARIABLE_COMPLEX32: retval = pughGH->PUGH_mpi_complex32; break;
#endif

    default: CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                         "Unsupported CCTK variable type %d", cctk_type);
             retval = MPI_DATATYPE_NULL; break;
  }

  return (retval);
}
#endif /* CCTK_MPI */


 /*@@
   @routine    PUGH_PrintStorageReport
   @author     Gabrielle Allen
   @date       16th Sept 2001
   @desc
               Print a report about the use of storage
   @enddesc
@@*/
void PUGH_PrintStorageReport (cGH *GH)
{
  DECLARE_CCTK_PARAMETERS


  if (storage_report_every > 0 &&
      GH->cctk_iteration % storage_report_every == 0)
  {
    PUGHi_PrintStorageReport();
  }
}


 /*@@
   @routine    PUGH_PrintFinalStorageReport
   @author     Gabrielle Allen
   @date       16th Sept 2001
   @desc
               Print a report about the use of storage
   @enddesc
@@*/
void PUGH_PrintFinalStorageReport (cGH *GH)
{
  GH = GH;
  PUGHi_PrintStorageReport ();
}


 /*@@
   @routine    PUGH_PrintStorage
   @author     Gabrielle Allen
   @date       13th Oct 2001
   @desc
               Print grid variables with storage assigned
   @enddesc
@@*/
void PUGH_PrintStorage (const cGH *GH)
{
  int i;
  int countgf;
  int countarray;
  int countscalar;
  char *messgf=NULL;
  char *messarray=NULL;
  char *messscalar=NULL;

  countgf = 0;
  countarray = 0;
  countscalar = 0;

  /* Work out how long the strings are */

  for (i=0;i<CCTK_NumGroups();i++)
  {
    if (PUGH_QueryGroupStorage(GH,i,NULL))
    {
      if (CCTK_GroupTypeI(i) == CCTK_GF)
      {
        countgf += strlen(CCTK_GroupName(i))+1;
      }
      else if (CCTK_GroupTypeI(i) == CCTK_ARRAY)
      {
        countarray += strlen(CCTK_GroupName(i))+1;
      }
      else
      {
        countscalar += strlen(CCTK_GroupName(i))+1;
      }
    }
  }

  /* Now assemble the strings */
  if (countgf)
  {
    messgf = (char *)malloc((countgf+1)*sizeof(char));
    strcpy(messgf,"");
  }
  if (countarray)
  {
    messarray = (char *)malloc((countarray+100)*sizeof(char));
    strcpy(messarray,"");
  }
  if (countscalar)
  {
    messscalar = (char *)malloc((countscalar+100)*sizeof(char));
    strcpy(messscalar,"");
  }

  for (i=0;i<CCTK_NumGroups();i++)
  {
    if (PUGH_QueryGroupStorage(GH,i,NULL))
    {
      if (CCTK_GroupTypeI(i) == CCTK_GF)
      {
        sprintf(messgf,"%s%s ",messgf,CCTK_GroupName(i));
      }
      else if (CCTK_GroupTypeI(i) == CCTK_ARRAY)
      {
        sprintf(messarray,"%s%s ",messarray,CCTK_GroupName(i));
      }
      else
      {
        sprintf(messscalar,"%s%s ",messscalar,CCTK_GroupName(i));
      }
    }
  }

  if (messgf||messarray||messscalar)
  {
    CCTK_VInfo(CCTK_THORNSTRING,"Grid Variables with storage enabled:");
    if (messgf)
    {
      CCTK_VInfo(CCTK_THORNSTRING," Grid Functions: %s",messgf);
    }
    if (messarray)
    {
      CCTK_VInfo(CCTK_THORNSTRING," Grid Arrays: %s",messarray);
    }
    if (messscalar)
    {
      CCTK_VInfo(CCTK_THORNSTRING," Grid Scalars: %s",messscalar);
    }
  }

  if (messgf) free(messgf);
  if (messarray) free(messarray);
  if (messscalar) free(messscalar);

}

void CCTK_FCALL CCTK_FNAME (PUGH_PrintStorage)
                           (const cGH **GH)
{
  PUGH_PrintStorage (*GH);
}