1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <petsclandau.h> /*I "petsclandau.h" I*/ 3 #include <petscts.h> 4 #include <petscdmforest.h> 5 6 /* Landau collision operator */ 7 #define PETSC_THREAD_SYNC 8 #define PETSC_DEVICE_FUNC_DECL static 9 #include "land_tensors.h" 10 11 /* vector padding not supported */ 12 #define LANDAU_VL 1 13 14 static PetscErrorCode LandauGPUMapsDestroy(void *ptr) 15 { 16 P4estVertexMaps *maps = (P4estVertexMaps *)ptr; 17 PetscErrorCode ierr; 18 PetscFunctionBegin; 19 if (maps->deviceType != LANDAU_CPU) { 20 #if defined(PETSC_HAVE_KOKKOS_KERNELS) 21 if (maps->deviceType == LANDAU_KOKKOS) { 22 ierr = LandauKokkosDestroyMatMaps(maps);CHKERRQ(ierr); // imples Kokkos does 23 } // else could be CUDA 24 #elif defined(PETSC_HAVE_CUDA) 25 if (maps->deviceType == LANDAU_CUDA){ 26 ierr = LandauCUDADestroyMatMaps(maps);CHKERRQ(ierr); 27 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "maps->deviceType %D ?????",maps->deviceType); 28 #endif 29 } 30 ierr = PetscFree(maps->c_maps);CHKERRQ(ierr); 31 ierr = PetscFree(maps->gIdx);CHKERRQ(ierr); 32 ierr = PetscFree(maps);CHKERRQ(ierr); 33 PetscFunctionReturn(0); 34 } 35 36 /* ------------------------------------------------------------------- */ 37 /* 38 LandauFormJacobian_Internal - Evaluates Jacobian matrix. 39 40 Input Parameters: 41 . globX - input vector 42 . actx - optional user-defined context 43 . dim - dimension 44 45 Output Parameters: 46 . J0acP - Jacobian matrix filled, not created 47 */ 48 static PetscErrorCode LandauFormJacobian_Internal(Vec a_X, Mat JacP, const PetscInt dim, PetscReal shift, void *a_ctx) 49 { 50 LandauCtx *ctx = (LandauCtx*)a_ctx; 51 PetscErrorCode ierr; 52 PetscInt cStart, cEnd, elemMatSize; 53 PetscDS prob; 54 PetscSection section,globsection; 55 PetscInt numCells,totDim,ej,Nq,*Nbf,*Ncf,Nb,Ncx,Nf,d,f,fieldA,qj,N; 56 PetscQuadrature quad; 57 const PetscReal *quadWeights; 58 PetscTabulation *Tf; // used for CPU and print info 59 PetscReal Eq_m[LANDAU_MAX_SPECIES], m_0=ctx->m_0; /* normalize mass -- not needed! */ 60 PetscScalar *IPf=NULL; 61 const PetscScalar *xdata=NULL; 62 PetscLogDouble flops; 63 PetscContainer container; 64 P4estVertexMaps *maps=NULL; 65 66 PetscFunctionBegin; 67 PetscValidHeaderSpecific(a_X,VEC_CLASSID,1); 68 PetscValidHeaderSpecific(JacP,MAT_CLASSID,2); 69 PetscValidPointer(ctx,5); 70 /* check for matrix container for GPU assembly */ 71 ierr = PetscLogEventBegin(ctx->events[10],0,0,0,0);CHKERRQ(ierr); 72 ierr = PetscObjectQuery((PetscObject) JacP, "assembly_maps", (PetscObject *) &container);CHKERRQ(ierr); 73 if (container /* && ctx->deviceType != LANDAU_CPU */) { 74 if (!ctx->gpu_assembly) SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"GPU matrix container but no GPU assembly"); 75 ierr = PetscContainerGetPointer(container, (void **) &maps);CHKERRQ(ierr); 76 if (!maps) SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"empty GPU matrix container"); 77 } 78 if (ctx->plex == NULL) { 79 ierr = DMConvert(ctx->dmv, DMPLEX, &ctx->plex);CHKERRQ(ierr); 80 } 81 ierr = DMPlexGetHeightStratum(ctx->plex, 0, &cStart, &cEnd);CHKERRQ(ierr); 82 ierr = DMGetLocalSection(ctx->plex, §ion);CHKERRQ(ierr); 83 ierr = DMGetGlobalSection(ctx->plex, &globsection);CHKERRQ(ierr); 84 ierr = DMGetDS(ctx->plex, &prob);CHKERRQ(ierr); 85 ierr = PetscDSGetTabulation(prob, &Tf);CHKERRQ(ierr); // Bf, &Df 86 ierr = PetscDSGetDimensions(prob, &Nbf);CHKERRQ(ierr); Nb = Nbf[0]; /* number of vertices*S */ 87 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); if (Nf!=ctx->num_species) SETERRQ1(ctx->comm, PETSC_ERR_PLIB, "Nf %D != S",Nf); 88 ierr = PetscDSGetComponents(prob, &Ncf);CHKERRQ(ierr); Ncx = Ncf[0]; if (Ncx!=1) SETERRQ1(ctx->comm, PETSC_ERR_PLIB, "Nc %D != 1",Ncx); 89 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 90 numCells = cEnd - cStart; 91 ierr = PetscFEGetQuadrature(ctx->fe[0], &quad);CHKERRQ(ierr); 92 ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, &quadWeights);CHKERRQ(ierr); 93 if (Nb!=Nq) SETERRQ4(ctx->comm, PETSC_ERR_PLIB, "Nb!=Nq %D %D over integration or simplices? Tf[0]->Nb=%D dim=%D",Nb,Nq,Tf[0]->Nb,dim); 94 if (Nq >LANDAU_MAX_NQ) SETERRQ2(ctx->comm,PETSC_ERR_ARG_WRONG,"Order too high. Nq = %D > LANDAU_MAX_NQ (%D)",Nq,LANDAU_MAX_NQ); 95 if (LANDAU_DIM != dim) SETERRQ2(ctx->comm, PETSC_ERR_PLIB, "dim %D != LANDAU_DIM %d",dim,LANDAU_DIM); 96 if (!ctx->SData_d) { 97 ierr = PetscNew(&ctx->SData_d);CHKERRQ(ierr); 98 } 99 elemMatSize = totDim*totDim; // used for CPU and print info 100 ierr = PetscLogEventEnd(ctx->events[10],0,0,0,0);CHKERRQ(ierr); 101 ierr = VecGetSize(a_X,&N);CHKERRQ(ierr); 102 if (!ctx->init) { /* create static point data, Jacobian called first */ 103 PetscReal *invJ,*ww,*xx,*yy,*zz=NULL,*mass_w,*invJ_a; 104 const PetscInt nip = Nq*numCells; 105 106 ierr = PetscLogEventBegin(ctx->events[7],0,0,0,0);CHKERRQ(ierr); 107 ctx->init = PETSC_TRUE; 108 ierr = PetscInfo(ctx->plex, "Initialize static data\n");CHKERRQ(ierr); 109 /* collect f data, first time is for Jacobian, but make mass now */ 110 if (ctx->verbose > 1 || ctx->verbose > 0) { 111 ierr = PetscPrintf(ctx->comm,"[%D]%s: %D IPs, %D cells, totDim=%D, Nb=%D, Nq=%D, elemMatSize=%D, dim=%D, Tab: Nb=%D Nf=%D Np=%D cdim=%D N=%D shift=%g\n", 112 0,"FormLandau",Nq*numCells,numCells, totDim, Nb, Nq, elemMatSize, dim, Tf[0]->Nb, Nf, Tf[0]->Np, Tf[0]->cdim, N, shift);CHKERRQ(ierr); 113 } 114 ierr = PetscMalloc5(nip,&mass_w,nip,&ww,nip,&xx,nip,&yy,nip*dim*dim,&invJ_a);CHKERRQ(ierr); 115 if (dim==3) { 116 ierr = PetscMalloc1(nip,&zz);CHKERRQ(ierr); 117 } 118 for (ej = 0 ; ej < numCells; ++ej) { 119 PetscReal vj[LANDAU_MAX_NQ*LANDAU_DIM],detJj[LANDAU_MAX_NQ], Jdummy[LANDAU_MAX_NQ*LANDAU_DIM*LANDAU_DIM]; 120 invJ = invJ_a ? invJ_a + ej * Nq*dim*dim : NULL; 121 ierr = DMPlexComputeCellGeometryFEM(ctx->plex, cStart+ej, quad, vj, Jdummy, invJ, detJj);CHKERRQ(ierr); 122 /* create dynamic point data */ 123 for (qj = 0; qj < Nq; ++qj) { 124 PetscInt gidx = (ej*Nq + qj); 125 mass_w[gidx] = detJj[qj] * quadWeights[qj]; 126 if (dim==2) mass_w[gidx] *= 2.*PETSC_PI*vj[qj * dim + 0]; /* cylindrical coordinate, w/o 2pi */ 127 xx[gidx] = vj[qj * dim + 0]; /* coordinate */ 128 yy[gidx] = vj[qj * dim + 1]; 129 if (dim==3) zz[gidx] = vj[qj * dim + 2]; 130 ww[gidx] = detJj[qj] * quadWeights[qj]; 131 if (dim==2) ww[gidx] *= xx[gidx]; /* cylindrical coordinate, w/o 2pi */ 132 } /* q */ 133 } /* ej */ 134 /* cache static data */ 135 if (ctx->deviceType == LANDAU_CUDA || ctx->deviceType == LANDAU_KOKKOS) { 136 #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_KOKKOS) 137 PetscReal invMass[LANDAU_MAX_SPECIES],nu_alpha[LANDAU_MAX_SPECIES], nu_beta[LANDAU_MAX_SPECIES]; 138 for (fieldA=0;fieldA<Nf;fieldA++) { 139 invMass[fieldA] = m_0/ctx->masses[fieldA]; 140 nu_alpha[fieldA] = PetscSqr(ctx->charges[fieldA]/m_0)*m_0/ctx->masses[fieldA]; 141 nu_beta[fieldA] = PetscSqr(ctx->charges[fieldA]/ctx->epsilon0)*ctx->lnLam / (8*PETSC_PI) * ctx->t_0*ctx->n_0/PetscPowReal(ctx->v_0,3); 142 } 143 if (ctx->deviceType == LANDAU_CUDA) { 144 #if defined(PETSC_HAVE_CUDA) 145 ierr = LandauCUDAStaticDataSet(ctx->plex,Nq,nu_alpha,nu_beta,invMass,invJ_a,mass_w,xx,yy,zz,ww,ctx->SData_d);CHKERRQ(ierr); 146 #else 147 SETERRQ1(ctx->comm,PETSC_ERR_ARG_WRONG,"-landau_device_type %s not built","cuda"); 148 #endif 149 } else if (ctx->deviceType == LANDAU_KOKKOS) { 150 #if defined(PETSC_HAVE_KOKKOS) 151 ierr = LandauKokkosStaticDataSet(ctx->plex,Nq,nu_alpha,nu_beta,invMass,invJ_a,mass_w,xx,yy,zz,ww,ctx->SData_d);CHKERRQ(ierr); 152 #else 153 SETERRQ1(ctx->comm,PETSC_ERR_ARG_WRONG,"-landau_device_type %s not built","kokkos"); 154 #endif 155 } 156 #endif 157 /* free */ 158 ierr = PetscFree5(mass_w,ww,xx,yy,invJ_a);CHKERRQ(ierr); 159 if (dim==3) { 160 ierr = PetscFree(zz);CHKERRQ(ierr); 161 } 162 } else { /* CPU version, just copy in, only use part */ 163 ctx->SData_d->w = (void*)ww; 164 ctx->SData_d->x = (void*)xx; 165 ctx->SData_d->y = (void*)yy; 166 ctx->SData_d->z = (void*)zz; 167 ctx->SData_d->invJ = (void*)invJ_a; 168 ctx->SData_d->mass_w = (void*)mass_w; 169 } 170 ierr = PetscLogEventEnd(ctx->events[7],0,0,0,0);CHKERRQ(ierr); 171 } 172 if (shift==0) { /* create dynamic point data */ 173 ierr = PetscLogEventBegin(ctx->events[1],0,0,0,0);CHKERRQ(ierr); 174 ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 175 flops = (PetscLogDouble)numCells*(PetscLogDouble)Nq*(PetscLogDouble)(5*dim*dim*Nf*Nf + 165); 176 for (fieldA=0;fieldA<Nf;fieldA++) { 177 Eq_m[fieldA] = ctx->Ez * ctx->t_0 * ctx->charges[fieldA] / (ctx->v_0 * ctx->masses[fieldA]); /* normalize dimensionless */ 178 if (dim==2) Eq_m[fieldA] *= 2 * PETSC_PI; /* add the 2pi term that is not in Landau */ 179 } 180 if (!ctx->gpu_assembly || !container) { 181 Vec locX; 182 ierr = DMGetLocalVector(ctx->plex, &locX);CHKERRQ(ierr); 183 ierr = VecZeroEntries(locX);CHKERRQ(ierr); /* zero BCs so don't set */ 184 ierr = DMGlobalToLocalBegin(ctx->plex, a_X, INSERT_VALUES, locX);CHKERRQ(ierr); 185 ierr = DMGlobalToLocalEnd (ctx->plex, a_X, INSERT_VALUES, locX);CHKERRQ(ierr); 186 ierr = PetscMalloc1(Nq*numCells*Nf,&IPf);CHKERRQ(ierr); 187 for (ej = 0 ; ej < numCells; ++ej) { 188 PetscScalar *coef = NULL; 189 ierr = DMPlexVecGetClosure(ctx->plex, section, locX, cStart+ej, NULL, &coef);CHKERRQ(ierr); 190 ierr = PetscMemcpy(&IPf[ej*Nb*Nf],coef,Nb*Nf*sizeof(PetscScalar));CHKERRQ(ierr); /* change if LandauIPReal != PetscScalar */ 191 ierr = DMPlexVecRestoreClosure(ctx->plex, section, locX, cStart+ej, NULL, &coef);CHKERRQ(ierr); 192 } /* ej */ 193 ierr = DMRestoreLocalVector(ctx->plex, &locX);CHKERRQ(ierr); 194 } else { 195 PetscMemType mtype; 196 ierr = VecGetArrayReadAndMemType(a_X,&xdata,&mtype);CHKERRQ(ierr); 197 } 198 ierr = PetscLogEventEnd(ctx->events[1],0,0,0,0);CHKERRQ(ierr); 199 } else { 200 flops = (PetscLogDouble)numCells*(PetscLogDouble)Nq*(PetscLogDouble)(5*dim*dim*Nf*Nf); 201 } 202 /* do it */ 203 if (ctx->deviceType == LANDAU_CUDA || ctx->deviceType == LANDAU_KOKKOS) { 204 if (ctx->deviceType == LANDAU_CUDA) { 205 #if defined(PETSC_HAVE_CUDA) 206 ierr = LandauCUDAJacobian(ctx->plex,Nq,Eq_m,IPf,N,xdata,ctx->SData_d,ctx->subThreadBlockSize,shift,ctx->events,JacP);CHKERRQ(ierr); 207 #else 208 SETERRQ1(ctx->comm,PETSC_ERR_ARG_WRONG,"-landau_device_type %s not built","cuda"); 209 #endif 210 } else if (ctx->deviceType == LANDAU_KOKKOS) { 211 #if defined(PETSC_HAVE_KOKKOS) 212 ierr = LandauKokkosJacobian(ctx->plex,Nq,Eq_m,IPf,N,xdata,ctx->SData_d,ctx->subThreadBlockSize,shift,ctx->events,JacP);CHKERRQ(ierr); 213 #else 214 SETERRQ1(ctx->comm,PETSC_ERR_ARG_WRONG,"-landau_device_type %s not built","kokkos"); 215 #endif 216 } 217 } else { /* CPU version */ 218 PetscInt ei, qi; 219 PetscScalar *elemMat,coef_buff[LANDAU_MAX_SPECIES*LANDAU_MAX_NQ]; 220 PetscReal *ff, *dudx, *dudy, *dudz, *invJ, *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, *mass_w = (PetscReal*)ctx->SData_d->mass_w; 221 const PetscInt nip = Nq*numCells; 222 const PetscReal *const BB = Tf[0]->T[0], * const DD = Tf[0]->T[1]; 223 PetscReal Eq_m[LANDAU_MAX_SPECIES], invMass[LANDAU_MAX_SPECIES], nu_alpha[LANDAU_MAX_SPECIES], nu_beta[LANDAU_MAX_SPECIES]; 224 if (shift!=0.0) { // mass 225 ierr = PetscMalloc1(elemMatSize, &elemMat);CHKERRQ(ierr); 226 } else { /* compute f and df and init data for Jacobian */ 227 ierr = PetscLogEventBegin(ctx->events[8],0,0,0,0);CHKERRQ(ierr); 228 for (fieldA=0;fieldA<Nf;fieldA++) { 229 invMass[fieldA] = m_0/ctx->masses[fieldA]; 230 Eq_m[fieldA] = ctx->Ez * ctx->t_0 * ctx->charges[fieldA] / (ctx->v_0 * ctx->masses[fieldA]); /* normalize dimensionless */ 231 if (dim==2) Eq_m[fieldA] *= 2 * PETSC_PI; /* add the 2pi term that is not in Landau */ 232 nu_alpha[fieldA] = PetscSqr(ctx->charges[fieldA]/m_0)*m_0/ctx->masses[fieldA]; 233 nu_beta[fieldA] = PetscSqr(ctx->charges[fieldA]/ctx->epsilon0)*ctx->lnLam / (8*PETSC_PI) * ctx->t_0*ctx->n_0/PetscPowReal(ctx->v_0,3); 234 } 235 ierr = PetscMalloc5(elemMatSize, &elemMat, nip*Nf, &ff, nip*Nf, &dudx, nip*Nf, &dudy, dim==3 ? nip*Nf : 0, &dudz);CHKERRQ(ierr); 236 for (ei = cStart, invJ = invJ_a; ei < cEnd; ++ei, invJ += Nq*dim*dim) { 237 PetscScalar *coef; 238 PetscInt b,f,idx,q; 239 PetscReal u_x[LANDAU_MAX_SPECIES][LANDAU_DIM]; 240 if (IPf) { 241 coef = &IPf[ei*Nb*Nf]; // this is const 242 } else { 243 if (!maps) SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"!maps"); 244 coef = coef_buff; 245 for (f = 0; f < Nf; ++f) { 246 LandauIdx *const Idxs = &maps->gIdx[ei-cStart][f][0]; 247 for (b = 0; b < Nb; ++b) { 248 idx = Idxs[b]; 249 if (idx >= 0) { 250 coef[f*Nb+b] = xdata[idx]; 251 } else { 252 idx = -idx - 1; 253 coef[f*Nb+b] = 0; 254 for (q = 0; q < maps->num_face; q++) { 255 PetscInt id = maps->c_maps[idx][q].gid; 256 PetscScalar scale = maps->c_maps[idx][q].scale; 257 coef[f*Nb+b] += scale*xdata[id]; 258 } 259 } 260 } 261 } 262 } 263 /* get f and df */ 264 for (qi = 0; qi < Nq; ++qi) { 265 const PetscReal *Bq = &BB[qi*Nb]; 266 const PetscReal *Dq = &DD[qi*Nb*dim]; 267 const PetscInt gidx = ei*Nq + qi; 268 /* get f & df */ 269 for (f = 0; f < Nf; ++f) { 270 PetscInt b, e; 271 PetscReal refSpaceDer[LANDAU_DIM]; 272 ff[gidx + f*nip] = 0.0; 273 for (d = 0; d < LANDAU_DIM; ++d) refSpaceDer[d] = 0.0; 274 for (b = 0; b < Nb; ++b) { 275 const PetscInt cidx = b; 276 ff[gidx + f*nip] += Bq[cidx]*PetscRealPart(coef[f*Nb+cidx]); 277 for (d = 0; d < dim; ++d) refSpaceDer[d] += Dq[cidx*dim+d]*PetscRealPart(coef[f*Nb+cidx]); 278 } 279 for (d = 0; d < dim; ++d) { 280 for (e = 0, u_x[f][d] = 0.0; e < dim; ++e) { 281 u_x[f][d] += invJ[qi * dim * dim + e*dim+d]*refSpaceDer[e]; 282 } 283 } 284 } 285 for (f=0;f<Nf;f++) { 286 dudx[gidx + f*nip] = u_x[f][0]; 287 dudy[gidx + f*nip] = u_x[f][1]; 288 #if LANDAU_DIM==3 289 dudz[gidx + f*nip] = u_x[f][2]; 290 #endif 291 } 292 } 293 } 294 ierr = PetscLogEventEnd(ctx->events[8],0,0,0,0);CHKERRQ(ierr); 295 } 296 for (ej = cStart, invJ = invJ_a; ej < cEnd; ++ej, invJ += Nq*dim*dim) { 297 ierr = PetscMemzero(elemMat, totDim *totDim * sizeof(PetscScalar));CHKERRQ(ierr); 298 ierr = PetscLogEventBegin(ctx->events[4],0,0,0,0);CHKERRQ(ierr); 299 ierr = PetscLogFlops((PetscLogDouble)Nq*flops);CHKERRQ(ierr); 300 //printf("\t:%d.%d) Invj[0] = %e (%d)\n",ej,qj,invJ[0],(int)(invJ-invJ_a)); 301 for (qj = 0; qj < Nq; ++qj) { 302 const PetscReal * const BB = Tf[0]->T[0], * const DD = Tf[0]->T[1]; 303 PetscReal g0[LANDAU_MAX_SPECIES], g2[LANDAU_MAX_SPECIES][LANDAU_DIM], g3[LANDAU_MAX_SPECIES][LANDAU_DIM][LANDAU_DIM]; 304 PetscInt d,d2,dp,d3,ipidx,fieldA; 305 const PetscInt jpidx = Nq*(ej-cStart) + qj; 306 if (shift==0.0) { 307 const PetscReal * const invJj = &invJ[qj*dim*dim]; 308 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]; 309 const PetscReal vj[3] = {xx[jpidx], yy[jpidx], zz ? zz[jpidx] : 0}, wj = ww[jpidx]; 310 311 // create g2 & g3 312 for (d=0;d<dim;d++) { // clear accumulation data D & K 313 gg2_temp[d] = 0; 314 for (d2=0;d2<dim;d2++) gg3_temp[d][d2] = 0; 315 } 316 for (ipidx = 0; ipidx < nip; ipidx++) { 317 const PetscReal wi = ww[ipidx], x = xx[ipidx], y = yy[ipidx]; 318 PetscReal temp1[3] = {0, 0, 0}, temp2 = 0; 319 #if LANDAU_DIM==2 320 PetscReal Ud[2][2], Uk[2][2]; 321 LandauTensor2D(vj, x, y, Ud, Uk, (ipidx==jpidx) ? 0. : 1.); 322 #else 323 PetscReal U[3][3], z = zz[ipidx]; 324 LandauTensor3D(vj, x, y, z, U, (ipidx==jpidx) ? 0. : 1.); 325 #endif 326 for (fieldA = 0; fieldA < Nf; ++fieldA) { 327 temp1[0] += dudx[ipidx + fieldA*nip]*nu_beta[fieldA]*invMass[fieldA]; 328 temp1[1] += dudy[ipidx + fieldA*nip]*nu_beta[fieldA]*invMass[fieldA]; 329 #if LANDAU_DIM==3 330 temp1[2] += dudz[ipidx + fieldA*nip]*nu_beta[fieldA]*invMass[fieldA]; 331 #endif 332 temp2 += ff[ipidx + fieldA*nip]*nu_beta[fieldA]; 333 } 334 temp1[0] *= wi; 335 temp1[1] *= wi; 336 #if LANDAU_DIM==3 337 temp1[2] *= wi; 338 #endif 339 temp2 *= wi; 340 #if LANDAU_DIM==2 341 for (d2 = 0; d2 < 2; d2++) { 342 for (d3 = 0; d3 < 2; ++d3) { 343 /* K = U * grad(f): g2=e: i,A */ 344 gg2_temp[d2] += Uk[d2][d3]*temp1[d3]; 345 /* D = -U * (I \kron (fx)): g3=f: i,j,A */ 346 gg3_temp[d2][d3] += Ud[d2][d3]*temp2; 347 } 348 } 349 #else 350 for (d2 = 0; d2 < 3; ++d2) { 351 for (d3 = 0; d3 < 3; ++d3) { 352 /* K = U * grad(f): g2 = e: i,A */ 353 gg2_temp[d2] += U[d2][d3]*temp1[d3]; 354 /* D = -U * (I \kron (fx)): g3 = f: i,j,A */ 355 gg3_temp[d2][d3] += U[d2][d3]*temp2; 356 } 357 } 358 #endif 359 } /* IPs */ 360 // add alpha and put in gg2/3 361 for (fieldA = 0; fieldA < Nf; ++fieldA) { 362 for (d2 = 0; d2 < dim; d2++) { 363 gg2[fieldA][d2] = gg2_temp[d2]*nu_alpha[fieldA]; 364 for (d3 = 0; d3 < dim; d3++) { 365 gg3[fieldA][d2][d3] = -gg3_temp[d2][d3]*nu_alpha[fieldA]*invMass[fieldA]; 366 } 367 } 368 } 369 /* add electric field term once per IP */ 370 for (fieldA = 0; fieldA < Nf; ++fieldA) { 371 gg2[fieldA][dim-1] += Eq_m[fieldA]; 372 } 373 /* Jacobian transform - g2, g3 */ 374 for (fieldA = 0; fieldA < Nf; ++fieldA) { 375 for (d = 0; d < dim; ++d) { 376 g2[fieldA][d] = 0.0; 377 for (d2 = 0; d2 < dim; ++d2) { 378 g2[fieldA][d] += invJj[d*dim+d2]*gg2[fieldA][d2]; 379 g3[fieldA][d][d2] = 0.0; 380 for (d3 = 0; d3 < dim; ++d3) { 381 for (dp = 0; dp < dim; ++dp) { 382 g3[fieldA][d][d2] += invJj[d*dim + d3]*gg3[fieldA][d3][dp]*invJj[d2*dim + dp]; 383 } 384 } 385 g3[fieldA][d][d2] *= wj; 386 } 387 g2[fieldA][d] *= wj; 388 } 389 } 390 } else { // mass 391 /* Jacobian transform - g0 */ 392 for (fieldA = 0; fieldA < Nf; ++fieldA) { 393 g0[fieldA] = mass_w[jpidx] * shift; // move this to below and remove g0 394 } 395 } 396 /* FE matrix construction */ 397 { 398 PetscInt fieldA,d,f,d2,g; 399 const PetscReal *BJq = &BB[qj*Nb], *DIq = &DD[qj*Nb*dim]; 400 /* assemble - on the diagonal (I,I) */ 401 for (fieldA = 0; fieldA < Nf ; fieldA++) { 402 for (f = 0; f < Nb ; f++) { 403 const PetscInt i = fieldA*Nb + f; /* Element matrix row */ 404 for (g = 0; g < Nb; ++g) { 405 const PetscInt j = fieldA*Nb + g; /* Element matrix column */ 406 const PetscInt fOff = i*totDim + j; 407 if (shift==0.0) { 408 for (d = 0; d < dim; ++d) { 409 elemMat[fOff] += DIq[f*dim+d]*g2[fieldA][d]*BJq[g]; 410 //printf("\t:%d.%d.%d.%d.%d.%d) elemMat=%e += %e %e %e\n",ej,qj,fieldA,f,g,d,elemMat[fOff],DIq[f*dim+d],g2[fieldA][d],BJq[g]); 411 for (d2 = 0; d2 < dim; ++d2) { 412 elemMat[fOff] += DIq[f*dim + d]*g3[fieldA][d][d2]*DIq[g*dim + d2]; 413 } 414 } 415 } else { // mass 416 elemMat[fOff] += BJq[f]*g0[fieldA]*BJq[g]; 417 } 418 } 419 } 420 } 421 } 422 } /* qj loop */ 423 ierr = PetscLogEventEnd(ctx->events[4],0,0,0,0);CHKERRQ(ierr); 424 /* assemble matrix */ 425 ierr = PetscLogEventBegin(ctx->events[6],0,0,0,0);CHKERRQ(ierr); 426 if (!maps) { 427 ierr = DMPlexMatSetClosure(ctx->plex, section, globsection, JacP, ej, elemMat, ADD_VALUES);CHKERRQ(ierr); 428 } else { // GPU like assembly for debugging 429 PetscInt fieldA,idx,q,f,g,d,nr,nc,rows0[LANDAU_MAX_Q_FACE],cols0[LANDAU_MAX_Q_FACE]={0},rows[LANDAU_MAX_Q_FACE],cols[LANDAU_MAX_Q_FACE]; 430 PetscScalar vals[LANDAU_MAX_Q_FACE*LANDAU_MAX_Q_FACE],row_scale[LANDAU_MAX_Q_FACE],col_scale[LANDAU_MAX_Q_FACE]={0}; 431 /* assemble - from the diagonal (I,I) in this format for DMPlexMatSetClosure */ 432 for (fieldA = 0; fieldA < Nf ; fieldA++) { 433 LandauIdx *const Idxs = &maps->gIdx[ej-cStart][fieldA][0]; 434 for (f = 0; f < Nb ; f++) { 435 idx = Idxs[f]; 436 if (idx >= 0) { 437 nr = 1; 438 rows0[0] = idx; 439 row_scale[0] = 1.; 440 } else { 441 idx = -idx - 1; 442 nr = maps->num_face; 443 for (q = 0; q < maps->num_face; q++) { 444 rows0[q] = maps->c_maps[idx][q].gid; 445 row_scale[q] = maps->c_maps[idx][q].scale; 446 } 447 } 448 for (g = 0; g < Nb; ++g) { 449 idx = Idxs[g]; 450 if (idx >= 0) { 451 nc = 1; 452 cols0[0] = idx; 453 col_scale[0] = 1.; 454 } else { 455 idx = -idx - 1; 456 nc = maps->num_face; 457 for (q = 0; q < maps->num_face; q++) { 458 cols0[q] = maps->c_maps[idx][q].gid; 459 col_scale[q] = maps->c_maps[idx][q].scale; 460 } 461 } 462 const PetscInt i = fieldA*Nb + f; /* Element matrix row */ 463 const PetscInt j = fieldA*Nb + g; /* Element matrix column */ 464 const PetscScalar Aij = elemMat[i*totDim + j]; 465 for (q = 0; q < nr; q++) rows[q] = rows0[q]; 466 for (q = 0; q < nc; q++) cols[q] = cols0[q]; 467 for (q = 0; q < nr; q++) { 468 for (d = 0; d < nc; d++) { 469 vals[q*nc + d] = row_scale[q]*col_scale[d]*Aij; 470 } 471 } 472 ierr = MatSetValues(JacP,nr,rows,nc,cols,vals,ADD_VALUES);CHKERRQ(ierr); 473 } 474 } 475 } 476 } 477 if (ej==-1) { 478 PetscErrorCode ierr2; 479 ierr2 = PetscPrintf(ctx->comm,"CPU Element matrix\n");CHKERRQ(ierr2); 480 for (d = 0; d < totDim; ++d){ 481 for (f = 0; f < totDim; ++f) {ierr2 = PetscPrintf(ctx->comm," %12.5e", PetscRealPart(elemMat[d*totDim + f]));CHKERRQ(ierr2);} 482 ierr2 = PetscPrintf(ctx->comm,"\n");CHKERRQ(ierr2); 483 } 484 exit(12); 485 } 486 ierr = PetscLogEventEnd(ctx->events[6],0,0,0,0);CHKERRQ(ierr); 487 } /* ej cells loop, not cuda */ 488 if (shift!=0.0) { // mass 489 ierr = PetscFree(elemMat);CHKERRQ(ierr); 490 } else { 491 ierr = PetscFree5(elemMat, ff, dudx, dudy, dudz);CHKERRQ(ierr); 492 } 493 } /* CPU version */ 494 495 /* assemble matrix or vector */ 496 ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 497 ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 498 #define MAP_BF_SIZE (128*LANDAU_DIM*LANDAU_MAX_Q_FACE*LANDAU_MAX_SPECIES) 499 if (ctx->gpu_assembly && !container) { 500 PetscScalar elemMatrix[LANDAU_MAX_NQ*LANDAU_MAX_NQ*LANDAU_MAX_SPECIES*LANDAU_MAX_SPECIES], *elMat; 501 pointInterpolationP4est pointMaps[MAP_BF_SIZE][LANDAU_MAX_Q_FACE]; 502 PetscInt q,eidx,fieldA; 503 MatType type; 504 ierr = PetscInfo1(ctx->plex, "Make GPU maps %D\n",1);CHKERRQ(ierr); 505 ierr = MatGetType(JacP,&type);CHKERRQ(ierr); 506 ierr = PetscLogEventBegin(ctx->events[2],0,0,0,0);CHKERRQ(ierr); 507 ierr = PetscMalloc(sizeof(P4estVertexMaps), &maps);CHKERRQ(ierr); 508 ierr = PetscContainerCreate(PETSC_COMM_SELF, &container);CHKERRQ(ierr); 509 ierr = PetscContainerSetPointer(container, (void *)maps);CHKERRQ(ierr); 510 ierr = PetscContainerSetUserDestroy(container, LandauGPUMapsDestroy);CHKERRQ(ierr); 511 ierr = PetscObjectCompose((PetscObject) JacP, "assembly_maps", (PetscObject) container);CHKERRQ(ierr); 512 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 513 // make maps 514 maps->data = NULL; 515 maps->num_elements = numCells; 516 maps->num_face = (PetscInt)(pow(Nq,1./((double)dim))+.001); // Q 517 maps->num_face = (PetscInt)(pow(maps->num_face,(double)(dim-1))+.001); // Q^2 518 maps->num_reduced = 0; 519 maps->deviceType = ctx->deviceType; 520 521 // count reduced and get 522 ierr = PetscMalloc(maps->num_elements * sizeof *maps->gIdx, &maps->gIdx);CHKERRQ(ierr); 523 for (fieldA=0;fieldA<Nf;fieldA++) { 524 for (ej = cStart, eidx = 0 ; ej < cEnd; ++ej, ++eidx) { 525 for (q = 0; q < Nb; ++q) { 526 PetscInt numindices,*indices; 527 PetscScalar *valuesOrig = elMat = elemMatrix; 528 ierr = PetscMemzero(elMat, totDim*totDim*sizeof(PetscScalar));CHKERRQ(ierr); 529 elMat[ (fieldA*Nb + q)*totDim + fieldA*Nb + q] = 1; 530 ierr = DMPlexGetClosureIndices(ctx->plex, section, globsection, ej, PETSC_TRUE, &numindices, &indices, NULL, (PetscScalar **) &elMat);CHKERRQ(ierr); 531 for (f = 0 ; f < numindices ; ++f) { // look for a non-zero on the diagonal 532 if (PetscAbs(PetscRealPart(elMat[f*numindices + f])) > PETSC_MACHINE_EPSILON) { 533 // found it 534 if (PetscAbs(PetscRealPart(elMat[f*numindices + f] - 1.)) < PETSC_MACHINE_EPSILON) { 535 maps->gIdx[eidx][fieldA][q] = (LandauIdx)indices[f]; // normal vertex 1.0 536 } else { //found a constraint 537 int jj = 0; 538 PetscReal sum = 0; 539 const PetscInt ff = f; 540 maps->gIdx[eidx][fieldA][q] = -maps->num_reduced - 1; // gid = -(idx+1): idx = -gid - 1 541 do { // constraints are continous in Plex - exploit that here 542 int ii; 543 for (ii = 0, pointMaps[maps->num_reduced][jj].scale = 0; ii < maps->num_face; ii++) { // DMPlex puts them all together 544 if (ff + ii < numindices) { 545 pointMaps[maps->num_reduced][jj].scale += PetscRealPart(elMat[f*numindices + ff + ii]); 546 } 547 } 548 sum += pointMaps[maps->num_reduced][jj].scale; 549 if (pointMaps[maps->num_reduced][jj].scale == 0) pointMaps[maps->num_reduced][jj].gid = -1; // 3D has Q and Q^2 interps -- all contiguous??? 550 else pointMaps[maps->num_reduced][jj].gid = indices[f]; 551 } while (++jj < maps->num_face && ++f < numindices); // jj is incremented if we hit the end 552 while (jj++ < maps->num_face) { 553 pointMaps[maps->num_reduced][jj].scale = 0; 554 pointMaps[maps->num_reduced][jj].gid = -1; 555 } 556 if (PetscAbs(sum-1.0) > 10*PETSC_MACHINE_EPSILON) { // debug 557 int d,f; 558 PetscReal tmp = 0; 559 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,sum,LANDAU_MAX_Q_FACE,maps->num_face); 560 for (d = 0, tmp = 0; d < numindices; ++d){ 561 if (tmp!=0 && PetscAbs(tmp-1.0) > 10*PETSC_MACHINE_EPSILON) ierr = PetscPrintf(PETSC_COMM_WORLD,"%3D) %3D: ",d,indices[d]);CHKERRQ(ierr); 562 for (f = 0; f < numindices; ++f) { 563 tmp += PetscRealPart(elMat[d*numindices + f]); 564 } 565 if (tmp!=0) ierr = PetscPrintf(ctx->comm," | %22.16e\n",tmp);CHKERRQ(ierr); 566 } 567 } 568 maps->num_reduced++; 569 if (maps->num_reduced>=MAP_BF_SIZE) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "maps->num_reduced %d > %d",maps->num_reduced,MAP_BF_SIZE); 570 } 571 break; 572 } 573 } 574 // cleanup 575 ierr = DMPlexRestoreClosureIndices(ctx->plex, section, globsection, ej, PETSC_TRUE, &numindices, &indices, NULL, (PetscScalar **) &elMat);CHKERRQ(ierr); 576 if (elMat != valuesOrig) {ierr = DMRestoreWorkArray(ctx->plex, numindices*numindices, MPIU_SCALAR, &elMat);} 577 } 578 } 579 } 580 // allocate and copy point datamaps->gIdx[eidx][field][q] -- for CPU version of this code, for debugging 581 ierr = PetscMalloc(maps->num_reduced * sizeof *maps->c_maps, &maps->c_maps);CHKERRQ(ierr); 582 for (ej = 0; ej < maps->num_reduced; ++ej) { 583 for (q = 0; q < maps->num_face; ++q) { 584 maps->c_maps[ej][q].scale = pointMaps[ej][q].scale; 585 maps->c_maps[ej][q].gid = pointMaps[ej][q].gid; 586 } 587 } 588 #if defined(PETSC_HAVE_KOKKOS_KERNELS) 589 if (ctx->deviceType == LANDAU_KOKKOS) { 590 ierr = LandauKokkosCreateMatMaps(maps, pointMaps,Nf,Nq);CHKERRQ(ierr); // imples Kokkos does 591 } // else could be CUDA 592 #endif 593 #if defined(PETSC_HAVE_CUDA) 594 if (ctx->deviceType == LANDAU_CUDA){ 595 ierr = LandauCUDACreateMatMaps(maps, pointMaps,Nf,Nq);CHKERRQ(ierr); 596 } 597 #endif 598 599 ierr = PetscLogEventEnd(ctx->events[2],0,0,0,0);CHKERRQ(ierr); 600 } 601 /* clean up */ 602 if (IPf) { 603 ierr = PetscFree(IPf);CHKERRQ(ierr); 604 } 605 if (xdata) { 606 ierr = VecRestoreArrayReadAndMemType(a_X,&xdata);CHKERRQ(ierr); 607 } 608 609 PetscFunctionReturn(0); 610 } 611 612 #if defined(LANDAU_ADD_BCS) 613 static void zero_bc(PetscInt dim, PetscInt Nf, PetscInt NfAux, 614 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 615 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 616 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uexact[]) 617 { 618 uexact[0] = 0; 619 } 620 #endif 621 622 #define MATVEC2(__a,__x,__p) {int i,j; for (i=0.; i<2; i++) {__p[i] = 0; for (j=0.; j<2; j++) __p[i] += __a[i][j]*__x[j]; }} 623 static void CircleInflate(PetscReal r1, PetscReal r2, PetscReal r0, PetscInt num_sections, PetscReal x, PetscReal y, 624 PetscReal *outX, PetscReal *outY) 625 { 626 PetscReal rr = PetscSqrtReal(x*x + y*y), outfact, efact; 627 if (rr < r1 + PETSC_SQRT_MACHINE_EPSILON) { 628 *outX = x; *outY = y; 629 } else { 630 const PetscReal xy[2] = {x,y}, sinphi=y/rr, cosphi=x/rr; 631 PetscReal cth,sth,xyprime[2],Rth[2][2],rotcos,newrr; 632 if (num_sections==2) { 633 rotcos = 0.70710678118654; 634 outfact = 1.5; efact = 2.5; 635 /* rotate normalized vector into [-pi/4,pi/4) */ 636 if (sinphi >= 0.) { /* top cell, -pi/2 */ 637 cth = 0.707106781186548; sth = -0.707106781186548; 638 } else { /* bottom cell -pi/8 */ 639 cth = 0.707106781186548; sth = .707106781186548; 640 } 641 } else if (num_sections==3) { 642 rotcos = 0.86602540378443; 643 outfact = 1.5; efact = 2.5; 644 /* rotate normalized vector into [-pi/6,pi/6) */ 645 if (sinphi >= 0.5) { /* top cell, -pi/3 */ 646 cth = 0.5; sth = -0.866025403784439; 647 } else if (sinphi >= -.5) { /* mid cell 0 */ 648 cth = 1.; sth = .0; 649 } else { /* bottom cell +pi/3 */ 650 cth = 0.5; sth = 0.866025403784439; 651 } 652 } else if (num_sections==4) { 653 rotcos = 0.9238795325112; 654 outfact = 1.5; efact = 3; 655 /* rotate normalized vector into [-pi/8,pi/8) */ 656 if (sinphi >= 0.707106781186548) { /* top cell, -3pi/8 */ 657 cth = 0.38268343236509; sth = -0.923879532511287; 658 } else if (sinphi >= 0.) { /* mid top cell -pi/8 */ 659 cth = 0.923879532511287; sth = -.38268343236509; 660 } else if (sinphi >= -0.707106781186548) { /* mid bottom cell + pi/8 */ 661 cth = 0.923879532511287; sth = 0.38268343236509; 662 } else { /* bottom cell + 3pi/8 */ 663 cth = 0.38268343236509; sth = .923879532511287; 664 } 665 } else { 666 cth = 0.; sth = 0.; rotcos = 0; efact = 0; 667 } 668 Rth[0][0] = cth; Rth[0][1] =-sth; 669 Rth[1][0] = sth; Rth[1][1] = cth; 670 MATVEC2(Rth,xy,xyprime); 671 if (num_sections==2) { 672 newrr = xyprime[0]/rotcos; 673 } else { 674 PetscReal newcosphi=xyprime[0]/rr, rin = r1, rout = rr - rin; 675 PetscReal routmax = r0*rotcos/newcosphi - rin, nroutmax = r0 - rin, routfrac = rout/routmax; 676 newrr = rin + routfrac*nroutmax; 677 } 678 *outX = cosphi*newrr; *outY = sinphi*newrr; 679 /* grade */ 680 PetscReal fact,tt,rs,re, rr = PetscSqrtReal(PetscSqr(*outX) + PetscSqr(*outY)); 681 if (rr > r2) { rs = r2; re = r0; fact = outfact;} /* outer zone */ 682 else { rs = r1; re = r2; fact = efact;} /* electron zone */ 683 tt = (rs + PetscPowReal((rr - rs)/(re - rs),fact) * (re-rs)) / rr; 684 *outX *= tt; 685 *outY *= tt; 686 } 687 } 688 689 static PetscErrorCode GeometryDMLandau(DM base, PetscInt point, PetscInt dim, const PetscReal abc[], PetscReal xyz[], void *a_ctx) 690 { 691 LandauCtx *ctx = (LandauCtx*)a_ctx; 692 PetscReal r = abc[0], z = abc[1]; 693 if (ctx->inflate) { 694 PetscReal absR, absZ; 695 absR = PetscAbs(r); 696 absZ = PetscAbs(z); 697 CircleInflate(ctx->i_radius,ctx->e_radius,ctx->radius,ctx->num_sections,absR,absZ,&absR,&absZ); 698 r = (r > 0) ? absR : -absR; 699 z = (z > 0) ? absZ : -absZ; 700 } 701 xyz[0] = r; 702 xyz[1] = z; 703 if (dim==3) xyz[2] = abc[2]; 704 705 PetscFunctionReturn(0); 706 } 707 708 static PetscErrorCode ErrorIndicator_Simple(PetscInt dim, PetscReal volume, PetscReal x[], PetscInt Nc, const PetscInt Nf[], const PetscScalar u[], const PetscScalar u_x[], PetscReal *error, void *actx) 709 { 710 PetscReal err = 0.0; 711 PetscInt f = *(PetscInt*)actx, j; 712 PetscFunctionBegin; 713 for (j = 0; j < dim; ++j) { 714 err += PetscSqr(PetscRealPart(u_x[f*dim+j])); 715 } 716 err = PetscRealPart(u[f]); /* just use rho */ 717 *error = volume * err; /* * (ctx->axisymmetric ? 2.*PETSC_PI * r : 1); */ 718 PetscFunctionReturn(0); 719 } 720 721 static PetscErrorCode LandauDMCreateVMesh(MPI_Comm comm, const PetscInt dim, const char prefix[], LandauCtx *ctx, DM *dm) 722 { 723 PetscErrorCode ierr; 724 PetscReal radius = ctx->radius; 725 size_t len; 726 char fname[128] = ""; /* we can add a file if we want */ 727 728 PetscFunctionBegin; 729 /* create DM */ 730 ierr = PetscStrlen(fname, &len);CHKERRQ(ierr); 731 if (len) { 732 PetscInt dim2; 733 ierr = DMPlexCreateFromFile(comm, fname, ctx->interpolate, dm);CHKERRQ(ierr); 734 ierr = DMGetDimension(*dm, &dim2);CHKERRQ(ierr); 735 if (LANDAU_DIM != dim2) SETERRQ2(comm, PETSC_ERR_PLIB, "dim %D != LANDAU_DIM %d",dim2,LANDAU_DIM); 736 } else { /* p4est, quads */ 737 /* Create plex mesh of Landau domain */ 738 if (!ctx->sphere) { 739 PetscInt cells[] = {2,2,2}; 740 PetscReal lo[] = {-radius,-radius,-radius}, hi[] = {radius,radius,radius}; 741 DMBoundaryType periodicity[3] = {DM_BOUNDARY_NONE, dim==2 ? DM_BOUNDARY_NONE : DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 742 if (dim==2) { lo[0] = 0; cells[0] = 1; } 743 ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, cells, lo, hi, periodicity, PETSC_TRUE, dm);CHKERRQ(ierr); 744 ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); /* needed for periodic */ 745 if (dim==3) {ierr = PetscObjectSetName((PetscObject) *dm, "cube");CHKERRQ(ierr);} 746 else {ierr = PetscObjectSetName((PetscObject) *dm, "half-plane");CHKERRQ(ierr);} 747 } else if (dim==2) { 748 PetscInt numCells,cells[16][4],i,j; 749 PetscInt numVerts; 750 PetscReal inner_radius1 = ctx->i_radius, inner_radius2 = ctx->e_radius; 751 PetscReal *flatCoords = NULL; 752 PetscInt *flatCells = NULL, *pcell; 753 if (ctx->num_sections==2) { 754 #if 1 755 numCells = 5; 756 numVerts = 10; 757 int cells2[][4] = { {0,1,4,3}, 758 {1,2,5,4}, 759 {3,4,7,6}, 760 {4,5,8,7}, 761 {6,7,8,9} }; 762 for (i = 0; i < numCells; i++) for (j = 0; j < 4; j++) cells[i][j] = cells2[i][j]; 763 ierr = PetscMalloc2(numVerts * 2, &flatCoords, numCells * 4, &flatCells);CHKERRQ(ierr); 764 { 765 PetscReal (*coords)[2] = (PetscReal (*) [2]) flatCoords; 766 for (j = 0; j < numVerts-1; j++) { 767 PetscReal z, r, theta = -PETSC_PI/2 + (j%3) * PETSC_PI/2; 768 PetscReal rad = (j >= 6) ? inner_radius1 : (j >= 3) ? inner_radius2 : ctx->radius; 769 z = rad * PetscSinReal(theta); 770 coords[j][1] = z; 771 r = rad * PetscCosReal(theta); 772 coords[j][0] = r; 773 } 774 coords[numVerts-1][0] = coords[numVerts-1][1] = 0; 775 } 776 #else 777 numCells = 4; 778 numVerts = 8; 779 static int cells2[][4] = {{0,1,2,3}, 780 {4,5,1,0}, 781 {5,6,2,1}, 782 {6,7,3,2}}; 783 for (i = 0; i < numCells; i++) for (j = 0; j < 4; j++) cells[i][j] = cells2[i][j]; 784 ierr = loc2(numVerts * 2, &flatCoords, numCells * 4, &flatCells);CHKERRQ(ierr); 785 { 786 PetscReal (*coords)[2] = (PetscReal (*) [2]) flatCoords; 787 PetscInt j; 788 for (j = 0; j < 8; j++) { 789 PetscReal z, r; 790 PetscReal theta = -PETSC_PI/2 + (j%4) * PETSC_PI/3.; 791 PetscReal rad = ctx->radius * ((j < 4) ? 0.5 : 1.0); 792 z = rad * PetscSinReal(theta); 793 coords[j][1] = z; 794 r = rad * PetscCosReal(theta); 795 coords[j][0] = r; 796 } 797 } 798 #endif 799 } else if (ctx->num_sections==3) { 800 numCells = 7; 801 numVerts = 12; 802 int cells2[][4] = { {0,1,5,4}, 803 {1,2,6,5}, 804 {2,3,7,6}, 805 {4,5,9,8}, 806 {5,6,10,9}, 807 {6,7,11,10}, 808 {8,9,10,11} }; 809 for (i = 0; i < numCells; i++) for (j = 0; j < 4; j++) cells[i][j] = cells2[i][j]; 810 ierr = PetscMalloc2(numVerts * 2, &flatCoords, numCells * 4, &flatCells);CHKERRQ(ierr); 811 { 812 PetscReal (*coords)[2] = (PetscReal (*) [2]) flatCoords; 813 for (j = 0; j < numVerts; j++) { 814 PetscReal z, r, theta = -PETSC_PI/2 + (j%4) * PETSC_PI/3; 815 PetscReal rad = (j >= 8) ? inner_radius1 : (j >= 4) ? inner_radius2 : ctx->radius; 816 z = rad * PetscSinReal(theta); 817 coords[j][1] = z; 818 r = rad * PetscCosReal(theta); 819 coords[j][0] = r; 820 } 821 } 822 } else if (ctx->num_sections==4) { 823 numCells = 10; 824 numVerts = 16; 825 int cells2[][4] = { {0,1,6,5}, 826 {1,2,7,6}, 827 {2,3,8,7}, 828 {3,4,9,8}, 829 {5,6,11,10}, 830 {6,7,12,11}, 831 {7,8,13,12}, 832 {8,9,14,13}, 833 {10,11,12,15}, 834 {12,13,14,15}}; 835 for (i = 0; i < numCells; i++) for (j = 0; j < 4; j++) cells[i][j] = cells2[i][j]; 836 ierr = PetscMalloc2(numVerts * 2, &flatCoords, numCells * 4, &flatCells);CHKERRQ(ierr); 837 { 838 PetscReal (*coords)[2] = (PetscReal (*) [2]) flatCoords; 839 for (j = 0; j < numVerts-1; j++) { 840 PetscReal z, r, theta = -PETSC_PI/2 + (j%5) * PETSC_PI/4; 841 PetscReal rad = (j >= 10) ? inner_radius1 : (j >= 5) ? inner_radius2 : ctx->radius; 842 z = rad * PetscSinReal(theta); 843 coords[j][1] = z; 844 r = rad * PetscCosReal(theta); 845 coords[j][0] = r; 846 } 847 coords[numVerts-1][0] = coords[numVerts-1][1] = 0; 848 } 849 } else { 850 numCells = 0; 851 numVerts = 0; 852 } 853 for (j = 0, pcell = flatCells; j < numCells; j++, pcell += 4) { 854 pcell[0] = cells[j][0]; pcell[1] = cells[j][1]; 855 pcell[2] = cells[j][2]; pcell[3] = cells[j][3]; 856 } 857 ierr = DMPlexCreateFromCellListPetsc(comm,2,numCells,numVerts,4,ctx->interpolate,flatCells,2,flatCoords,dm);CHKERRQ(ierr); 858 ierr = PetscFree2(flatCoords,flatCells);CHKERRQ(ierr); 859 ierr = PetscObjectSetName((PetscObject) *dm, "semi-circle");CHKERRQ(ierr); 860 } else SETERRQ(ctx->comm, PETSC_ERR_PLIB, "Velocity space meshes does not support cubed sphere"); 861 } 862 ierr = PetscObjectSetOptionsPrefix((PetscObject)*dm,prefix);CHKERRQ(ierr); 863 864 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); /* Plex refine */ 865 866 { /* p4est? */ 867 char convType[256]; 868 PetscBool flg; 869 ierr = PetscOptionsBegin(ctx->comm, prefix, "Mesh conversion options", "DMPLEX");CHKERRQ(ierr); 870 ierr = PetscOptionsFList("-dm_landau_type","Convert DMPlex to another format (should not be Plex!)","plexland.c",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr); 871 ierr = PetscOptionsEnd(); 872 if (flg) { 873 DM dmforest; 874 ierr = DMConvert(*dm,convType,&dmforest);CHKERRQ(ierr); 875 if (dmforest) { 876 PetscBool isForest; 877 if (dmforest->prealloc_only != (*dm)->prealloc_only) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"plex->prealloc_only != dm->prealloc_only"); 878 ierr = PetscObjectSetOptionsPrefix((PetscObject)dmforest,prefix);CHKERRQ(ierr); 879 ierr = DMIsForest(dmforest,&isForest);CHKERRQ(ierr); 880 if (isForest) { 881 if (ctx->sphere && ctx->inflate) { 882 ierr = DMForestSetBaseCoordinateMapping(dmforest,GeometryDMLandau,ctx);CHKERRQ(ierr); 883 } 884 if (dmforest->prealloc_only != (*dm)->prealloc_only) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"plex->prealloc_only != dm->prealloc_only"); 885 ierr = DMDestroy(dm);CHKERRQ(ierr); 886 *dm = dmforest; 887 ctx->errorIndicator = ErrorIndicator_Simple; /* flag for Forest */ 888 } else SETERRQ(ctx->comm, PETSC_ERR_USER, "Converted to non Forest?"); 889 } else SETERRQ(ctx->comm, PETSC_ERR_USER, "Convert failed?"); 890 } 891 } 892 ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr); 893 PetscFunctionReturn(0); 894 } 895 896 static PetscErrorCode SetupDS(DM dm, PetscInt dim, LandauCtx *ctx) 897 { 898 PetscErrorCode ierr; 899 PetscInt ii; 900 PetscFunctionBegin; 901 for (ii=0;ii<ctx->num_species;ii++) { 902 char buf[256]; 903 if (ii==0) ierr = PetscSNPrintf(buf, 256, "e"); 904 else {ierr = PetscSNPrintf(buf, 256, "i%D", ii);CHKERRQ(ierr);} 905 /* Setup Discretization - FEM */ 906 ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, PETSC_FALSE, NULL, PETSC_DECIDE, &ctx->fe[ii]);CHKERRQ(ierr); 907 ierr = PetscObjectSetName((PetscObject) ctx->fe[ii], buf);CHKERRQ(ierr); 908 ierr = DMSetField(dm, ii, NULL, (PetscObject) ctx->fe[ii]);CHKERRQ(ierr); 909 } 910 ierr = DMCreateDS(dm);CHKERRQ(ierr); 911 if (1) { 912 PetscInt ii; 913 PetscSection section; 914 ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 915 for (ii=0;ii<ctx->num_species;ii++){ 916 char buf[256]; 917 if (ii==0) ierr = PetscSNPrintf(buf, 256, "se"); 918 else ierr = PetscSNPrintf(buf, 256, "si%D", ii); 919 ierr = PetscSectionSetComponentName(section, ii, 0, buf);CHKERRQ(ierr); 920 } 921 } 922 PetscFunctionReturn(0); 923 } 924 925 /* Define a Maxwellian function for testing out the operator. */ 926 927 /* Using cartesian velocity space coordinates, the particle */ 928 /* density, [1/m^3], is defined according to */ 929 930 /* $$ n=\int_{R^3} dv^3 \left(\frac{m}{2\pi T}\right)^{3/2}\exp [- mv^2/(2T)] $$ */ 931 932 /* Using some constant, c, we normalize the velocity vector into a */ 933 /* dimensionless variable according to v=c*x. Thus the density, $n$, becomes */ 934 935 /* $$ n=\int_{R^3} dx^3 \left(\frac{mc^2}{2\pi T}\right)^{3/2}\exp [- mc^2/(2T)*x^2] $$ */ 936 937 /* Defining $\theta=2T/mc^2$, we thus find that the probability density */ 938 /* for finding the particle within the interval in a box dx^3 around x is */ 939 940 /* f(x;\theta)=\left(\frac{1}{\pi\theta}\right)^{3/2} \exp [ -x^2/\theta ] */ 941 942 typedef struct { 943 LandauCtx *ctx; 944 PetscReal kT_m; 945 PetscReal n; 946 PetscReal shift; 947 } MaxwellianCtx; 948 949 static PetscErrorCode maxwellian(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx) 950 { 951 MaxwellianCtx *mctx = (MaxwellianCtx*)actx; 952 LandauCtx *ctx = mctx->ctx; 953 PetscInt i; 954 PetscReal v2 = 0, theta = 2*mctx->kT_m/(ctx->v_0*ctx->v_0); /* theta = 2kT/mc^2 */ 955 PetscFunctionBegin; 956 /* compute the exponents, v^2 */ 957 for (i = 0; i < dim; ++i) v2 += x[i]*x[i]; 958 /* evaluate the Maxwellian */ 959 u[0] = mctx->n*PetscPowReal(PETSC_PI*theta,-1.5)*(PetscExpReal(-v2/theta)); 960 if (mctx->shift!=0.) { 961 v2 = 0; 962 for (i = 0; i < dim-1; ++i) v2 += x[i]*x[i]; 963 v2 += (x[dim-1]-mctx->shift)*(x[dim-1]-mctx->shift); 964 /* evaluate the shifted Maxwellian */ 965 u[0] += mctx->n*PetscPowReal(PETSC_PI*theta,-1.5)*(PetscExpReal(-v2/theta)); 966 } 967 PetscFunctionReturn(0); 968 } 969 970 /*@ 971 LandauAddMaxwellians - Add a Maxwellian distribution to a state 972 973 Collective on X 974 975 Input Parameters: 976 . dm - The mesh 977 + time - Current time 978 - temps - Temperatures of each species 979 . ns - Number density of each species 980 + actx - Landau context 981 982 Output Parameter: 983 . X - The state 984 985 Level: beginner 986 987 .keywords: mesh 988 .seealso: LandauCreateVelocitySpace() 989 @*/ 990 PetscErrorCode LandauAddMaxwellians(DM dm, Vec X, PetscReal time, PetscReal temps[], PetscReal ns[], void *actx) 991 { 992 LandauCtx *ctx = (LandauCtx*)actx; 993 PetscErrorCode (*initu[LANDAU_MAX_SPECIES])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar [], void *); 994 PetscErrorCode ierr,ii; 995 PetscInt dim; 996 MaxwellianCtx *mctxs[LANDAU_MAX_SPECIES], data[LANDAU_MAX_SPECIES]; 997 998 PetscFunctionBegin; 999 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1000 if (!ctx) { ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr); } 1001 for (ii=0;ii<ctx->num_species;ii++) { 1002 mctxs[ii] = &data[ii]; 1003 data[ii].ctx = ctx; 1004 data[ii].kT_m = ctx->k*temps[ii]/ctx->masses[ii]; /* kT/m */ 1005 data[ii].n = ns[ii]; 1006 initu[ii] = maxwellian; 1007 data[ii].shift = 0; 1008 } 1009 data[0].shift = ctx->electronShift; 1010 /* need to make ADD_ALL_VALUES work - TODO */ 1011 ierr = DMProjectFunction(dm, time, initu, (void**)mctxs, INSERT_ALL_VALUES, X);CHKERRQ(ierr); 1012 PetscFunctionReturn(0); 1013 } 1014 1015 /* 1016 LandauSetInitialCondition - Addes Maxwellians with context 1017 1018 Collective on X 1019 1020 Input Parameters: 1021 . dm - The mesh 1022 + actx - Landau context with T and n 1023 1024 Output Parameter: 1025 . X - The state 1026 1027 Level: beginner 1028 1029 .keywords: mesh 1030 .seealso: LandauCreateVelocitySpace(), LandauAddMaxwellians() 1031 */ 1032 static PetscErrorCode LandauSetInitialCondition(DM dm, Vec X, void *actx) 1033 { 1034 LandauCtx *ctx = (LandauCtx*)actx; 1035 PetscErrorCode ierr; 1036 PetscFunctionBegin; 1037 if (!ctx) { ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr); } 1038 ierr = VecZeroEntries(X);CHKERRQ(ierr); 1039 ierr = LandauAddMaxwellians(dm, X, 0.0, ctx->thermal_temps, ctx->n, ctx);CHKERRQ(ierr); 1040 PetscFunctionReturn(0); 1041 } 1042 1043 static PetscErrorCode adaptToleranceFEM(PetscFE fem, Vec sol, PetscReal refineTol[], PetscReal coarsenTol[], PetscInt type, LandauCtx *ctx, DM *newDM) 1044 { 1045 DM dm, plex, adaptedDM = NULL; 1046 PetscDS prob; 1047 PetscBool isForest; 1048 PetscQuadrature quad; 1049 PetscInt Nq, *Nb, cStart, cEnd, c, dim, qj, k; 1050 DMLabel adaptLabel = NULL; 1051 PetscErrorCode ierr; 1052 1053 PetscFunctionBegin; 1054 ierr = VecGetDM(sol, &dm);CHKERRQ(ierr); 1055 ierr = DMCreateDS(dm);CHKERRQ(ierr); 1056 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1057 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1058 ierr = DMIsForest(dm, &isForest);CHKERRQ(ierr); 1059 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 1060 ierr = DMPlexGetHeightStratum(plex,0,&cStart,&cEnd);CHKERRQ(ierr); 1061 ierr = DMLabelCreate(PETSC_COMM_SELF,"adapt",&adaptLabel);CHKERRQ(ierr); 1062 ierr = PetscFEGetQuadrature(fem, &quad);CHKERRQ(ierr); 1063 ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 1064 if (Nq >LANDAU_MAX_NQ) SETERRQ2(ctx->comm,PETSC_ERR_ARG_WRONG,"Order too high. Nq = %D > LANDAU_MAX_NQ (%D)",Nq,LANDAU_MAX_NQ); 1065 ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 1066 if (type==4) { 1067 for (c = cStart; c < cEnd; c++) { 1068 ierr = DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE);CHKERRQ(ierr); 1069 } 1070 ierr = PetscInfo1(sol, "Phase:%s: Uniform refinement\n","adaptToleranceFEM");CHKERRQ(ierr); 1071 } else if (type==2) { 1072 PetscInt rCellIdx[8], eCellIdx[64], iCellIdx[64], eMaxIdx = -1, iMaxIdx = -1, nr = 0, nrmax = (dim==3) ? 8 : 2; 1073 PetscReal minRad = PETSC_INFINITY, r, eMinRad = PETSC_INFINITY, iMinRad = PETSC_INFINITY; 1074 for (c = 0; c < 64; c++) { eCellIdx[c] = iCellIdx[c] = -1; } 1075 for (c = cStart; c < cEnd; c++) { 1076 PetscReal tt, v0[LANDAU_MAX_NQ*3], detJ[LANDAU_MAX_NQ]; 1077 ierr = DMPlexComputeCellGeometryFEM(plex, c, quad, v0, NULL, NULL, detJ);CHKERRQ(ierr); 1078 for (qj = 0; qj < Nq; ++qj) { 1079 tt = PetscSqr(v0[dim*qj+0]) + PetscSqr(v0[dim*qj+1]) + PetscSqr(((dim==3) ? v0[dim*qj+2] : 0)); 1080 r = PetscSqrtReal(tt); 1081 if (r < minRad - PETSC_SQRT_MACHINE_EPSILON*10.) { 1082 minRad = r; 1083 nr = 0; 1084 rCellIdx[nr++]= c; 1085 ierr = PetscInfo4(sol, "\t\tPhase: adaptToleranceFEM Found first inner r=%e, cell %D, qp %D/%D\n", r, c, qj+1, Nq);CHKERRQ(ierr); 1086 } else if ((r-minRad) < PETSC_SQRT_MACHINE_EPSILON*100. && nr < nrmax) { 1087 for (k=0;k<nr;k++) if (c == rCellIdx[k]) break; 1088 if (k==nr) { 1089 rCellIdx[nr++]= c; 1090 ierr = PetscInfo5(sol, "\t\t\tPhase: adaptToleranceFEM Found another inner r=%e, cell %D, qp %D/%D, d=%e\n", r, c, qj+1, Nq, r-minRad);CHKERRQ(ierr); 1091 } 1092 } 1093 if (ctx->sphere) { 1094 if ((tt=r-ctx->e_radius) > 0) { 1095 PetscInfo2(sol, "\t\t\t %D cell r=%g\n",c,tt); 1096 if (tt < eMinRad - PETSC_SQRT_MACHINE_EPSILON*100.) { 1097 eMinRad = tt; 1098 eMaxIdx = 0; 1099 eCellIdx[eMaxIdx++] = c; 1100 } else if (eMaxIdx > 0 && (tt-eMinRad) <= PETSC_SQRT_MACHINE_EPSILON && c != eCellIdx[eMaxIdx-1]) { 1101 eCellIdx[eMaxIdx++] = c; 1102 } 1103 } 1104 if ((tt=r-ctx->i_radius) > 0) { 1105 if (tt < iMinRad - 1.e-5) { 1106 iMinRad = tt; 1107 iMaxIdx = 0; 1108 iCellIdx[iMaxIdx++] = c; 1109 } else if (iMaxIdx > 0 && (tt-iMinRad) <= PETSC_SQRT_MACHINE_EPSILON && c != iCellIdx[iMaxIdx-1]) { 1110 iCellIdx[iMaxIdx++] = c; 1111 } 1112 } 1113 } 1114 } 1115 } 1116 for (k=0;k<nr;k++) { 1117 ierr = DMLabelSetValue(adaptLabel, rCellIdx[k], DM_ADAPT_REFINE);CHKERRQ(ierr); 1118 } 1119 if (ctx->sphere) { 1120 for (c = 0; c < eMaxIdx; c++) { 1121 ierr = DMLabelSetValue(adaptLabel, eCellIdx[c], DM_ADAPT_REFINE);CHKERRQ(ierr); 1122 ierr = PetscInfo3(sol, "\t\tPhase:%s: refine sphere e cell %D r=%g\n","adaptToleranceFEM",eCellIdx[c],eMinRad); 1123 } 1124 for (c = 0; c < iMaxIdx; c++) { 1125 ierr = DMLabelSetValue(adaptLabel, iCellIdx[c], DM_ADAPT_REFINE);CHKERRQ(ierr); 1126 ierr = PetscInfo3(sol, "\t\tPhase:%s: refine sphere i cell %D r=%g\n","adaptToleranceFEM",iCellIdx[c],iMinRad); 1127 } 1128 } 1129 ierr = PetscInfo4(sol, "Phase:%s: Adaptive refine origin cells %D,%D r=%g\n","adaptToleranceFEM",rCellIdx[0],rCellIdx[1],minRad); 1130 } else if (type==0 || type==1 || type==3) { /* refine along r=0 axis */ 1131 PetscScalar *coef = NULL; 1132 Vec coords; 1133 PetscInt csize,Nv,d,nz; 1134 DM cdm; 1135 PetscSection cs; 1136 ierr = DMGetCoordinatesLocal(dm, &coords);CHKERRQ(ierr); 1137 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 1138 ierr = DMGetLocalSection(cdm, &cs);CHKERRQ(ierr); 1139 for (c = cStart; c < cEnd; c++) { 1140 PetscInt doit = 0, outside = 0; 1141 ierr = DMPlexVecGetClosure(cdm, cs, coords, c, &csize, &coef);CHKERRQ(ierr); 1142 Nv = csize/dim; 1143 for (nz = d = 0; d < Nv; d++) { 1144 PetscReal z = PetscRealPart(coef[d*dim + (dim-1)]), x = PetscSqr(PetscRealPart(coef[d*dim + 0])) + ((dim==3) ? PetscSqr(PetscRealPart(coef[d*dim + 1])) : 0); 1145 x = PetscSqrtReal(x); 1146 if (x < PETSC_MACHINE_EPSILON*10. && PetscAbs(z)<PETSC_MACHINE_EPSILON*10.) doit = 1; /* refine origin */ 1147 else if (type==0 && (z < -PETSC_MACHINE_EPSILON*10. || z > ctx->re_radius+PETSC_MACHINE_EPSILON*10.)) outside++; /* first pass don't refine bottom */ 1148 else if (type==1 && (z > ctx->vperp0_radius1 || z < -ctx->vperp0_radius1)) outside++; /* don't refine outside electron refine radius */ 1149 else if (type==3 && (z > ctx->vperp0_radius2 || z < -ctx->vperp0_radius2)) outside++; /* don't refine outside ion refine radius */ 1150 if (x < PETSC_MACHINE_EPSILON*10.) nz++; 1151 } 1152 ierr = DMPlexVecRestoreClosure(cdm, cs, coords, c, &csize, &coef);CHKERRQ(ierr); 1153 if (doit || (outside<Nv && nz)) { 1154 ierr = DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE);CHKERRQ(ierr); 1155 } 1156 } 1157 ierr = PetscInfo1(sol, "Phase:%s: RE refinement\n","adaptToleranceFEM"); 1158 } 1159 ierr = DMDestroy(&plex);CHKERRQ(ierr); 1160 ierr = DMAdaptLabel(dm, adaptLabel, &adaptedDM);CHKERRQ(ierr); 1161 ierr = DMLabelDestroy(&adaptLabel);CHKERRQ(ierr); 1162 *newDM = adaptedDM; 1163 if (adaptedDM) { 1164 if (isForest) { 1165 ierr = DMForestSetAdaptivityForest(adaptedDM,NULL);CHKERRQ(ierr); 1166 } 1167 ierr = DMConvert(adaptedDM, DMPLEX, &plex);CHKERRQ(ierr); 1168 ierr = DMPlexGetHeightStratum(plex,0,&cStart,&cEnd);CHKERRQ(ierr); 1169 ierr = PetscInfo2(sol, "\tPhase: adaptToleranceFEM: %D cells, %d total quadrature points\n",cEnd-cStart,Nq*(cEnd-cStart));CHKERRQ(ierr); 1170 ierr = DMDestroy(&plex);CHKERRQ(ierr); 1171 } 1172 PetscFunctionReturn(0); 1173 } 1174 1175 static PetscErrorCode adapt(DM *dm, LandauCtx *ctx, Vec *uu) 1176 { 1177 PetscErrorCode ierr; 1178 PetscInt type, limits[5] = {ctx->numRERefine,ctx->nZRefine1,ctx->maxRefIts,ctx->nZRefine2,ctx->postAMRRefine}; 1179 PetscInt adaptIter; 1180 1181 PetscFunctionBegin; 1182 for (type=0;type<5;type++) { 1183 for (adaptIter = 0; adaptIter<limits[type];adaptIter++) { 1184 DM dmNew = NULL; 1185 ierr = adaptToleranceFEM(ctx->fe[0], *uu, ctx->refineTol, ctx->coarsenTol, type, ctx, &dmNew);CHKERRQ(ierr); 1186 if (!dmNew) SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"should not happen"); 1187 else { 1188 ierr = DMDestroy(dm);CHKERRQ(ierr); 1189 ierr = VecDestroy(uu);CHKERRQ(ierr); 1190 ierr = DMCreateGlobalVector(dmNew,uu);CHKERRQ(ierr); 1191 ierr = PetscObjectSetName((PetscObject) *uu, "u");CHKERRQ(ierr); 1192 ierr = LandauSetInitialCondition(dmNew, *uu, ctx);CHKERRQ(ierr); 1193 *dm = dmNew; 1194 } 1195 } 1196 } 1197 PetscFunctionReturn(0); 1198 } 1199 1200 static PetscErrorCode ProcessOptions(LandauCtx *ctx, const char prefix[]) 1201 { 1202 PetscErrorCode ierr; 1203 PetscBool flg, sph_flg; 1204 PetscInt ii,nt,nm,nc; 1205 DM dummy; 1206 1207 PetscFunctionBegin; 1208 ierr = DMCreate(ctx->comm,&dummy);CHKERRQ(ierr); 1209 /* get options - initialize context */ 1210 ctx->verbose = 1; 1211 ctx->interpolate = PETSC_TRUE; 1212 ctx->gpu_assembly = PETSC_TRUE; 1213 ctx->sphere = PETSC_FALSE; 1214 ctx->inflate = PETSC_FALSE; 1215 ctx->electronShift = 0; 1216 ctx->errorIndicator = NULL; 1217 ctx->radius = 5.; /* electron thermal radius (velocity) */ 1218 ctx->re_radius = 0.; 1219 ctx->vperp0_radius1 = 0; 1220 ctx->vperp0_radius2 = 0; 1221 ctx->e_radius = .1; 1222 ctx->i_radius = .01; 1223 ctx->maxRefIts = 5; 1224 ctx->postAMRRefine = 0; 1225 ctx->nZRefine1 = 0; 1226 ctx->nZRefine2 = 0; 1227 ctx->numRERefine = 0; 1228 ctx->aux_bool = PETSC_FALSE; 1229 ctx->num_sections = 3; /* 2, 3 or 4 */ 1230 /* species - [0] electrons, [1] one ion species eg, duetarium, [2] heavy impurity ion, ... */ 1231 ctx->charges[0] = -1; /* electron charge (MKS) */ 1232 ctx->masses[0] = 1/1835.5; /* temporary value in proton mass */ 1233 ctx->n[0] = 1; 1234 ctx->thermal_temps[0] = 1; 1235 /* constants, etc. */ 1236 ctx->epsilon0 = 8.8542e-12; /* permittivity of free space (MKS) F/m */ 1237 ctx->k = 1.38064852e-23; /* Boltzmann constant (MKS) J/K */ 1238 ctx->lnLam = 10; /* cross section ratio large - small angle collisions */ 1239 ctx->n_0 = 1.e20; /* typical plasma n, but could set it to 1 */ 1240 ctx->Ez = 0; 1241 ctx->v_0 = 1; /* in electron thermal velocity */ 1242 ctx->subThreadBlockSize = 1; /* for device and maybe OMP */ 1243 ctx->numConcurrency = 1; /* for device */ 1244 ctx->SData_d = NULL; /* for device */ 1245 ctx->times[0] = 0; 1246 ctx->init = PETSC_FALSE; // doit first time 1247 ctx->use_matrix_mass = PETSC_FALSE; /* fast but slightly fragile */ 1248 ctx->plex = NULL; /* cache as expensive to Convert */ 1249 ierr = PetscOptionsBegin(ctx->comm, prefix, "Options for Fokker-Plank-Landau collision operator", "none");CHKERRQ(ierr); 1250 { 1251 char opstring[256]; 1252 #if defined(PETSC_HAVE_KOKKOS) 1253 ctx->deviceType = LANDAU_KOKKOS; 1254 ierr = PetscStrcpy(opstring,"kokkos");CHKERRQ(ierr); 1255 #if defined(PETSC_HAVE_CUDA) 1256 ctx->subThreadBlockSize = 16; 1257 #endif 1258 #elif defined(PETSC_HAVE_CUDA) 1259 ctx->deviceType = LANDAU_CUDA; 1260 ierr = PetscStrcpy(opstring,"cuda");CHKERRQ(ierr); 1261 #else 1262 ctx->deviceType = LANDAU_CPU; 1263 ierr = PetscStrcpy(opstring,"cpu");CHKERRQ(ierr); 1264 ctx->subThreadBlockSize = 0; 1265 #endif 1266 ierr = PetscOptionsString("-dm_landau_device_type","Use kernels on 'cpu', 'cuda', or 'kokkos'","plexland.c",opstring,opstring,256,NULL);CHKERRQ(ierr); 1267 ierr = PetscStrcmp("cpu",opstring,&flg);CHKERRQ(ierr); 1268 if (flg) { 1269 ctx->deviceType = LANDAU_CPU; 1270 ctx->subThreadBlockSize = 0; 1271 } else { 1272 ierr = PetscStrcmp("cuda",opstring,&flg);CHKERRQ(ierr); 1273 if (flg) { 1274 ctx->deviceType = LANDAU_CUDA; 1275 ctx->subThreadBlockSize = 0; 1276 } else { 1277 ierr = PetscStrcmp("kokkos",opstring,&flg);CHKERRQ(ierr); 1278 if (flg) ctx->deviceType = LANDAU_KOKKOS; 1279 else SETERRQ1(ctx->comm,PETSC_ERR_ARG_WRONG,"-dm_landau_device_type %s",opstring); 1280 } 1281 } 1282 } 1283 ierr = PetscOptionsBool("-dm_landau_gpu_assembly", "Assemble Jacobian on GPU", "plexland.c", ctx->gpu_assembly, &ctx->gpu_assembly, NULL);CHKERRQ(ierr); 1284 ierr = PetscOptionsReal("-dm_landau_electron_shift","Shift in thermal velocity of electrons","none",ctx->electronShift,&ctx->electronShift, NULL);CHKERRQ(ierr); 1285 ierr = PetscOptionsBool("-dm_landau_sphere", "use sphere/semi-circle domain instead of rectangle", "plexland.c", ctx->sphere, &ctx->sphere, &sph_flg);CHKERRQ(ierr); 1286 ierr = PetscOptionsBool("-dm_landau_inflate", "With sphere, inflate for curved edges (no AMR)", "plexland.c", ctx->inflate, &ctx->inflate, NULL);CHKERRQ(ierr); 1287 ierr = PetscOptionsInt("-dm_landau_amr_re_levels", "Number of levels to refine along v_perp=0, z>0", "plexland.c", ctx->numRERefine, &ctx->numRERefine, NULL);CHKERRQ(ierr); 1288 ierr = PetscOptionsInt("-dm_landau_amr_z_refine1", "Number of levels to refine along v_perp=0", "plexland.c", ctx->nZRefine1, &ctx->nZRefine1, NULL);CHKERRQ(ierr); 1289 ierr = PetscOptionsInt("-dm_landau_amr_z_refine2", "Number of levels to refine along v_perp=0", "plexland.c", ctx->nZRefine2, &ctx->nZRefine2, NULL);CHKERRQ(ierr); 1290 ierr = PetscOptionsInt("-dm_landau_amr_levels_max", "Number of AMR levels of refinement around origin after r=0 refinements", "plexland.c", ctx->maxRefIts, &ctx->maxRefIts, NULL);CHKERRQ(ierr); 1291 ierr = PetscOptionsInt("-dm_landau_amr_post_refine", "Number of levels to uniformly refine after AMR", "plexland.c", ctx->postAMRRefine, &ctx->postAMRRefine, NULL);CHKERRQ(ierr); 1292 ierr = PetscOptionsInt("-dm_landau_verbose", "", "plexland.c", ctx->verbose, &ctx->verbose, NULL);CHKERRQ(ierr); 1293 ierr = 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);CHKERRQ(ierr); 1294 ierr = PetscOptionsReal("-dm_landau_z_radius1","velocity range to refine r=0 axis (for electrons)","plexland.c",ctx->vperp0_radius1,&ctx->vperp0_radius1, &flg);CHKERRQ(ierr); 1295 ierr = PetscOptionsReal("-dm_landau_z_radius2","velocity range to refine r=0 axis (for ions) after origin AMR","plexland.c",ctx->vperp0_radius2,&ctx->vperp0_radius2, &flg);CHKERRQ(ierr); 1296 ierr = PetscOptionsReal("-dm_landau_Ez","Initial parallel electric field in unites of Conner-Hastie criticle field","plexland.c",ctx->Ez,&ctx->Ez, NULL);CHKERRQ(ierr); 1297 ierr = PetscOptionsReal("-dm_landau_n_0","Normalization constant for number density","plexland.c",ctx->n_0,&ctx->n_0, NULL);CHKERRQ(ierr); 1298 ierr = PetscOptionsReal("-dm_landau_ln_lambda","Cross section parameter","plexland.c",ctx->lnLam,&ctx->lnLam, NULL);CHKERRQ(ierr); 1299 ierr = PetscOptionsInt("-dm_landau_num_sections", "Number of tangential section in (2D) grid, 2, 3, of 4", "plexland.c", ctx->num_sections, &ctx->num_sections, NULL);CHKERRQ(ierr); 1300 ierr = PetscOptionsInt("-dm_landau_num_thread_teams", "The number of other concurrent runs to make room for", "plexland.c", ctx->numConcurrency, &ctx->numConcurrency, NULL);CHKERRQ(ierr); 1301 ierr = 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);CHKERRQ(ierr); 1302 1303 /* get num species with tempurature*/ 1304 { 1305 PetscReal arr[100]; 1306 nt = 100; 1307 ierr = PetscOptionsRealArray("-dm_landau_thermal_temps", "Temperature of each species [e,i_0,i_1,...] in keV", "plexland.c", arr, &nt, &flg);CHKERRQ(ierr); 1308 if (flg && nt > LANDAU_MAX_SPECIES) SETERRQ2(ctx->comm,PETSC_ERR_ARG_WRONG,"-thermal_temps ,t1,t2,.. number of species %D > MAX %D",nt,LANDAU_MAX_SPECIES); 1309 } 1310 nt = LANDAU_MAX_SPECIES; 1311 for (ii=1;ii<LANDAU_MAX_SPECIES;ii++) { 1312 ctx->thermal_temps[ii] = 1.; 1313 ctx->charges[ii] = 1; 1314 ctx->masses[ii] = 1; 1315 ctx->n[ii] = (ii==1) ? 1 : 0; 1316 } 1317 ierr = 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);CHKERRQ(ierr); 1318 if (flg) { 1319 PetscInfo1(dummy, "num_species set to number of thermal temps provided (%D)\n",nt); 1320 ctx->num_species = nt; 1321 } else SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"-dm_landau_thermal_temps ,t1,t2,.. must be provided to set the number of species"); 1322 for (ii=0;ii<ctx->num_species;ii++) ctx->thermal_temps[ii] *= 1.1604525e7; /* convert to Kelvin */ 1323 nm = LANDAU_MAX_SPECIES-1; 1324 ierr = 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);CHKERRQ(ierr); 1325 if (flg && nm != ctx->num_species-1) { 1326 SETERRQ2(ctx->comm,PETSC_ERR_ARG_WRONG,"num ion masses %D != num species %D",nm,ctx->num_species-1); 1327 } 1328 nm = LANDAU_MAX_SPECIES; 1329 ierr = PetscOptionsRealArray("-dm_landau_n", "Normalized (by -n_0) number density of each species", "plexland.c", ctx->n, &nm, &flg);CHKERRQ(ierr); 1330 if (flg && nm != ctx->num_species) SETERRQ2(ctx->comm,PETSC_ERR_ARG_WRONG,"wrong num n: %D != num species %D",nm,ctx->num_species); 1331 ctx->n_0 *= ctx->n[0]; /* normalized number density */ 1332 for (ii=1;ii<ctx->num_species;ii++) ctx->n[ii] = ctx->n[ii]/ctx->n[0]; 1333 ctx->n[0] = 1; 1334 for (ii=0;ii<LANDAU_MAX_SPECIES;ii++) ctx->masses[ii] *= 1.6720e-27; /* scale by proton mass kg */ 1335 ctx->masses[0] = 9.10938356e-31; /* electron mass kg (should be about right already) */ 1336 ctx->m_0 = ctx->masses[0]; /* arbitrary reference mass, electrons */ 1337 ierr = PetscOptionsReal("-dm_landau_v_0","Velocity to normalize with in units of initial electrons thermal velocity (not recommended to change default)","plexland.c",ctx->v_0,&ctx->v_0, NULL);CHKERRQ(ierr); 1338 ctx->v_0 *= PetscSqrtReal(ctx->k*ctx->thermal_temps[0]/(ctx->masses[0])); /* electron mean velocity in 1D (need 3D form in computing T from FE integral) */ 1339 nc = LANDAU_MAX_SPECIES-1; 1340 ierr = 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);CHKERRQ(ierr); 1341 if (flg && nc != ctx->num_species-1) SETERRQ2(ctx->comm,PETSC_ERR_ARG_WRONG,"num charges %D != num species %D",nc,ctx->num_species-1); 1342 for (ii=0;ii<LANDAU_MAX_SPECIES;ii++) ctx->charges[ii] *= 1.6022e-19; /* electron/proton charge (MKS) */ 1343 ctx->t_0 = 8*PETSC_PI*PetscSqr(ctx->epsilon0*ctx->m_0/PetscSqr(ctx->charges[0]))/ctx->lnLam/ctx->n_0*PetscPowReal(ctx->v_0,3); /* note, this t_0 makes nu[0,0]=1 */ 1344 /* geometry */ 1345 for (ii=0;ii<ctx->num_species;ii++) ctx->refineTol[ii] = PETSC_MAX_REAL; 1346 for (ii=0;ii<ctx->num_species;ii++) ctx->coarsenTol[ii] = 0.; 1347 ii = LANDAU_MAX_SPECIES; 1348 ierr = PetscOptionsRealArray("-dm_landau_refine_tol","tolerance for refining cells in AMR","plexland.c",ctx->refineTol, &ii, &flg);CHKERRQ(ierr); 1349 if (flg && ii != ctx->num_species) ierr = PetscInfo2(dummy, "Phase: Warning, #refine_tol %D != num_species %D\n",ii,ctx->num_species);CHKERRQ(ierr); 1350 ii = LANDAU_MAX_SPECIES; 1351 ierr = PetscOptionsRealArray("-dm_landau_coarsen_tol","tolerance for coarsening cells in AMR","plexland.c",ctx->coarsenTol, &ii, &flg);CHKERRQ(ierr); 1352 if (flg && ii != ctx->num_species) ierr = PetscInfo2(dummy, "Phase: Warning, #coarsen_tol %D != num_species %D\n",ii,ctx->num_species);CHKERRQ(ierr); 1353 ierr = PetscOptionsReal("-dm_landau_domain_radius","Phase space size in units of electron thermal velocity","plexland.c",ctx->radius,&ctx->radius, &flg);CHKERRQ(ierr); 1354 if (flg && ctx->radius <= 0) { /* negative is ratio of c */ 1355 if (ctx->radius == 0) ctx->radius = 0.75; 1356 else ctx->radius = -ctx->radius; 1357 ctx->radius = ctx->radius*299792458.0/ctx->v_0; 1358 ierr = PetscInfo1(dummy, "Change domain radius to %e\n",ctx->radius);CHKERRQ(ierr); 1359 } 1360 ierr = PetscOptionsReal("-dm_landau_i_radius","Ion thermal velocity, used for circular meshes","plexland.c",ctx->i_radius,&ctx->i_radius, &flg);CHKERRQ(ierr); 1361 if (flg && !sph_flg) ctx->sphere = PETSC_TRUE; /* you gave me an ion radius but did not set sphere, user error really */ 1362 if (!flg) { 1363 ctx->i_radius = 1.5*PetscSqrtReal(8*ctx->k*ctx->thermal_temps[1]/ctx->masses[1]/PETSC_PI)/ctx->v_0; /* normalized radius with thermal velocity of first ion */ 1364 } 1365 ierr = PetscOptionsReal("-dm_landau_e_radius","Electron thermal velocity, used for circular meshes","plexland.c",ctx->e_radius,&ctx->e_radius, &flg);CHKERRQ(ierr); 1366 if (flg && !sph_flg) ctx->sphere = PETSC_TRUE; /* you gave me an e radius but did not set sphere, user error really */ 1367 if (!flg) { 1368 ctx->e_radius = 1.5*PetscSqrtReal(8*ctx->k*ctx->thermal_temps[0]/ctx->masses[0]/PETSC_PI)/ctx->v_0; /* normalized radius with thermal velocity of electrons */ 1369 } 1370 if (ctx->sphere && (ctx->e_radius <= ctx->i_radius || ctx->radius <= ctx->e_radius)) SETERRQ3(ctx->comm,PETSC_ERR_ARG_WRONG,"bad radii: %g < %g < %g",ctx->i_radius,ctx->e_radius,ctx->radius); 1371 ierr = PetscOptionsInt("-dm_landau_sub_thread_block_size", "Number of threads in Kokkos integration point subblock", "plexland.c", ctx->subThreadBlockSize, &ctx->subThreadBlockSize, NULL);CHKERRQ(ierr); 1372 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1373 for (ii=ctx->num_species;ii<LANDAU_MAX_SPECIES;ii++) ctx->masses[ii] = ctx->thermal_temps[ii] = ctx->charges[ii] = 0; 1374 if (ctx->verbose > 0) { 1375 ierr = PetscPrintf(ctx->comm, "masses: e=%10.3e; ions in proton mass units: %10.3e %10.3e ...\n",ctx->masses[0],ctx->masses[1]/1.6720e-27,ctx->num_species>2 ? ctx->masses[2]/1.6720e-27 : 0);CHKERRQ(ierr); 1376 ierr = PetscPrintf(ctx->comm, "charges: e=%10.3e; charges in elementary units: %10.3e %10.3e\n", ctx->charges[0],-ctx->charges[1]/ctx->charges[0],ctx->num_species>2 ? -ctx->charges[2]/ctx->charges[0] : 0);CHKERRQ(ierr); 1377 ierr = PetscPrintf(ctx->comm, "thermal T (K): e=%10.3e i=%10.3e imp=%10.3e. v_0=%10.3e n_0=%10.3e t_0=%10.3e domain=%10.3e\n",ctx->thermal_temps[0],ctx->thermal_temps[1],ctx->num_species>2 ? ctx->thermal_temps[2] : 0,ctx->v_0,ctx->n_0,ctx->t_0,ctx->radius);CHKERRQ(ierr); 1378 } 1379 ierr = DMDestroy(&dummy);CHKERRQ(ierr); 1380 { 1381 PetscMPIInt rank; 1382 ierr = MPI_Comm_rank(ctx->comm, &rank);CHKERRMPI(ierr); 1383 /* PetscLogStage setup_stage; */ 1384 ierr = PetscLogEventRegister("Landau Operator", DM_CLASSID, &ctx->events[11]);CHKERRQ(ierr); /* 11 */ 1385 ierr = PetscLogEventRegister("Landau Jacobian", DM_CLASSID, &ctx->events[0]);CHKERRQ(ierr); /* 0 */ 1386 ierr = PetscLogEventRegister("Landau Mass", DM_CLASSID, &ctx->events[9]);CHKERRQ(ierr); /* 9 */ 1387 ierr = PetscLogEventRegister(" Preamble", DM_CLASSID, &ctx->events[10]);CHKERRQ(ierr); /* 10 */ 1388 ierr = PetscLogEventRegister(" static IP Data", DM_CLASSID, &ctx->events[7]);CHKERRQ(ierr); /* 7 */ 1389 ierr = PetscLogEventRegister(" dynamic IP-Jac", DM_CLASSID, &ctx->events[1]);CHKERRQ(ierr); /* 1 */ 1390 ierr = PetscLogEventRegister(" Kernel-init", DM_CLASSID, &ctx->events[3]);CHKERRQ(ierr); /* 3 */ 1391 ierr = PetscLogEventRegister(" Jac-f-df (GPU)", DM_CLASSID, &ctx->events[8]);CHKERRQ(ierr); /* 8 */ 1392 ierr = PetscLogEventRegister(" Kernel (GPU)", DM_CLASSID, &ctx->events[4]);CHKERRQ(ierr); /* 4 */ 1393 ierr = PetscLogEventRegister(" Copy to CPU", DM_CLASSID, &ctx->events[5]);CHKERRQ(ierr); /* 5 */ 1394 ierr = PetscLogEventRegister(" Jac-assemble", DM_CLASSID, &ctx->events[6]);CHKERRQ(ierr); /* 6 */ 1395 ierr = PetscLogEventRegister(" Jac asmbl setup", DM_CLASSID, &ctx->events[2]);CHKERRQ(ierr); /* 2 */ 1396 ierr = PetscLogEventRegister(" Other", DM_CLASSID, &ctx->events[13]);CHKERRQ(ierr); /* 13 */ 1397 1398 if (rank) { /* turn off output stuff for duplicate runs - do we need to add the prefix to all this? */ 1399 ierr = PetscOptionsClearValue(NULL,"-snes_converged_reason");CHKERRQ(ierr); 1400 ierr = PetscOptionsClearValue(NULL,"-ksp_converged_reason");CHKERRQ(ierr); 1401 ierr = PetscOptionsClearValue(NULL,"-snes_monitor");CHKERRQ(ierr); 1402 ierr = PetscOptionsClearValue(NULL,"-ksp_monitor");CHKERRQ(ierr); 1403 ierr = PetscOptionsClearValue(NULL,"-ts_monitor");CHKERRQ(ierr); 1404 ierr = PetscOptionsClearValue(NULL,"-ts_adapt_monitor");CHKERRQ(ierr); 1405 ierr = PetscOptionsClearValue(NULL,"-dm_landau_amr_dm_view");CHKERRQ(ierr); 1406 ierr = PetscOptionsClearValue(NULL,"-dm_landau_amr_vec_view");CHKERRQ(ierr); 1407 ierr = PetscOptionsClearValue(NULL,"-dm_landau_pre_dm_view");CHKERRQ(ierr); 1408 ierr = PetscOptionsClearValue(NULL,"-dm_landau_pre_vec_view");CHKERRQ(ierr); 1409 ierr = PetscOptionsClearValue(NULL,"-info");CHKERRQ(ierr); 1410 } 1411 } 1412 PetscFunctionReturn(0); 1413 } 1414 1415 /*@C 1416 LandauCreateVelocitySpace - Create a DMPlex velocity space mesh 1417 1418 Collective on comm 1419 1420 Input Parameters: 1421 + comm - The MPI communicator 1422 . dim - velocity space dimension (2 for axisymmetric, 3 for full 3X + 3V solver) 1423 - prefix - prefix for options 1424 1425 Output Parameter: 1426 . dm - The DM object representing the mesh 1427 + X - A vector (user destroys) 1428 - J - Optional matrix (object destroys) 1429 1430 Level: beginner 1431 1432 .keywords: mesh 1433 .seealso: DMPlexCreate(), LandauDestroyVelocitySpace() 1434 @*/ 1435 PetscErrorCode LandauCreateVelocitySpace(MPI_Comm comm, PetscInt dim, const char prefix[], Vec *X, Mat *J, DM *dm) 1436 { 1437 PetscErrorCode ierr; 1438 LandauCtx *ctx; 1439 PetscBool prealloc_only,flg; 1440 1441 PetscFunctionBegin; 1442 if (dim!=2 && dim!=3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Only 2D and 3D supported"); 1443 ierr = PetscNew(&ctx);CHKERRQ(ierr); 1444 ctx->comm = comm; /* used for diagnostics and global errors */ 1445 /* process options */ 1446 ierr = ProcessOptions(ctx,prefix);CHKERRQ(ierr); 1447 /* Create Mesh */ 1448 ierr = LandauDMCreateVMesh(PETSC_COMM_SELF, dim, prefix, ctx, dm);CHKERRQ(ierr); 1449 prealloc_only = (*dm)->prealloc_only; 1450 ierr = DMViewFromOptions(*dm,NULL,"-dm_landau_pre_dm_view");CHKERRQ(ierr); 1451 ierr = DMSetApplicationContext(*dm, ctx);CHKERRQ(ierr); 1452 /* create FEM */ 1453 ierr = SetupDS(*dm,dim,ctx);CHKERRQ(ierr); 1454 /* set initial state */ 1455 ierr = DMCreateGlobalVector(*dm,X);CHKERRQ(ierr); 1456 ierr = PetscObjectSetName((PetscObject) *X, "u");CHKERRQ(ierr); 1457 /* initial static refinement, no solve */ 1458 ierr = LandauSetInitialCondition(*dm, *X, ctx);CHKERRQ(ierr); 1459 ierr = VecViewFromOptions(*X, NULL, "-dm_landau_pre_vec_view");CHKERRQ(ierr); 1460 /* forest refinement */ 1461 if (ctx->errorIndicator) { 1462 /* AMR */ 1463 ierr = adapt(dm,ctx,X);CHKERRQ(ierr); 1464 if ((*dm)->prealloc_only != prealloc_only) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"(*dm)->prealloc_only != prealloc_only"); 1465 ierr = DMViewFromOptions(*dm,NULL,"-dm_landau_amr_dm_view");CHKERRQ(ierr); 1466 ierr = VecViewFromOptions(*X, NULL, "-dm_landau_amr_vec_view");CHKERRQ(ierr); 1467 } 1468 ierr = DMSetApplicationContext(*dm, ctx);CHKERRQ(ierr); 1469 ctx->dmv = *dm; 1470 if (ctx->dmv->prealloc_only != prealloc_only) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"ctx->dmv->prealloc_only != prealloc_only"); 1471 ierr = DMCreateMatrix(ctx->dmv, &ctx->J);CHKERRQ(ierr); 1472 ierr = MatSetOption(ctx->J, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);CHKERRQ(ierr); 1473 ierr = MatSetOption(ctx->J, MAT_STRUCTURALLY_SYMMETRIC, PETSC_TRUE);CHKERRQ(ierr); 1474 if (J) *J = ctx->J; 1475 /* check for types that we need */ 1476 #if defined(PETSC_HAVE_KOKKOS) 1477 if (ctx->deviceType == LANDAU_CPU) { 1478 ierr = PetscObjectTypeCompareAny((PetscObject)ctx->J,&flg,MATSEQAIJKOKKOS,MATMPIAIJKOKKOS,MATAIJKOKKOS,"");CHKERRQ(ierr); 1479 } 1480 #elif defined(PETSC_HAVE_CUDA) 1481 if (ctx->deviceType == LANDAU_CPU) { 1482 ierr = PetscObjectTypeCompareAny((PetscObject)ctx->J,&flg,MATSEQAIJCUSPARSE,MATMPIAIJCUSPARSE,MATAIJCUSPARSE,"");CHKERRQ(ierr); 1483 } 1484 #endif 1485 if (ctx->gpu_assembly) { /* we need GPU object with GPU assembly */ 1486 if (ctx->deviceType == LANDAU_CUDA) { 1487 ierr = PetscObjectTypeCompareAny((PetscObject)ctx->J,&flg,MATSEQAIJCUSPARSE,MATMPIAIJCUSPARSE,MATAIJCUSPARSE,"");CHKERRQ(ierr); 1488 if (!flg) SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"must use '-dm_mat_type aijcusparse -dm_vec_type cuda' for GPU assembly and Cuda"); 1489 } else if (ctx->deviceType == LANDAU_KOKKOS) { 1490 ierr = PetscObjectTypeCompareAny((PetscObject)ctx->J,&flg,MATSEQAIJKOKKOS,MATMPIAIJKOKKOS,MATAIJKOKKOS,"");CHKERRQ(ierr); 1491 #if defined(PETSC_HAVE_KOKKOS_KERNELS) 1492 if (!flg) SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"must use '-dm_mat_type aijkokkos -dm_vec_type kokkos' for GPU assembly and Kokkos"); 1493 #else 1494 if (!flg) SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"must configure with '--download-kokkos-kernels=1' for GPU assembly and Kokkos"); 1495 #endif 1496 } 1497 } 1498 PetscFunctionReturn(0); 1499 } 1500 1501 /*@ 1502 LandauDestroyVelocitySpace - Destroy a DMPlex velocity space mesh 1503 1504 Collective on dm 1505 1506 Input/Output Parameters: 1507 . dm - the dm to destroy 1508 1509 Level: beginner 1510 1511 .keywords: mesh 1512 .seealso: LandauCreateVelocitySpace() 1513 @*/ 1514 PetscErrorCode LandauDestroyVelocitySpace(DM *dm) 1515 { 1516 PetscErrorCode ierr,ii; 1517 LandauCtx *ctx; 1518 PetscContainer container = NULL; 1519 PetscFunctionBegin; 1520 ierr = DMGetApplicationContext(*dm, &ctx);CHKERRQ(ierr); 1521 ierr = PetscObjectQuery((PetscObject)ctx->J,"coloring", (PetscObject*)&container);CHKERRQ(ierr); 1522 if (container) { 1523 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 1524 } 1525 ierr = MatDestroy(&ctx->M);CHKERRQ(ierr); 1526 ierr = MatDestroy(&ctx->J);CHKERRQ(ierr); 1527 for (ii=0;ii<ctx->num_species;ii++) { 1528 ierr = PetscFEDestroy(&ctx->fe[ii]);CHKERRQ(ierr); 1529 } 1530 if (ctx->deviceType == LANDAU_CUDA) { 1531 #if defined(PETSC_HAVE_CUDA) 1532 ierr = LandauCUDAStaticDataClear(ctx->SData_d);CHKERRQ(ierr); 1533 #else 1534 SETERRQ1(ctx->comm,PETSC_ERR_ARG_WRONG,"-landau_device_type %s not built","cuda"); 1535 #endif 1536 } else if (ctx->deviceType == LANDAU_KOKKOS) { 1537 #if defined(PETSC_HAVE_KOKKOS) 1538 ierr = LandauKokkosStaticDataClear(ctx->SData_d);CHKERRQ(ierr); 1539 #else 1540 SETERRQ1(ctx->comm,PETSC_ERR_ARG_WRONG,"-landau_device_type %s not built","kokkos"); 1541 #endif 1542 } else { 1543 if (ctx->SData_d) { /* in a CPU run */ 1544 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, *mass_w = (PetscReal*)ctx->SData_d->mass_w; 1545 ierr = PetscFree5(mass_w,ww,xx,yy,invJ);CHKERRQ(ierr); 1546 if (zz) { 1547 ierr = PetscFree(zz);CHKERRQ(ierr); 1548 } 1549 } 1550 } 1551 ierr = PetscFree(ctx->SData_d);CHKERRQ(ierr); 1552 if (ctx->times[0] > 0) { 1553 ierr = PetscPrintf(ctx->comm, "Landau Operator %d 1.0 %10.3e ....\n",10000,ctx->times[0]);CHKERRQ(ierr); 1554 } 1555 if (ctx->plex != NULL) { 1556 ierr = DMDestroy(&ctx->plex);CHKERRQ(ierr); 1557 } 1558 PetscFree(ctx); 1559 ierr = DMDestroy(dm);CHKERRQ(ierr); 1560 PetscFunctionReturn(0); 1561 } 1562 1563 /* < v, ru > */ 1564 static void f0_s_den(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1565 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1566 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1567 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 1568 { 1569 PetscInt ii = (PetscInt)PetscRealPart(constants[0]); 1570 f0[0] = u[ii]; 1571 } 1572 1573 /* < v, ru > */ 1574 static void f0_s_mom(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1575 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1576 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1577 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 1578 { 1579 PetscInt ii = (PetscInt)PetscRealPart(constants[0]), jj = (PetscInt)PetscRealPart(constants[1]); 1580 f0[0] = x[jj]*u[ii]; /* x momentum */ 1581 } 1582 1583 static void f0_s_v2(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1584 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1585 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1586 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 1587 { 1588 PetscInt i, ii = (PetscInt)PetscRealPart(constants[0]); 1589 double tmp1 = 0.; 1590 for (i = 0; i < dim; ++i) tmp1 += x[i]*x[i]; 1591 f0[0] = tmp1*u[ii]; 1592 } 1593 1594 /* < v, ru > */ 1595 static void f0_s_rden(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1596 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1597 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1598 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 1599 { 1600 PetscInt ii = (PetscInt)PetscRealPart(constants[0]); 1601 f0[0] = 2.*PETSC_PI*x[0]*u[ii]; 1602 } 1603 1604 /* < v, ru > */ 1605 static void f0_s_rmom(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1606 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1607 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1608 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 1609 { 1610 PetscInt ii = (PetscInt)PetscRealPart(constants[0]); 1611 f0[0] = 2.*PETSC_PI*x[0]*x[1]*u[ii]; 1612 } 1613 1614 static void f0_s_rv2(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1615 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1616 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1617 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 1618 { 1619 PetscInt ii = (PetscInt)PetscRealPart(constants[0]); 1620 f0[0] = 2.*PETSC_PI*x[0]*(x[0]*x[0] + x[1]*x[1])*u[ii]; 1621 } 1622 1623 /*@ 1624 LandauPrintNorms - collects moments and prints them 1625 1626 Collective on dm 1627 1628 Input Parameters: 1629 + X - the state 1630 - stepi - current step to print 1631 1632 Level: beginner 1633 1634 .keywords: mesh 1635 .seealso: LandauCreateVelocitySpace() 1636 @*/ 1637 PetscErrorCode LandauPrintNorms(Vec X, PetscInt stepi) 1638 { 1639 PetscErrorCode ierr; 1640 LandauCtx *ctx; 1641 PetscDS prob; 1642 DM dm; 1643 PetscInt cStart, cEnd, dim, ii; 1644 PetscScalar xmomentumtot=0, ymomentumtot=0, zmomentumtot=0, energytot=0, densitytot=0, tt[LANDAU_MAX_SPECIES]; 1645 PetscScalar xmomentum[LANDAU_MAX_SPECIES], ymomentum[LANDAU_MAX_SPECIES], zmomentum[LANDAU_MAX_SPECIES], energy[LANDAU_MAX_SPECIES], density[LANDAU_MAX_SPECIES]; 1646 1647 PetscFunctionBegin; 1648 ierr = VecGetDM(X, &dm);CHKERRQ(ierr); 1649 if (!dm) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no DM"); 1650 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1651 ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr); 1652 if (!ctx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context"); 1653 if (!ctx->plex) { 1654 ierr = DMConvert(ctx->dmv, DMPLEX, &ctx->plex);CHKERRQ(ierr); 1655 } 1656 ierr = DMCreateDS(ctx->plex);CHKERRQ(ierr); 1657 ierr = DMGetDS(ctx->plex, &prob);CHKERRQ(ierr); 1658 /* print momentum and energy */ 1659 for (ii=0;ii<ctx->num_species;ii++) { 1660 PetscScalar user[2] = { (PetscScalar)ii, (PetscScalar)ctx->charges[ii]}; 1661 ierr = PetscDSSetConstants(prob, 2, user);CHKERRQ(ierr); 1662 if (dim==2) { /* 2/3X + 3V (cylindrical coordinates) */ 1663 ierr = PetscDSSetObjective(prob, 0, &f0_s_rden);CHKERRQ(ierr); 1664 ierr = DMPlexComputeIntegralFEM(ctx->plex,X,tt,ctx);CHKERRQ(ierr); 1665 density[ii] = tt[0]*ctx->n_0*ctx->charges[ii]; 1666 ierr = PetscDSSetObjective(prob, 0, &f0_s_rmom);CHKERRQ(ierr); 1667 ierr = DMPlexComputeIntegralFEM(ctx->plex,X,tt,ctx);CHKERRQ(ierr); 1668 zmomentum[ii] = tt[0]*ctx->n_0*ctx->v_0*ctx->masses[ii]; 1669 ierr = PetscDSSetObjective(prob, 0, &f0_s_rv2);CHKERRQ(ierr); 1670 ierr = DMPlexComputeIntegralFEM(ctx->plex,X,tt,ctx);CHKERRQ(ierr); 1671 energy[ii] = tt[0]*0.5*ctx->n_0*ctx->v_0*ctx->v_0*ctx->masses[ii]; 1672 zmomentumtot += zmomentum[ii]; 1673 energytot += energy[ii]; 1674 densitytot += density[ii]; 1675 ierr = PetscPrintf(ctx->comm, "%3D) species-%D: charge density= %20.13e z-momentum= %20.13e energy= %20.13e",stepi,ii,PetscRealPart(density[ii]),PetscRealPart(zmomentum[ii]),PetscRealPart(energy[ii]));CHKERRQ(ierr); 1676 } else { /* 2/3X + 3V */ 1677 ierr = PetscDSSetObjective(prob, 0, &f0_s_den);CHKERRQ(ierr); 1678 ierr = DMPlexComputeIntegralFEM(ctx->plex,X,tt,ctx);CHKERRQ(ierr); 1679 density[ii] = tt[0]*ctx->n_0*ctx->charges[ii]; 1680 ierr = PetscDSSetObjective(prob, 0, &f0_s_mom);CHKERRQ(ierr); 1681 user[1] = 0; 1682 ierr = DMPlexComputeIntegralFEM(ctx->plex,X,tt,ctx);CHKERRQ(ierr); 1683 xmomentum[ii] = tt[0]*ctx->n_0*ctx->v_0*ctx->masses[ii]; 1684 user[1] = 1; 1685 ierr = DMPlexComputeIntegralFEM(ctx->plex,X,tt,ctx);CHKERRQ(ierr); 1686 ymomentum[ii] = tt[0]*ctx->n_0*ctx->v_0*ctx->masses[ii]; 1687 user[1] = 2; 1688 ierr = DMPlexComputeIntegralFEM(ctx->plex,X,tt,ctx);CHKERRQ(ierr); 1689 zmomentum[ii] = tt[0]*ctx->n_0*ctx->v_0*ctx->masses[ii]; 1690 ierr = PetscDSSetObjective(prob, 0, &f0_s_v2);CHKERRQ(ierr); 1691 ierr = DMPlexComputeIntegralFEM(ctx->plex,X,tt,ctx);CHKERRQ(ierr); 1692 energy[ii] = 0.5*tt[0]*ctx->n_0*ctx->v_0*ctx->v_0*ctx->masses[ii]; 1693 ierr = PetscPrintf(ctx->comm, "%3D) species %D: density=%20.13e, x-momentum=%20.13e, y-momentum=%20.13e, z-momentum=%20.13e, energy=%21.13e", 1694 stepi,ii,PetscRealPart(density[ii]),PetscRealPart(xmomentum[ii]),PetscRealPart(ymomentum[ii]),PetscRealPart(zmomentum[ii]),PetscRealPart(energy[ii]));CHKERRQ(ierr); 1695 xmomentumtot += xmomentum[ii]; 1696 ymomentumtot += ymomentum[ii]; 1697 zmomentumtot += zmomentum[ii]; 1698 energytot += energy[ii]; 1699 densitytot += density[ii]; 1700 } 1701 if (ctx->num_species>1) PetscPrintf(ctx->comm, "\n"); 1702 } 1703 /* totals */ 1704 ierr = DMPlexGetHeightStratum(ctx->plex,0,&cStart,&cEnd);CHKERRQ(ierr); 1705 if (ctx->num_species>1) { 1706 if (dim==2) { 1707 ierr = PetscPrintf(ctx->comm, "\t%3D) Total: charge density=%21.13e, momentum=%21.13e, energy=%21.13e (m_i[0]/m_e = %g, %D cells)", 1708 stepi,(double)PetscRealPart(densitytot),(double)PetscRealPart(zmomentumtot),(double)PetscRealPart(energytot),(double)(ctx->masses[1]/ctx->masses[0]),cEnd-cStart);CHKERRQ(ierr); 1709 } else { 1710 ierr = PetscPrintf(ctx->comm, "\t%3D) 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, %D cells)", 1711 stepi,(double)PetscRealPart(densitytot),(double)PetscRealPart(xmomentumtot),(double)PetscRealPart(ymomentumtot),(double)PetscRealPart(zmomentumtot),(double)PetscRealPart(energytot),(double)(ctx->masses[1]/ctx->masses[0]),cEnd-cStart);CHKERRQ(ierr); 1712 } 1713 } else { 1714 ierr = PetscPrintf(ctx->comm, " -- %D cells",cEnd-cStart);CHKERRQ(ierr); 1715 } 1716 if (ctx->verbose > 1) {ierr = PetscPrintf(ctx->comm,", %D sub (vector) threads\n",ctx->subThreadBlockSize);CHKERRQ(ierr);} 1717 else {ierr = PetscPrintf(ctx->comm,"\n");CHKERRQ(ierr);} 1718 PetscFunctionReturn(0); 1719 } 1720 1721 static PetscErrorCode destroy_coloring (void *is) 1722 { 1723 ISColoring tmp = (ISColoring)is; 1724 return ISColoringDestroy(&tmp); 1725 } 1726 1727 /*@ 1728 LandauCreateColoring - create a coloring and add to matrix (Landau context used just for 'print' flag, should be in DMPlex) 1729 1730 Collective on JacP 1731 1732 Input Parameters: 1733 + JacP - matrix to add coloring to 1734 - plex - The DM 1735 1736 Output Parameter: 1737 . container - Container with coloring 1738 1739 Level: beginner 1740 1741 .keywords: mesh 1742 .seealso: LandauCreateVelocitySpace() 1743 @*/ 1744 PetscErrorCode LandauCreateColoring(Mat JacP, DM plex, PetscContainer *container) 1745 { 1746 PetscErrorCode ierr; 1747 PetscInt dim,cell,i,ej,nc,Nv,totDim,numGCells,cStart,cEnd; 1748 ISColoring iscoloring = NULL; 1749 Mat G,Q; 1750 PetscScalar ones[128]; 1751 MatColoring mc; 1752 IS *is; 1753 PetscInt csize,colour,j,k; 1754 const PetscInt *indices; 1755 PetscInt numComp[1]; 1756 PetscInt numDof[4]; 1757 PetscFE fe; 1758 DM colordm; 1759 PetscSection csection, section, globalSection; 1760 PetscDS prob; 1761 LandauCtx *ctx; 1762 1763 PetscFunctionBegin; 1764 ierr = DMGetApplicationContext(plex, &ctx);CHKERRQ(ierr); 1765 ierr = DMGetLocalSection(plex, §ion);CHKERRQ(ierr); 1766 ierr = DMGetGlobalSection(plex, &globalSection);CHKERRQ(ierr); 1767 ierr = DMGetDimension(plex, &dim);CHKERRQ(ierr); 1768 ierr = DMGetDS(plex, &prob);CHKERRQ(ierr); 1769 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 1770 ierr = DMPlexGetHeightStratum(plex,0,&cStart,&cEnd);CHKERRQ(ierr); 1771 numGCells = cEnd - cStart; 1772 /* create cell centered DM */ 1773 ierr = DMClone(plex, &colordm);CHKERRQ(ierr); 1774 ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) plex), dim, 1, PETSC_FALSE, "color_", PETSC_DECIDE, &fe);CHKERRQ(ierr); 1775 ierr = PetscObjectSetName((PetscObject) fe, "color");CHKERRQ(ierr); 1776 ierr = DMSetField(colordm, 0, NULL, (PetscObject)fe);CHKERRQ(ierr); 1777 ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 1778 for (i = 0; i < (dim+1); ++i) numDof[i] = 0; 1779 numDof[dim] = 1; 1780 numComp[0] = 1; 1781 ierr = DMPlexCreateSection(colordm, NULL, numComp, numDof, 0, NULL, NULL, NULL, NULL, &csection);CHKERRQ(ierr); 1782 ierr = PetscSectionSetFieldName(csection, 0, "color");CHKERRQ(ierr); 1783 ierr = DMSetLocalSection(colordm, csection);CHKERRQ(ierr); 1784 ierr = DMViewFromOptions(colordm,NULL,"-color_dm_view");CHKERRQ(ierr); 1785 /* get vertex to element map Q and colroing graph G */ 1786 ierr = MatGetSize(JacP,NULL,&Nv);CHKERRQ(ierr); 1787 ierr = MatCreateAIJ(PETSC_COMM_SELF,PETSC_DECIDE,PETSC_DECIDE,numGCells,Nv,totDim,NULL,0,NULL,&Q);CHKERRQ(ierr); 1788 for (i=0;i<128;i++) ones[i] = 1.0; 1789 for (cell = cStart, ej = 0 ; cell < cEnd; ++cell, ++ej) { 1790 PetscInt numindices,*indices; 1791 ierr = DMPlexGetClosureIndices(plex, section, globalSection, cell, PETSC_TRUE, &numindices, &indices, NULL, NULL);CHKERRQ(ierr); 1792 if (numindices>128) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "too many indices. %D > %D",numindices,128); 1793 ierr = MatSetValues(Q,1,&ej,numindices,indices,ones,ADD_VALUES);CHKERRQ(ierr); 1794 ierr = DMPlexRestoreClosureIndices(plex, section, globalSection, cell, PETSC_TRUE, &numindices, &indices, NULL, NULL);CHKERRQ(ierr); 1795 } 1796 ierr = MatAssemblyBegin(Q, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1797 ierr = MatAssemblyEnd(Q, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1798 ierr = MatMatTransposeMult(Q,Q,MAT_INITIAL_MATRIX,4.0,&G);CHKERRQ(ierr); 1799 ierr = PetscObjectSetName((PetscObject) Q, "Q");CHKERRQ(ierr); 1800 ierr = PetscObjectSetName((PetscObject) G, "coloring graph");CHKERRQ(ierr); 1801 ierr = MatViewFromOptions(G,NULL,"-coloring_mat_view");CHKERRQ(ierr); 1802 ierr = MatViewFromOptions(Q,NULL,"-coloring_mat_view");CHKERRQ(ierr); 1803 ierr = MatDestroy(&Q);CHKERRQ(ierr); 1804 /* coloring */ 1805 ierr = MatColoringCreate(G,&mc);CHKERRQ(ierr); 1806 ierr = MatColoringSetDistance(mc,1);CHKERRQ(ierr); 1807 ierr = MatColoringSetType(mc,MATCOLORINGJP);CHKERRQ(ierr); 1808 ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr); 1809 ierr = MatColoringApply(mc,&iscoloring);CHKERRQ(ierr); 1810 ierr = MatColoringDestroy(&mc);CHKERRQ(ierr); 1811 /* view */ 1812 ierr = ISColoringViewFromOptions(iscoloring,NULL,"-coloring_is_view");CHKERRQ(ierr); 1813 ierr = ISColoringGetIS(iscoloring,PETSC_USE_POINTER,&nc,&is);CHKERRQ(ierr); 1814 if (ctx && ctx->verbose > 2) { 1815 PetscViewer viewer; 1816 Vec color_vec, eidx_vec; 1817 ierr = DMGetGlobalVector(colordm, &color_vec);CHKERRQ(ierr); 1818 ierr = DMGetGlobalVector(colordm, &eidx_vec);CHKERRQ(ierr); 1819 for (colour=0; colour<nc; colour++) { 1820 ierr = ISGetLocalSize(is[colour],&csize);CHKERRQ(ierr); 1821 ierr = ISGetIndices(is[colour],&indices);CHKERRQ(ierr); 1822 for (j=0; j<csize; j++) { 1823 PetscScalar v = (PetscScalar)colour; 1824 k = indices[j]; 1825 ierr = VecSetValues(color_vec,1,&k,&v,INSERT_VALUES); 1826 v = (PetscScalar)k; 1827 ierr = VecSetValues(eidx_vec,1,&k,&v,INSERT_VALUES); 1828 } 1829 ierr = ISRestoreIndices(is[colour],&indices);CHKERRQ(ierr); 1830 } 1831 /* view */ 1832 ierr = PetscViewerVTKOpen(ctx->comm, "color.vtu", FILE_MODE_WRITE, &viewer);CHKERRQ(ierr); 1833 ierr = PetscObjectSetName((PetscObject) color_vec, "color");CHKERRQ(ierr); 1834 ierr = VecView(color_vec, viewer);CHKERRQ(ierr); 1835 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1836 ierr = PetscViewerVTKOpen(ctx->comm, "eidx.vtu", FILE_MODE_WRITE, &viewer);CHKERRQ(ierr); 1837 ierr = PetscObjectSetName((PetscObject) eidx_vec, "element-idx");CHKERRQ(ierr); 1838 ierr = VecView(eidx_vec, viewer);CHKERRQ(ierr); 1839 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1840 ierr = DMRestoreGlobalVector(colordm, &color_vec);CHKERRQ(ierr); 1841 ierr = DMRestoreGlobalVector(colordm, &eidx_vec);CHKERRQ(ierr); 1842 } 1843 ierr = PetscSectionDestroy(&csection);CHKERRQ(ierr); 1844 ierr = DMDestroy(&colordm);CHKERRQ(ierr); 1845 ierr = ISColoringRestoreIS(iscoloring,PETSC_USE_POINTER,&is);CHKERRQ(ierr); 1846 ierr = MatDestroy(&G);CHKERRQ(ierr); 1847 /* stash coloring */ 1848 ierr = PetscContainerCreate(PETSC_COMM_SELF, container);CHKERRQ(ierr); 1849 ierr = PetscContainerSetPointer(*container,(void*)iscoloring);CHKERRQ(ierr); 1850 ierr = PetscContainerSetUserDestroy(*container, destroy_coloring);CHKERRQ(ierr); 1851 ierr = PetscObjectCompose((PetscObject)JacP,"coloring",(PetscObject)*container);CHKERRQ(ierr); 1852 if (ctx && ctx->verbose > 0) { 1853 ierr = PetscPrintf(ctx->comm, "Made coloring with %D colors\n", nc);CHKERRQ(ierr); 1854 } 1855 PetscFunctionReturn(0); 1856 } 1857 1858 PetscErrorCode LandauAssembleOpenMP(PetscInt cStart, PetscInt cEnd, PetscInt totDim, DM plex, PetscSection section, PetscSection globalSection, Mat JacP, PetscScalar elemMats[], PetscContainer container) 1859 { 1860 PetscErrorCode ierr; 1861 IS *is; 1862 PetscInt nc,colour,j; 1863 const PetscInt *clr_idxs; 1864 ISColoring iscoloring; 1865 PetscFunctionBegin; 1866 ierr = PetscContainerGetPointer(container,(void**)&iscoloring);CHKERRQ(ierr); 1867 ierr = ISColoringGetIS(iscoloring,PETSC_USE_POINTER,&nc,&is);CHKERRQ(ierr); 1868 for (colour=0; colour<nc; colour++) { 1869 PetscInt *idx_arr[1024]; /* need to make dynamic for general use */ 1870 PetscScalar *new_el_mats[1024]; 1871 PetscInt idx_size[1024],csize; 1872 ierr = ISGetLocalSize(is[colour],&csize);CHKERRQ(ierr); 1873 if (csize>1024) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "too many elements in color. %D > %D",csize,1024); 1874 ierr = ISGetIndices(is[colour],&clr_idxs);CHKERRQ(ierr); 1875 /* get indices and mats */ 1876 for (j=0; j<csize; j++) { 1877 PetscInt cell = cStart + clr_idxs[j]; 1878 PetscInt numindices,*indices; 1879 PetscScalar *elMat = &elemMats[clr_idxs[j]*totDim*totDim]; 1880 PetscScalar *valuesOrig = elMat; 1881 ierr = DMPlexGetClosureIndices(plex, section, globalSection, cell, PETSC_TRUE, &numindices, &indices, NULL, (PetscScalar **) &elMat);CHKERRQ(ierr); 1882 idx_size[j] = numindices; 1883 ierr = PetscMalloc2(numindices,&idx_arr[j],numindices*numindices,&new_el_mats[j]);CHKERRQ(ierr); 1884 ierr = PetscMemcpy(idx_arr[j],indices,numindices*sizeof(PetscInt));CHKERRQ(ierr); 1885 ierr = PetscMemcpy(new_el_mats[j],elMat,numindices*numindices*sizeof(PetscScalar));CHKERRQ(ierr); 1886 ierr = DMPlexRestoreClosureIndices(plex, section, globalSection, cell, PETSC_TRUE, &numindices, &indices, NULL, (PetscScalar **) &elMat);CHKERRQ(ierr); 1887 if (elMat != valuesOrig) {ierr = DMRestoreWorkArray(plex, numindices*numindices, MPIU_SCALAR, &elMat);} 1888 } 1889 /* assemble matrix */ 1890 for (j=0; j<csize; j++) { 1891 PetscInt numindices = idx_size[j], *indices = idx_arr[j]; 1892 PetscScalar *elMat = new_el_mats[j]; 1893 MatSetValues(JacP,numindices,indices,numindices,indices,elMat,ADD_VALUES); 1894 } 1895 /* free */ 1896 ierr = ISRestoreIndices(is[colour],&clr_idxs);CHKERRQ(ierr); 1897 for (j=0; j<csize; j++) { 1898 ierr = PetscFree2(idx_arr[j],new_el_mats[j]);CHKERRQ(ierr); 1899 } 1900 } 1901 ierr = ISColoringRestoreIS(iscoloring,PETSC_USE_POINTER,&is);CHKERRQ(ierr); 1902 PetscFunctionReturn(0); 1903 } 1904 1905 /* < v, u > */ 1906 static void g0_1(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1907 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1908 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1909 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 1910 { 1911 g0[0] = 1.; 1912 } 1913 1914 /* < v, u > */ 1915 static void g0_r(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1916 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1917 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1918 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 1919 { 1920 g0[0] = 2.*PETSC_PI*x[0]; 1921 } 1922 1923 /*@ 1924 LandauCreateMassMatrix - Create mass matrix for Landau 1925 1926 Collective on dm 1927 1928 Input Parameters: 1929 . dm - the DM object 1930 1931 Output Parameters: 1932 . Amat - The mass matrix (optional), mass matrix is added to the DM context 1933 1934 Level: beginner 1935 1936 .keywords: mesh 1937 .seealso: LandauCreateVelocitySpace() 1938 @*/ 1939 PetscErrorCode LandauCreateMassMatrix(DM dm, Mat *Amat) 1940 { 1941 DM massDM; 1942 PetscDS prob; 1943 PetscInt ii,dim,N1=1,N2; 1944 PetscErrorCode ierr; 1945 LandauCtx *ctx; 1946 Mat M; 1947 1948 PetscFunctionBegin; 1949 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1950 if (Amat) PetscValidPointer(Amat,2); 1951 ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr); 1952 if (!ctx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context"); 1953 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1954 ierr = DMClone(dm, &massDM);CHKERRQ(ierr); 1955 ierr = DMCopyFields(dm, massDM);CHKERRQ(ierr); 1956 ierr = DMCreateDS(massDM);CHKERRQ(ierr); 1957 ierr = DMGetDS(massDM, &prob);CHKERRQ(ierr); 1958 for (ii=0;ii<ctx->num_species;ii++) { 1959 if (dim==3) {ierr = PetscDSSetJacobian(prob, ii, ii, g0_1, NULL, NULL, NULL);CHKERRQ(ierr);} 1960 else {ierr = PetscDSSetJacobian(prob, ii, ii, g0_r, NULL, NULL, NULL);CHKERRQ(ierr);} 1961 } 1962 ierr = DMViewFromOptions(massDM,NULL,"-dm_landau_mass_dm_view");CHKERRQ(ierr); 1963 ierr = DMCreateMatrix(massDM, &M);CHKERRQ(ierr); 1964 ierr = MatSetOption(M, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);CHKERRQ(ierr); 1965 { 1966 Vec locX; 1967 DM plex; 1968 ierr = DMConvert(massDM, DMPLEX, &plex);CHKERRQ(ierr); 1969 ierr = DMGetLocalVector(massDM, &locX);CHKERRQ(ierr); 1970 /* Mass matrix is independent of the input, so no need to fill locX */ 1971 if (plex->prealloc_only != dm->prealloc_only) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "plex->prealloc_only = massDM->prealloc_only %D, =%D",plex->prealloc_only,massDM->prealloc_only); 1972 ierr = DMPlexSNESComputeJacobianFEM(plex, locX, M, M, ctx);CHKERRQ(ierr); 1973 ierr = DMRestoreLocalVector(massDM, &locX);CHKERRQ(ierr); 1974 ierr = DMDestroy(&plex);CHKERRQ(ierr); 1975 } 1976 ierr = DMDestroy(&massDM);CHKERRQ(ierr); 1977 ierr = MatGetSize(ctx->J, &N1, NULL);CHKERRQ(ierr); 1978 ierr = MatGetSize(M, &N2, NULL);CHKERRQ(ierr); 1979 if (N1 != N2) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Incorrect matrix sizes: |Jacobian| = %D, |Mass|=%D",N1,N2); 1980 ierr = PetscObjectSetName((PetscObject)M, "mass");CHKERRQ(ierr); 1981 ierr = MatViewFromOptions(M,NULL,"-dm_landau_mass_mat_view");CHKERRQ(ierr); 1982 ctx->M = M; /* this could be a noop, a = a */ 1983 if (Amat) *Amat = M; 1984 PetscFunctionReturn(0); 1985 } 1986 1987 /*@ 1988 LandauIFunction - TS residual calculation 1989 1990 Collective on ts 1991 1992 Input Parameters: 1993 + TS - The time stepping context 1994 . time_dummy - current time (not used) 1995 - X - Current state 1996 + X_t - Time derivative of current state 1997 . actx - Landau context 1998 1999 Output Parameter: 2000 . F - The residual 2001 2002 Level: beginner 2003 2004 .keywords: mesh 2005 .seealso: LandauCreateVelocitySpace(), LandauIJacobian() 2006 @*/ 2007 PetscErrorCode LandauIFunction(TS ts, PetscReal time_dummy, Vec X, Vec X_t, Vec F, void *actx) 2008 { 2009 PetscErrorCode ierr; 2010 LandauCtx *ctx=(LandauCtx*)actx; 2011 PetscInt dim; 2012 DM dm; 2013 #if defined(PETSC_HAVE_THREADSAFETY) 2014 double starttime, endtime; 2015 #endif 2016 2017 PetscFunctionBegin; 2018 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 2019 ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr); 2020 if (!ctx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context"); 2021 //ierr = MPI_Barrier(MPI_COMM_WORLD);CHKERRMPI(ierr); 2022 ierr = PetscLogEventBegin(ctx->events[11],0,0,0,0);CHKERRQ(ierr); 2023 ierr = PetscLogEventBegin(ctx->events[0],0,0,0,0);CHKERRQ(ierr); 2024 #if defined(PETSC_HAVE_THREADSAFETY) 2025 starttime = MPI_Wtime(); 2026 #endif 2027 ierr = DMGetDimension(ctx->dmv, &dim);CHKERRQ(ierr); 2028 if (!ctx->aux_bool) { 2029 ierr = PetscInfo3(ts, "Create Landau Jacobian t=%g X=%p %s\n",time_dummy,X_t,ctx->aux_bool ? " -- seems to be in line search" : "");CHKERRQ(ierr); 2030 ierr = LandauFormJacobian_Internal(X,ctx->J,dim,0.0,(void*)ctx);CHKERRQ(ierr); 2031 ctx->aux_bool = PETSC_TRUE; 2032 } else { 2033 ierr = PetscInfo(ts, "Skip forming Jacobian, has not changed (should check norm)\n");CHKERRQ(ierr); 2034 } 2035 ierr = MatViewFromOptions(ctx->J,NULL,"-landau_jacobian_mat_view");CHKERRQ(ierr); 2036 /* mat vec for op */ 2037 ierr = MatMult(ctx->J,X,F);CHKERRQ(ierr);CHKERRQ(ierr); /* C*f */ 2038 /* add time term */ 2039 if (X_t) { 2040 ierr = MatMultAdd(ctx->M,X_t,F,F);CHKERRQ(ierr); 2041 } 2042 #if defined(PETSC_HAVE_THREADSAFETY) 2043 endtime = MPI_Wtime(); 2044 ctx->times[0] += (endtime - starttime); 2045 #endif 2046 ierr = PetscLogEventEnd(ctx->events[0],0,0,0,0);CHKERRQ(ierr); 2047 ierr = PetscLogEventEnd(ctx->events[11],0,0,0,0);CHKERRQ(ierr); 2048 PetscFunctionReturn(0); 2049 } 2050 static PetscErrorCode MatrixNfDestroy(void *ptr) 2051 { 2052 PetscInt *nf = (PetscInt *)ptr; 2053 PetscErrorCode ierr; 2054 PetscFunctionBegin; 2055 ierr = PetscFree(nf);CHKERRQ(ierr); 2056 PetscFunctionReturn(0); 2057 } 2058 /*@ 2059 LandauIJacobian - TS Jacobian construction 2060 2061 Collective on ts 2062 2063 Input Parameters: 2064 + TS - The time stepping context 2065 . time_dummy - current time (not used) 2066 - X - Current state 2067 + U_tdummy - Time derivative of current state (not used) 2068 . shift - shift for du/dt term 2069 - actx - Landau context 2070 2071 Output Parameter: 2072 . Amat - Jacobian 2073 + Pmat - same as Amat 2074 2075 Level: beginner 2076 2077 .keywords: mesh 2078 .seealso: LandauCreateVelocitySpace(), LandauIFunction() 2079 @*/ 2080 PetscErrorCode LandauIJacobian(TS ts, PetscReal time_dummy, Vec X, Vec U_tdummy, PetscReal shift, Mat Amat, Mat Pmat, void *actx) 2081 { 2082 PetscErrorCode ierr; 2083 LandauCtx *ctx=(LandauCtx*)actx; 2084 PetscInt dim; 2085 DM dm; 2086 PetscContainer container; 2087 #if defined(PETSC_HAVE_THREADSAFETY) 2088 double starttime, endtime; 2089 #endif 2090 2091 PetscFunctionBegin; 2092 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 2093 ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr); 2094 if (!ctx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context"); 2095 if (Amat!=Pmat || Amat!=ctx->J) SETERRQ(ctx->comm, PETSC_ERR_PLIB, "Amat!=Pmat || Amat!=ctx->J"); 2096 ierr = DMGetDimension(ctx->dmv, &dim);CHKERRQ(ierr); 2097 /* get collision Jacobian into A */ 2098 ierr = PetscLogEventBegin(ctx->events[11],0,0,0,0);CHKERRQ(ierr); 2099 ierr = PetscLogEventBegin(ctx->events[9],0,0,0,0);CHKERRQ(ierr); 2100 #if defined(PETSC_HAVE_THREADSAFETY) 2101 starttime = MPI_Wtime(); 2102 #endif 2103 ierr = PetscInfo2(ts, "Adding just mass to Jacobian t=%g, shift=%g\n",(double)time_dummy,(double)shift);CHKERRQ(ierr); 2104 if (shift==0.0) SETERRQ(ctx->comm, PETSC_ERR_PLIB, "zero shift"); 2105 if (!ctx->aux_bool) SETERRQ(ctx->comm, PETSC_ERR_PLIB, "wrong state"); 2106 if (!ctx->use_matrix_mass) { 2107 ierr = LandauFormJacobian_Internal(X,ctx->J,dim,shift,(void*)ctx);CHKERRQ(ierr); 2108 } else { /* add mass */ 2109 ierr = MatAXPY(Pmat,shift,ctx->M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 2110 } 2111 ctx->aux_bool = PETSC_FALSE; 2112 ierr = MatViewFromOptions(Pmat,NULL,"-landau_mat_view");CHKERRQ(ierr); 2113 /* set number species in Jacobian */ 2114 ierr = PetscObjectQuery((PetscObject) ctx->J, "Nf", (PetscObject *) &container);CHKERRQ(ierr); 2115 if (!container) { 2116 PetscInt *pNf; 2117 ierr = PetscContainerCreate(PETSC_COMM_SELF, &container);CHKERRQ(ierr); 2118 ierr = PetscMalloc(sizeof(PetscInt), &pNf);CHKERRQ(ierr); 2119 *pNf = ctx->num_species + 1000*ctx->numConcurrency; 2120 ierr = PetscContainerSetPointer(container, (void *)pNf);CHKERRQ(ierr); 2121 ierr = PetscContainerSetUserDestroy(container, MatrixNfDestroy);CHKERRQ(ierr); 2122 ierr = PetscObjectCompose((PetscObject)ctx->J, "Nf", (PetscObject) container);CHKERRQ(ierr); 2123 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 2124 } 2125 #if defined(PETSC_HAVE_THREADSAFETY) 2126 endtime = MPI_Wtime(); 2127 ctx->times[0] += (endtime - starttime); 2128 #endif 2129 ierr = PetscLogEventEnd(ctx->events[9],0,0,0,0);CHKERRQ(ierr); 2130 ierr = PetscLogEventEnd(ctx->events[11],0,0,0,0);CHKERRQ(ierr); 2131 PetscFunctionReturn(0); 2132 } 2133