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