1 /* 2 Implements the Kokkos kernel 3 */ 4 #include <petscconf.h> 5 #include <petscvec_kokkos.hpp> 6 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 7 #include <petsc/private/deviceimpl.h> 8 #include <petsclandau.h> 9 #include <petscts.h> 10 11 #include <Kokkos_Core.hpp> 12 #include <cstdio> 13 typedef Kokkos::TeamPolicy<>::member_type team_member; 14 #include "../land_tensors.h" 15 #include <petscaijdevice.h> 16 17 namespace landau_inner_red { // namespace helps with name resolution in reduction identity 18 template <class ScalarType> 19 struct array_type { 20 ScalarType gg2[LANDAU_DIM]; 21 ScalarType gg3[LANDAU_DIM][LANDAU_DIM]; 22 23 KOKKOS_INLINE_FUNCTION // Default constructor - Initialize to 0's 24 array_type() { 25 for (int j = 0; j < LANDAU_DIM; j++) { 26 gg2[j] = 0; 27 for (int k = 0; k < LANDAU_DIM; k++) gg3[j][k] = 0; 28 } 29 } 30 KOKKOS_INLINE_FUNCTION // Copy Constructor 31 array_type(const array_type &rhs) { 32 for (int j = 0; j < LANDAU_DIM; j++) { 33 gg2[j] = rhs.gg2[j]; 34 for (int k = 0; k < LANDAU_DIM; k++) gg3[j][k] = rhs.gg3[j][k]; 35 } 36 } 37 KOKKOS_INLINE_FUNCTION // add operator 38 array_type & 39 operator+=(const array_type &src) { 40 for (int j = 0; j < LANDAU_DIM; j++) { 41 gg2[j] += src.gg2[j]; 42 for (int k = 0; k < LANDAU_DIM; k++) gg3[j][k] += src.gg3[j][k]; 43 } 44 return *this; 45 } 46 KOKKOS_INLINE_FUNCTION // volatile add operator 47 void 48 operator+=(const volatile array_type &src) volatile { 49 for (int j = 0; j < LANDAU_DIM; j++) { 50 gg2[j] += src.gg2[j]; 51 for (int k = 0; k < LANDAU_DIM; k++) gg3[j][k] += src.gg3[j][k]; 52 } 53 } 54 }; 55 typedef array_type<PetscReal> TensorValueType; // used to simplify code below 56 } // namespace landau_inner_red 57 58 namespace Kokkos { //reduction identity must be defined in Kokkos namespace 59 template <> 60 struct reduction_identity<landau_inner_red::TensorValueType> { 61 KOKKOS_FORCEINLINE_FUNCTION static landau_inner_red::TensorValueType sum() { return landau_inner_red::TensorValueType(); } 62 }; 63 } // namespace Kokkos 64 65 extern "C" { 66 PetscErrorCode LandauKokkosCreateMatMaps(P4estVertexMaps maps[], pointInterpolationP4est (*pointMaps)[LANDAU_MAX_Q_FACE], PetscInt Nf[], PetscInt Nq, PetscInt grid) { 67 P4estVertexMaps h_maps; /* host container */ 68 const Kokkos::View<pointInterpolationP4est *[LANDAU_MAX_Q_FACE], Kokkos::LayoutRight, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_points((pointInterpolationP4est *)pointMaps, maps[grid].num_reduced); 69 const Kokkos::View<LandauIdx *[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ], Kokkos::LayoutRight, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_gidx((LandauIdx *)maps[grid].gIdx, maps[grid].num_elements); 70 Kokkos::View<pointInterpolationP4est *[LANDAU_MAX_Q_FACE], Kokkos::LayoutRight> *d_points = new Kokkos::View<pointInterpolationP4est *[LANDAU_MAX_Q_FACE], Kokkos::LayoutRight>("points", maps[grid].num_reduced); 71 Kokkos::View<LandauIdx *[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ], Kokkos::LayoutRight> *d_gidx = new Kokkos::View<LandauIdx *[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ], Kokkos::LayoutRight>("gIdx", maps[grid].num_elements); 72 73 PetscFunctionBegin; 74 Kokkos::deep_copy(*d_gidx, h_gidx); 75 Kokkos::deep_copy(*d_points, h_points); 76 h_maps.num_elements = maps[grid].num_elements; 77 h_maps.num_face = maps[grid].num_face; 78 h_maps.num_reduced = maps[grid].num_reduced; 79 h_maps.deviceType = maps[grid].deviceType; 80 h_maps.numgrids = maps[grid].numgrids; 81 h_maps.Nf = Nf[grid]; 82 h_maps.c_maps = (pointInterpolationP4est(*)[LANDAU_MAX_Q_FACE])d_points->data(); 83 maps[grid].vp1 = (void *)d_points; 84 h_maps.gIdx = (LandauIdx(*)[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ])d_gidx->data(); 85 maps[grid].vp2 = (void *)d_gidx; 86 { 87 Kokkos::View<P4estVertexMaps, Kokkos::HostSpace> h_maps_k(&h_maps); 88 Kokkos::View<P4estVertexMaps> *d_maps_k = new Kokkos::View<P4estVertexMaps>(Kokkos::create_mirror(Kokkos::DefaultExecutionSpace::memory_space(), h_maps_k)); 89 Kokkos::deep_copy(*d_maps_k, h_maps_k); 90 maps[grid].d_self = d_maps_k->data(); 91 maps[grid].vp3 = (void *)d_maps_k; 92 } 93 PetscFunctionReturn(0); 94 } 95 PetscErrorCode LandauKokkosDestroyMatMaps(P4estVertexMaps maps[], PetscInt num_grids) { 96 PetscFunctionBegin; 97 for (PetscInt grid = 0; grid < num_grids; grid++) { 98 auto a = static_cast<Kokkos::View<pointInterpolationP4est *[LANDAU_MAX_Q_FACE], Kokkos::LayoutRight> *>(maps[grid].vp1); 99 auto b = static_cast<Kokkos::View<LandauIdx *[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ], Kokkos::LayoutRight> *>(maps[grid].vp2); 100 auto c = static_cast<Kokkos::View<P4estVertexMaps *> *>(maps[grid].vp3); 101 delete a; 102 delete b; 103 delete c; 104 } 105 PetscFunctionReturn(0); 106 } 107 108 PetscErrorCode LandauKokkosStaticDataSet(DM plex, const PetscInt Nq, const PetscInt batch_sz, const PetscInt num_grids, PetscInt a_numCells[], PetscInt a_species_offset[], PetscInt a_mat_offset[], PetscReal a_nu_alpha[], PetscReal a_nu_beta[], PetscReal a_invMass[], PetscReal a_invJ[], PetscReal a_x[], PetscReal a_y[], PetscReal a_z[], PetscReal a_w[], LandauStaticData *SData_d) { 109 PetscReal *BB, *DD; 110 PetscTabulation *Tf; 111 PetscInt dim; 112 PetscInt Nb = Nq, ip_offset[LANDAU_MAX_GRIDS + 1], ipf_offset[LANDAU_MAX_GRIDS + 1], elem_offset[LANDAU_MAX_GRIDS + 1], nip, IPf_sz, Nftot; 113 PetscDS prob; 114 115 PetscFunctionBegin; 116 PetscCall(DMGetDimension(plex, &dim)); 117 PetscCall(DMGetDS(plex, &prob)); 118 PetscCheck(LANDAU_DIM == dim, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dim %" PetscInt_FMT " != LANDAU_DIM %d", dim, LANDAU_DIM); 119 PetscCall(PetscDSGetTabulation(prob, &Tf)); 120 BB = Tf[0]->T[0]; 121 DD = Tf[0]->T[1]; 122 ip_offset[0] = ipf_offset[0] = elem_offset[0] = 0; 123 nip = 0; 124 IPf_sz = 0; 125 for (PetscInt grid = 0; grid < num_grids; grid++) { 126 PetscInt nfloc = a_species_offset[grid + 1] - a_species_offset[grid]; 127 elem_offset[grid + 1] = elem_offset[grid] + a_numCells[grid]; 128 nip += a_numCells[grid] * Nq; 129 ip_offset[grid + 1] = nip; 130 IPf_sz += Nq * nfloc * a_numCells[grid]; 131 ipf_offset[grid + 1] = IPf_sz; 132 } 133 Nftot = a_species_offset[num_grids]; 134 PetscCall(PetscKokkosInitializeCheck()); 135 { 136 const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_alpha(a_nu_alpha, Nftot); 137 auto alpha = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("alpha", Nftot); 138 SData_d->alpha = static_cast<void *>(alpha); 139 const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_beta(a_nu_beta, Nftot); 140 auto beta = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("beta", Nftot); 141 SData_d->beta = static_cast<void *>(beta); 142 const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_invMass(a_invMass, Nftot); 143 auto invMass = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("invMass", Nftot); 144 SData_d->invMass = static_cast<void *>(invMass); 145 const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_BB(BB, Nq * Nb); 146 auto B = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("B", Nq * Nb); 147 SData_d->B = static_cast<void *>(B); 148 const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_DD(DD, Nq * Nb * dim); 149 auto D = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("D", Nq * Nb * dim); 150 SData_d->D = static_cast<void *>(D); 151 const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_invJ(a_invJ, nip * dim * dim); 152 auto invJ = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("invJ", nip * dim * dim); 153 SData_d->invJ = static_cast<void *>(invJ); 154 const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_x(a_x, nip); 155 auto x = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("x", nip); 156 SData_d->x = static_cast<void *>(x); 157 const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_y(a_y, nip); 158 auto y = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("y", nip); 159 SData_d->y = static_cast<void *>(y); 160 const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_w(a_w, nip); 161 auto w = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("w", nip); 162 SData_d->w = static_cast<void *>(w); 163 164 Kokkos::deep_copy(*alpha, h_alpha); 165 Kokkos::deep_copy(*beta, h_beta); 166 Kokkos::deep_copy(*invMass, h_invMass); 167 Kokkos::deep_copy(*B, h_BB); 168 Kokkos::deep_copy(*D, h_DD); 169 Kokkos::deep_copy(*invJ, h_invJ); 170 Kokkos::deep_copy(*x, h_x); 171 Kokkos::deep_copy(*y, h_y); 172 Kokkos::deep_copy(*w, h_w); 173 174 if (dim == 3) { 175 const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_z(a_z, nip); 176 auto z = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("z", nip); 177 SData_d->z = static_cast<void *>(z); 178 Kokkos::deep_copy(*z, h_z); 179 } else SData_d->z = NULL; 180 181 // 182 const Kokkos::View<PetscInt *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_NCells(a_numCells, num_grids); 183 auto nc = new Kokkos::View<PetscInt *, Kokkos::LayoutLeft>("NCells", num_grids); 184 SData_d->NCells = static_cast<void *>(nc); 185 Kokkos::deep_copy(*nc, h_NCells); 186 187 const Kokkos::View<PetscInt *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_species_offset(a_species_offset, num_grids + 1); 188 auto soff = new Kokkos::View<PetscInt *, Kokkos::LayoutLeft>("species_offset", num_grids + 1); 189 SData_d->species_offset = static_cast<void *>(soff); 190 Kokkos::deep_copy(*soff, h_species_offset); 191 192 const Kokkos::View<PetscInt *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_mat_offset(a_mat_offset, num_grids + 1); 193 auto moff = new Kokkos::View<PetscInt *, Kokkos::LayoutLeft>("mat_offset", num_grids + 1); 194 SData_d->mat_offset = static_cast<void *>(moff); 195 Kokkos::deep_copy(*moff, h_mat_offset); 196 197 const Kokkos::View<PetscInt *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_ip_offset(ip_offset, num_grids + 1); 198 auto ipoff = new Kokkos::View<PetscInt *, Kokkos::LayoutLeft>("ip_offset", num_grids + 1); 199 SData_d->ip_offset = static_cast<void *>(ipoff); 200 Kokkos::deep_copy(*ipoff, h_ip_offset); 201 202 const Kokkos::View<PetscInt *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_elem_offset(elem_offset, num_grids + 1); 203 auto eoff = new Kokkos::View<PetscInt *, Kokkos::LayoutLeft>("elem_offset", num_grids + 1); 204 SData_d->elem_offset = static_cast<void *>(eoff); 205 Kokkos::deep_copy(*eoff, h_elem_offset); 206 207 const Kokkos::View<PetscInt *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_ipf_offset(ipf_offset, num_grids + 1); 208 auto ipfoff = new Kokkos::View<PetscInt *, Kokkos::LayoutLeft>("ipf_offset", num_grids + 1); 209 SData_d->ipf_offset = static_cast<void *>(ipfoff); 210 Kokkos::deep_copy(*ipfoff, h_ipf_offset); 211 #if defined(LANDAU_LAYOUT_LEFT) // preallocate dynamic data, no copy 212 auto ipfdf_data = new Kokkos::View<PetscReal ***, Kokkos::LayoutLeft>("fdf", batch_sz, dim + 1, IPf_sz); 213 #else 214 auto ipfdf_data = new Kokkos::View<PetscReal ***, Kokkos::LayoutRight>("fdf", batch_sz, dim + 1, IPf_sz); 215 #endif 216 SData_d->ipfdf_data = static_cast<void *>(ipfdf_data); 217 auto Eq_m = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("Eq_m", Nftot); // allocate but do not set, same for whole batch 218 SData_d->Eq_m = static_cast<void *>(Eq_m); 219 const Kokkos::View<LandauIdx *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_coo_elem_offsets((LandauIdx *)SData_d->coo_elem_offsets, SData_d->coo_n_cellsTot + 1); 220 auto coo_elem_offsets = new Kokkos::View<LandauIdx *, Kokkos::LayoutLeft>("coo_elem_offsets", SData_d->coo_n_cellsTot + 1); 221 const Kokkos::View<LandauIdx *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_coo_elem_fullNb((LandauIdx *)SData_d->coo_elem_fullNb, SData_d->coo_n_cellsTot); 222 auto coo_elem_fullNb = new Kokkos::View<LandauIdx *, Kokkos::LayoutLeft>("coo_elem_offsets", SData_d->coo_n_cellsTot); 223 const Kokkos::View<LandauIdx *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_coo_elem_point_offsets((LandauIdx *)SData_d->coo_elem_point_offsets, SData_d->coo_n_cellsTot * (LANDAU_MAX_NQ + 1)); 224 auto coo_elem_point_offsets = new Kokkos::View<LandauIdx *, Kokkos::LayoutLeft>("coo_elem_point_offsets", SData_d->coo_n_cellsTot * (LANDAU_MAX_NQ + 1)); 225 // assign & copy 226 Kokkos::deep_copy(*coo_elem_offsets, h_coo_elem_offsets); 227 Kokkos::deep_copy(*coo_elem_fullNb, h_coo_elem_fullNb); 228 Kokkos::deep_copy(*coo_elem_point_offsets, h_coo_elem_point_offsets); 229 // need to free this now and use pointer space 230 PetscCall(PetscFree3(SData_d->coo_elem_offsets, SData_d->coo_elem_fullNb, SData_d->coo_elem_point_offsets)); 231 SData_d->coo_elem_offsets = static_cast<void *>(coo_elem_offsets); 232 SData_d->coo_elem_fullNb = static_cast<void *>(coo_elem_fullNb); 233 SData_d->coo_elem_point_offsets = static_cast<void *>(coo_elem_point_offsets); 234 auto coo_vals = new Kokkos::View<PetscScalar *, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace>("coo_vals", SData_d->coo_size); 235 SData_d->coo_vals = static_cast<void *>(coo_vals); 236 } 237 SData_d->maps = NULL; // not used 238 PetscFunctionReturn(0); 239 } 240 241 PetscErrorCode LandauKokkosStaticDataClear(LandauStaticData *SData_d) { 242 PetscFunctionBegin; 243 if (SData_d->alpha) { 244 auto alpha = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->alpha); 245 delete alpha; 246 SData_d->alpha = NULL; 247 auto beta = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->beta); 248 delete beta; 249 auto invMass = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->invMass); 250 delete invMass; 251 auto B = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->B); 252 delete B; 253 auto D = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->D); 254 delete D; 255 auto invJ = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->invJ); 256 delete invJ; 257 auto x = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->x); 258 delete x; 259 auto y = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->y); 260 delete y; 261 if (SData_d->z) { 262 auto z = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->z); 263 delete z; 264 } 265 #if defined(LANDAU_LAYOUT_LEFT) // preallocate dynamic data, no copy 266 auto z = static_cast<Kokkos::View<PetscReal ***, Kokkos::LayoutLeft> *>(SData_d->ipfdf_data); 267 #else 268 auto z = static_cast<Kokkos::View<PetscReal ***, Kokkos::LayoutRight> *>(SData_d->ipfdf_data); 269 #endif 270 delete z; 271 auto w = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->w); 272 delete w; 273 auto Eq_m = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->Eq_m); 274 delete Eq_m; 275 // offset 276 auto nc = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->NCells); 277 delete nc; 278 auto soff = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->species_offset); 279 delete soff; 280 auto moff = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->mat_offset); 281 delete moff; 282 auto ipoff = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->ip_offset); 283 delete ipoff; 284 auto eoff = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->elem_offset); 285 delete eoff; 286 auto ipfoff = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->ipf_offset); 287 delete ipfoff; 288 auto coo_elem_offsets = static_cast<Kokkos::View<LandauIdx *, Kokkos::LayoutLeft> *>((void *)SData_d->coo_elem_offsets); 289 delete coo_elem_offsets; 290 auto coo_elem_fullNb = static_cast<Kokkos::View<LandauIdx *, Kokkos::LayoutLeft> *>((void *)SData_d->coo_elem_fullNb); 291 delete coo_elem_fullNb; 292 auto coo_elem_point_offsets = static_cast<Kokkos::View<LandauIdx *, Kokkos::LayoutLeft> *>((void *)SData_d->coo_elem_point_offsets); 293 delete coo_elem_point_offsets; 294 SData_d->coo_elem_offsets = NULL; 295 SData_d->coo_elem_point_offsets = NULL; 296 SData_d->coo_elem_fullNb = NULL; 297 auto coo_vals = static_cast<Kokkos::View<PetscScalar *, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace> *>((void *)SData_d->coo_vals); 298 delete coo_vals; 299 } 300 PetscFunctionReturn(0); 301 } 302 303 #define KOKKOS_SHARED_LEVEL 0 // 0 is shared, 1 is global 304 305 KOKKOS_INLINE_FUNCTION 306 PetscErrorCode landau_mat_assemble(PetscSplitCSRDataStructure d_mat, PetscScalar *coo_vals, const PetscScalar Aij, const PetscInt f, const PetscInt g, const PetscInt Nb, PetscInt moffset, const PetscInt elem, const PetscInt fieldA, const P4estVertexMaps *d_maps, const LandauIdx coo_elem_offsets[], const LandauIdx coo_elem_fullNb[], const LandauIdx (*coo_elem_point_offsets)[LANDAU_MAX_NQ + 1], const PetscInt glb_elem_idx, const PetscInt bid_coo_sz_batch) { 307 PetscInt idx, q, nr, nc, d; 308 const LandauIdx *const Idxs = &d_maps->gIdx[elem][fieldA][0]; 309 PetscScalar row_scale[LANDAU_MAX_Q_FACE] = {0}, col_scale[LANDAU_MAX_Q_FACE] = {0}; 310 if (coo_vals) { // mirror (i,j) in CreateStaticGPUData 311 const int fullNb = coo_elem_fullNb[glb_elem_idx], fullNb2 = fullNb * fullNb; 312 idx = Idxs[f]; 313 if (idx >= 0) { 314 nr = 1; 315 row_scale[0] = 1.; 316 } else { 317 idx = -idx - 1; 318 for (q = 0, nr = 0; q < d_maps->num_face; q++, nr++) { 319 if (d_maps->c_maps[idx][q].gid < 0) break; 320 row_scale[q] = d_maps->c_maps[idx][q].scale; 321 } 322 } 323 idx = Idxs[g]; 324 if (idx >= 0) { 325 nc = 1; 326 col_scale[0] = 1.; 327 } else { 328 idx = -idx - 1; 329 nc = d_maps->num_face; 330 for (q = 0, nc = 0; q < d_maps->num_face; q++, nc++) { 331 if (d_maps->c_maps[idx][q].gid < 0) break; 332 col_scale[q] = d_maps->c_maps[idx][q].scale; 333 } 334 } 335 const int idx0 = bid_coo_sz_batch + coo_elem_offsets[glb_elem_idx] + fieldA * fullNb2 + fullNb * coo_elem_point_offsets[glb_elem_idx][f] + nr * coo_elem_point_offsets[glb_elem_idx][g]; 336 for (int q = 0, idx2 = idx0; q < nr; q++) { 337 for (int d = 0; d < nc; d++, idx2++) coo_vals[idx2] = row_scale[q] * col_scale[d] * Aij; 338 } 339 } else { 340 PetscScalar vals[LANDAU_MAX_Q_FACE * LANDAU_MAX_Q_FACE] = {0}; 341 PetscInt rows[LANDAU_MAX_Q_FACE], cols[LANDAU_MAX_Q_FACE]; 342 idx = Idxs[f]; 343 if (idx >= 0) { 344 nr = 1; 345 rows[0] = idx; 346 row_scale[0] = 1.; 347 } else { 348 idx = -idx - 1; 349 for (q = 0, nr = 0; q < d_maps->num_face; q++, nr++) { 350 if (d_maps->c_maps[idx][q].gid < 0) break; 351 rows[q] = d_maps->c_maps[idx][q].gid; 352 row_scale[q] = d_maps->c_maps[idx][q].scale; 353 } 354 } 355 idx = Idxs[g]; 356 if (idx >= 0) { 357 nc = 1; 358 cols[0] = idx; 359 col_scale[0] = 1.; 360 } else { 361 idx = -idx - 1; 362 for (q = 0, nc = 0; q < d_maps->num_face; q++, nc++) { 363 if (d_maps->c_maps[idx][q].gid < 0) break; 364 cols[q] = d_maps->c_maps[idx][q].gid; 365 col_scale[q] = d_maps->c_maps[idx][q].scale; 366 } 367 } 368 369 for (q = 0; q < nr; q++) rows[q] = rows[q] + moffset; 370 for (q = 0; q < nc; q++) cols[q] = cols[q] + moffset; 371 for (q = 0; q < nr; q++) { 372 for (d = 0; d < nc; d++) vals[q * nc + d] = row_scale[q] * col_scale[d] * Aij; 373 } 374 MatSetValuesDevice(d_mat, nr, rows, nc, cols, vals, ADD_VALUES); 375 } 376 return 0; 377 } 378 379 PetscErrorCode LandauKokkosJacobian(DM plex[], const PetscInt Nq, const PetscInt batch_sz, const PetscInt num_grids, const PetscInt a_numCells[], PetscReal a_Eq_m[], PetscScalar a_elem_closure[], const PetscScalar a_xarray[], const LandauStaticData *SData_d, const PetscReal shift, const PetscLogEvent events[], const PetscInt a_mat_offset[], const PetscInt a_species_offset[], Mat subJ[], Mat JacP) { 380 using scr_mem_t = Kokkos::DefaultExecutionSpace::scratch_memory_space; 381 using real2_scr_t = Kokkos::View<PetscScalar **, Kokkos::LayoutRight, scr_mem_t>; 382 using g2_scr_t = Kokkos::View<PetscReal ***, Kokkos::LayoutRight, scr_mem_t>; 383 using g3_scr_t = Kokkos::View<PetscReal ****, Kokkos::LayoutRight, scr_mem_t>; 384 PetscInt Nb = Nq, dim, num_cells_max, Nf_max, num_cells_batch; 385 int nfaces = 0, vector_size = 512 / Nq; 386 LandauCtx *ctx; 387 PetscReal *d_Eq_m = NULL; 388 PetscScalar *d_vertex_f = NULL; 389 P4estVertexMaps *maps[LANDAU_MAX_GRIDS]; // this gets captured 390 PetscSplitCSRDataStructure d_mat; 391 PetscContainer container; 392 const int conc = Kokkos::DefaultExecutionSpace().concurrency(), openmp = !!(conc < 1000), team_size = (openmp == 0) ? Nq : 1; 393 const PetscInt coo_sz_batch = SData_d->coo_size / batch_sz; // capture 394 auto d_alpha_k = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->alpha); //static data 395 const PetscReal *d_alpha = d_alpha_k->data(); 396 const PetscInt Nftot = d_alpha_k->size(); // total number of species 397 auto d_beta_k = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->beta); 398 const PetscReal *d_beta = d_beta_k->data(); 399 auto d_invMass_k = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->invMass); 400 const PetscReal *d_invMass = d_invMass_k->data(); 401 auto d_B_k = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->B); 402 const PetscReal *d_BB = d_B_k->data(); 403 auto d_D_k = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->D); 404 const PetscReal *d_DD = d_D_k->data(); 405 auto d_invJ_k = *static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->invJ); // use Kokkos vector in kernels 406 auto d_x_k = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->x); //static data 407 const PetscReal *d_x = d_x_k->data(); 408 auto d_y_k = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->y); //static data 409 const PetscReal *d_y = d_y_k->data(); 410 auto d_z_k = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->z); //static data 411 const PetscReal *d_z = (LANDAU_DIM == 3) ? d_z_k->data() : NULL; 412 auto d_w_k = *static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->w); //static data 413 const PetscReal *d_w = d_w_k.data(); 414 // grid offsets - single vertex grid data 415 auto d_numCells_k = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->NCells); 416 const PetscInt *d_numCells = d_numCells_k->data(); 417 auto d_species_offset_k = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->species_offset); 418 const PetscInt *d_species_offset = d_species_offset_k->data(); 419 auto d_mat_offset_k = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->mat_offset); 420 const PetscInt *d_mat_offset = d_mat_offset_k->data(); 421 auto d_ip_offset_k = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->ip_offset); 422 const PetscInt *d_ip_offset = d_ip_offset_k->data(); 423 auto d_ipf_offset_k = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->ipf_offset); 424 const PetscInt *d_ipf_offset = d_ipf_offset_k->data(); 425 auto d_elem_offset_k = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->elem_offset); 426 const PetscInt *d_elem_offset = d_elem_offset_k->data(); 427 #if defined(LANDAU_LAYOUT_LEFT) // preallocate dynamic data, including batched vertices 428 Kokkos::View<PetscReal ***, Kokkos::LayoutLeft> d_fdf_k = *static_cast<Kokkos::View<PetscReal ***, Kokkos::LayoutLeft> *>(SData_d->ipfdf_data); 429 #else 430 Kokkos::View<PetscReal ***, Kokkos::LayoutRight> d_fdf_k = *static_cast<Kokkos::View<PetscReal ***, Kokkos::LayoutRight> *>(SData_d->ipfdf_data); 431 #endif 432 auto d_Eq_m_k = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->Eq_m); // static storage, dynamci data - E(t), copy later, single vertex 433 // COO 434 auto d_coo_elem_offsets_k = static_cast<Kokkos::View<LandauIdx *, Kokkos::LayoutLeft> *>(SData_d->coo_elem_offsets); 435 LandauIdx *d_coo_elem_offsets = (SData_d->coo_size == 0) ? NULL : d_coo_elem_offsets_k->data(); 436 auto d_coo_elem_fullNb_k = static_cast<Kokkos::View<LandauIdx *, Kokkos::LayoutLeft> *>(SData_d->coo_elem_fullNb); 437 LandauIdx *d_coo_elem_fullNb = (SData_d->coo_size == 0) ? NULL : d_coo_elem_fullNb_k->data(); 438 auto d_coo_elem_point_offsets_k = static_cast<Kokkos::View<LandauIdx *, Kokkos::LayoutLeft> *>(SData_d->coo_elem_point_offsets); 439 LandauIdx(*d_coo_elem_point_offsets)[LANDAU_MAX_NQ + 1] = (SData_d->coo_size == 0) ? NULL : (LandauIdx(*)[LANDAU_MAX_NQ + 1]) d_coo_elem_point_offsets_k->data(); 440 auto d_coo_vals_k = static_cast<Kokkos::View<PetscScalar *, Kokkos::LayoutRight> *>(SData_d->coo_vals); 441 PetscScalar *d_coo_vals = (SData_d->coo_size == 0) ? NULL : d_coo_vals_k->data(); 442 443 PetscFunctionBegin; 444 while (vector_size & (vector_size - 1)) vector_size = vector_size & (vector_size - 1); 445 if (vector_size > 16) vector_size = 16; // printf("DEBUG\n"); 446 PetscCall(PetscLogEventBegin(events[3], 0, 0, 0, 0)); 447 PetscCall(DMGetApplicationContext(plex[0], &ctx)); 448 PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context"); 449 PetscCall(DMGetDimension(plex[0], &dim)); 450 PetscCheck(LANDAU_DIM == dim, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dim %" PetscInt_FMT " != LANDAU_DIM %d", dim, LANDAU_DIM); 451 if (ctx->gpu_assembly) { 452 PetscCall(PetscObjectQuery((PetscObject)JacP, "assembly_maps", (PetscObject *)&container)); 453 if (container) { 454 P4estVertexMaps *h_maps = NULL; 455 PetscCall(PetscContainerGetPointer(container, (void **)&h_maps)); 456 for (PetscInt grid = 0; grid < num_grids; grid++) { 457 if (h_maps[grid].d_self) { 458 maps[grid] = h_maps[grid].d_self; 459 nfaces = h_maps[grid].num_face; // nface=0 for CPU assembly 460 } else { 461 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "GPU assembly but no metadata in container"); 462 } 463 } 464 if (!d_coo_vals) { 465 // this does the setup the first time called 466 PetscCall(MatKokkosGetDeviceMatWrite(JacP, &d_mat)); 467 } else { 468 d_mat = NULL; 469 } 470 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "GPU assembly but no metadata in container"); 471 } else { 472 for (PetscInt grid = 0; grid < num_grids; grid++) maps[grid] = NULL; 473 nfaces = 0; 474 d_mat = NULL; 475 container = NULL; 476 } 477 num_cells_batch = Nf_max = num_cells_max = 0; 478 for (PetscInt grid = 0; grid < num_grids; grid++) { 479 int Nfloc = a_species_offset[grid + 1] - a_species_offset[grid]; 480 if (Nfloc > Nf_max) Nf_max = Nfloc; 481 if (a_numCells[grid] > num_cells_max) num_cells_max = a_numCells[grid]; 482 num_cells_batch += a_numCells[grid]; // we don't have a host element offset here (but in ctx) 483 } 484 const int elem_mat_num_cells_max_grid = container ? 0 : num_cells_max; 485 #ifdef LAND_SUPPORT_CPU_ASS 486 const int totDim_max = Nf_max * Nq, elem_mat_size_max = totDim_max * totDim_max; 487 Kokkos::View<PetscScalar ****, Kokkos::LayoutRight> d_elem_mats("element matrices", batch_sz, num_grids, elem_mat_num_cells_max_grid, elem_mat_size_max); // not used (cpu assembly) 488 #endif 489 const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_Eq_m_k(a_Eq_m, Nftot); 490 if (a_elem_closure || a_xarray) { 491 Kokkos::deep_copy(*d_Eq_m_k, h_Eq_m_k); 492 d_Eq_m = d_Eq_m_k->data(); 493 } else d_Eq_m = NULL; 494 PetscCall(PetscKokkosInitializeCheck()); 495 PetscCall(PetscLogEventEnd(events[3], 0, 0, 0, 0)); 496 if (a_elem_closure || a_xarray) { // Jacobian, create f & df 497 Kokkos::View<PetscScalar *, Kokkos::LayoutLeft> *d_vertex_f_k = NULL; 498 PetscCall(PetscLogEventBegin(events[1], 0, 0, 0, 0)); 499 if (a_elem_closure) { 500 PetscInt closure_sz = 0; // argh, don't have this on the host!!! 501 for (PetscInt grid = 0; grid < num_grids; grid++) { 502 PetscInt nfloc = a_species_offset[grid + 1] - a_species_offset[grid]; 503 closure_sz += Nq * nfloc * a_numCells[grid]; 504 } 505 closure_sz *= batch_sz; 506 d_vertex_f_k = new Kokkos::View<PetscScalar *, Kokkos::LayoutLeft>("closure", closure_sz); 507 const Kokkos::View<PetscScalar *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_closure_k(a_elem_closure, closure_sz); // Vertex data for each element 508 Kokkos::deep_copy(*d_vertex_f_k, h_closure_k); 509 d_vertex_f = d_vertex_f_k->data(); 510 } else { 511 d_vertex_f = (PetscScalar *)a_xarray; 512 } 513 PetscCall(PetscLogEventEnd(events[1], 0, 0, 0, 0)); 514 PetscCall(PetscLogEventBegin(events[8], 0, 0, 0, 0)); 515 PetscCall(PetscLogGpuTimeBegin()); 516 517 const int scr_bytes_fdf = real2_scr_t::shmem_size(Nf_max, Nb); 518 PetscCall(PetscInfo(plex[0], "df & f shared memory size: %d bytes in level, %d num cells total=%d, team size=%d, vector size=%d, #face=%d, Nf_max=%d\n", scr_bytes_fdf, KOKKOS_SHARED_LEVEL, num_cells_batch * batch_sz, team_size, vector_size, nfaces, Nf_max)); 519 Kokkos::parallel_for( 520 "f, df", Kokkos::TeamPolicy<>(num_cells_batch * batch_sz, team_size, /* Kokkos::AUTO */ vector_size).set_scratch_size(KOKKOS_SHARED_LEVEL, Kokkos::PerTeam(scr_bytes_fdf)), KOKKOS_LAMBDA(const team_member team) { 521 const PetscInt b_Nelem = d_elem_offset[num_grids], b_elem_idx = team.league_rank() % b_Nelem, b_id = team.league_rank() / b_Nelem, IPf_sz_glb = d_ipf_offset[num_grids]; 522 // find my grid 523 PetscInt grid = 0; 524 while (b_elem_idx >= d_elem_offset[grid + 1]) grid++; 525 { 526 const PetscInt loc_nip = d_numCells[grid] * Nq, loc_Nf = d_species_offset[grid + 1] - d_species_offset[grid], loc_elem = b_elem_idx - d_elem_offset[grid]; 527 const PetscInt moffset = LAND_MOFFSET(b_id, grid, batch_sz, num_grids, d_mat_offset); 528 { 529 real2_scr_t s_coef_k(team.team_scratch(KOKKOS_SHARED_LEVEL), Nf_max, Nb); 530 PetscScalar *coef; 531 const PetscReal *invJe = &d_invJ_k((d_ip_offset[grid] + loc_elem * Nq) * dim * dim); 532 // un pack IPData 533 if (!maps[grid]) { 534 coef = &d_vertex_f[b_id * IPf_sz_glb + d_ipf_offset[grid] + loc_elem * Nb * loc_Nf]; // closure and IP indexing are the same 535 } else { 536 coef = s_coef_k.data(); 537 Kokkos::parallel_for(Kokkos::TeamThreadRange(team, 0, (int)loc_Nf), [=](int f) { 538 //for (int f = 0; f < loc_Nf; ++f) { 539 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, 0, (int)Nb), [=](int b) { 540 //for (int b = 0; b < Nb; ++b) { 541 LandauIdx idx = maps[grid]->gIdx[loc_elem][f][b]; 542 if (idx >= 0) { 543 coef[f * Nb + b] = d_vertex_f[idx + moffset]; // xarray data, not IP, need raw array to deal with CPU assemble case (not used) 544 } else { 545 idx = -idx - 1; 546 coef[f * Nb + b] = 0; 547 for (int q = 0; q < maps[grid]->num_face; q++) { 548 PetscInt id = maps[grid]->c_maps[idx][q].gid; 549 PetscScalar scale = maps[grid]->c_maps[idx][q].scale; 550 if (id >= 0) coef[f * Nb + b] += scale * d_vertex_f[id + moffset]; 551 } 552 } 553 }); 554 }); 555 } 556 team.team_barrier(); 557 Kokkos::parallel_for(Kokkos::TeamThreadRange(team, 0, Nq), [=](int myQi) { 558 const PetscReal *const invJ = &invJe[myQi * dim * dim]; // b_elem_idx: batch element index 559 const PetscReal *Bq = &d_BB[myQi * Nb], *Dq = &d_DD[myQi * Nb * dim]; 560 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, 0, (int)loc_Nf), [=](int f) { 561 PetscInt b, e, d; 562 PetscReal refSpaceDer[LANDAU_DIM]; 563 const PetscInt idx = d_ipf_offset[grid] + f * loc_nip + loc_elem * Nq + myQi; 564 d_fdf_k(b_id, 0, idx) = 0.0; 565 for (d = 0; d < LANDAU_DIM; ++d) refSpaceDer[d] = 0.0; 566 for (b = 0; b < Nb; ++b) { 567 d_fdf_k(b_id, 0, idx) += Bq[b] * PetscRealPart(coef[f * Nb + b]); 568 for (d = 0; d < dim; ++d) refSpaceDer[d] += Dq[b * dim + d] * PetscRealPart(coef[f * Nb + b]); 569 } 570 for (d = 0; d < dim; ++d) { 571 for (e = 0, d_fdf_k(b_id, d + 1, idx) = 0.0; e < dim; ++e) d_fdf_k(b_id, d + 1, idx) += invJ[e * dim + d] * refSpaceDer[e]; 572 } 573 }); // f 574 }); // myQi 575 } 576 team.team_barrier(); 577 } // 'grid' scope 578 }); // global elems - fdf 579 Kokkos::fence(); 580 PetscCall(PetscLogGpuTimeEnd()); 581 PetscCall(PetscLogEventEnd(events[8], 0, 0, 0, 0)); 582 // Jacobian 583 int maximum_shared_mem_size = 64000; 584 585 if (PetscDeviceConfiguredFor_Internal(PETSC_DEVICE_DEFAULT)) { 586 // FIXME: remove the configuredfor check once the PETSC_DEVICE_HOST MR is merged in 587 PetscDevice device; 588 589 PetscCall(PetscDeviceGetDefault_Internal(&device)); 590 PetscCall(PetscDeviceGetAttribute(device, PETSC_DEVICE_ATTR_SIZE_T_SHARED_MEM_PER_BLOCK, &maximum_shared_mem_size)); 591 } 592 const int jac_scr_bytes = 2 * (g2_scr_t::shmem_size(dim, Nf_max, Nq) + g3_scr_t::shmem_size(dim, dim, Nf_max, Nq)); 593 const int jac_shared_level = (jac_scr_bytes > maximum_shared_mem_size) ? 1 : KOKKOS_SHARED_LEVEL; 594 auto jac_lambda = KOKKOS_LAMBDA(const team_member team) { 595 const PetscInt b_Nelem = d_elem_offset[num_grids], b_elem_idx = team.league_rank() % b_Nelem, b_id = team.league_rank() / b_Nelem; 596 // find my grid 597 PetscInt grid = 0; 598 while (b_elem_idx >= d_elem_offset[grid + 1]) grid++; 599 { 600 const PetscInt loc_Nf = d_species_offset[grid + 1] - d_species_offset[grid], loc_elem = b_elem_idx - d_elem_offset[grid]; 601 const PetscInt moffset = LAND_MOFFSET(b_id, grid, batch_sz, num_grids, d_mat_offset); 602 const PetscInt f_off = d_species_offset[grid]; 603 #ifdef LAND_SUPPORT_CPU_ASS 604 const PetscInt totDim = loc_Nf * Nq; 605 #endif 606 g2_scr_t g2(team.team_scratch(jac_shared_level), dim, loc_Nf, Nq); 607 g3_scr_t g3(team.team_scratch(jac_shared_level), dim, dim, loc_Nf, Nq); 608 g2_scr_t gg2(team.team_scratch(jac_shared_level), dim, loc_Nf, Nq); 609 g3_scr_t gg3(team.team_scratch(jac_shared_level), dim, dim, loc_Nf, Nq); 610 // get g2[] & g3[] 611 Kokkos::parallel_for(Kokkos::TeamThreadRange(team, 0, Nq), [=](int myQi) { 612 using Kokkos::parallel_reduce; 613 const PetscInt jpidx_glb = d_ip_offset[grid] + loc_elem * Nq + myQi; 614 const PetscReal *invJ = &d_invJ_k(jpidx_glb * dim * dim); 615 const PetscReal vj[3] = {d_x[jpidx_glb], d_y[jpidx_glb], d_z ? d_z[jpidx_glb] : 0}, wj = d_w[jpidx_glb]; 616 landau_inner_red::TensorValueType gg_temp; // reduce on part of gg2 and g33 for IP jpidx_g 617 Kokkos::parallel_reduce( 618 Kokkos::ThreadVectorRange(team, (int)d_ip_offset[num_grids]), 619 [=](const int &ipidx, landau_inner_red::TensorValueType &ggg) { 620 const PetscReal wi = d_w[ipidx], x = d_x[ipidx], y = d_y[ipidx]; 621 PetscReal temp1[3] = {0, 0, 0}, temp2 = 0; 622 PetscInt fieldA, d2, d3, f_off_r, grid_r, ipidx_g, nip_loc_r, loc_Nf_r; 623 #if LANDAU_DIM == 2 624 PetscReal Ud[2][2], Uk[2][2], mask = (PetscAbs(vj[0] - x) < 100 * PETSC_SQRT_MACHINE_EPSILON && PetscAbs(vj[1] - y) < 100 * PETSC_SQRT_MACHINE_EPSILON) ? 0. : 1.; 625 LandauTensor2D(vj, x, y, Ud, Uk, mask); 626 #else 627 PetscReal U[3][3], z = d_z[ipidx], mask = (PetscAbs(vj[0]-x) < 100*PETSC_SQRT_MACHINE_EPSILON && PetscAbs(vj[1]-y) < 100*PETSC_SQRT_MACHINE_EPSILON && PetscAbs(vj[2]-z) < 100*PETSC_SQRT_MACHINE_EPSILON) ? 0. : 1.; 628 LandauTensor3D(vj, x, y, z, U, mask); 629 #endif 630 grid_r = 0; 631 while (ipidx >= d_ip_offset[grid_r + 1]) grid_r++; // yuck search for grid 632 f_off_r = d_species_offset[grid_r]; 633 ipidx_g = ipidx - d_ip_offset[grid_r]; 634 nip_loc_r = d_numCells[grid_r] * Nq; 635 loc_Nf_r = d_species_offset[grid_r + 1] - d_species_offset[grid_r]; 636 for (fieldA = 0; fieldA < loc_Nf_r; ++fieldA) { 637 const PetscInt idx = d_ipf_offset[grid_r] + fieldA * nip_loc_r + ipidx_g; 638 temp1[0] += d_fdf_k(b_id, 1, idx) * d_beta[fieldA + f_off_r] * d_invMass[fieldA + f_off_r]; 639 temp1[1] += d_fdf_k(b_id, 2, idx) * d_beta[fieldA + f_off_r] * d_invMass[fieldA + f_off_r]; 640 #if LANDAU_DIM == 3 641 temp1[2] += d_fdf_k(b_id, 3, idx) * d_beta[fieldA + f_off_r] * d_invMass[fieldA + f_off_r]; 642 #endif 643 temp2 += d_fdf_k(b_id, 0, idx) * d_beta[fieldA + f_off_r]; 644 } 645 temp1[0] *= wi; 646 temp1[1] *= wi; 647 #if LANDAU_DIM == 3 648 temp1[2] *= wi; 649 #endif 650 temp2 *= wi; 651 #if LANDAU_DIM == 2 652 for (d2 = 0; d2 < 2; d2++) { 653 for (d3 = 0; d3 < 2; ++d3) { 654 /* K = U * grad(f): g2=e: i,A */ 655 ggg.gg2[d2] += Uk[d2][d3] * temp1[d3]; 656 /* D = -U * (I \kron (fx)): g3=f: i,j,A */ 657 ggg.gg3[d2][d3] += Ud[d2][d3] * temp2; 658 } 659 } 660 #else 661 for (d2 = 0; d2 < 3; ++d2) { 662 for (d3 = 0; d3 < 3; ++d3) { 663 /* K = U * grad(f): g2 = e: i,A */ 664 ggg.gg2[d2] += U[d2][d3]*temp1[d3]; 665 /* D = -U * (I \kron (fx)): g3 = f: i,j,A */ 666 ggg.gg3[d2][d3] += U[d2][d3]*temp2; 667 } 668 } 669 #endif 670 }, 671 Kokkos::Sum<landau_inner_red::TensorValueType>(gg_temp)); 672 // add alpha and put in gg2/3 673 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, (int)loc_Nf), [&](const int &fieldA) { 674 PetscInt d2, d3; 675 for (d2 = 0; d2 < dim; d2++) { 676 gg2(d2, fieldA, myQi) = gg_temp.gg2[d2] * d_alpha[fieldA + f_off]; 677 for (d3 = 0; d3 < dim; d3++) gg3(d2, d3, fieldA, myQi) = -gg_temp.gg3[d2][d3] * d_alpha[fieldA + f_off] * d_invMass[fieldA + f_off]; 678 } 679 }); 680 /* add electric field term once per IP */ 681 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, (int)loc_Nf), [&](const int &fieldA) { gg2(dim - 1, fieldA, myQi) += d_Eq_m[fieldA + f_off]; }); 682 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, (int)loc_Nf), [=](const int &fieldA) { 683 int d, d2, d3, dp; 684 /* Jacobian transform - g2, g3 - per thread (2D) */ 685 for (d = 0; d < dim; ++d) { 686 g2(d, fieldA, myQi) = 0; 687 for (d2 = 0; d2 < dim; ++d2) { 688 g2(d, fieldA, myQi) += invJ[d * dim + d2] * gg2(d2, fieldA, myQi); 689 g3(d, d2, fieldA, myQi) = 0; 690 for (d3 = 0; d3 < dim; ++d3) { 691 for (dp = 0; dp < dim; ++dp) g3(d, d2, fieldA, myQi) += invJ[d * dim + d3] * gg3(d3, dp, fieldA, myQi) * invJ[d2 * dim + dp]; 692 } 693 g3(d, d2, fieldA, myQi) *= wj; 694 } 695 g2(d, fieldA, myQi) *= wj; 696 } 697 }); 698 }); // Nq 699 team.team_barrier(); 700 { /* assemble */ 701 for (PetscInt fieldA = 0; fieldA < loc_Nf; fieldA++) { 702 /* assemble */ 703 Kokkos::parallel_for(Kokkos::TeamThreadRange(team, 0, Nb), [=](int f) { 704 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, 0, Nb), [=](int g) { 705 PetscScalar t = 0; 706 for (int qj = 0; qj < Nq; qj++) { // look at others integration points 707 const PetscReal *BJq = &d_BB[qj * Nb], *DIq = &d_DD[qj * Nb * dim]; 708 for (int d = 0; d < dim; ++d) { 709 t += DIq[f * dim + d] * g2(d, fieldA, qj) * BJq[g]; 710 for (int d2 = 0; d2 < dim; ++d2) t += DIq[f * dim + d] * g3(d, d2, fieldA, qj) * DIq[g * dim + d2]; 711 } 712 } 713 if (elem_mat_num_cells_max_grid) { // CPU assembly 714 #ifdef LAND_SUPPORT_CPU_ASS 715 const PetscInt fOff = (fieldA * Nb + f) * totDim + fieldA * Nb + g; 716 d_elem_mats(b_id, grid, loc_elem, fOff) = t; 717 #endif 718 } else { 719 landau_mat_assemble(d_mat, d_coo_vals, t, f, g, Nb, moffset, loc_elem, fieldA, maps[grid], d_coo_elem_offsets, d_coo_elem_fullNb, d_coo_elem_point_offsets, b_elem_idx, b_id * coo_sz_batch); 720 } 721 }); 722 }); 723 } 724 } 725 } // scope with 'grid' 726 }; 727 PetscCall(PetscLogEventBegin(events[4], 0, 0, 0, 0)); 728 PetscCall(PetscLogGpuTimeBegin()); 729 PetscCall(PetscInfo(plex[0], "Jacobian shared memory size: %d bytes in level %s (max shared=%d), num cells total=%d, team size=%d, vector size=%d, #face=%d, Nf_max=%d\n", jac_scr_bytes, jac_shared_level == 0 ? "local" : "global", maximum_shared_mem_size, num_cells_batch * batch_sz, team_size, vector_size, nfaces, Nf_max)); 730 Kokkos::parallel_for("Jacobian", Kokkos::TeamPolicy<>(num_cells_batch * batch_sz, team_size, vector_size).set_scratch_size(jac_shared_level, Kokkos::PerTeam(jac_scr_bytes)), jac_lambda); // Kokkos::LaunchBounds<512,2> 731 Kokkos::fence(); 732 PetscCall(PetscLogGpuTimeEnd()); 733 PetscCall(PetscLogEventEnd(events[4], 0, 0, 0, 0)); 734 if (d_vertex_f_k) delete d_vertex_f_k; 735 } else { // mass 736 PetscCall(PetscLogEventBegin(events[16], 0, 0, 0, 0)); 737 PetscCall(PetscLogGpuTimeBegin()); 738 PetscCall(PetscInfo(plex[0], "Mass team size=%d vector size=%d #face=%d Nb=%" PetscInt_FMT ", %s assembly\n", team_size, vector_size, nfaces, Nb, d_coo_vals ? "COO" : "CSR")); 739 Kokkos::parallel_for( 740 "Mass", Kokkos::TeamPolicy<>(num_cells_batch * batch_sz, team_size, vector_size), KOKKOS_LAMBDA(const team_member team) { // Kokkos::LaunchBounds<512,4> 741 const PetscInt b_Nelem = d_elem_offset[num_grids], b_elem_idx = team.league_rank() % b_Nelem, b_id = team.league_rank() / b_Nelem; 742 // find my grid 743 PetscInt grid = 0; 744 while (b_elem_idx >= d_elem_offset[grid + 1]) grid++; 745 { 746 const PetscInt loc_Nf = d_species_offset[grid + 1] - d_species_offset[grid], loc_elem = b_elem_idx - d_elem_offset[grid], jpidx_0 = d_ip_offset[grid] + loc_elem * Nq; 747 #ifdef LAND_SUPPORT_CPU_ASS 748 const PetscInt totDim = loc_Nf * Nq; 749 #endif 750 const PetscInt moffset = LAND_MOFFSET(b_id, grid, batch_sz, num_grids, d_mat_offset); 751 for (int fieldA = 0; fieldA < loc_Nf; fieldA++) { 752 /* assemble */ 753 Kokkos::parallel_for(Kokkos::TeamThreadRange(team, 0, Nb), [=](int f) { 754 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, 0, Nb), [=](int g) { 755 PetscScalar t = 0; 756 for (int qj = 0; qj < Nq; qj++) { // look at others integration points 757 const PetscReal *BJq = &d_BB[qj * Nb]; 758 const PetscInt jpidx_glb = jpidx_0 + qj; 759 if (dim == 2) { 760 t += BJq[f] * d_w_k(jpidx_glb) * shift * BJq[g] * 2. * PETSC_PI; 761 } else { 762 t += BJq[f] * d_w_k(jpidx_glb) * shift * BJq[g]; 763 } 764 } 765 if (elem_mat_num_cells_max_grid) { 766 #ifdef LAND_SUPPORT_CPU_ASS 767 const PetscInt fOff = (fieldA * Nb + f) * totDim + fieldA * Nb + g; 768 d_elem_mats(b_id, grid, loc_elem, fOff) = t; 769 #endif 770 } else { 771 landau_mat_assemble(d_mat, d_coo_vals, t, f, g, Nb, moffset, loc_elem, fieldA, maps[grid], d_coo_elem_offsets, d_coo_elem_fullNb, d_coo_elem_point_offsets, b_elem_idx, b_id * coo_sz_batch); 772 } 773 }); 774 }); 775 } // field 776 } // grid 777 }); 778 Kokkos::fence(); 779 PetscCall(PetscLogGpuTimeEnd()); 780 PetscCall(PetscLogEventEnd(events[16], 0, 0, 0, 0)); 781 } 782 if (d_coo_vals) { 783 PetscCall(MatSetValuesCOO(JacP, d_coo_vals, ADD_VALUES)); 784 #ifndef LAND_SUPPORT_CPU_ASS 785 Kokkos::fence(); // for timers 786 #endif 787 } else if (elem_mat_num_cells_max_grid) { // CPU assembly 788 #ifdef LAND_SUPPORT_CPU_ASS 789 Kokkos::View<PetscScalar ****, Kokkos::LayoutRight>::HostMirror h_elem_mats = Kokkos::create_mirror_view(d_elem_mats); 790 Kokkos::deep_copy(h_elem_mats, d_elem_mats); 791 PetscCheck(!container, PETSC_COMM_SELF, PETSC_ERR_PLIB, "?????"); 792 for (PetscInt b_id = 0; b_id < batch_sz; b_id++) { // OpenMP (once) 793 for (PetscInt grid = 0; grid < num_grids; grid++) { 794 PetscSection section, globalSection; 795 PetscInt cStart, cEnd; 796 Mat B = subJ[LAND_PACK_IDX(b_id, grid)]; 797 PetscInt moffset = LAND_MOFFSET(b_id, grid, batch_sz, num_grids, a_mat_offset), nloc, nzl, colbuf[1024], row; 798 const PetscInt *cols; 799 const PetscScalar *vals; 800 PetscCall(PetscLogEventBegin(events[5], 0, 0, 0, 0)); 801 PetscCall(DMPlexGetHeightStratum(plex[grid], 0, &cStart, &cEnd)); 802 PetscCall(DMGetLocalSection(plex[grid], §ion)); 803 PetscCall(DMGetGlobalSection(plex[grid], &globalSection)); 804 PetscCall(PetscLogEventEnd(events[5], 0, 0, 0, 0)); 805 PetscCall(PetscLogEventBegin(events[6], 0, 0, 0, 0)); 806 for (PetscInt ej = cStart; ej < cEnd; ++ej) { 807 const PetscScalar *elMat = &h_elem_mats(b_id, grid, ej - cStart, 0); 808 PetscCall(DMPlexMatSetClosure(plex[grid], section, globalSection, B, ej, elMat, ADD_VALUES)); 809 if (grid == 0 && ej == -1) { 810 const PetscInt loc_Nf = a_species_offset[grid + 1] - a_species_offset[grid], totDim = loc_Nf * Nq; 811 int d, f; 812 PetscPrintf(PETSC_COMM_SELF, "Kokkos Element matrix %d/%d\n", 1, (int)a_numCells[grid]); 813 for (d = 0; d < totDim; ++d) { 814 for (f = 0; f < totDim; ++f) PetscPrintf(PETSC_COMM_SELF, " %12.5e", PetscRealPart(elMat[d * totDim + f])); 815 PetscPrintf(PETSC_COMM_SELF, "\n"); 816 } 817 exit(14); 818 } 819 } 820 // move nest matrix to global JacP 821 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 822 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 823 PetscCall(MatGetSize(B, &nloc, NULL)); 824 for (int i = 0; i < nloc; i++) { 825 PetscCall(MatGetRow(B, i, &nzl, &cols, &vals)); 826 PetscCheck(nzl <= 1024, PetscObjectComm((PetscObject)B), PETSC_ERR_PLIB, "Row too big: %" PetscInt_FMT, nzl); 827 for (int j = 0; j < nzl; j++) colbuf[j] = cols[j] + moffset; 828 row = i + moffset; 829 PetscCall(MatSetValues(JacP, 1, &row, nzl, colbuf, vals, ADD_VALUES)); 830 PetscCall(MatRestoreRow(B, i, &nzl, &cols, &vals)); 831 } 832 PetscCall(MatDestroy(&subJ[LAND_PACK_IDX(b_id, grid)])); 833 PetscCall(PetscLogEventEnd(events[6], 0, 0, 0, 0)); 834 } // grids 835 } 836 #else 837 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "CPU assembly not supported"); 838 #endif 839 } 840 PetscFunctionReturn(0); 841 } 842 } // extern "C" 843