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