aboutsummaryrefslogtreecommitdiff
path: root/Carpet/Carpet/src/Restrict.cc
blob: 6dc040c2383719461d39fc8d22ddcb3bdcb0d411 (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
#include <cassert>
#include <cmath>
#include <cstdlib>

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

#include <Requirements.hh>

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

#include <carpet.hh>
#include <Timers.hh>



namespace Carpet {
  
  using namespace std;
  
  // restricts a set of groups
  static void RestrictGroups (const cGH* cctkGH, const vector<int>& groups);
  
  
  void Restrict (const cGH* cctkGH)
  {
    DECLARE_CCTK_PARAMETERS;
    
    assert (is_level_mode());
    
    if (suppress_restriction) {
      Checkpoint ("Restriction suppressed");
      return;
    }
    
    Checkpoint ("Restrict");
    
    // all but the finest level are restricted
    if (reflevel == reflevels-1) {
      return;
    }

    // remove all groups that are non-GFs, empty, or have no storage assigned
    vector<int> groups;
    groups.reserve (CCTK_NumGroups());

    for (int group = 0; group < CCTK_NumGroups(); ++group) {
      operator_type const op = groupdata.AT(group).transport_operator;
      bool const do_restrict = op != op_none and op != op_sync;
      if (do_restrict
          and CCTK_NumVarsInGroupI(group) > 0
          and CCTK_QueryGroupStorageI(cctkGH, group))
      {
        groups.push_back (group);
      }
    }
    if (groups.empty()) return;
    
#ifdef REQUIREMENTS_HH
    Requirements::Restrict(groups, cctkGH->cctk_iteration, reflevel);
#endif
    
    // Restrict
    {
      static Timer timer ("Restrict");
      timer.start();
      RestrictGroups (cctkGH, groups);
      timer.stop();
    }

    // Synchronise
    // TODO: Why? Restricted grid functions need to have boundary
    // conditions applied in postrestrict anyway, so this should not
    // be necessary.
#if 0
    {
      static Timer timer ("RestrictSync");
      timer.start();
      SyncGroups (cctkGH, groups);
      timer.stop();
    }
#endif
  }
  

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

    static vector<Timer*> timers;
    if (timers.empty()) {
      timers.push_back(new 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 Timer(name1.str()));
        ostringstream name2;
        name2 << "comm_state[" << timers.size() << "]"
              << "." << tostring(state) << ".step";
        timers.push_back(new Timer(name2.str()));
      }
    }

    vector<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 active_tl = CCTK_ActiveTimeLevelsGI (cctkGH, g);
        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_restrict_all (state, tl, reflevel, mglevel);
          }
        }
      } // loop over groups
      (*ti)->stop(); ++ti; (*ti)->start();
    } // for state
    (*ti)->stop();
    ++ti; assert(ti == timers.end());
  }

} // namespace Carpet