aboutsummaryrefslogtreecommitdiff
path: root/src/qlm_broadcast.c
blob: 32661773e38b7fd52030c7bb3bf1d2383760decf (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
#include <assert.h>
#include <stdio.h>
#include <stdlib.h>

#include <cctk.h>
#include <cctk_Arguments.h>
#include <cctk_Parameters.h>

#ifdef CCTK_MPI
#  include <mpi.h>
#endif



/* Broadcast a vector element of a grid group */
static void
bcast (cGH const * restrict const cctkGH,
       char const * restrict const group,
       int const vecind,
       int const root)
{
#ifdef CCTK_MPI
  
  int ierr;
  
  assert (group);
  
  MPI_Comm comm;
  if (CCTK_IsFunctionAliased ("GetMPICommWorld"))
  {
    comm = * (MPI_Comm const *) GetMPICommWorld (cctkGH);
  }
  else
  {
    comm = MPI_COMM_WORLD;
  }
  
  int num_procs;
  MPI_Comm_size (comm, & num_procs);
  
  assert (root >= 0 && root < num_procs);
  
  int const gi = CCTK_GroupIndex (group);
  assert (gi >= 0);
  
  cGroup data;
  ierr = CCTK_GroupData (gi, & data);
  assert (! ierr);
  
  assert (vecind >= 0 && vecind < data.vectorlength);
  
  if (data.numvars == 0) return;
  
  MPI_Datatype mpitype;
  int items;
  switch (data.vartype)
  {
  case CCTK_VARIABLE_INT:
    assert (sizeof (CCTK_INT) == sizeof (int));
    mpitype = MPI_INT;
    items = 1;
    break;
  case CCTK_VARIABLE_REAL:
    assert (sizeof (CCTK_REAL) == sizeof (double));
    mpitype = MPI_DOUBLE;
    items = 1;
    break;
  case CCTK_VARIABLE_COMPLEX:
    assert (sizeof (CCTK_COMPLEX) == 2 * sizeof (double));
    mpitype = MPI_DOUBLE;
    items = 2;
    break;
  default:
    assert (0);
  }
  assert (items > 0);
  
  cGroupDynamicData dyndata;
  ierr = CCTK_GroupDynamicData (cctkGH, gi, & dyndata);
  assert (! ierr);
  
  int elems = 1;
  {
    int d;
    for (d = 0; d < dyndata.dim; ++ d)
    {
      elems *= dyndata.lsh[d];
    }
  }
  assert (elems >= 0);
  
  int const v0 = CCTK_FirstVarIndexI (gi);
  assert (v0 >= 0);
  
  int size;
  MPI_Type_size (mpitype, & size);
  assert (items * size == CCTK_VarTypeSize (data.vartype));
  
  int const numvectors = data.numvars / data.vectorlength;
  
  {
    int vec;
    for (vec = 0; vec < numvectors; ++ vec)
    {
      int const vi = v0 + vecind + data.vectorlength * vec;
      void * restrict const ptr = CCTK_VarDataPtrI (cctkGH, 0, vi);
      assert (ptr);
      
      MPI_Bcast (ptr, elems * items, mpitype, root, comm);
    }
  }

#endif  /* ifdef CCTK_MPI */
}


void CCTK_FCALL
CCTK_FNAME(qlm_broadcast) (CCTK_POINTER_TO_CONST * restrict cctkGH_);

void CCTK_FCALL
CCTK_FNAME(qlm_broadcast) (CCTK_POINTER_TO_CONST * restrict const cctkGH_)
{
  cGH const * restrict const cctkGH = * (cGH const * const *) cctkGH_;
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  
  if (veryverbose)
  {
    CCTK_INFO ("Broadcasting quasi-local measures");
  }
  
  int const num_procs = CCTK_nProcs (cctkGH);
  int hn; 
  for (hn = 0; hn < num_surfaces; ++ hn)
  {
    int const root = hn % num_procs;
    
    bcast (cctkGH, "QuasiLocalMeasures::qlm_state"      , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_state_p"    , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_grid_int"   , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_grid_real"  , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_grid_real_p", hn, root);
    
    bcast (cctkGH, "QuasiLocalMeasures::qlm_shapes"               , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_shapes_p"             , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_coordinates"          , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_tetrad_l"             , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_tetrad_n"             , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_tetrad_m"             , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_newman_penrose"       , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_weyl_scalars"         , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_ricci_scalars"        , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_twometric"            , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_killing_vector"       , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_killed_twometric"     , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_invariant_coordinates", hn, root);
    
    bcast (cctkGH, "QuasiLocalMeasures::qlm_multipole_moments", hn, root);
    
    bcast (cctkGH, "QuasiLocalMeasures::qlm_3determinant", hn, root);
    
    bcast (cctkGH, "QuasiLocalMeasures::qlm_scalars"  , hn, root);
    bcast (cctkGH, "QuasiLocalMeasures::qlm_scalars_p", hn, root);
    
  } /* for hn */
}