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