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