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