aboutsummaryrefslogtreecommitdiff
path: root/src/include/pGH.h
blob: f21db097a4449397aff1fc21c0568553142f5d2f (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
 /*@@
   @header    pGH.h
   @date      Fri Feb 21 10:12:58 1997
   @author    Tom Goodale, Paul Walker
   @desc
              The Pugh Grid Hierarchy structure. 
   @enddesc
   @version   $Header$
 @@*/

#include "pugh_constants.h"

typedef struct PGH 
{

  /* pGH identifier */
  void     *callerid;
  int       dim;                /* The dimension of the GH */

  /* Size of the processor group */
  int       nprocs;             /* Number of processors */
  int       myproc;             /* My processor */
  int       commmodel;          /* Comm model is PUGH_ALLOCATEDBUFFERS */ 
                                /* or PUGH_DERIVEDTYPES. Currently unused */
  
  /* Size of the problems */
  int       nvariables;         /* Number of Grid variables */
  int       narrays;            /* Number of Grid arrays (inc. timelevels) */
  void   ***variables;          /* Pointers to them */

  /* What time level we're on */
  int       timelevel;

  int      *perme;              /* Periodic in each direction? */
  int       periodic;           /* Is the system periodic? */

  int       forceSync;          /* force synchronisation of GFs with storage */

  /* Coordinate information */
  /* FIXME: this is only needed for BAM */
  CCTK_REAL dx0, dy0, dz0, dt0; /* Delta of our system */

  /* Identification for a pGH */
  int       level;              /* level in Berger-Oliger stack */
  int       mglevel;            /* level in multigrid stack */
  int       convlevel;          /* level in convergence stack */

  /* pGH flags */
  int       active;             /* 1 -> GH active */
                                /* 0 -> GH is skipped for conv.test */
                                /* is used in NanCheckGH.c */

  int       comm_time;          /* time spent in communication */

#ifdef CCTK_MPI
  MPI_Datatype PUGH_mpi_complex;/* MPI datatype for generic CCTK_COMPLEX type */
  MPI_Datatype PUGH_mpi_complex8;/* same for fixed-precision CCTK_COMPLEX */
  MPI_Datatype PUGH_mpi_complex16;
  MPI_Datatype PUGH_mpi_complex32;
  MPI_Comm PUGH_COMM_WORLD;     /* the MPI communicator */

#ifdef PUGH_WITH_DERIVED_DATATYPES
  /* Derived data types for communication */
  MPI_Datatype send_char_dt [PUGH_NSTAGGER][6];
  MPI_Datatype recv_char_dt [PUGH_NSTAGGER][6];
  MPI_Datatype send_int_dt  [PUGH_NSTAGGER][6];
  MPI_Datatype recv_int_dt  [PUGH_NSTAGGER][6];
  MPI_Datatype send_real_dt [PUGH_NSTAGGER][6];
  MPI_Datatype recv_real_dt [PUGH_NSTAGGER][6];
  MPI_Datatype send_complex_dt [PUGH_NSTAGGER][6];
  MPI_Datatype recv_complex_dt [PUGH_NSTAGGER][6];
#endif
#endif

  /* Used for all grid functions */
  pGExtras **GFExtras;          /* [dim] stagger ? */
  pConnectivity **Connectivity; /* [dim] */

  char *identity_string;        /* identifier for this pGH */

} pGH;