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