xref: /petsc/include/petsclandau.h (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 *);
1024ded41bSmarkadams4 PETSC_EXTERN PetscErrorCode DMPlexLandauAddToFunction(DM, Vec, PetscErrorCode (*)(DM, Vec, PetscInt, PetscInt, PetscInt, void *), void *);
11c3e4dd79SMark Adams PETSC_EXTERN PetscErrorCode DMPlexLandauAddMaxwellians(DM, Vec, PetscReal, PetscReal[], PetscReal[], PetscInt, PetscInt, PetscInt, void *);
128594ddcfSMark Adams PETSC_EXTERN PetscErrorCode DMPlexLandauCreateMassMatrix(DM dm, Mat *Amat);
138594ddcfSMark Adams PETSC_EXTERN PetscErrorCode DMPlexLandauIFunction(TS, PetscReal, Vec, Vec, Vec, void *);
148594ddcfSMark Adams PETSC_EXTERN PetscErrorCode DMPlexLandauIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
15e0eea495SMark 
16bfc784b7SMark Adams typedef int LandauIdx;
17bfc784b7SMark Adams 
18e0eea495SMark /* the Fokker-Planck-Landau context */
19e0eea495SMark #if !defined(LANDAU_MAX_SPECIES)
20984ed092SMark Adams #if defined(PETSC_USE_DMLANDAU_2D)
21e0eea495SMark #define LANDAU_MAX_SPECIES 10
228a6f2e61SMark Adams #define LANDAU_MAX_GRIDS   3
2352cdd6eaSMark #else
24763ae2f8SMark Adams #define LANDAU_MAX_SPECIES 10
25c3e4dd79SMark Adams #define LANDAU_MAX_GRIDS   3
2652cdd6eaSMark #endif
278a6f2e61SMark Adams #else
288a6f2e61SMark Adams #define LANDAU_MAX_GRIDS 3
29e0eea495SMark #endif
30a587d139SMark 
31a587d139SMark #if !defined(LANDAU_MAX_Q)
32a587d139SMark #if defined(LANDAU_MAX_NQ)
33a587d139SMark #error "LANDAU_MAX_NQ but not LANDAU_MAX_Q. Use -DLANDAU_MAX_Q=4 for Q3 elements"
34a587d139SMark #endif
35984ed092SMark Adams #if defined(PETSC_USE_DMLANDAU_2D)
363b65d667Smarkadams4 #define LANDAU_MAX_Q 5
37c3e4dd79SMark Adams #else
383b65d667Smarkadams4 // 3D CUDA fails with > 3 (KK-CUDA is OK)
39c3e4dd79SMark Adams #define LANDAU_MAX_Q 3
40a587d139SMark #endif
41a587d139SMark #else
42a587d139SMark #undef LANDAU_MAX_NQ
43e0eea495SMark #endif
44a587d139SMark 
45984ed092SMark Adams #if defined(PETSC_USE_DMLANDAU_2D)
46930e68a5SMark Adams #define LANDAU_MAX_Q_FACE   LANDAU_MAX_Q
47a587d139SMark #define LANDAU_MAX_NQ       (LANDAU_MAX_Q * LANDAU_MAX_Q)
483b65d667Smarkadams4 #define LANDAU_MAX_BATCH_SZ 1024
49984ed092SMark Adams #define LANDAU_DIM          2
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)
533b65d667Smarkadams4 #define LANDAU_MAX_BATCH_SZ 64
54984ed092SMark Adams #define LANDAU_DIM          3
55a587d139SMark #endif
56a587d139SMark 
57*9371c9d4SSatish Balay typedef enum {
58*9371c9d4SSatish Balay   LANDAU_CUDA,
59*9371c9d4SSatish Balay   LANDAU_KOKKOS,
60*9371c9d4SSatish Balay   LANDAU_CPU
61*9371c9d4SSatish Balay } LandauDeviceType;
62e8d2b73aSMark Adams 
638fdabdddSMark Adams // static data - will be "device" data
64e8d2b73aSMark Adams typedef struct {
65e8d2b73aSMark Adams   void     *invJ;    // nip*dim*dim
66e8d2b73aSMark Adams   void     *D;       // nq*nb*dim
67e8d2b73aSMark Adams   void     *B;       // nq*nb
68e8d2b73aSMark Adams   void     *alpha;   // ns
69e8d2b73aSMark Adams   void     *beta;    // ns
70e8d2b73aSMark Adams   void     *invMass; // ns
71e8d2b73aSMark Adams   void     *w;       // nip
72e8d2b73aSMark Adams   void     *x;       // nip
73e8d2b73aSMark Adams   void     *y;       // nip
74e8d2b73aSMark Adams   void     *z;       // nip
75e8d2b73aSMark Adams   void     *Eq_m;    // ns - dynamic
76e8d2b73aSMark Adams   void     *f;       //  nip*Nf - dynamic (IP)
77e8d2b73aSMark Adams   void     *dfdx;    // nip*Nf - dynamic (IP)
78e8d2b73aSMark Adams   void     *dfdy;    // nip*Nf - dynamic (IP)
79e8d2b73aSMark Adams   void     *dfdz;    // nip*Nf - dynamic (IP)
80e8d2b73aSMark Adams   int       dim_, ns_, nip_, nq_, nb_;
818fdabdddSMark Adams   void     *NCells;         // remove and ise elem_offset - TODO
828fdabdddSMark Adams   void     *species_offset; // for each grid, but same for all batched vertices
838fdabdddSMark Adams   void     *mat_offset;     // for each grid, but same for all batched vertices
848fdabdddSMark Adams   void     *elem_offset;    // for each grid, but same for all batched vertices
858fdabdddSMark Adams   void     *ip_offset;      // for each grid, but same for all batched vertices
868fdabdddSMark Adams   void     *ipf_offset;     // for each grid, but same for all batched vertices
878fdabdddSMark Adams   void     *ipfdf_data;     // for each grid, but same for all batched vertices
888fdabdddSMark Adams   void     *maps;           // for each grid, but same for all batched vertices
89bfc784b7SMark Adams   // COO
90cada7fc7SMark Adams   void     *coo_elem_offsets;
91cada7fc7SMark Adams   void     *coo_elem_point_offsets;
92cada7fc7SMark Adams   void     *coo_elem_fullNb;
93a31f6053SMark Adams   void     *coo_vals;
94cada7fc7SMark Adams   LandauIdx coo_n_cellsTot;
95bfc784b7SMark Adams   LandauIdx coo_size;
9613241b68SMark Adams   LandauIdx coo_max_fullnb;
978a6f2e61SMark Adams } LandauStaticData;
98e8d2b73aSMark Adams 
99*9371c9d4SSatish Balay typedef enum {
100*9371c9d4SSatish Balay   LANDAU_EX2_TSSOLVE,
101*9371c9d4SSatish Balay   LANDAU_MATRIX_TOTAL,
102*9371c9d4SSatish Balay   LANDAU_OPERATOR,
103*9371c9d4SSatish Balay   LANDAU_JACOBIAN_COUNT,
104*9371c9d4SSatish Balay   LANDAU_JACOBIAN,
105*9371c9d4SSatish Balay   LANDAU_MASS,
106*9371c9d4SSatish Balay   LANDAU_F_DF,
107*9371c9d4SSatish Balay   LANDAU_KERNEL,
108*9371c9d4SSatish Balay   KSP_FACTOR,
109*9371c9d4SSatish Balay   KSP_SOLVE,
110*9371c9d4SSatish Balay   LANDAU_NUM_TIMERS
111*9371c9d4SSatish Balay } LandauOMPTimers;
1128fdabdddSMark Adams 
113e0eea495SMark typedef struct {
114e0eea495SMark   PetscBool  interpolate; /* Generate intermediate mesh elements */
115a587d139SMark   PetscBool  gpu_assembly;
1168a6f2e61SMark Adams   MPI_Comm   comm; /* global communicator to use for errors and diagnostics */
1178fdabdddSMark Adams   double     times[LANDAU_NUM_TIMERS];
1188a6f2e61SMark Adams   PetscBool  use_matrix_mass;
1198a6f2e61SMark Adams   /* FE */
120e0eea495SMark   PetscFE    fe[LANDAU_MAX_SPECIES];
121e0eea495SMark   /* geometry  */
1228a6f2e61SMark Adams   PetscReal  radius[LANDAU_MAX_GRIDS];
1238a6f2e61SMark Adams   PetscReal  re_radius;                  /* RE: radius of refinement along v_perp=0, z>0 */
1248a6f2e61SMark Adams   PetscReal  vperp0_radius1;             /* RE: radius of refinement along v_perp=0 */
1258a6f2e61SMark Adams   PetscReal  vperp0_radius2;             /* RE: radius of refinement along v_perp=0 after origin AMR refinement */
1268a6f2e61SMark Adams   PetscBool  sphere;                     // not used
1278a6f2e61SMark Adams   PetscBool  inflate;                    // not used
1288a6f2e61SMark Adams   PetscReal  i_radius[LANDAU_MAX_GRIDS]; // not used
1298a6f2e61SMark Adams   PetscReal  e_radius;                   // not used
1308a6f2e61SMark Adams   PetscInt   num_sections;               // not used
1318a6f2e61SMark Adams   /* AMR */
1328a6f2e61SMark Adams   PetscBool  use_p4est;
1338a6f2e61SMark Adams   PetscInt   numRERefine;                     /* RE: refinement along v_perp=0, z > 0 */
1348a6f2e61SMark Adams   PetscInt   nZRefine1;                       /* RE: origin refinement after v_perp=0 refinement */
1358a6f2e61SMark Adams   PetscInt   nZRefine2;                       /* RE: origin refinement after origin AMR refinement */
1368a6f2e61SMark Adams   PetscInt   numAMRRefine[LANDAU_MAX_GRIDS];  /* normal AMR - refine from origin */
1378a6f2e61SMark Adams   PetscInt   postAMRRefine[LANDAU_MAX_GRIDS]; /* uniform refinement of AMR */
1388a6f2e61SMark Adams   /* relativistic */
139cefb98e8SMark Adams   PetscBool  use_energy_tensor_trick;
140cefb98e8SMark Adams   PetscBool  use_relativistic_corrections;
141e0eea495SMark   /* physics */
142e0eea495SMark   PetscReal  thermal_temps[LANDAU_MAX_SPECIES];
143e0eea495SMark   PetscReal  masses[LANDAU_MAX_SPECIES];  /* mass of each species  */
144e0eea495SMark   PetscReal  charges[LANDAU_MAX_SPECIES]; /* charge of each species  */
145e0eea495SMark   PetscReal  n[LANDAU_MAX_SPECIES];       /* number density of each species  */
146e0eea495SMark   PetscReal  m_0;                         /* reference mass */
147e0eea495SMark   PetscReal  v_0;                         /* reference velocity */
148e0eea495SMark   PetscReal  n_0;                         /* reference number density */
149e0eea495SMark   PetscReal  t_0;                         /* reference time */
150e0eea495SMark   PetscReal  Ez;
151e0eea495SMark   PetscReal  epsilon0;
152e0eea495SMark   PetscReal  k;
153e0eea495SMark   PetscReal  lnLam;
1548a6f2e61SMark Adams   PetscReal  electronShift;
155e0eea495SMark   PetscInt   num_species;
1568a6f2e61SMark Adams   PetscInt   num_grids;
1578fdabdddSMark Adams   PetscInt   species_offset[LANDAU_MAX_GRIDS + 1]; // for each grid, but same for all batched vertices
1588fdabdddSMark Adams   PetscInt   mat_offset[LANDAU_MAX_GRIDS + 1];     // for each grid, but same for all batched vertices
159cb25d741SMark Adams   // batching
160cb25d741SMark Adams   PetscBool  jacobian_field_major_order; // this could be a type but lets not get pedantic
161cb25d741SMark Adams   VecScatter plex_batch;
162cb25d741SMark Adams   Vec        work_vec;
163cb25d741SMark Adams   IS         batch_is;
164cb25d741SMark Adams   PetscErrorCode (*seqaij_mult)(Mat, Vec, Vec);
165cb25d741SMark Adams   PetscErrorCode (*seqaij_multtranspose)(Mat, Vec, Vec);
166cb25d741SMark Adams   PetscErrorCode (*seqaij_solve)(Mat, Vec, Vec);
167cb25d741SMark Adams   PetscErrorCode (*seqaij_getdiagonal)(Mat, Vec);
168bfc784b7SMark Adams   /* COO */
169aab769cbSMark Adams   PetscBool        coo_assembly;
170e0eea495SMark   /* cache */
171e0eea495SMark   Mat              J;
172e0eea495SMark   Mat              M;
173e0eea495SMark   Vec              X;
174e0eea495SMark   /* derived type */
175e0eea495SMark   void            *data;
176e0eea495SMark   /* computing */
177e0eea495SMark   LandauDeviceType deviceType;
1788a6f2e61SMark Adams   DM               pack;
1798a6f2e61SMark Adams   DM               plex[LANDAU_MAX_GRIDS];
1808fdabdddSMark Adams   LandauStaticData SData_d; /* static geometric data on device */
18154545eeeSMark Adams   /* diagnostics */
18254545eeeSMark Adams   PetscInt         verbose;
18354545eeeSMark Adams   PetscLogEvent    events[20];
184c751c0a2SMark Adams   PetscLogStage    stage;
185984ed092SMark Adams   PetscObjectState norm_state;
1868fdabdddSMark Adams   PetscInt         batch_sz;
1878fdabdddSMark Adams   PetscInt         batch_view_idx;
188e0eea495SMark } LandauCtx;
189e0eea495SMark 
1908fdabdddSMark Adams #define LANDAU_SPECIES_MAJOR
1918fdabdddSMark Adams #if !defined(LANDAU_SPECIES_MAJOR)
1928fdabdddSMark Adams #define LAND_PACK_IDX(_b, _g)                         (_b * ctx->num_grids + _g)
1938fdabdddSMark Adams #define LAND_MOFFSET(_b, _g, _nbch, _ngrid, _mat_off) (_b * _mat_off[_ngrid] + _mat_off[_g])
1948fdabdddSMark Adams #else
1958fdabdddSMark Adams #define LAND_PACK_IDX(_b, _g)                         (_g * ctx->batch_sz + _b)
1968fdabdddSMark Adams #define LAND_MOFFSET(_b, _g, _nbch, _ngrid, _mat_off) (_nbch * _mat_off[_g] + _b * (_mat_off[_g + 1] - _mat_off[_g]))
1978fdabdddSMark Adams #endif
1988fdabdddSMark Adams 
199a587d139SMark typedef struct {
200a587d139SMark   PetscReal scale;
201bfc784b7SMark Adams   LandauIdx gid; // Landau matrix index (<10,000)
202a587d139SMark } pointInterpolationP4est;
203a587d139SMark typedef struct _lP4estVertexMaps {
204a587d139SMark   LandauIdx (*gIdx)[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ]; // #elems *  LANDAU_MAX_NQ (spoof for max , Nb) on device,
205a587d139SMark   LandauIdx        num_elements;
206a587d139SMark   LandauIdx        num_reduced;
207a587d139SMark   LandauIdx        num_face; // (Q or Q^2 for 3D)
208a587d139SMark   LandauDeviceType deviceType;
209a587d139SMark   PetscInt         Nf;
210a587d139SMark   pointInterpolationP4est (*c_maps)[LANDAU_MAX_Q_FACE];
2118a6f2e61SMark Adams   struct _lP4estVertexMaps *d_self;
212a587d139SMark   void                     *vp1, *vp2, *vp3;
2138a6f2e61SMark Adams   PetscInt                  numgrids;
214a587d139SMark } P4estVertexMaps;
215a587d139SMark 
216e0eea495SMark PETSC_EXTERN PetscErrorCode LandauCreateColoring(Mat, DM, PetscContainer *);
217e0eea495SMark #if defined(PETSC_HAVE_CUDA)
218*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode LandauCUDAJacobian(DM[], const PetscInt, const PetscInt, const PetscInt, const PetscInt[], PetscReal[], PetscScalar[], const PetscScalar[], const LandauStaticData *, const PetscReal, const PetscLogEvent[], const PetscInt[], const PetscInt[], Mat[], Mat);
2198a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauCUDACreateMatMaps(P4estVertexMaps *, pointInterpolationP4est (*)[LANDAU_MAX_Q_FACE], PetscInt[], PetscInt, PetscInt);
2208a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauCUDADestroyMatMaps(P4estVertexMaps *, PetscInt);
221*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode LandauCUDAStaticDataSet(DM, const PetscInt, const PetscInt, const PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], LandauStaticData *);
2228a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauCUDAStaticDataClear(LandauStaticData *);
223e0eea495SMark #endif
224e0eea495SMark #if defined(PETSC_HAVE_KOKKOS)
225*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode LandauKokkosJacobian(DM[], const PetscInt, const PetscInt, const PetscInt, const PetscInt[], PetscReal[], PetscScalar[], const PetscScalar[], const LandauStaticData *, const PetscReal, const PetscLogEvent[], const PetscInt[], const PetscInt[], Mat[], Mat);
2268a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauKokkosCreateMatMaps(P4estVertexMaps *, pointInterpolationP4est (*)[LANDAU_MAX_Q_FACE], PetscInt[], PetscInt, PetscInt);
2278a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauKokkosDestroyMatMaps(P4estVertexMaps *, PetscInt);
228*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode LandauKokkosStaticDataSet(DM, const PetscInt, const PetscInt, const PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], LandauStaticData *);
2298a6f2e61SMark Adams PETSC_EXTERN PetscErrorCode LandauKokkosStaticDataClear(LandauStaticData *);
230e0eea495SMark #endif
231e0eea495SMark 
232e0eea495SMark #endif /* PETSCLANDAU_H */
233