xref: /petsc/include/petsclandau.h (revision cefb98e87952b36bfd83b9ff5af63d155f2c71df)
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*);
10e0eea495SMark PETSC_EXTERN PetscErrorCode LandauAddMaxwellians(DM, Vec, PetscReal, PetscReal[], PetscReal[], 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
2352cdd6eaSMark #else
2452cdd6eaSMark #define LANDAU_MAX_SPECIES 3
2552cdd6eaSMark #endif
26e0eea495SMark #endif
27a587d139SMark 
28a587d139SMark #if !defined(LANDAU_MAX_Q)
29a587d139SMark #if defined(LANDAU_MAX_NQ)
30a587d139SMark #error"LANDAU_MAX_NQ but not LANDAU_MAX_Q. Use -DLANDAU_MAX_Q=4 for Q3 elements"
31a587d139SMark #endif
3252cdd6eaSMark #if LANDAU_DIM==2
33a587d139SMark #define LANDAU_MAX_Q 5
34e0eea495SMark #else
35a587d139SMark #define LANDAU_MAX_Q 3
36a587d139SMark #endif
37a587d139SMark #else
38a587d139SMark #undef LANDAU_MAX_NQ
39e0eea495SMark #endif
40a587d139SMark 
41a587d139SMark #if LANDAU_DIM==2
42930e68a5SMark Adams #define LANDAU_MAX_Q_FACE LANDAU_MAX_Q
43a587d139SMark #define LANDAU_MAX_NQ (LANDAU_MAX_Q*LANDAU_MAX_Q)
44a587d139SMark #else
45930e68a5SMark Adams #define LANDAU_MAX_Q_FACE (LANDAU_MAX_Q*LANDAU_MAX_Q)
46a587d139SMark #define LANDAU_MAX_NQ (LANDAU_MAX_Q*LANDAU_MAX_Q*LANDAU_MAX_Q)
47a587d139SMark #endif
48a587d139SMark 
49e0eea495SMark typedef enum {LANDAU_CUDA, LANDAU_KOKKOS, LANDAU_CPU} LandauDeviceType;
50e8d2b73aSMark Adams 
51e8d2b73aSMark Adams /* typedef PetscReal LandauIPReal; */
52e8d2b73aSMark Adams /* typedef struct { */
53e8d2b73aSMark Adams /*   LandauIPReal  *coefs; */
54e8d2b73aSMark Adams /*   int           dim_,ns_,nip_; */
55e8d2b73aSMark Adams /* } LandauIPFdF; */
56e8d2b73aSMark Adams 
57e8d2b73aSMark Adams typedef struct {
58e8d2b73aSMark Adams   void  *invJ;  // nip*dim*dim
59e8d2b73aSMark Adams   void  *D;     // nq*nb*dim
60e8d2b73aSMark Adams   void  *B;     // nq*nb
61e8d2b73aSMark Adams   void  *alpha; // ns
62e8d2b73aSMark Adams   void  *beta;  // ns
63e8d2b73aSMark Adams   void  *invMass; // ns
64e8d2b73aSMark Adams   void  *mass_w;  // nip
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   void  *IPf;  // Ncells*Nb*Nf - dynamic (vertex in cells)
75e8d2b73aSMark Adams   int   dim_,ns_,nip_,nq_,nb_;
76e8d2b73aSMark Adams } LandauGeomData;
77e8d2b73aSMark Adams 
78e0eea495SMark typedef struct {
79e0eea495SMark   PetscBool      interpolate;                  /* Generate intermediate mesh elements */
80a587d139SMark   PetscBool      gpu_assembly;
81e0eea495SMark   PetscFE        fe[LANDAU_MAX_SPECIES];
82e0eea495SMark   /* geometry  */
83e0eea495SMark   PetscReal      i_radius;
84e0eea495SMark   PetscReal      e_radius;
85e0eea495SMark   PetscInt       num_sections;
86e0eea495SMark   PetscReal      radius;
87e0eea495SMark   PetscReal      re_radius;           /* radius of refinement along v_perp=0, z>0 */
88e0eea495SMark   PetscReal      vperp0_radius1;      /* radius of refinement along v_perp=0 */
89e0eea495SMark   PetscReal      vperp0_radius2;      /* radius of refinement along v_perp=0 after origin AMR refinement */
90e0eea495SMark   PetscBool      sphere;
91e0eea495SMark   PetscBool      inflate;
92e0eea495SMark   PetscInt       numRERefine;       /* refinement along v_perp=0, z > 0 */
93e0eea495SMark   PetscInt       nZRefine1;          /* origin refinement after v_perp=0 refinement */
94e0eea495SMark   PetscInt       nZRefine2;          /* origin refinement after origin AMR refinement */
95e0eea495SMark   PetscInt       maxRefIts;         /* normal AMR - refine from origin */
96e0eea495SMark   PetscInt       postAMRRefine;     /* uniform refinement of AMR */
97*cefb98e8SMark Adams   PetscInt       preAMRRefine;     /* uniform refinement of AMR */
98e0eea495SMark   /* discretization - AMR */
99e0eea495SMark   PetscErrorCode (*errorIndicator)(PetscInt, PetscReal, PetscReal [], PetscInt, const PetscInt[], const PetscScalar[], const PetscScalar[], PetscReal *, void *);
100e0eea495SMark   PetscReal      refineTol[LANDAU_MAX_SPECIES];
101e0eea495SMark   PetscReal      coarsenTol[LANDAU_MAX_SPECIES];
102*cefb98e8SMark Adams   PetscBool      use_energy_tensor_trick;
103*cefb98e8SMark Adams   PetscBool      use_relativistic_corrections;
104e0eea495SMark   /* physics */
105e0eea495SMark   PetscReal      thermal_temps[LANDAU_MAX_SPECIES];
106e0eea495SMark   PetscReal      masses[LANDAU_MAX_SPECIES];  /* mass of each species  */
107e0eea495SMark   PetscReal      charges[LANDAU_MAX_SPECIES]; /* charge of each species  */
108e0eea495SMark   PetscReal      n[LANDAU_MAX_SPECIES];       /* number density of each species  */
109e0eea495SMark   PetscReal      m_0;      /* reference mass */
110e0eea495SMark   PetscReal      v_0;      /* reference velocity */
111e0eea495SMark   PetscReal      n_0;      /* reference number density */
112e0eea495SMark   PetscReal      t_0;      /* reference time */
113e0eea495SMark   PetscReal      Ez;
114e0eea495SMark   PetscReal      epsilon0;
115e0eea495SMark   PetscReal      k;
116e0eea495SMark   PetscReal      lnLam;
117e0eea495SMark   PetscReal      electronShift; /* for tests */
118e0eea495SMark   PetscInt       num_species;
119e0eea495SMark   /* cache */
120e0eea495SMark   Mat            J;
121e0eea495SMark   Mat            M;
122e0eea495SMark   Vec            X;
123e0eea495SMark   /* derived type */
124e0eea495SMark   void          *data;
125e0eea495SMark   PetscBool      aux_bool;  /* helper */
126e0eea495SMark   /* computing */
127e0eea495SMark   LandauDeviceType deviceType;
128bddcd29dSMark Adams   PetscInt       subThreadBlockSize;
129bddcd29dSMark Adams   PetscInt       numConcurrency; /* number of SMs in Cuda to use */
130930e68a5SMark Adams   MPI_Comm       comm; /* global communicator to use for errors and diagnostics */
131e8d2b73aSMark Adams   LandauGeomData *SData_d; /* static geometric data on device, but this pointer is a host pointer */
132e8d2b73aSMark Adams   double         times[1];
133e8d2b73aSMark Adams   PetscBool      init;
134e8d2b73aSMark Adams   PetscBool      use_matrix_mass;
13554545eeeSMark Adams   DM             dmv;
13654545eeeSMark Adams   DM             plex;
13754545eeeSMark Adams   /* diagnostics */
13854545eeeSMark Adams   PetscInt       verbose;
13954545eeeSMark Adams   PetscLogEvent  events[20];
140e0eea495SMark } LandauCtx;
141e0eea495SMark 
142a587d139SMark typedef int LandauIdx;
143a587d139SMark typedef struct {
144a587d139SMark   PetscReal scale;
145a587d139SMark   LandauIdx gid;   // Lanadu matrix index (<10,000)
146a587d139SMark } pointInterpolationP4est;
147a587d139SMark typedef struct _lP4estVertexMaps {
148a587d139SMark   LandauIdx                (*gIdx)[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ]; // #elems *  LANDAU_MAX_NQ (spoof for max , Nb) on device,
149a587d139SMark   LandauIdx                num_elements;
150a587d139SMark   LandauIdx                num_reduced;
151a587d139SMark   LandauIdx                num_face;  // (Q or Q^2 for 3D)
152a587d139SMark   LandauDeviceType         deviceType;
153a587d139SMark   PetscInt                 Nf;
154a587d139SMark   PetscInt                 Nq;
155a587d139SMark   pointInterpolationP4est (*c_maps)[LANDAU_MAX_Q_FACE];
156a587d139SMark   struct _lP4estVertexMaps*data;
157a587d139SMark   void                    *vp1,*vp2,*vp3;
158a587d139SMark } P4estVertexMaps;
159a587d139SMark 
160e0eea495SMark PETSC_EXTERN PetscErrorCode LandauCreateColoring(Mat, DM, PetscContainer *);
161e0eea495SMark #if defined(PETSC_HAVE_CUDA)
162b8d122ceSMark Adams PETSC_EXTERN PetscErrorCode LandauCUDAJacobian(DM, const PetscInt, PetscReal[], PetscScalar[], const PetscInt, const PetscScalar[], LandauGeomData *, const PetscInt, PetscReal, const PetscLogEvent[], Mat);
163a587d139SMark PETSC_EXTERN PetscErrorCode LandauCUDACreateMatMaps(P4estVertexMaps *, pointInterpolationP4est (*)[LANDAU_MAX_Q_FACE], PetscInt, PetscInt);
164a587d139SMark PETSC_EXTERN PetscErrorCode LandauCUDADestroyMatMaps(P4estVertexMaps *);
165b8d122ceSMark Adams PETSC_EXTERN PetscErrorCode LandauCUDAStaticDataSet(DM, const PetscInt, PetscReal [], PetscReal [], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], LandauGeomData *);
166b8d122ceSMark Adams PETSC_EXTERN PetscErrorCode LandauCUDAStaticDataClear(LandauGeomData *);
167e0eea495SMark #endif
168e0eea495SMark #if defined(PETSC_HAVE_KOKKOS)
169af0767daSMark Adams PETSC_EXTERN PetscErrorCode LandauKokkosJacobian(DM, const PetscInt, PetscReal[], PetscScalar[],  const PetscInt, const PetscScalar[], LandauGeomData *, const PetscInt, PetscReal, const PetscLogEvent[], Mat);
170a587d139SMark PETSC_EXTERN PetscErrorCode LandauKokkosCreateMatMaps(P4estVertexMaps *, pointInterpolationP4est (*)[LANDAU_MAX_Q_FACE], PetscInt, PetscInt);
171a587d139SMark PETSC_EXTERN PetscErrorCode LandauKokkosDestroyMatMaps(P4estVertexMaps *);
172e8d2b73aSMark Adams PETSC_EXTERN PetscErrorCode LandauKokkosStaticDataSet(DM, const PetscInt, PetscReal [], PetscReal [], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], LandauGeomData *);
173e8d2b73aSMark Adams PETSC_EXTERN PetscErrorCode LandauKokkosStaticDataClear(LandauGeomData *);
174e0eea495SMark #endif
175e0eea495SMark 
176e0eea495SMark #endif /* PETSCLANDAU_H */
177