1 static char help[] = "Grid based Landau collision operator with PIC interface with OpenMP setup. (one species per grid)\n"; 2 3 /* 4 Support 2.5V with axisymmetric coordinates 5 - r,z coordinates 6 - Domain and species data input by Landau operator 7 - "radius" for each grid, normalized with electron thermal velocity 8 - Domain: (0,radius) x (-radius,radius), thus first coordinate x[0] is perpendicular velocity and 2pi*x[0] term is added for axisymmetric 9 Supports full 3V 10 11 */ 12 13 #include <petscdmplex.h> 14 #include <petscds.h> 15 #include <petscdmswarm.h> 16 #include <petscksp.h> 17 #include <petsc/private/petscimpl.h> 18 #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY) 19 #include <omp.h> 20 #endif 21 #include <petsclandau.h> 22 #include <petscdmcomposite.h> 23 24 typedef struct { 25 Mat MpTrans; 26 Mat Mp; 27 Vec ff; 28 Vec uu; 29 } MatShellCtx; 30 31 typedef struct { 32 PetscInt v_target; 33 PetscInt g_target; 34 PetscInt global_vertex_id_0; 35 DM *globSwarmArray; 36 LandauCtx *ctx; 37 DM *grid_dm; 38 Mat *g_Mass; 39 Mat *globMpArray; 40 Vec *globXArray; 41 PetscBool print; 42 PetscBool print_entropy; 43 } PrintCtx; 44 45 PetscErrorCode MatMultMtM_SeqAIJ(Mat MtM, Vec xx, Vec yy) 46 { 47 MatShellCtx *matshellctx; 48 49 PetscFunctionBeginUser; 50 PetscCall(MatShellGetContext(MtM, &matshellctx)); 51 PetscCheck(matshellctx, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context"); 52 PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff)); 53 PetscCall(MatMult(matshellctx->MpTrans, matshellctx->ff, yy)); 54 PetscFunctionReturn(PETSC_SUCCESS); 55 } 56 57 PetscErrorCode MatMultAddMtM_SeqAIJ(Mat MtM, Vec xx, Vec yy, Vec zz) 58 { 59 MatShellCtx *matshellctx; 60 61 PetscFunctionBeginUser; 62 PetscCall(MatShellGetContext(MtM, &matshellctx)); 63 PetscCheck(matshellctx, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context"); 64 PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff)); 65 PetscCall(MatMultAdd(matshellctx->MpTrans, matshellctx->ff, yy, zz)); 66 PetscFunctionReturn(PETSC_SUCCESS); 67 } 68 69 PetscErrorCode createSwarm(const DM dm, PetscInt dim, DM *sw) 70 { 71 PetscInt Nc = 1; 72 73 PetscFunctionBeginUser; 74 PetscCall(DMCreate(PETSC_COMM_SELF, sw)); 75 PetscCall(DMSetType(*sw, DMSWARM)); 76 PetscCall(DMSetDimension(*sw, dim)); 77 PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC)); 78 PetscCall(DMSwarmSetCellDM(*sw, dm)); 79 PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", Nc, PETSC_REAL)); 80 PetscCall(DMSwarmFinalizeFieldRegister(*sw)); 81 PetscCall(DMSetFromOptions(*sw)); 82 PetscCall(PetscObjectSetName((PetscObject)*sw, "Particle Grid")); 83 PetscFunctionReturn(PETSC_SUCCESS); 84 } 85 86 static PetscErrorCode makeSwarm(DM sw, const PetscInt dim, const PetscInt Np, const PetscReal xx[], const PetscReal yy[], const PetscReal zz[]) 87 { 88 PetscReal *coords; 89 PetscDataType dtype; 90 PetscInt bs, p, zero = 0; 91 92 PetscFunctionBeginUser; 93 PetscCall(DMSwarmSetLocalSizes(sw, Np, zero)); 94 PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords)); 95 for (p = 0; p < Np; p++) { 96 coords[p * dim + 0] = xx[p]; 97 coords[p * dim + 1] = yy[p]; 98 if (dim == 3) coords[p * dim + 2] = zz[p]; 99 } 100 PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords)); 101 PetscFunctionReturn(PETSC_SUCCESS); 102 } 103 104 static PetscErrorCode createMp(const DM dm, DM sw, Mat *Mp_out) 105 { 106 PetscBool removePoints = PETSC_TRUE; 107 Mat M_p; 108 109 PetscFunctionBeginUser; 110 // migrate after coords are set 111 PetscCall(DMSwarmMigrate(sw, removePoints)); 112 // 113 PetscCall(PetscObjectSetName((PetscObject)sw, "Particle Grid")); 114 115 /* PetscInt N,*count,nmin=10000,nmax=0,ntot=0; */ 116 /* // count */ 117 /* PetscCall(DMSwarmCreatePointPerCellCount(sw, &N, &count)); */ 118 /* for (int i=0, n; i< N ; i++) { */ 119 /* if ((n=count[i]) > nmax) nmax = n; */ 120 /* if (n < nmin) nmin = n; */ 121 /* PetscCall(PetscInfo(dm, " %d) %d particles\n", i, n)); */ 122 /* ntot += n; */ 123 /* } */ 124 /* PetscCall(PetscFree(count)); */ 125 /* PetscCall(PetscInfo(dm, " %" PetscInt_FMT " max particle / cell, and %" PetscInt_FMT " min, ratio = %g, %" PetscInt_FMT " total\n", nmax, nmin, (double)nmax/(double)nmin,ntot)); */ 126 127 /* This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE */ 128 PetscCall(DMCreateMassMatrix(sw, dm, &M_p)); 129 PetscCall(DMViewFromOptions(sw, NULL, "-ex30_sw_view")); 130 // output 131 *Mp_out = M_p; 132 PetscFunctionReturn(PETSC_SUCCESS); 133 } 134 135 static PetscErrorCode particlesToGrid(const DM dm, DM sw, const PetscInt a_tid, const PetscInt dim, const PetscReal a_wp[], Vec rho, Mat M_p) 136 { 137 PetscReal *wq; 138 PetscDataType dtype; 139 Vec ff; 140 PetscInt bs, p, Np; 141 142 PetscFunctionBeginUser; 143 PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wq)); 144 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 145 for (p = 0; p < Np; p++) wq[p] = a_wp[p]; 146 PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wq)); 147 PetscCall(PetscObjectSetName((PetscObject)rho, "rho")); 148 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff)); 149 PetscCall(PetscObjectSetName((PetscObject)ff, "weights")); 150 PetscCall(MatMultTranspose(M_p, ff, rho)); 151 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff)); 152 PetscFunctionReturn(PETSC_SUCCESS); 153 } 154 155 // 156 // add grid to arg 'sw.w_q' 157 // 158 PetscErrorCode gridToParticles(const DM dm, DM sw, const Vec rhs, Vec work_ferhs, Mat M_p, Mat Mass) 159 { 160 PetscBool is_lsqr; 161 KSP ksp; 162 Mat PM_p = NULL, MtM, D = NULL; 163 Vec ff; 164 PetscInt N, M, nzl; 165 MatShellCtx *matshellctx = NULL; 166 PC pc; 167 168 PetscFunctionBeginUser; 169 // (Mp' Mp)^-1 M 170 PetscCall(MatMult(Mass, rhs, work_ferhs)); 171 // pseudo-inverse 172 PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp)); 173 PetscCall(KSPSetType(ksp, KSPCG)); 174 PetscCall(KSPGetPC(ksp, &pc)); 175 PetscCall(PCSetType(pc, PCJACOBI)); 176 PetscCall(KSPSetOptionsPrefix(ksp, "ftop_")); 177 PetscCall(KSPSetFromOptions(ksp)); 178 PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPLSQR, &is_lsqr)); 179 if (!is_lsqr) { 180 PetscCall(MatGetLocalSize(M_p, &M, &N)); 181 if (N > M) { 182 PetscCall(PetscInfo(ksp, " M (%" PetscInt_FMT ") < M (%" PetscInt_FMT ") more vertices than particles: revert to lsqr\n", M, N)); 183 is_lsqr = PETSC_TRUE; 184 PetscCall(KSPSetType(ksp, KSPLSQR)); 185 PetscCall(PCSetType(pc, PCNONE)); // should not happen, but could solve stable (Mp^T Mp), move projection Mp before solve 186 } else { 187 PetscCall(PetscNew(&matshellctx)); 188 PetscCall(MatCreateVecs(M_p, &matshellctx->uu, &matshellctx->ff)); 189 if (0) { 190 PetscCall(MatTransposeMatMult(M_p, M_p, MAT_INITIAL_MATRIX, 4, &MtM)); 191 PetscCall(KSPSetOperators(ksp, MtM, MtM)); 192 PetscCall(PetscInfo(M_p, "createMtM KSP with explicit Mp'Mp\n")); 193 PetscCall(MatViewFromOptions(MtM, NULL, "-ftop2_MtM_mat_view")); 194 } else { 195 PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dm), N, N, PETSC_DECIDE, PETSC_DECIDE, matshellctx, &MtM)); 196 PetscCall(MatTranspose(M_p, MAT_INITIAL_MATRIX, &matshellctx->MpTrans)); 197 matshellctx->Mp = M_p; 198 PetscCall(MatShellSetOperation(MtM, MATOP_MULT, (void (*)(void))MatMultMtM_SeqAIJ)); 199 PetscCall(MatShellSetOperation(MtM, MATOP_MULT_ADD, (void (*)(void))MatMultAddMtM_SeqAIJ)); 200 PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, N, N, 1, NULL, &D)); 201 PetscCall(MatViewFromOptions(matshellctx->MpTrans, NULL, "-ftop2_MpT_mat_view")); 202 for (int i = 0; i < N; i++) { 203 const PetscScalar *vals; 204 const PetscInt *cols; 205 PetscScalar dot = 0; 206 PetscCall(MatGetRow(matshellctx->MpTrans, i, &nzl, &cols, &vals)); 207 for (int ii = 0; ii < nzl; ii++) dot += PetscSqr(vals[ii]); 208 if (dot < PETSC_MACHINE_EPSILON) { 209 PetscCall(PetscInfo(ksp, "empty row in pseudo-inverse %d\n", i)); 210 is_lsqr = PETSC_TRUE; // empty rows 211 PetscCall(KSPSetType(ksp, KSPLSQR)); 212 PetscCall(PCSetType(pc, PCNONE)); // should not happen, but could solve stable (Mp Mp^T), move projection Mp before solve 213 // clean up 214 PetscCall(MatDestroy(&matshellctx->MpTrans)); 215 PetscCall(VecDestroy(&matshellctx->ff)); 216 PetscCall(VecDestroy(&matshellctx->uu)); 217 PetscCall(MatDestroy(&D)); 218 PetscCall(MatDestroy(&MtM)); 219 PetscCall(PetscFree(matshellctx)); 220 D = NULL; 221 break; 222 } 223 PetscCall(MatSetValue(D, i, i, dot, INSERT_VALUES)); 224 } 225 if (D) { 226 PetscCall(MatAssemblyBegin(D, MAT_FINAL_ASSEMBLY)); 227 PetscCall(MatAssemblyEnd(D, MAT_FINAL_ASSEMBLY)); 228 PetscCall(PetscInfo(M_p, "createMtMKSP Have %" PetscInt_FMT " eqs, nzl = %" PetscInt_FMT "\n", N, nzl)); 229 PetscCall(KSPSetOperators(ksp, MtM, D)); 230 PetscCall(MatViewFromOptions(D, NULL, "-ftop2_D_mat_view")); 231 PetscCall(MatViewFromOptions(M_p, NULL, "-ftop2_Mp_mat_view")); 232 PetscCall(MatViewFromOptions(matshellctx->MpTrans, NULL, "-ftop2_MpTranspose_mat_view")); 233 PetscCall(MatViewFromOptions(MtM, NULL, "-ftop2_MtM_mat_view")); 234 } 235 } 236 } 237 } 238 if (is_lsqr) { 239 PC pc2; 240 PetscBool is_bjac; 241 PetscCall(KSPGetPC(ksp, &pc2)); 242 PetscCall(PetscObjectTypeCompare((PetscObject)pc2, PCBJACOBI, &is_bjac)); 243 if (is_bjac) { 244 PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p)); 245 PetscCall(KSPSetOperators(ksp, M_p, PM_p)); 246 } else { 247 PetscCall(KSPSetOperators(ksp, M_p, M_p)); 248 } 249 PetscCall(MatViewFromOptions(M_p, NULL, "-ftop2_Mp_mat_view")); 250 } 251 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff)); // this grabs access 252 if (!is_lsqr) { 253 PetscErrorCode ierr; 254 ierr = KSPSolve(ksp, work_ferhs, matshellctx->uu); 255 if (!ierr) { 256 PetscCall(MatMult(M_p, matshellctx->uu, ff)); 257 } else { // failed 258 PC pc2; 259 PetscBool is_bjac; 260 PetscCall(PetscInfo(ksp, "Solver failed, probably singular, try lsqr\n")); 261 PetscCall(KSPReset(ksp)); 262 PetscCall(KSPSetType(ksp, KSPLSQR)); 263 PetscCall(KSPGetPC(ksp, &pc2)); 264 PetscCall(PCSetType(pc2, PCNONE)); // should not happen, but could solve stable (Mp Mp^T), move projection Mp before solve 265 PetscCall(KSPSetOptionsPrefix(ksp, "ftop_")); 266 PetscCall(KSPSetFromOptions(ksp)); 267 PetscCall(PetscObjectTypeCompare((PetscObject)pc2, PCBJACOBI, &is_bjac)); 268 if (is_bjac) { 269 PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p)); 270 PetscCall(KSPSetOperators(ksp, M_p, PM_p)); 271 } else { 272 PetscCall(KSPSetOperators(ksp, M_p, M_p)); 273 } 274 ierr = KSPSolveTranspose(ksp, work_ferhs, ff); 275 if (ierr) { PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "backup LSQR solver failed - need to add N_v > N_p Moore-Penrose pseudo-inverse"); } 276 } 277 if (D) PetscCall(MatDestroy(&D)); 278 PetscCall(MatDestroy(&MtM)); 279 if (matshellctx->MpTrans) PetscCall(MatDestroy(&matshellctx->MpTrans)); 280 PetscCall(VecDestroy(&matshellctx->ff)); 281 PetscCall(VecDestroy(&matshellctx->uu)); 282 PetscCall(PetscFree(matshellctx)); 283 } else { 284 PetscErrorCode ierr; 285 ierr = KSPSolveTranspose(ksp, work_ferhs, ff); 286 if (ierr) { PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "backup LSQR solver failed - need to add N_v > N_p Moore-Penrose pseudo-inverse"); } 287 } 288 PetscCall(KSPDestroy(&ksp)); 289 PetscCall(MatDestroy(&PM_p)); 290 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff)); 291 PetscFunctionReturn(PETSC_SUCCESS); 292 } 293 294 #define EX30_MAX_NUM_THRDS 12 295 #define EX30_MAX_BATCH_SZ 1024 296 // 297 // add grid to arg 'globSwarmArray[].w_q' 298 // 299 PetscErrorCode gridToParticles_private(DM grid_dm[], DM globSwarmArray[], const PetscInt dim, const PetscInt v_target, const PetscInt numthreads, const PetscInt num_vertices, const PetscInt global_vertex_id, Mat globMpArray[], Mat g_Mass[], Vec t_fhat[][EX30_MAX_NUM_THRDS], PetscReal moments[], Vec globXArray[], LandauCtx *ctx) 300 { 301 PetscErrorCode ierr = (PetscErrorCode)0; // used for inside thread loops 302 303 PetscFunctionBeginUser; 304 // map back to particles 305 for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) { 306 PetscCall(PetscInfo(grid_dm[0], "g2p: global batch %" PetscInt_FMT " of %" PetscInt_FMT ", Landau batch %" PetscInt_FMT " of %" PetscInt_FMT ": map back to particles\n", global_vertex_id + 1, num_vertices, v_id_0 + 1, ctx->batch_sz)); 307 //PetscPragmaOMP(parallel for) 308 for (int tid = 0; tid < numthreads; tid++) { 309 const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id + v_id; 310 if (glb_v_id < num_vertices) { 311 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids 312 PetscErrorCode ierr_t; 313 ierr_t = PetscInfo(grid_dm[0], "gridToParticles: global batch %" PetscInt_FMT ", local batch b=%" PetscInt_FMT ", grid g=%" PetscInt_FMT ", index(b,g) %" PetscInt_FMT "\n", global_vertex_id, v_id, grid, LAND_PACK_IDX(v_id, grid)); 314 ierr_t = gridToParticles(grid_dm[grid], globSwarmArray[LAND_PACK_IDX(v_id, grid)], globXArray[LAND_PACK_IDX(v_id, grid)], t_fhat[grid][tid], globMpArray[LAND_PACK_IDX(v_id, grid)], g_Mass[grid]); 315 if (ierr_t) ierr = ierr_t; 316 } 317 } 318 } 319 PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in OMP loop. ierr = %d", (int)ierr); 320 /* Get moments */ 321 PetscCall(PetscInfo(grid_dm[0], "Cleanup batches %" PetscInt_FMT " to %" PetscInt_FMT "\n", v_id_0, v_id_0 + numthreads)); 322 for (int tid = 0; tid < numthreads; tid++) { 323 const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id + v_id; 324 if (glb_v_id == v_target) { 325 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 326 PetscDataType dtype; 327 PetscReal *wp, *coords; 328 DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)]; 329 PetscInt npoints, bs = 1; 330 PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); // take data out here 331 PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords)); 332 PetscCall(DMSwarmGetLocalSize(sw, &npoints)); 333 for (int p = 0; p < npoints; p++) { 334 PetscReal v2 = 0, fact = (dim == 2) ? 2.0 * PETSC_PI * coords[p * dim + 0] : 1, w = fact * wp[p] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]]; 335 for (int i = 0; i < dim; ++i) v2 += PetscSqr(coords[p * dim + i]); 336 moments[0] += w; 337 moments[1] += w * ctx->v_0 * coords[p * dim + 1]; // z-momentum 338 moments[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2; 339 } 340 PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp)); 341 PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords)); 342 } 343 const PetscReal N_inv = 1 / moments[0]; 344 PetscCall(PetscInfo(grid_dm[0], "gridToParticles_private [%" PetscInt_FMT "], n = %g\n", v_id, (double)moments[0])); 345 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 346 PetscDataType dtype; 347 PetscReal *wp, *coords; 348 DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)]; 349 PetscInt npoints, bs = 1; 350 PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); // take data out here 351 PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords)); 352 PetscCall(DMSwarmGetLocalSize(sw, &npoints)); 353 for (int p = 0; p < npoints; p++) { 354 const PetscReal fact = dim == 2 ? 2.0 * PETSC_PI * coords[p * dim + 0] : 1, w = fact * wp[p] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]], ww = w * N_inv; 355 if (w > PETSC_REAL_MIN) { 356 moments[3] -= ww * PetscLogReal(ww); 357 PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "ww (%g) > 1", (double)ww); 358 } else moments[4] -= w; // keep track of density that is lost 359 } 360 PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp)); 361 PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords)); 362 } 363 } 364 } // thread batch 365 } // batch 366 PetscFunctionReturn(PETSC_SUCCESS); 367 } 368 369 static void maxwellian(PetscInt dim, const PetscReal x[], PetscReal kt_m, PetscReal n, PetscReal shift, PetscScalar *u) 370 { 371 PetscInt i; 372 PetscReal v2 = 0, theta = 2.0 * kt_m; /* theta = 2kT/mc^2 */ 373 374 if (shift != 0.) { 375 v2 = 0; 376 for (i = 0; i < dim - 1; ++i) v2 += x[i] * x[i]; 377 v2 += (x[dim - 1] - shift) * (x[dim - 1] - shift); 378 /* evaluate the shifted Maxwellian */ 379 u[0] += n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta)); 380 } else { 381 /* compute the exponents, v^2 */ 382 for (i = 0; i < dim; ++i) v2 += x[i] * x[i]; 383 /* evaluate the Maxwellian */ 384 u[0] += n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta)); 385 } 386 } 387 388 static PetscErrorCode PostStep(TS ts) 389 { 390 PetscInt n, dim, nDMs, v_id; 391 PetscReal t; 392 LandauCtx *ctx; 393 Vec X; 394 PrintCtx *printCtx; 395 DM pack; 396 PetscReal moments[5], e_grid[LANDAU_MAX_GRIDS]; 397 398 PetscFunctionBeginUser; 399 PetscCall(TSGetApplicationContext(ts, &printCtx)); 400 if (!printCtx->print && !printCtx->print_entropy) PetscFunctionReturn(PETSC_SUCCESS); 401 ctx = printCtx->ctx; 402 if (printCtx->v_target < printCtx->global_vertex_id_0 || printCtx->v_target >= printCtx->global_vertex_id_0 + ctx->batch_sz) PetscFunctionReturn(PETSC_SUCCESS); 403 for (int i = 0; i < 5; i++) moments[i] = 0; 404 for (int i = 0; i < LANDAU_MAX_GRIDS; i++) e_grid[i] = 0; 405 v_id = printCtx->v_target % ctx->batch_sz; 406 PetscCall(TSGetDM(ts, &pack)); 407 PetscCall(DMGetDimension(pack, &dim)); 408 PetscCall(DMCompositeGetNumberDM(pack, &nDMs)); // number of vertices * number of grids 409 PetscCall(TSGetSolution(ts, &X)); 410 PetscCall(TSGetStepNumber(ts, &n)); 411 PetscCall(TSGetTime(ts, &t)); 412 PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, printCtx->globXArray)); 413 if (printCtx->print_entropy && printCtx->v_target >= 0 && 0) { 414 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 415 PetscDataType dtype; 416 PetscReal *wp, *coords; 417 DM sw = printCtx->globSwarmArray[LAND_PACK_IDX(v_id, grid)]; 418 Vec work, subX = printCtx->globXArray[LAND_PACK_IDX(v_id, grid)]; 419 PetscInt bs, NN; 420 // C-G moments 421 PetscCall(VecDuplicate(subX, &work)); 422 PetscCall(gridToParticles(printCtx->grid_dm[grid], sw, subX, work, printCtx->globMpArray[LAND_PACK_IDX(v_id, grid)], printCtx->g_Mass[grid])); 423 PetscCall(VecDestroy(&work)); 424 // moments 425 PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords)); 426 PetscCall(DMSwarmGetLocalSize(sw, &NN)); 427 PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); 428 for (int pp = 0; pp < NN; pp++) { 429 PetscReal v2 = 0, fact = (dim == 2) ? 2.0 * PETSC_PI * coords[pp * dim + 0] : 1, w = fact * wp[pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]]; 430 for (int i = 0; i < dim; ++i) v2 += PetscSqr(coords[pp * dim + i]); 431 moments[0] += w; 432 moments[1] += w * ctx->v_0 * coords[pp * dim + 1]; // z-momentum 433 moments[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2; 434 e_grid[grid] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2; 435 } 436 PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords)); 437 PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp)); 438 } 439 // entropy 440 const PetscReal N_inv = 1 / moments[0]; 441 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 442 PetscDataType dtype; 443 PetscReal *wp, *coords; 444 DM sw = printCtx->globSwarmArray[LAND_PACK_IDX(v_id, grid)]; 445 PetscInt bs, NN; 446 PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords)); 447 PetscCall(DMSwarmGetLocalSize(sw, &NN)); 448 PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); 449 for (int pp = 0; pp < NN; pp++) { 450 PetscReal fact = (dim == 2) ? 2.0 * PETSC_PI * coords[pp * dim + 0] : 1, w = fact * wp[pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]], ww = w * N_inv; 451 if (w > PETSC_REAL_MIN) { 452 moments[3] -= ww * PetscLogReal(ww); 453 } else moments[4] -= w; 454 } 455 PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords)); 456 PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp)); 457 } 458 PetscCall(PetscInfo(X, "%4d) time %e, Landau particle moments: 0: %18.12e 1: %19.12e 2: %18.12e entropy: %e loss %e. energy = %e + %e + %e\n", (int)n, (double)t, (double)moments[0], (double)moments[1], (double)moments[2], (double)moments[3], (double)(moments[4] / moments[0]), (double)e_grid[0], (double)e_grid[1], (double)e_grid[2])); 459 } 460 if (printCtx->print && printCtx->g_target >= 0) { 461 PetscInt grid = printCtx->g_target, id; 462 static PetscReal last_t = -100000, period = .5; 463 if (last_t == -100000) last_t = -period + t; 464 if (t >= last_t + period) { 465 last_t = t; 466 PetscCall(DMGetOutputSequenceNumber(ctx->plex[grid], &id, NULL)); 467 PetscCall(DMSetOutputSequenceNumber(ctx->plex[grid], id + 1, t)); 468 PetscCall(VecViewFromOptions(printCtx->globXArray[LAND_PACK_IDX(v_id % ctx->batch_sz, grid)], NULL, "-ex30_vec_view")); 469 if (ctx->num_grids > grid + 1) { 470 PetscCall(DMSetOutputSequenceNumber(ctx->plex[grid + 1], id + 1, t)); 471 PetscCall(VecViewFromOptions(printCtx->globXArray[LAND_PACK_IDX(v_id % ctx->batch_sz, grid + 1)], NULL, "-ex30_vec_view2")); 472 } 473 PetscCall(PetscInfo(X, "%4d) time %e View\n", (int)n, (double)t)); 474 } 475 } 476 PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, printCtx->globXArray)); 477 PetscFunctionReturn(PETSC_SUCCESS); 478 } 479 480 PetscErrorCode go(TS ts, Vec X, const PetscInt num_vertices, const PetscInt a_Np, const PetscInt dim, const PetscInt v_target, const PetscInt g_target, PetscReal shift, PetscBool use_uniform_particle_grid) 481 { 482 DM pack, *globSwarmArray, grid_dm[LANDAU_MAX_GRIDS]; 483 Mat *globMpArray, g_Mass[LANDAU_MAX_GRIDS]; 484 KSP t_ksp[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS]; 485 Vec t_fhat[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS]; 486 PetscInt nDMs; 487 PetscErrorCode ierr = (PetscErrorCode)0; // used for inside thread loops 488 #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY) 489 PetscInt numthreads = PetscNumOMPThreads; 490 #else 491 PetscInt numthreads = 1; 492 #endif 493 LandauCtx *ctx; 494 Vec *globXArray; 495 PetscReal moments_0[5], moments_1a[5], moments_1b[5], dt_init; 496 PrintCtx *printCtx; 497 498 PetscFunctionBeginUser; 499 PetscCheck(numthreads <= EX30_MAX_NUM_THRDS, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Too many threads %" PetscInt_FMT " > %d", numthreads, EX30_MAX_NUM_THRDS); 500 PetscCheck(numthreads > 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number threads %" PetscInt_FMT " > %d", numthreads, EX30_MAX_NUM_THRDS); 501 PetscCall(TSGetDM(ts, &pack)); 502 PetscCall(DMGetApplicationContext(pack, &ctx)); 503 PetscCheck(ctx->batch_sz % numthreads == 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "batch size (-dm_landau_batch_size) %" PetscInt_FMT " mod #threads %" PetscInt_FMT " must equal zero", ctx->batch_sz, numthreads); 504 PetscCall(DMCompositeGetNumberDM(pack, &nDMs)); // number of vertices * number of grids 505 PetscCall(PetscInfo(pack, "Have %" PetscInt_FMT " total grids, with %" PetscInt_FMT " Landau local batched and %" PetscInt_FMT " global items (vertices) %d DMs\n", ctx->num_grids, ctx->batch_sz, num_vertices, (int)nDMs)); 506 PetscCall(PetscMalloc(sizeof(*globXArray) * nDMs, &globXArray)); 507 PetscCall(PetscMalloc(sizeof(*globMpArray) * nDMs, &globMpArray)); 508 PetscCall(PetscMalloc(sizeof(*globSwarmArray) * nDMs, &globSwarmArray)); 509 // print ctx 510 PetscCall(PetscNew(&printCtx)); 511 PetscCall(TSSetApplicationContext(ts, printCtx)); 512 printCtx->v_target = v_target; 513 printCtx->g_target = g_target; 514 printCtx->ctx = ctx; 515 printCtx->globSwarmArray = globSwarmArray; 516 printCtx->grid_dm = grid_dm; 517 printCtx->globMpArray = globMpArray; 518 printCtx->g_Mass = g_Mass; 519 printCtx->globXArray = globXArray; 520 printCtx->print_entropy = PETSC_FALSE; 521 PetscOptionsBegin(PETSC_COMM_SELF, "", "Print Options", "DMPLEX"); 522 PetscCall(PetscOptionsBool("-print_entropy", "Print entropy and moments at each time step", "ex30.c", printCtx->print_entropy, &printCtx->print_entropy, NULL)); 523 PetscOptionsEnd(); 524 // view 525 PetscCall(DMViewFromOptions(ctx->plex[g_target], NULL, "-ex30_dm_view")); 526 if (ctx->num_grids > g_target + 1) { PetscCall(DMViewFromOptions(ctx->plex[g_target + 1], NULL, "-ex30_dm_view2")); } 527 // create mesh mass matrices 528 PetscCall(VecZeroEntries(X)); 529 PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray)); // just to duplicate 530 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids 531 Vec subX = globXArray[LAND_PACK_IDX(0, grid)]; 532 DM dm = ctx->plex[grid]; 533 PetscSection s; 534 grid_dm[grid] = dm; 535 PetscCall(DMCreateMassMatrix(dm, dm, &g_Mass[grid])); 536 // 537 PetscCall(DMGetLocalSection(dm, &s)); 538 PetscCall(DMPlexCreateClosureIndex(dm, s)); 539 for (int tid = 0; tid < numthreads; tid++) { 540 PC pc; 541 PetscCall(VecDuplicate(subX, &t_fhat[grid][tid])); 542 PetscCall(KSPCreate(PETSC_COMM_SELF, &t_ksp[grid][tid])); 543 PetscCall(KSPSetType(t_ksp[grid][tid], KSPCG)); 544 PetscCall(KSPGetPC(t_ksp[grid][tid], &pc)); 545 PetscCall(PCSetType(pc, PCJACOBI)); 546 PetscCall(KSPSetOptionsPrefix(t_ksp[grid][tid], "ptof_")); 547 PetscCall(KSPSetOperators(t_ksp[grid][tid], g_Mass[grid], g_Mass[grid])); 548 PetscCall(KSPSetFromOptions(t_ksp[grid][tid])); 549 } 550 } 551 PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray)); 552 PetscCall(TSGetTimeStep(ts, &dt_init)); // we could have an adaptive time stepper 553 // loop over all vertices in chucks that are batched for TSSolve 554 for (int i = 0; i < 5; i++) moments_0[i] = moments_1a[i] = moments_1b[i] = 0; 555 for (PetscInt global_vertex_id_0 = 0; global_vertex_id_0 < num_vertices; global_vertex_id_0 += ctx->batch_sz, shift /= 2) { // outer vertex loop 556 PetscCall(TSSetTime(ts, 0)); 557 PetscCall(TSSetStepNumber(ts, 0)); 558 PetscCall(TSSetTimeStep(ts, dt_init)); 559 PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray)); 560 printCtx->global_vertex_id_0 = global_vertex_id_0; 561 if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) { 562 PetscCall(PetscObjectSetName((PetscObject)globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], "rho")); 563 printCtx->print = PETSC_TRUE; 564 } else printCtx->print = PETSC_FALSE; 565 // create fake particles in batches with threads 566 for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) { 567 PetscReal *xx_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *yy_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *zz_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *wp_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS] /* , radiuses[80000] */; 568 PetscInt Np_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS]; 569 // make particles 570 for (int tid = 0; tid < numthreads; tid++) { 571 const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id; 572 if (glb_v_id < num_vertices) { // the ragged edge (in last batch) 573 PetscInt Npp0 = a_Np + (glb_v_id % (a_Np / 10 + 1)), nTargetP[LANDAU_MAX_GRIDS]; // n of particels in each dim with load imbalance 574 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids 575 // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) { 576 const PetscReal kT_m = ctx->k * ctx->thermal_temps[ctx->species_offset[grid]] / ctx->masses[ctx->species_offset[grid]] / (ctx->v_0 * ctx->v_0); /* theta = 2kT/mc^2 per species */ 577 PetscReal lo[3] = {-ctx->radius[grid], -ctx->radius[grid], -ctx->radius[grid]}, hi[3] = {ctx->radius[grid], ctx->radius[grid], ctx->radius[grid]}, hp[3], vole; // would be nice to get box from DM 578 PetscInt Npi = Npp0, Npj = 2 * Npp0, Npk = 1; 579 PetscRandom rand; 580 PetscReal sigma = ctx->thermal_speed[grid] / ctx->thermal_speed[0], p2_shift = grid == 0 ? shift : -shift; // symmetric shift of e vs ions 581 PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand)); 582 PetscCall(PetscRandomSetInterval(rand, 0., 1.)); 583 PetscCall(PetscRandomSetFromOptions(rand)); 584 if (dim == 2) lo[0] = 0; // Landau coordinate (r,z) 585 else Npi = Npj = Npk = Npp0; 586 // User: use glb_v_id to index into your data 587 const PetscInt NNreal = Npi * Npj * Npk, NN = NNreal + (dim == 2 ? 3 : 6); // make room for bounding box 588 Np_t[grid][tid] = NN; 589 if (glb_v_id == v_target) nTargetP[grid] = NN; 590 PetscCall(PetscMalloc4(NN, &xx_t[grid][tid], NN, &yy_t[grid][tid], NN, &wp_t[grid][tid], dim == 2 ? 1 : NN, &zz_t[grid][tid])); 591 hp[0] = (hi[0] - lo[0]) / Npi; 592 hp[1] = (hi[1] - lo[1]) / Npj; 593 hp[2] = (hi[2] - lo[2]) / Npk; 594 if (dim == 2) hp[2] = 1; 595 PetscCall(PetscInfo(pack, " lo = %14.7e, hi = %14.7e; hp = %14.7e, %14.7e; kT_m = %g; \n", (double)lo[1], (double)hi[1], (double)hp[0], (double)hp[1], (double)kT_m)); // temp 596 vole = hp[0] * hp[1] * hp[2] * ctx->n[grid]; // fix for multi-species 597 PetscCall(PetscInfo(pack, "Vertex %" PetscInt_FMT ", grid %" PetscInt_FMT " with %" PetscInt_FMT " particles (diagnostic target = %" PetscInt_FMT ")\n", glb_v_id, grid, NN, v_target)); 598 for (int pj = 0, pp = 0; pj < Npj; pj++) { 599 for (int pk = 0; pk < Npk; pk++) { 600 for (int pi = 0; pi < Npi; pi++, pp++) { 601 PetscReal p_shift = p2_shift; 602 wp_t[grid][tid][pp] = 0; 603 if (use_uniform_particle_grid) { 604 xx_t[grid][tid][pp] = lo[0] + hp[0] / 2.0 + pi * hp[0]; 605 yy_t[grid][tid][pp] = lo[1] + hp[1] / 2.0 + pj * hp[1]; 606 if (dim == 3) zz_t[grid][tid][pp] = lo[2] + hp[2] / 2.0 + pk * hp[2]; 607 PetscReal x[] = {xx_t[grid][tid][pp], yy_t[grid][tid][pp], dim == 2 ? 0 : zz_t[grid][tid][pp]}; 608 p_shift *= ctx->thermal_speed[grid] / ctx->v_0; 609 if (ctx->sphere && PetscSqrtReal(PetscSqr(xx_t[grid][tid][pp]) + PetscSqr(yy_t[grid][tid][pp])) > 0.92 * hi[0]) { 610 wp_t[grid][tid][pp] = 0; 611 } else { 612 maxwellian(dim, x, kT_m, vole, p_shift, &wp_t[grid][tid][pp]); 613 if (ctx->num_grids == 1 && shift != 0) { // bi-maxwellian, electron plasma 614 maxwellian(dim, x, kT_m, vole, -p_shift, &wp_t[grid][tid][pp]); // symmetric shift of electron plasma 615 } 616 } 617 } else { 618 PetscReal u1, u2; 619 do { 620 do { 621 PetscCall(PetscRandomGetValueReal(rand, &u1)); 622 } while (u1 == 0); 623 PetscCall(PetscRandomGetValueReal(rand, &u2)); 624 //compute z0 and z1 625 PetscReal mag = sigma * PetscSqrtReal(-2.0 * PetscLogReal(u1)); // is this the same scale grid Maxwellian? t_therm = sigma 626 xx_t[grid][tid][pp] = mag * PetscCosReal(2.0 * PETSC_PI * u2); 627 yy_t[grid][tid][pp] = mag * PetscSinReal(2.0 * PETSC_PI * u2); 628 if (dim == 2 && xx_t[grid][tid][pp] < lo[0]) xx_t[grid][tid][pp] = -xx_t[grid][tid][pp]; 629 if (dim == 3) zz_t[grid][tid][pp] = lo[2] + hp[2] / 2.0 + pk * hp[2]; 630 if (!ctx->sphere) { 631 if (dim == 2 && xx_t[grid][tid][pp] < 0) xx_t[grid][tid][pp] = -xx_t[grid][tid][pp]; // ??? 632 else if (dim == 3) { 633 while (zz_t[grid][tid][pp] >= hi[2] || zz_t[grid][tid][pp] <= lo[2]) zz_t[grid][tid][pp] *= .9; 634 } 635 while (xx_t[grid][tid][pp] >= hi[0] || xx_t[grid][tid][pp] <= lo[0]) xx_t[grid][tid][pp] *= .9; 636 while (yy_t[grid][tid][pp] >= hi[1] || yy_t[grid][tid][pp] <= lo[1]) yy_t[grid][tid][pp] *= .9; 637 } else { // 2D 638 //if (glb_v_id == v_target && pp < 80000) radiuses[pp] = PetscSqrtReal(PetscSqr(xx_t[grid][tid][pp]) + PetscSqr(yy_t[grid][tid][pp])); 639 while (PetscSqrtReal(PetscSqr(xx_t[grid][tid][pp]) + PetscSqr(yy_t[grid][tid][pp])) > 0.92 * hi[0]) { // safety factor for facets of sphere 640 xx_t[grid][tid][pp] *= .9; 641 yy_t[grid][tid][pp] *= .9; 642 } 643 } 644 if (ctx->num_grids == 1 && pp % 2 == 0) p_shift = 0; // one species, split bi-max 645 p_shift *= ctx->thermal_speed[grid] / ctx->v_0; 646 if (dim == 3) zz_t[grid][tid][pp] += p_shift; 647 else yy_t[grid][tid][pp] += p_shift; 648 wp_t[grid][tid][pp] += ctx->n[grid] / NNreal * PetscSqrtReal(ctx->masses[ctx->species_offset[grid]] / ctx->masses[0]); 649 if (p_shift <= 0) break; // add bi-max for electron plasma only 650 p_shift = -p_shift; 651 } while (ctx->num_grids == 1); // add bi-max for electron plasma only 652 } 653 { 654 if (glb_v_id == v_target) { 655 PetscReal x[] = {xx_t[grid][tid][pp], yy_t[grid][tid][pp], dim == 2 ? 0 : zz_t[grid][tid][pp]}; 656 PetscReal v2 = 0, fact = dim == 2 ? 2.0 * PETSC_PI * x[0] : 1, w = fact * wp_t[grid][tid][pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]]; 657 for (int i = 0; i < dim; ++i) v2 += PetscSqr(x[i]); 658 moments_0[0] += w; // not thread safe 659 moments_0[1] += w * ctx->v_0 * x[1]; // z-momentum 660 moments_0[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2; 661 } 662 } 663 } 664 } 665 } 666 if (dim == 2) { // fix bounding box 667 int pp = NNreal; 668 wp_t[grid][tid][pp] = 0; 669 xx_t[grid][tid][pp] = 1.e-7; 670 yy_t[grid][tid][pp++] = hi[1] - 5.e-7; 671 wp_t[grid][tid][pp] = 0; 672 xx_t[grid][tid][pp] = hi[0] - 5.e-7; 673 yy_t[grid][tid][pp++] = 0; 674 wp_t[grid][tid][pp] = 0; 675 xx_t[grid][tid][pp] = 1.e-7; 676 yy_t[grid][tid][pp++] = lo[1] + 5.e-7; 677 } else { 678 const int p0 = NNreal; 679 for (int pj = 0; pj < 6; pj++) { xx_t[grid][tid][p0 + pj] = yy_t[grid][tid][p0 + pj] = zz_t[grid][tid][p0 + pj] = wp_t[grid][tid][p0 + pj] = 0; } 680 xx_t[grid][tid][p0 + 0] = lo[0]; 681 xx_t[grid][tid][p0 + 1] = hi[0]; 682 yy_t[grid][tid][p0 + 2] = lo[1]; 683 yy_t[grid][tid][p0 + 3] = hi[1]; 684 zz_t[grid][tid][p0 + 4] = lo[2]; 685 zz_t[grid][tid][p0 + 5] = hi[2]; 686 } 687 PetscCall(PetscRandomDestroy(&rand)); 688 } 689 // entropy init, need global n 690 if (glb_v_id == v_target) { 691 const PetscReal N_inv = 1 / moments_0[0]; 692 PetscCall(PetscInfo(pack, "Target %" PetscInt_FMT " with %" PetscInt_FMT " particels\n", glb_v_id, nTargetP[0])); 693 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 694 const PetscInt NN = nTargetP[grid]; 695 for (int pp = 0; pp < NN; pp++) { 696 const PetscReal fact = dim == 2 ? 2.0 * PETSC_PI * xx_t[grid][tid][pp] : 1, w = fact * ctx->n_0 * ctx->masses[ctx->species_offset[grid]] * wp_t[grid][tid][pp], ww = w * N_inv; 697 if (w > PETSC_REAL_MIN) { 698 moments_0[3] -= ww * PetscLogReal(ww); 699 PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "ww (%g) > 1", (double)ww); 700 } else moments_0[4] -= w; 701 } 702 } // grid 703 } // target 704 } // active 705 } // threads 706 /* Create particle swarm */ 707 for (int tid = 0; tid < numthreads; tid++) { 708 const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id; 709 if (glb_v_id < num_vertices) { // the ragged edge of the last batch 710 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids 711 PetscSection section; 712 PetscInt Nf; 713 DM dm = grid_dm[grid]; 714 PetscCall(DMGetLocalSection(dm, §ion)); 715 PetscCall(PetscSectionGetNumFields(section, &Nf)); 716 PetscCheck(Nf == 1, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Only on species per grid supported -- todo"); 717 PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 718 PetscCall(PetscInfo(pack, "call createSwarm [%" PetscInt_FMT ".%" PetscInt_FMT "] local block index %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid))); 719 PetscCall(createSwarm(dm, dim, &globSwarmArray[LAND_PACK_IDX(v_id, grid)])); 720 } 721 } // active 722 } // threads 723 PetscCheck(ierr != 9999, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Only support one species per grid"); 724 // make globMpArray 725 PetscPragmaOMP(parallel for) 726 for (int tid = 0; tid < numthreads; tid++) { 727 const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id; 728 if (glb_v_id < num_vertices) { 729 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids 730 // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) -- loop over species for Nf > 1 -- TODO 731 PetscErrorCode ierr_t; 732 DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)]; 733 ierr_t = PetscInfo(pack, "makeSwarm %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid)); 734 ierr_t = makeSwarm(sw, dim, Np_t[grid][tid], xx_t[grid][tid], yy_t[grid][tid], zz_t[grid][tid]); 735 if (ierr_t) ierr = ierr_t; 736 } 737 } 738 } 739 for (int tid = 0; tid < numthreads; tid++) { 740 const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id; 741 if (glb_v_id < num_vertices) { 742 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids 743 DM dm = grid_dm[grid]; 744 DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)]; 745 PetscCall(PetscInfo(pack, "createMp %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid))); 746 PetscCall(createMp(dm, sw, &globMpArray[LAND_PACK_IDX(v_id, grid)])); 747 PetscCall(MatViewFromOptions(globMpArray[LAND_PACK_IDX(v_id, grid)], NULL, "-mp_mat_view")); 748 } 749 } 750 } 751 // p --> g: set X 752 // PetscPragmaOMP(parallel for) 753 for (int tid = 0; tid < numthreads; tid++) { 754 const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id; 755 if (glb_v_id < num_vertices) { 756 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids 757 PetscErrorCode ierr_t; 758 DM dm = grid_dm[grid]; 759 DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)]; 760 Vec subX = globXArray[LAND_PACK_IDX(v_id, grid)], work = t_fhat[grid][tid]; 761 ierr_t = PetscInfo(pack, "particlesToGrid %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid)); 762 ierr_t = particlesToGrid(dm, sw, tid, dim, wp_t[grid][tid], subX, globMpArray[LAND_PACK_IDX(v_id, grid)]); 763 if (ierr_t) ierr = ierr_t; 764 // u = M^_1 f_w 765 ierr_t = VecCopy(subX, work); 766 ierr_t = KSPSolve(t_ksp[grid][tid], work, subX); 767 if (ierr_t) ierr = ierr_t; 768 } 769 } 770 } 771 PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in OMP loop. ierr = %d", (int)ierr); 772 /* Cleanup */ 773 for (int tid = 0; tid < numthreads; tid++) { 774 const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id; 775 if (glb_v_id < num_vertices) { 776 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids 777 PetscCall(PetscFree4(xx_t[grid][tid], yy_t[grid][tid], wp_t[grid][tid], zz_t[grid][tid])); 778 } 779 } // active 780 } // threads 781 } // (fake) particle loop 782 // standard view of initial conditions 783 if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) { 784 PetscCall(DMSetOutputSequenceNumber(ctx->plex[g_target], 0, 0.0)); 785 PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], NULL, "-ex30_vec_view")); 786 if (ctx->num_grids > g_target + 1) { 787 PetscCall(DMSetOutputSequenceNumber(ctx->plex[g_target + 1], 0, 0.0)); 788 PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target + 1)], NULL, "-ex30_vec_view2")); 789 } 790 PetscCall(MatViewFromOptions(globMpArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], NULL, "-ex30_mass_mat_view")); 791 PetscCall(DMViewFromOptions(globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], NULL, "-ex30_sw_view")); 792 PetscCall(DMSwarmViewXDMF(globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], "initial_swarm.xmf")); // writes a file by default!!! 793 } 794 // coarse graining moments_1a, bring f back from grid before advance 795 if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz && printCtx->print_entropy) { 796 const PetscInt v_id = v_target % ctx->batch_sz; 797 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 798 PetscDataType dtype; 799 PetscReal *wp, *coords; 800 DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)]; 801 Vec work, subX = globXArray[LAND_PACK_IDX(v_id, grid)]; 802 PetscInt bs, NN; 803 // C-G moments 804 PetscCall(VecDuplicate(subX, &work)); 805 PetscCall(gridToParticles(grid_dm[grid], sw, subX, work, globMpArray[LAND_PACK_IDX(v_id, grid)], g_Mass[grid])); 806 PetscCall(VecDestroy(&work)); 807 // moments 808 PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords)); 809 PetscCall(DMSwarmGetLocalSize(sw, &NN)); 810 PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); 811 for (int pp = 0; pp < NN; pp++) { 812 PetscReal v2 = 0, fact = (dim == 2) ? 2.0 * PETSC_PI * coords[pp * dim + 0] : 1, w = fact * wp[pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]]; 813 for (int i = 0; i < dim; ++i) v2 += PetscSqr(coords[pp * dim + i]); 814 moments_1a[0] += w; 815 moments_1a[1] += w * ctx->v_0 * coords[pp * dim + 1]; // z-momentum 816 moments_1a[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2; 817 } 818 PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords)); 819 PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp)); 820 } 821 // entropy 822 const PetscReal N_inv = 1 / moments_1a[0]; 823 PetscCall(PetscInfo(pack, "Entropy batch %" PetscInt_FMT " of %" PetscInt_FMT ", n = %g\n", v_target, num_vertices, (double)(1 / N_inv))); 824 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 825 PetscDataType dtype; 826 PetscReal *wp, *coords; 827 DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)]; 828 PetscInt bs, NN; 829 PetscCall(DMSwarmGetLocalSize(sw, &NN)); 830 PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); 831 PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords)); 832 for (int pp = 0; pp < NN; pp++) { 833 PetscReal fact = (dim == 2) ? 2.0 * PETSC_PI * coords[pp * dim + 0] : 1, w = fact * wp[pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]], ww = w * N_inv; 834 if (w > PETSC_REAL_MIN) { 835 moments_1a[3] -= ww * PetscLogReal(ww); 836 PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "ww (%g) > 1", (double)ww); 837 } else moments_1a[4] -= w; 838 } 839 PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp)); 840 PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords)); 841 } 842 } 843 // restore vector 844 PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray)); 845 // view initial grid 846 if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) { PetscCall(DMPlexLandauPrintNorms(X, 0)); } 847 // advance 848 PetscCall(TSSetSolution(ts, X)); 849 PetscCall(PetscInfo(pack, "Advance vertex %" PetscInt_FMT " to %" PetscInt_FMT "\n", global_vertex_id_0, global_vertex_id_0 + ctx->batch_sz)); 850 PetscCall(TSSetPostStep(ts, PostStep)); 851 PetscCall(PostStep(ts)); 852 PetscCall(TSSolve(ts, X)); 853 // view 854 PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray)); 855 if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) { 856 /* Visualize original particle field */ 857 DM sw = globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)]; 858 Vec f; 859 PetscCall(DMSetOutputSequenceNumber(sw, 0, 0.0)); 860 PetscCall(DMViewFromOptions(grid_dm[g_target], NULL, "-weights_dm_view")); 861 PetscCall(DMViewFromOptions(sw, NULL, "-weights_sw_view")); 862 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 863 PetscCall(PetscObjectSetName((PetscObject)f, "weights")); 864 PetscCall(VecViewFromOptions(f, NULL, "-weights_vec_view")); 865 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 866 // 867 PetscCall(DMPlexLandauPrintNorms(X, 1)); 868 } 869 if (!use_uniform_particle_grid) { // resample to uniform grid 870 for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) { 871 PetscReal *xx_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *yy_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *zz_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *wp_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS]; 872 PetscInt Np_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS]; 873 for (int tid = 0; tid < numthreads; tid++) { 874 const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id; 875 if (glb_v_id < num_vertices) { 876 // create uniform grid w/o weights & smaller 877 PetscInt Npp0 = (a_Np + (glb_v_id % (a_Np / 10 + 1))) / 2, Nv; // 1/2 of uniform particle grid size 878 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 879 // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) 880 PetscReal lo[3] = {-ctx->radius[grid], -ctx->radius[grid], -ctx->radius[grid]}, hi[3] = {ctx->radius[grid], ctx->radius[grid], ctx->radius[grid]}, hp[3]; 881 PetscInt Npi = Npp0, Npj = 2 * Npp0, Npk = 1, NN; 882 // delete old particles and particle mass matrix 883 PetscCall(DMDestroy(&globSwarmArray[LAND_PACK_IDX(v_id, grid)])); 884 PetscCall(MatDestroy(&globMpArray[LAND_PACK_IDX(v_id, grid)])); 885 // create fake particles in batches with threads 886 PetscCall(MatGetLocalSize(g_Mass[grid], &Nv, NULL)); 887 if (dim == 2) lo[0] = 0; 888 else Npi = Npj = Npk = Npp0; 889 NN = Npi * Npj * Npk + (dim == 2 ? 3 : 6); // make a regular grid of particles Npp x Npp 890 while (Npi * Npj * Npk < Nv) { // make stable - no LS 891 Npi++; 892 Npj++; 893 Npk++; 894 NN = Npi * Npj * Npk + (dim == 2 ? 3 : 6); 895 } 896 Np_t[grid][tid] = NN; 897 PetscCall(PetscMalloc4(NN, &xx_t[grid][tid], NN, &yy_t[grid][tid], NN, &wp_t[grid][tid], dim == 2 ? 1 : NN, &zz_t[grid][tid])); 898 hp[0] = (hi[0] - lo[0]) / Npi; 899 hp[1] = (hi[1] - lo[1]) / Npj; 900 hp[2] = (hi[2] - lo[2]) / Npk; 901 if (dim == 2) hp[2] = 1; 902 PetscCall(PetscInfo(pack, "Resampling %d particles, %d vertices\n", (int)NN, (int)Nv)); // temp 903 for (int pj = 0, pp = 0; pj < Npj; pj++) { 904 for (int pk = 0; pk < Npk; pk++) { 905 for (int pi = 0; pi < Npi; pi++, pp++) { 906 wp_t[grid][tid][pp] = 0; 907 xx_t[grid][tid][pp] = lo[0] + hp[0] / 2.0 + pi * hp[0]; 908 yy_t[grid][tid][pp] = lo[1] + hp[1] / 2.0 + pj * hp[1]; 909 if (dim == 3) zz_t[grid][tid][pp] = lo[2] + hp[2] / 2.0 + pk * hp[2]; 910 } 911 } 912 } 913 if (dim == 2) { // fix bounding box 914 int pp = NN - 3; 915 wp_t[grid][tid][pp] = 0; 916 xx_t[grid][tid][pp] = 1.e-7; 917 yy_t[grid][tid][pp++] = hi[1] - 5.e-7; 918 wp_t[grid][tid][pp] = 0; 919 xx_t[grid][tid][pp] = hi[0] - 5.e-7; 920 yy_t[grid][tid][pp++] = 0; 921 wp_t[grid][tid][pp] = 0; 922 xx_t[grid][tid][pp] = 1.e-7; 923 yy_t[grid][tid][pp++] = lo[1] + 5.e-7; 924 } else { 925 const int p0 = NN - 6; 926 for (int pj = 0; pj < 6; pj++) { xx_t[grid][tid][p0 + pj] = yy_t[grid][tid][p0 + pj] = zz_t[grid][tid][p0 + pj] = wp_t[grid][tid][p0 + pj] = 0; } 927 xx_t[grid][tid][p0 + 0] = lo[0]; 928 xx_t[grid][tid][p0 + 1] = hi[0]; 929 yy_t[grid][tid][p0 + 2] = lo[1]; 930 yy_t[grid][tid][p0 + 3] = hi[1]; 931 zz_t[grid][tid][p0 + 4] = lo[2]; 932 zz_t[grid][tid][p0 + 5] = hi[2]; 933 } 934 } 935 } // active 936 } // threads 937 /* Create particle swarm */ 938 for (int tid = 0; tid < numthreads; tid++) { 939 const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id; 940 if (glb_v_id < num_vertices) { // the ragged edge of the last batch 941 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids 942 // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) -- loop over species for Nf > 1 -- TODO 943 PetscErrorCode ierr_t; 944 PetscSection section; 945 PetscInt Nf; 946 DM dm = grid_dm[grid]; 947 ierr_t = DMGetLocalSection(dm, §ion); 948 ierr_t = PetscSectionGetNumFields(section, &Nf); 949 if (Nf != 1) ierr_t = (PetscErrorCode)9999; 950 else { 951 ierr_t = DMViewFromOptions(dm, NULL, "-dm_view"); 952 ierr_t = PetscInfo(pack, "call createSwarm [%" PetscInt_FMT ".%" PetscInt_FMT "] local block index %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid)); 953 ierr_t = createSwarm(dm, dim, &globSwarmArray[LAND_PACK_IDX(v_id, grid)]); 954 } 955 if (ierr_t) ierr = ierr_t; 956 } 957 } // active 958 } // threads 959 PetscCheck(ierr != 9999, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Only support one species per grid"); 960 PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Error in OMP loop. ierr = %d", (int)ierr); 961 // make globMpArray 962 PetscPragmaOMP(parallel for) 963 for (int tid = 0; tid < numthreads; tid++) { 964 const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id; 965 if (glb_v_id < num_vertices) { 966 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids 967 // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) -- loop over species for Nf > 1 -- TODO 968 PetscErrorCode ierr_t; 969 DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)]; 970 ierr_t = PetscInfo(pack, "makeSwarm %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid)); 971 ierr_t = makeSwarm(sw, dim, Np_t[grid][tid], xx_t[grid][tid], yy_t[grid][tid], zz_t[grid][tid]); 972 if (ierr_t) ierr = ierr_t; 973 } 974 } // active 975 } // threads 976 // create particle mass matrices 977 //PetscPragmaOMP(parallel for) 978 for (int tid = 0; tid < numthreads; tid++) { 979 const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id; 980 if (glb_v_id < num_vertices) { 981 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids 982 PetscErrorCode ierr_t; 983 DM dm = grid_dm[grid]; 984 DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)]; 985 ierr_t = PetscInfo(pack, "createMp %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid)); 986 ierr_t = createMp(dm, sw, &globMpArray[LAND_PACK_IDX(v_id, grid)]); 987 if (ierr_t) ierr = ierr_t; 988 } 989 } // active 990 } // threads 991 PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in OMP loop. ierr = %d", (int)ierr); 992 /* Cleanup */ 993 for (int tid = 0; tid < numthreads; tid++) { 994 const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id; 995 if (glb_v_id < num_vertices) { 996 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids 997 PetscCall(PetscFree4(xx_t[grid][tid], yy_t[grid][tid], wp_t[grid][tid], zz_t[grid][tid])); 998 } 999 } // active 1000 } // threads 1001 } // batch 1002 // view 1003 if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) { 1004 /* Visualize particle field */ 1005 DM sw = globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)]; 1006 Vec f; 1007 PetscCall(DMSetOutputSequenceNumber(sw, 0, 0.0)); 1008 PetscCall(DMViewFromOptions(sw, NULL, "-resampled_weights_sw_view")); 1009 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 1010 PetscCall(PetscObjectSetName((PetscObject)f, "resampled_weights")); 1011 PetscCall(VecViewFromOptions(f, NULL, "-resampled_weights_vec_view")); 1012 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 1013 PetscCall(DMSwarmViewXDMF(sw, "resampled.xmf")); 1014 } 1015 } // !uniform 1016 // particles to grid, compute moments and entropy, for target vertex only 1017 if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz && printCtx->print_entropy) { 1018 PetscReal energy_error_rel; 1019 PetscCall(gridToParticles_private(grid_dm, globSwarmArray, dim, v_target, numthreads, num_vertices, global_vertex_id_0, globMpArray, g_Mass, t_fhat, moments_1b, globXArray, ctx)); 1020 energy_error_rel = PetscAbsReal((moments_1b[2] - moments_0[2])) / moments_0[2]; 1021 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Particle Moments:\t number density momentum (par) energy entropy negative weights : # OMP threads %g\n", (double)numthreads)); 1022 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\tInitial: %18.12e %19.12e %18.12e %18.12e %g %%\n", (double)moments_0[0], (double)moments_0[1], (double)moments_0[2], (double)moments_0[3], 100 * (double)(moments_0[4] / moments_0[0]))); 1023 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\tCoarse-graining: %18.12e %19.12e %18.12e %18.12e %g %%\n", (double)moments_1a[0], (double)moments_1a[1], (double)moments_1a[2], (double)moments_1a[3], 100 * (double)(moments_1a[4] / moments_0[0]))); 1024 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\tLandau: %18.12e %19.12e %18.12e %18.12e %g %%\n", (double)moments_1b[0], (double)moments_1b[1], (double)moments_1b[2], (double)moments_1b[3], 100 * (double)(moments_1b[4] / moments_0[0]))); 1025 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Coarse-graining entropy generation = %e ; Landau entropy generation = %e\n", (double)(moments_1a[3] - moments_0[3]), (double)(moments_1b[3] - moments_0[3]))); 1026 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "(relative) energy conservation: Coarse-graining = %e, Landau = %e (%g %d)\n", (double)(moments_1a[2] - moments_0[2]) / (double)moments_0[2], (double)energy_error_rel, (double)PetscLog10Real(energy_error_rel), (int)(PetscLog10Real(energy_error_rel) + .5))); 1027 } 1028 // restore vector 1029 PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray)); 1030 // cleanup 1031 for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) { 1032 for (int tid = 0; tid < numthreads; tid++) { 1033 const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id; 1034 if (glb_v_id < num_vertices) { 1035 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 1036 PetscCall(DMDestroy(&globSwarmArray[LAND_PACK_IDX(v_id, grid)])); 1037 PetscCall(MatDestroy(&globMpArray[LAND_PACK_IDX(v_id, grid)])); 1038 } 1039 } 1040 } 1041 } 1042 } // user batch, not used 1043 /* Cleanup */ 1044 PetscCall(PetscFree(globXArray)); 1045 PetscCall(PetscFree(globSwarmArray)); 1046 PetscCall(PetscFree(globMpArray)); 1047 PetscCall(PetscFree(printCtx)); 1048 // clean up mass matrices 1049 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids 1050 PetscCall(MatDestroy(&g_Mass[grid])); 1051 for (int tid = 0; tid < numthreads; tid++) { 1052 PetscCall(VecDestroy(&t_fhat[grid][tid])); 1053 PetscCall(KSPDestroy(&t_ksp[grid][tid])); 1054 } 1055 } 1056 PetscFunctionReturn(PETSC_SUCCESS); 1057 } 1058 1059 int main(int argc, char **argv) 1060 { 1061 DM pack; 1062 Vec X; 1063 PetscInt dim = 2, num_vertices = 1, Np = 10, v_target = 0, g_target = 0; 1064 TS ts; 1065 Mat J; 1066 LandauCtx *ctx; 1067 PetscReal shift = 0; 1068 PetscBool use_uniform_particle_grid = PETSC_TRUE; 1069 1070 PetscFunctionBeginUser; 1071 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1072 // process args 1073 PetscOptionsBegin(PETSC_COMM_SELF, "", "Collision Options", "DMPLEX"); 1074 PetscCall(PetscOptionsInt("-dim", "Velocity space dimension", "ex30.c", dim, &dim, NULL)); 1075 PetscCall(PetscOptionsInt("-number_spatial_vertices", "Number of user spatial vertices to be batched for Landau", "ex30.c", num_vertices, &num_vertices, NULL)); 1076 PetscCall(PetscOptionsInt("-number_particles_per_dimension", "Number of particles per grid, with slight modification per spatial vertex, in each dimension of base Cartesian grid", "ex30.c", Np, &Np, NULL)); 1077 PetscCall(PetscOptionsBool("-use_uniform_particle_grid", "Use uniform particle grid", "ex30.c", use_uniform_particle_grid, &use_uniform_particle_grid, NULL)); 1078 PetscCall(PetscOptionsInt("-vertex_view_target", "Global vertex for diagnostics", "ex30.c", v_target, &v_target, NULL)); 1079 PetscCall(PetscOptionsReal("-e_shift", "Bi-Maxwellian shift", "ex30.c", shift, &shift, NULL)); 1080 PetscCheck(v_target < num_vertices, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Batch to view %" PetscInt_FMT " should be < number of vertices %" PetscInt_FMT, v_target, num_vertices); 1081 PetscCall(PetscOptionsInt("-grid_view_target", "Grid to view with diagnostics", "ex30.c", g_target, &g_target, NULL)); 1082 PetscOptionsEnd(); 1083 /* Create a mesh */ 1084 PetscCall(DMPlexLandauCreateVelocitySpace(PETSC_COMM_SELF, dim, "", &X, &J, &pack)); 1085 PetscCall(DMGetApplicationContext(pack, &ctx)); 1086 PetscCall(DMSetUp(pack)); 1087 PetscCall(DMSetOutputSequenceNumber(pack, 0, 0.0)); 1088 PetscCheck(g_target < ctx->num_grids, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Grid to view %" PetscInt_FMT " should be < number of grids %" PetscInt_FMT, g_target, ctx->num_grids); 1089 PetscCheck(ctx->batch_view_idx == v_target % ctx->batch_sz, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Global view index %" PetscInt_FMT " mode batch size %" PetscInt_FMT " != ctx->batch_view_idx %" PetscInt_FMT, v_target, ctx->batch_sz, ctx->batch_view_idx); 1090 // PetscCheck(!use_uniform_particle_grid || !ctx->sphere, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Can not use -use_uniform_particle_grid and -dm_landau_sphere"); 1091 /* Create timestepping solver context */ 1092 PetscCall(TSCreate(PETSC_COMM_SELF, &ts)); 1093 PetscCall(TSSetDM(ts, pack)); 1094 PetscCall(TSSetIFunction(ts, NULL, DMPlexLandauIFunction, NULL)); 1095 PetscCall(TSSetIJacobian(ts, J, J, DMPlexLandauIJacobian, NULL)); 1096 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 1097 PetscCall(TSSetFromOptions(ts)); 1098 PetscCall(PetscObjectSetName((PetscObject)X, "X")); 1099 // do particle advance 1100 PetscCall(go(ts, X, num_vertices, Np, dim, v_target, g_target, shift, use_uniform_particle_grid)); 1101 PetscCall(MatZeroEntries(J)); // need to zero out so as to not reuse it in Landau's logic 1102 /* clean up */ 1103 PetscCall(DMPlexLandauDestroyVelocitySpace(&pack)); 1104 PetscCall(TSDestroy(&ts)); 1105 PetscCall(VecDestroy(&X)); 1106 PetscCall(PetscFinalize()); 1107 return 0; 1108 } 1109 1110 /*TEST 1111 1112 build: 1113 requires: !complex 1114 1115 testset: 1116 requires: double defined(PETSC_USE_DMLANDAU_2D) 1117 output_file: output/ex30_0.out 1118 args: -dim 2 -petscspace_degree 3 -dm_landau_num_species_grid 1,1,1 -dm_refine 1 -number_particles_per_dimension 20 \ 1119 -dm_landau_batch_size 4 -number_spatial_vertices 6 -vertex_view_target 5 -grid_view_target 1 -dm_landau_batch_view_idx 1 \ 1120 -dm_landau_n 1.000018,1,1e-6 -dm_landau_thermal_temps 2,1,1 -dm_landau_ion_masses 2,180 -dm_landau_ion_charges 1,18 \ 1121 -ftop_ksp_rtol 1e-10 -ftop_ksp_type lsqr -ftop_pc_type bjacobi -ftop_sub_pc_factor_shift_type nonzero -ftop_sub_pc_type lu -ftop_ksp_error_if_not_converged \ 1122 -ksp_type gmres -ksp_error_if_not_converged -dm_landau_verbose 4 -print_entropy \ 1123 -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_rtol 1e-12 -ptof_ksp_error_if_not_converged\ 1124 -snes_converged_reason -snes_monitor -snes_rtol 1e-12 -snes_stol 1e-12 \ 1125 -ts_dt 0.01 -ts_rtol 1e-1 -ts_exact_final_time stepover -ts_max_snes_failures -1 -ts_max_steps 1 -ts_monitor -ts_type beuler 1126 test: 1127 suffix: cpu 1128 args: -dm_landau_device_type cpu -pc_type jacobi 1129 test: 1130 suffix: kokkos 1131 # failed on Sunspot@ALCF with sycl 1132 requires: kokkos_kernels !openmp !sycl 1133 args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -pc_type bjkokkos -pc_bjkokkos_ksp_type tfqmr -pc_bjkokkos_pc_type jacobi 1134 1135 testset: 1136 requires: double !defined(PETSC_USE_DMLANDAU_2D) 1137 output_file: output/ex30_3d.out 1138 args: -dim 3 -petscspace_degree 2 -dm_landau_num_species_grid 1,1 -dm_refine 0 -number_particles_per_dimension 10 -dm_plex_hash_location \ 1139 -dm_landau_batch_size 1 -number_spatial_vertices 1 -vertex_view_target 0 -grid_view_target 0 -dm_landau_batch_view_idx 0 \ 1140 -dm_landau_n 1.000018,1 -dm_landau_thermal_temps 2,1 -dm_landau_ion_masses 2 -dm_landau_ion_charges 1 \ 1141 -ftop_ksp_type cg -ftop_pc_type jacobi -ftop_ksp_rtol 1e-12 -ftop_ksp_error_if_not_converged -ksp_type preonly -pc_type lu -ksp_error_if_not_converged \ 1142 -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_rtol 1e-12 -ptof_ksp_error_if_not_converged \ 1143 -snes_converged_reason -snes_monitor -snes_rtol 1e-12 -snes_stol 1e-12 \ 1144 -ts_dt 0.1 -ts_exact_final_time stepover -ts_max_snes_failures -1 -ts_max_steps 1 -ts_monitor -ts_type beuler -print_entropy 1145 test: 1146 suffix: cpu_3d 1147 args: -dm_landau_device_type cpu 1148 test: 1149 suffix: kokkos_3d 1150 requires: kokkos_kernels !openmp 1151 args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -pc_type bjkokkos -pc_bjkokkos_ksp_type tfqmr -pc_bjkokkos_pc_type jacobi 1152 1153 test: 1154 suffix: conserve 1155 requires: !complex double defined(PETSC_USE_DMLANDAU_2D) !cuda 1156 args: -dm_landau_batch_size 4 -dm_refine 0 -dm_landau_num_species_grid 1 -dm_landau_thermal_temps 1 -petscspace_degree 3 -snes_converged_reason -ts_type beuler -ts_dt .1 \ 1157 -ts_max_steps 1 -ksp_type preonly -ksp_error_if_not_converged -snes_rtol 1e-14 -snes_stol 1e-14 -dm_landau_device_type cpu -number_particles_per_dimension 20 \ 1158 -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_rtol 1e-14 -ptof_ksp_error_if_not_converged -pc_type lu -dm_landau_simplex 1 -use_uniform_particle_grid false -dm_landau_sphere -print_entropy -number_particles_per_dimension 50 -ftop_ksp_type cg -ftop_pc_type jacobi -ftop_ksp_rtol 1e-14 1159 1160 TEST*/ 1161