aboutsummaryrefslogtreecommitdiff
path: root/Carpet/Carpet/src/Comm.cc
blob: 4d23b675455b10a7eea96d8fc6f2bf64a700e459 (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
319
320
321
322
323
324
325
#include <algorithm>
#include <cassert>
#include <cmath>
#include <cstdlib>
#include <iostream>
#include <vector>

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

#include <Requirements.hh>

#include <Timer.hh>

#include <ggf.hh>
#include <gh.hh>

#include <carpet.hh>



namespace Carpet {

  using namespace std;



  static void ProlongateGroupBoundaries (const cGH* cctkGH,
                                         const vector<int>& groups);


  // Carpet's overload function for CCTK_SyncGroupsByDirI()
  // which synchronises a set of groups given by their indices.
  // Synchronisation of individual directions is not (yet?) implemented.
  //
  // returns the number of groups successfully synchronised
  //        -1 if a group in the set doesn't have storage assigned
  int SyncGroupsByDirI (const cGH* cctkGH,
                        int num_groups,
                        const int *groups,
                        const int *directions)
  {
    int group, retval = 0;
    vector<int> groups_set;
    groups_set.reserve(num_groups);

    // individual directions aren't supported (yet?)
    if (directions != NULL) {
      CCTK_WARN (0, "Carpet doesn't support synchronisation of individual "
                    "directions");
    }

    for (group = 0; group < num_groups; group++) {
      if (CCTK_NumVarsInGroupI (groups[group]) > 0) {
        groups_set.push_back (groups[group]);
      }
    }

    if (groups_set.size() > 0) {
      retval = SyncProlongateGroups (cctkGH, groups_set);
      if (retval == 0) {
        retval = groups_set.size();
      }
    }

    return retval;
  }


  // synchronises ghostzones and prolongates boundaries of a set of groups
  //
  // returns 0 for success and -1 if the set contains a group with no storage
  int SyncProlongateGroups (const cGH* cctkGH, const vector<int>& groups,
                            cFunctionData const* function_data)
  {
    int retval = 0;
    DECLARE_CCTK_PARAMETERS;

    assert (groups.size() > 0);

    // check consistency of all groups:
    // create a new set with empty and no-storage groups removed
    vector<int> goodgroups;
    goodgroups.reserve (groups.size());
    for (size_t group = 0; group < groups.size(); group++) {
      const int g = groups.AT(group);
      const int grouptype = CCTK_GroupTypeI (g);
      char* const groupname = CCTK_GroupName (g);
      Checkpoint ("SyncGroup \"%s\" iteration=%d time=%g",
                  groupname,
                  cctkGH->cctk_iteration, (double) cctkGH->cctk_time);

      if (grouptype == CCTK_GF) {
        if (reflevel == -1) {
          CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
                      "Cannot synchronise in global mode "
                      "(Tried to synchronise group \"%s\")",
                      groupname);
        }
        if (map != -1 and component == -1) {
          if (maps == 1) {
            CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                        "Synchronising group \"%s\" in singlemap mode",
                        groupname);
          } else {
            CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
                        "Cannot synchronise in singlemap mode "
                        "(Tried to synchronise group \"%s\")",
                        groupname);
          }
        }
        if (component != -1) {
          if (maps == 1 and vhh.AT(map)->local_components(reflevel) == 1) {
            CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                        "Synchronising group \"%s\" in local mode",
                        groupname);
          } else {
            CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
                        "Cannot synchronise in local mode "
                        "(Tried to synchronise group \"%s\")",
                        groupname);
          }
        }
      }

      if (not CCTK_QueryGroupStorageI (cctkGH, g)) {
        CCTK_VWarn (4, __LINE__, __FILE__, CCTK_THORNSTRING,
                    "Cannot synchronise group \"%s\" because it has no storage",
                    groupname);
        retval = -1;
      }
      else if (CCTK_NumVarsInGroupI (g) > 0) {
        goodgroups.push_back(g);
      }

      free (groupname);
    } // for g

    if (goodgroups.size() > 0) {
      
#ifdef REQUIREMENTS_HH
      Requirements::Sync(function_data, goodgroups,
                         cctkGH->cctk_iteration, reflevel, timelevel);
#endif
      
      // prolongate boundaries
      bool const local_do_prolongate = do_prolongate and not do_taper;
      if (local_do_prolongate) {
        static Timers::Timer timer ("Prolongate");
        timer.start();
        ProlongateGroupBoundaries (cctkGH, goodgroups);
        timer.stop();
      }

      // This was found to be necessary on Hopper, otherwise memory
      // seems to be overwritten while prolongating/syncronizing. It
      // looks liks this might be an MPI implementation issue, but
      // this is not clear. A barrier at this point seems to be a
      // sufficient workaround, and is now used on Hopper. For more
      // information about this ask Frank Loeffler
      // <knarf@cct.lsu.edu>.
#ifdef CARPET_MPI_BARRIER_PROLONGATE_SYNC
      Carpet::NamedBarrier(cctkGH,
                           8472211063, "CARPET_MPI_BARRIER_PROLONGATE_SYNC");
#endif

      if (sync_barriers)
      {
        static Timers::Timer barrier_timer ("ProlongateSyncBarrier",0,true);
        barrier_timer.start();
        CCTK_Barrier(cctkGH);
        barrier_timer.stop();
      }

      // synchronise ghostzones
      if (sync_during_time_integration or local_do_prolongate) {
        static Timers::Timer timer ("Sync");
        timer.start();
        SyncGroups (cctkGH, goodgroups);
        timer.stop();
      }
      
    }
    
    return retval;
  }

  // Prolongate the boundaries of all CCTK_GF groups in the given set
  static void ProlongateGroupBoundaries (const cGH* cctkGH,
                                         const vector<int>& groups)
  {
    DECLARE_CCTK_PARAMETERS;
    
    if (reflevel == 0) return;
    
    Checkpoint ("ProlongateGroups");

    assert (groups.size() > 0);

    // use the current time here (which may be modified by the user)
    const CCTK_REAL time = cctkGH->cctk_time;

    static vector<Timers::Timer*> timers;
    if (timers.empty()) {
      timers.push_back(new Timers::Timer("comm_state[0].create"));
      for (astate state = static_cast<astate>(0);
           state != state_done;
           state = static_cast<astate>(static_cast<int>(state)+1))
      {
        ostringstream name1;
        name1 << "comm_state[" << timers.size() << "]"
              << "." << tostring(state) << ".user";
        timers.push_back(new Timers::Timer(name1.str()));
        ostringstream name2;
        name2 << "comm_state[" << timers.size() << "]"
              << "." << tostring(state) << ".step";
        timers.push_back(new Timers::Timer(name2.str()));
      }
    }
    
    vector<Timers::Timer*>::iterator ti = timers.begin();
    (*ti)->start();
    for (comm_state state; not state.done(); state.step()) {
      (*ti)->stop(); ++ti; (*ti)->start();
      for (int group = 0; group < (int)groups.size(); ++group) {
        const int g = groups.AT(group);
        const int grouptype = CCTK_GroupTypeI (g);
        if (grouptype != CCTK_GF) {
          continue;
        }
        assert (reflevel>=0 and reflevel<reflevels);
        const int active_tl =
          groupdata.AT(g).activetimelevels.AT(mglevel).AT(reflevel);
        assert (active_tl>=0);
        const int tl = active_tl > 1 ? timelevel : 0;
        
        for (int m = 0; m < (int)arrdata.AT(g).size(); ++m) {
          for (int v = 0; v < (int)arrdata.AT(g).AT(m).data.size(); ++v) {
            ggf *const gv = arrdata.AT(g).AT(m).data.AT(v);
            gv->ref_bnd_prolongate_all
              (state, tl, reflevel, mglevel, time);
          }
        }
      }
      (*ti)->stop(); ++ti; (*ti)->start();
    }
    (*ti)->stop();
    ++ti; assert(ti == timers.end());
  }


  // synchronises a set of groups
  void SyncGroups (const cGH* cctkGH, const vector<int>& groups)
  {
    DECLARE_CCTK_PARAMETERS;

    Checkpoint ("SyncGroups");

    assert (groups.size() > 0);
    
    if (CCTK_IsFunctionAliased("Accelerator_PreSync")) {
      Accelerator_PreSync(cctkGH, &groups.front(), groups.size());
    }

    static vector<Timers::Timer*> timers;
    if (timers.empty()) {
      timers.push_back(new Timers::Timer("comm_state[0].create"));
      for (astate state = static_cast<astate>(0);
           state != state_done;
           state = static_cast<astate>(static_cast<int>(state)+1))
      {
        ostringstream name1;
        name1 << "comm_state[" << timers.size() << "]"
              << "." << tostring(state) << ".user";
        timers.push_back(new Timers::Timer(name1.str()));
        ostringstream name2;
        name2 << "comm_state[" << timers.size() << "]"
              << "." << tostring(state) << ".step";
        timers.push_back(new Timers::Timer(name2.str()));
      }
    }
    
    vector<Timers::Timer*>::iterator ti = timers.begin();
    (*ti)->start();
    for (comm_state state; not state.done(); state.step()) {
      (*ti)->stop(); ++ti; (*ti)->start();
      for (int group = 0; group < (int)groups.size(); ++group) {
        const int g = groups.AT(group);
        const int grouptype = CCTK_GroupTypeI (g);
        const int ml = grouptype == CCTK_GF ? mglevel : 0;
        const int rl = grouptype == CCTK_GF ? reflevel : 0;
        const int active_tl = groupdata.AT(g).activetimelevels.AT(ml).AT(rl);
        assert (active_tl>=0);
        const int tl = active_tl > 1 ? timelevel : 0;
        for (int m = 0; m < (int)arrdata.AT(g).size(); ++m) {
          for (int v = 0; v < (int)arrdata.AT(g).AT(m).data.size(); ++v) {
            arrdesc& array = arrdata.AT(g).AT(m);
            array.data.AT(v)->sync_all (state, tl, rl, ml);
          }
        }
      }
      (*ti)->stop(); ++ti; (*ti)->start();
    }
    (*ti)->stop();
    ++ti; assert(ti == timers.end());
    
    if (CCTK_IsFunctionAliased("Accelerator_PostSync")) {
      Accelerator_PostSync(cctkGH, &groups.front(), groups.size());
    }
  }


  int EnableGroupComm (const cGH* cctkGH, const char* groupname)
  {
    // Communication is always enabled
    return 0;
  }

  int DisableGroupComm (const cGH* cctkGH, const char* groupname)
  {
    // Communication is always enabled
    return -1;
  }

} // namespace Carpet