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

#include <stdlib.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)

void PUGH_Report(CCTK_ARGUMENTS);

 /*@@
   @routine    PUGH_Report
   @date       Sunday 12 September 1999
   @author     Gabrielle Allen
   @desc 
   Report on PUGH set up
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

void PUGH_Report(CCTK_ARGUMENTS)
{
  DECLARE_CCTK_PARAMETERS                                                   

  pGH *pughGH;           
  int i,gi,dim;
  int *havedims;
  char *mess;
#ifdef CCTK_MPI
  int in;
#endif

  pughGH = PUGH_pGH(cctkGH);

  havedims = (int *)calloc(cctkGH->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);
      havedims[dim-1]=1;
    }
  }

  for (dim=0;dim<cctkGH->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);
      } 
      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     
   @calledby   
   @history 
 
   @endhistory 

@@*/
pGH *PUGH_pGH(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      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 (pGH *pughGH, int cctk_type)
{
  MPI_Datatype retval;


  switch (cctk_type)
  {
    case CCTK_VARIABLE_CHAR:      retval = PUGH_MPI_CHAR; 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_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 */