XGCa
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
DM_wrapper.hpp
Go to the documentation of this file.
1 #ifndef DM_WRAPPER_HPP
2 #define DM_WRAPPER_HPP
3 
4 #define VELOCITY_SPACE_DIMENSION 2
5 
6 #ifdef NO_PETSC
7 // To avoid having to check for PETSc elsewhere define an empty struct if no PETSc
8 struct DMWrapper{
9 };
10 #else
11 
12 #include <string>
13 #include <Kokkos_Core.hpp>
14 #include "space_settings.hpp"
15 #include <petscdm.h>
16 #include <petscversion.h>
17 #include <petscsys.h>
18 #include <petscdmplex.h>
19 #if PETSC_VERSION_GE(3,14,0)
20 #include <petscds.h>
21 #include "petscdmswarm.h"
22 #include "petscksp.h"
23 #endif
24 
25 // For older PETSc versions
26 #if !defined(PETSC_USE_STRICT_PETSCERRORCODE) && !defined(PETSC_SUCCESS)
27 #define PETSC_SUCCESS 0
28 #endif
29 
30 // For some reason Kokkos cannot handle Views of PETSc objects, but it works when putting them in a wrapper
31 // Pseudo-inverse mesh PETSc objects
32 struct DMObjects{
33  DM dm;
34  Mat MM;
35 };
36 
37 // Wrapper for PETSc DM object to facilitate cleanup in destructor
38 struct DMWrapper{
39  bool is_owner;
40  int dm_len;
41  Kokkos::View<DMObjects*, Kokkos::LayoutRight, HostType> obj;
42 
44  : is_owner(false),
45  dm_len(0){
46  Kokkos::resize(obj, 0);
47  }
48 
49  ~DMWrapper(){ // destructor
50  if(is_owner){
51  for (int id = 0; id < dm_len; id++){
52  DMDestroy(&obj(id).dm);
53  MatDestroy(&obj(id).MM);
54  }
55  }
56  }
57 
58  DMWrapper(const DMWrapper& original) // copy constructor
59  : is_owner(false),
60  dm_len(original.dm_len){
61  Kokkos::resize(obj, dm_len);
62  for (int id = 0; id < dm_len; id++){
63  obj(id).dm = original.obj(id).dm;
64  obj(id).MM = original.obj(id).MM;
65  }
66  }
67 
68  DMWrapper& operator=(const DMWrapper& original) { // copy assignment
69  if (this == &original){
70  return *this;
71  }
72  is_owner = false;
73  dm_len = original.dm_len;
74  Kokkos::resize(obj, dm_len);
75  for (int id = 0; id < dm_len; id++){
76  obj(id).dm = original.obj(id).dm;
77  obj(id).MM = original.obj(id).MM;
78  }
79  return *this;
80  }
81 
99  PetscErrorCode initialize_pseudo_inverse_interpolation(int nthreads,
100  PetscInt order, const PetscInt Nspecies,
101  const PetscReal vPara_min, const PetscReal vPara_max,
102  const PetscInt num_vPara, const PetscReal vPerp_min,
103  const PetscReal vPerp_max, const PetscInt num_vPerp)
104  {
105 #if PETSC_VERSION_LT(3,14,0)
106  // To be able to compile with older PETSc versions
107  printf("\n\n\nERROR: initialize_pseudo_inverse_interpolation (should never reach here)\n");
108  return 0;
109 #else
110  is_owner = true;
111  dm_len = nthreads;
112  Kokkos::resize(obj, dm_len);
113 
114  PetscErrorCode ierr;
115 
116  PetscInt cells[2];
117  PetscReal lo[2], hi[2];
118  cells[0] = num_vPara - 1; //2 grid points give 1 cell, etc
119  cells[1] = num_vPerp - 1;
120  lo[0] = vPara_min;
121  lo[1] = vPerp_min;
122  hi[0] = vPara_max;
123  hi[1] = vPerp_max;
124 
125  PetscFunctionBeginUser;
126 
127  std::string order_str = std::to_string(order);
128  ierr = PetscOptionsSetValue(NULL, "-petscspace_degree", order_str.c_str());CHKERRQ(ierr);
129 
130  /* create mesh */
131  PetscFE fe;
132  PetscInt field = 0, dim = VELOCITY_SPACE_DIMENSION;
133  DMBoundaryType periodicity[2] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
134  PetscBool interpolate = PETSC_TRUE;
135 
136  for (int id = 0; id < dm_len; id++){
137  ierr = DMPlexCreateBoxMesh(PETSC_COMM_SELF, dim, PETSC_FALSE, cells, lo, hi, periodicity, interpolate, &obj(id).dm);CHKERRQ(ierr);
138  ierr = PetscObjectSetName((PetscObject)obj(id).dm, "Potential Grid");CHKERRQ(ierr);
139  if (id == dm_len - 1) {ierr = DMViewFromOptions(obj(id).dm, NULL, "-dm_view");CHKERRQ(ierr);}
140  ierr = DMSetFromOptions(obj(id).dm);CHKERRQ(ierr);
141  ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, Nspecies, PETSC_FALSE, "", PETSC_DECIDE, &fe);CHKERRQ(ierr);
142  ierr = PetscFESetFromOptions(fe);CHKERRQ(ierr);
143  ierr = PetscObjectSetName((PetscObject)fe, "fe");CHKERRQ(ierr);
144  ierr = DMSetField(obj(id).dm, field, NULL, (PetscObject)fe);CHKERRQ(ierr);
145  ierr = DMCreateDS(obj(id).dm);CHKERRQ(ierr);
146  ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
147 
148  ierr = DMCreateMassMatrix(obj(id).dm, obj(id).dm, &obj(id).MM);CHKERRQ(ierr);
149  }
150 
151  PetscFunctionReturn(PETSC_SUCCESS);
152 #endif
153 //PETSC_VERSION_LT(3,14,0)
154  } //end initialize_pseudo_inverse_interpolation()
155 
162  PetscErrorCode destroy_DM(){
163  PetscErrorCode ierr;
164  PetscFunctionBeginUser;
165  if(is_owner){
166  for (int id = 0; id < dm_len; id++){
167  if(obj(id).dm != NULL){
168  ierr = DMDestroy(&obj(id).dm);CHKERRQ(ierr);
169  obj(id).dm = NULL;
170  }
171  if(obj(id).MM != NULL){
172  ierr = MatDestroy(&obj(id).MM);CHKERRQ(ierr);
173  obj(id).MM = NULL;
174  }
175  }
176  }
177  PetscFunctionReturn(PETSC_SUCCESS);
178  }
179 };
180 #endif
181 //NO_PETSC
182 
183 #endif
PetscErrorCode initialize_pseudo_inverse_interpolation(int nthreads, PetscInt order, const PetscInt Nspecies, const PetscReal vPara_min, const PetscReal vPara_max, const PetscInt num_vPara, const PetscReal vPerp_min, const PetscReal vPerp_max, const PetscInt num_vPerp)
Definition: DM_wrapper.hpp:99
DMWrapper()
Definition: DM_wrapper.hpp:43
DM dm
PETSc DM mesh object representing the velocity grid when pseudo-inverse is used.
Definition: DM_wrapper.hpp:33
#define VELOCITY_SPACE_DIMENSION
Definition: DM_wrapper.hpp:4
bool is_owner
Definition: DM_wrapper.hpp:39
~DMWrapper()
Definition: DM_wrapper.hpp:49
Definition: DM_wrapper.hpp:38
DMWrapper(const DMWrapper &original)
Definition: DM_wrapper.hpp:58
DMWrapper & operator=(const DMWrapper &original)
Definition: DM_wrapper.hpp:68
Kokkos::View< DMObjects *, Kokkos::LayoutRight, HostType > obj
Stores number of threads DMObjects.
Definition: DM_wrapper.hpp:41
PetscErrorCode destroy_DM()
Definition: DM_wrapper.hpp:162
#define PETSC_SUCCESS
Definition: DM_wrapper.hpp:27
int dm_len
Number of PETSc mesh objects to store.
Definition: DM_wrapper.hpp:40
Definition: DM_wrapper.hpp:32
Mat MM
PETSc mass matrix object for the velocity grid when pseudo-inverse is used.
Definition: DM_wrapper.hpp:34