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