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 <petsclandau.h> 8 #include <petscts.h> 9 10 #include <Kokkos_Core.hpp> 11 #include <cstdio> 12 typedef Kokkos::TeamPolicy<>::member_type team_member; 13 #include "../land_tensors.h" 14 #include <petscaijdevice.h> 15 16 namespace landau_inner_red { // namespace helps with name resolution in reduction identity 17 template< class ScalarType > 18 struct array_type { 19 ScalarType gg2[LANDAU_DIM]; 20 ScalarType gg3[LANDAU_DIM][LANDAU_DIM]; 21 22 KOKKOS_INLINE_FUNCTION // Default constructor - Initialize to 0's 23 array_type() { 24 for (int j = 0; j < LANDAU_DIM; j++) { 25 gg2[j] = 0; 26 for (int k = 0; k < LANDAU_DIM; k++) { 27 gg3[j][k] = 0; 28 } 29 } 30 } 31 KOKKOS_INLINE_FUNCTION // Copy Constructor 32 array_type(const array_type & rhs) { 33 for (int j = 0; j < LANDAU_DIM; j++) { 34 gg2[j] = rhs.gg2[j]; 35 for (int k = 0; k < LANDAU_DIM; k++) { 36 gg3[j][k] = rhs.gg3[j][k]; 37 } 38 } 39 } 40 KOKKOS_INLINE_FUNCTION // add operator 41 array_type& operator += (const array_type& src) 42 { 43 for (int j = 0; j < LANDAU_DIM; j++) { 44 gg2[j] += src.gg2[j]; 45 for (int k = 0; k < LANDAU_DIM; k++) { 46 gg3[j][k] += src.gg3[j][k]; 47 } 48 } 49 return *this; 50 } 51 KOKKOS_INLINE_FUNCTION // volatile add operator 52 void operator += (const volatile array_type& src) volatile 53 { 54 for (int j = 0; j < LANDAU_DIM; j++) { 55 gg2[j] += src.gg2[j]; 56 for (int k = 0; k < LANDAU_DIM; k++) { 57 gg3[j][k] += src.gg3[j][k]; 58 } 59 } 60 } 61 }; 62 typedef array_type<PetscReal> TensorValueType; // used to simplify code below 63 } 64 65 namespace Kokkos { //reduction identity must be defined in Kokkos namespace 66 template<> 67 struct reduction_identity< landau_inner_red::TensorValueType > { 68 KOKKOS_FORCEINLINE_FUNCTION static landau_inner_red::TensorValueType sum() { 69 return landau_inner_red::TensorValueType(); 70 } 71 }; 72 } 73 74 extern "C" { 75 PetscErrorCode LandauKokkosCreateMatMaps(P4estVertexMaps *maps, pointInterpolationP4est (*pointMaps)[LANDAU_MAX_Q_FACE], PetscInt Nf, PetscInt Nq) 76 { 77 P4estVertexMaps h_maps; /* host container */ 78 const Kokkos::View<pointInterpolationP4est*[LANDAU_MAX_Q_FACE], Kokkos::LayoutRight, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_points ((pointInterpolationP4est*)pointMaps, maps->num_reduced); 79 const Kokkos::View< LandauIdx*[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ], Kokkos::LayoutRight, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_gidx ((LandauIdx*)maps->gIdx, maps->num_elements); 80 Kokkos::View<pointInterpolationP4est*[LANDAU_MAX_Q_FACE], Kokkos::LayoutRight> *d_points = new Kokkos::View<pointInterpolationP4est*[LANDAU_MAX_Q_FACE], Kokkos::LayoutRight>("points", maps->num_reduced); 81 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->num_elements); 82 83 PetscFunctionBegin; 84 Kokkos::deep_copy (*d_gidx, h_gidx); 85 Kokkos::deep_copy (*d_points, h_points); 86 h_maps.num_elements = maps->num_elements; 87 h_maps.num_face = maps->num_face; 88 h_maps.num_reduced = maps->num_reduced; 89 h_maps.deviceType = maps->deviceType; 90 h_maps.Nf = Nf; 91 h_maps.Nq = Nq; 92 h_maps.c_maps = (pointInterpolationP4est (*)[LANDAU_MAX_Q_FACE]) d_points->data(); 93 maps->vp1 = (void*)d_points; 94 h_maps.gIdx = (LandauIdx (*)[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ]) d_gidx->data(); 95 maps->vp2 = (void*)d_gidx; 96 { 97 Kokkos::View<P4estVertexMaps, Kokkos::HostSpace> h_maps_k(&h_maps); 98 Kokkos::View<P4estVertexMaps> *d_maps_k = new Kokkos::View<P4estVertexMaps>(Kokkos::create_mirror(Kokkos::DefaultExecutionSpace::memory_space(),h_maps_k)); 99 Kokkos::deep_copy (*d_maps_k, h_maps_k); 100 maps->data = d_maps_k->data(); 101 maps->vp3 = (void*)d_maps_k; 102 } 103 PetscFunctionReturn(0); 104 } 105 PetscErrorCode LandauKokkosDestroyMatMaps(P4estVertexMaps *maps) 106 { 107 Kokkos::View<pointInterpolationP4est*[LANDAU_MAX_Q_FACE], Kokkos::LayoutRight> *a = (Kokkos::View<pointInterpolationP4est*[LANDAU_MAX_Q_FACE], Kokkos::LayoutRight>*)maps->vp1; 108 Kokkos::View<LandauIdx*[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ], Kokkos::LayoutRight>*b = (Kokkos::View<LandauIdx*[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ], Kokkos::LayoutRight>*)maps->vp2; 109 Kokkos::View<P4estVertexMaps*> *c = (Kokkos::View<P4estVertexMaps*>*)maps->vp3; 110 delete a; delete b; delete c; 111 PetscFunctionReturn(0); 112 } 113 114 PetscErrorCode LandauKokkosStaticDataSet(DM plex, const PetscInt Nq, PetscReal nu_alpha[], PetscReal nu_beta[], PetscReal a_invMass[], PetscReal a_invJ[], PetscReal a_mass_w[], 115 PetscReal a_x[], PetscReal a_y[], PetscReal a_z[], PetscReal a_w[], LandauGeomData *SData_d) 116 { 117 PetscReal *BB,*DD; 118 PetscErrorCode ierr; 119 PetscTabulation *Tf; 120 LandauCtx *ctx; 121 PetscInt *Nbf,dim,Nf,Nb,nip,cStart,cEnd; 122 PetscDS prob; 123 124 PetscFunctionBegin; 125 ierr = DMGetApplicationContext(plex, &ctx);CHKERRQ(ierr); 126 if (!ctx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context"); 127 ierr = DMGetDimension(plex, &dim);CHKERRQ(ierr); 128 ierr = DMPlexGetHeightStratum(plex,0,&cStart,&cEnd);CHKERRQ(ierr); 129 nip = (cEnd - cStart)*Nq; 130 ierr = DMGetDS(plex, &prob);CHKERRQ(ierr); 131 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 132 ierr = PetscDSGetDimensions(prob, &Nbf);CHKERRQ(ierr); Nb = Nbf[0]; 133 if (Nq != Nb) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Nq != Nb. %D %D",Nq,Nb); 134 if (LANDAU_DIM != dim) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dim %D != LANDAU_DIM %d",dim,LANDAU_DIM); 135 ierr = PetscDSGetTabulation(prob, &Tf);CHKERRQ(ierr); 136 BB = Tf[0]->T[0]; DD = Tf[0]->T[1]; 137 ierr = PetscKokkosInitializeCheck();CHKERRQ(ierr); 138 { 139 const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_alpha (nu_alpha, Nf); 140 auto alpha = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("alpha", Nf); 141 SData_d->alpha = static_cast<void*>(alpha); 142 const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_beta (nu_beta, Nf); 143 auto beta = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("beta", Nf); 144 SData_d->beta = static_cast<void*>(beta); 145 const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_invMass (a_invMass,Nf); 146 auto invMass = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("invMass", Nf); 147 SData_d->invMass = static_cast<void*>(invMass); 148 const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_BB (BB,Nq*Nb); 149 auto B = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("B", Nq*Nb); 150 SData_d->B = static_cast<void*>(B); 151 const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_DD (DD,Nq*Nb*dim); 152 auto D = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("D", Nq*Nb*dim); 153 SData_d->D = static_cast<void*>(D); 154 const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_mass_w (a_mass_w, nip); 155 auto mass_w = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("mass_w", nip); 156 SData_d->mass_w = static_cast<void*>(mass_w); 157 const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_invJ (a_invJ, nip*dim*dim); 158 auto invJ = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("invJ", nip*dim*dim); 159 SData_d->invJ = static_cast<void*>(invJ); 160 const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_x (a_x, nip); 161 auto x = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("x", nip); 162 SData_d->x = static_cast<void*>(x); 163 const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_y (a_y, nip); 164 auto y = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("y", nip); 165 SData_d->y = static_cast<void*>(y); 166 const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_w (a_w, nip); 167 auto w = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("w", nip); 168 SData_d->w = static_cast<void*>(w); 169 const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_z (a_z , dim==3 ? nip : 0); 170 auto z = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("z", dim==3 ? nip : 0); 171 SData_d->z = static_cast<void*>(z); 172 173 Kokkos::deep_copy (*mass_w, h_mass_w); 174 Kokkos::deep_copy (*alpha, h_alpha); 175 Kokkos::deep_copy (*beta, h_beta); 176 Kokkos::deep_copy (*invMass, h_invMass); 177 Kokkos::deep_copy (*B, h_BB); 178 Kokkos::deep_copy (*D, h_DD); 179 Kokkos::deep_copy (*invJ, h_invJ); 180 Kokkos::deep_copy (*x, h_x); 181 Kokkos::deep_copy (*y, h_y); 182 Kokkos::deep_copy (*z, h_z); 183 Kokkos::deep_copy (*w, h_w); 184 185 auto Eq_m = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("Eq_m",Nf); 186 SData_d->Eq_m = static_cast<void*>(Eq_m); 187 auto IPf = new Kokkos::View<PetscScalar*, Kokkos::LayoutLeft> ("IPf",nip*Nf); // Nq==Nb 188 SData_d->IPf = static_cast<void*>(IPf); 189 } 190 PetscFunctionReturn(0); 191 } 192 193 PetscErrorCode LandauKokkosStaticDataClear(LandauGeomData *SData_d) 194 { 195 PetscFunctionBegin; 196 if (SData_d->alpha) { 197 auto alpha = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->alpha); 198 delete alpha; 199 auto beta = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->beta); 200 delete beta; 201 auto invMass = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->invMass); 202 delete invMass; 203 auto B = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->B); 204 delete B; 205 auto D = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->D); 206 delete D; 207 auto mass_w = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->mass_w); 208 delete mass_w; 209 auto invJ = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->invJ); 210 delete invJ; 211 auto x = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->x); 212 delete x; 213 auto y = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->y); 214 delete y; 215 auto z = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->z); 216 delete z; 217 auto w = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->w); 218 delete w; 219 auto Eq_m = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->Eq_m); 220 delete Eq_m; 221 if (SData_d->IPf) { 222 auto IPf = static_cast<Kokkos::View<PetscScalar*, Kokkos::LayoutLeft>*>(SData_d->IPf); 223 delete IPf; 224 } 225 } 226 PetscFunctionReturn(0); 227 } 228 229 PetscErrorCode LandauKokkosJacobian(DM plex, const PetscInt Nq, PetscReal a_Eq_m[], PetscScalar a_IPf[], const PetscInt N, const PetscScalar a_xarray[], LandauGeomData *SData_d, 230 const PetscInt num_sub_blocks, PetscReal shift, const PetscLogEvent events[], Mat JacP) 231 { 232 using scr_mem_t = Kokkos::DefaultExecutionSpace::scratch_memory_space; 233 using fieldMats_scr_t = Kokkos::View<PetscScalar**, Kokkos::LayoutRight, scr_mem_t>; 234 using idx_scr_t = Kokkos::View<PetscInt**, Kokkos::LayoutRight, scr_mem_t>; 235 using scale_scr_t = Kokkos::View<PetscReal**, Kokkos::LayoutRight, scr_mem_t>; 236 using g2_scr_t = Kokkos::View<PetscReal***, Kokkos::LayoutRight, scr_mem_t>; 237 using g3_scr_t = Kokkos::View<PetscReal****, Kokkos::LayoutRight, scr_mem_t>; 238 PetscErrorCode ierr; 239 PetscInt *Nbf,Nb,cStart,cEnd,Nf,dim,numCells,totDim,global_elem_mat_sz,nip,nfaces=0; 240 PetscDS prob; 241 LandauCtx *ctx; 242 PetscReal *d_Eq_m=NULL; 243 PetscScalar *d_IPf=NULL; 244 P4estVertexMaps *d_maps=NULL; 245 PetscSplitCSRDataStructure d_mat=NULL; 246 const int conc = Kokkos::DefaultExecutionSpace().concurrency(), openmp = !!(conc < 1000), team_size = (openmp==0) ? Nq : 1; 247 auto d_alpha_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->alpha); //static data 248 const PetscReal *d_alpha = d_alpha_k->data(); 249 auto d_beta_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->beta); 250 const PetscReal *d_beta = d_beta_k->data(); 251 auto d_invMass_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->invMass); 252 const PetscReal *d_invMass = d_invMass_k->data(); 253 auto d_B_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->B); 254 const PetscReal *d_BB = d_B_k->data(); 255 auto d_D_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->D); 256 const PetscReal *d_DD = d_D_k->data(); 257 auto d_mass_w_k = *static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->mass_w); 258 auto d_invJ_k = *static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->invJ); // use Kokkos vector in kernels 259 auto d_x_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->x); //static data 260 const PetscReal *d_x = d_x_k->data(); 261 auto d_y_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->y); //static data 262 const PetscReal *d_y = d_y_k->data(); 263 auto d_z_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->z); //static data 264 const PetscReal *d_z = d_z_k->data(); 265 auto d_w_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->w); //static data 266 const PetscReal *d_w = d_w_k->data(); 267 // dynamic data pointers 268 auto d_Eq_m_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->Eq_m); //static data 269 auto d_IPf_k = static_cast<Kokkos::View<PetscScalar*, Kokkos::LayoutLeft>*>(SData_d->IPf); //static data 270 271 PetscFunctionBegin; 272 ierr = PetscLogEventBegin(events[3],0,0,0,0);CHKERRQ(ierr); 273 ierr = DMGetApplicationContext(plex, &ctx);CHKERRQ(ierr); 274 if (!ctx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context"); 275 ierr = DMGetDimension(plex, &dim);CHKERRQ(ierr); 276 ierr = DMPlexGetHeightStratum(plex,0,&cStart,&cEnd);CHKERRQ(ierr); 277 numCells = cEnd - cStart; 278 nip = numCells*Nq; 279 ierr = DMGetDS(plex, &prob);CHKERRQ(ierr); 280 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 281 ierr = PetscDSGetDimensions(prob, &Nbf);CHKERRQ(ierr); Nb = Nbf[0]; 282 if (Nq != Nb) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Nq != Nb. %D %D",Nq,Nb); 283 if (LANDAU_DIM != dim) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dim %D != LANDAU_DIM %d",dim,LANDAU_DIM); 284 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 285 if (ctx->gpu_assembly) { 286 PetscContainer container; 287 ierr = PetscObjectQuery((PetscObject) JacP, "assembly_maps", (PetscObject *) &container);CHKERRQ(ierr); 288 if (container) { // not here first call 289 #if defined(PETSC_HAVE_KOKKOS_KERNELS) 290 P4estVertexMaps *h_maps=NULL; 291 ierr = PetscContainerGetPointer(container, (void **) &h_maps);CHKERRQ(ierr); 292 if (h_maps->data) { 293 d_maps = h_maps->data; 294 nfaces = h_maps->num_face; 295 } else { 296 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "GPU assembly but no metadata in container"); 297 } 298 // this does the setup the first time called 299 ierr = MatKokkosGetDeviceMatWrite(JacP,&d_mat);CHKERRQ(ierr); 300 global_elem_mat_sz = 0; 301 #else 302 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "GPU assembly w/o kokkos kernels -- should not be here"); 303 #endif 304 } else { // kernel output - first call assembled on device 305 global_elem_mat_sz = numCells; 306 } 307 } else { 308 global_elem_mat_sz = numCells; // no device assembly 309 } 310 ierr = PetscKokkosInitializeCheck();CHKERRQ(ierr); 311 ierr = PetscLogEventEnd(events[3],0,0,0,0);CHKERRQ(ierr); 312 313 Kokkos::View<PetscScalar**, Kokkos::LayoutRight> d_elem_mats("element matrices", global_elem_mat_sz, totDim*totDim); 314 if (a_IPf || a_xarray) { // Jacobian 315 const int num_face = nfaces; 316 const int scr_bytes = 2*(g2_scr_t::shmem_size(dim,Nf,Nq) + g3_scr_t::shmem_size(dim,dim,Nf,Nq)) + fieldMats_scr_t::shmem_size(Nb,Nb) + idx_scr_t::shmem_size(Nb,num_face) + scale_scr_t::shmem_size(Nb,num_face); 317 #if defined(LANDAU_LAYOUT_LEFT) 318 Kokkos::View<PetscReal***, Kokkos::LayoutLeft > d_fdf_k("df", dim+1, Nf, (a_IPf || a_xarray) ? nip : 0); 319 #else 320 Kokkos::View<PetscReal***, Kokkos::LayoutRight > d_fdf_k("df", dim+1, Nf, (a_IPf || a_xarray) ? nip : 0); 321 #endif 322 const Kokkos::View<PetscScalar*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> >h_IPf_k (a_IPf, (a_IPf || a_xarray) ? nip*Nf : 0); 323 const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_Eq_m_k (a_Eq_m, (a_IPf || a_xarray) ? Nf : 0); 324 Kokkos::deep_copy (*d_Eq_m_k, h_Eq_m_k); 325 d_Eq_m = d_Eq_m_k->data(); 326 // copy dynamic data to device 327 ierr = PetscLogEventBegin(events[1],0,0,0,0);CHKERRQ(ierr); 328 if (a_IPf) { 329 Kokkos::deep_copy (*d_IPf_k, h_IPf_k); 330 d_IPf = d_IPf_k->data(); 331 } else { 332 d_IPf = (PetscScalar*)a_xarray; 333 } 334 ierr = PetscLogEventEnd(events[1],0,0,0,0);CHKERRQ(ierr); 335 #define KOKKOS_SHARED_LEVEL 1 336 // get f and df 337 ierr = PetscLogEventBegin(events[8],0,0,0,0);CHKERRQ(ierr); 338 ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr); 339 Kokkos::parallel_for("f, df", Kokkos::TeamPolicy<>(numCells, team_size, /* Kokkos::AUTO */ 16), KOKKOS_LAMBDA (const team_member team) { 340 const PetscInt elem = team.league_rank(); 341 const PetscScalar *coef; 342 PetscScalar coef_buff[LANDAU_MAX_SPECIES*LANDAU_MAX_NQ]; 343 // un pack IPData 344 if (!d_maps) { 345 coef = &d_IPf[elem*Nb*Nf]; 346 } else { 347 coef = coef_buff; 348 for (int f = 0; f < Nf; ++f) { 349 LandauIdx *const Idxs = &d_maps->gIdx[elem][f][0]; 350 for (int b = 0; b < Nb; ++b) { 351 PetscInt idx = Idxs[b]; 352 if (idx >= 0) { 353 coef_buff[f*Nb+b] = d_IPf[idx]; 354 } else { 355 idx = -idx - 1; 356 coef_buff[f*Nb+b] = 0; 357 for (int q = 0; q < d_maps->num_face; q++) { 358 PetscInt id = d_maps->c_maps[idx][q].gid; 359 PetscScalar scale = d_maps->c_maps[idx][q].scale; 360 coef_buff[f*Nb+b] += scale*d_IPf[id]; 361 } 362 } 363 } 364 } 365 } 366 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,Nq), [=] (int myQi) { 367 const PetscInt ipidx = myQi + elem * Nq; 368 const PetscReal *const invJj = &d_invJ_k(ipidx*dim*dim); 369 const PetscReal *Bq = &d_BB[myQi*Nb], *Dq = &d_DD[myQi*Nb*dim]; 370 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team,0,(int)Nf), [=] (int f) { 371 PetscInt b, e, d; 372 PetscReal refSpaceDer[LANDAU_DIM]; 373 d_fdf_k(0,f,ipidx) = 0.0; 374 for (d = 0; d < LANDAU_DIM; ++d) refSpaceDer[d] = 0.0; 375 for (b = 0; b < Nb; ++b) { 376 const PetscInt cidx = b; 377 d_fdf_k(0,f,ipidx) += Bq[cidx]*PetscRealPart(coef[f*Nb+cidx]); 378 for (d = 0; d < dim; ++d) refSpaceDer[d] += Dq[cidx*dim+d]*PetscRealPart(coef[f*Nb+cidx]); 379 } 380 for (d = 0; d < dim; ++d) { 381 for (e = 0, d_fdf_k(d+1,f,ipidx) = 0.0; e < dim; ++e) { 382 d_fdf_k(d+1,f,ipidx) += invJj[e*dim+d]*refSpaceDer[e]; 383 } 384 } 385 }); // Nf 386 }); // Nq 387 }); // elems 388 #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP) 389 ierr = PetscLogGpuFlops(nip*(PetscLogDouble)(2*Nb*(1+dim)));CHKERRQ(ierr); 390 #else 391 ierr = PetscLogFlops(nip*(PetscLogDouble)(2*Nb*(1+dim)));CHKERRQ(ierr); 392 #endif 393 ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr); 394 ierr = PetscLogEventEnd(events[8],0,0,0,0);CHKERRQ(ierr); 395 // Jacobian 396 ierr = PetscLogEventBegin(events[4],0,0,0,0);CHKERRQ(ierr); 397 ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr); 398 #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP) 399 ierr = PetscLogGpuFlops(nip*(PetscLogDouble)((nip*(11*Nf+ 4*dim*dim) + 6*Nf*dim*dim*dim + 10*Nf*dim*dim + 4*Nf*dim + Nb*Nf*Nb*Nq*dim*dim*5)));CHKERRQ(ierr); 400 if (ctx->deviceType == LANDAU_CPU) PetscInfo(plex, "Warning: Landau selected CPU but no support for Kokkos using CPU\n"); 401 #else 402 ierr = PetscLogFlops(nip*(PetscLogDouble)((nip*(11*Nf+ 4*dim*dim) + 6*Nf*dim*dim*dim + 10*Nf*dim*dim + 4*Nf*dim + Nb*Nf*Nb*Nq*dim*dim*5)));CHKERRQ(ierr); 403 #endif 404 ierr = PetscInfo5(plex, "Jacobian shared memory size: %d bytes in level %d conc=%D team size=%D #face=%D\n",scr_bytes,KOKKOS_SHARED_LEVEL,conc,team_size,nfaces);CHKERRQ(ierr); 405 Kokkos::parallel_for("Landau Jacobian", Kokkos::TeamPolicy<>(numCells, team_size, /*Kokkos::AUTO*/ 16).set_scratch_size(KOKKOS_SHARED_LEVEL, Kokkos::PerTeam(scr_bytes)), KOKKOS_LAMBDA (const team_member team) { 406 const PetscInt elem = team.league_rank(); 407 g2_scr_t g2(team.team_scratch(KOKKOS_SHARED_LEVEL),dim,Nf,Nq); 408 g3_scr_t g3(team.team_scratch(KOKKOS_SHARED_LEVEL),dim,dim,Nf,Nq); 409 g2_scr_t gg2(team.team_scratch(KOKKOS_SHARED_LEVEL),dim,Nf,Nq); 410 g3_scr_t gg3(team.team_scratch(KOKKOS_SHARED_LEVEL),dim,dim,Nf,Nq); 411 // get g2[] & g3[] 412 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,Nq), [=] (int myQi) { 413 using Kokkos::parallel_reduce; 414 const PetscInt jpidx = myQi + elem * Nq; 415 const PetscReal* const invJj = &d_invJ_k(jpidx*dim*dim); 416 const PetscReal vj[3] = {d_x[jpidx], d_y[jpidx], d_z ? d_z[jpidx] : 0}, wj = d_w[jpidx]; 417 landau_inner_red::TensorValueType gg_temp; // reduce on part of gg2 and g33 for IP jpidx 418 Kokkos::parallel_reduce(Kokkos::ThreadVectorRange (team, (int)nip), [=] (const int& ipidx, landau_inner_red::TensorValueType & ggg) { 419 const PetscReal wi = d_w[ipidx], x = d_x[ipidx], y = d_y[ipidx]; 420 PetscReal temp1[3] = {0, 0, 0}, temp2 = 0; 421 PetscInt fieldA,d2,d3; 422 #if LANDAU_DIM==2 423 PetscReal Ud[2][2], Uk[2][2]; 424 LandauTensor2D(vj, x, y, Ud, Uk, (ipidx==jpidx) ? 0. : 1.); 425 #else 426 PetscReal U[3][3], z = d_z[ipidx]; 427 LandauTensor3D(vj, x, y, z, U, (ipidx==jpidx) ? 0. : 1.); 428 #endif 429 for (fieldA = 0; fieldA < Nf; ++fieldA) { 430 temp1[0] += d_fdf_k(1,fieldA,ipidx)*d_beta[fieldA]*d_invMass[fieldA]; 431 temp1[1] += d_fdf_k(2,fieldA,ipidx)*d_beta[fieldA]*d_invMass[fieldA]; 432 #if LANDAU_DIM==3 433 temp1[2] += d_fdf_k(3,fieldA,ipidx)*d_beta[fieldA]*d_invMass[fieldA]; 434 #endif 435 temp2 += d_fdf_k(0,fieldA,ipidx)*d_beta[fieldA]; 436 } 437 temp1[0] *= wi; 438 temp1[1] *= wi; 439 #if LANDAU_DIM==3 440 temp1[2] *= wi; 441 #endif 442 temp2 *= wi; 443 #if LANDAU_DIM==2 444 for (d2 = 0; d2 < 2; d2++) { 445 for (d3 = 0; d3 < 2; ++d3) { 446 /* K = U * grad(f): g2=e: i,A */ 447 ggg.gg2[d2] += Uk[d2][d3]*temp1[d3]; 448 /* D = -U * (I \kron (fx)): g3=f: i,j,A */ 449 ggg.gg3[d2][d3] += Ud[d2][d3]*temp2; 450 } 451 } 452 #else 453 for (d2 = 0; d2 < 3; ++d2) { 454 for (d3 = 0; d3 < 3; ++d3) { 455 /* K = U * grad(f): g2 = e: i,A */ 456 ggg.gg2[d2] += U[d2][d3]*temp1[d3]; 457 /* D = -U * (I \kron (fx)): g3 = f: i,j,A */ 458 ggg.gg3[d2][d3] += U[d2][d3]*temp2; 459 } 460 } 461 #endif 462 }, Kokkos::Sum<landau_inner_red::TensorValueType>(gg_temp)); 463 // add alpha and put in gg2/3 464 Kokkos::parallel_for(Kokkos::ThreadVectorRange (team, (int)Nf), [&] (const int& fieldA) { 465 PetscInt d2,d3; 466 for (d2 = 0; d2 < dim; d2++) { 467 gg2(d2,fieldA,myQi) = gg_temp.gg2[d2]*d_alpha[fieldA]; 468 for (d3 = 0; d3 < dim; d3++) { 469 gg3(d2,d3,fieldA,myQi) = -gg_temp.gg3[d2][d3]*d_alpha[fieldA]*d_invMass[fieldA]; 470 } 471 } 472 }); 473 /* add electric field term once per IP */ 474 Kokkos::parallel_for(Kokkos::ThreadVectorRange (team, (int)Nf), [&] (const int& fieldA) { 475 gg2(dim-1,fieldA,myQi) += d_Eq_m[fieldA]; 476 }); 477 Kokkos::parallel_for(Kokkos::ThreadVectorRange (team, (int)Nf), [=] (const int& fieldA) { 478 int d,d2,d3,dp; 479 /* Jacobian transform - g2, g3 - per thread (2D) */ 480 for (d = 0; d < dim; ++d) { 481 g2(d,fieldA,myQi) = 0; 482 for (d2 = 0; d2 < dim; ++d2) { 483 g2(d,fieldA,myQi) += invJj[d*dim+d2]*gg2(d2,fieldA,myQi); 484 g3(d,d2,fieldA,myQi) = 0; 485 for (d3 = 0; d3 < dim; ++d3) { 486 for (dp = 0; dp < dim; ++dp) { 487 g3(d,d2,fieldA,myQi) += invJj[d*dim + d3]*gg3(d3,dp,fieldA,myQi)*invJj[d2*dim + dp]; 488 } 489 } 490 g3(d,d2,fieldA,myQi) *= wj; 491 } 492 g2(d,fieldA,myQi) *= wj; 493 } 494 }); 495 }); // Nq 496 team.team_barrier(); 497 { /* assemble */ 498 fieldMats_scr_t s_fieldMats(team.team_scratch(KOKKOS_SHARED_LEVEL),Nq,Nq); // Only used for GPU assembly (ie, not first pass) 499 idx_scr_t s_idx(team.team_scratch(KOKKOS_SHARED_LEVEL),Nq,num_face); 500 scale_scr_t s_scale(team.team_scratch(KOKKOS_SHARED_LEVEL),Nq,num_face); 501 PetscInt fieldA; 502 for (fieldA = 0; fieldA < Nf; fieldA++) { 503 /* assemble */ 504 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,Nb), [=] (int f) { 505 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team,0,Nb), [=] (int g) { 506 PetscScalar t = 0; 507 for (int qj = 0 ; qj < Nq ; qj++) { // look at others integration points 508 const PetscReal *BJq = &d_BB[qj*Nb], *DIq = &d_DD[qj*Nb*dim]; 509 for (int d = 0; d < dim; ++d) { 510 t += DIq[f*dim+d]*g2(d,fieldA,qj)*BJq[g]; 511 for (int d2 = 0; d2 < dim; ++d2) { 512 t += DIq[f*dim + d]*g3(d,d2,fieldA,qj)*DIq[g*dim + d2]; 513 } 514 } 515 } 516 if (global_elem_mat_sz) { 517 const PetscInt fOff = (fieldA*Nb + f)*totDim + fieldA*Nb + g; 518 d_elem_mats(elem,fOff) = t; 519 } else s_fieldMats(f,g) = t; 520 }); 521 }); 522 if (!global_elem_mat_sz) { 523 const LandauIdx *const Idxs = &d_maps->gIdx[elem][fieldA][0]; 524 team.team_barrier(); 525 Kokkos::parallel_for(Kokkos::TeamVectorRange(team,0,Nb), [=] (int f) { 526 PetscInt q,idx = Idxs[f]; 527 if (idx >= 0) { 528 s_idx(f,0) = idx; 529 s_scale(f,0) = 1.; 530 } else { 531 idx = -idx - 1; 532 for (q = 0; q < nfaces; q++) { 533 s_idx(f,q) = d_maps->c_maps[idx][q].gid; 534 s_scale(f,q) = d_maps->c_maps[idx][q].scale; 535 } 536 } 537 }); 538 team.team_barrier(); 539 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,Nb), [=] (int f) { 540 PetscInt nr,idx = Idxs[f]; 541 if (idx >= 0) { 542 nr = 1; 543 } else { 544 nr = nfaces; 545 } 546 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team,0,Nb), [=] (int g) { 547 PetscScalar vals[LANDAU_MAX_Q_FACE*LANDAU_MAX_Q_FACE]; 548 PetscInt q,d,nc,idx = Idxs[g]; 549 if (idx >= 0) { 550 nc = 1; 551 } else { 552 nc = nfaces; 553 } 554 for (q = 0; q < nr; q++) { 555 for (d = 0; d < nc; d++) { 556 vals[q*nc + d] = s_scale(f,q)*s_scale(g,d)*s_fieldMats(f,g); 557 } 558 } 559 MatSetValuesDevice(d_mat,nr,&s_idx(f,0),nc,&s_idx(g,0),vals,ADD_VALUES); 560 }); 561 }); 562 } 563 } 564 } 565 }); 566 ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr); 567 ierr = PetscLogEventEnd(events[4],0,0,0,0);CHKERRQ(ierr); 568 } else { // mass 569 int scr_bytes = fieldMats_scr_t::shmem_size(Nq,Nq) + idx_scr_t::shmem_size(Nb,nfaces) + scale_scr_t::shmem_size(Nb,nfaces); 570 ierr = PetscLogEventBegin(events[4],0,0,0,0);CHKERRQ(ierr); 571 ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr); 572 #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP) 573 ierr = PetscLogGpuFlops(nip*(PetscLogDouble)(Nb*Nf*Nb*Nq*4));CHKERRQ(ierr); 574 if (ctx->deviceType == LANDAU_CPU) PetscInfo(plex, "Warning: Landau selected CPU but no support for Kokkos using CPU\n"); 575 #else 576 ierr = PetscLogFlops(nip*(PetscLogDouble)(Nb*Nf*Nb*Nq*4));CHKERRQ(ierr); 577 #endif 578 ierr = PetscInfo5(plex, "Mass shared memory size: %d bytes in level %d conc=%D team size=%D #face=%D\n",scr_bytes,KOKKOS_SHARED_LEVEL,conc,team_size,nfaces);CHKERRQ(ierr); 579 Kokkos::parallel_for("Landau mass", Kokkos::TeamPolicy<>(numCells, team_size, /*Kokkos::AUTO*/ 16).set_scratch_size(KOKKOS_SHARED_LEVEL, Kokkos::PerTeam(scr_bytes)), KOKKOS_LAMBDA (const team_member team) { 580 const PetscInt elem = team.league_rank(); 581 fieldMats_scr_t s_fieldMats(team.team_scratch(KOKKOS_SHARED_LEVEL),Nb,Nb); 582 idx_scr_t s_idx(team.team_scratch(KOKKOS_SHARED_LEVEL),Nb,nfaces); 583 scale_scr_t s_scale(team.team_scratch(KOKKOS_SHARED_LEVEL),Nb,nfaces); 584 for (int fieldA = 0; fieldA < Nf; fieldA++) { 585 /* assemble */ 586 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,Nb), [=] (int f) { 587 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team,0,Nb), [=] (int g) { 588 PetscScalar t = 0; 589 for (int qj = 0 ; qj < Nq ; qj++) { // look at others integration points 590 const PetscReal *BJq = &d_BB[qj*Nb]; 591 const PetscInt jpidx = qj + elem * Nq; 592 t += BJq[f] * d_mass_w_k(jpidx) * shift * BJq[g]; 593 } 594 if (global_elem_mat_sz) { 595 const PetscInt fOff = (fieldA*Nb + f)*totDim + fieldA*Nb + g; 596 d_elem_mats(elem,fOff) = t; 597 } else s_fieldMats(f,g) = t; 598 }); 599 }); 600 if (!global_elem_mat_sz) { 601 const LandauIdx *const Idxs = &d_maps->gIdx[elem][fieldA][0]; 602 team.team_barrier(); // needed ???? 603 Kokkos::parallel_for(Kokkos::TeamVectorRange(team,0,Nb), [=] (int f) { 604 PetscInt q,idx = Idxs[f]; 605 if (idx >= 0) { 606 s_idx(f,0) = idx; 607 s_scale(f,0) = 1.; 608 } else { 609 idx = -idx - 1; 610 for (q = 0; q < nfaces; q++) { 611 s_idx(f,q) = d_maps->c_maps[idx][q].gid; 612 s_scale(f,q) = d_maps->c_maps[idx][q].scale; 613 } 614 } 615 }); 616 team.team_barrier(); 617 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,Nb), [=] (int f) { 618 PetscInt nr,idx = Idxs[f]; 619 if (idx >= 0) { 620 nr = 1; 621 } else { 622 nr = nfaces; 623 } 624 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team,0,Nb), [=] (int g) { 625 PetscScalar vals[LANDAU_MAX_Q_FACE*LANDAU_MAX_Q_FACE]; 626 PetscInt q,d,nc,idx = Idxs[g]; 627 if (idx >= 0) { 628 nc = 1; 629 } else { 630 nc = nfaces; 631 } 632 for (q = 0; q < nr; q++) { 633 for (d = 0; d < nc; d++) { 634 vals[q*nc + d] = s_scale(f,q)*s_scale(g,d)*s_fieldMats(f,g); 635 } 636 } 637 MatSetValuesDevice(d_mat,nr,&s_idx(f,0),nc,&s_idx(g,0),vals,ADD_VALUES); 638 }); 639 }); 640 } 641 } 642 }); 643 ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr); 644 ierr = PetscLogEventEnd(events[4],0,0,0,0);CHKERRQ(ierr); 645 } 646 Kokkos::fence(); 647 648 if (global_elem_mat_sz) { 649 PetscSection section, globalSection; 650 Kokkos::View<PetscScalar**, Kokkos::LayoutRight>::HostMirror h_elem_mats = Kokkos::create_mirror_view(d_elem_mats); 651 ierr = PetscLogEventBegin(events[5],0,0,0,0);CHKERRQ(ierr); 652 ierr = DMGetLocalSection(plex, §ion);CHKERRQ(ierr); 653 ierr = DMGetGlobalSection(plex, &globalSection);CHKERRQ(ierr); 654 Kokkos::deep_copy (h_elem_mats, d_elem_mats); 655 ierr = PetscLogEventEnd(events[5],0,0,0,0);CHKERRQ(ierr); 656 ierr = PetscLogEventBegin(events[6],0,0,0,0);CHKERRQ(ierr); 657 PetscInt ej; 658 for (ej = cStart ; ej < cEnd; ++ej) { 659 const PetscScalar *elMat = &h_elem_mats(ej-cStart,0); 660 ierr = DMPlexMatSetClosure(plex, section, globalSection, JacP, ej, elMat, ADD_VALUES);CHKERRQ(ierr); 661 if (ej==-1) { 662 int d,f; 663 PetscPrintf(PETSC_COMM_SELF,"Kokkos Element matrix %d/%d\n",1,(int)numCells); 664 for (d = 0; d < totDim; ++d){ 665 for (f = 0; f < totDim; ++f) PetscPrintf(PETSC_COMM_SELF," %12.5e", PetscRealPart(elMat[d*totDim + f])); 666 PetscPrintf(PETSC_COMM_SELF,"\n"); 667 } 668 exit(14); 669 } 670 } 671 ierr = PetscLogEventEnd(events[6],0,0,0,0);CHKERRQ(ierr); 672 // transition to use of maps for VecGetClosure 673 if (ctx->gpu_assembly) { 674 auto IPf = static_cast<Kokkos::View<PetscScalar*, Kokkos::LayoutLeft>*>(SData_d->IPf); 675 delete IPf; 676 SData_d->IPf = NULL; 677 if (!(a_IPf || a_xarray)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "transition without Jacobian"); 678 } 679 } 680 PetscFunctionReturn(0); 681 } 682 } // extern "C" 683