InterpolatorPUPMPI.h (8640B)
1 // File : InterpolatorPUPMPI.h 2 // Date : Tue Nov 29 21:00:01 2016 3 // Author : Fabian Wermelinger 4 // Description: pack/unpack kernels 5 // Copyright 2016 ETH Zurich. All Rights Reserved. 6 #ifndef INTERPOLATORPUPMPI_H_EPMRY0VT 7 #define INTERPOLATORPUPMPI_H_EPMRY0VT 8 9 /////////////////////////////////////////////////////////////////////////////// 10 // FACES 11 /////////////////////////////////////////////////////////////////////////////// 12 template <typename Tinterp> 13 void InterpolatorMPI<Tinterp>::_process_faceX(func_t& eval) 14 { 15 for (int iz=0; iz<this->m_data.Nz(); ++iz) 16 for (int iy=0; iy<this->m_data.Ny(); ++iy) 17 for (int ix=0; ix<Face_t::_P; ++ix) 18 eval(ix,iy,iz); 19 } 20 21 template <typename Tinterp> 22 void InterpolatorMPI<Tinterp>::_process_faceY(func_t& eval) 23 { 24 for (int iz=0; iz<this->m_data.Nz(); ++iz) 25 for (int iy=0; iy<Face_t::_P; ++iy) 26 for (int ix=0; ix<this->m_data.Nx(); ++ix) 27 eval(ix,iy,iz); 28 } 29 30 template <typename Tinterp> 31 void InterpolatorMPI<Tinterp>::_process_faceZ(func_t& eval) 32 { 33 for (int iz=0; iz<Face_t::_P; ++iz) 34 for (int iy=0; iy<this->m_data.Ny(); ++iy) 35 for (int ix=0; ix<this->m_data.Nx(); ++ix) 36 eval(ix,iy,iz); 37 } 38 39 template <typename Tinterp> 40 void InterpolatorMPI<Tinterp>::_pack_faceX(Face_t& f, const int _s0) 41 { 42 func_t pack = [this, &f, _s0](const int ix, const int iy, const int iz) { f(iy,iz,ix) = this->m_data(ix+_s0,iy,iz); }; 43 _process_faceX(pack); 44 } 45 46 template <typename Tinterp> 47 void InterpolatorMPI<Tinterp>::_pack_faceY(Face_t& f, const int _s0) 48 { 49 func_t pack = [this, &f, _s0](const int ix, const int iy, const int iz) { f(ix,iz,iy) = this->m_data(ix,iy+_s0,iz); }; 50 _process_faceY(pack); 51 } 52 53 template <typename Tinterp> 54 void InterpolatorMPI<Tinterp>::_pack_faceZ(Face_t& f, const int _s0) 55 { 56 func_t pack = [this, &f, _s0](const int ix, const int iy, const int iz) { f(ix,iy,iz) = this->m_data(ix,iy,iz+_s0); }; 57 _process_faceZ(pack); 58 } 59 60 template <typename Tinterp> 61 void InterpolatorMPI<Tinterp>::_upack_faceX(const Face_t& f, const int _d0) 62 { 63 func_t unpack = [this, &f, _d0](const int ix, const int iy, const int iz) { this->m_data(ix+_d0,iy,iz) = f(iy,iz,ix); }; 64 _process_faceX(unpack); 65 } 66 67 template <typename Tinterp> 68 void InterpolatorMPI<Tinterp>::_upack_faceY(const Face_t& f, const int _d0) 69 { 70 func_t unpack = [this, &f, _d0](const int ix, const int iy, const int iz) { this->m_data(ix,iy+_d0,iz) = f(ix,iz,iy); }; 71 _process_faceY(unpack); 72 } 73 74 template <typename Tinterp> 75 void InterpolatorMPI<Tinterp>::_upack_faceZ(const Face_t& f, const int _d0) 76 { 77 func_t unpack = [this, &f, _d0](const int ix, const int iy, const int iz) { this->m_data(ix,iy,iz+_d0) = f(ix,iy,iz); }; 78 _process_faceZ(unpack); 79 } 80 81 template <typename Tinterp> 82 void InterpolatorMPI<Tinterp>::_halo_faceX(const int _s0, const int _d0) 83 { 84 func_t halo = [this, _s0, _d0](const int ix, const int iy, const int iz) { this->m_data(ix+_d0,iy,iz) = this->m_data(_s0,iy,iz); }; 85 _process_faceX(halo); 86 } 87 88 template <typename Tinterp> 89 void InterpolatorMPI<Tinterp>::_halo_faceY(const int _s0, const int _d0) 90 { 91 func_t halo = [this, _s0, _d0](const int ix, const int iy, const int iz) { this->m_data(ix,iy+_d0,iz) = this->m_data(ix,_s0,iz); }; 92 _process_faceY(halo); 93 } 94 95 template <typename Tinterp> 96 void InterpolatorMPI<Tinterp>::_halo_faceZ(const int _s0, const int _d0) 97 { 98 func_t halo = [this, _s0, _d0](const int ix, const int iy, const int iz) { this->m_data(ix,iy,iz+_d0) = this->m_data(ix,iy,_s0); }; 99 _process_faceZ(halo); 100 } 101 102 /////////////////////////////////////////////////////////////////////////////// 103 // EDGES 104 /////////////////////////////////////////////////////////////////////////////// 105 template <typename Tinterp> 106 void InterpolatorMPI<Tinterp>::_process_edgeX(func_t& eval) 107 { 108 for (int iz=0; iz<Edge_t::_Q; ++iz) 109 for (int iy=0; iy<Edge_t::_P; ++iy) 110 for (int ix=0; ix<this->m_data.Nx(); ++ix) 111 eval(ix,iy,iz); 112 } 113 114 template <typename Tinterp> 115 void InterpolatorMPI<Tinterp>::_process_edgeY(func_t& eval) 116 { 117 for (int iz=0; iz<Edge_t::_Q; ++iz) 118 for (int iy=0; iy<this->m_data.Ny(); ++iy) 119 for (int ix=0; ix<Edge_t::_P; ++ix) 120 eval(ix,iy,iz); 121 } 122 123 template <typename Tinterp> 124 void InterpolatorMPI<Tinterp>::_process_edgeZ(func_t& eval) 125 { 126 for (int iz=0; iz<this->m_data.Nz(); ++iz) 127 for (int iy=0; iy<Edge_t::_Q; ++iy) 128 for (int ix=0; ix<Edge_t::_P; ++ix) 129 eval(ix,iy,iz); 130 } 131 132 template <typename Tinterp> 133 void InterpolatorMPI<Tinterp>::_pack_edgeX(Edge_t& e, const int _s0, const int _s1) 134 { 135 func_t pack = [this, &e, _s0, _s1](const int ix, const int iy, const int iz) { e(ix,iy,iz) = this->m_data(ix,iy+_s0,iz+_s1); }; 136 _process_edgeX(pack); 137 } 138 139 template <typename Tinterp> 140 void InterpolatorMPI<Tinterp>::_pack_edgeY(Edge_t& e, const int _s0, const int _s1) 141 { 142 func_t pack = [this, &e, _s0, _s1](const int ix, const int iy, const int iz) { e(iy,ix,iz) = this->m_data(ix+_s0,iy,iz+_s1); }; 143 _process_edgeY(pack); 144 } 145 146 template <typename Tinterp> 147 void InterpolatorMPI<Tinterp>::_pack_edgeZ(Edge_t& e, const int _s0, const int _s1) 148 { 149 func_t pack = [this, &e, _s0, _s1](const int ix, const int iy, const int iz) { e(iz,ix,iy) = this->m_data(ix+_s0,iy+_s1,iz); }; 150 _process_edgeZ(pack); 151 } 152 153 template <typename Tinterp> 154 void InterpolatorMPI<Tinterp>::_upack_edgeX(const Edge_t& e, const int _d0, const int _d1) 155 { 156 func_t unpack = [this, &e, _d0, _d1](const int ix, const int iy, const int iz) { this->m_data(ix,iy+_d0,iz+_d1) = e(ix,iy,iz); }; 157 _process_edgeX(unpack); 158 } 159 160 template <typename Tinterp> 161 void InterpolatorMPI<Tinterp>::_upack_edgeY(const Edge_t& e, const int _d0, const int _d1) 162 { 163 func_t unpack = [this, &e, _d0, _d1](const int ix, const int iy, const int iz) { this->m_data(ix+_d0,iy,iz+_d1) = e(iy,ix,iz); }; 164 _process_edgeY(unpack); 165 } 166 167 template <typename Tinterp> 168 void InterpolatorMPI<Tinterp>::_upack_edgeZ(const Edge_t& e, const int _d0, const int _d1) 169 { 170 func_t unpack = [this, &e, _d0, _d1](const int ix, const int iy, const int iz) { this->m_data(ix+_d0,iy+_d1,iz) = e(iz,ix,iy); }; 171 _process_edgeZ(unpack); 172 } 173 174 template <typename Tinterp> 175 void InterpolatorMPI<Tinterp>::_halo_edgeX(const int _s0, const int _s1, const int _d0, const int _d1) 176 { 177 func_t halo = [this, _s0, _s1, _d0, _d1](const int ix, const int iy, const int iz) { this->m_data(ix,iy+_d0,iz+_d1) = this->m_data(ix,_s0,_s1); }; 178 _process_edgeX(halo); 179 } 180 181 template <typename Tinterp> 182 void InterpolatorMPI<Tinterp>::_halo_edgeY(const int _s0, const int _s1, const int _d0, const int _d1) 183 { 184 func_t halo = [this, _s0, _s1, _d0, _d1](const int ix, const int iy, const int iz) { this->m_data(ix+_d0,iy,iz+_d1) = this->m_data(_s0,iy,_s1); }; 185 _process_edgeY(halo); 186 } 187 188 template <typename Tinterp> 189 void InterpolatorMPI<Tinterp>::_halo_edgeZ(const int _s0, const int _s1, const int _d0, const int _d1) 190 { 191 func_t halo = [this, _s0, _s1, _d0, _d1](const int ix, const int iy, const int iz) { this->m_data(ix+_d0,iy+_d1,iz) = this->m_data(_s0,_s1,iz); }; 192 _process_edgeZ(halo); 193 } 194 195 /////////////////////////////////////////////////////////////////////////////// 196 // CORNERS 197 /////////////////////////////////////////////////////////////////////////////// 198 template <typename Tinterp> 199 void InterpolatorMPI<Tinterp>::_process_corner(func_t& eval) 200 { 201 for (int iz=0; iz<Corner_t::_R; ++iz) 202 for (int iy=0; iy<Corner_t::_Q; ++iy) 203 for (int ix=0; ix<Corner_t::_P; ++ix) 204 eval(ix,iy,iz); 205 } 206 207 template <typename Tinterp> 208 void InterpolatorMPI<Tinterp>::_pack_corner(Corner_t& c, const int _s0, const int _s1, const int _s2) 209 { 210 func_t pack = [this, &c, _s0, _s1, _s2](const int ix, const int iy, const int iz) { c(ix,iy,iz) = this->m_data(ix+_s0,iy+_s1,iz+_s2); }; 211 _process_corner(pack); 212 } 213 214 template <typename Tinterp> 215 void InterpolatorMPI<Tinterp>::_upack_corner(const Corner_t& c, const int _d0, const int _d1, const int _d2) 216 { 217 func_t unpack = [this, &c, _d0, _d1, _d2](const int ix, const int iy, const int iz) { this->m_data(ix+_d0,iy+_d1,iz+_d2) = c(ix,iy,iz); }; 218 _process_corner(unpack); 219 } 220 221 template <typename Tinterp> 222 void InterpolatorMPI<Tinterp>::_halo_corner(const int _s0, const int _s1, const int _s2, const int _d0, const int _d1, const int _d2) 223 { 224 func_t halo = [this, _s0, _s1, _s2, _d0, _d1, _d2](const int ix, const int iy, const int iz) { this->m_data(ix+_d0,iy+_d1,iz+_d2) = this->m_data(_s0,_s1,_s2); }; 225 _process_corner(halo); 226 } 227 228 #endif /* INTERPOLATORPUPMPI_H_EPMRY0VT */