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