aboutsummaryrefslogtreecommitdiff
path: root/src/include/pugh.h
blob: 8c528af4357262f765bcad4ed67c1f9881255e95 (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
/*@@
  @header pugh.h
  @author Paul Walker
  @date March 1997
  @desc
  This is just a bunch of includes and some simple definitions.
  You should really see @seefile pGF.h @seefile pGH.h and
  @seefile pughProtos.h for the pugh bit, and @seefile
  pughProblem.h for the cactus-specific bits.
  @enddesc
  @version $Header$
@@*/

#ifndef _PUGH_H_
#define _PUGH_H_ 1

#include "cctk.h"

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

/***
 Define the different datatypes used for MPI communication
 NOTE: the complex datatype is defined dynamically at runtime in SetupPGH.c
 ***/
/* byte type is easy */
#define PUGH_MPI_BYTE      MPI_UNSIGNED_CHAR

/* floating point types are architecture-independent,
   ie. a float has always 4 bytes, and a double has 8 bytes

   PUGH_MPI_REAL  is used for communicating reals of the generic CCTK_REAL type
   PUGH_MPI_REALn is used to explicitely communicate n-byte reals */
#ifdef  CCTK_REAL4
#define PUGH_MPI_REAL4  MPI_FLOAT
#endif
#ifdef  CCTK_REAL8
#define PUGH_MPI_REAL8  MPI_DOUBLE
#endif
#ifdef  CCTK_REAL16
#define PUGH_MPI_REAL16  (sizeof (CCTK_REAL16) == sizeof (long double) ?      \
                          MPI_LONG_DOUBLE : MPI_DATATYPE_NULL)
#endif


#ifdef  CCTK_REAL_PRECISION_16
#define PUGH_MPI_REAL   PUGH_MPI_REAL16
#elif   CCTK_REAL_PRECISION_8
#define PUGH_MPI_REAL   PUGH_MPI_REAL8
#elif   CCTK_REAL_PRECISION_4
#define PUGH_MPI_REAL   PUGH_MPI_REAL4
#endif


/* char type is easy */
#define PUGH_MPI_CHAR      MPI_CHAR


/* integer types are architecture-dependent:
   PUGH_MPI_INT  is used for communicating integers of the generic CCTK_INT type
   PUGH_MPI_INTn is used to explicitely communicate n-byte integers */
#ifdef  CCTK_INT8
#define PUGH_MPI_INT8  (sizeof (CCTK_INT8) == sizeof (int) ? MPI_INT :        \
                        sizeof (CCTK_INT8) == sizeof (long) ? MPI_LONG :      \
                        sizeof (CCTK_INT8) == sizeof (long long) ?            \
                        MPI_LONG_LONG_INT : MPI_DATATYPE_NULL)
#endif

#ifdef  CCTK_INT4
#define PUGH_MPI_INT4  (sizeof (CCTK_INT4) == sizeof (int) ? MPI_INT :        \
                        sizeof (CCTK_INT4) == sizeof (short) ? MPI_SHORT :    \
                        MPI_DATATYPE_NULL)
#endif

#ifdef  CCTK_INT2
#define PUGH_MPI_INT2  (sizeof (CCTK_INT2) == sizeof (short) ? MPI_SHORT :    \
                        MPI_DATATYPE_NULL)
#endif

#ifdef  CCTK_INT1
#define PUGH_MPI_INT1  MPI_CHAR
#endif

#ifdef  CCTK_INTEGER_PRECISION_8
#define PUGH_MPI_INT    PUGH_MPI_INT8
#elif   CCTK_INTEGER_PRECISION_4
#define PUGH_MPI_INT    PUGH_MPI_INT4
#elif   CCTK_INTEGER_PRECISION_2
#define PUGH_MPI_INT    PUGH_MPI_INT2
#endif


#define HEREINPUGH printf("I'm in %s at line %d\n",__FILE__,__LINE__);

#include "pugh_constants.h"

#include "pGV.h"
#include "pGH.h"

#ifdef CCTK_MPI

#define CACTUS_MPI_ERROR(fn_call)                                             \
          do {                                                                \
            int errcode;                                                      \
                                                                              \
            if ((errcode = fn_call) != 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", \
                               #fn_call, errcode, mpi_error_string);          \
              fprintf(stderr, "At line %d of file %s\n", __LINE__, __FILE__); \
            }                                                                 \
          } while (0)
#endif

#ifdef __cplusplus
extern "C"
{
#endif


int PUGH_SetupGroup (pGH *newGH,
                     int *nsize,
                     int *nghostsize,
                     int  gtype,
                     int  vtype,
                     int  dim,
                     int  n_variables,
                     int  vectorlength,
                     int  n_timelevels,
                     int  vectorgroup);

pGH *PUGH_SetupPGH(void *callerid,
                   int dim,
                   int *nsize,
                   int *nghostzones,
                   int *perme);

void PUGH_GFSize(int dim, int *nsize);

void PUGH_GFGhostsize(int dim, int *ghostsize);
void PUGH_GFPeriodic(int dim, int *perme);

pGH *PUGH_pGH(const cGH *GH);

int PUGH_Evolve(tFleshConfig *config);

int PUGH_GetBounds(int is_gf,
                   int dim,
                   int **bounds,
                   int *nprocs,
                   int *nsize);

void PUGH_InitializeMemory (const char *do_initialize_memory,
                            int vtype,
                            int bytes,
                            void *data);

void PUGH_Terminate (cGH *GH);

int PUGH_ParallelInit(cGH *GH);

int PUGH_Abort(cGH *GH, int retval);

int PUGH_MyProc(const cGH *GH);

int PUGH_nProcs(const cGH *GH);

int PUGH_Exit(cGH *GH, int retval);

const int *PUGH_Topology(const cGH *GH, int dim);

int PUGH_SetTopology (int dim, const int topology[]);

extern int (*PUGH_GenerateTopology)(int dim, 
                                    int total_procs, 
                                    const int *nsize,
                                    const int *nghostzones, 
                                    int *nprocs);

#ifdef CCTK_MPI
MPI_Datatype PUGH_MPIDataType (const pGH *pughGH, int cctk_type);
#endif

#ifdef __cplusplus
}
#endif

#endif /* defined _PUGH_H_ */