1e0eea495SMark #if !defined(PETSCLANDAU_H) 2e0eea495SMark #define PETSCLANDAU_H 3e0eea495SMark 4e0eea495SMark #include <petscdmplex.h> /*I "petscdmplex.h" I*/ 5e0eea495SMark #include <petscts.h> 6e0eea495SMark 7e0eea495SMark PETSC_EXTERN PetscErrorCode LandauPrintNorms(Vec, PetscInt); 8e0eea495SMark PETSC_EXTERN PetscErrorCode LandauCreateVelocitySpace(MPI_Comm,PetscInt,const char[],Vec*,Mat*,DM*); 9e0eea495SMark PETSC_EXTERN PetscErrorCode LandauDestroyVelocitySpace(DM*); 108fdabdddSMark Adams PETSC_EXTERN PetscErrorCode LandauAddMaxwellians(DM, Vec, PetscReal, PetscReal[], PetscReal[], PetscInt, PetscInt, void *); 11e0eea495SMark PETSC_EXTERN PetscErrorCode LandauCreateMassMatrix(DM dm, Mat *Amat); 12e0eea495SMark PETSC_EXTERN PetscErrorCode LandauIFunction(TS, PetscReal,Vec,Vec,Vec,void *); 13e0eea495SMark PETSC_EXTERN PetscErrorCode LandauIJacobian(TS, PetscReal,Vec,Vec,PetscReal,Mat,Mat,void *); 14e0eea495SMark 15e0eea495SMark /* the Fokker-Planck-Landau context */ 16e0eea495SMark #if !defined(LANDAU_DIM) 17e0eea495SMark #define LANDAU_DIM 2 18e0eea495SMark #endif 19a587d139SMark 20e0eea495SMark #if !defined(LANDAU_MAX_SPECIES) 2152cdd6eaSMark #if LANDAU_DIM==2 22e0eea495SMark #define LANDAU_MAX_SPECIES 10 238a6f2e61SMark Adams #define LANDAU_MAX_GRIDS 3 2452cdd6eaSMark #else 2552cdd6eaSMark #define LANDAU_MAX_SPECIES 3 268a6f2e61SMark Adams #define LANDAU_MAX_GRIDS 2 2752cdd6eaSMark #endif 288a6f2e61SMark Adams #else 298a6f2e61SMark Adams #define LANDAU_MAX_GRIDS 3 30e0eea495SMark #endif 31a587d139SMark 32f37ccb5fSMark Adams #define LANDAU_MAX_BATCH_SZ 320 338fdabdddSMark Adams 34a587d139SMark #if !defined(LANDAU_MAX_Q) 35a587d139SMark #if defined(LANDAU_MAX_NQ) 36a587d139SMark #error"LANDAU_MAX_NQ but not LANDAU_MAX_Q. Use -DLANDAU_MAX_Q=4 for Q3 elements" 37a587d139SMark #endif 3852cdd6eaSMark #if LANDAU_DIM==2 39a587d139SMark #define LANDAU_MAX_Q 5 40e0eea495SMark #else 41a587d139SMark #define LANDAU_MAX_Q 3 42a587d139SMark #endif 43a587d139SMark #else 44a587d139SMark #undef LANDAU_MAX_NQ 45e0eea495SMark #endif 46a587d139SMark 47a587d139SMark #if LANDAU_DIM==2 48930e68a5SMark Adams #define LANDAU_MAX_Q_FACE LANDAU_MAX_Q 49a587d139SMark #define LANDAU_MAX_NQ (LANDAU_MAX_Q*LANDAU_MAX_Q) 50a587d139SMark #else 51930e68a5SMark Adams #define LANDAU_MAX_Q_FACE (LANDAU_MAX_Q*LANDAU_MAX_Q) 52a587d139SMark #define LANDAU_MAX_NQ (LANDAU_MAX_Q*LANDAU_MAX_Q*LANDAU_MAX_Q) 53a587d139SMark #endif 54a587d139SMark 55e0eea495SMark typedef enum {LANDAU_CUDA, LANDAU_KOKKOS, LANDAU_CPU} LandauDeviceType; 56e8d2b73aSMark Adams 578fdabdddSMark Adams // static data - will be "device" data 58e8d2b73aSMark Adams typedef struct { 59e8d2b73aSMark Adams void *invJ; // nip*dim*dim 60e8d2b73aSMark Adams void *D; // nq*nb*dim 61e8d2b73aSMark Adams void *B; // nq*nb 62e8d2b73aSMark Adams void *alpha; // ns 63e8d2b73aSMark Adams void *beta; // ns 64e8d2b73aSMark Adams void *invMass; // ns 65e8d2b73aSMark Adams void *w; // nip 66e8d2b73aSMark Adams void *x; // nip 67e8d2b73aSMark Adams void *y; // nip 68e8d2b73aSMark Adams void *z; // nip 69e8d2b73aSMark Adams void *Eq_m; // ns - dynamic 70e8d2b73aSMark Adams void *f; // nip*Nf - dynamic (IP) 71e8d2b73aSMark Adams void *dfdx; // nip*Nf - dynamic (IP) 72e8d2b73aSMark Adams void *dfdy; // nip*Nf - dynamic (IP) 73e8d2b73aSMark Adams void *dfdz; // nip*Nf - dynamic (IP) 74e8d2b73aSMark Adams int dim_,ns_,nip_,nq_,nb_; 758fdabdddSMark Adams void *NCells; // remove and ise elem_offset - TODO 768fdabdddSMark Adams void *species_offset; // for each grid, but same for all batched vertices 778fdabdddSMark Adams void *mat_offset; // for each grid, but same for all batched vertices 788fdabdddSMark Adams void *elem_offset; // for each grid, but same for all batched vertices 798fdabdddSMark Adams void *ip_offset; // for each grid, but same for all batched vertices 808fdabdddSMark Adams void *ipf_offset; // for each grid, but same for all batched vertices 818fdabdddSMark Adams void *ipfdf_data; // for each grid, but same for all batched vertices 828fdabdddSMark Adams void *maps; // for each grid, but same for all batched vertices 838a6f2e61SMark Adams } LandauStaticData; 84e8d2b73aSMark Adams 858fdabdddSMark Adams typedef enum {LANDAU_EX2_TSSOLVE, LANDAU_MATRIX_TOTAL, LANDAU_OPERATOR, LANDAU_JACOBIAN_COUNT, LANDAU_JACOBIAN, LANDAU_MASS, LANDAU_F_DF, LANDAU_KERNEL, KSP_FACTOR, KSP_SOLVE, LANDAU_NUM_TIMERS} LandauOMPTimers; 868fdabdddSMark Adams 87e0eea495SMark typedef struct { 88e0eea495SMark PetscBool interpolate; /* Generate intermediate mesh elements */ 89a587d139SMark PetscBool gpu_assembly; 908a6f2e61SMark Adams MPI_Comm comm; /* global communicator to use for errors and diagnostics */ 918fdabdddSMark Adams double times[LANDAU_NUM_TIMERS]; 928a6f2e61SMark Adams PetscBool use_matrix_mass; 938a6f2e61SMark Adams /* FE */ 94e0eea495SMark PetscFE fe[LANDAU_MAX_SPECIES]; 95e0eea495SMark /* geometry */ 968a6f2e61SMark Adams PetscReal radius[LANDAU_MAX_GRIDS]; 978a6f2e61SMark Adams PetscReal re_radius; /* RE: radius of refinement along v_perp=0, z>0 */ 988a6f2e61SMark Adams PetscReal vperp0_radius1; /* RE: radius of refinement along v_perp=0 */ 998a6f2e61SMark Adams PetscReal vperp0_radius2; /* RE: radius of refinement along v_perp=0 after origin AMR refinement */ 1008a6f2e61SMark Adams PetscBool sphere; // not used 1018a6f2e61SMark Adams PetscBool inflate; // not used 1028a6f2e61SMark Adams PetscReal i_radius[LANDAU_MAX_GRIDS]; // not used 1038a6f2e61SMark Adams PetscReal e_radius; // not used 1048a6f2e61SMark Adams PetscInt num_sections; // not used 1058a6f2e61SMark Adams /* AMR */ 1068a6f2e61SMark Adams PetscBool use_p4est; 1078a6f2e61SMark Adams PetscInt numRERefine; /* RE: refinement along v_perp=0, z > 0 */ 1088a6f2e61SMark Adams PetscInt nZRefine1; /* RE: origin refinement after v_perp=0 refinement */ 1098a6f2e61SMark Adams PetscInt nZRefine2; /* RE: origin refinement after origin AMR refinement */ 1108a6f2e61SMark Adams PetscInt numAMRRefine[LANDAU_MAX_GRIDS]; /* normal AMR - refine from origin */ 1118a6f2e61SMark Adams PetscInt postAMRRefine[LANDAU_MAX_GRIDS]; /* uniform refinement of AMR */ 1128a6f2e61SMark Adams /* relativistic */ 113cefb98e8SMark Adams PetscBool use_energy_tensor_trick; 114cefb98e8SMark Adams PetscBool use_relativistic_corrections; 115e0eea495SMark /* physics */ 116e0eea495SMark PetscReal thermal_temps[LANDAU_MAX_SPECIES]; 117e0eea495SMark PetscReal masses[LANDAU_MAX_SPECIES]; /* mass of each species */ 118e0eea495SMark PetscReal charges[LANDAU_MAX_SPECIES]; /* charge of each species */ 119e0eea495SMark PetscReal n[LANDAU_MAX_SPECIES]; /* number density of each species */ 120e0eea495SMark PetscReal m_0; /* reference mass */ 121e0eea495SMark PetscReal v_0; /* reference velocity */ 122e0eea495SMark PetscReal n_0; /* reference number density */ 123e0eea495SMark PetscReal t_0; /* reference time */ 124e0eea495SMark PetscReal Ez; 125e0eea495SMark PetscReal epsilon0; 126e0eea495SMark PetscReal k; 127e0eea495SMark PetscReal lnLam; 1288a6f2e61SMark Adams PetscReal electronShift; 129e0eea495SMark PetscInt num_species; 1308a6f2e61SMark Adams PetscInt num_grids; 1318fdabdddSMark Adams PetscInt species_offset[LANDAU_MAX_GRIDS+1]; // for each grid, but same for all batched vertices 1328fdabdddSMark Adams PetscInt mat_offset[LANDAU_MAX_GRIDS+1]; // for each grid, but same for all batched vertices 133*cb25d741SMark Adams // batching 134*cb25d741SMark Adams PetscBool jacobian_field_major_order; // this could be a type but lets not get pedantic 135*cb25d741SMark Adams VecScatter plex_batch; 136*cb25d741SMark Adams Vec work_vec; 137*cb25d741SMark Adams IS batch_is; 138*cb25d741SMark Adams PetscErrorCode (*seqaij_mult)(Mat,Vec,Vec); 139*cb25d741SMark Adams PetscErrorCode (*seqaij_multtranspose)(Mat,Vec,Vec); 140*cb25d741SMark Adams PetscErrorCode (*seqaij_solve)(Mat,Vec,Vec); 141*cb25d741SMark Adams PetscErrorCode (*seqaij_getdiagonal)(Mat,Vec); 142*cb25d741SMark Adams // 143e0eea495SMark /* cache */ 144e0eea495SMark Mat J; 145e0eea495SMark Mat M; 146e0eea495SMark Vec X; 147e0eea495SMark /* derived type */ 148e0eea495SMark void *data; 149e0eea495SMark /* computing */ 150e0eea495SMark LandauDeviceType deviceType; 151bddcd29dSMark Adams PetscInt subThreadBlockSize; 1528a6f2e61SMark Adams DM pack; 1538a6f2e61SMark Adams DM plex[LANDAU_MAX_GRIDS]; 1548fdabdddSMark Adams LandauStaticData SData_d; /* static geometric data on device */ 15554545eeeSMark Adams /* diagnostics */ 15654545eeeSMark Adams PetscInt verbose; 15754545eeeSMark Adams PetscLogEvent events[20]; 158c751c0a2SMark Adams PetscLogStage stage; 1598fdabdddSMark Adams PetscBool aux_bool; /* helper */ 1608fdabdddSMark Adams PetscInt batch_sz; 1618fdabdddSMark Adams PetscInt batch_view_idx; 162e0eea495SMark } LandauCtx; 163e0eea495SMark 1648fdabdddSMark Adams #define LANDAU_SPECIES_MAJOR 1658fdabdddSMark Adams #if !defined(LANDAU_SPECIES_MAJOR) 1668fdabdddSMark Adams #define LAND_PACK_IDX(_b,_g) (_b*ctx->num_grids + _g) 1678fdabdddSMark Adams #define LAND_MOFFSET(_b,_g,_nbch,_ngrid,_mat_off) (_b*_mat_off[_ngrid] + _mat_off[_g]) 1688fdabdddSMark Adams #else 1698fdabdddSMark Adams #define LAND_PACK_IDX(_b,_g) (_g*ctx->batch_sz + _b) 1708fdabdddSMark Adams #define LAND_MOFFSET(_b,_g,_nbch,_ngrid,_mat_off) (_nbch*_mat_off[_g] + _b*(_mat_off[_g+1] - _mat_off[_g])) 1718fdabdddSMark Adams #endif 1728fdabdddSMark Adams 173a587d139SMark typedef int LandauIdx; 174a587d139SMark typedef struct { 175a587d139SMark PetscReal scale; 176a587d139SMark LandauIdx gid; // Lanadu matrix index (<10,000) 177a587d139SMark } pointInterpolationP4est; 178a587d139SMark typedef struct _lP4estVertexMaps { 179a587d139SMark LandauIdx (*gIdx)[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ]; // #elems * LANDAU_MAX_NQ (spoof for max , Nb) on device, 180a587d139SMark LandauIdx num_elements; 181a587d139SMark LandauIdx num_reduced; 182a587d139SMark LandauIdx num_face; // (Q or Q^2 for 3D) 183a587d139SMark LandauDeviceType deviceType; 184a587d139SMark PetscInt Nf; 185a587d139SMark PetscInt Nq; 186a587d139SMark pointInterpolationP4est (*c_maps)[LANDAU_MAX_Q_FACE]; 1878a6f2e61SMark Adams struct _lP4estVertexMaps*d_self; 188a587d139SMark void *vp1,*vp2,*vp3; 1898a6f2e61SMark Adams PetscInt numgrids; 190a587d139SMark } P4estVertexMaps; 191a587d139SMark 192e0eea495SMark PETSC_EXTERN PetscErrorCode LandauCreateColoring(Mat, DM, PetscContainer *); 193e0eea495SMark #if defined(PETSC_HAVE_CUDA) 194a82411e9SMark Adams PETSC_EXTERN PetscErrorCode LandauCUDAJacobian(DM[], const PetscInt, const PetscInt, const PetscInt, const PetscInt[], PetscReal[], PetscScalar[], const PetscScalar[], 1958a6f2e61SMark Adams const LandauStaticData *, const PetscInt, const PetscReal, const PetscLogEvent[],const PetscInt[], const PetscInt[], Mat[], Mat); 1968a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauCUDACreateMatMaps(P4estVertexMaps *, pointInterpolationP4est (*)[LANDAU_MAX_Q_FACE], PetscInt[], PetscInt, PetscInt); 1978a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauCUDADestroyMatMaps(P4estVertexMaps *, PetscInt); 1988fdabdddSMark Adams PETSC_EXTERN PetscErrorCode LandauCUDAStaticDataSet(DM, const PetscInt, const PetscInt, const PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscReal [], PetscReal [], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], 1998a6f2e61SMark Adams PetscReal[], LandauStaticData *); 2008a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauCUDAStaticDataClear(LandauStaticData *); 201e0eea495SMark #endif 202e0eea495SMark #if defined(PETSC_HAVE_KOKKOS) 203a82411e9SMark Adams PETSC_EXTERN PetscErrorCode LandauKokkosJacobian(DM[], const PetscInt, const PetscInt, const PetscInt, const PetscInt[], PetscReal[], PetscScalar[], const PetscScalar[], 2048a6f2e61SMark Adams const LandauStaticData *, const PetscInt, const PetscReal, const PetscLogEvent[], const PetscInt[], const PetscInt[], Mat[], Mat); 2058a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauKokkosCreateMatMaps(P4estVertexMaps *, pointInterpolationP4est (*)[LANDAU_MAX_Q_FACE], PetscInt[], PetscInt, PetscInt); 2068a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauKokkosDestroyMatMaps(P4estVertexMaps *, PetscInt); 2078fdabdddSMark Adams PETSC_EXTERN PetscErrorCode LandauKokkosStaticDataSet(DM, const PetscInt, const PetscInt, const PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscReal [], PetscReal [], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], 2088a6f2e61SMark Adams PetscReal[], LandauStaticData *); 2098a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauKokkosStaticDataClear(LandauStaticData*); 210e0eea495SMark #endif 211e0eea495SMark 212e0eea495SMark #endif /* PETSCLANDAU_H */ 213