aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetLib/src/dist.cc
blob: 89acacfad9a9d37c73a9b428d90862f977a6d454 (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
#include <cassert>
#include <typeinfo>

#include <mpi.h>
#ifdef _OPENMP
#  include <omp.h>
#endif

#include "cctk.h"
#include "cctk_Parameters.h"

#include "defs.hh"
#include "limits.hh"
#include "startup_time.hh"

#include "dist.hh"

using namespace std;



namespace dist {
  
  MPI_Comm comm_ = MPI_COMM_NULL;
  
  MPI_Datatype mpi_complex8  = MPI_DATATYPE_NULL;
  MPI_Datatype mpi_complex16 = MPI_DATATYPE_NULL;
  MPI_Datatype mpi_complex32 = MPI_DATATYPE_NULL;
  
  int total_num_threads_ = -1;
  
  void init (int& argc, char**& argv) {
    MPI_Init (&argc, &argv);
    pseudoinit (MPI_COMM_WORLD);
  }
  
  void pseudoinit (MPI_Comm const c) {
    comm_ = c;
    
#ifdef HAVE_CCTK_REAL4
    CCTK_REAL4 dummy4;
    MPI_Type_contiguous (2, mpi_datatype(dummy4), &mpi_complex8);
    MPI_Type_commit (&mpi_complex8);
#endif
#ifdef HAVE_CCTK_REAL8
    CCTK_REAL8 dummy8;
    MPI_Type_contiguous (2, mpi_datatype(dummy8), &mpi_complex16);
    MPI_Type_commit (&mpi_complex16);
#endif
#ifdef HAVE_CCTK_REAL16
    CCTK_REAL16 dummy16;
    MPI_Type_contiguous (2, mpi_datatype(dummy16), &mpi_complex32);
    MPI_Type_commit (&mpi_complex32);
#endif
    
    // Output startup time
    // cerr << "QQQ: pseudoinit[1]" << endl;
    CarpetLib::output_startup_time ();
    // cerr << "QQQ: pseudoinit[2]" << endl;
    
    // Check and/or modify system limits
    CarpetLib::set_system_limits ();
    // cerr << "QQQ: pseudoinit[3]" << endl;
    
    // cerr << "QQQ: pseudoinit[4]" << endl;
    collect_total_num_threads ();
    // cerr << "QQQ: pseudoinit[5]" << endl;
  }
  
  void finalize () {
    MPI_Finalize ();
  }
  
  
  
  // Create an MPI datatype from a C datatype description
  
  ostream& operator<< (ostream& os, mpi_struct_descr_t const& descr)
  {
    int type_size;
    MPI_Type_size (descr.type, &type_size);
    os << "{"
       << "blocklength:" << descr.blocklength << ","
       << "displacement:" << descr.displacement << ","
       << "type:" << descr.type << ","
       << "type_size:" << type_size << ","
       << "field_name:" << descr.field_name << ","
       << "type_name:" << descr.type_name
       << "}";
    return os;
  }
  
  MPI_Datatype create_mpi_datatype (size_t const count,
                                    mpi_struct_descr_t const descr[],
                                    char const * const name, size_t const size)
  {
    DECLARE_CCTK_PARAMETERS;
    int blocklengths[count];
    MPI_Aint displacements[count];
    MPI_Datatype types[count];
    for (size_t n=0; n<count; ++n) {
      blocklengths [n] = descr[n].blocklength;
      displacements[n] = descr[n].displacement;
      types        [n] = descr[n].type;
    }
    MPI_Datatype newtype;
    MPI_Type_struct (count, blocklengths, displacements, types, &newtype);
    MPI_Type_commit (&newtype);
    if (verbose) {
      CCTK_VInfo (CCTK_THORNSTRING,
                  "Creating new MPI type for C type %s:", name);
      cout << "   Type has " << count << " components" << endl;
      for (size_t n=0; n<count; ++n) {
        cout << "   [" << n << "]: " << descr[n] << endl;
      }
      cout << "   New MPI type ID is " << newtype << endl;
      int datatypesize;
      MPI_Type_size (newtype, &datatypesize);
      cout << "   C type size is " << size << endl;
      cout << "   MPI type size is " << datatypesize << endl;
    }
    return newtype;
  }
  
#if 0
  
  ostream&
  generic_mpi_datatype_t::field_t::output (ostream& os) const
  {
    int type_size;
    MPI_Type_size (mpi_datatype, &type_size);
    os << "{"
       << "offset:" << offset << ","
       << "count:" << count << ","
       << "mpi_datatype:" << mpi_datatype << ","
       << "type_size:" << type_size << ","
       << "field_name:" << field_name << ","
       << "type_name:" << type_name
       << "}";
    return os;
  }
  
  generic_mpi_datatype_t::generic_mpi_datatype_t (string const type_name_)
    : type_name (type_name_), type_is_committed (false)
  {
  }
  
  template <typename U>
  void
  generic_mpi_datatype_t::add_field (size_t const offset, size_t const count,
                                     string const field_name)
  {
    assert (not type_is_committed);
    U u;
    entries.push_back (field_t (offset, count, mpi_datatype(u),
                                field_name, typeid(U).name()));
  }
  
  void
  generic_mpi_datatype_t::commit ()
  {
    DECLARE_CCTK_PARAMETERS;
    
    // Debug output
    if (verbose) {
      CCTK_VInfo (CCTK_THORNSTRING,
                  "Creating new MPI type for C type %s:", type_name.c_str());
      cout << *this;
    }
    
    assert (not type_is_committed);
    type_is_committed = true;
    
    // Out of caution -- this could be allowed
    assert (not entries.empty());
    
    // Create MPI type
    size_t const count = entries.size();
    int          blocklengths [count+1];
    MPI_Aint     displacements[count+1];
    MPI_Datatype types        [count+1];
    {
      size_t n = 0;
      for (list<field_t>::const_iterator ifield =
             entries.begin(); ifield!=entries.end(); ++ifield, ++n)
      {
        blocklengths [n] = ifield->count;
        displacements[n] = ifield->offset;
        types        [n] = ifield->mpi_datatype;
      }
      assert (n == count);
      // Add MPI_UB
      blocklengths [n] = 1;
      displacements[n] = type_size();
      types        [n] = MPI_UB;
    }
    
    MPI_Type_struct
      (count+1, blocklengths, displacements, types, &mpi_datatype);
    MPI_Type_commit (&mpi_datatype);
  }
  
  ostream&
  generic_mpi_datatype_t::output (ostream& os) const
  {
    cout << "Datatype: " << type_name << endl;
    size_t const count = entries.size();
    cout << "   Type has " << count << " components" << endl;
    {
      size_t n = 0;
      for (list<field_t>::const_iterator ifield =
             entries.begin(); ifield!=entries.end(); ++ifield, ++n)
      {
        cout << "   [" << n << "]: " << *ifield << endl;
      }
      assert (n == count);
    }
    cout << "   MPI type ID: " << mpi_datatype << endl;
    int datatypesize;
    MPI_Type_size (mpi_datatype, &datatypesize);
    cout << "   C   type size: " << size << endl;
    cout << "   MPI type size: " << datatypesize << endl;
    return os;
  }
  
#endif
  
  
  
  void checkpoint (const char* file, int line) {
    DECLARE_CCTK_PARAMETERS;
    if (verbose) {
      int rank;
      MPI_Comm_rank (comm(), &rank);
      printf ("CHECKPOINT: processor %d, file %s, line %d\n",
	      rank, file, line);
    }
    if (barriers) {
      MPI_Barrier (comm());
    }
  }
  
  // Set number of threads
  void set_num_threads (int const num_threads)
  {
#ifdef _OPENMP
    if (num_threads > 0) {
      // Set number of threads which should be used
      // TODO: do this at startup, not in this routine
      omp_set_num_threads (num_threads);
    }
#else
    if (num_threads > 0 and num_threads != 1) {
      CCTK_WARN (CCTK_WARN_ABORT,
                 "OpenMP is not enabled.  Cannot set the number of threads.");
    }
#endif
  }
  
  // Global number of threads
  void collect_total_num_threads ()
  {
    int const mynthreads = num_threads();
    // cerr << "QQQ: collect_total_num_threads[1]" << endl;
    MPI_Allreduce
      (const_cast <int *> (& mynthreads), & total_num_threads_, 1, MPI_INT,
       MPI_SUM, comm());
    // cerr << "QQQ: collect_total_num_threads[2]" << endl;
    assert (total_num_threads_ >= size());
  }
  
  
  
  char const * c_datatype_name (unsigned const type)
  {
    switch (type) {
    case  0: return "char";
    case  1: return "signed char";
    case  2: return "unsigned char";
    case  3: return "short";
    case  4: return "unsigned short";
    case  5: return "int";
    case  6: return "unsigned int";
    case  7: return "long";
    case  8: return "unsigned long";
    case  9: return "long long";
    case 10: return "unsigned long long";
    case 11: return "float";
    case 12: return "double";
    case 13: return "long double";
#ifdef HAVE_CCTK_COMPLEX8
    case 14: return "CCTK_COMPLEX8";
#endif
#ifdef HAVE_CCTK_COMPLEX16
    case 15: return "CCTK_COMPLEX16";
#endif
#ifdef HAVE_CCTK_COMPLEX32
    case 16: return "CCTK_COMPLEX32";
#endif
    }
    assert (0); abort();
    return NULL;
  }
  
} // namespace dist