4 #define VELOCITY_SPACE_DIMENSION 2
13 #include <Kokkos_Core.hpp>
16 #include <petscversion.h>
18 #include <petscdmplex.h>
19 #if PETSC_VERSION_GE(3,14,0)
21 #include "petscdmswarm.h"
26 #if !defined(PETSC_USE_STRICT_PETSCERRORCODE) && !defined(PETSC_SUCCESS)
27 #define PETSC_SUCCESS 0
41 Kokkos::View<DMObjects*, Kokkos::LayoutRight, HostType>
obj;
46 Kokkos::resize(
obj, 0);
51 for (
int id = 0;
id <
dm_len;
id++){
52 DMDestroy(&
obj(
id).dm);
53 MatDestroy(&
obj(
id).MM);
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;
69 if (
this == &original){
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;
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)
105 #if PETSC_VERSION_LT(3,14,0)
107 printf(
"\n\n\nERROR: initialize_pseudo_inverse_interpolation (should never reach here)\n");
117 PetscReal lo[2], hi[2];
118 cells[0] = num_vPara - 1;
119 cells[1] = num_vPerp - 1;
125 PetscFunctionBeginUser;
127 std::string order_str = std::to_string(order);
128 ierr = PetscOptionsSetValue(NULL,
"-petscspace_degree", order_str.c_str());CHKERRQ(ierr);
133 DMBoundaryType periodicity[2] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
134 PetscBool interpolate = PETSC_TRUE;
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);
148 ierr = DMCreateMassMatrix(
obj(
id).dm,
obj(
id).dm, &
obj(
id).MM);CHKERRQ(ierr);
164 PetscFunctionBeginUser;
166 for (
int id = 0;
id <
dm_len;
id++){
167 if(
obj(
id).dm != NULL){
168 ierr = DMDestroy(&
obj(
id).dm);CHKERRQ(ierr);
171 if(
obj(
id).MM != NULL){
172 ierr = MatDestroy(&
obj(
id).MM);CHKERRQ(ierr);
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