xref: /petsc/include/petsclandau.h (revision e0eea49528ddcce32a775a862994a11c74fcac87)
1*e0eea495SMark #if !defined(PETSCLANDAU_H)
2*e0eea495SMark #define PETSCLANDAU_H
3*e0eea495SMark 
4*e0eea495SMark #include <petscdmplex.h> /*I      "petscdmplex.h"    I*/
5*e0eea495SMark #include <petscts.h>
6*e0eea495SMark 
7*e0eea495SMark PETSC_EXTERN PetscErrorCode LandauPrintNorms(Vec, PetscInt);
8*e0eea495SMark PETSC_EXTERN PetscErrorCode LandauCreateVelocitySpace(MPI_Comm,PetscInt,const char[],Vec*,Mat*,DM*);
9*e0eea495SMark PETSC_EXTERN PetscErrorCode LandauDestroyVelocitySpace(DM*);
10*e0eea495SMark PETSC_EXTERN PetscErrorCode LandauAddMaxwellians(DM, Vec, PetscReal, PetscReal[], PetscReal[], void *);
11*e0eea495SMark PETSC_EXTERN PetscErrorCode LandauCreateMassMatrix(DM dm, Mat *Amat);
12*e0eea495SMark PETSC_EXTERN PetscErrorCode LandauIFunction(TS, PetscReal,Vec,Vec,Vec,void *);
13*e0eea495SMark PETSC_EXTERN PetscErrorCode LandauIJacobian(TS, PetscReal,Vec,Vec,PetscReal,Mat,Mat,void *);
14*e0eea495SMark 
15*e0eea495SMark /* the Fokker-Planck-Landau context */
16*e0eea495SMark #if !defined(LANDAU_DIM)
17*e0eea495SMark #define LANDAU_DIM 2
18*e0eea495SMark #endif
19*e0eea495SMark #if !defined(LANDAU_MAX_SPECIES)
20*e0eea495SMark #define LANDAU_MAX_SPECIES 10
21*e0eea495SMark #endif
22*e0eea495SMark #if !defined(LANDAU_MAX_NQ)
23*e0eea495SMark #define LANDAU_MAX_NQ 25
24*e0eea495SMark #endif
25*e0eea495SMark #if !defined(LANDAU_MAX_SUB_THREAD_BLOCKS)
26*e0eea495SMark #if defined(PETSC_HAVE_CUDA)
27*e0eea495SMark #define LANDAU_MAX_SUB_THREAD_BLOCKS 4
28*e0eea495SMark #else
29*e0eea495SMark #define LANDAU_MAX_SUB_THREAD_BLOCKS 1
30*e0eea495SMark #endif
31*e0eea495SMark #endif
32*e0eea495SMark typedef enum {LANDAU_CUDA, LANDAU_KOKKOS, LANDAU_CPU} LandauDeviceType;
33*e0eea495SMark typedef struct {
34*e0eea495SMark   PetscBool      interpolate;                  /* Generate intermediate mesh elements */
35*e0eea495SMark   PetscBool      simplex;
36*e0eea495SMark   PetscFE        fe[LANDAU_MAX_SPECIES];
37*e0eea495SMark   /* geometry  */
38*e0eea495SMark   PetscReal      i_radius;
39*e0eea495SMark   PetscReal      e_radius;
40*e0eea495SMark   PetscInt       num_sections;
41*e0eea495SMark   PetscReal      radius;
42*e0eea495SMark   PetscReal      re_radius;           /* radius of refinement along v_perp=0, z>0 */
43*e0eea495SMark   PetscReal      vperp0_radius1;      /* radius of refinement along v_perp=0 */
44*e0eea495SMark   PetscReal      vperp0_radius2;      /* radius of refinement along v_perp=0 after origin AMR refinement */
45*e0eea495SMark   PetscBool      sphere;
46*e0eea495SMark   PetscBool      inflate;
47*e0eea495SMark   PetscInt       numRERefine;       /* refinement along v_perp=0, z > 0 */
48*e0eea495SMark   PetscInt       nZRefine1;          /* origin refinement after v_perp=0 refinement */
49*e0eea495SMark   PetscInt       nZRefine2;          /* origin refinement after origin AMR refinement */
50*e0eea495SMark   PetscInt       maxRefIts;         /* normal AMR - refine from origin */
51*e0eea495SMark   PetscInt       postAMRRefine;     /* uniform refinement of AMR */
52*e0eea495SMark   PetscBool      quarter3DDomain;   /* bilateral symetry, 1/4 x-y domain */
53*e0eea495SMark   /* discretization - AMR */
54*e0eea495SMark   PetscErrorCode (*errorIndicator)(PetscInt, PetscReal, PetscReal [], PetscInt, const PetscInt[], const PetscScalar[], const PetscScalar[], PetscReal *, void *);
55*e0eea495SMark   PetscReal      refineTol[LANDAU_MAX_SPECIES];
56*e0eea495SMark   PetscReal      coarsenTol[LANDAU_MAX_SPECIES];
57*e0eea495SMark   /* physics */
58*e0eea495SMark   PetscReal      thermal_temps[LANDAU_MAX_SPECIES];
59*e0eea495SMark   PetscReal      masses[LANDAU_MAX_SPECIES];  /* mass of each species  */
60*e0eea495SMark   PetscReal      charges[LANDAU_MAX_SPECIES]; /* charge of each species  */
61*e0eea495SMark   PetscReal      n[LANDAU_MAX_SPECIES];       /* number density of each species  */
62*e0eea495SMark   PetscReal      m_0;      /* reference mass */
63*e0eea495SMark   PetscReal      v_0;      /* reference velocity */
64*e0eea495SMark   PetscReal      n_0;      /* reference number density */
65*e0eea495SMark   PetscReal      t_0;      /* reference time */
66*e0eea495SMark   PetscReal      Ez;
67*e0eea495SMark   PetscReal      epsilon0;
68*e0eea495SMark   PetscReal      k;
69*e0eea495SMark   PetscReal      lnLam;
70*e0eea495SMark   PetscReal      electronShift; /* for tests */
71*e0eea495SMark   PetscInt       num_species;
72*e0eea495SMark   /* diagnostics */
73*e0eea495SMark   PetscInt       verbose;
74*e0eea495SMark   PetscLogEvent  events[20];
75*e0eea495SMark   DM             dmv;
76*e0eea495SMark   /* cache */
77*e0eea495SMark   Mat            J;
78*e0eea495SMark   Mat            M;
79*e0eea495SMark   Vec            X;
80*e0eea495SMark   PetscReal      normJ; /* used to see if function changed */
81*e0eea495SMark   /* derived type */
82*e0eea495SMark   void          *data;
83*e0eea495SMark   PetscBool      aux_bool;  /* helper */
84*e0eea495SMark   /* computing */
85*e0eea495SMark   LandauDeviceType deviceType;
86*e0eea495SMark   PetscInt       subThreadBlockSize;
87*e0eea495SMark } LandauCtx;
88*e0eea495SMark 
89*e0eea495SMark typedef struct {
90*e0eea495SMark   PetscReal     f;
91*e0eea495SMark   PetscReal     df[LANDAU_DIM];
92*e0eea495SMark } LandauFDF;
93*e0eea495SMark 
94*e0eea495SMark typedef struct {
95*e0eea495SMark   /* coordinate */
96*e0eea495SMark   PetscReal   crd[LANDAU_DIM];
97*e0eea495SMark   /* f; df data [Nc] */
98*e0eea495SMark   LandauFDF     fdf[LANDAU_MAX_SPECIES];
99*e0eea495SMark } LandauPointData;
100*e0eea495SMark 
101*e0eea495SMark PETSC_EXTERN PetscErrorCode LandauAssembleOpenMP(PetscInt cStart, PetscInt cEnd, PetscInt totDim, DM plex, PetscSection section, PetscSection globalSection, Mat JacP, PetscScalar elemMats[], PetscContainer container);
102*e0eea495SMark PETSC_EXTERN PetscErrorCode LandauCreateColoring(Mat, DM, PetscContainer *);
103*e0eea495SMark PETSC_EXTERN PetscErrorCode LandauFormJacobian_Internal(Vec, Mat, const PetscInt, void *);
104*e0eea495SMark #if defined(PETSC_HAVE_CUDA)
105*e0eea495SMark PETSC_EXTERN PetscErrorCode LandauCUDAJacobian(DM, const PetscInt, const PetscReal [], const PetscReal [], const PetscReal[], const PetscReal[],
106*e0eea495SMark                                              const PetscReal * const, const PetscReal[], const PetscReal [],const PetscInt, const PetscLogEvent[], PetscBool, Mat);
107*e0eea495SMark #endif
108*e0eea495SMark #if defined(PETSC_HAVE_KOKKOS)
109*e0eea495SMark /* TODO: this won't work if PETSc is built with C++ */
110*e0eea495SMark #if !defined(__cplusplus)
111*e0eea495SMark PETSC_EXTERN PetscErrorCode LandauKokkosJacobian(DM, const PetscInt, const PetscReal [], const PetscReal [], const PetscReal[], const PetscReal[],
112*e0eea495SMark                                                const PetscReal * const, const PetscReal[], const PetscReal [],const PetscInt, const PetscLogEvent[], PetscBool, Mat);
113*e0eea495SMark #endif
114*e0eea495SMark #endif
115*e0eea495SMark 
116*e0eea495SMark #endif /* PETSCLANDAU_H */
117