xref: /petsc/include/petsclandau.h (revision 984ed092ecc7313e63e6aa733432dfb60d929fa3)
1e0eea495SMark #if !defined(PETSCLANDAU_H)
2e0eea495SMark #define PETSCLANDAU_H
3e0eea495SMark 
4e0eea495SMark #include <petscdmplex.h> /*I      "petscdmplex.h"    I*/
5e0eea495SMark #include <petscts.h>
6e0eea495SMark 
78594ddcfSMark Adams PETSC_EXTERN PetscErrorCode DMPlexLandauPrintNorms(Vec, PetscInt);
88594ddcfSMark Adams PETSC_EXTERN PetscErrorCode DMPlexLandauCreateVelocitySpace(MPI_Comm,PetscInt,const char[],Vec*,Mat*,DM*);
98594ddcfSMark Adams PETSC_EXTERN PetscErrorCode DMPlexLandauDestroyVelocitySpace(DM*);
108594ddcfSMark Adams PETSC_EXTERN PetscErrorCode DMPlexLandauAddMaxwellians(DM, Vec, PetscReal, PetscReal[], PetscReal[], PetscInt, PetscInt, void *);
118594ddcfSMark Adams PETSC_EXTERN PetscErrorCode DMPlexLandauCreateMassMatrix(DM dm, Mat *Amat);
128594ddcfSMark Adams PETSC_EXTERN PetscErrorCode DMPlexLandauIFunction(TS, PetscReal,Vec,Vec,Vec,void *);
138594ddcfSMark Adams PETSC_EXTERN PetscErrorCode DMPlexLandauIJacobian(TS, PetscReal,Vec,Vec,PetscReal,Mat,Mat,void *);
14e0eea495SMark 
15bfc784b7SMark Adams typedef int LandauIdx;
16bfc784b7SMark Adams 
17e0eea495SMark /* the Fokker-Planck-Landau context */
18e0eea495SMark #if !defined(LANDAU_MAX_SPECIES)
19*984ed092SMark Adams #if defined(PETSC_USE_DMLANDAU_2D)
20e0eea495SMark #define LANDAU_MAX_SPECIES 10
218a6f2e61SMark Adams #define LANDAU_MAX_GRIDS 3
2252cdd6eaSMark #else
2352cdd6eaSMark #define LANDAU_MAX_SPECIES 3
24*984ed092SMark Adams #define LANDAU_MAX_GRIDS 3
2552cdd6eaSMark #endif
268a6f2e61SMark Adams #else
278a6f2e61SMark Adams #define LANDAU_MAX_GRIDS 3
28e0eea495SMark #endif
29a587d139SMark 
30a587d139SMark #if !defined(LANDAU_MAX_Q)
31a587d139SMark #if defined(LANDAU_MAX_NQ)
32a587d139SMark #error"LANDAU_MAX_NQ but not LANDAU_MAX_Q. Use -DLANDAU_MAX_Q=4 for Q3 elements"
33a587d139SMark #endif
34*984ed092SMark Adams #if defined(PETSC_USE_DMLANDAU_2D)
35a587d139SMark #define LANDAU_MAX_Q 5
36e0eea495SMark #else
37a587d139SMark #define LANDAU_MAX_Q 3
38a587d139SMark #endif
39a587d139SMark #else
40a587d139SMark #undef LANDAU_MAX_NQ
41e0eea495SMark #endif
42a587d139SMark 
43*984ed092SMark Adams #if defined(PETSC_USE_DMLANDAU_2D)
44930e68a5SMark Adams #define LANDAU_MAX_Q_FACE LANDAU_MAX_Q
45a587d139SMark #define LANDAU_MAX_NQ (LANDAU_MAX_Q*LANDAU_MAX_Q)
4628a07af2SMark Adams #define LANDAU_MAX_BATCH_SZ 320
47*984ed092SMark Adams #define LANDAU_DIM 2
48a587d139SMark #else
49930e68a5SMark Adams #define LANDAU_MAX_Q_FACE (LANDAU_MAX_Q*LANDAU_MAX_Q)
50a587d139SMark #define LANDAU_MAX_NQ (LANDAU_MAX_Q*LANDAU_MAX_Q*LANDAU_MAX_Q)
5128a07af2SMark Adams #define LANDAU_MAX_BATCH_SZ 32
52*984ed092SMark Adams #define LANDAU_DIM 3
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
83bfc784b7SMark Adams   // COO
84cada7fc7SMark Adams   void             *coo_elem_offsets;
85cada7fc7SMark Adams   void             *coo_elem_point_offsets;
86cada7fc7SMark Adams   void             *coo_elem_fullNb;
87cada7fc7SMark Adams   LandauIdx        coo_n_cellsTot;
88bfc784b7SMark Adams   LandauIdx        coo_size;
8913241b68SMark Adams   LandauIdx        coo_max_fullnb;
908a6f2e61SMark Adams } LandauStaticData;
91e8d2b73aSMark Adams 
928fdabdddSMark 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;
938fdabdddSMark Adams 
94e0eea495SMark typedef struct {
95e0eea495SMark   PetscBool      interpolate;                  /* Generate intermediate mesh elements */
96a587d139SMark   PetscBool      gpu_assembly;
978a6f2e61SMark Adams   MPI_Comm       comm; /* global communicator to use for errors and diagnostics */
988fdabdddSMark Adams   double         times[LANDAU_NUM_TIMERS];
998a6f2e61SMark Adams   PetscBool      use_matrix_mass;
1008a6f2e61SMark Adams   /* FE */
101e0eea495SMark   PetscFE        fe[LANDAU_MAX_SPECIES];
102e0eea495SMark   /* geometry  */
1038a6f2e61SMark Adams   PetscReal      radius[LANDAU_MAX_GRIDS];
1048a6f2e61SMark Adams   PetscReal      re_radius;           /* RE: radius of refinement along v_perp=0, z>0 */
1058a6f2e61SMark Adams   PetscReal      vperp0_radius1;      /* RE: radius of refinement along v_perp=0 */
1068a6f2e61SMark Adams   PetscReal      vperp0_radius2;      /* RE: radius of refinement along v_perp=0 after origin AMR refinement */
1078a6f2e61SMark Adams   PetscBool      sphere; // not used
1088a6f2e61SMark Adams   PetscBool      inflate; // not used
1098a6f2e61SMark Adams   PetscReal      i_radius[LANDAU_MAX_GRIDS]; // not used
1108a6f2e61SMark Adams   PetscReal      e_radius; // not used
1118a6f2e61SMark Adams   PetscInt       num_sections; // not used
1128a6f2e61SMark Adams   /* AMR */
1138a6f2e61SMark Adams   PetscBool      use_p4est;
1148a6f2e61SMark Adams   PetscInt       numRERefine;        /* RE: refinement along v_perp=0, z > 0 */
1158a6f2e61SMark Adams   PetscInt       nZRefine1;          /* RE: origin refinement after v_perp=0 refinement */
1168a6f2e61SMark Adams   PetscInt       nZRefine2;          /* RE: origin refinement after origin AMR refinement */
1178a6f2e61SMark Adams   PetscInt       numAMRRefine[LANDAU_MAX_GRIDS];  /* normal AMR - refine from origin */
1188a6f2e61SMark Adams   PetscInt       postAMRRefine[LANDAU_MAX_GRIDS]; /* uniform refinement of AMR */
1198a6f2e61SMark Adams   /* relativistic */
120cefb98e8SMark Adams   PetscBool      use_energy_tensor_trick;
121cefb98e8SMark Adams   PetscBool      use_relativistic_corrections;
122e0eea495SMark   /* physics */
123e0eea495SMark   PetscReal      thermal_temps[LANDAU_MAX_SPECIES];
124e0eea495SMark   PetscReal      masses[LANDAU_MAX_SPECIES];  /* mass of each species  */
125e0eea495SMark   PetscReal      charges[LANDAU_MAX_SPECIES]; /* charge of each species  */
126e0eea495SMark   PetscReal      n[LANDAU_MAX_SPECIES];       /* number density of each species  */
127e0eea495SMark   PetscReal      m_0;      /* reference mass */
128e0eea495SMark   PetscReal      v_0;      /* reference velocity */
129e0eea495SMark   PetscReal      n_0;      /* reference number density */
130e0eea495SMark   PetscReal      t_0;      /* reference time */
131e0eea495SMark   PetscReal      Ez;
132e0eea495SMark   PetscReal      epsilon0;
133e0eea495SMark   PetscReal      k;
134e0eea495SMark   PetscReal      lnLam;
1358a6f2e61SMark Adams   PetscReal      electronShift;
136e0eea495SMark   PetscInt       num_species;
1378a6f2e61SMark Adams   PetscInt       num_grids;
1388fdabdddSMark Adams   PetscInt       species_offset[LANDAU_MAX_GRIDS+1]; // for each grid, but same for all batched vertices
1398fdabdddSMark Adams   PetscInt       mat_offset[LANDAU_MAX_GRIDS+1]; // for each grid, but same for all batched vertices
140cb25d741SMark Adams   // batching
141cb25d741SMark Adams   PetscBool      jacobian_field_major_order; // this could be a type but lets not get pedantic
142cb25d741SMark Adams   VecScatter     plex_batch;
143cb25d741SMark Adams   Vec            work_vec;
144cb25d741SMark Adams   IS             batch_is;
145cb25d741SMark Adams   PetscErrorCode (*seqaij_mult)(Mat,Vec,Vec);
146cb25d741SMark Adams   PetscErrorCode (*seqaij_multtranspose)(Mat,Vec,Vec);
147cb25d741SMark Adams   PetscErrorCode (*seqaij_solve)(Mat,Vec,Vec);
148cb25d741SMark Adams   PetscErrorCode (*seqaij_getdiagonal)(Mat,Vec);
149bfc784b7SMark Adams   /* COO */
150aab769cbSMark Adams   PetscBool        coo_assembly;
151e0eea495SMark   /* cache */
152e0eea495SMark   Mat              J;
153e0eea495SMark   Mat              M;
154e0eea495SMark   Vec              X;
155e0eea495SMark   /* derived type */
156e0eea495SMark   void             *data;
157e0eea495SMark   /* computing */
158e0eea495SMark   LandauDeviceType deviceType;
1598a6f2e61SMark Adams   DM               pack;
1608a6f2e61SMark Adams   DM               plex[LANDAU_MAX_GRIDS];
1618fdabdddSMark Adams   LandauStaticData SData_d; /* static geometric data on device */
16254545eeeSMark Adams   /* diagnostics */
16354545eeeSMark Adams   PetscInt         verbose;
16454545eeeSMark Adams   PetscLogEvent    events[20];
165c751c0a2SMark Adams   PetscLogStage    stage;
166*984ed092SMark Adams   PetscObjectState norm_state;
1678fdabdddSMark Adams   PetscInt         batch_sz;
1688fdabdddSMark Adams   PetscInt         batch_view_idx;
169e0eea495SMark } LandauCtx;
170e0eea495SMark 
1718fdabdddSMark Adams #define LANDAU_SPECIES_MAJOR
1728fdabdddSMark Adams #if !defined(LANDAU_SPECIES_MAJOR)
1738fdabdddSMark Adams #define LAND_PACK_IDX(_b,_g) (_b*ctx->num_grids + _g)
1748fdabdddSMark Adams #define LAND_MOFFSET(_b,_g,_nbch,_ngrid,_mat_off) (_b*_mat_off[_ngrid] + _mat_off[_g])
1758fdabdddSMark Adams #else
1768fdabdddSMark Adams #define LAND_PACK_IDX(_b,_g) (_g*ctx->batch_sz + _b)
1778fdabdddSMark Adams #define LAND_MOFFSET(_b,_g,_nbch,_ngrid,_mat_off) (_nbch*_mat_off[_g] + _b*(_mat_off[_g+1] - _mat_off[_g]))
1788fdabdddSMark Adams #endif
1798fdabdddSMark Adams 
180a587d139SMark typedef struct {
181a587d139SMark   PetscReal scale;
182bfc784b7SMark Adams   LandauIdx gid;   // Landau matrix index (<10,000)
183a587d139SMark } pointInterpolationP4est;
184a587d139SMark typedef struct _lP4estVertexMaps {
185a587d139SMark   LandauIdx                (*gIdx)[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ]; // #elems *  LANDAU_MAX_NQ (spoof for max , Nb) on device,
186a587d139SMark   LandauIdx                num_elements;
187a587d139SMark   LandauIdx                num_reduced;
188a587d139SMark   LandauIdx                num_face;  // (Q or Q^2 for 3D)
189a587d139SMark   LandauDeviceType         deviceType;
190a587d139SMark   PetscInt                 Nf;
191a587d139SMark   pointInterpolationP4est (*c_maps)[LANDAU_MAX_Q_FACE];
1928a6f2e61SMark Adams   struct _lP4estVertexMaps*d_self;
193a587d139SMark   void                    *vp1,*vp2,*vp3;
1948a6f2e61SMark Adams   PetscInt                numgrids;
195a587d139SMark } P4estVertexMaps;
196a587d139SMark 
197e0eea495SMark PETSC_EXTERN PetscErrorCode LandauCreateColoring(Mat, DM, PetscContainer *);
198e0eea495SMark #if defined(PETSC_HAVE_CUDA)
199a82411e9SMark Adams PETSC_EXTERN PetscErrorCode LandauCUDAJacobian(DM[], const PetscInt, const PetscInt, const PetscInt, const PetscInt[], PetscReal[], PetscScalar[], const PetscScalar[],
200d2319daeSMark Adams                                                const LandauStaticData *, const PetscReal, const PetscLogEvent[],const PetscInt[], const PetscInt[], Mat[], Mat);
2018a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauCUDACreateMatMaps(P4estVertexMaps *, pointInterpolationP4est (*)[LANDAU_MAX_Q_FACE], PetscInt[], PetscInt, PetscInt);
2028a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauCUDADestroyMatMaps(P4estVertexMaps *, PetscInt);
2038fdabdddSMark Adams PETSC_EXTERN PetscErrorCode LandauCUDAStaticDataSet(DM, const PetscInt, const PetscInt, const PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscReal [], PetscReal [], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[],
2048a6f2e61SMark Adams                                                     PetscReal[], LandauStaticData *);
2058a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauCUDAStaticDataClear(LandauStaticData *);
206e0eea495SMark #endif
207e0eea495SMark #if defined(PETSC_HAVE_KOKKOS)
208a82411e9SMark Adams PETSC_EXTERN PetscErrorCode LandauKokkosJacobian(DM[], const PetscInt, const PetscInt, const PetscInt, const PetscInt[], PetscReal[], PetscScalar[],  const PetscScalar[],
209d2319daeSMark Adams                                                  const LandauStaticData *, const PetscReal, const PetscLogEvent[], const PetscInt[], const PetscInt[], Mat[], Mat);
2108a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauKokkosCreateMatMaps(P4estVertexMaps *, pointInterpolationP4est (*)[LANDAU_MAX_Q_FACE], PetscInt[], PetscInt, PetscInt);
2118a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauKokkosDestroyMatMaps(P4estVertexMaps *, PetscInt);
2128fdabdddSMark Adams PETSC_EXTERN PetscErrorCode LandauKokkosStaticDataSet(DM, const PetscInt, const PetscInt, const PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscReal [], PetscReal [], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[],
2138a6f2e61SMark Adams                                                       PetscReal[], LandauStaticData *);
2148a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauKokkosStaticDataClear(LandauStaticData*);
215e0eea495SMark #endif
216e0eea495SMark 
217e0eea495SMark #endif /* PETSCLANDAU_H */
218