aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetLib/src/restrict_3d_rf2.cc
blob: 0885479d811f6a7dee045791a5e59313425d94e5 (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
#include <cctk.h>
#include <cctk_Parameters.h>

#include <loopcontrol.h>

#include <algorithm>
#include <cassert>
#include <cmath>
#include <cstdlib>
#include <iostream>

#include "operator_prototypes_3d.hh"
#include "typeprops.hh"

using namespace std;



namespace CarpetLib {


  
#define SRCIND3(i,j,k)                                  \
  index3 (srcioff + (i), srcjoff + (j), srckoff + (k),  \
          srcipadext, srcjpadext, srckpadext,           \
          srciext, srcjext, srckext)
#define DSTIND3(i,j,k)                                  \
  index3 (dstioff + (i), dstjoff + (j), dstkoff + (k),  \
          dstipadext, dstjpadext, dstkpadext,           \
          dstiext, dstjext, dstkext)
  
  
  
  template <typename T>
  void
  restrict_3d_rf2 (T const * restrict const src,
                   ivect3 const & restrict srcpadext,
                   ivect3 const & restrict srcext,
                   T * restrict const dst,
                   ivect3 const & restrict dstpadext,
                   ivect3 const & restrict dstext,
                   ibbox3 const & restrict srcbbox,
                   ibbox3 const & restrict dstbbox,
                   ibbox3 const & restrict,
                   ibbox3 const & restrict regbbox,
                   void * extraargs)
  {
    DECLARE_CCTK_PARAMETERS;
    
    assert (not extraargs);
    
    if (any (srcbbox.stride() >= regbbox.stride() or
             dstbbox.stride() != regbbox.stride()))
    {
      CCTK_WARN (0, "Internal error: strides disagree");
    }
    
    if (any (reffact2 * srcbbox.stride() != dstbbox.stride())) {
      CCTK_WARN (0, "Internal error: destination strides are not twice the source strides");
    }
    
    // This could be handled, but is likely to point to an error
    // elsewhere
    if (regbbox.empty()) {
      CCTK_WARN (0, "Internal error: region extent is empty");
    }
    
    if (not regbbox.is_contained_in(srcbbox) or
        not regbbox.is_contained_in(dstbbox))
    {
      cerr << "srcbbox: " << srcbbox << endl
           << "dstbbox: " << dstbbox << endl
           << "regbbox: " << regbbox << endl;
      CCTK_WARN (0, "Internal error: region extent is not contained in array extent");
    }
    
    
    
    ivect3 const regext = regbbox.shape() / regbbox.stride();
    assert (all ((regbbox.lower() - srcbbox.lower()) % srcbbox.stride() == 0));
    ivect3 const srcoff =
      (regbbox.lower() - srcbbox.lower()) / srcbbox.stride();
    assert (all ((regbbox.lower() - dstbbox.lower()) % dstbbox.stride() == 0));
    ivect3 const dstoff =
      (regbbox.lower() - dstbbox.lower()) / dstbbox.stride();
    
    
    
    int const srcipadext = srcpadext[0];
    int const srcjpadext = srcpadext[1];
    int const srckpadext = srcpadext[2];
    
    int const dstipadext = dstpadext[0];
    int const dstjpadext = dstpadext[1];
    int const dstkpadext = dstpadext[2];
    
    int const srciext = srcext[0];
    int const srcjext = srcext[1];
    int const srckext = srcext[2];
    
    int const dstiext = dstext[0];
    int const dstjext = dstext[1];
    int const dstkext = dstext[2];
    
    int const regiext = regext[0];
    int const regjext = regext[1];
    int const regkext = regext[2];
    
    int const srcioff = srcoff[0];
    int const srcjoff = srcoff[1];
    int const srckoff = srcoff[2];
    
    int const dstioff = dstoff[0];
    int const dstjoff = dstoff[1];
    int const dstkoff = dstoff[2];
    
    
    
    if (not use_loopcontrol_in_operators) {
      
      // Loop over coarse region
      for (int k=0; k<regkext; ++k) {
        for (int j=0; j<regjext; ++j) {
          for (int i=0; i<regiext; ++i) {
            
            dst [DSTIND3(i, j, k)] = src [SRCIND3(2*i, 2*j, 2*k)];
            
          }
        }
      }
      
    } else {
      
      // Loop over coarse region
#pragma omp parallel
      CCTK_LOOP3(restrict_3d_rf2,
                 i,j,k, 0,0,0, regiext,regjext,regkext,
                 dstipadext,dstjpadext,dstkpadext)
      {
        
        dst [DSTIND3(i, j, k)] = src [SRCIND3(2*i, 2*j, 2*k)];
        
      } CCTK_ENDLOOP3(restrict_3d_rf2);
      
    }
    
  }
  
  
  
#define TYPECASE(N,T)                                   \
  template                                              \
  void                                                  \
  restrict_3d_rf2 (T const * restrict const src,        \
                   ivect3 const & restrict srcpadext,   \
                   ivect3 const & restrict srcext,      \
                   T * restrict const dst,              \
                   ivect3 const & restrict dstpadext,   \
                   ivect3 const & restrict dstext,      \
                   ibbox3 const & restrict srcbbox,     \
                   ibbox3 const & restrict dstbbox,     \
                   ibbox3 const & restrict,             \
                   ibbox3 const & restrict regbbox,     \
                   void * extraargs);
#include "typecase.hh"
#undef TYPECASE
  
  
  
} // namespace CarpetLib