xref: /petsc/include/petsclandau.h (revision 984ed092ecc7313e63e6aa733432dfb60d929fa3)
1 #if !defined(PETSCLANDAU_H)
2 #define PETSCLANDAU_H
3 
4 #include <petscdmplex.h> /*I      "petscdmplex.h"    I*/
5 #include <petscts.h>
6 
7 PETSC_EXTERN PetscErrorCode DMPlexLandauPrintNorms(Vec, PetscInt);
8 PETSC_EXTERN PetscErrorCode DMPlexLandauCreateVelocitySpace(MPI_Comm,PetscInt,const char[],Vec*,Mat*,DM*);
9 PETSC_EXTERN PetscErrorCode DMPlexLandauDestroyVelocitySpace(DM*);
10 PETSC_EXTERN PetscErrorCode DMPlexLandauAddMaxwellians(DM, Vec, PetscReal, PetscReal[], PetscReal[], PetscInt, PetscInt, void *);
11 PETSC_EXTERN PetscErrorCode DMPlexLandauCreateMassMatrix(DM dm, Mat *Amat);
12 PETSC_EXTERN PetscErrorCode DMPlexLandauIFunction(TS, PetscReal,Vec,Vec,Vec,void *);
13 PETSC_EXTERN PetscErrorCode DMPlexLandauIJacobian(TS, PetscReal,Vec,Vec,PetscReal,Mat,Mat,void *);
14 
15 typedef int LandauIdx;
16 
17 /* the Fokker-Planck-Landau context */
18 #if !defined(LANDAU_MAX_SPECIES)
19 #if defined(PETSC_USE_DMLANDAU_2D)
20 #define LANDAU_MAX_SPECIES 10
21 #define LANDAU_MAX_GRIDS 3
22 #else
23 #define LANDAU_MAX_SPECIES 3
24 #define LANDAU_MAX_GRIDS 3
25 #endif
26 #else
27 #define LANDAU_MAX_GRIDS 3
28 #endif
29 
30 #if !defined(LANDAU_MAX_Q)
31 #if defined(LANDAU_MAX_NQ)
32 #error"LANDAU_MAX_NQ but not LANDAU_MAX_Q. Use -DLANDAU_MAX_Q=4 for Q3 elements"
33 #endif
34 #if defined(PETSC_USE_DMLANDAU_2D)
35 #define LANDAU_MAX_Q 5
36 #else
37 #define LANDAU_MAX_Q 3
38 #endif
39 #else
40 #undef LANDAU_MAX_NQ
41 #endif
42 
43 #if defined(PETSC_USE_DMLANDAU_2D)
44 #define LANDAU_MAX_Q_FACE LANDAU_MAX_Q
45 #define LANDAU_MAX_NQ (LANDAU_MAX_Q*LANDAU_MAX_Q)
46 #define LANDAU_MAX_BATCH_SZ 320
47 #define LANDAU_DIM 2
48 #else
49 #define LANDAU_MAX_Q_FACE (LANDAU_MAX_Q*LANDAU_MAX_Q)
50 #define LANDAU_MAX_NQ (LANDAU_MAX_Q*LANDAU_MAX_Q*LANDAU_MAX_Q)
51 #define LANDAU_MAX_BATCH_SZ 32
52 #define LANDAU_DIM 3
53 #endif
54 
55 typedef enum {LANDAU_CUDA, LANDAU_KOKKOS, LANDAU_CPU} LandauDeviceType;
56 
57 // static data - will be "device" data
58 typedef struct {
59   void  *invJ;  // nip*dim*dim
60   void  *D;     // nq*nb*dim
61   void  *B;     // nq*nb
62   void  *alpha; // ns
63   void  *beta;  // ns
64   void  *invMass; // ns
65   void  *w; // nip
66   void  *x; // nip
67   void  *y; // nip
68   void  *z; // nip
69   void  *Eq_m; // ns - dynamic
70   void  *f; //  nip*Nf - dynamic (IP)
71   void  *dfdx; // nip*Nf - dynamic (IP)
72   void  *dfdy; // nip*Nf - dynamic (IP)
73   void  *dfdz; // nip*Nf - dynamic (IP)
74   int   dim_,ns_,nip_,nq_,nb_;
75   void  *NCells; // remove and ise elem_offset - TODO
76   void  *species_offset; // for each grid, but same for all batched vertices
77   void  *mat_offset; // for each grid, but same for all batched vertices
78   void  *elem_offset; // for each grid, but same for all batched vertices
79   void  *ip_offset; // for each grid, but same for all batched vertices
80   void  *ipf_offset; // for each grid, but same for all batched vertices
81   void  *ipfdf_data; // for each grid, but same for all batched vertices
82   void  *maps; // for each grid, but same for all batched vertices
83   // COO
84   void             *coo_elem_offsets;
85   void             *coo_elem_point_offsets;
86   void             *coo_elem_fullNb;
87   LandauIdx        coo_n_cellsTot;
88   LandauIdx        coo_size;
89   LandauIdx        coo_max_fullnb;
90 } LandauStaticData;
91 
92 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;
93 
94 typedef struct {
95   PetscBool      interpolate;                  /* Generate intermediate mesh elements */
96   PetscBool      gpu_assembly;
97   MPI_Comm       comm; /* global communicator to use for errors and diagnostics */
98   double         times[LANDAU_NUM_TIMERS];
99   PetscBool      use_matrix_mass;
100   /* FE */
101   PetscFE        fe[LANDAU_MAX_SPECIES];
102   /* geometry  */
103   PetscReal      radius[LANDAU_MAX_GRIDS];
104   PetscReal      re_radius;           /* RE: radius of refinement along v_perp=0, z>0 */
105   PetscReal      vperp0_radius1;      /* RE: radius of refinement along v_perp=0 */
106   PetscReal      vperp0_radius2;      /* RE: radius of refinement along v_perp=0 after origin AMR refinement */
107   PetscBool      sphere; // not used
108   PetscBool      inflate; // not used
109   PetscReal      i_radius[LANDAU_MAX_GRIDS]; // not used
110   PetscReal      e_radius; // not used
111   PetscInt       num_sections; // not used
112   /* AMR */
113   PetscBool      use_p4est;
114   PetscInt       numRERefine;        /* RE: refinement along v_perp=0, z > 0 */
115   PetscInt       nZRefine1;          /* RE: origin refinement after v_perp=0 refinement */
116   PetscInt       nZRefine2;          /* RE: origin refinement after origin AMR refinement */
117   PetscInt       numAMRRefine[LANDAU_MAX_GRIDS];  /* normal AMR - refine from origin */
118   PetscInt       postAMRRefine[LANDAU_MAX_GRIDS]; /* uniform refinement of AMR */
119   /* relativistic */
120   PetscBool      use_energy_tensor_trick;
121   PetscBool      use_relativistic_corrections;
122   /* physics */
123   PetscReal      thermal_temps[LANDAU_MAX_SPECIES];
124   PetscReal      masses[LANDAU_MAX_SPECIES];  /* mass of each species  */
125   PetscReal      charges[LANDAU_MAX_SPECIES]; /* charge of each species  */
126   PetscReal      n[LANDAU_MAX_SPECIES];       /* number density of each species  */
127   PetscReal      m_0;      /* reference mass */
128   PetscReal      v_0;      /* reference velocity */
129   PetscReal      n_0;      /* reference number density */
130   PetscReal      t_0;      /* reference time */
131   PetscReal      Ez;
132   PetscReal      epsilon0;
133   PetscReal      k;
134   PetscReal      lnLam;
135   PetscReal      electronShift;
136   PetscInt       num_species;
137   PetscInt       num_grids;
138   PetscInt       species_offset[LANDAU_MAX_GRIDS+1]; // for each grid, but same for all batched vertices
139   PetscInt       mat_offset[LANDAU_MAX_GRIDS+1]; // for each grid, but same for all batched vertices
140   // batching
141   PetscBool      jacobian_field_major_order; // this could be a type but lets not get pedantic
142   VecScatter     plex_batch;
143   Vec            work_vec;
144   IS             batch_is;
145   PetscErrorCode (*seqaij_mult)(Mat,Vec,Vec);
146   PetscErrorCode (*seqaij_multtranspose)(Mat,Vec,Vec);
147   PetscErrorCode (*seqaij_solve)(Mat,Vec,Vec);
148   PetscErrorCode (*seqaij_getdiagonal)(Mat,Vec);
149   /* COO */
150   PetscBool        coo_assembly;
151   /* cache */
152   Mat              J;
153   Mat              M;
154   Vec              X;
155   /* derived type */
156   void             *data;
157   /* computing */
158   LandauDeviceType deviceType;
159   DM               pack;
160   DM               plex[LANDAU_MAX_GRIDS];
161   LandauStaticData SData_d; /* static geometric data on device */
162   /* diagnostics */
163   PetscInt         verbose;
164   PetscLogEvent    events[20];
165   PetscLogStage    stage;
166   PetscObjectState norm_state;
167   PetscInt         batch_sz;
168   PetscInt         batch_view_idx;
169 } LandauCtx;
170 
171 #define LANDAU_SPECIES_MAJOR
172 #if !defined(LANDAU_SPECIES_MAJOR)
173 #define LAND_PACK_IDX(_b,_g) (_b*ctx->num_grids + _g)
174 #define LAND_MOFFSET(_b,_g,_nbch,_ngrid,_mat_off) (_b*_mat_off[_ngrid] + _mat_off[_g])
175 #else
176 #define LAND_PACK_IDX(_b,_g) (_g*ctx->batch_sz + _b)
177 #define LAND_MOFFSET(_b,_g,_nbch,_ngrid,_mat_off) (_nbch*_mat_off[_g] + _b*(_mat_off[_g+1] - _mat_off[_g]))
178 #endif
179 
180 typedef struct {
181   PetscReal scale;
182   LandauIdx gid;   // Landau matrix index (<10,000)
183 } pointInterpolationP4est;
184 typedef struct _lP4estVertexMaps {
185   LandauIdx                (*gIdx)[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ]; // #elems *  LANDAU_MAX_NQ (spoof for max , Nb) on device,
186   LandauIdx                num_elements;
187   LandauIdx                num_reduced;
188   LandauIdx                num_face;  // (Q or Q^2 for 3D)
189   LandauDeviceType         deviceType;
190   PetscInt                 Nf;
191   pointInterpolationP4est (*c_maps)[LANDAU_MAX_Q_FACE];
192   struct _lP4estVertexMaps*d_self;
193   void                    *vp1,*vp2,*vp3;
194   PetscInt                numgrids;
195 } P4estVertexMaps;
196 
197 PETSC_EXTERN PetscErrorCode LandauCreateColoring(Mat, DM, PetscContainer *);
198 #if defined(PETSC_HAVE_CUDA)
199 PETSC_EXTERN PetscErrorCode LandauCUDAJacobian(DM[], const PetscInt, const PetscInt, const PetscInt, const PetscInt[], PetscReal[], PetscScalar[], const PetscScalar[],
200                                                const LandauStaticData *, const PetscReal, const PetscLogEvent[],const PetscInt[], const PetscInt[], Mat[], Mat);
201 PETSC_EXTERN PetscErrorCode LandauCUDACreateMatMaps(P4estVertexMaps *, pointInterpolationP4est (*)[LANDAU_MAX_Q_FACE], PetscInt[], PetscInt, PetscInt);
202 PETSC_EXTERN PetscErrorCode LandauCUDADestroyMatMaps(P4estVertexMaps *, PetscInt);
203 PETSC_EXTERN PetscErrorCode LandauCUDAStaticDataSet(DM, const PetscInt, const PetscInt, const PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscReal [], PetscReal [], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[],
204                                                     PetscReal[], LandauStaticData *);
205 PETSC_EXTERN PetscErrorCode LandauCUDAStaticDataClear(LandauStaticData *);
206 #endif
207 #if defined(PETSC_HAVE_KOKKOS)
208 PETSC_EXTERN PetscErrorCode LandauKokkosJacobian(DM[], const PetscInt, const PetscInt, const PetscInt, const PetscInt[], PetscReal[], PetscScalar[],  const PetscScalar[],
209                                                  const LandauStaticData *, const PetscReal, const PetscLogEvent[], const PetscInt[], const PetscInt[], Mat[], Mat);
210 PETSC_EXTERN PetscErrorCode LandauKokkosCreateMatMaps(P4estVertexMaps *, pointInterpolationP4est (*)[LANDAU_MAX_Q_FACE], PetscInt[], PetscInt, PetscInt);
211 PETSC_EXTERN PetscErrorCode LandauKokkosDestroyMatMaps(P4estVertexMaps *, PetscInt);
212 PETSC_EXTERN PetscErrorCode LandauKokkosStaticDataSet(DM, const PetscInt, const PetscInt, const PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscReal [], PetscReal [], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[],
213                                                       PetscReal[], LandauStaticData *);
214 PETSC_EXTERN PetscErrorCode LandauKokkosStaticDataClear(LandauStaticData*);
215 #endif
216 
217 #endif /* PETSCLANDAU_H */
218