1 #include <../src/mat/impls/aij/seq/aij.h> 2 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3 #include <petsclandau.h> /*I "petsclandau.h" I*/ 4 #include <petscts.h> 5 #include <petscdmforest.h> 6 #include <petscdmcomposite.h> 7 8 /* Landau collision operator */ 9 10 /* relativistic terms */ 11 #if defined(PETSC_USE_REAL_SINGLE) 12 #define SPEED_OF_LIGHT 2.99792458e8F 13 #define C_0(v0) (SPEED_OF_LIGHT / v0) /* needed for relativistic tensor on all architectures */ 14 #else 15 #define SPEED_OF_LIGHT 2.99792458e8 16 #define C_0(v0) (SPEED_OF_LIGHT / v0) /* needed for relativistic tensor on all architectures */ 17 #endif 18 19 #include "land_tensors.h" 20 21 #if defined(PETSC_HAVE_OPENMP) 22 #include <omp.h> 23 #endif 24 25 static PetscErrorCode LandauGPUMapsDestroy(void *ptr) 26 { 27 P4estVertexMaps *maps = (P4estVertexMaps *)ptr; 28 PetscFunctionBegin; 29 // free device data 30 if (maps[0].deviceType != LANDAU_CPU) { 31 #if defined(PETSC_HAVE_KOKKOS_KERNELS) 32 if (maps[0].deviceType == LANDAU_KOKKOS) { 33 PetscCall(LandauKokkosDestroyMatMaps(maps, maps[0].numgrids)); // implies Kokkos does 34 } // else could be CUDA 35 #elif defined(PETSC_HAVE_CUDA) 36 if (maps[0].deviceType == LANDAU_CUDA) { 37 PetscCall(LandauCUDADestroyMatMaps(maps, maps[0].numgrids)); 38 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "maps->deviceType %d ?????", maps->deviceType); 39 #endif 40 } 41 // free host data 42 for (PetscInt grid = 0; grid < maps[0].numgrids; grid++) { 43 PetscCall(PetscFree(maps[grid].c_maps)); 44 PetscCall(PetscFree(maps[grid].gIdx)); 45 } 46 PetscCall(PetscFree(maps)); 47 48 PetscFunctionReturn(PETSC_SUCCESS); 49 } 50 static PetscErrorCode energy_f(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx) 51 { 52 PetscReal v2 = 0; 53 PetscFunctionBegin; 54 /* compute v^2 / 2 */ 55 for (int i = 0; i < dim; ++i) v2 += x[i] * x[i]; 56 /* evaluate the Maxwellian */ 57 u[0] = v2 / 2; 58 PetscFunctionReturn(PETSC_SUCCESS); 59 } 60 61 /* needs double */ 62 static PetscErrorCode gamma_m1_f(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx) 63 { 64 PetscReal *c2_0_arr = ((PetscReal *)actx); 65 double u2 = 0, c02 = (double)*c2_0_arr, xx; 66 67 PetscFunctionBegin; 68 /* compute u^2 / 2 */ 69 for (int i = 0; i < dim; ++i) u2 += x[i] * x[i]; 70 /* gamma - 1 = g_eps, for conditioning and we only take derivatives */ 71 xx = u2 / c02; 72 #if defined(PETSC_USE_DEBUG) 73 u[0] = PetscSqrtReal(1. + xx); 74 #else 75 u[0] = xx / (PetscSqrtReal(1. + xx) + 1.) - 1.; // better conditioned. -1 might help condition and only used for derivative 76 #endif 77 PetscFunctionReturn(PETSC_SUCCESS); 78 } 79 80 /* 81 LandauFormJacobian_Internal - Evaluates Jacobian matrix. 82 83 Input Parameters: 84 . globX - input vector 85 . actx - optional user-defined context 86 . dim - dimension 87 88 Output Parameter: 89 . J0acP - Jacobian matrix filled, not created 90 */ 91 static PetscErrorCode LandauFormJacobian_Internal(Vec a_X, Mat JacP, const PetscInt dim, PetscReal shift, void *a_ctx) 92 { 93 LandauCtx *ctx = (LandauCtx *)a_ctx; 94 PetscInt numCells[LANDAU_MAX_GRIDS], Nq, Nb; 95 PetscQuadrature quad; 96 PetscReal Eq_m[LANDAU_MAX_SPECIES]; // could be static data w/o quench (ex2) 97 PetscScalar *cellClosure = NULL; 98 const PetscScalar *xdata = NULL; 99 PetscDS prob; 100 PetscContainer container; 101 P4estVertexMaps *maps; 102 Mat subJ[LANDAU_MAX_GRIDS * LANDAU_MAX_BATCH_SZ]; 103 104 PetscFunctionBegin; 105 PetscValidHeaderSpecific(a_X, VEC_CLASSID, 1); 106 PetscValidHeaderSpecific(JacP, MAT_CLASSID, 2); 107 PetscValidPointer(ctx, 5); 108 /* check for matrix container for GPU assembly. Support CPU assembly for debugging */ 109 PetscCheck(ctx->plex[0] != NULL, ctx->comm, PETSC_ERR_ARG_WRONG, "Plex not created"); 110 PetscCall(PetscLogEventBegin(ctx->events[10], 0, 0, 0, 0)); 111 PetscCall(DMGetDS(ctx->plex[0], &prob)); // same DS for all grids 112 PetscCall(PetscObjectQuery((PetscObject)JacP, "assembly_maps", (PetscObject *)&container)); 113 if (container) { 114 PetscCheck(ctx->gpu_assembly, ctx->comm, PETSC_ERR_ARG_WRONG, "maps but no GPU assembly"); 115 PetscCall(PetscContainerGetPointer(container, (void **)&maps)); 116 PetscCheck(maps, ctx->comm, PETSC_ERR_ARG_WRONG, "empty GPU matrix container"); 117 for (PetscInt i = 0; i < ctx->num_grids * ctx->batch_sz; i++) subJ[i] = NULL; 118 } else { 119 PetscCheck(!ctx->gpu_assembly, ctx->comm, PETSC_ERR_ARG_WRONG, "No maps but GPU assembly"); 120 for (PetscInt tid = 0; tid < ctx->batch_sz; tid++) { 121 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(DMCreateMatrix(ctx->plex[grid], &subJ[LAND_PACK_IDX(tid, grid)])); 122 } 123 maps = NULL; 124 } 125 // get dynamic data (Eq is odd, for quench and Spitzer test) for CPU assembly and raw data for Jacobian GPU assembly. Get host numCells[], Nq (yuck) 126 PetscCall(PetscFEGetQuadrature(ctx->fe[0], &quad)); 127 PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL)); 128 Nb = Nq; 129 PetscCheck(Nq <= LANDAU_MAX_NQ, ctx->comm, PETSC_ERR_ARG_WRONG, "Order too high. Nq = %" PetscInt_FMT " > LANDAU_MAX_NQ (%d)", Nq, LANDAU_MAX_NQ); 130 // get metadata for collecting dynamic data 131 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 132 PetscInt cStart, cEnd; 133 PetscCheck(ctx->plex[grid] != NULL, ctx->comm, PETSC_ERR_ARG_WRONG, "Plex not created"); 134 PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, &cEnd)); 135 numCells[grid] = cEnd - cStart; // grids can have different topology 136 } 137 PetscCall(PetscLogEventEnd(ctx->events[10], 0, 0, 0, 0)); 138 if (shift == 0) { /* create dynamic point data: f_alpha for closure of each cell (cellClosure[nbatch,ngrids,ncells[g],f[Nb,ns[g]]]) or xdata */ 139 DM pack; 140 PetscCall(VecGetDM(a_X, &pack)); 141 PetscCheck(pack, PETSC_COMM_SELF, PETSC_ERR_PLIB, "pack has no DM"); 142 PetscCall(PetscLogEventBegin(ctx->events[1], 0, 0, 0, 0)); 143 for (PetscInt fieldA = 0; fieldA < ctx->num_species; fieldA++) { 144 Eq_m[fieldA] = ctx->Ez * ctx->t_0 * ctx->charges[fieldA] / (ctx->v_0 * ctx->masses[fieldA]); /* normalize dimensionless */ 145 if (dim == 2) Eq_m[fieldA] *= 2 * PETSC_PI; /* add the 2pi term that is not in Landau */ 146 } 147 if (!ctx->gpu_assembly) { 148 Vec *locXArray, *globXArray; 149 PetscScalar *cellClosure_it; 150 PetscInt cellClosure_sz = 0, nDMs, Nf[LANDAU_MAX_GRIDS]; 151 PetscSection section[LANDAU_MAX_GRIDS], globsection[LANDAU_MAX_GRIDS]; 152 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 153 PetscCall(DMGetLocalSection(ctx->plex[grid], §ion[grid])); 154 PetscCall(DMGetGlobalSection(ctx->plex[grid], &globsection[grid])); 155 PetscCall(PetscSectionGetNumFields(section[grid], &Nf[grid])); 156 } 157 /* count cellClosure size */ 158 PetscCall(DMCompositeGetNumberDM(pack, &nDMs)); 159 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) cellClosure_sz += Nb * Nf[grid] * numCells[grid]; 160 PetscCall(PetscMalloc1(cellClosure_sz * ctx->batch_sz, &cellClosure)); 161 cellClosure_it = cellClosure; 162 PetscCall(PetscMalloc(sizeof(*locXArray) * nDMs, &locXArray)); 163 PetscCall(PetscMalloc(sizeof(*globXArray) * nDMs, &globXArray)); 164 PetscCall(DMCompositeGetLocalAccessArray(pack, a_X, nDMs, NULL, locXArray)); 165 PetscCall(DMCompositeGetAccessArray(pack, a_X, nDMs, NULL, globXArray)); 166 for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { // OpenMP (once) 167 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 168 Vec locX = locXArray[LAND_PACK_IDX(b_id, grid)], globX = globXArray[LAND_PACK_IDX(b_id, grid)], locX2; 169 PetscInt cStart, cEnd, ei; 170 PetscCall(VecDuplicate(locX, &locX2)); 171 PetscCall(DMGlobalToLocalBegin(ctx->plex[grid], globX, INSERT_VALUES, locX2)); 172 PetscCall(DMGlobalToLocalEnd(ctx->plex[grid], globX, INSERT_VALUES, locX2)); 173 PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, &cEnd)); 174 for (ei = cStart; ei < cEnd; ++ei) { 175 PetscScalar *coef = NULL; 176 PetscCall(DMPlexVecGetClosure(ctx->plex[grid], section[grid], locX2, ei, NULL, &coef)); 177 PetscCall(PetscMemcpy(cellClosure_it, coef, Nb * Nf[grid] * sizeof(*cellClosure_it))); /* change if LandauIPReal != PetscScalar */ 178 PetscCall(DMPlexVecRestoreClosure(ctx->plex[grid], section[grid], locX2, ei, NULL, &coef)); 179 cellClosure_it += Nb * Nf[grid]; 180 } 181 PetscCall(VecDestroy(&locX2)); 182 } 183 } 184 PetscCheck(cellClosure_it - cellClosure == cellClosure_sz * ctx->batch_sz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "iteration wrong %" PetscCount_FMT " != cellClosure_sz = %" PetscInt_FMT, (PetscCount)(cellClosure_it - cellClosure), 185 cellClosure_sz * ctx->batch_sz); 186 PetscCall(DMCompositeRestoreLocalAccessArray(pack, a_X, nDMs, NULL, locXArray)); 187 PetscCall(DMCompositeRestoreAccessArray(pack, a_X, nDMs, NULL, globXArray)); 188 PetscCall(PetscFree(locXArray)); 189 PetscCall(PetscFree(globXArray)); 190 xdata = NULL; 191 } else { 192 PetscMemType mtype; 193 if (ctx->jacobian_field_major_order) { // get data in batch ordering 194 PetscCall(VecScatterBegin(ctx->plex_batch, a_X, ctx->work_vec, INSERT_VALUES, SCATTER_FORWARD)); 195 PetscCall(VecScatterEnd(ctx->plex_batch, a_X, ctx->work_vec, INSERT_VALUES, SCATTER_FORWARD)); 196 PetscCall(VecGetArrayReadAndMemType(ctx->work_vec, &xdata, &mtype)); 197 } else { 198 PetscCall(VecGetArrayReadAndMemType(a_X, &xdata, &mtype)); 199 } 200 PetscCheck(mtype == PETSC_MEMTYPE_HOST || ctx->deviceType != LANDAU_CPU, ctx->comm, PETSC_ERR_ARG_WRONG, "CPU run with device data: use -mat_type aij"); 201 cellClosure = NULL; 202 } 203 PetscCall(PetscLogEventEnd(ctx->events[1], 0, 0, 0, 0)); 204 } else xdata = cellClosure = NULL; 205 206 /* do it */ 207 if (ctx->deviceType == LANDAU_CUDA || ctx->deviceType == LANDAU_KOKKOS) { 208 if (ctx->deviceType == LANDAU_CUDA) { 209 #if defined(PETSC_HAVE_CUDA) 210 PetscCall(LandauCUDAJacobian(ctx->plex, Nq, ctx->batch_sz, ctx->num_grids, numCells, Eq_m, cellClosure, xdata, &ctx->SData_d, shift, ctx->events, ctx->mat_offset, ctx->species_offset, subJ, JacP)); 211 #else 212 SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-landau_device_type %s not built", "cuda"); 213 #endif 214 } else if (ctx->deviceType == LANDAU_KOKKOS) { 215 #if defined(PETSC_HAVE_KOKKOS_KERNELS) 216 PetscCall(LandauKokkosJacobian(ctx->plex, Nq, ctx->batch_sz, ctx->num_grids, numCells, Eq_m, cellClosure, xdata, &ctx->SData_d, shift, ctx->events, ctx->mat_offset, ctx->species_offset, subJ, JacP)); 217 #else 218 SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-landau_device_type %s not built", "kokkos"); 219 #endif 220 } 221 } else { /* CPU version */ 222 PetscTabulation *Tf; // used for CPU and print info. Same on all grids and all species 223 PetscInt ip_offset[LANDAU_MAX_GRIDS + 1], ipf_offset[LANDAU_MAX_GRIDS + 1], elem_offset[LANDAU_MAX_GRIDS + 1], IPf_sz_glb, IPf_sz_tot, num_grids = ctx->num_grids, Nf[LANDAU_MAX_GRIDS]; 224 PetscReal *ff, *dudx, *dudy, *dudz, *invJ_a = (PetscReal *)ctx->SData_d.invJ, *xx = (PetscReal *)ctx->SData_d.x, *yy = (PetscReal *)ctx->SData_d.y, *zz = (PetscReal *)ctx->SData_d.z, *ww = (PetscReal *)ctx->SData_d.w; 225 PetscReal *nu_alpha = (PetscReal *)ctx->SData_d.alpha, *nu_beta = (PetscReal *)ctx->SData_d.beta, *invMass = (PetscReal *)ctx->SData_d.invMass; 226 PetscReal(*lambdas)[LANDAU_MAX_GRIDS][LANDAU_MAX_GRIDS] = (PetscReal(*)[LANDAU_MAX_GRIDS][LANDAU_MAX_GRIDS])ctx->SData_d.lambdas; 227 PetscSection section[LANDAU_MAX_GRIDS], globsection[LANDAU_MAX_GRIDS]; 228 PetscScalar *coo_vals = NULL; 229 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 230 PetscCall(DMGetLocalSection(ctx->plex[grid], §ion[grid])); 231 PetscCall(DMGetGlobalSection(ctx->plex[grid], &globsection[grid])); 232 PetscCall(PetscSectionGetNumFields(section[grid], &Nf[grid])); 233 } 234 /* count IPf size, etc */ 235 PetscCall(PetscDSGetTabulation(prob, &Tf)); // Bf, &Df same for all grids 236 const PetscReal *const BB = Tf[0]->T[0], *const DD = Tf[0]->T[1]; 237 ip_offset[0] = ipf_offset[0] = elem_offset[0] = 0; 238 for (PetscInt grid = 0; grid < num_grids; grid++) { 239 PetscInt nfloc = ctx->species_offset[grid + 1] - ctx->species_offset[grid]; 240 elem_offset[grid + 1] = elem_offset[grid] + numCells[grid]; 241 ip_offset[grid + 1] = ip_offset[grid] + numCells[grid] * Nq; 242 ipf_offset[grid + 1] = ipf_offset[grid] + Nq * nfloc * numCells[grid]; 243 } 244 IPf_sz_glb = ipf_offset[num_grids]; 245 IPf_sz_tot = IPf_sz_glb * ctx->batch_sz; 246 // prep COO 247 if (ctx->coo_assembly) { 248 PetscCall(PetscMalloc1(ctx->SData_d.coo_size, &coo_vals)); // allocate every time? 249 PetscCall(PetscInfo(ctx->plex[0], "COO Allocate %" PetscInt_FMT " values\n", (PetscInt)ctx->SData_d.coo_size)); 250 } 251 if (shift == 0.0) { /* compute dynamic data f and df and init data for Jacobian */ 252 #if defined(PETSC_HAVE_THREADSAFETY) 253 double starttime, endtime; 254 starttime = MPI_Wtime(); 255 #endif 256 PetscCall(PetscLogEventBegin(ctx->events[8], 0, 0, 0, 0)); 257 PetscCall(PetscMalloc4(IPf_sz_tot, &ff, IPf_sz_tot, &dudx, IPf_sz_tot, &dudy, dim == 3 ? IPf_sz_tot : 0, &dudz)); 258 // F df/dx 259 for (PetscInt tid = 0; tid < ctx->batch_sz * elem_offset[num_grids]; tid++) { // for each element 260 const PetscInt b_Nelem = elem_offset[num_grids], b_elem_idx = tid % b_Nelem, b_id = tid / b_Nelem; // b_id == OMP thd_id in batch 261 // find my grid: 262 PetscInt grid = 0; 263 while (b_elem_idx >= elem_offset[grid + 1]) grid++; // yuck search for grid 264 { 265 const PetscInt loc_nip = numCells[grid] * Nq, loc_Nf = ctx->species_offset[grid + 1] - ctx->species_offset[grid], loc_elem = b_elem_idx - elem_offset[grid]; 266 const PetscInt moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset); //b_id*b_N + ctx->mat_offset[grid]; 267 PetscScalar *coef, coef_buff[LANDAU_MAX_SPECIES * LANDAU_MAX_NQ]; 268 PetscReal *invJe = &invJ_a[(ip_offset[grid] + loc_elem * Nq) * dim * dim]; // ingJ is static data on batch 0 269 PetscInt b, f, q; 270 if (cellClosure) { 271 coef = &cellClosure[b_id * IPf_sz_glb + ipf_offset[grid] + loc_elem * Nb * loc_Nf]; // this is const 272 } else { 273 coef = coef_buff; 274 for (f = 0; f < loc_Nf; ++f) { 275 LandauIdx *const Idxs = &maps[grid].gIdx[loc_elem][f][0]; 276 for (b = 0; b < Nb; ++b) { 277 PetscInt idx = Idxs[b]; 278 if (idx >= 0) { 279 coef[f * Nb + b] = xdata[idx + moffset]; 280 } else { 281 idx = -idx - 1; 282 coef[f * Nb + b] = 0; 283 for (q = 0; q < maps[grid].num_face; q++) { 284 PetscInt id = maps[grid].c_maps[idx][q].gid; 285 PetscScalar scale = maps[grid].c_maps[idx][q].scale; 286 coef[f * Nb + b] += scale * xdata[id + moffset]; 287 } 288 } 289 } 290 } 291 } 292 /* get f and df */ 293 for (PetscInt qi = 0; qi < Nq; qi++) { 294 const PetscReal *invJ = &invJe[qi * dim * dim]; 295 const PetscReal *Bq = &BB[qi * Nb]; 296 const PetscReal *Dq = &DD[qi * Nb * dim]; 297 PetscReal u_x[LANDAU_DIM]; 298 /* get f & df */ 299 for (f = 0; f < loc_Nf; ++f) { 300 const PetscInt idx = b_id * IPf_sz_glb + ipf_offset[grid] + f * loc_nip + loc_elem * Nq + qi; 301 PetscInt b, e; 302 PetscReal refSpaceDer[LANDAU_DIM]; 303 ff[idx] = 0.0; 304 for (int d = 0; d < LANDAU_DIM; ++d) refSpaceDer[d] = 0.0; 305 for (b = 0; b < Nb; ++b) { 306 const PetscInt cidx = b; 307 ff[idx] += Bq[cidx] * PetscRealPart(coef[f * Nb + cidx]); 308 for (int d = 0; d < dim; ++d) refSpaceDer[d] += Dq[cidx * dim + d] * PetscRealPart(coef[f * Nb + cidx]); 309 } 310 for (int d = 0; d < LANDAU_DIM; ++d) { 311 for (e = 0, u_x[d] = 0.0; e < LANDAU_DIM; ++e) u_x[d] += invJ[e * dim + d] * refSpaceDer[e]; 312 } 313 dudx[idx] = u_x[0]; 314 dudy[idx] = u_x[1]; 315 #if LANDAU_DIM == 3 316 dudz[idx] = u_x[2]; 317 #endif 318 } 319 } // q 320 } // grid 321 } // grid*batch 322 PetscCall(PetscLogEventEnd(ctx->events[8], 0, 0, 0, 0)); 323 #if defined(PETSC_HAVE_THREADSAFETY) 324 endtime = MPI_Wtime(); 325 if (ctx->stage) ctx->times[LANDAU_F_DF] += (endtime - starttime); 326 #endif 327 } // Jacobian setup 328 // assemble Jacobian (or mass) 329 for (PetscInt tid = 0; tid < ctx->batch_sz * elem_offset[num_grids]; tid++) { // for each element 330 const PetscInt b_Nelem = elem_offset[num_grids]; 331 const PetscInt glb_elem_idx = tid % b_Nelem, b_id = tid / b_Nelem; 332 PetscInt grid = 0; 333 #if defined(PETSC_HAVE_THREADSAFETY) 334 double starttime, endtime; 335 starttime = MPI_Wtime(); 336 #endif 337 while (glb_elem_idx >= elem_offset[grid + 1]) grid++; 338 { 339 const PetscInt loc_Nf = ctx->species_offset[grid + 1] - ctx->species_offset[grid], loc_elem = glb_elem_idx - elem_offset[grid]; 340 const PetscInt moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset), totDim = loc_Nf * Nq, elemMatSize = totDim * totDim; 341 PetscScalar *elemMat; 342 const PetscReal *invJe = &invJ_a[(ip_offset[grid] + loc_elem * Nq) * dim * dim]; 343 PetscCall(PetscMalloc1(elemMatSize, &elemMat)); 344 PetscCall(PetscMemzero(elemMat, elemMatSize * sizeof(*elemMat))); 345 if (shift == 0.0) { // Jacobian 346 PetscCall(PetscLogEventBegin(ctx->events[4], 0, 0, 0, 0)); 347 } else { // mass 348 PetscCall(PetscLogEventBegin(ctx->events[16], 0, 0, 0, 0)); 349 } 350 for (PetscInt qj = 0; qj < Nq; ++qj) { 351 const PetscInt jpidx_glb = ip_offset[grid] + qj + loc_elem * Nq; 352 PetscReal g0[LANDAU_MAX_SPECIES], g2[LANDAU_MAX_SPECIES][LANDAU_DIM], g3[LANDAU_MAX_SPECIES][LANDAU_DIM][LANDAU_DIM]; // could make a LANDAU_MAX_SPECIES_GRID ~ number of ions - 1 353 PetscInt d, d2, dp, d3, IPf_idx; 354 if (shift == 0.0) { // Jacobian 355 const PetscReal *const invJj = &invJe[qj * dim * dim]; 356 PetscReal gg2[LANDAU_MAX_SPECIES][LANDAU_DIM], gg3[LANDAU_MAX_SPECIES][LANDAU_DIM][LANDAU_DIM], gg2_temp[LANDAU_DIM], gg3_temp[LANDAU_DIM][LANDAU_DIM]; 357 const PetscReal vj[3] = {xx[jpidx_glb], yy[jpidx_glb], zz ? zz[jpidx_glb] : 0}, wj = ww[jpidx_glb]; 358 // create g2 & g3 359 for (d = 0; d < LANDAU_DIM; d++) { // clear accumulation data D & K 360 gg2_temp[d] = 0; 361 for (d2 = 0; d2 < LANDAU_DIM; d2++) gg3_temp[d][d2] = 0; 362 } 363 /* inner beta reduction */ 364 IPf_idx = 0; 365 for (PetscInt grid_r = 0, f_off = 0, ipidx = 0; grid_r < ctx->num_grids; grid_r++, f_off = ctx->species_offset[grid_r]) { // IPf_idx += nip_loc_r*Nfloc_r 366 PetscInt nip_loc_r = numCells[grid_r] * Nq, Nfloc_r = Nf[grid_r]; 367 for (PetscInt ei_r = 0, loc_fdf_idx = 0; ei_r < numCells[grid_r]; ++ei_r) { 368 for (PetscInt qi = 0; qi < Nq; qi++, ipidx++, loc_fdf_idx++) { 369 const PetscReal wi = ww[ipidx], x = xx[ipidx], y = yy[ipidx]; 370 PetscReal temp1[3] = {0, 0, 0}, temp2 = 0; 371 #if LANDAU_DIM == 2 372 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.; 373 LandauTensor2D(vj, x, y, Ud, Uk, mask); 374 #else 375 PetscReal U[3][3], z = zz[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.; 376 if (ctx->use_relativistic_corrections) { 377 LandauTensor3DRelativistic(vj, x, y, z, U, mask, C_0(ctx->v_0)); 378 } else { 379 LandauTensor3D(vj, x, y, z, U, mask); 380 } 381 #endif 382 for (int f = 0; f < Nfloc_r; ++f) { 383 const PetscInt idx = b_id * IPf_sz_glb + ipf_offset[grid_r] + f * nip_loc_r + ei_r * Nq + qi; // IPf_idx + f*nip_loc_r + loc_fdf_idx; 384 temp1[0] += dudx[idx] * nu_beta[f + f_off] * invMass[f + f_off] * (*lambdas)[grid][grid_r]; 385 temp1[1] += dudy[idx] * nu_beta[f + f_off] * invMass[f + f_off] * (*lambdas)[grid][grid_r]; 386 #if LANDAU_DIM == 3 387 temp1[2] += dudz[idx] * nu_beta[f + f_off] * invMass[f + f_off] * (*lambdas)[grid][grid_r]; 388 #endif 389 temp2 += ff[idx] * nu_beta[f + f_off] * (*lambdas)[grid][grid_r]; 390 } 391 temp1[0] *= wi; 392 temp1[1] *= wi; 393 #if LANDAU_DIM == 3 394 temp1[2] *= wi; 395 #endif 396 temp2 *= wi; 397 #if LANDAU_DIM == 2 398 for (d2 = 0; d2 < 2; d2++) { 399 for (d3 = 0; d3 < 2; ++d3) { 400 /* K = U * grad(f): g2=e: i,A */ 401 gg2_temp[d2] += Uk[d2][d3] * temp1[d3]; 402 /* D = -U * (I \kron (fx)): g3=f: i,j,A */ 403 gg3_temp[d2][d3] += Ud[d2][d3] * temp2; 404 } 405 } 406 #else 407 for (d2 = 0; d2 < 3; ++d2) { 408 for (d3 = 0; d3 < 3; ++d3) { 409 /* K = U * grad(f): g2 = e: i,A */ 410 gg2_temp[d2] += U[d2][d3] * temp1[d3]; 411 /* D = -U * (I \kron (fx)): g3 = f: i,j,A */ 412 gg3_temp[d2][d3] += U[d2][d3] * temp2; 413 } 414 } 415 #endif 416 } // qi 417 } // ei_r 418 IPf_idx += nip_loc_r * Nfloc_r; 419 } /* grid_r - IPs */ 420 PetscCheck(IPf_idx == IPf_sz_glb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "IPf_idx != IPf_sz %" PetscInt_FMT " %" PetscInt_FMT, IPf_idx, IPf_sz_glb); 421 // add alpha and put in gg2/3 422 for (PetscInt fieldA = 0, f_off = ctx->species_offset[grid]; fieldA < loc_Nf; ++fieldA) { 423 for (d2 = 0; d2 < LANDAU_DIM; d2++) { 424 gg2[fieldA][d2] = gg2_temp[d2] * nu_alpha[fieldA + f_off]; 425 for (d3 = 0; d3 < LANDAU_DIM; d3++) gg3[fieldA][d2][d3] = -gg3_temp[d2][d3] * nu_alpha[fieldA + f_off] * invMass[fieldA + f_off]; 426 } 427 } 428 /* add electric field term once per IP */ 429 for (PetscInt fieldA = 0, f_off = ctx->species_offset[grid]; fieldA < loc_Nf; ++fieldA) gg2[fieldA][LANDAU_DIM - 1] += Eq_m[fieldA + f_off]; 430 /* Jacobian transform - g2, g3 */ 431 for (PetscInt fieldA = 0; fieldA < loc_Nf; ++fieldA) { 432 for (d = 0; d < dim; ++d) { 433 g2[fieldA][d] = 0.0; 434 for (d2 = 0; d2 < dim; ++d2) { 435 g2[fieldA][d] += invJj[d * dim + d2] * gg2[fieldA][d2]; 436 g3[fieldA][d][d2] = 0.0; 437 for (d3 = 0; d3 < dim; ++d3) { 438 for (dp = 0; dp < dim; ++dp) g3[fieldA][d][d2] += invJj[d * dim + d3] * gg3[fieldA][d3][dp] * invJj[d2 * dim + dp]; 439 } 440 g3[fieldA][d][d2] *= wj; 441 } 442 g2[fieldA][d] *= wj; 443 } 444 } 445 } else { // mass 446 PetscReal wj = ww[jpidx_glb]; 447 /* Jacobian transform - g0 */ 448 for (PetscInt fieldA = 0; fieldA < loc_Nf; ++fieldA) { 449 if (dim == 2) { 450 g0[fieldA] = wj * shift * 2. * PETSC_PI; // move this to below and remove g0 451 } else { 452 g0[fieldA] = wj * shift; // move this to below and remove g0 453 } 454 } 455 } 456 /* FE matrix construction */ 457 { 458 PetscInt fieldA, d, f, d2, g; 459 const PetscReal *BJq = &BB[qj * Nb], *DIq = &DD[qj * Nb * dim]; 460 /* assemble - on the diagonal (I,I) */ 461 for (fieldA = 0; fieldA < loc_Nf; fieldA++) { 462 for (f = 0; f < Nb; f++) { 463 const PetscInt i = fieldA * Nb + f; /* Element matrix row */ 464 for (g = 0; g < Nb; ++g) { 465 const PetscInt j = fieldA * Nb + g; /* Element matrix column */ 466 const PetscInt fOff = i * totDim + j; 467 if (shift == 0.0) { 468 for (d = 0; d < dim; ++d) { 469 elemMat[fOff] += DIq[f * dim + d] * g2[fieldA][d] * BJq[g]; 470 for (d2 = 0; d2 < dim; ++d2) elemMat[fOff] += DIq[f * dim + d] * g3[fieldA][d][d2] * DIq[g * dim + d2]; 471 } 472 } else { // mass 473 elemMat[fOff] += BJq[f] * g0[fieldA] * BJq[g]; 474 } 475 } 476 } 477 } 478 } 479 } /* qj loop */ 480 if (shift == 0.0) { // Jacobian 481 PetscCall(PetscLogEventEnd(ctx->events[4], 0, 0, 0, 0)); 482 } else { 483 PetscCall(PetscLogEventEnd(ctx->events[16], 0, 0, 0, 0)); 484 } 485 #if defined(PETSC_HAVE_THREADSAFETY) 486 endtime = MPI_Wtime(); 487 if (ctx->stage) ctx->times[LANDAU_KERNEL] += (endtime - starttime); 488 #endif 489 /* assemble matrix */ 490 if (!container) { 491 PetscInt cStart; 492 PetscCall(PetscLogEventBegin(ctx->events[6], 0, 0, 0, 0)); 493 PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, NULL)); 494 PetscCall(DMPlexMatSetClosure(ctx->plex[grid], section[grid], globsection[grid], subJ[LAND_PACK_IDX(b_id, grid)], loc_elem + cStart, elemMat, ADD_VALUES)); 495 PetscCall(PetscLogEventEnd(ctx->events[6], 0, 0, 0, 0)); 496 } else { // GPU like assembly for debugging 497 PetscInt fieldA, q, f, g, d, nr, nc, rows0[LANDAU_MAX_Q_FACE] = {0}, cols0[LANDAU_MAX_Q_FACE] = {0}, rows[LANDAU_MAX_Q_FACE], cols[LANDAU_MAX_Q_FACE]; 498 PetscScalar vals[LANDAU_MAX_Q_FACE * LANDAU_MAX_Q_FACE] = {0}, row_scale[LANDAU_MAX_Q_FACE] = {0}, col_scale[LANDAU_MAX_Q_FACE] = {0}; 499 LandauIdx *coo_elem_offsets = (LandauIdx *)ctx->SData_d.coo_elem_offsets, *coo_elem_fullNb = (LandauIdx *)ctx->SData_d.coo_elem_fullNb, (*coo_elem_point_offsets)[LANDAU_MAX_NQ + 1] = (LandauIdx(*)[LANDAU_MAX_NQ + 1]) ctx->SData_d.coo_elem_point_offsets; 500 /* assemble - from the diagonal (I,I) in this format for DMPlexMatSetClosure */ 501 for (fieldA = 0; fieldA < loc_Nf; fieldA++) { 502 LandauIdx *const Idxs = &maps[grid].gIdx[loc_elem][fieldA][0]; 503 for (f = 0; f < Nb; f++) { 504 PetscInt idx = Idxs[f]; 505 if (idx >= 0) { 506 nr = 1; 507 rows0[0] = idx; 508 row_scale[0] = 1.; 509 } else { 510 idx = -idx - 1; 511 for (q = 0, nr = 0; q < maps[grid].num_face; q++, nr++) { 512 if (maps[grid].c_maps[idx][q].gid < 0) break; 513 rows0[q] = maps[grid].c_maps[idx][q].gid; 514 row_scale[q] = maps[grid].c_maps[idx][q].scale; 515 } 516 } 517 for (g = 0; g < Nb; ++g) { 518 idx = Idxs[g]; 519 if (idx >= 0) { 520 nc = 1; 521 cols0[0] = idx; 522 col_scale[0] = 1.; 523 } else { 524 idx = -idx - 1; 525 nc = maps[grid].num_face; 526 for (q = 0, nc = 0; q < maps[grid].num_face; q++, nc++) { 527 if (maps[grid].c_maps[idx][q].gid < 0) break; 528 cols0[q] = maps[grid].c_maps[idx][q].gid; 529 col_scale[q] = maps[grid].c_maps[idx][q].scale; 530 } 531 } 532 const PetscInt i = fieldA * Nb + f; /* Element matrix row */ 533 const PetscInt j = fieldA * Nb + g; /* Element matrix column */ 534 const PetscScalar Aij = elemMat[i * totDim + j]; 535 if (coo_vals) { // mirror (i,j) in CreateStaticGPUData 536 const int fullNb = coo_elem_fullNb[glb_elem_idx], fullNb2 = fullNb * fullNb; 537 const int idx0 = b_id * coo_elem_offsets[elem_offset[num_grids]] + 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]; 538 for (int q = 0, idx2 = idx0; q < nr; q++) { 539 for (int d = 0; d < nc; d++, idx2++) coo_vals[idx2] = row_scale[q] * col_scale[d] * Aij; 540 } 541 } else { 542 for (q = 0; q < nr; q++) rows[q] = rows0[q] + moffset; 543 for (d = 0; d < nc; d++) cols[d] = cols0[d] + moffset; 544 for (q = 0; q < nr; q++) { 545 for (d = 0; d < nc; d++) vals[q * nc + d] = row_scale[q] * col_scale[d] * Aij; 546 } 547 PetscCall(MatSetValues(JacP, nr, rows, nc, cols, vals, ADD_VALUES)); 548 } 549 } 550 } 551 } 552 } 553 if (loc_elem == -1) { 554 PetscCall(PetscPrintf(ctx->comm, "CPU Element matrix\n")); 555 for (int d = 0; d < totDim; ++d) { 556 for (int f = 0; f < totDim; ++f) PetscCall(PetscPrintf(ctx->comm, " %12.5e", (double)PetscRealPart(elemMat[d * totDim + f]))); 557 PetscCall(PetscPrintf(ctx->comm, "\n")); 558 } 559 exit(12); 560 } 561 PetscCall(PetscFree(elemMat)); 562 } /* grid */ 563 } /* outer element & batch loop */ 564 if (shift == 0.0) { // mass 565 PetscCall(PetscFree4(ff, dudx, dudy, dudz)); 566 } 567 if (!container) { // 'CPU' assembly move nest matrix to global JacP 568 for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { // OpenMP 569 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 570 const PetscInt moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset); // b_id*b_N + ctx->mat_offset[grid]; 571 PetscInt nloc, nzl, colbuf[1024], row; 572 const PetscInt *cols; 573 const PetscScalar *vals; 574 Mat B = subJ[LAND_PACK_IDX(b_id, grid)]; 575 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 576 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 577 PetscCall(MatGetSize(B, &nloc, NULL)); 578 for (int i = 0; i < nloc; i++) { 579 PetscCall(MatGetRow(B, i, &nzl, &cols, &vals)); 580 PetscCheck(nzl <= 1024, PetscObjectComm((PetscObject)B), PETSC_ERR_PLIB, "Row too big: %" PetscInt_FMT, nzl); 581 for (int j = 0; j < nzl; j++) colbuf[j] = moffset + cols[j]; 582 row = moffset + i; 583 PetscCall(MatSetValues(JacP, 1, &row, nzl, colbuf, vals, ADD_VALUES)); 584 PetscCall(MatRestoreRow(B, i, &nzl, &cols, &vals)); 585 } 586 PetscCall(MatDestroy(&B)); 587 } 588 } 589 } 590 if (coo_vals) { 591 PetscCall(MatSetValuesCOO(JacP, coo_vals, ADD_VALUES)); 592 PetscCall(PetscFree(coo_vals)); 593 } 594 } /* CPU version */ 595 PetscCall(MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY)); 596 PetscCall(MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY)); 597 /* clean up */ 598 if (cellClosure) PetscCall(PetscFree(cellClosure)); 599 if (xdata) PetscCall(VecRestoreArrayReadAndMemType(a_X, &xdata)); 600 PetscFunctionReturn(PETSC_SUCCESS); 601 } 602 603 static PetscErrorCode GeometryDMLandau(DM base, PetscInt point, PetscInt dim, const PetscReal abc[], PetscReal xyz[], void *a_ctx) 604 { 605 PetscReal r = abc[0], z = abc[1]; 606 607 PetscFunctionBegin; 608 xyz[0] = r; 609 xyz[1] = z; 610 if (dim == 3) xyz[2] = abc[2]; 611 612 PetscFunctionReturn(PETSC_SUCCESS); 613 } 614 615 /* create DMComposite of meshes for each species group */ 616 static PetscErrorCode LandauDMCreateVMeshes(MPI_Comm comm_self, const PetscInt dim, const char prefix[], LandauCtx *ctx, DM pack) 617 { 618 PetscFunctionBegin; 619 { /* p4est, quads */ 620 /* Create plex mesh of Landau domain */ 621 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 622 PetscReal par_radius = ctx->radius_par[grid], perp_radius = ctx->radius_perp[grid]; 623 if (!ctx->sphere) { 624 PetscReal lo[] = {-perp_radius, -par_radius, -par_radius}, hi[] = {perp_radius, par_radius, par_radius}; 625 DMBoundaryType periodicity[3] = {DM_BOUNDARY_NONE, dim == 2 ? DM_BOUNDARY_NONE : DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 626 if (dim == 2) lo[0] = 0; 627 else { 628 lo[1] = -perp_radius; 629 hi[1] = perp_radius; // 3D y is a perp 630 } 631 PetscCall(DMPlexCreateBoxMesh(comm_self, dim, PETSC_FALSE, ctx->cells0, lo, hi, periodicity, PETSC_TRUE, &ctx->plex[grid])); // todo: make composite and create dm[grid] here 632 PetscCall(DMLocalizeCoordinates(ctx->plex[grid])); /* needed for periodic */ 633 if (dim == 3) PetscCall(PetscObjectSetName((PetscObject)ctx->plex[grid], "cube")); 634 else PetscCall(PetscObjectSetName((PetscObject)ctx->plex[grid], "half-plane")); 635 } else if (dim == 2) { 636 PetscInt numCells = 6, cells[11][4], i, j; 637 PetscInt numVerts = 11; 638 PetscReal *flatCoords = NULL; 639 PetscInt *flatCells = NULL, *pcell; 640 int cells2[][4] = { 641 {0, 1, 6, 5 }, 642 {1, 2, 7, 6 }, 643 {2, 3, 8, 7 }, 644 {3, 4, 9, 8 }, 645 {5, 6, 7, 10}, 646 {10, 7, 8, 9 }, 647 }; 648 for (i = 0; i < numCells; i++) 649 for (j = 0; j < 4; j++) cells[i][j] = cells2[i][j]; 650 PetscCall(PetscMalloc2(numVerts * 2, &flatCoords, numCells * 4, &flatCells)); 651 { 652 PetscReal(*coords)[2] = (PetscReal(*)[2])flatCoords; 653 PetscReal rad = ctx->radius[grid]; 654 for (j = 0; j < 5; j++) { 655 PetscReal z, r, theta = -PETSC_PI / 2 + (j % 5) * PETSC_PI / 4; 656 r = rad * PetscCosReal(theta); 657 coords[j][0] = r; 658 z = rad * PetscSinReal(theta); 659 coords[j][1] = z; 660 } 661 coords[j][0] = 0; 662 coords[j++][1] = -rad * ctx->sphere_inner_radius_90degree; 663 coords[j][0] = rad * ctx->sphere_inner_radius_45degree; 664 coords[j++][1] = -rad * ctx->sphere_inner_radius_45degree; 665 coords[j][0] = rad * ctx->sphere_inner_radius_90degree; 666 coords[j++][1] = 0; 667 coords[j][0] = rad * ctx->sphere_inner_radius_45degree; 668 coords[j++][1] = rad * ctx->sphere_inner_radius_45degree; 669 coords[j][0] = 0; 670 coords[j++][1] = rad * ctx->sphere_inner_radius_90degree; 671 coords[j][0] = 0; 672 coords[j++][1] = 0; 673 } 674 for (j = 0, pcell = flatCells; j < numCells; j++, pcell += 4) { 675 for (int jj = 0; jj < 4; jj++) pcell[jj] = cells[j][jj]; 676 } 677 PetscCall(DMPlexCreateFromCellListPetsc(comm_self, 2, numCells, numVerts, 4, ctx->interpolate, flatCells, 2, flatCoords, &ctx->plex[grid])); 678 PetscCall(PetscFree2(flatCoords, flatCells)); 679 PetscCall(PetscObjectSetName((PetscObject)ctx->plex[grid], "semi-circle")); 680 } else SETERRQ(ctx->comm, PETSC_ERR_PLIB, "Velocity space meshes does not support 3V cubed sphere"); 681 PetscCall(DMSetFromOptions(ctx->plex[grid])); 682 } // grid loop 683 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pack, prefix)); 684 685 { /* convert to p4est (or whatever), wait for discretization to create pack */ 686 char convType[256]; 687 PetscBool flg; 688 689 PetscOptionsBegin(ctx->comm, prefix, "Mesh conversion options", "DMPLEX"); 690 PetscCall(PetscOptionsFList("-dm_landau_type", "Convert DMPlex to another format (p4est)", "plexland.c", DMList, DMPLEX, convType, 256, &flg)); 691 PetscOptionsEnd(); 692 if (flg) { 693 ctx->use_p4est = PETSC_TRUE; /* flag for Forest */ 694 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 695 DM dmforest; 696 PetscCall(DMConvert(ctx->plex[grid], convType, &dmforest)); 697 if (dmforest) { 698 PetscBool isForest; 699 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dmforest, prefix)); 700 PetscCall(DMIsForest(dmforest, &isForest)); 701 if (isForest) { 702 if (ctx->sphere) PetscCall(DMForestSetBaseCoordinateMapping(dmforest, GeometryDMLandau, ctx)); 703 PetscCall(DMDestroy(&ctx->plex[grid])); 704 ctx->plex[grid] = dmforest; // Forest for adaptivity 705 } else SETERRQ(ctx->comm, PETSC_ERR_PLIB, "Converted to non Forest?"); 706 } else SETERRQ(ctx->comm, PETSC_ERR_PLIB, "Convert failed?"); 707 } 708 } else ctx->use_p4est = PETSC_FALSE; /* flag for Forest */ 709 } 710 } /* non-file */ 711 PetscCall(DMSetDimension(pack, dim)); 712 PetscCall(PetscObjectSetName((PetscObject)pack, "Mesh")); 713 PetscCall(DMSetApplicationContext(pack, ctx)); 714 715 PetscFunctionReturn(PETSC_SUCCESS); 716 } 717 718 static PetscErrorCode SetupDS(DM pack, PetscInt dim, PetscInt grid, LandauCtx *ctx) 719 { 720 PetscInt ii, i0; 721 char buf[256]; 722 PetscSection section; 723 724 PetscFunctionBegin; 725 for (ii = ctx->species_offset[grid], i0 = 0; ii < ctx->species_offset[grid + 1]; ii++, i0++) { 726 if (ii == 0) PetscCall(PetscSNPrintf(buf, sizeof(buf), "e")); 727 else PetscCall(PetscSNPrintf(buf, sizeof(buf), "i%" PetscInt_FMT, ii)); 728 /* Setup Discretization - FEM */ 729 PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, NULL, PETSC_DECIDE, &ctx->fe[ii])); 730 PetscCall(PetscObjectSetName((PetscObject)ctx->fe[ii], buf)); 731 PetscCall(DMSetField(ctx->plex[grid], i0, NULL, (PetscObject)ctx->fe[ii])); 732 } 733 PetscCall(DMCreateDS(ctx->plex[grid])); 734 PetscCall(DMGetSection(ctx->plex[grid], §ion)); 735 for (PetscInt ii = ctx->species_offset[grid], i0 = 0; ii < ctx->species_offset[grid + 1]; ii++, i0++) { 736 if (ii == 0) PetscCall(PetscSNPrintf(buf, sizeof(buf), "se")); 737 else PetscCall(PetscSNPrintf(buf, sizeof(buf), "si%" PetscInt_FMT, ii)); 738 PetscCall(PetscSectionSetComponentName(section, i0, 0, buf)); 739 } 740 PetscFunctionReturn(PETSC_SUCCESS); 741 } 742 743 /* Define a Maxwellian function for testing out the operator. */ 744 745 /* Using cartesian velocity space coordinates, the particle */ 746 /* density, [1/m^3], is defined according to */ 747 748 /* $$ n=\int_{R^3} dv^3 \left(\frac{m}{2\pi T}\right)^{3/2}\exp [- mv^2/(2T)] $$ */ 749 750 /* Using some constant, c, we normalize the velocity vector into a */ 751 /* dimensionless variable according to v=c*x. Thus the density, $n$, becomes */ 752 753 /* $$ n=\int_{R^3} dx^3 \left(\frac{mc^2}{2\pi T}\right)^{3/2}\exp [- mc^2/(2T)*x^2] $$ */ 754 755 /* Defining $\theta=2T/mc^2$, we thus find that the probability density */ 756 /* for finding the particle within the interval in a box dx^3 around x is */ 757 758 /* f(x;\theta)=\left(\frac{1}{\pi\theta}\right)^{3/2} \exp [ -x^2/\theta ] */ 759 760 typedef struct { 761 PetscReal v_0; 762 PetscReal kT_m; 763 PetscReal n; 764 PetscReal shift; 765 } MaxwellianCtx; 766 767 static PetscErrorCode maxwellian(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx) 768 { 769 MaxwellianCtx *mctx = (MaxwellianCtx *)actx; 770 PetscInt i; 771 PetscReal v2 = 0, theta = 2 * mctx->kT_m / (mctx->v_0 * mctx->v_0), shift; /* theta = 2kT/mc^2 */ 772 PetscFunctionBegin; 773 /* compute the exponents, v^2 */ 774 for (i = 0; i < dim; ++i) v2 += x[i] * x[i]; 775 /* evaluate the Maxwellian */ 776 if (mctx->shift < 0) shift = -mctx->shift; 777 else { 778 u[0] = mctx->n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta)); 779 shift = mctx->shift; 780 } 781 if (shift != 0.) { 782 v2 = 0; 783 for (i = 0; i < dim - 1; ++i) v2 += x[i] * x[i]; 784 v2 += (x[dim - 1] - shift) * (x[dim - 1] - shift); 785 /* evaluate the shifted Maxwellian */ 786 u[0] += mctx->n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta)); 787 } 788 PetscFunctionReturn(PETSC_SUCCESS); 789 } 790 791 /*@ 792 DMPlexLandauAddMaxwellians - Add a Maxwellian distribution to a state 793 794 Collective 795 796 Input Parameters: 797 . dm - The mesh (local) 798 + time - Current time 799 - temps - Temperatures of each species (global) 800 . ns - Number density of each species (global) 801 - grid - index into current grid - just used for offset into temp and ns 802 . b_id - batch index 803 - n_batch - number of batches 804 + actx - Landau context 805 806 Output Parameter: 807 . X - The state (local to this grid) 808 809 Level: beginner 810 811 .keywords: mesh 812 .seealso: `DMPlexLandauCreateVelocitySpace()` 813 @*/ 814 PetscErrorCode DMPlexLandauAddMaxwellians(DM dm, Vec X, PetscReal time, PetscReal temps[], PetscReal ns[], PetscInt grid, PetscInt b_id, PetscInt n_batch, void *actx) 815 { 816 LandauCtx *ctx = (LandauCtx *)actx; 817 PetscErrorCode (*initu[LANDAU_MAX_SPECIES])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 818 PetscInt dim; 819 MaxwellianCtx *mctxs[LANDAU_MAX_SPECIES], data[LANDAU_MAX_SPECIES]; 820 821 PetscFunctionBegin; 822 PetscCall(DMGetDimension(dm, &dim)); 823 if (!ctx) PetscCall(DMGetApplicationContext(dm, &ctx)); 824 for (PetscInt ii = ctx->species_offset[grid], i0 = 0; ii < ctx->species_offset[grid + 1]; ii++, i0++) { 825 mctxs[i0] = &data[i0]; 826 data[i0].v_0 = ctx->v_0; // v_0 same for all grids 827 data[i0].kT_m = ctx->k * temps[ii] / ctx->masses[ii]; /* kT/m */ 828 data[i0].n = ns[ii]; 829 initu[i0] = maxwellian; 830 data[i0].shift = 0; 831 } 832 data[0].shift = ctx->electronShift; 833 /* need to make ADD_ALL_VALUES work - TODO */ 834 PetscCall(DMProjectFunction(dm, time, initu, (void **)mctxs, INSERT_ALL_VALUES, X)); 835 PetscFunctionReturn(PETSC_SUCCESS); 836 } 837 838 /* 839 LandauSetInitialCondition - Adds Maxwellians with context 840 841 Collective 842 843 Input Parameters: 844 . dm - The mesh 845 - grid - index into current grid - just used for offset into temp and ns 846 . b_id - batch index 847 - n_batch - number of batches 848 + actx - Landau context with T and n 849 850 Output Parameter: 851 . X - The state 852 853 Level: beginner 854 855 .keywords: mesh 856 .seealso: `DMPlexLandauCreateVelocitySpace()`, `DMPlexLandauAddMaxwellians()` 857 */ 858 static PetscErrorCode LandauSetInitialCondition(DM dm, Vec X, PetscInt grid, PetscInt b_id, PetscInt n_batch, void *actx) 859 { 860 LandauCtx *ctx = (LandauCtx *)actx; 861 PetscFunctionBegin; 862 if (!ctx) PetscCall(DMGetApplicationContext(dm, &ctx)); 863 PetscCall(VecZeroEntries(X)); 864 PetscCall(DMPlexLandauAddMaxwellians(dm, X, 0.0, ctx->thermal_temps, ctx->n, grid, b_id, n_batch, ctx)); 865 PetscFunctionReturn(PETSC_SUCCESS); 866 } 867 868 // adapt a level once. Forest in/out 869 #if defined(PETSC_USE_INFO) 870 static const char *s_refine_names[] = {"RE", "Z1", "Origin", "Z2", "Uniform"}; 871 #endif 872 static PetscErrorCode adaptToleranceFEM(PetscFE fem, Vec sol, PetscInt type, PetscInt grid, LandauCtx *ctx, DM *newForest) 873 { 874 DM forest, plex, adaptedDM = NULL; 875 PetscDS prob; 876 PetscBool isForest; 877 PetscQuadrature quad; 878 PetscInt Nq, *Nb, cStart, cEnd, c, dim, qj, k; 879 DMLabel adaptLabel = NULL; 880 881 PetscFunctionBegin; 882 forest = ctx->plex[grid]; 883 PetscCall(DMCreateDS(forest)); 884 PetscCall(DMGetDS(forest, &prob)); 885 PetscCall(DMGetDimension(forest, &dim)); 886 PetscCall(DMIsForest(forest, &isForest)); 887 PetscCheck(isForest, ctx->comm, PETSC_ERR_ARG_WRONG, "! Forest"); 888 PetscCall(DMConvert(forest, DMPLEX, &plex)); 889 PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd)); 890 PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 891 PetscCall(PetscFEGetQuadrature(fem, &quad)); 892 PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL)); 893 PetscCheck(Nq <= LANDAU_MAX_NQ, ctx->comm, PETSC_ERR_ARG_WRONG, "Order too high. Nq = %" PetscInt_FMT " > LANDAU_MAX_NQ (%d)", Nq, LANDAU_MAX_NQ); 894 PetscCall(PetscDSGetDimensions(prob, &Nb)); 895 PetscCall(PetscInfo(sol, "%" PetscInt_FMT ") Refine phase: %s\n", grid, s_refine_names[type])); 896 if (type == 4) { 897 for (c = cStart; c < cEnd; c++) PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE)); 898 } else if (type == 2) { 899 PetscInt rCellIdx[8], nr = 0, nrmax = (dim == 3) ? 8 : 2; 900 PetscReal minRad = PETSC_INFINITY, r; 901 for (c = cStart; c < cEnd; c++) { 902 PetscReal tt, v0[LANDAU_MAX_NQ * 3], detJ[LANDAU_MAX_NQ]; 903 PetscCall(DMPlexComputeCellGeometryFEM(plex, c, quad, v0, NULL, NULL, detJ)); 904 for (qj = 0; qj < Nq; ++qj) { 905 tt = PetscSqr(v0[dim * qj + 0]) + PetscSqr(v0[dim * qj + 1]) + PetscSqr(((dim == 3) ? v0[dim * qj + 2] : 0)); 906 r = PetscSqrtReal(tt); 907 if (r < minRad - PETSC_SQRT_MACHINE_EPSILON * 10.) { 908 minRad = r; 909 nr = 0; 910 rCellIdx[nr++] = c; 911 PetscCall(PetscInfo(sol, "\t\t%" PetscInt_FMT ") Found first inner r=%e, cell %" PetscInt_FMT ", qp %" PetscInt_FMT "/%" PetscInt_FMT "\n", grid, (double)r, c, qj + 1, Nq)); 912 } else if ((r - minRad) < PETSC_SQRT_MACHINE_EPSILON * 100. && nr < nrmax) { 913 for (k = 0; k < nr; k++) 914 if (c == rCellIdx[k]) break; 915 if (k == nr) { 916 rCellIdx[nr++] = c; 917 PetscCall(PetscInfo(sol, "\t\t\t%" PetscInt_FMT ") Found another inner r=%e, cell %" PetscInt_FMT ", qp %" PetscInt_FMT "/%" PetscInt_FMT ", d=%e\n", grid, (double)r, c, qj + 1, Nq, (double)(r - minRad))); 918 } 919 } 920 } 921 } 922 for (k = 0; k < nr; k++) PetscCall(DMLabelSetValue(adaptLabel, rCellIdx[k], DM_ADAPT_REFINE)); 923 PetscCall(PetscInfo(sol, "\t\t\t%" PetscInt_FMT ") Refined %" PetscInt_FMT " origin cells %" PetscInt_FMT ",%" PetscInt_FMT " r=%g\n", grid, nr, rCellIdx[0], rCellIdx[1], (double)minRad)); 924 } else if (type == 0 || type == 1 || type == 3) { /* refine along r=0 axis */ 925 PetscScalar *coef = NULL; 926 Vec coords; 927 PetscInt csize, Nv, d, nz, nrefined = 0; 928 DM cdm; 929 PetscSection cs; 930 PetscCall(DMGetCoordinatesLocal(forest, &coords)); 931 PetscCall(DMGetCoordinateDM(forest, &cdm)); 932 PetscCall(DMGetLocalSection(cdm, &cs)); 933 for (c = cStart; c < cEnd; c++) { 934 PetscInt doit = 0, outside = 0; 935 PetscCall(DMPlexVecGetClosure(cdm, cs, coords, c, &csize, &coef)); 936 Nv = csize / dim; 937 for (nz = d = 0; d < Nv; d++) { 938 PetscReal z = PetscRealPart(coef[d * dim + (dim - 1)]), x = PetscSqr(PetscRealPart(coef[d * dim + 0])) + ((dim == 3) ? PetscSqr(PetscRealPart(coef[d * dim + 1])) : 0); 939 x = PetscSqrtReal(x); 940 if (type == 0) { 941 if (ctx->re_radius > PETSC_SQRT_MACHINE_EPSILON && (z < -PETSC_MACHINE_EPSILON * 10. || z > ctx->re_radius + PETSC_MACHINE_EPSILON * 10.)) outside++; /* first pass don't refine bottom */ 942 } else if (type == 1 && (z > ctx->vperp0_radius1 || z < -ctx->vperp0_radius1)) { 943 outside++; /* don't refine outside electron refine radius */ 944 PetscCall(PetscInfo(sol, "\t%" PetscInt_FMT ") (debug) found %s cells\n", grid, s_refine_names[type])); 945 } else if (type == 3 && (z > ctx->vperp0_radius2 || z < -ctx->vperp0_radius2)) { 946 outside++; /* refine r=0 cells on refinement front */ 947 PetscCall(PetscInfo(sol, "\t%" PetscInt_FMT ") (debug) found %s cells\n", grid, s_refine_names[type])); 948 } 949 if (x < PETSC_MACHINE_EPSILON * 10. && (type != 0 || ctx->re_radius > PETSC_SQRT_MACHINE_EPSILON)) nz++; 950 } 951 PetscCall(DMPlexVecRestoreClosure(cdm, cs, coords, c, &csize, &coef)); 952 if (doit || (outside < Nv && nz)) { 953 PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE)); 954 nrefined++; 955 } 956 } 957 PetscCall(PetscInfo(sol, "\t%" PetscInt_FMT ") Refined %" PetscInt_FMT " cells\n", grid, nrefined)); 958 } 959 PetscCall(DMDestroy(&plex)); 960 PetscCall(DMAdaptLabel(forest, adaptLabel, &adaptedDM)); 961 PetscCall(DMLabelDestroy(&adaptLabel)); 962 *newForest = adaptedDM; 963 if (adaptedDM) { 964 if (isForest) { 965 PetscCall(DMForestSetAdaptivityForest(adaptedDM, NULL)); // ???? 966 } 967 PetscCall(DMConvert(adaptedDM, DMPLEX, &plex)); 968 PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd)); 969 PetscCall(PetscInfo(sol, "\t\t\t\t%" PetscInt_FMT ") %" PetscInt_FMT " cells, %" PetscInt_FMT " total quadrature points\n", grid, cEnd - cStart, Nq * (cEnd - cStart))); 970 PetscCall(DMDestroy(&plex)); 971 } else *newForest = NULL; 972 PetscFunctionReturn(PETSC_SUCCESS); 973 } 974 975 // forest goes in (ctx->plex[grid]), plex comes out 976 static PetscErrorCode adapt(PetscInt grid, LandauCtx *ctx, Vec *uu) 977 { 978 PetscInt adaptIter; 979 980 PetscFunctionBegin; 981 PetscInt type, limits[5] = {(grid == 0) ? ctx->numRERefine : 0, (grid == 0) ? ctx->nZRefine1 : 0, ctx->numAMRRefine[grid], (grid == 0) ? ctx->nZRefine2 : 0, ctx->postAMRRefine[grid]}; 982 for (type = 0; type < 5; type++) { 983 for (adaptIter = 0; adaptIter < limits[type]; adaptIter++) { 984 DM newForest = NULL; 985 PetscCall(adaptToleranceFEM(ctx->fe[0], *uu, type, grid, ctx, &newForest)); 986 if (newForest) { 987 PetscCall(DMDestroy(&ctx->plex[grid])); 988 PetscCall(VecDestroy(uu)); 989 PetscCall(DMCreateGlobalVector(newForest, uu)); 990 PetscCall(PetscObjectSetName((PetscObject)*uu, "uAMR")); 991 PetscCall(LandauSetInitialCondition(newForest, *uu, grid, 0, 1, ctx)); 992 ctx->plex[grid] = newForest; 993 } else { 994 PetscCall(PetscInfo(*uu, "No refinement\n")); 995 } 996 } 997 } 998 PetscFunctionReturn(PETSC_SUCCESS); 999 } 1000 1001 // make log(Lambdas) from NRL Plasma formulary 1002 static PetscErrorCode makeLambdas(LandauCtx *ctx) 1003 { 1004 PetscFunctionBegin; 1005 for (PetscInt gridi = 0; gridi < ctx->num_grids; gridi++) { 1006 int iii = ctx->species_offset[gridi]; 1007 PetscReal Ti_ev = (ctx->thermal_temps[iii] / 1.1604525e7) * 1000; // convert (back) to eV 1008 PetscReal ni = ctx->n[iii] * ctx->n_0; 1009 for (PetscInt gridj = gridi; gridj < ctx->num_grids; gridj++) { 1010 PetscInt jjj = ctx->species_offset[gridj]; 1011 PetscReal Zj = ctx->charges[jjj] / 1.6022e-19; 1012 if (gridi == 0) { 1013 if (gridj == 0) { // lam_ee 1014 ctx->lambdas[gridi][gridj] = 23.5 - PetscLogReal(PetscSqrtReal(ni) * PetscPowReal(Ti_ev, -1.25)) - PetscSqrtReal(1e-5 + PetscSqr(PetscLogReal(Ti_ev) - 2) / 16); 1015 } else { // lam_ei == lam_ie 1016 if (10 * Zj * Zj > Ti_ev) { 1017 ctx->lambdas[gridi][gridj] = ctx->lambdas[gridj][gridi] = 23 - PetscLogReal(PetscSqrtReal(ni) * Zj * PetscPowReal(Ti_ev, -1.5)); 1018 } else { 1019 ctx->lambdas[gridi][gridj] = ctx->lambdas[gridj][gridi] = 24 - PetscLogReal(PetscSqrtReal(ni) / Ti_ev); 1020 } 1021 } 1022 } else { // lam_ii' 1023 PetscReal mui = ctx->masses[iii] / 1.6720e-27, Zi = ctx->charges[iii] / 1.6022e-19; 1024 PetscReal Tj_ev = (ctx->thermal_temps[jjj] / 1.1604525e7) * 1000; // convert (back) to eV 1025 PetscReal muj = ctx->masses[jjj] / 1.6720e-27; 1026 PetscReal nj = ctx->n[jjj] * ctx->n_0; 1027 ctx->lambdas[gridi][gridj] = ctx->lambdas[gridj][gridi] = 23 - PetscLogReal(Zi * Zj * (mui + muj) / (mui * Tj_ev + muj * Ti_ev) * PetscSqrtReal(ni * Zi * Zi / Ti_ev + nj * Zj * Zj / Tj_ev)); 1028 } 1029 } 1030 } 1031 //PetscReal v0 = PetscSqrtReal(ctx->k * ctx->thermal_temps[iii] / ctx->masses[iii]); /* arbitrary units for non-dimensionalization: plasma formulary def */ 1032 PetscFunctionReturn(PETSC_SUCCESS); 1033 } 1034 1035 static PetscErrorCode ProcessOptions(LandauCtx *ctx, const char prefix[]) 1036 { 1037 PetscBool flg; 1038 PetscInt ii, nt, nm, nc, num_species_grid[LANDAU_MAX_GRIDS], non_dim_grid; 1039 PetscReal v0_grid[LANDAU_MAX_GRIDS], lnLam = 10; 1040 DM dummy; 1041 1042 PetscFunctionBegin; 1043 PetscCall(DMCreate(ctx->comm, &dummy)); 1044 /* get options - initialize context */ 1045 ctx->verbose = 1; // should be 0 for silent compliance 1046 ctx->batch_sz = 1; 1047 ctx->batch_view_idx = 0; 1048 ctx->interpolate = PETSC_TRUE; 1049 ctx->gpu_assembly = PETSC_TRUE; 1050 ctx->norm_state = 0; 1051 ctx->electronShift = 0; 1052 ctx->M = NULL; 1053 ctx->J = NULL; 1054 /* geometry and grids */ 1055 ctx->sphere = PETSC_FALSE; 1056 ctx->use_p4est = PETSC_FALSE; 1057 for (PetscInt grid = 0; grid < LANDAU_MAX_GRIDS; grid++) { 1058 ctx->radius[grid] = 5.; /* thermal radius (velocity) */ 1059 ctx->radius_perp[grid] = 5.; /* thermal radius (velocity) */ 1060 ctx->radius_par[grid] = 5.; /* thermal radius (velocity) */ 1061 ctx->numAMRRefine[grid] = 0; 1062 ctx->postAMRRefine[grid] = 0; 1063 ctx->species_offset[grid + 1] = 1; // one species default 1064 num_species_grid[grid] = 0; 1065 ctx->plex[grid] = NULL; /* cache as expensive to Convert */ 1066 } 1067 ctx->species_offset[0] = 0; 1068 ctx->re_radius = 0.; 1069 ctx->vperp0_radius1 = 0; 1070 ctx->vperp0_radius2 = 0; 1071 ctx->nZRefine1 = 0; 1072 ctx->nZRefine2 = 0; 1073 ctx->numRERefine = 0; 1074 num_species_grid[0] = 1; // one species default 1075 /* species - [0] electrons, [1] one ion species eg, duetarium, [2] heavy impurity ion, ... */ 1076 ctx->charges[0] = -1; /* electron charge (MKS) */ 1077 ctx->masses[0] = 1 / 1835.469965278441013; /* temporary value in proton mass */ 1078 ctx->n[0] = 1; 1079 ctx->v_0 = 1; /* thermal velocity, we could start with a scale != 1 */ 1080 ctx->thermal_temps[0] = 1; 1081 /* constants, etc. */ 1082 ctx->epsilon0 = 8.8542e-12; /* permittivity of free space (MKS) F/m */ 1083 ctx->k = 1.38064852e-23; /* Boltzmann constant (MKS) J/K */ 1084 ctx->n_0 = 1.e20; /* typical plasma n, but could set it to 1 */ 1085 ctx->Ez = 0; 1086 for (PetscInt grid = 0; grid < LANDAU_NUM_TIMERS; grid++) ctx->times[grid] = 0; 1087 for (PetscInt ii = 0; ii < LANDAU_DIM; ii++) ctx->cells0[ii] = 2; 1088 if (LANDAU_DIM == 2) ctx->cells0[0] = 1; 1089 ctx->use_matrix_mass = PETSC_FALSE; 1090 ctx->use_relativistic_corrections = PETSC_FALSE; 1091 ctx->use_energy_tensor_trick = PETSC_FALSE; /* Use Eero's trick for energy conservation v --> grad(v^2/2) */ 1092 ctx->SData_d.w = NULL; 1093 ctx->SData_d.x = NULL; 1094 ctx->SData_d.y = NULL; 1095 ctx->SData_d.z = NULL; 1096 ctx->SData_d.invJ = NULL; 1097 ctx->jacobian_field_major_order = PETSC_FALSE; 1098 ctx->SData_d.coo_elem_offsets = NULL; 1099 ctx->SData_d.coo_elem_point_offsets = NULL; 1100 ctx->coo_assembly = PETSC_FALSE; 1101 ctx->SData_d.coo_elem_fullNb = NULL; 1102 ctx->SData_d.coo_size = 0; 1103 PetscOptionsBegin(ctx->comm, prefix, "Options for Fokker-Plank-Landau collision operator", "none"); 1104 { 1105 char opstring[256]; 1106 #if defined(PETSC_HAVE_KOKKOS_KERNELS) 1107 ctx->deviceType = LANDAU_KOKKOS; 1108 PetscCall(PetscStrncpy(opstring, "kokkos", sizeof(opstring))); 1109 #elif defined(PETSC_HAVE_CUDA) 1110 ctx->deviceType = LANDAU_CUDA; 1111 PetscCall(PetscStrncpy(opstring, "cuda", sizeof(opstring))); 1112 #else 1113 ctx->deviceType = LANDAU_CPU; 1114 PetscCall(PetscStrncpy(opstring, "cpu", sizeof(opstring))); 1115 #endif 1116 PetscCall(PetscOptionsString("-dm_landau_device_type", "Use kernels on 'cpu', 'cuda', or 'kokkos'", "plexland.c", opstring, opstring, sizeof(opstring), NULL)); 1117 PetscCall(PetscStrcmp("cpu", opstring, &flg)); 1118 if (flg) { 1119 ctx->deviceType = LANDAU_CPU; 1120 } else { 1121 PetscCall(PetscStrcmp("cuda", opstring, &flg)); 1122 if (flg) { 1123 ctx->deviceType = LANDAU_CUDA; 1124 } else { 1125 PetscCall(PetscStrcmp("kokkos", opstring, &flg)); 1126 if (flg) ctx->deviceType = LANDAU_KOKKOS; 1127 else SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_device_type %s", opstring); 1128 } 1129 } 1130 } 1131 PetscCall(PetscOptionsReal("-dm_landau_electron_shift", "Shift in thermal velocity of electrons", "none", ctx->electronShift, &ctx->electronShift, NULL)); 1132 PetscCall(PetscOptionsInt("-dm_landau_verbose", "Level of verbosity output", "plexland.c", ctx->verbose, &ctx->verbose, NULL)); 1133 PetscCall(PetscOptionsInt("-dm_landau_batch_size", "Number of 'vertices' to batch", "ex2.c", ctx->batch_sz, &ctx->batch_sz, NULL)); 1134 PetscCheck(LANDAU_MAX_BATCH_SZ >= ctx->batch_sz, ctx->comm, PETSC_ERR_ARG_WRONG, "LANDAU_MAX_BATCH_SZ %" PetscInt_FMT " < ctx->batch_sz %" PetscInt_FMT, (PetscInt)LANDAU_MAX_BATCH_SZ, ctx->batch_sz); 1135 PetscCall(PetscOptionsInt("-dm_landau_batch_view_idx", "Index of batch for diagnostics like plotting", "ex2.c", ctx->batch_view_idx, &ctx->batch_view_idx, NULL)); 1136 PetscCheck(ctx->batch_view_idx < ctx->batch_sz, ctx->comm, PETSC_ERR_ARG_WRONG, "-ctx->batch_view_idx %" PetscInt_FMT " > ctx->batch_sz %" PetscInt_FMT, ctx->batch_view_idx, ctx->batch_sz); 1137 PetscCall(PetscOptionsReal("-dm_landau_Ez", "Initial parallel electric field in unites of Conner-Hastie critical field", "plexland.c", ctx->Ez, &ctx->Ez, NULL)); 1138 PetscCall(PetscOptionsReal("-dm_landau_n_0", "Normalization constant for number density", "plexland.c", ctx->n_0, &ctx->n_0, NULL)); 1139 PetscCall(PetscOptionsBool("-dm_landau_use_mataxpy_mass", "Use fast but slightly fragile MATAXPY to add mass term", "plexland.c", ctx->use_matrix_mass, &ctx->use_matrix_mass, NULL)); 1140 PetscCall(PetscOptionsBool("-dm_landau_use_relativistic_corrections", "Use relativistic corrections", "plexland.c", ctx->use_relativistic_corrections, &ctx->use_relativistic_corrections, NULL)); 1141 if (LANDAU_DIM == 2 && ctx->use_relativistic_corrections) ctx->use_relativistic_corrections = PETSC_FALSE; // should warn 1142 PetscCall(PetscOptionsBool("-dm_landau_use_energy_tensor_trick", "Use Eero's trick of using grad(v^2/2) instead of v as args to Landau tensor to conserve energy with relativistic corrections and Q1 elements", "plexland.c", ctx->use_energy_tensor_trick, 1143 &ctx->use_energy_tensor_trick, NULL)); 1144 1145 /* get num species with temperature, set defaults */ 1146 for (ii = 1; ii < LANDAU_MAX_SPECIES; ii++) { 1147 ctx->thermal_temps[ii] = 1; 1148 ctx->charges[ii] = 1; 1149 ctx->masses[ii] = 1; 1150 ctx->n[ii] = 1; 1151 } 1152 nt = LANDAU_MAX_SPECIES; 1153 PetscCall(PetscOptionsRealArray("-dm_landau_thermal_temps", "Temperature of each species [e,i_0,i_1,...] in keV (must be set to set number of species)", "plexland.c", ctx->thermal_temps, &nt, &flg)); 1154 if (flg) { 1155 PetscCall(PetscInfo(dummy, "num_species set to number of thermal temps provided (%" PetscInt_FMT ")\n", nt)); 1156 ctx->num_species = nt; 1157 } else SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_thermal_temps ,t1,t2,.. must be provided to set the number of species"); 1158 for (ii = 0; ii < ctx->num_species; ii++) ctx->thermal_temps[ii] *= 1.1604525e7; /* convert to Kelvin */ 1159 nm = LANDAU_MAX_SPECIES - 1; 1160 PetscCall(PetscOptionsRealArray("-dm_landau_ion_masses", "Mass of each species in units of proton mass [i_0=2,i_1=40...]", "plexland.c", &ctx->masses[1], &nm, &flg)); 1161 PetscCheck(!flg || nm == ctx->num_species - 1, ctx->comm, PETSC_ERR_ARG_WRONG, "num ion masses %" PetscInt_FMT " != num species %" PetscInt_FMT, nm, ctx->num_species - 1); 1162 nm = LANDAU_MAX_SPECIES; 1163 PetscCall(PetscOptionsRealArray("-dm_landau_n", "Number density of each species = n_s * n_0", "plexland.c", ctx->n, &nm, &flg)); 1164 PetscCheck(!flg || nm == ctx->num_species, ctx->comm, PETSC_ERR_ARG_WRONG, "wrong num n: %" PetscInt_FMT " != num species %" PetscInt_FMT, nm, ctx->num_species); 1165 for (ii = 0; ii < LANDAU_MAX_SPECIES; ii++) ctx->masses[ii] *= 1.6720e-27; /* scale by proton mass kg */ 1166 ctx->masses[0] = 9.10938356e-31; /* electron mass kg (should be about right already) */ 1167 nc = LANDAU_MAX_SPECIES - 1; 1168 PetscCall(PetscOptionsRealArray("-dm_landau_ion_charges", "Charge of each species in units of proton charge [i_0=2,i_1=18,...]", "plexland.c", &ctx->charges[1], &nc, &flg)); 1169 if (flg) PetscCheck(nc == ctx->num_species - 1, ctx->comm, PETSC_ERR_ARG_WRONG, "num charges %" PetscInt_FMT " != num species %" PetscInt_FMT, nc, ctx->num_species - 1); 1170 for (ii = 0; ii < LANDAU_MAX_SPECIES; ii++) ctx->charges[ii] *= 1.6022e-19; /* electron/proton charge (MKS) */ 1171 /* geometry and grids */ 1172 nt = LANDAU_MAX_GRIDS; 1173 PetscCall(PetscOptionsIntArray("-dm_landau_num_species_grid", "Number of species on each grid: [ 1, ....] or [S, 0 ....] for single grid", "plexland.c", num_species_grid, &nt, &flg)); 1174 if (flg) { 1175 ctx->num_grids = nt; 1176 for (ii = nt = 0; ii < ctx->num_grids; ii++) nt += num_species_grid[ii]; 1177 PetscCheck(ctx->num_species == nt, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_num_species_grid: sum %" PetscInt_FMT " != num_species = %" PetscInt_FMT ". %" PetscInt_FMT " grids (check that number of grids <= LANDAU_MAX_GRIDS = %d)", nt, ctx->num_species, 1178 ctx->num_grids, LANDAU_MAX_GRIDS); 1179 } else { 1180 ctx->num_grids = 1; // go back to a single grid run 1181 num_species_grid[0] = ctx->num_species; 1182 } 1183 for (ctx->species_offset[0] = ii = 0; ii < ctx->num_grids; ii++) ctx->species_offset[ii + 1] = ctx->species_offset[ii] + num_species_grid[ii]; 1184 PetscCheck(ctx->species_offset[ctx->num_grids] == ctx->num_species, ctx->comm, PETSC_ERR_ARG_WRONG, "ctx->species_offset[ctx->num_grids] %" PetscInt_FMT " != ctx->num_species = %" PetscInt_FMT " ???????????", ctx->species_offset[ctx->num_grids], 1185 ctx->num_species); 1186 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 1187 int iii = ctx->species_offset[grid]; // normalize with first (arbitrary) species on grid 1188 v0_grid[grid] = PetscSqrtReal(ctx->k * ctx->thermal_temps[iii] / ctx->masses[iii]); /* arbitrary units for non-dimensionalization: plasma formulary def */ 1189 } 1190 // get lambdas here because we need them for t_0 etc 1191 PetscCall(PetscOptionsReal("-dm_landau_ln_lambda", "Universal cross section parameter. Default uses NRL formulas", "plexland.c", lnLam, &lnLam, &flg)); 1192 if (flg) { 1193 for (PetscInt grid = 0; grid < LANDAU_MAX_GRIDS; grid++) { 1194 for (PetscInt gridj = 0; gridj < LANDAU_MAX_GRIDS; gridj++) ctx->lambdas[gridj][grid] = lnLam; /* cross section ratio large - small angle collisions */ 1195 } 1196 } else { 1197 PetscCall(makeLambdas(ctx)); 1198 } 1199 non_dim_grid = 0; 1200 PetscCall(PetscOptionsInt("-dm_landau_normalization_grid", "Index of grid to use for setting v_0, m_0, t_0. (Not recommended)", "plexland.c", non_dim_grid, &non_dim_grid, &flg)); 1201 if (non_dim_grid != 0) PetscCall(PetscInfo(dummy, "Normalization grid set to %" PetscInt_FMT ", but non-default not well verified\n", non_dim_grid)); 1202 PetscCheck(non_dim_grid >= 0 && non_dim_grid < ctx->num_species, ctx->comm, PETSC_ERR_ARG_WRONG, "Normalization grid wrong: %" PetscInt_FMT, non_dim_grid); 1203 ctx->v_0 = v0_grid[non_dim_grid]; /* arbitrary units for non dimensionalization: global mean velocity in 1D of electrons */ 1204 ctx->m_0 = ctx->masses[non_dim_grid]; /* arbitrary reference mass, electrons */ 1205 ctx->t_0 = 8 * PETSC_PI * PetscSqr(ctx->epsilon0 * ctx->m_0 / PetscSqr(ctx->charges[non_dim_grid])) / ctx->lambdas[non_dim_grid][non_dim_grid] / ctx->n_0 * PetscPowReal(ctx->v_0, 3); /* note, this t_0 makes nu[non_dim_grid,non_dim_grid]=1 */ 1206 /* domain */ 1207 nt = LANDAU_MAX_GRIDS; 1208 PetscCall(PetscOptionsRealArray("-dm_landau_domain_radius", "Phase space size in units of thermal velocity of grid", "plexland.c", ctx->radius, &nt, &flg)); 1209 if (flg) { 1210 PetscCheck(nt >= ctx->num_grids, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_domain_radius: given %" PetscInt_FMT " radius != number grids %" PetscInt_FMT, nt, ctx->num_grids); 1211 while (nt--) ctx->radius_par[nt] = ctx->radius_perp[nt] = ctx->radius[nt]; 1212 } else { 1213 nt = LANDAU_MAX_GRIDS; 1214 PetscCall(PetscOptionsRealArray("-dm_landau_domain_max_par", "Parallel velocity domain size in units of thermal velocity of grid", "plexland.c", ctx->radius_par, &nt, &flg)); 1215 if (flg) PetscCheck(nt >= ctx->num_grids, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_domain_max_par: given %" PetscInt_FMT " radius != number grids %" PetscInt_FMT, nt, ctx->num_grids); 1216 PetscCall(PetscOptionsRealArray("-dm_landau_domain_max_perp", "Perpendicular velocity domain size in units of thermal velocity of grid", "plexland.c", ctx->radius_perp, &nt, &flg)); 1217 if (flg) PetscCheck(nt >= ctx->num_grids, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_domain_max_perp: given %" PetscInt_FMT " radius != number grids %" PetscInt_FMT, nt, ctx->num_grids); 1218 } 1219 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 1220 if (flg && ctx->radius[grid] <= 0) { /* negative is ratio of c - need to set par and perp with this -- todo */ 1221 if (ctx->radius[grid] == 0) ctx->radius[grid] = 0.75; 1222 else ctx->radius[grid] = -ctx->radius[grid]; 1223 ctx->radius[grid] = ctx->radius[grid] * SPEED_OF_LIGHT / ctx->v_0; // use any species on grid to normalize (v_0 same for all on grid) 1224 PetscCall(PetscInfo(dummy, "Change domain radius to %g for grid %" PetscInt_FMT "\n", (double)ctx->radius[grid], grid)); 1225 } 1226 ctx->radius[grid] *= v0_grid[grid] / ctx->v_0; // scale domain by thermal radius relative to v_0 1227 ctx->radius_perp[grid] *= v0_grid[grid] / ctx->v_0; // scale domain by thermal radius relative to v_0 1228 ctx->radius_par[grid] *= v0_grid[grid] / ctx->v_0; // scale domain by thermal radius relative to v_0 1229 } 1230 /* amr parameters */ 1231 nt = LANDAU_MAX_GRIDS; 1232 PetscCall(PetscOptionsIntArray("-dm_landau_amr_levels_max", "Number of AMR levels of refinement around origin, after (RE) refinements along z", "plexland.c", ctx->numAMRRefine, &nt, &flg)); 1233 PetscCheck(!flg || nt >= ctx->num_grids, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_amr_levels_max: given %" PetscInt_FMT " != number grids %" PetscInt_FMT, nt, ctx->num_grids); 1234 nt = LANDAU_MAX_GRIDS; 1235 PetscCall(PetscOptionsIntArray("-dm_landau_amr_post_refine", "Number of levels to uniformly refine after AMR", "plexland.c", ctx->postAMRRefine, &nt, &flg)); 1236 for (ii = 1; ii < ctx->num_grids; ii++) ctx->postAMRRefine[ii] = ctx->postAMRRefine[0]; // all grids the same now 1237 PetscCall(PetscOptionsInt("-dm_landau_amr_re_levels", "Number of levels to refine along v_perp=0, z>0", "plexland.c", ctx->numRERefine, &ctx->numRERefine, &flg)); 1238 PetscCall(PetscOptionsInt("-dm_landau_amr_z_refine_pre", "Number of levels to refine along v_perp=0 before origin refine", "plexland.c", ctx->nZRefine1, &ctx->nZRefine1, &flg)); 1239 PetscCall(PetscOptionsInt("-dm_landau_amr_z_refine_post", "Number of levels to refine along v_perp=0 after origin refine", "plexland.c", ctx->nZRefine2, &ctx->nZRefine2, &flg)); 1240 PetscCall(PetscOptionsReal("-dm_landau_re_radius", "velocity range to refine on positive (z>0) r=0 axis for runaways", "plexland.c", ctx->re_radius, &ctx->re_radius, &flg)); 1241 PetscCall(PetscOptionsReal("-dm_landau_z_radius_pre", "velocity range to refine r=0 axis (for electrons)", "plexland.c", ctx->vperp0_radius1, &ctx->vperp0_radius1, &flg)); 1242 PetscCall(PetscOptionsReal("-dm_landau_z_radius_post", "velocity range to refine r=0 axis (for electrons) after origin AMR", "plexland.c", ctx->vperp0_radius2, &ctx->vperp0_radius2, &flg)); 1243 /* spherical domain (not used) */ 1244 PetscCall(PetscOptionsBool("-dm_landau_sphere", "use sphere/semi-circle domain instead of rectangle", "plexland.c", ctx->sphere, &ctx->sphere, NULL)); 1245 if (ctx->sphere) { 1246 ctx->sphere_inner_radius_90degree = 0.50; 1247 ctx->sphere_inner_radius_45degree = 0.45; 1248 PetscCall(PetscOptionsReal("-dm_landau_sphere_inner_radius_90degree_scale", "Scaling of radius for inner circle on 90 degree grid", "plexland.c", ctx->sphere_inner_radius_90degree, &ctx->sphere_inner_radius_90degree, NULL)); 1249 PetscCall(PetscOptionsReal("-dm_landau_sphere_inner_radius_45degree_scale", "Scaling of radius for inner circle on 45 degree grid", "plexland.c", ctx->sphere_inner_radius_45degree, &ctx->sphere_inner_radius_45degree, NULL)); 1250 } else { 1251 nt = LANDAU_DIM; 1252 PetscCall(PetscOptionsIntArray("-dm_landau_num_cells", "Number of cells in each dimension of base grid", "plexland.c", ctx->cells0, &nt, &flg)); 1253 } 1254 /* processing options */ 1255 PetscCall(PetscOptionsBool("-dm_landau_gpu_assembly", "Assemble Jacobian on GPU", "plexland.c", ctx->gpu_assembly, &ctx->gpu_assembly, NULL)); 1256 if (ctx->deviceType == LANDAU_CPU || ctx->deviceType == LANDAU_KOKKOS) { // make Kokkos 1257 PetscCall(PetscOptionsBool("-dm_landau_coo_assembly", "Assemble Jacobian with Kokkos on 'device'", "plexland.c", ctx->coo_assembly, &ctx->coo_assembly, NULL)); 1258 if (ctx->coo_assembly) PetscCheck(ctx->gpu_assembly, ctx->comm, PETSC_ERR_ARG_WRONG, "COO assembly requires 'gpu assembly' even if Kokkos 'CPU' back-end %d", ctx->coo_assembly); 1259 } 1260 PetscCall(PetscOptionsBool("-dm_landau_jacobian_field_major_order", "Reorder Jacobian for GPU assembly with field major, or block diagonal, ordering (DEPRECATED)", "plexland.c", ctx->jacobian_field_major_order, &ctx->jacobian_field_major_order, NULL)); 1261 if (ctx->jacobian_field_major_order) PetscCheck(ctx->gpu_assembly, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_jacobian_field_major_order requires -dm_landau_gpu_assembly"); 1262 PetscCheck(!ctx->jacobian_field_major_order, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_jacobian_field_major_order DEPRECATED"); 1263 PetscOptionsEnd(); 1264 1265 for (ii = ctx->num_species; ii < LANDAU_MAX_SPECIES; ii++) ctx->masses[ii] = ctx->thermal_temps[ii] = ctx->charges[ii] = 0; 1266 if (ctx->verbose != 0) { 1267 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "masses: e=%10.3e; ions in proton mass units: %10.3e %10.3e ...\n", (double)ctx->masses[0], (double)(ctx->masses[1] / 1.6720e-27), (double)(ctx->num_species > 2 ? ctx->masses[2] / 1.6720e-27 : 0))); 1268 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "charges: e=%10.3e; charges in elementary units: %10.3e %10.3e\n", (double)ctx->charges[0], (double)(-ctx->charges[1] / ctx->charges[0]), (double)(ctx->num_species > 2 ? -ctx->charges[2] / ctx->charges[0] : 0))); 1269 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "n: e: %10.3e i: %10.3e %10.3e\n", (double)ctx->n[0], (double)ctx->n[1], (double)(ctx->num_species > 2 ? ctx->n[2] : 0))); 1270 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "thermal T (K): e=%10.3e i=%10.3e %10.3e. Normalization grid %d: v_0=%10.3e (%10.3ec) n_0=%10.3e t_0=%10.3e %" PetscInt_FMT " batched, view batch %" PetscInt_FMT "\n", (double)ctx->thermal_temps[0], 1271 (double)ctx->thermal_temps[1], (double)((ctx->num_species > 2) ? ctx->thermal_temps[2] : 0), (int)non_dim_grid, (double)ctx->v_0, (double)(ctx->v_0 / SPEED_OF_LIGHT), (double)ctx->n_0, (double)ctx->t_0, ctx->batch_sz, ctx->batch_view_idx)); 1272 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Domain radius (AMR levels) grid %d: par=%10.3e perp=%10.3e (%" PetscInt_FMT ") ", 0, (double)ctx->radius_par[0], (double)ctx->radius_perp[0], ctx->numAMRRefine[0])); 1273 for (ii = 1; ii < ctx->num_grids; ii++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", %" PetscInt_FMT ": par=%10.3e perp=%10.3e (%" PetscInt_FMT ") ", ii, (double)ctx->radius_par[ii], (double)ctx->radius_perp[ii], ctx->numAMRRefine[ii])); 1274 if (ctx->use_relativistic_corrections) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nUse relativistic corrections\n")); 1275 else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n")); 1276 } 1277 PetscCall(DMDestroy(&dummy)); 1278 { 1279 PetscMPIInt rank; 1280 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1281 ctx->stage = 0; 1282 PetscCall(PetscLogEventRegister("Landau Create", DM_CLASSID, &ctx->events[13])); /* 13 */ 1283 PetscCall(PetscLogEventRegister(" GPU ass. setup", DM_CLASSID, &ctx->events[2])); /* 2 */ 1284 PetscCall(PetscLogEventRegister(" Build matrix", DM_CLASSID, &ctx->events[12])); /* 12 */ 1285 PetscCall(PetscLogEventRegister(" Assembly maps", DM_CLASSID, &ctx->events[15])); /* 15 */ 1286 PetscCall(PetscLogEventRegister("Landau Mass mat", DM_CLASSID, &ctx->events[14])); /* 14 */ 1287 PetscCall(PetscLogEventRegister("Landau Operator", DM_CLASSID, &ctx->events[11])); /* 11 */ 1288 PetscCall(PetscLogEventRegister("Landau Jacobian", DM_CLASSID, &ctx->events[0])); /* 0 */ 1289 PetscCall(PetscLogEventRegister("Landau Mass", DM_CLASSID, &ctx->events[9])); /* 9 */ 1290 PetscCall(PetscLogEventRegister(" Preamble", DM_CLASSID, &ctx->events[10])); /* 10 */ 1291 PetscCall(PetscLogEventRegister(" static IP Data", DM_CLASSID, &ctx->events[7])); /* 7 */ 1292 PetscCall(PetscLogEventRegister(" dynamic IP-Jac", DM_CLASSID, &ctx->events[1])); /* 1 */ 1293 PetscCall(PetscLogEventRegister(" Kernel-init", DM_CLASSID, &ctx->events[3])); /* 3 */ 1294 PetscCall(PetscLogEventRegister(" Jac-f-df (GPU)", DM_CLASSID, &ctx->events[8])); /* 8 */ 1295 PetscCall(PetscLogEventRegister(" J Kernel (GPU)", DM_CLASSID, &ctx->events[4])); /* 4 */ 1296 PetscCall(PetscLogEventRegister(" M Kernel (GPU)", DM_CLASSID, &ctx->events[16])); /* 16 */ 1297 PetscCall(PetscLogEventRegister(" Copy to CPU", DM_CLASSID, &ctx->events[5])); /* 5 */ 1298 PetscCall(PetscLogEventRegister(" CPU assemble", DM_CLASSID, &ctx->events[6])); /* 6 */ 1299 1300 if (rank) { /* turn off output stuff for duplicate runs - do we need to add the prefix to all this? */ 1301 PetscCall(PetscOptionsClearValue(NULL, "-snes_converged_reason")); 1302 PetscCall(PetscOptionsClearValue(NULL, "-ksp_converged_reason")); 1303 PetscCall(PetscOptionsClearValue(NULL, "-snes_monitor")); 1304 PetscCall(PetscOptionsClearValue(NULL, "-ksp_monitor")); 1305 PetscCall(PetscOptionsClearValue(NULL, "-ts_monitor")); 1306 PetscCall(PetscOptionsClearValue(NULL, "-ts_view")); 1307 PetscCall(PetscOptionsClearValue(NULL, "-ts_adapt_monitor")); 1308 PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_amr_dm_view")); 1309 PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_amr_vec_view")); 1310 PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_mass_dm_view")); 1311 PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_mass_view")); 1312 PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_jacobian_view")); 1313 PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_mat_view")); 1314 PetscCall(PetscOptionsClearValue(NULL, "-pc_bjkokkos_ksp_converged_reason")); 1315 PetscCall(PetscOptionsClearValue(NULL, "-pc_bjkokkos_ksp_monitor")); 1316 PetscCall(PetscOptionsClearValue(NULL, "-")); 1317 PetscCall(PetscOptionsClearValue(NULL, "-info")); 1318 } 1319 } 1320 PetscFunctionReturn(PETSC_SUCCESS); 1321 } 1322 1323 static PetscErrorCode CreateStaticData(PetscInt dim, IS grid_batch_is_inv[], LandauCtx *ctx) 1324 { 1325 PetscSection section[LANDAU_MAX_GRIDS], globsection[LANDAU_MAX_GRIDS]; 1326 PetscQuadrature quad; 1327 const PetscReal *quadWeights; 1328 PetscReal invMass[LANDAU_MAX_SPECIES], nu_alpha[LANDAU_MAX_SPECIES], nu_beta[LANDAU_MAX_SPECIES]; 1329 PetscInt numCells[LANDAU_MAX_GRIDS], Nq, Nf[LANDAU_MAX_GRIDS], ncellsTot = 0, MAP_BF_SIZE = 64 * LANDAU_DIM * LANDAU_DIM * LANDAU_MAX_Q_FACE * LANDAU_MAX_SPECIES; 1330 PetscTabulation *Tf; 1331 PetscDS prob; 1332 1333 PetscFunctionBegin; 1334 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 1335 for (PetscInt ii = ctx->species_offset[grid]; ii < ctx->species_offset[grid + 1]; ii++) { 1336 invMass[ii] = ctx->m_0 / ctx->masses[ii]; 1337 nu_alpha[ii] = PetscSqr(ctx->charges[ii] / ctx->m_0) * ctx->m_0 / ctx->masses[ii]; 1338 nu_beta[ii] = PetscSqr(ctx->charges[ii] / ctx->epsilon0) / (8 * PETSC_PI) * ctx->t_0 * ctx->n_0 / PetscPowReal(ctx->v_0, 3); 1339 } 1340 } 1341 if (ctx->verbose == 4) { 1342 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "nu_alpha: ")); 1343 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 1344 int iii = ctx->species_offset[grid]; 1345 for (PetscInt ii = iii; ii < ctx->species_offset[grid + 1]; ii++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %e", (double)nu_alpha[ii])); 1346 } 1347 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nnu_beta: ")); 1348 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 1349 int iii = ctx->species_offset[grid]; 1350 for (PetscInt ii = iii; ii < ctx->species_offset[grid + 1]; ii++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %e", (double)nu_beta[ii])); 1351 } 1352 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nnu_alpha[i]*nu_beta[j]*lambda[i][j]:\n")); 1353 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 1354 int iii = ctx->species_offset[grid]; 1355 for (PetscInt ii = iii; ii < ctx->species_offset[grid + 1]; ii++) { 1356 for (PetscInt gridj = 0; gridj < ctx->num_grids; gridj++) { 1357 int jjj = ctx->species_offset[gridj]; 1358 for (PetscInt jj = jjj; jj < ctx->species_offset[gridj + 1]; jj++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %14.9e", (double)(nu_alpha[ii] * nu_beta[jj] * ctx->lambdas[grid][gridj]))); 1359 } 1360 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n")); 1361 } 1362 } 1363 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "lambda[i][j]:\n")); 1364 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 1365 int iii = ctx->species_offset[grid]; 1366 for (PetscInt ii = iii; ii < ctx->species_offset[grid + 1]; ii++) { 1367 for (PetscInt gridj = 0; gridj < ctx->num_grids; gridj++) { 1368 int jjj = ctx->species_offset[gridj]; 1369 for (PetscInt jj = jjj; jj < ctx->species_offset[gridj + 1]; jj++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %14.9e", (double)ctx->lambdas[grid][gridj])); 1370 } 1371 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n")); 1372 } 1373 } 1374 } 1375 PetscCall(DMGetDS(ctx->plex[0], &prob)); // same DS for all grids 1376 PetscCall(PetscDSGetTabulation(prob, &Tf)); // Bf, &Df same for all grids 1377 /* DS, Tab and quad is same on all grids */ 1378 PetscCheck(ctx->plex[0], ctx->comm, PETSC_ERR_ARG_WRONG, "Plex not created"); 1379 PetscCall(PetscFEGetQuadrature(ctx->fe[0], &quad)); 1380 PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, &quadWeights)); 1381 PetscCheck(Nq <= LANDAU_MAX_NQ, ctx->comm, PETSC_ERR_ARG_WRONG, "Order too high. Nq = %" PetscInt_FMT " > LANDAU_MAX_NQ (%d)", Nq, LANDAU_MAX_NQ); 1382 /* setup each grid */ 1383 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 1384 PetscInt cStart, cEnd; 1385 PetscCheck(ctx->plex[grid] != NULL, ctx->comm, PETSC_ERR_ARG_WRONG, "Plex not created"); 1386 PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, &cEnd)); 1387 numCells[grid] = cEnd - cStart; // grids can have different topology 1388 PetscCall(DMGetLocalSection(ctx->plex[grid], §ion[grid])); 1389 PetscCall(DMGetGlobalSection(ctx->plex[grid], &globsection[grid])); 1390 PetscCall(PetscSectionGetNumFields(section[grid], &Nf[grid])); 1391 ncellsTot += numCells[grid]; 1392 } 1393 /* create GPU assembly data */ 1394 if (ctx->gpu_assembly) { /* we need GPU object with GPU assembly */ 1395 PetscContainer container; 1396 PetscScalar *elemMatrix, *elMat; 1397 pointInterpolationP4est(*pointMaps)[LANDAU_MAX_Q_FACE]; 1398 P4estVertexMaps *maps; 1399 const PetscInt *plex_batch = NULL, Nb = Nq, elMatSz = Nq * Nq * ctx->num_species * ctx->num_species; // tensor elements; 1400 LandauIdx *coo_elem_offsets = NULL, *coo_elem_fullNb = NULL, (*coo_elem_point_offsets)[LANDAU_MAX_NQ + 1] = NULL; 1401 /* create GPU assembly data */ 1402 PetscCall(PetscInfo(ctx->plex[0], "Make GPU maps %d\n", 1)); 1403 PetscCall(PetscLogEventBegin(ctx->events[2], 0, 0, 0, 0)); 1404 PetscCall(PetscMalloc(sizeof(*maps) * ctx->num_grids, &maps)); 1405 PetscCall(PetscMalloc(sizeof(*pointMaps) * MAP_BF_SIZE, &pointMaps)); 1406 PetscCall(PetscMalloc(sizeof(*elemMatrix) * elMatSz, &elemMatrix)); 1407 1408 if (ctx->coo_assembly) { // setup COO assembly -- put COO metadata directly in ctx->SData_d 1409 PetscCall(PetscMalloc3(ncellsTot + 1, &coo_elem_offsets, ncellsTot, &coo_elem_fullNb, ncellsTot, &coo_elem_point_offsets)); // array of integer pointers 1410 coo_elem_offsets[0] = 0; // finish later 1411 PetscCall(PetscInfo(ctx->plex[0], "COO initialization, %" PetscInt_FMT " cells\n", ncellsTot)); 1412 ctx->SData_d.coo_n_cellsTot = ncellsTot; 1413 ctx->SData_d.coo_elem_offsets = (void *)coo_elem_offsets; 1414 ctx->SData_d.coo_elem_fullNb = (void *)coo_elem_fullNb; 1415 ctx->SData_d.coo_elem_point_offsets = (void *)coo_elem_point_offsets; 1416 } else { 1417 ctx->SData_d.coo_elem_offsets = ctx->SData_d.coo_elem_fullNb = NULL; 1418 ctx->SData_d.coo_elem_point_offsets = NULL; 1419 ctx->SData_d.coo_n_cellsTot = 0; 1420 } 1421 1422 ctx->SData_d.coo_max_fullnb = 0; 1423 for (PetscInt grid = 0, glb_elem_idx = 0; grid < ctx->num_grids; grid++) { 1424 PetscInt cStart, cEnd, Nfloc = Nf[grid], totDim = Nfloc * Nq; 1425 if (grid_batch_is_inv[grid]) PetscCall(ISGetIndices(grid_batch_is_inv[grid], &plex_batch)); 1426 PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, &cEnd)); 1427 // make maps 1428 maps[grid].d_self = NULL; 1429 maps[grid].num_elements = numCells[grid]; 1430 maps[grid].num_face = (PetscInt)(pow(Nq, 1. / ((double)dim)) + .001); // Q 1431 maps[grid].num_face = (PetscInt)(pow(maps[grid].num_face, (double)(dim - 1)) + .001); // Q^2 1432 maps[grid].num_reduced = 0; 1433 maps[grid].deviceType = ctx->deviceType; 1434 maps[grid].numgrids = ctx->num_grids; 1435 // count reduced and get 1436 PetscCall(PetscMalloc(maps[grid].num_elements * sizeof(*maps[grid].gIdx), &maps[grid].gIdx)); 1437 for (int ej = cStart, eidx = 0; ej < cEnd; ++ej, ++eidx, glb_elem_idx++) { 1438 if (coo_elem_offsets) coo_elem_offsets[glb_elem_idx + 1] = coo_elem_offsets[glb_elem_idx]; // start with last one, then add 1439 for (int fieldA = 0; fieldA < Nf[grid]; fieldA++) { 1440 int fullNb = 0; 1441 for (int q = 0; q < Nb; ++q) { 1442 PetscInt numindices, *indices; 1443 PetscScalar *valuesOrig = elMat = elemMatrix; 1444 PetscCall(PetscArrayzero(elMat, totDim * totDim)); 1445 elMat[(fieldA * Nb + q) * totDim + fieldA * Nb + q] = 1; 1446 PetscCall(DMPlexGetClosureIndices(ctx->plex[grid], section[grid], globsection[grid], ej, PETSC_TRUE, &numindices, &indices, NULL, (PetscScalar **)&elMat)); 1447 for (PetscInt f = 0; f < numindices; ++f) { // look for a non-zero on the diagonal 1448 if (PetscAbs(PetscRealPart(elMat[f * numindices + f])) > PETSC_MACHINE_EPSILON) { 1449 // found it 1450 if (PetscAbs(PetscRealPart(elMat[f * numindices + f] - 1.)) < PETSC_MACHINE_EPSILON) { // normal vertex 1.0 1451 if (plex_batch) { 1452 maps[grid].gIdx[eidx][fieldA][q] = (LandauIdx)plex_batch[indices[f]]; 1453 } else { 1454 maps[grid].gIdx[eidx][fieldA][q] = (LandauIdx)indices[f]; 1455 } 1456 fullNb++; 1457 } else { //found a constraint 1458 int jj = 0; 1459 PetscReal sum = 0; 1460 const PetscInt ff = f; 1461 maps[grid].gIdx[eidx][fieldA][q] = -maps[grid].num_reduced - 1; // store (-)index: id = -(idx+1): idx = -id - 1 1462 1463 do { // constraints are continuous in Plex - exploit that here 1464 int ii; // get 'scale' 1465 for (ii = 0, pointMaps[maps[grid].num_reduced][jj].scale = 0; ii < maps[grid].num_face; ii++) { // sum row of outer product to recover vector value 1466 if (ff + ii < numindices) { // 3D has Q and Q^2 interps so might run off end. We could test that elMat[f*numindices + ff + ii] > 0, and break if not 1467 pointMaps[maps[grid].num_reduced][jj].scale += PetscRealPart(elMat[f * numindices + ff + ii]); 1468 } 1469 } 1470 sum += pointMaps[maps[grid].num_reduced][jj].scale; // diagnostic 1471 // get 'gid' 1472 if (pointMaps[maps[grid].num_reduced][jj].scale == 0) pointMaps[maps[grid].num_reduced][jj].gid = -1; // 3D has Q and Q^2 interps 1473 else { 1474 if (plex_batch) { 1475 pointMaps[maps[grid].num_reduced][jj].gid = plex_batch[indices[f]]; 1476 } else { 1477 pointMaps[maps[grid].num_reduced][jj].gid = indices[f]; 1478 } 1479 fullNb++; 1480 } 1481 } while (++jj < maps[grid].num_face && ++f < numindices); // jj is incremented if we hit the end 1482 while (jj < maps[grid].num_face) { 1483 pointMaps[maps[grid].num_reduced][jj].scale = 0; 1484 pointMaps[maps[grid].num_reduced][jj].gid = -1; 1485 jj++; 1486 } 1487 if (PetscAbs(sum - 1.0) > 10 * PETSC_MACHINE_EPSILON) { // debug 1488 int d, f; 1489 PetscReal tmp = 0; 1490 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\t\t%d.%d.%d) ERROR total I = %22.16e (LANDAU_MAX_Q_FACE=%d, #face=%d)\n", eidx, q, fieldA, (double)sum, LANDAU_MAX_Q_FACE, maps[grid].num_face)); 1491 for (d = 0, tmp = 0; d < numindices; ++d) { 1492 if (tmp != 0 && PetscAbs(tmp - 1.0) > 10 * PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%3d) %3" PetscInt_FMT ": ", d, indices[d])); 1493 for (f = 0; f < numindices; ++f) tmp += PetscRealPart(elMat[d * numindices + f]); 1494 if (tmp != 0) PetscCall(PetscPrintf(ctx->comm, " | %22.16e\n", (double)tmp)); 1495 } 1496 } 1497 maps[grid].num_reduced++; 1498 PetscCheck(maps[grid].num_reduced < MAP_BF_SIZE, PETSC_COMM_SELF, PETSC_ERR_PLIB, "maps[grid].num_reduced %d > %" PetscInt_FMT, maps[grid].num_reduced, MAP_BF_SIZE); 1499 } 1500 break; 1501 } 1502 } 1503 // cleanup 1504 PetscCall(DMPlexRestoreClosureIndices(ctx->plex[grid], section[grid], globsection[grid], ej, PETSC_TRUE, &numindices, &indices, NULL, (PetscScalar **)&elMat)); 1505 if (elMat != valuesOrig) PetscCall(DMRestoreWorkArray(ctx->plex[grid], numindices * numindices, MPIU_SCALAR, &elMat)); 1506 } 1507 if (ctx->coo_assembly) { // setup COO assembly 1508 coo_elem_offsets[glb_elem_idx + 1] += fullNb * fullNb; // one species block, adds a block for each species, on this element in this grid 1509 if (fieldA == 0) { // cache full Nb for this element, on this grid per species 1510 coo_elem_fullNb[glb_elem_idx] = fullNb; 1511 if (fullNb > ctx->SData_d.coo_max_fullnb) ctx->SData_d.coo_max_fullnb = fullNb; 1512 } else PetscCheck(coo_elem_fullNb[glb_elem_idx] == fullNb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "full element size change with species %d %d", coo_elem_fullNb[glb_elem_idx], fullNb); 1513 } 1514 } // field 1515 } // cell 1516 // allocate and copy point data maps[grid].gIdx[eidx][field][q] 1517 PetscCall(PetscMalloc(maps[grid].num_reduced * sizeof(*maps[grid].c_maps), &maps[grid].c_maps)); 1518 for (int ej = 0; ej < maps[grid].num_reduced; ++ej) { 1519 for (int q = 0; q < maps[grid].num_face; ++q) { 1520 maps[grid].c_maps[ej][q].scale = pointMaps[ej][q].scale; 1521 maps[grid].c_maps[ej][q].gid = pointMaps[ej][q].gid; 1522 } 1523 } 1524 #if defined(PETSC_HAVE_KOKKOS_KERNELS) 1525 if (ctx->deviceType == LANDAU_KOKKOS) { 1526 PetscCall(LandauKokkosCreateMatMaps(maps, pointMaps, Nf, Nq, grid)); // implies Kokkos does 1527 } // else could be CUDA 1528 #endif 1529 #if defined(PETSC_HAVE_CUDA) 1530 if (ctx->deviceType == LANDAU_CUDA) PetscCall(LandauCUDACreateMatMaps(maps, pointMaps, Nf, Nq, grid)); 1531 #endif 1532 if (plex_batch) { 1533 PetscCall(ISRestoreIndices(grid_batch_is_inv[grid], &plex_batch)); 1534 PetscCall(ISDestroy(&grid_batch_is_inv[grid])); // we are done with this 1535 } 1536 } /* grids */ 1537 // finish COO 1538 if (ctx->coo_assembly) { // setup COO assembly 1539 PetscInt *oor, *ooc; 1540 ctx->SData_d.coo_size = coo_elem_offsets[ncellsTot] * ctx->batch_sz; 1541 PetscCall(PetscMalloc2(ctx->SData_d.coo_size, &oor, ctx->SData_d.coo_size, &ooc)); 1542 for (int i = 0; i < ctx->SData_d.coo_size; i++) oor[i] = ooc[i] = -1; 1543 // get 1544 for (int grid = 0, glb_elem_idx = 0; grid < ctx->num_grids; grid++) { 1545 for (int ej = 0; ej < numCells[grid]; ++ej, glb_elem_idx++) { 1546 const int fullNb = coo_elem_fullNb[glb_elem_idx]; 1547 const LandauIdx *const Idxs = &maps[grid].gIdx[ej][0][0]; // just use field-0 maps, They should be the same but this is just for COO storage 1548 coo_elem_point_offsets[glb_elem_idx][0] = 0; 1549 for (int f = 0, cnt2 = 0; f < Nb; f++) { 1550 int idx = Idxs[f]; 1551 coo_elem_point_offsets[glb_elem_idx][f + 1] = coo_elem_point_offsets[glb_elem_idx][f]; // start at last 1552 if (idx >= 0) { 1553 cnt2++; 1554 coo_elem_point_offsets[glb_elem_idx][f + 1]++; // inc 1555 } else { 1556 idx = -idx - 1; 1557 for (int q = 0; q < maps[grid].num_face; q++) { 1558 if (maps[grid].c_maps[idx][q].gid < 0) break; 1559 cnt2++; 1560 coo_elem_point_offsets[glb_elem_idx][f + 1]++; // inc 1561 } 1562 } 1563 PetscCheck(cnt2 <= fullNb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "wrong count %d < %d", fullNb, cnt2); 1564 } 1565 PetscCheck(coo_elem_point_offsets[glb_elem_idx][Nb] == fullNb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "coo_elem_point_offsets size %d != fullNb=%d", coo_elem_point_offsets[glb_elem_idx][Nb], fullNb); 1566 } 1567 } 1568 // set 1569 for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { 1570 for (int grid = 0, glb_elem_idx = 0; grid < ctx->num_grids; grid++) { 1571 const PetscInt moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset); 1572 for (int ej = 0; ej < numCells[grid]; ++ej, glb_elem_idx++) { 1573 const int fullNb = coo_elem_fullNb[glb_elem_idx], fullNb2 = fullNb * fullNb; 1574 // set (i,j) 1575 for (int fieldA = 0; fieldA < Nf[grid]; fieldA++) { 1576 const LandauIdx *const Idxs = &maps[grid].gIdx[ej][fieldA][0]; 1577 int rows[LANDAU_MAX_Q_FACE], cols[LANDAU_MAX_Q_FACE]; 1578 for (int f = 0; f < Nb; ++f) { 1579 const int nr = coo_elem_point_offsets[glb_elem_idx][f + 1] - coo_elem_point_offsets[glb_elem_idx][f]; 1580 if (nr == 1) rows[0] = Idxs[f]; 1581 else { 1582 const int idx = -Idxs[f] - 1; 1583 for (int q = 0; q < nr; q++) rows[q] = maps[grid].c_maps[idx][q].gid; 1584 } 1585 for (int g = 0; g < Nb; ++g) { 1586 const int nc = coo_elem_point_offsets[glb_elem_idx][g + 1] - coo_elem_point_offsets[glb_elem_idx][g]; 1587 if (nc == 1) cols[0] = Idxs[g]; 1588 else { 1589 const int idx = -Idxs[g] - 1; 1590 for (int q = 0; q < nc; q++) cols[q] = maps[grid].c_maps[idx][q].gid; 1591 } 1592 const int idx0 = b_id * coo_elem_offsets[ncellsTot] + 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]; 1593 for (int q = 0, idx = idx0; q < nr; q++) { 1594 for (int d = 0; d < nc; d++, idx++) { 1595 oor[idx] = rows[q] + moffset; 1596 ooc[idx] = cols[d] + moffset; 1597 } 1598 } 1599 } 1600 } 1601 } 1602 } // cell 1603 } // grid 1604 } // batch 1605 PetscCall(MatSetPreallocationCOO(ctx->J, ctx->SData_d.coo_size, oor, ooc)); 1606 PetscCall(PetscFree2(oor, ooc)); 1607 } 1608 PetscCall(PetscFree(pointMaps)); 1609 PetscCall(PetscFree(elemMatrix)); 1610 PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container)); 1611 PetscCall(PetscContainerSetPointer(container, (void *)maps)); 1612 PetscCall(PetscContainerSetUserDestroy(container, LandauGPUMapsDestroy)); 1613 PetscCall(PetscObjectCompose((PetscObject)ctx->J, "assembly_maps", (PetscObject)container)); 1614 PetscCall(PetscContainerDestroy(&container)); 1615 PetscCall(PetscLogEventEnd(ctx->events[2], 0, 0, 0, 0)); 1616 } // end GPU assembly 1617 { /* create static point data, Jacobian called first, only one vertex copy */ 1618 PetscReal *invJe, *ww, *xx, *yy, *zz = NULL, *invJ_a; 1619 PetscInt outer_ipidx, outer_ej, grid, nip_glb = 0; 1620 PetscFE fe; 1621 const PetscInt Nb = Nq; 1622 PetscCall(PetscLogEventBegin(ctx->events[7], 0, 0, 0, 0)); 1623 PetscCall(PetscInfo(ctx->plex[0], "Initialize static data\n")); 1624 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) nip_glb += Nq * numCells[grid]; 1625 /* collect f data, first time is for Jacobian, but make mass now */ 1626 if (ctx->verbose != 0) { 1627 PetscInt ncells = 0, N; 1628 PetscCall(MatGetSize(ctx->J, &N, NULL)); 1629 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) ncells += numCells[grid]; 1630 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d) %s %" PetscInt_FMT " IPs, %" PetscInt_FMT " cells total, Nb=%" PetscInt_FMT ", Nq=%" PetscInt_FMT ", dim=%" PetscInt_FMT ", Tab: Nb=%" PetscInt_FMT " Nf=%" PetscInt_FMT " Np=%" PetscInt_FMT " cdim=%" PetscInt_FMT " N=%" PetscInt_FMT "\n", 0, "FormLandau", nip_glb, ncells, Nb, Nq, dim, Nb, 1631 ctx->num_species, Nb, dim, N)); 1632 } 1633 PetscCall(PetscMalloc4(nip_glb, &ww, nip_glb, &xx, nip_glb, &yy, nip_glb * dim * dim, &invJ_a)); 1634 if (dim == 3) PetscCall(PetscMalloc1(nip_glb, &zz)); 1635 if (ctx->use_energy_tensor_trick) { 1636 PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, NULL, PETSC_DECIDE, &fe)); 1637 PetscCall(PetscObjectSetName((PetscObject)fe, "energy")); 1638 } 1639 /* init each grids static data - no batch */ 1640 for (grid = 0, outer_ipidx = 0, outer_ej = 0; grid < ctx->num_grids; grid++) { // OpenMP (once) 1641 Vec v2_2 = NULL; // projected function: v^2/2 for non-relativistic, gamma... for relativistic 1642 PetscSection e_section; 1643 DM dmEnergy; 1644 PetscInt cStart, cEnd, ej; 1645 1646 PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, &cEnd)); 1647 // prep energy trick, get v^2 / 2 vector 1648 if (ctx->use_energy_tensor_trick) { 1649 PetscErrorCode (*energyf[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *) = {ctx->use_relativistic_corrections ? gamma_m1_f : energy_f}; 1650 Vec glob_v2; 1651 PetscReal *c2_0[1], data[1] = {PetscSqr(C_0(ctx->v_0))}; 1652 1653 PetscCall(DMClone(ctx->plex[grid], &dmEnergy)); 1654 PetscCall(PetscObjectSetName((PetscObject)dmEnergy, "energy")); 1655 PetscCall(DMSetField(dmEnergy, 0, NULL, (PetscObject)fe)); 1656 PetscCall(DMCreateDS(dmEnergy)); 1657 PetscCall(DMGetSection(dmEnergy, &e_section)); 1658 PetscCall(DMGetGlobalVector(dmEnergy, &glob_v2)); 1659 PetscCall(PetscObjectSetName((PetscObject)glob_v2, "trick")); 1660 c2_0[0] = &data[0]; 1661 PetscCall(DMProjectFunction(dmEnergy, 0., energyf, (void **)c2_0, INSERT_ALL_VALUES, glob_v2)); 1662 PetscCall(DMGetLocalVector(dmEnergy, &v2_2)); 1663 PetscCall(VecZeroEntries(v2_2)); /* zero BCs so don't set */ 1664 PetscCall(DMGlobalToLocalBegin(dmEnergy, glob_v2, INSERT_VALUES, v2_2)); 1665 PetscCall(DMGlobalToLocalEnd(dmEnergy, glob_v2, INSERT_VALUES, v2_2)); 1666 PetscCall(DMViewFromOptions(dmEnergy, NULL, "-energy_dm_view")); 1667 PetscCall(VecViewFromOptions(glob_v2, NULL, "-energy_vec_view")); 1668 PetscCall(DMRestoreGlobalVector(dmEnergy, &glob_v2)); 1669 } 1670 /* append part of the IP data for each grid */ 1671 for (ej = 0; ej < numCells[grid]; ++ej, ++outer_ej) { 1672 PetscScalar *coefs = NULL; 1673 PetscReal vj[LANDAU_MAX_NQ * LANDAU_DIM], detJj[LANDAU_MAX_NQ], Jdummy[LANDAU_MAX_NQ * LANDAU_DIM * LANDAU_DIM], c0 = C_0(ctx->v_0), c02 = PetscSqr(c0); 1674 invJe = invJ_a + outer_ej * Nq * dim * dim; 1675 PetscCall(DMPlexComputeCellGeometryFEM(ctx->plex[grid], ej + cStart, quad, vj, Jdummy, invJe, detJj)); 1676 if (ctx->use_energy_tensor_trick) PetscCall(DMPlexVecGetClosure(dmEnergy, e_section, v2_2, ej + cStart, NULL, &coefs)); 1677 /* create static point data */ 1678 for (PetscInt qj = 0; qj < Nq; qj++, outer_ipidx++) { 1679 const PetscInt gidx = outer_ipidx; 1680 const PetscReal *invJ = &invJe[qj * dim * dim]; 1681 ww[gidx] = detJj[qj] * quadWeights[qj]; 1682 if (dim == 2) ww[gidx] *= vj[qj * dim + 0]; /* cylindrical coordinate, w/o 2pi */ 1683 // get xx, yy, zz 1684 if (ctx->use_energy_tensor_trick) { 1685 double refSpaceDer[3], eGradPhi[3]; 1686 const PetscReal *const DD = Tf[0]->T[1]; 1687 const PetscReal *Dq = &DD[qj * Nb * dim]; 1688 for (int d = 0; d < 3; ++d) refSpaceDer[d] = eGradPhi[d] = 0.0; 1689 for (int b = 0; b < Nb; ++b) { 1690 for (int d = 0; d < dim; ++d) refSpaceDer[d] += Dq[b * dim + d] * PetscRealPart(coefs[b]); 1691 } 1692 xx[gidx] = 1e10; 1693 if (ctx->use_relativistic_corrections) { 1694 double dg2_c2 = 0; 1695 //for (int d = 0; d < dim; ++d) refSpaceDer[d] *= c02; 1696 for (int d = 0; d < dim; ++d) dg2_c2 += PetscSqr(refSpaceDer[d]); 1697 dg2_c2 *= (double)c02; 1698 if (dg2_c2 >= .999) { 1699 xx[gidx] = vj[qj * dim + 0]; /* coordinate */ 1700 yy[gidx] = vj[qj * dim + 1]; 1701 if (dim == 3) zz[gidx] = vj[qj * dim + 2]; 1702 PetscCall(PetscPrintf(ctx->comm, "Error: %12.5e %" PetscInt_FMT ".%" PetscInt_FMT ") dg2/c02 = %12.5e x= %12.5e %12.5e %12.5e\n", (double)PetscSqrtReal(xx[gidx] * xx[gidx] + yy[gidx] * yy[gidx] + zz[gidx] * zz[gidx]), ej, qj, dg2_c2, (double)xx[gidx], (double)yy[gidx], (double)zz[gidx])); 1703 } else { 1704 PetscReal fact = c02 / PetscSqrtReal(1. - dg2_c2); 1705 for (int d = 0; d < dim; ++d) refSpaceDer[d] *= fact; 1706 // could test with other point u' that (grad - grad') * U (refSpaceDer, refSpaceDer') == 0 1707 } 1708 } 1709 if (xx[gidx] == 1e10) { 1710 for (int d = 0; d < dim; ++d) { 1711 for (int e = 0; e < dim; ++e) eGradPhi[d] += invJ[e * dim + d] * refSpaceDer[e]; 1712 } 1713 xx[gidx] = eGradPhi[0]; 1714 yy[gidx] = eGradPhi[1]; 1715 if (dim == 3) zz[gidx] = eGradPhi[2]; 1716 } 1717 } else { 1718 xx[gidx] = vj[qj * dim + 0]; /* coordinate */ 1719 yy[gidx] = vj[qj * dim + 1]; 1720 if (dim == 3) zz[gidx] = vj[qj * dim + 2]; 1721 } 1722 } /* q */ 1723 if (ctx->use_energy_tensor_trick) PetscCall(DMPlexVecRestoreClosure(dmEnergy, e_section, v2_2, ej + cStart, NULL, &coefs)); 1724 } /* ej */ 1725 if (ctx->use_energy_tensor_trick) { 1726 PetscCall(DMRestoreLocalVector(dmEnergy, &v2_2)); 1727 PetscCall(DMDestroy(&dmEnergy)); 1728 } 1729 } /* grid */ 1730 if (ctx->use_energy_tensor_trick) PetscCall(PetscFEDestroy(&fe)); 1731 /* cache static data */ 1732 if (ctx->deviceType == LANDAU_CUDA || ctx->deviceType == LANDAU_KOKKOS) { 1733 #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_KOKKOS_KERNELS) 1734 if (ctx->deviceType == LANDAU_CUDA) { 1735 #if defined(PETSC_HAVE_CUDA) 1736 PetscCall(LandauCUDAStaticDataSet(ctx->plex[0], Nq, ctx->batch_sz, ctx->num_grids, numCells, ctx->species_offset, ctx->mat_offset, nu_alpha, nu_beta, invMass, (PetscReal *)ctx->lambdas, invJ_a, xx, yy, zz, ww, &ctx->SData_d)); 1737 #else 1738 SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-landau_device_type cuda not built"); 1739 #endif 1740 } else if (ctx->deviceType == LANDAU_KOKKOS) { 1741 #if defined(PETSC_HAVE_KOKKOS_KERNELS) 1742 PetscCall(LandauKokkosStaticDataSet(ctx->plex[0], Nq, ctx->batch_sz, ctx->num_grids, numCells, ctx->species_offset, ctx->mat_offset, nu_alpha, nu_beta, invMass, (PetscReal *)ctx->lambdas, invJ_a, xx, yy, zz, ww, &ctx->SData_d)); 1743 #else 1744 SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-landau_device_type kokkos not built"); 1745 #endif 1746 } 1747 #endif 1748 /* free */ 1749 PetscCall(PetscFree4(ww, xx, yy, invJ_a)); 1750 if (dim == 3) PetscCall(PetscFree(zz)); 1751 } else { /* CPU version, just copy in, only use part */ 1752 PetscReal *nu_alpha_p = (PetscReal *)ctx->SData_d.alpha, *nu_beta_p = (PetscReal *)ctx->SData_d.beta, *invMass_p = (PetscReal *)ctx->SData_d.invMass, *lambdas_p = NULL; // why set these ? 1753 ctx->SData_d.w = (void *)ww; 1754 ctx->SData_d.x = (void *)xx; 1755 ctx->SData_d.y = (void *)yy; 1756 ctx->SData_d.z = (void *)zz; 1757 ctx->SData_d.invJ = (void *)invJ_a; 1758 PetscCall(PetscMalloc4(ctx->num_species, &nu_alpha_p, ctx->num_species, &nu_beta_p, ctx->num_species, &invMass_p, LANDAU_MAX_GRIDS * LANDAU_MAX_GRIDS, &lambdas_p)); 1759 for (PetscInt ii = 0; ii < ctx->num_species; ii++) { 1760 nu_alpha_p[ii] = nu_alpha[ii]; 1761 nu_beta_p[ii] = nu_beta[ii]; 1762 invMass_p[ii] = invMass[ii]; 1763 } 1764 ctx->SData_d.alpha = (void *)nu_alpha_p; 1765 ctx->SData_d.beta = (void *)nu_beta_p; 1766 ctx->SData_d.invMass = (void *)invMass_p; 1767 ctx->SData_d.lambdas = (void *)lambdas_p; 1768 for (PetscInt grid = 0; grid < LANDAU_MAX_GRIDS; grid++) { 1769 PetscReal(*lambdas)[LANDAU_MAX_GRIDS][LANDAU_MAX_GRIDS] = (PetscReal(*)[LANDAU_MAX_GRIDS][LANDAU_MAX_GRIDS])ctx->SData_d.lambdas; 1770 for (PetscInt gridj = 0; gridj < LANDAU_MAX_GRIDS; gridj++) { (*lambdas)[grid][gridj] = ctx->lambdas[grid][gridj]; } 1771 } 1772 } 1773 PetscCall(PetscLogEventEnd(ctx->events[7], 0, 0, 0, 0)); 1774 } // initialize 1775 PetscFunctionReturn(PETSC_SUCCESS); 1776 } 1777 1778 /* < v, u > */ 1779 static void g0_1(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 1780 { 1781 g0[0] = 1.; 1782 } 1783 1784 /* < v, u > */ 1785 static void g0_fake(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 1786 { 1787 static double ttt = 1e-12; 1788 g0[0] = ttt++; 1789 } 1790 1791 /* < v, u > */ 1792 static void g0_r(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 1793 { 1794 g0[0] = 2. * PETSC_PI * x[0]; 1795 } 1796 1797 static PetscErrorCode MatrixNfDestroy(void *ptr) 1798 { 1799 PetscInt *nf = (PetscInt *)ptr; 1800 PetscFunctionBegin; 1801 PetscCall(PetscFree(nf)); 1802 PetscFunctionReturn(PETSC_SUCCESS); 1803 } 1804 1805 /* 1806 LandauCreateJacobianMatrix - creates ctx->J with without real data. Hard to keep sparse. 1807 - Like DMPlexLandauCreateMassMatrix. Should remove one and combine 1808 - has old support for field major ordering 1809 */ 1810 static PetscErrorCode LandauCreateJacobianMatrix(MPI_Comm comm, Vec X, IS grid_batch_is_inv[LANDAU_MAX_GRIDS], LandauCtx *ctx) 1811 { 1812 PetscInt *idxs = NULL; 1813 Mat subM[LANDAU_MAX_GRIDS]; 1814 1815 PetscFunctionBegin; 1816 if (!ctx->gpu_assembly) { /* we need GPU object with GPU assembly */ 1817 PetscFunctionReturn(PETSC_SUCCESS); 1818 } 1819 // get the RCM for this grid to separate out species into blocks -- create 'idxs' & 'ctx->batch_is' -- not used 1820 if (ctx->gpu_assembly && ctx->jacobian_field_major_order) PetscCall(PetscMalloc1(ctx->mat_offset[ctx->num_grids] * ctx->batch_sz, &idxs)); 1821 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 1822 const PetscInt *values, n = ctx->mat_offset[grid + 1] - ctx->mat_offset[grid]; 1823 Mat gMat; 1824 DM massDM; 1825 PetscDS prob; 1826 Vec tvec; 1827 // get "mass" matrix for reordering 1828 PetscCall(DMClone(ctx->plex[grid], &massDM)); 1829 PetscCall(DMCopyFields(ctx->plex[grid], massDM)); 1830 PetscCall(DMCreateDS(massDM)); 1831 PetscCall(DMGetDS(massDM, &prob)); 1832 for (int ix = 0, ii = ctx->species_offset[grid]; ii < ctx->species_offset[grid + 1]; ii++, ix++) PetscCall(PetscDSSetJacobian(prob, ix, ix, g0_fake, NULL, NULL, NULL)); 1833 PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only")); // this trick is need to both sparsify the matrix and avoid runtime error 1834 PetscCall(DMCreateMatrix(massDM, &gMat)); 1835 PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only false")); 1836 PetscCall(MatSetOption(gMat, MAT_STRUCTURALLY_SYMMETRIC, PETSC_TRUE)); 1837 PetscCall(MatSetOption(gMat, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE)); 1838 PetscCall(DMCreateLocalVector(ctx->plex[grid], &tvec)); 1839 PetscCall(DMPlexSNESComputeJacobianFEM(massDM, tvec, gMat, gMat, ctx)); 1840 PetscCall(MatViewFromOptions(gMat, NULL, "-dm_landau_reorder_mat_view")); 1841 PetscCall(DMDestroy(&massDM)); 1842 PetscCall(VecDestroy(&tvec)); 1843 subM[grid] = gMat; 1844 if (ctx->gpu_assembly && ctx->jacobian_field_major_order) { 1845 MatOrderingType rtype = MATORDERINGRCM; 1846 IS isrow, isicol; 1847 PetscCall(MatGetOrdering(gMat, rtype, &isrow, &isicol)); 1848 PetscCall(ISInvertPermutation(isrow, PETSC_DECIDE, &grid_batch_is_inv[grid])); 1849 PetscCall(ISGetIndices(isrow, &values)); 1850 for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { // add batch size DMs for this species grid 1851 #if !defined(LANDAU_SPECIES_MAJOR) 1852 PetscInt N = ctx->mat_offset[ctx->num_grids], n0 = ctx->mat_offset[grid] + b_id * N; 1853 for (int ii = 0; ii < n; ++ii) idxs[n0 + ii] = values[ii] + n0; 1854 #else 1855 PetscInt n0 = ctx->mat_offset[grid] * ctx->batch_sz + b_id * n; 1856 for (int ii = 0; ii < n; ++ii) idxs[n0 + ii] = values[ii] + n0; 1857 #endif 1858 } 1859 PetscCall(ISRestoreIndices(isrow, &values)); 1860 PetscCall(ISDestroy(&isrow)); 1861 PetscCall(ISDestroy(&isicol)); 1862 } 1863 } 1864 if (ctx->gpu_assembly && ctx->jacobian_field_major_order) PetscCall(ISCreateGeneral(comm, ctx->mat_offset[ctx->num_grids] * ctx->batch_sz, idxs, PETSC_OWN_POINTER, &ctx->batch_is)); 1865 // get a block matrix 1866 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 1867 Mat B = subM[grid]; 1868 PetscInt nloc, nzl, *colbuf, row, COL_BF_SIZE = 1024; 1869 PetscCall(PetscMalloc(sizeof(*colbuf) * COL_BF_SIZE, &colbuf)); 1870 PetscCall(MatGetSize(B, &nloc, NULL)); 1871 for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { 1872 const PetscInt moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset); 1873 const PetscInt *cols; 1874 const PetscScalar *vals; 1875 for (int i = 0; i < nloc; i++) { 1876 PetscCall(MatGetRow(B, i, &nzl, NULL, NULL)); 1877 if (nzl > COL_BF_SIZE) { 1878 PetscCall(PetscFree(colbuf)); 1879 PetscCall(PetscInfo(ctx->plex[grid], "Realloc buffer %" PetscInt_FMT " to %" PetscInt_FMT " (row size %" PetscInt_FMT ") \n", COL_BF_SIZE, 2 * COL_BF_SIZE, nzl)); 1880 COL_BF_SIZE = nzl; 1881 PetscCall(PetscMalloc(sizeof(*colbuf) * COL_BF_SIZE, &colbuf)); 1882 } 1883 PetscCall(MatGetRow(B, i, &nzl, &cols, &vals)); 1884 for (int j = 0; j < nzl; j++) colbuf[j] = cols[j] + moffset; 1885 row = i + moffset; 1886 PetscCall(MatSetValues(ctx->J, 1, &row, nzl, colbuf, vals, INSERT_VALUES)); 1887 PetscCall(MatRestoreRow(B, i, &nzl, &cols, &vals)); 1888 } 1889 } 1890 PetscCall(PetscFree(colbuf)); 1891 } 1892 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(MatDestroy(&subM[grid])); 1893 PetscCall(MatAssemblyBegin(ctx->J, MAT_FINAL_ASSEMBLY)); 1894 PetscCall(MatAssemblyEnd(ctx->J, MAT_FINAL_ASSEMBLY)); 1895 1896 // debug 1897 PetscCall(MatViewFromOptions(ctx->J, NULL, "-dm_landau_mat_view")); 1898 if (ctx->gpu_assembly && ctx->jacobian_field_major_order) { 1899 Mat mat_block_order; 1900 PetscCall(MatCreateSubMatrix(ctx->J, ctx->batch_is, ctx->batch_is, MAT_INITIAL_MATRIX, &mat_block_order)); // use MatPermute 1901 PetscCall(MatViewFromOptions(mat_block_order, NULL, "-dm_landau_mat_view")); 1902 PetscCall(MatDestroy(&mat_block_order)); 1903 PetscCall(VecScatterCreate(X, ctx->batch_is, X, NULL, &ctx->plex_batch)); 1904 PetscCall(VecDuplicate(X, &ctx->work_vec)); 1905 } 1906 1907 PetscFunctionReturn(PETSC_SUCCESS); 1908 } 1909 1910 PetscErrorCode DMPlexLandauCreateMassMatrix(DM pack, Mat *Amat); 1911 /*@C 1912 DMPlexLandauCreateVelocitySpace - Create a DMPlex velocity space mesh 1913 1914 Collective 1915 1916 Input Parameters: 1917 + comm - The MPI communicator 1918 . dim - velocity space dimension (2 for axisymmetric, 3 for full 3X + 3V solver) 1919 - prefix - prefix for options (not tested) 1920 1921 Output Parameter: 1922 . pack - The DM object representing the mesh 1923 + X - A vector (user destroys) 1924 - J - Optional matrix (object destroys) 1925 1926 Level: beginner 1927 1928 .keywords: mesh 1929 .seealso: `DMPlexCreate()`, `DMPlexLandauDestroyVelocitySpace()` 1930 @*/ 1931 PetscErrorCode DMPlexLandauCreateVelocitySpace(MPI_Comm comm, PetscInt dim, const char prefix[], Vec *X, Mat *J, DM *pack) 1932 { 1933 LandauCtx *ctx; 1934 Vec Xsub[LANDAU_MAX_GRIDS]; 1935 IS grid_batch_is_inv[LANDAU_MAX_GRIDS]; 1936 1937 PetscFunctionBegin; 1938 PetscCheck(dim == 2 || dim == 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Only 2D and 3D supported"); 1939 PetscCheck(LANDAU_DIM == dim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "dim %" PetscInt_FMT " != LANDAU_DIM %d", dim, LANDAU_DIM); 1940 PetscCall(PetscNew(&ctx)); 1941 ctx->comm = comm; /* used for diagnostics and global errors */ 1942 /* process options */ 1943 PetscCall(ProcessOptions(ctx, prefix)); 1944 if (dim == 2) ctx->use_relativistic_corrections = PETSC_FALSE; 1945 /* Create Mesh */ 1946 PetscCall(DMCompositeCreate(PETSC_COMM_SELF, pack)); 1947 PetscCall(PetscLogEventBegin(ctx->events[13], 0, 0, 0, 0)); 1948 PetscCall(PetscLogEventBegin(ctx->events[15], 0, 0, 0, 0)); 1949 PetscCall(LandauDMCreateVMeshes(PETSC_COMM_SELF, dim, prefix, ctx, *pack)); // creates grids (Forest of AMR) 1950 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 1951 /* create FEM */ 1952 PetscCall(SetupDS(ctx->plex[grid], dim, grid, ctx)); 1953 /* set initial state */ 1954 PetscCall(DMCreateGlobalVector(ctx->plex[grid], &Xsub[grid])); 1955 PetscCall(PetscObjectSetName((PetscObject)Xsub[grid], "u_orig")); 1956 /* initial static refinement, no solve */ 1957 PetscCall(LandauSetInitialCondition(ctx->plex[grid], Xsub[grid], grid, 0, 1, ctx)); 1958 /* forest refinement - forest goes in (if forest), plex comes out */ 1959 if (ctx->use_p4est) { 1960 DM plex; 1961 PetscCall(adapt(grid, ctx, &Xsub[grid])); // forest goes in, plex comes out 1962 PetscCall(DMViewFromOptions(ctx->plex[grid], NULL, "-dm_landau_amr_dm_view")); // need to differentiate - todo 1963 PetscCall(VecViewFromOptions(Xsub[grid], NULL, "-dm_landau_amr_vec_view")); 1964 // convert to plex, all done with this level 1965 PetscCall(DMConvert(ctx->plex[grid], DMPLEX, &plex)); 1966 PetscCall(DMDestroy(&ctx->plex[grid])); 1967 ctx->plex[grid] = plex; 1968 } 1969 #if !defined(LANDAU_SPECIES_MAJOR) 1970 PetscCall(DMCompositeAddDM(*pack, ctx->plex[grid])); 1971 #else 1972 for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { // add batch size DMs for this species grid 1973 PetscCall(DMCompositeAddDM(*pack, ctx->plex[grid])); 1974 } 1975 #endif 1976 PetscCall(DMSetApplicationContext(ctx->plex[grid], ctx)); 1977 } 1978 #if !defined(LANDAU_SPECIES_MAJOR) 1979 // stack the batched DMs, could do it all here!!! b_id=0 1980 for (PetscInt b_id = 1; b_id < ctx->batch_sz; b_id++) { 1981 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(DMCompositeAddDM(*pack, ctx->plex[grid])); 1982 } 1983 #endif 1984 // create ctx->mat_offset 1985 ctx->mat_offset[0] = 0; 1986 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 1987 PetscInt n; 1988 PetscCall(VecGetLocalSize(Xsub[grid], &n)); 1989 ctx->mat_offset[grid + 1] = ctx->mat_offset[grid] + n; 1990 } 1991 // creat DM & Jac 1992 PetscCall(DMSetApplicationContext(*pack, ctx)); 1993 PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only")); 1994 PetscCall(DMCreateMatrix(*pack, &ctx->J)); 1995 PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only false")); 1996 PetscCall(MatSetOption(ctx->J, MAT_STRUCTURALLY_SYMMETRIC, PETSC_TRUE)); 1997 PetscCall(MatSetOption(ctx->J, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE)); 1998 PetscCall(PetscObjectSetName((PetscObject)ctx->J, "Jac")); 1999 // construct initial conditions in X 2000 PetscCall(DMCreateGlobalVector(*pack, X)); 2001 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 2002 PetscInt n; 2003 PetscCall(VecGetLocalSize(Xsub[grid], &n)); 2004 for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { 2005 PetscScalar const *values; 2006 const PetscInt moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset); 2007 PetscCall(LandauSetInitialCondition(ctx->plex[grid], Xsub[grid], grid, b_id, ctx->batch_sz, ctx)); 2008 PetscCall(VecGetArrayRead(Xsub[grid], &values)); // Drop whole grid in Plex ordering 2009 for (int i = 0, idx = moffset; i < n; i++, idx++) PetscCall(VecSetValue(*X, idx, values[i], INSERT_VALUES)); 2010 PetscCall(VecRestoreArrayRead(Xsub[grid], &values)); 2011 } 2012 } 2013 // cleanup 2014 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(VecDestroy(&Xsub[grid])); 2015 /* check for correct matrix type */ 2016 if (ctx->gpu_assembly) { /* we need GPU object with GPU assembly */ 2017 PetscBool flg; 2018 if (ctx->deviceType == LANDAU_CUDA) { 2019 PetscCall(PetscObjectTypeCompareAny((PetscObject)ctx->J, &flg, MATSEQAIJCUSPARSE, MATMPIAIJCUSPARSE, MATAIJCUSPARSE, "")); 2020 PetscCheck(flg, ctx->comm, PETSC_ERR_ARG_WRONG, "must use '-dm_mat_type aijcusparse -dm_vec_type cuda' for GPU assembly and Cuda or use '-dm_landau_device_type cpu'"); 2021 } else if (ctx->deviceType == LANDAU_KOKKOS) { 2022 PetscCall(PetscObjectTypeCompareAny((PetscObject)ctx->J, &flg, MATSEQAIJKOKKOS, MATMPIAIJKOKKOS, MATAIJKOKKOS, "")); 2023 #if defined(PETSC_HAVE_KOKKOS_KERNELS) 2024 PetscCheck(flg, ctx->comm, PETSC_ERR_ARG_WRONG, "must use '-dm_mat_type aijkokkos -dm_vec_type kokkos' for GPU assembly and Kokkos or use '-dm_landau_device_type cpu'"); 2025 #else 2026 PetscCheck(flg, ctx->comm, PETSC_ERR_ARG_WRONG, "must configure with '--download-kokkos-kernels' for GPU assembly and Kokkos or use '-dm_landau_device_type cpu'"); 2027 #endif 2028 } 2029 } 2030 PetscCall(PetscLogEventEnd(ctx->events[15], 0, 0, 0, 0)); 2031 2032 // create field major ordering 2033 ctx->work_vec = NULL; 2034 ctx->plex_batch = NULL; 2035 ctx->batch_is = NULL; 2036 for (int i = 0; i < LANDAU_MAX_GRIDS; i++) grid_batch_is_inv[i] = NULL; 2037 PetscCall(PetscLogEventBegin(ctx->events[12], 0, 0, 0, 0)); 2038 PetscCall(LandauCreateJacobianMatrix(comm, *X, grid_batch_is_inv, ctx)); 2039 PetscCall(PetscLogEventEnd(ctx->events[12], 0, 0, 0, 0)); 2040 2041 // create AMR GPU assembly maps and static GPU data 2042 PetscCall(CreateStaticData(dim, grid_batch_is_inv, ctx)); 2043 2044 PetscCall(PetscLogEventEnd(ctx->events[13], 0, 0, 0, 0)); 2045 2046 // create mass matrix 2047 PetscCall(DMPlexLandauCreateMassMatrix(*pack, NULL)); 2048 2049 if (J) *J = ctx->J; 2050 2051 if (ctx->gpu_assembly && ctx->jacobian_field_major_order) { 2052 PetscContainer container; 2053 // cache ctx for KSP with batch/field major Jacobian ordering -ksp_type gmres/etc -dm_landau_jacobian_field_major_order 2054 PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container)); 2055 PetscCall(PetscContainerSetPointer(container, (void *)ctx)); 2056 PetscCall(PetscObjectCompose((PetscObject)ctx->J, "LandauCtx", (PetscObject)container)); 2057 PetscCall(PetscContainerDestroy(&container)); 2058 // batch solvers need to map -- can batch solvers work 2059 PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container)); 2060 PetscCall(PetscContainerSetPointer(container, (void *)ctx->plex_batch)); 2061 PetscCall(PetscObjectCompose((PetscObject)ctx->J, "plex_batch_is", (PetscObject)container)); 2062 PetscCall(PetscContainerDestroy(&container)); 2063 } 2064 // for batch solvers 2065 { 2066 PetscContainer container; 2067 PetscInt *pNf; 2068 PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container)); 2069 PetscCall(PetscMalloc1(sizeof(*pNf), &pNf)); 2070 *pNf = ctx->batch_sz; 2071 PetscCall(PetscContainerSetPointer(container, (void *)pNf)); 2072 PetscCall(PetscContainerSetUserDestroy(container, MatrixNfDestroy)); 2073 PetscCall(PetscObjectCompose((PetscObject)ctx->J, "batch size", (PetscObject)container)); 2074 PetscCall(PetscContainerDestroy(&container)); 2075 } 2076 2077 PetscFunctionReturn(PETSC_SUCCESS); 2078 } 2079 2080 /*@ 2081 DMPlexLandauAccess - Access to the distribution function with user callback 2082 2083 Collective 2084 2085 Input Parameters: 2086 . pack - the DMComposite 2087 + func - call back function 2088 . user_ctx - user context 2089 2090 Input/Output Parameter: 2091 . X - Vector to data to 2092 2093 Level: advanced 2094 2095 .keywords: mesh 2096 .seealso: `DMPlexLandauCreateVelocitySpace()` 2097 @*/ 2098 PetscErrorCode DMPlexLandauAccess(DM pack, Vec X, PetscErrorCode (*func)(DM, Vec, PetscInt, PetscInt, PetscInt, void *), void *user_ctx) 2099 { 2100 LandauCtx *ctx; 2101 PetscFunctionBegin; 2102 PetscCall(DMGetApplicationContext(pack, &ctx)); // uses ctx->num_grids; ctx->plex[grid]; ctx->batch_sz; ctx->mat_offset 2103 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 2104 PetscInt dim, n; 2105 PetscCall(DMGetDimension(pack, &dim)); 2106 for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) { 2107 Vec vec; 2108 PetscInt vf[1] = {i0}; 2109 IS vis; 2110 DM vdm; 2111 PetscCall(DMCreateSubDM(ctx->plex[grid], 1, vf, &vis, &vdm)); 2112 PetscCall(DMSetApplicationContext(vdm, ctx)); // the user might want this 2113 PetscCall(DMCreateGlobalVector(vdm, &vec)); 2114 PetscCall(VecGetSize(vec, &n)); 2115 for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { 2116 const PetscInt moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset); 2117 PetscCall(VecZeroEntries(vec)); 2118 /* Add your data with 'dm' for species 'sp' to 'vec' */ 2119 PetscCall(func(vdm, vec, i0, grid, b_id, user_ctx)); 2120 /* add to global */ 2121 PetscScalar const *values; 2122 const PetscInt *offsets; 2123 PetscCall(VecGetArrayRead(vec, &values)); 2124 PetscCall(ISGetIndices(vis, &offsets)); 2125 for (int i = 0; i < n; i++) PetscCall(VecSetValue(X, moffset + offsets[i], values[i], ADD_VALUES)); 2126 PetscCall(VecRestoreArrayRead(vec, &values)); 2127 PetscCall(ISRestoreIndices(vis, &offsets)); 2128 } // batch 2129 PetscCall(VecDestroy(&vec)); 2130 PetscCall(ISDestroy(&vis)); 2131 PetscCall(DMDestroy(&vdm)); 2132 } 2133 } // grid 2134 PetscFunctionReturn(PETSC_SUCCESS); 2135 } 2136 2137 /*@ 2138 DMPlexLandauDestroyVelocitySpace - Destroy a DMPlex velocity space mesh 2139 2140 Collective 2141 2142 Input/Output Parameters: 2143 . dm - the dm to destroy 2144 2145 Level: beginner 2146 2147 .keywords: mesh 2148 .seealso: `DMPlexLandauCreateVelocitySpace()` 2149 @*/ 2150 PetscErrorCode DMPlexLandauDestroyVelocitySpace(DM *dm) 2151 { 2152 LandauCtx *ctx; 2153 PetscFunctionBegin; 2154 PetscCall(DMGetApplicationContext(*dm, &ctx)); 2155 PetscCall(MatDestroy(&ctx->M)); 2156 PetscCall(MatDestroy(&ctx->J)); 2157 for (PetscInt ii = 0; ii < ctx->num_species; ii++) PetscCall(PetscFEDestroy(&ctx->fe[ii])); 2158 PetscCall(ISDestroy(&ctx->batch_is)); 2159 PetscCall(VecDestroy(&ctx->work_vec)); 2160 PetscCall(VecScatterDestroy(&ctx->plex_batch)); 2161 if (ctx->deviceType == LANDAU_CUDA) { 2162 #if defined(PETSC_HAVE_CUDA) 2163 PetscCall(LandauCUDAStaticDataClear(&ctx->SData_d)); 2164 #else 2165 SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-landau_device_type %s not built", "cuda"); 2166 #endif 2167 } else if (ctx->deviceType == LANDAU_KOKKOS) { 2168 #if defined(PETSC_HAVE_KOKKOS_KERNELS) 2169 PetscCall(LandauKokkosStaticDataClear(&ctx->SData_d)); 2170 #else 2171 SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-landau_device_type %s not built", "kokkos"); 2172 #endif 2173 } else { 2174 if (ctx->SData_d.x) { /* in a CPU run */ 2175 PetscReal *invJ = (PetscReal *)ctx->SData_d.invJ, *xx = (PetscReal *)ctx->SData_d.x, *yy = (PetscReal *)ctx->SData_d.y, *zz = (PetscReal *)ctx->SData_d.z, *ww = (PetscReal *)ctx->SData_d.w; 2176 LandauIdx *coo_elem_offsets = (LandauIdx *)ctx->SData_d.coo_elem_offsets, *coo_elem_fullNb = (LandauIdx *)ctx->SData_d.coo_elem_fullNb, (*coo_elem_point_offsets)[LANDAU_MAX_NQ + 1] = (LandauIdx(*)[LANDAU_MAX_NQ + 1]) ctx->SData_d.coo_elem_point_offsets; 2177 PetscCall(PetscFree4(ww, xx, yy, invJ)); 2178 if (zz) PetscCall(PetscFree(zz)); 2179 if (coo_elem_offsets) { 2180 PetscCall(PetscFree3(coo_elem_offsets, coo_elem_fullNb, coo_elem_point_offsets)); // could be NULL 2181 } 2182 PetscCall(PetscFree4(ctx->SData_d.alpha, ctx->SData_d.beta, ctx->SData_d.invMass, ctx->SData_d.lambdas)); 2183 } 2184 } 2185 2186 if (ctx->times[LANDAU_MATRIX_TOTAL] > 0) { // OMP timings 2187 PetscCall(PetscPrintf(ctx->comm, "TSStep N 1.0 %10.3e\n", ctx->times[LANDAU_EX2_TSSOLVE])); 2188 PetscCall(PetscPrintf(ctx->comm, "2: Solve: %10.3e with %" PetscInt_FMT " threads\n", ctx->times[LANDAU_EX2_TSSOLVE] - ctx->times[LANDAU_MATRIX_TOTAL], ctx->batch_sz)); 2189 PetscCall(PetscPrintf(ctx->comm, "3: Landau: %10.3e\n", ctx->times[LANDAU_MATRIX_TOTAL])); 2190 PetscCall(PetscPrintf(ctx->comm, "Landau Jacobian %" PetscInt_FMT " 1.0 %10.3e\n", (PetscInt)ctx->times[LANDAU_JACOBIAN_COUNT], ctx->times[LANDAU_JACOBIAN])); 2191 PetscCall(PetscPrintf(ctx->comm, "Landau Operator N 1.0 %10.3e\n", ctx->times[LANDAU_OPERATOR])); 2192 PetscCall(PetscPrintf(ctx->comm, "Landau Mass N 1.0 %10.3e\n", ctx->times[LANDAU_MASS])); 2193 PetscCall(PetscPrintf(ctx->comm, " Jac-f-df (GPU) N 1.0 %10.3e\n", ctx->times[LANDAU_F_DF])); 2194 PetscCall(PetscPrintf(ctx->comm, " Kernel (GPU) N 1.0 %10.3e\n", ctx->times[LANDAU_KERNEL])); 2195 PetscCall(PetscPrintf(ctx->comm, "MatLUFactorNum X 1.0 %10.3e\n", ctx->times[KSP_FACTOR])); 2196 PetscCall(PetscPrintf(ctx->comm, "MatSolve X 1.0 %10.3e\n", ctx->times[KSP_SOLVE])); 2197 } 2198 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(DMDestroy(&ctx->plex[grid])); 2199 PetscCall(PetscFree(ctx)); 2200 PetscCall(DMDestroy(dm)); 2201 PetscFunctionReturn(PETSC_SUCCESS); 2202 } 2203 2204 /* < v, ru > */ 2205 static void f0_s_den(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 2206 { 2207 PetscInt ii = (PetscInt)PetscRealPart(constants[0]); 2208 f0[0] = u[ii]; 2209 } 2210 2211 /* < v, ru > */ 2212 static void f0_s_mom(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 2213 { 2214 PetscInt ii = (PetscInt)PetscRealPart(constants[0]), jj = (PetscInt)PetscRealPart(constants[1]); 2215 f0[0] = x[jj] * u[ii]; /* x momentum */ 2216 } 2217 2218 static void f0_s_v2(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 2219 { 2220 PetscInt i, ii = (PetscInt)PetscRealPart(constants[0]); 2221 double tmp1 = 0.; 2222 for (i = 0; i < dim; ++i) tmp1 += x[i] * x[i]; 2223 f0[0] = tmp1 * u[ii]; 2224 } 2225 2226 static PetscErrorCode gamma_n_f(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *actx) 2227 { 2228 const PetscReal *c2_0_arr = ((PetscReal *)actx); 2229 const PetscReal c02 = c2_0_arr[0]; 2230 2231 PetscFunctionBegin; 2232 for (int s = 0; s < Nf; s++) { 2233 PetscReal tmp1 = 0.; 2234 for (int i = 0; i < dim; ++i) tmp1 += x[i] * x[i]; 2235 #if defined(PETSC_USE_DEBUG) 2236 u[s] = PetscSqrtReal(1. + tmp1 / c02); // u[0] = PetscSqrtReal(1. + xx); 2237 #else 2238 { 2239 PetscReal xx = tmp1 / c02; 2240 u[s] = xx / (PetscSqrtReal(1. + xx) + 1.); // better conditioned = xx/(PetscSqrtReal(1. + xx) + 1.) 2241 } 2242 #endif 2243 } 2244 PetscFunctionReturn(PETSC_SUCCESS); 2245 } 2246 2247 /* < v, ru > */ 2248 static void f0_s_rden(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 2249 { 2250 PetscInt ii = (PetscInt)PetscRealPart(constants[0]); 2251 f0[0] = 2. * PETSC_PI * x[0] * u[ii]; 2252 } 2253 2254 /* < v, ru > */ 2255 static void f0_s_rmom(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 2256 { 2257 PetscInt ii = (PetscInt)PetscRealPart(constants[0]); 2258 f0[0] = 2. * PETSC_PI * x[0] * x[1] * u[ii]; 2259 } 2260 2261 static void f0_s_rv2(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 2262 { 2263 PetscInt ii = (PetscInt)PetscRealPart(constants[0]); 2264 f0[0] = 2. * PETSC_PI * x[0] * (x[0] * x[0] + x[1] * x[1]) * u[ii]; 2265 } 2266 2267 /*@ 2268 DMPlexLandauPrintNorms - collects moments and prints them 2269 2270 Collective 2271 2272 Input Parameters: 2273 + X - the state 2274 - stepi - current step to print 2275 2276 Level: beginner 2277 2278 .keywords: mesh 2279 .seealso: `DMPlexLandauCreateVelocitySpace()` 2280 @*/ 2281 PetscErrorCode DMPlexLandauPrintNorms(Vec X, PetscInt stepi) 2282 { 2283 LandauCtx *ctx; 2284 PetscDS prob; 2285 DM pack; 2286 PetscInt cStart, cEnd, dim, ii, i0, nDMs; 2287 PetscScalar xmomentumtot = 0, ymomentumtot = 0, zmomentumtot = 0, energytot = 0, densitytot = 0, tt[LANDAU_MAX_SPECIES]; 2288 PetscScalar xmomentum[LANDAU_MAX_SPECIES], ymomentum[LANDAU_MAX_SPECIES], zmomentum[LANDAU_MAX_SPECIES], energy[LANDAU_MAX_SPECIES], density[LANDAU_MAX_SPECIES]; 2289 Vec *globXArray; 2290 2291 PetscFunctionBegin; 2292 PetscCall(VecGetDM(X, &pack)); 2293 PetscCheck(pack, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Vector has no DM"); 2294 PetscCall(DMGetDimension(pack, &dim)); 2295 PetscCheck(dim == 2 || dim == 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "dim %" PetscInt_FMT " not in [2,3]", dim); 2296 PetscCall(DMGetApplicationContext(pack, &ctx)); 2297 PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context"); 2298 /* print momentum and energy */ 2299 PetscCall(DMCompositeGetNumberDM(pack, &nDMs)); 2300 PetscCheck(nDMs == ctx->num_grids * ctx->batch_sz, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "#DM wrong %" PetscInt_FMT " %" PetscInt_FMT, nDMs, ctx->num_grids * ctx->batch_sz); 2301 PetscCall(PetscMalloc(sizeof(*globXArray) * nDMs, &globXArray)); 2302 PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray)); 2303 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 2304 Vec Xloc = globXArray[LAND_PACK_IDX(ctx->batch_view_idx, grid)]; 2305 PetscCall(DMGetDS(ctx->plex[grid], &prob)); 2306 for (ii = ctx->species_offset[grid], i0 = 0; ii < ctx->species_offset[grid + 1]; ii++, i0++) { 2307 PetscScalar user[2] = {(PetscScalar)i0, (PetscScalar)ctx->charges[ii]}; 2308 PetscCall(PetscDSSetConstants(prob, 2, user)); 2309 if (dim == 2) { /* 2/3X + 3V (cylindrical coordinates) */ 2310 PetscCall(PetscDSSetObjective(prob, 0, &f0_s_rden)); 2311 PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx)); 2312 density[ii] = tt[0] * ctx->n_0 * ctx->charges[ii]; 2313 PetscCall(PetscDSSetObjective(prob, 0, &f0_s_rmom)); 2314 PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx)); 2315 zmomentum[ii] = tt[0] * ctx->n_0 * ctx->v_0 * ctx->masses[ii]; 2316 PetscCall(PetscDSSetObjective(prob, 0, &f0_s_rv2)); 2317 PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx)); 2318 energy[ii] = tt[0] * 0.5 * ctx->n_0 * ctx->v_0 * ctx->v_0 * ctx->masses[ii]; 2319 zmomentumtot += zmomentum[ii]; 2320 energytot += energy[ii]; 2321 densitytot += density[ii]; 2322 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%3" PetscInt_FMT ") species-%" PetscInt_FMT ": charge density= %20.13e z-momentum= %20.13e energy= %20.13e", stepi, ii, (double)PetscRealPart(density[ii]), (double)PetscRealPart(zmomentum[ii]), (double)PetscRealPart(energy[ii]))); 2323 } else { /* 2/3Xloc + 3V */ 2324 PetscCall(PetscDSSetObjective(prob, 0, &f0_s_den)); 2325 PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx)); 2326 density[ii] = tt[0] * ctx->n_0 * ctx->charges[ii]; 2327 PetscCall(PetscDSSetObjective(prob, 0, &f0_s_mom)); 2328 user[1] = 0; 2329 PetscCall(PetscDSSetConstants(prob, 2, user)); 2330 PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx)); 2331 xmomentum[ii] = tt[0] * ctx->n_0 * ctx->v_0 * ctx->masses[ii]; 2332 user[1] = 1; 2333 PetscCall(PetscDSSetConstants(prob, 2, user)); 2334 PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx)); 2335 ymomentum[ii] = tt[0] * ctx->n_0 * ctx->v_0 * ctx->masses[ii]; 2336 user[1] = 2; 2337 PetscCall(PetscDSSetConstants(prob, 2, user)); 2338 PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx)); 2339 zmomentum[ii] = tt[0] * ctx->n_0 * ctx->v_0 * ctx->masses[ii]; 2340 if (ctx->use_relativistic_corrections) { 2341 /* gamma * M * f */ 2342 if (ii == 0 && grid == 0) { // do all at once 2343 Vec Mf, globGamma, *globMfArray, *globGammaArray; 2344 PetscErrorCode (*gammaf[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *) = {gamma_n_f}; 2345 PetscReal *c2_0[1], data[1]; 2346 2347 PetscCall(VecDuplicate(X, &globGamma)); 2348 PetscCall(VecDuplicate(X, &Mf)); 2349 PetscCall(PetscMalloc(sizeof(*globMfArray) * nDMs, &globMfArray)); 2350 PetscCall(PetscMalloc(sizeof(*globMfArray) * nDMs, &globGammaArray)); 2351 /* M * f */ 2352 PetscCall(MatMult(ctx->M, X, Mf)); 2353 /* gamma */ 2354 PetscCall(DMCompositeGetAccessArray(pack, globGamma, nDMs, NULL, globGammaArray)); 2355 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // yes a grid loop in a grid loop to print nice, need to fix for batching 2356 Vec v1 = globGammaArray[LAND_PACK_IDX(ctx->batch_view_idx, grid)]; 2357 data[0] = PetscSqr(C_0(ctx->v_0)); 2358 c2_0[0] = &data[0]; 2359 PetscCall(DMProjectFunction(ctx->plex[grid], 0., gammaf, (void **)c2_0, INSERT_ALL_VALUES, v1)); 2360 } 2361 PetscCall(DMCompositeRestoreAccessArray(pack, globGamma, nDMs, NULL, globGammaArray)); 2362 /* gamma * Mf */ 2363 PetscCall(DMCompositeGetAccessArray(pack, globGamma, nDMs, NULL, globGammaArray)); 2364 PetscCall(DMCompositeGetAccessArray(pack, Mf, nDMs, NULL, globMfArray)); 2365 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // yes a grid loop in a grid loop to print nice 2366 PetscInt Nf = ctx->species_offset[grid + 1] - ctx->species_offset[grid], N, bs; 2367 Vec Mfsub = globMfArray[LAND_PACK_IDX(ctx->batch_view_idx, grid)], Gsub = globGammaArray[LAND_PACK_IDX(ctx->batch_view_idx, grid)], v1, v2; 2368 // get each component 2369 PetscCall(VecGetSize(Mfsub, &N)); 2370 PetscCall(VecCreate(ctx->comm, &v1)); 2371 PetscCall(VecSetSizes(v1, PETSC_DECIDE, N / Nf)); 2372 PetscCall(VecCreate(ctx->comm, &v2)); 2373 PetscCall(VecSetSizes(v2, PETSC_DECIDE, N / Nf)); 2374 PetscCall(VecSetFromOptions(v1)); // ??? 2375 PetscCall(VecSetFromOptions(v2)); 2376 // get each component 2377 PetscCall(VecGetBlockSize(Gsub, &bs)); 2378 PetscCheck(bs == Nf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " != num_species %" PetscInt_FMT " in Gsub", bs, Nf); 2379 PetscCall(VecGetBlockSize(Mfsub, &bs)); 2380 PetscCheck(bs == Nf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " != num_species %" PetscInt_FMT, bs, Nf); 2381 for (int i = 0, ix = ctx->species_offset[grid]; i < Nf; i++, ix++) { 2382 PetscScalar val; 2383 PetscCall(VecStrideGather(Gsub, i, v1, INSERT_VALUES)); // this is not right -- TODO 2384 PetscCall(VecStrideGather(Mfsub, i, v2, INSERT_VALUES)); 2385 PetscCall(VecDot(v1, v2, &val)); 2386 energy[ix] = PetscRealPart(val) * ctx->n_0 * ctx->v_0 * ctx->v_0 * ctx->masses[ix]; 2387 } 2388 PetscCall(VecDestroy(&v1)); 2389 PetscCall(VecDestroy(&v2)); 2390 } /* grids */ 2391 PetscCall(DMCompositeRestoreAccessArray(pack, globGamma, nDMs, NULL, globGammaArray)); 2392 PetscCall(DMCompositeRestoreAccessArray(pack, Mf, nDMs, NULL, globMfArray)); 2393 PetscCall(PetscFree(globGammaArray)); 2394 PetscCall(PetscFree(globMfArray)); 2395 PetscCall(VecDestroy(&globGamma)); 2396 PetscCall(VecDestroy(&Mf)); 2397 } 2398 } else { 2399 PetscCall(PetscDSSetObjective(prob, 0, &f0_s_v2)); 2400 PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx)); 2401 energy[ii] = 0.5 * tt[0] * ctx->n_0 * ctx->v_0 * ctx->v_0 * ctx->masses[ii]; 2402 } 2403 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%3" PetscInt_FMT ") species %" PetscInt_FMT ": density=%20.13e, x-momentum=%20.13e, y-momentum=%20.13e, z-momentum=%20.13e, energy=%21.13e", stepi, ii, (double)PetscRealPart(density[ii]), (double)PetscRealPart(xmomentum[ii]), (double)PetscRealPart(ymomentum[ii]), (double)PetscRealPart(zmomentum[ii]), (double)PetscRealPart(energy[ii]))); 2404 xmomentumtot += xmomentum[ii]; 2405 ymomentumtot += ymomentum[ii]; 2406 zmomentumtot += zmomentum[ii]; 2407 energytot += energy[ii]; 2408 densitytot += density[ii]; 2409 } 2410 if (ctx->num_species > 1) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n")); 2411 } 2412 } 2413 PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray)); 2414 PetscCall(PetscFree(globXArray)); 2415 /* totals */ 2416 PetscCall(DMPlexGetHeightStratum(ctx->plex[0], 0, &cStart, &cEnd)); 2417 if (ctx->num_species > 1) { 2418 if (dim == 2) { 2419 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\t%3" PetscInt_FMT ") Total: charge density=%21.13e, momentum=%21.13e, energy=%21.13e (m_i[0]/m_e = %g, %" PetscInt_FMT " cells on electron grid)", stepi, (double)PetscRealPart(densitytot), (double)PetscRealPart(zmomentumtot), (double)PetscRealPart(energytot), 2420 (double)(ctx->masses[1] / ctx->masses[0]), cEnd - cStart)); 2421 } else { 2422 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\t%3" PetscInt_FMT ") Total: charge density=%21.13e, x-momentum=%21.13e, y-momentum=%21.13e, z-momentum=%21.13e, energy=%21.13e (m_i[0]/m_e = %g, %" PetscInt_FMT " cells)", stepi, (double)PetscRealPart(densitytot), (double)PetscRealPart(xmomentumtot), (double)PetscRealPart(ymomentumtot), (double)PetscRealPart(zmomentumtot), (double)PetscRealPart(energytot), 2423 (double)(ctx->masses[1] / ctx->masses[0]), cEnd - cStart)); 2424 } 2425 } else PetscCall(PetscPrintf(PETSC_COMM_WORLD, " -- %" PetscInt_FMT " cells", cEnd - cStart)); 2426 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n")); 2427 PetscFunctionReturn(PETSC_SUCCESS); 2428 } 2429 2430 /*@ 2431 DMPlexLandauCreateMassMatrix - Create mass matrix for Landau in Plex space (not field major order of Jacobian) 2432 - puts mass matrix into ctx->M 2433 2434 Collective 2435 2436 Input Parameter: 2437 . pack - the DM object. Puts matrix in Landau context M field 2438 2439 Output Parameter: 2440 . Amat - The mass matrix (optional), mass matrix is added to the DM context 2441 2442 Level: beginner 2443 2444 .keywords: mesh 2445 .seealso: `DMPlexLandauCreateVelocitySpace()` 2446 @*/ 2447 PetscErrorCode DMPlexLandauCreateMassMatrix(DM pack, Mat *Amat) 2448 { 2449 DM mass_pack, massDM[LANDAU_MAX_GRIDS]; 2450 PetscDS prob; 2451 PetscInt ii, dim, N1 = 1, N2; 2452 LandauCtx *ctx; 2453 Mat packM, subM[LANDAU_MAX_GRIDS]; 2454 2455 PetscFunctionBegin; 2456 PetscValidHeaderSpecific(pack, DM_CLASSID, 1); 2457 if (Amat) PetscValidPointer(Amat, 2); 2458 PetscCall(DMGetApplicationContext(pack, &ctx)); 2459 PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context"); 2460 PetscCall(PetscLogEventBegin(ctx->events[14], 0, 0, 0, 0)); 2461 PetscCall(DMGetDimension(pack, &dim)); 2462 PetscCall(DMCompositeCreate(PetscObjectComm((PetscObject)pack), &mass_pack)); 2463 /* create pack mass matrix */ 2464 for (PetscInt grid = 0, ix = 0; grid < ctx->num_grids; grid++) { 2465 PetscCall(DMClone(ctx->plex[grid], &massDM[grid])); 2466 PetscCall(DMCopyFields(ctx->plex[grid], massDM[grid])); 2467 PetscCall(DMCreateDS(massDM[grid])); 2468 PetscCall(DMGetDS(massDM[grid], &prob)); 2469 for (ix = 0, ii = ctx->species_offset[grid]; ii < ctx->species_offset[grid + 1]; ii++, ix++) { 2470 if (dim == 3) PetscCall(PetscDSSetJacobian(prob, ix, ix, g0_1, NULL, NULL, NULL)); 2471 else PetscCall(PetscDSSetJacobian(prob, ix, ix, g0_r, NULL, NULL, NULL)); 2472 } 2473 #if !defined(LANDAU_SPECIES_MAJOR) 2474 PetscCall(DMCompositeAddDM(mass_pack, massDM[grid])); 2475 #else 2476 for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { // add batch size DMs for this species grid 2477 PetscCall(DMCompositeAddDM(mass_pack, massDM[grid])); 2478 } 2479 #endif 2480 PetscCall(DMCreateMatrix(massDM[grid], &subM[grid])); 2481 } 2482 #if !defined(LANDAU_SPECIES_MAJOR) 2483 // stack the batched DMs 2484 for (PetscInt b_id = 1; b_id < ctx->batch_sz; b_id++) { 2485 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(DMCompositeAddDM(mass_pack, massDM[grid])); 2486 } 2487 #endif 2488 PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only")); 2489 PetscCall(DMCreateMatrix(mass_pack, &packM)); 2490 PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only false")); 2491 PetscCall(MatSetOption(packM, MAT_STRUCTURALLY_SYMMETRIC, PETSC_TRUE)); 2492 PetscCall(MatSetOption(packM, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE)); 2493 PetscCall(DMDestroy(&mass_pack)); 2494 /* make mass matrix for each block */ 2495 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 2496 Vec locX; 2497 DM plex = massDM[grid]; 2498 PetscCall(DMGetLocalVector(plex, &locX)); 2499 /* Mass matrix is independent of the input, so no need to fill locX */ 2500 PetscCall(DMPlexSNESComputeJacobianFEM(plex, locX, subM[grid], subM[grid], ctx)); 2501 PetscCall(DMRestoreLocalVector(plex, &locX)); 2502 PetscCall(DMDestroy(&massDM[grid])); 2503 } 2504 PetscCall(MatGetSize(ctx->J, &N1, NULL)); 2505 PetscCall(MatGetSize(packM, &N2, NULL)); 2506 PetscCheck(N1 == N2, PetscObjectComm((PetscObject)pack), PETSC_ERR_PLIB, "Incorrect matrix sizes: |Jacobian| = %" PetscInt_FMT ", |Mass|=%" PetscInt_FMT, N1, N2); 2507 /* assemble block diagonals */ 2508 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 2509 Mat B = subM[grid]; 2510 PetscInt nloc, nzl, *colbuf, COL_BF_SIZE = 1024, row; 2511 PetscCall(PetscMalloc(sizeof(*colbuf) * COL_BF_SIZE, &colbuf)); 2512 PetscCall(MatGetSize(B, &nloc, NULL)); 2513 for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { 2514 const PetscInt moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset); 2515 const PetscInt *cols; 2516 const PetscScalar *vals; 2517 for (int i = 0; i < nloc; i++) { 2518 PetscCall(MatGetRow(B, i, &nzl, NULL, NULL)); 2519 if (nzl > COL_BF_SIZE) { 2520 PetscCall(PetscFree(colbuf)); 2521 PetscCall(PetscInfo(pack, "Realloc buffer %" PetscInt_FMT " to %" PetscInt_FMT " (row size %" PetscInt_FMT ") \n", COL_BF_SIZE, 2 * COL_BF_SIZE, nzl)); 2522 COL_BF_SIZE = nzl; 2523 PetscCall(PetscMalloc(sizeof(*colbuf) * COL_BF_SIZE, &colbuf)); 2524 } 2525 PetscCall(MatGetRow(B, i, &nzl, &cols, &vals)); 2526 for (int j = 0; j < nzl; j++) colbuf[j] = cols[j] + moffset; 2527 row = i + moffset; 2528 PetscCall(MatSetValues(packM, 1, &row, nzl, colbuf, vals, INSERT_VALUES)); 2529 PetscCall(MatRestoreRow(B, i, &nzl, &cols, &vals)); 2530 } 2531 } 2532 PetscCall(PetscFree(colbuf)); 2533 } 2534 // cleanup 2535 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(MatDestroy(&subM[grid])); 2536 PetscCall(MatAssemblyBegin(packM, MAT_FINAL_ASSEMBLY)); 2537 PetscCall(MatAssemblyEnd(packM, MAT_FINAL_ASSEMBLY)); 2538 PetscCall(PetscObjectSetName((PetscObject)packM, "mass")); 2539 PetscCall(MatViewFromOptions(packM, NULL, "-dm_landau_mass_view")); 2540 ctx->M = packM; 2541 if (Amat) *Amat = packM; 2542 PetscCall(PetscLogEventEnd(ctx->events[14], 0, 0, 0, 0)); 2543 PetscFunctionReturn(PETSC_SUCCESS); 2544 } 2545 2546 /*@ 2547 DMPlexLandauIFunction - TS residual calculation, confusingly this computes the Jacobian w/o mass 2548 2549 Collective 2550 2551 Input Parameters: 2552 + TS - The time stepping context 2553 . time_dummy - current time (not used) 2554 . X - Current state 2555 . X_t - Time derivative of current state 2556 - actx - Landau context 2557 2558 Output Parameter: 2559 . F - The residual 2560 2561 Level: beginner 2562 2563 .keywords: mesh 2564 .seealso: `DMPlexLandauCreateVelocitySpace()`, `DMPlexLandauIJacobian()` 2565 @*/ 2566 PetscErrorCode DMPlexLandauIFunction(TS ts, PetscReal time_dummy, Vec X, Vec X_t, Vec F, void *actx) 2567 { 2568 LandauCtx *ctx = (LandauCtx *)actx; 2569 PetscInt dim; 2570 DM pack; 2571 #if defined(PETSC_HAVE_THREADSAFETY) 2572 double starttime, endtime; 2573 #endif 2574 PetscObjectState state; 2575 2576 PetscFunctionBegin; 2577 PetscCall(TSGetDM(ts, &pack)); 2578 PetscCall(DMGetApplicationContext(pack, &ctx)); 2579 PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context"); 2580 if (ctx->stage) PetscCall(PetscLogStagePush(ctx->stage)); 2581 PetscCall(PetscLogEventBegin(ctx->events[11], 0, 0, 0, 0)); 2582 PetscCall(PetscLogEventBegin(ctx->events[0], 0, 0, 0, 0)); 2583 #if defined(PETSC_HAVE_THREADSAFETY) 2584 starttime = MPI_Wtime(); 2585 #endif 2586 PetscCall(DMGetDimension(pack, &dim)); 2587 PetscCall(PetscObjectStateGet((PetscObject)ctx->J, &state)); 2588 if (state != ctx->norm_state) { 2589 PetscCall(PetscInfo(ts, "Create Landau Jacobian t=%g J.state %" PetscInt64_FMT " --> %" PetscInt64_FMT "\n", (double)time_dummy, ctx->norm_state, state)); 2590 PetscCall(MatZeroEntries(ctx->J)); 2591 PetscCall(LandauFormJacobian_Internal(X, ctx->J, dim, 0.0, (void *)ctx)); 2592 PetscCall(MatViewFromOptions(ctx->J, NULL, "-dm_landau_jacobian_view")); 2593 PetscCall(PetscObjectStateGet((PetscObject)ctx->J, &state)); 2594 ctx->norm_state = state; 2595 } else { 2596 PetscCall(PetscInfo(ts, "WARNING Skip forming Jacobian, has not changed %" PetscInt64_FMT "\n", state)); 2597 } 2598 /* mat vec for op */ 2599 PetscCall(MatMult(ctx->J, X, F)); /* C*f */ 2600 /* add time term */ 2601 if (X_t) PetscCall(MatMultAdd(ctx->M, X_t, F, F)); 2602 #if defined(PETSC_HAVE_THREADSAFETY) 2603 if (ctx->stage) { 2604 endtime = MPI_Wtime(); 2605 ctx->times[LANDAU_OPERATOR] += (endtime - starttime); 2606 ctx->times[LANDAU_JACOBIAN] += (endtime - starttime); 2607 ctx->times[LANDAU_MATRIX_TOTAL] += (endtime - starttime); 2608 ctx->times[LANDAU_JACOBIAN_COUNT] += 1; 2609 } 2610 #endif 2611 PetscCall(PetscLogEventEnd(ctx->events[0], 0, 0, 0, 0)); 2612 PetscCall(PetscLogEventEnd(ctx->events[11], 0, 0, 0, 0)); 2613 if (ctx->stage) PetscCall(PetscLogStagePop()); 2614 PetscFunctionReturn(PETSC_SUCCESS); 2615 } 2616 2617 /*@ 2618 DMPlexLandauIJacobian - TS Jacobian construction, confusingly this adds mass 2619 2620 Collective 2621 2622 Input Parameters: 2623 + TS - The time stepping context 2624 . time_dummy - current time (not used) 2625 . X - Current state 2626 . U_tdummy - Time derivative of current state (not used) 2627 . shift - shift for du/dt term 2628 - actx - Landau context 2629 2630 Output Parameters: 2631 + Amat - Jacobian 2632 - Pmat - same as Amat 2633 2634 Level: beginner 2635 2636 .keywords: mesh 2637 .seealso: `DMPlexLandauCreateVelocitySpace()`, `DMPlexLandauIFunction()` 2638 @*/ 2639 PetscErrorCode DMPlexLandauIJacobian(TS ts, PetscReal time_dummy, Vec X, Vec U_tdummy, PetscReal shift, Mat Amat, Mat Pmat, void *actx) 2640 { 2641 LandauCtx *ctx = NULL; 2642 PetscInt dim; 2643 DM pack; 2644 #if defined(PETSC_HAVE_THREADSAFETY) 2645 double starttime, endtime; 2646 #endif 2647 PetscObjectState state; 2648 2649 PetscFunctionBegin; 2650 PetscCall(TSGetDM(ts, &pack)); 2651 PetscCall(DMGetApplicationContext(pack, &ctx)); 2652 PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context"); 2653 PetscCheck(Amat == Pmat && Amat == ctx->J, ctx->comm, PETSC_ERR_PLIB, "Amat!=Pmat || Amat!=ctx->J"); 2654 PetscCall(DMGetDimension(pack, &dim)); 2655 /* get collision Jacobian into A */ 2656 if (ctx->stage) PetscCall(PetscLogStagePush(ctx->stage)); 2657 PetscCall(PetscLogEventBegin(ctx->events[11], 0, 0, 0, 0)); 2658 PetscCall(PetscLogEventBegin(ctx->events[9], 0, 0, 0, 0)); 2659 #if defined(PETSC_HAVE_THREADSAFETY) 2660 starttime = MPI_Wtime(); 2661 #endif 2662 PetscCall(PetscInfo(ts, "Adding mass to Jacobian t=%g, shift=%g\n", (double)time_dummy, (double)shift)); 2663 PetscCheck(shift != 0.0, ctx->comm, PETSC_ERR_PLIB, "zero shift"); 2664 PetscCall(PetscObjectStateGet((PetscObject)ctx->J, &state)); 2665 PetscCheck(state == ctx->norm_state, ctx->comm, PETSC_ERR_PLIB, "wrong state, %" PetscInt64_FMT " %" PetscInt64_FMT "", ctx->norm_state, state); 2666 if (!ctx->use_matrix_mass) { 2667 PetscCall(LandauFormJacobian_Internal(X, ctx->J, dim, shift, (void *)ctx)); 2668 } else { /* add mass */ 2669 PetscCall(MatAXPY(Pmat, shift, ctx->M, SAME_NONZERO_PATTERN)); 2670 } 2671 #if defined(PETSC_HAVE_THREADSAFETY) 2672 if (ctx->stage) { 2673 endtime = MPI_Wtime(); 2674 ctx->times[LANDAU_OPERATOR] += (endtime - starttime); 2675 ctx->times[LANDAU_MASS] += (endtime - starttime); 2676 ctx->times[LANDAU_MATRIX_TOTAL] += (endtime - starttime); 2677 } 2678 #endif 2679 PetscCall(PetscLogEventEnd(ctx->events[9], 0, 0, 0, 0)); 2680 PetscCall(PetscLogEventEnd(ctx->events[11], 0, 0, 0, 0)); 2681 if (ctx->stage) PetscCall(PetscLogStagePop()); 2682 PetscFunctionReturn(PETSC_SUCCESS); 2683 } 2684