xref: /petsc/include/petsclandau.h (revision a587d1398bd356a7db20a09283619c724e82d622)
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
19*a587d139SMark 
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
27*a587d139SMark 
28*a587d139SMark #if !defined(LANDAU_MAX_Q)
29*a587d139SMark #if defined(LANDAU_MAX_NQ)
30*a587d139SMark #error"LANDAU_MAX_NQ but not LANDAU_MAX_Q. Use -DLANDAU_MAX_Q=4 for Q3 elements"
31*a587d139SMark #endif
3252cdd6eaSMark #if LANDAU_DIM==2
33*a587d139SMark #define LANDAU_MAX_Q 5
34*a587d139SMark #define LANDAU_MAX_Q_FACE LANDAU_MAX_Q
35e0eea495SMark #else
36*a587d139SMark #define LANDAU_MAX_Q 3
37*a587d139SMark #define LANDAU_MAX_Q_FACE (LANDAU_MAX_Q*LANDAU_MAX_Q)
38*a587d139SMark #endif
39*a587d139SMark #else
40*a587d139SMark #undef LANDAU_MAX_NQ
41*a587d139SMark #if LANDAU_DIM==2
42*a587d139SMark #define LANDAU_MAX_NQ (LANDAU_MAX_Q*LANDAU_MAX_Q)
43*a587d139SMark #define LANDAU_MAX_Q_FACE LANDAU_MAX_Q
44*a587d139SMark #else
45*a587d139SMark #define LANDAU_MAX_NQ (LANDAU_MAX_Q*LANDAU_MAX_Q*LANDAU_MAX_Q)
46*a587d139SMark #define LANDAU_MAX_Q_FACE (LANDAU_MAX_Q*LANDAU_MAX_Q)
47e0eea495SMark #endif
48e0eea495SMark #endif
49*a587d139SMark 
50*a587d139SMark #if LANDAU_DIM==2
51*a587d139SMark #define LANDAU_MAX_NQ (LANDAU_MAX_Q*LANDAU_MAX_Q)
52*a587d139SMark #else
53*a587d139SMark #define LANDAU_MAX_NQ (LANDAU_MAX_Q*LANDAU_MAX_Q*LANDAU_MAX_Q)
54*a587d139SMark #endif
55*a587d139SMark 
56e0eea495SMark typedef enum {LANDAU_CUDA, LANDAU_KOKKOS, LANDAU_CPU} LandauDeviceType;
57e0eea495SMark typedef struct {
58e0eea495SMark   PetscBool      interpolate;                  /* Generate intermediate mesh elements */
59*a587d139SMark   PetscBool      gpu_assembly;
60e0eea495SMark   PetscFE        fe[LANDAU_MAX_SPECIES];
61e0eea495SMark   /* geometry  */
62e0eea495SMark   PetscReal      i_radius;
63e0eea495SMark   PetscReal      e_radius;
64e0eea495SMark   PetscInt       num_sections;
65e0eea495SMark   PetscReal      radius;
66e0eea495SMark   PetscReal      re_radius;           /* radius of refinement along v_perp=0, z>0 */
67e0eea495SMark   PetscReal      vperp0_radius1;      /* radius of refinement along v_perp=0 */
68e0eea495SMark   PetscReal      vperp0_radius2;      /* radius of refinement along v_perp=0 after origin AMR refinement */
69e0eea495SMark   PetscBool      sphere;
70e0eea495SMark   PetscBool      inflate;
71e0eea495SMark   PetscInt       numRERefine;       /* refinement along v_perp=0, z > 0 */
72e0eea495SMark   PetscInt       nZRefine1;          /* origin refinement after v_perp=0 refinement */
73e0eea495SMark   PetscInt       nZRefine2;          /* origin refinement after origin AMR refinement */
74e0eea495SMark   PetscInt       maxRefIts;         /* normal AMR - refine from origin */
75e0eea495SMark   PetscInt       postAMRRefine;     /* uniform refinement of AMR */
76e0eea495SMark   /* discretization - AMR */
77e0eea495SMark   PetscErrorCode (*errorIndicator)(PetscInt, PetscReal, PetscReal [], PetscInt, const PetscInt[], const PetscScalar[], const PetscScalar[], PetscReal *, void *);
78e0eea495SMark   PetscReal      refineTol[LANDAU_MAX_SPECIES];
79e0eea495SMark   PetscReal      coarsenTol[LANDAU_MAX_SPECIES];
80e0eea495SMark   /* physics */
81e0eea495SMark   PetscReal      thermal_temps[LANDAU_MAX_SPECIES];
82e0eea495SMark   PetscReal      masses[LANDAU_MAX_SPECIES];  /* mass of each species  */
83e0eea495SMark   PetscReal      charges[LANDAU_MAX_SPECIES]; /* charge of each species  */
84e0eea495SMark   PetscReal      n[LANDAU_MAX_SPECIES];       /* number density of each species  */
85e0eea495SMark   PetscReal      m_0;      /* reference mass */
86e0eea495SMark   PetscReal      v_0;      /* reference velocity */
87e0eea495SMark   PetscReal      n_0;      /* reference number density */
88e0eea495SMark   PetscReal      t_0;      /* reference time */
89e0eea495SMark   PetscReal      Ez;
90e0eea495SMark   PetscReal      epsilon0;
91e0eea495SMark   PetscReal      k;
92e0eea495SMark   PetscReal      lnLam;
93e0eea495SMark   PetscReal      electronShift; /* for tests */
94e0eea495SMark   PetscInt       num_species;
95e0eea495SMark   /* diagnostics */
96e0eea495SMark   PetscInt       verbose;
97e0eea495SMark   PetscLogEvent  events[20];
98e0eea495SMark   DM             dmv;
99e0eea495SMark   /* cache */
100e0eea495SMark   Mat            J;
101e0eea495SMark   Mat            M;
102e0eea495SMark   Vec            X;
103e0eea495SMark   PetscReal      normJ; /* used to see if function changed */
104e0eea495SMark   /* derived type */
105e0eea495SMark   void          *data;
106e0eea495SMark   PetscBool      aux_bool;  /* helper */
107e0eea495SMark   /* computing */
108e0eea495SMark   LandauDeviceType deviceType;
109e0eea495SMark   PetscInt       subThreadBlockSize;
110e0eea495SMark } LandauCtx;
111e0eea495SMark 
112*a587d139SMark typedef int LandauIdx;
113*a587d139SMark typedef struct {
114*a587d139SMark   PetscReal scale;
115*a587d139SMark   LandauIdx gid;   // Lanadu matrix index (<10,000)
116*a587d139SMark } pointInterpolationP4est;
117*a587d139SMark typedef struct _lP4estVertexMaps {
118*a587d139SMark   LandauIdx                (*gIdx)[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ]; // #elems *  LANDAU_MAX_NQ (spoof for max , Nb) on device,
119*a587d139SMark   LandauIdx                num_elements;
120*a587d139SMark   LandauIdx                num_reduced;
121*a587d139SMark   LandauIdx                num_face;  // (Q or Q^2 for 3D)
122*a587d139SMark   LandauDeviceType         deviceType;
123*a587d139SMark   PetscInt                 Nf;
124*a587d139SMark   PetscInt                 Nq;
125*a587d139SMark   pointInterpolationP4est (*c_maps)[LANDAU_MAX_Q_FACE];
126*a587d139SMark   struct _lP4estVertexMaps*data;
127*a587d139SMark   void                    *vp1,*vp2,*vp3;
128*a587d139SMark } P4estVertexMaps;
129*a587d139SMark 
13052cdd6eaSMark typedef PetscReal LandauIPReal;
131e0eea495SMark typedef struct {
132*a587d139SMark   LandauIPReal  *w;
13352cdd6eaSMark   LandauIPReal  *x;
13452cdd6eaSMark   LandauIPReal  *y;
13552cdd6eaSMark   LandauIPReal  *z;
136*a587d139SMark   LandauIPReal  *coefs;
13752cdd6eaSMark   int           dim_,ns_,nip_;
13852cdd6eaSMark } LandauIPData;
139e0eea495SMark 
140e0eea495SMark PETSC_EXTERN PetscErrorCode LandauCreateColoring(Mat, DM, PetscContainer *);
141e0eea495SMark PETSC_EXTERN PetscErrorCode LandauFormJacobian_Internal(Vec, Mat, const PetscInt, void *);
14252cdd6eaSMark PETSC_EXTERN int LandauGetIPDataSize(const LandauIPData *const);
143e0eea495SMark #if defined(PETSC_HAVE_CUDA)
144e0eea495SMark PETSC_EXTERN PetscErrorCode LandauCUDAJacobian(DM, const PetscInt, const PetscReal [], const PetscReal [], const PetscReal[], const PetscReal[],
145*a587d139SMark   const LandauIPData *const, const PetscReal [], const PetscLogEvent[], Mat);
146*a587d139SMark PETSC_EXTERN PetscErrorCode LandauCUDACreateMatMaps(P4estVertexMaps *, pointInterpolationP4est (*)[LANDAU_MAX_Q_FACE], PetscInt, PetscInt);
147*a587d139SMark PETSC_EXTERN PetscErrorCode LandauCUDADestroyMatMaps(P4estVertexMaps *);
148*a587d139SMark 
149e0eea495SMark #endif
150e0eea495SMark #if defined(PETSC_HAVE_KOKKOS)
151e0eea495SMark   /* TODO: this won't work if PETSc is built with C++ */
152e0eea495SMark #if !defined(__cplusplus)
153e0eea495SMark PETSC_EXTERN PetscErrorCode LandauKokkosJacobian(DM, const PetscInt, const PetscReal [], const PetscReal [], const PetscReal[], const PetscReal[],
15452cdd6eaSMark                                                  const LandauIPData *const, const PetscReal [],const PetscInt, const PetscLogEvent[], Mat);
155*a587d139SMark PETSC_EXTERN PetscErrorCode LandauKokkosCreateMatMaps(P4estVertexMaps *, pointInterpolationP4est (*)[LANDAU_MAX_Q_FACE], PetscInt, PetscInt);
156*a587d139SMark PETSC_EXTERN PetscErrorCode LandauKokkosDestroyMatMaps(P4estVertexMaps *);
157e0eea495SMark #endif
158e0eea495SMark #endif
159e0eea495SMark 
160e0eea495SMark #endif /* PETSCLANDAU_H */
161