summaryrefslogtreecommitdiff
path: root/src/comm/CactusDefaultComm.c
blob: ccbc203d30640c63245fe37bed40ecbbfcfc3c68 (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
 /*@@
   @file      CactusDefaultComm.c
   @date      Tue Sep 29 15:06:22 1998
   @author    Tom Goodale
   @desc 
   Default communication routines.
   @enddesc 
   @version $Header$
 @@*/


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

#include "cctk_Flesh.h"
#include "cctk_Groups.h"
#include "cctk_Constants.h"
#include "CactusMainDefaults.h"
#include "cctk_GHExtensions.h"
#include "cctki_GHExtensions.h"

#include "cctk_ParamCheck.h"

#ifdef MPI
#include "mpi.h"
#endif

static char *rcsid = "$Header$";

CCTK_FILEVERSION(comm_CactusDefaultComm_c)

/********************************************************************
 *********************     Global Data   *****************************
 ********************************************************************/

/* FIXME:  This should be in a header somewhere */ 
#ifdef MPI
extern char MPI_Active;
#endif

/********************************************************************
 *********************     Local Definitions   **********************
 ********************************************************************/

#ifdef MPI
#define CACTUS_MPI_ERROR(xf)  do {int errcode; \
                                    if((errcode = xf) != MPI_SUCCESS)                     \
                                    {                                                     \
                                      char mpi_error_string[MPI_MAX_ERROR_STRING+1];      \
                                      int resultlen;                                      \
                                      MPI_Error_string(errcode, mpi_error_string, &resultlen);\
                                      fprintf(stderr, "MPI Call %s returned error code %d (%s)\n", \
                                      #xf, errcode, mpi_error_string);                    \
                                      fprintf(stderr, "At line %d of file %s\n",                   \
                                             __LINE__, __FILE__);                         \
                                    }                                                     \
                                  } while (0)
#endif

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


 /*@@
   @routine    CactusDefaultSetupGH
   @date       Tue Sep 29 15:06:22 1998
   @author     Tom Goodale
   @desc 
   Default cactus SetupGH routine.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/
cGH *CactusDefaultSetupGH(tFleshConfig *config, int convergence_level)
{
  cGH *retval;

  cGH *thisGH;

  int n_groups;
  int n_variables;
  int variable;
  int ntimelevels;
  int level;
  int cctk_dim;

  retval = NULL;

  /* Put this in for the moment until parameter stuff is done. */
  if(convergence_level > 0)
  {
    return retval;
  }

  /* Initialise this since it is used later and in exceptional
   * circumstances might not be initialsed beforehand. 
   */

  variable = -1;

  /* Create a new Grid Hierarchy */
  thisGH = (cGH *)malloc(sizeof(cGH));

  if(thisGH)
  {
    thisGH->cctk_dim = CCTK_MaxDim();

    /* Need this to be at least one otherwise the memory allocation will fail. */
    cctk_dim = thisGH->cctk_dim;
    if(thisGH->cctk_dim == 0) cctk_dim = 1;
    thisGH->cctk_iteration    = 0;
    thisGH->cctk_gsh          = (int *)malloc(cctk_dim*sizeof(int));
    thisGH->cctk_lsh          = (int *)malloc(cctk_dim*sizeof(int));
    thisGH->cctk_lbnd         = (int *)malloc(cctk_dim*sizeof(int));
    thisGH->cctk_ubnd         = (int *)malloc(cctk_dim*sizeof(int));
    
    thisGH->cctk_lssh         = (int *)malloc(CCTK_NSTAGGER*cctk_dim*sizeof(int));
    thisGH->cctk_to           = (int *)malloc(cctk_dim*sizeof(int));
    thisGH->cctk_from         = (int *)malloc(cctk_dim*sizeof(int));
    thisGH->cctk_bbox         = (int *)malloc(2*cctk_dim*sizeof(int));
    thisGH->cctk_nghostzones  = (int *)malloc(2*cctk_dim*sizeof(int));
    thisGH->cctk_levfac       = (int *)malloc(cctk_dim*sizeof(int));
    thisGH->cctk_delta_space  = (CCTK_REAL *)malloc(cctk_dim*sizeof(CCTK_REAL));
    /* FIXME : Next line goes when coords are done properly */
    thisGH->cctk_origin_space = (CCTK_REAL *)malloc(cctk_dim*sizeof(CCTK_REAL));

    thisGH->cctk_delta_time = 1;
    thisGH->cctk_convlevel = 0;

    n_variables = CCTK_NumVars();

    /* Allocate memory for the variable data pointers.
     * Note we want at least one to prevent memory allocattion from failing !
     */
    thisGH->data = (void ***)malloc((n_variables ? n_variables:1)*sizeof(void **));

    if(thisGH->data)
    {
      for(variable = 0; variable < n_variables; variable++)
      {
        ntimelevels = CCTK_NumTimeLevelsFromVarI(variable);

        thisGH->data[variable] = (void **)malloc(ntimelevels*sizeof(void *));
        if(thisGH->data[variable])
        {
          for(level = 0; level < ntimelevels; level++)
          {
            thisGH->data[variable][level] = NULL;
          }
        }
        else
        {
          break;
        }
      }
    }

    thisGH->extensions = NULL;

    /* Allocate memory for the group data pointers.
     * Note we want at least one to prevent memory allocattion from failing !
     */
    n_groups = CCTK_NumGroups();
    thisGH->GroupData = (cGHGroupData *)malloc((n_groups ? n_groups:1)*sizeof(cGHGroupData));

  }
  
  if(thisGH && 
     thisGH->cctk_gsh &&
     thisGH->cctk_lsh &&
     thisGH->cctk_lbnd &&
     thisGH->cctk_ubnd &&
     thisGH->cctk_lssh &&
     thisGH->cctk_from &&
     thisGH->cctk_to &&
     thisGH->cctk_bbox &&
     thisGH->cctk_nghostzones &&
     thisGH->cctk_levfac &&
     thisGH->cctk_delta_space &&
     thisGH->cctk_origin_space &&
     thisGH->data &&
     variable == n_variables &&
     thisGH->GroupData)
  {
    /* Traverse list of GH setup routines. */
    CCTKi_SetupGHExtensions(config, convergence_level, thisGH);

    retval = thisGH;
  }

  return retval;
}

 /*@@
   @routine    CactusDefaultMyProc
   @date       Tue Jan 23 1999
   @author     Gabrielle Allen
   @desc 
   Default cactus MyProc routine.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

int CactusDefaultMyProc(cGH *GH)
{
  int myproc;

  if(CCTK_ParamChecking())
  {
    myproc = 0;
  }
  else
  {
#ifdef MPI
    if(MPI_Active)
    {
      CACTUS_MPI_ERROR(MPI_Comm_rank(MPI_COMM_WORLD, &myproc));
    }
#else
    myproc = 0;
#endif
  }

  return myproc;
}

 /*@@
   @routine    CactusDefaultnProcs
   @date       Tue Jan 23 1999
   @author     Gabrielle Allen
   @desc 
   Default cactus nProcs routine.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

int CactusDefaultnProcs(cGH *GH)
{
  int nprocs;

  if(CCTK_ParamChecking())
  {
    nprocs = CCTK_ParamCheckNProcs();
  }
  else
  {
#ifdef MPI
    if(MPI_Active)
    {
      CACTUS_MPI_ERROR(MPI_Comm_size(MPI_COMM_WORLD, &nprocs));
    }
#else
    nprocs = 1;
#endif
  }
  
  return nprocs;
}

 /*@@
   @routine    CactusDefaultExit
   @date       Tue Apr 18 15:21:15 2000
   @author     Gerd Lanfermann
   @desc 
   The default for when people call CCTK_Exit.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/
int CactusDefaultExit(cGH *GH, int retval)
{
#ifdef MPI  
  if(MPI_Active)
  {
    CACTUS_MPI_ERROR(MPI_Finalize());
  }
#endif
  exit(retval);
  return (0);
}

 /*@@
   @routine    CactusDefaultBarrier
   @date       Tue Apr 18 15:21:42 2000
   @author     Tom Goodale
   @desc 
   The default for when people call CCTK_Barrier
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/
int CactusDefaultBarrier(cGH *GH) 
{
  return 0;
}