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