15 #ifndef PSEUDO_INVERSE_HPP
16 #define PSEUDO_INVERSE_HPP
17 #include <Kokkos_Core.hpp>
24 #include <petscversion.h>
26 #include <petscdmplex.h>
27 #if PETSC_VERSION_GE(3,14,0)
29 #include "petscdmswarm.h"
45 template<
class Device>
52 Kokkos::View<PseudoInversePetscObjects*, Kokkos::LayoutRight, HostType>
obj;
56 Kokkos::View<double***, Kokkos::LayoutRight, Device>
ptl_coords;
57 Kokkos::View<double***, Kokkos::LayoutRight, Kokkos::HostSpace>
ptl_coords_h;
59 Kokkos::View<double***, Kokkos::LayoutRight, Device>
ptl_weights;
60 Kokkos::View<double***, Kokkos::LayoutRight, Kokkos::HostSpace>
ptl_weights_h;
65 Kokkos::View<int***, Kokkos::LayoutRight, Device>
ptl_indices;
66 Kokkos::View<int***, Kokkos::LayoutRight, Kokkos::HostSpace>
ptl_indices_h;
68 Kokkos::View<int**, Kokkos::LayoutRight, Device>
n_ptl;
69 Kokkos::View<int**, Kokkos::LayoutRight, Kokkos::HostSpace>
n_ptl_h;
71 Kokkos::View<double***, Kokkos::LayoutRight, Device>
grid_coords;
72 Kokkos::View<double***, Kokkos::LayoutRight, Kokkos::HostSpace>
grid_coords_h;
75 Kokkos::View<double***, Kokkos::LayoutRight, Kokkos::HostSpace>
grid_weights_h;
115 Kokkos::resize(
obj, 0);
124 Kokkos::resize(
obj, nthreads);
126 Kokkos::resize(
ptl_coords, n_non_ad, nnode, 2*max_n_ptl_on_node);
127 Kokkos::resize(
ptl_coords_h, n_non_ad, nnode, 2*max_n_ptl_on_node);
128 Kokkos::resize(
ptl_weights, n_non_ad, nnode, max_n_ptl_on_node);
129 Kokkos::resize(
ptl_weights_h, n_non_ad, nnode, max_n_ptl_on_node);
132 Kokkos::resize(
ptl_indices, n_non_ad, nnode, max_n_ptl_on_node);
133 Kokkos::resize(
ptl_indices_h, n_non_ad, nnode, max_n_ptl_on_node);
134 Kokkos::resize(
n_ptl, n_non_ad, nnode);
135 Kokkos::resize(
n_ptl_h, n_non_ad, nnode);
136 Kokkos::resize(
grid_coords, n_non_ad, nnode, 2*n_grid_points);
137 Kokkos::resize(
grid_coords_h, n_non_ad, nnode, 2*n_grid_points);
138 Kokkos::resize(
grid_weights, n_non_ad, nnode, n_grid_points);
150 Kokkos::resize(
obj, 0);
160 Kokkos::resize(
n_ptl, 0, 0);
221 void set_vgrid_filled(Kokkos::View<bool**,Kokkos::LayoutRight,DeviceType> vgrid_filled_in){
235 Kokkos::View<PetscReal*, Kokkos::LayoutRight, Kokkos::HostSpace> particle_coords,
236 Kokkos::View<PetscReal*, Kokkos::LayoutRight, Kokkos::HostSpace> particle_weights,
237 const PetscInt Nparticles,
const PetscInt Nspecies,
238 Kokkos::View<PetscReal*, Kokkos::LayoutRight, Kokkos::HostSpace> grid_coords,
239 Kokkos::View<PetscReal*, Kokkos::LayoutRight, Kokkos::HostSpace> grid_weights);
246 const PetscInt thread_id,
const PetscInt target_id,
247 Kokkos::View<PetscReal*, Kokkos::LayoutRight, Kokkos::HostSpace> particle_coords,
248 Kokkos::View<PetscReal*, Kokkos::LayoutRight, Kokkos::HostSpace> particle_weights,
249 const PetscInt Nparticles,
250 Vec rho, Mat *Mp_out);
254 Kokkos::View<PetscReal*, Kokkos::LayoutRight, Kokkos::HostSpace> grid_weights,
255 const PetscInt Ngrid, Vec rho, Vec rhs);
259 Vec rhs, Mat M_p, Mat PM_p,
260 const PetscInt Ngrid,
261 const PetscInt Nparticles,
262 const PetscInt order,
263 Kokkos::View<PetscReal*, Kokkos::LayoutRight, Kokkos::HostSpace> particle_weights,
265 bool &ksp_converged);
296 const std::string &msg,
297 Kokkos::View<bool**,Kokkos::LayoutRight,DeviceType> vgrid_filled);
300 return (pseudo_inv.
n_ptl(isp_non_ad, node) > 0 && (pseudo_inv.
vgrid_filled(isp_non_ad, node) || pseudo_inv.
force));
void send_host_to_device()
Definition: pseudo_inverse.hpp:194
Kokkos::View< double ***, Kokkos::LayoutRight, Device > grid_weights
Stores velocity grid weights: weight(non-adiabatic species, node, grid point index) ...
Definition: pseudo_inverse.hpp:74
int max_n_ptl_on_node_and_vgrid_filled(const Simulation< DeviceType > &sml, const Grid< DeviceType > &grid, const MagneticField< DeviceType > &magnetic_field, Plasma &plasma, const VelocityGrid &vgrid, const DomainDecomposition< DeviceType > &pol_decomp, const std::string &msg, Kokkos::View< bool **, Kokkos::LayoutRight, DeviceType > vgrid_filled)
Definition: pseudo_inverse.cpp:967
Kokkos::View< double ***, Kokkos::LayoutRight, Kokkos::HostSpace > ptl_weights_h
host version of ptl_weights
Definition: pseudo_inverse.hpp:60
For description see source file.
Definition: pseudo_inverse.hpp:46
Kokkos::View< int **, Kokkos::LayoutRight, Device > n_ptl
Stores number of particles (non-adiabatic species, node)
Definition: pseudo_inverse.hpp:68
int n_nonadiabatic_species
Number of nonadiabatic species.
Definition: pseudo_inverse.hpp:87
PetscErrorCode pseudo_inv_create_projection(const DM dm, DM sw, const PetscInt thread_id, const PetscInt target_id, Kokkos::View< PetscReal *, Kokkos::LayoutRight, Kokkos::HostSpace > particle_coords, Kokkos::View< PetscReal *, Kokkos::LayoutRight, Kokkos::HostSpace > particle_weights, const PetscInt Nparticles, Vec rho, Mat *Mp_out)
For description see source file.
Definition: pseudo_inverse.cpp:344
Definition: velocity_grid.hpp:8
subroutine plasma(grid, itr, p, dene_out, deni_out, Te_out, Ti_out, Vparai_out)
Calculate the plasma density, temperature, and parallel velocity for a point in triangle itr using pl...
Definition: neutral_totalf.F90:1224
Kokkos::View< int ***, Kokkos::LayoutRight, Kokkos::HostSpace > ptl_indices_h
host version of ptl_indices
Definition: pseudo_inverse.hpp:66
bool use_delta_grid
Whether to only map delta-grid weights.
Definition: pseudo_inverse.hpp:86
Definition: DM_wrapper.hpp:38
Kokkos::View< int **, Kokkos::LayoutRight, Kokkos::HostSpace > n_ptl_h
host version of n_ptl
Definition: pseudo_inverse.hpp:69
Definition: magnetic_field.hpp:12
void pseudo_inv_allocate(const Simulation< DeviceType > &sml, const Grid< DeviceType > &grid, const MagneticField< DeviceType > &magnetic_field, Plasma &plasma, const VelocityGrid &vgrid, const DomainDecomposition< DeviceType > &pol_decomp, Pseudo_inverse< DeviceType > &pseudo_inv)
< For description see source file
Definition: pseudo_inverse.cpp:719
Vec rho
PETSc helper vector.
Definition: pseudo_inverse.hpp:39
void resize_to_zero()
Definition: pseudo_inverse.hpp:148
Mat M_p
PETSc projection matrix.
Definition: pseudo_inverse.hpp:41
Kokkos::View< double ***, Kokkos::LayoutRight, Device > grid_alpha_weights
Stores velocity grid alpha (see sml_f0_grid_alpha) weights: weight(non-adiabatic species, node, grid point index)
Definition: pseudo_inverse.hpp:77
Kokkos::View< bool **, Kokkos::LayoutRight, Device > vgrid_filled
Stores whether a node has no empty velocity cells (non-adiabatic species, node)
Definition: pseudo_inverse.hpp:80
Kokkos::View< double ***, Kokkos::LayoutRight, Device > grid_coords
Stores velocity grid coordinates: v_para(non-adiabatic species, node, 2*(grid point index)); v_perp(n...
Definition: pseudo_inverse.hpp:71
Kokkos::View< double ***, Kokkos::LayoutRight, Device > ptl_alpha_weights
Stores particle alpha (see sml_f0_grid_alpha) weights: weight(non-adiabatic species, node, particle index on node)
Definition: pseudo_inverse.hpp:62
void resize_on_device_and_host(int nthreads, int n_non_ad, int nnode, int max_n_ptl_on_node, int n_grid_points)
Definition: pseudo_inverse.hpp:122
Kokkos::View< double ***, Kokkos::LayoutRight, Device > ptl_weights
Stores particle weights: weight(non-adiabatic species, node, particle index on node) ...
Definition: pseudo_inverse.hpp:59
Kokkos::View< int ***, Kokkos::LayoutRight, Device > ptl_indices
Stores particle indices: index(non-adiabatic species, node, particle index on node) ...
Definition: pseudo_inverse.hpp:65
bool force
Whether to force pseudo-inverse even if a node has empty velocity cells.
Definition: pseudo_inverse.hpp:84
bool do_particles_to_grid_petsc
Whether the particles -> velocity grid interpolation should be done without or with (probably slower)...
Definition: pseudo_inverse.hpp:85
void all_species_update_pseudo_inv_weights(const Simulation< DeviceType > &sml, const Grid< DeviceType > &grid, const MagneticField< DeviceType > &magnetic_field, Plasma &plasma, const VelocityGrid &vgrid, const DomainDecomposition< DeviceType > &pol_decomp, DMWrapper &pseudo_inv_dm, Pseudo_inverse< DeviceType > &pseudo_inv)
For description see source file.
Definition: pseudo_inverse.cpp:765
void resize_alpha_weights_to_zero()
Definition: pseudo_inverse.hpp:174
PetscInt pseudo_inverse_interpolate_particles_to_velocity_grid(const DM &dm, Kokkos::View< PetscReal *, Kokkos::LayoutRight, Kokkos::HostSpace > particle_coords, Kokkos::View< PetscReal *, Kokkos::LayoutRight, Kokkos::HostSpace > particle_weights, const PetscInt Nparticles, const PetscInt Nspecies, Kokkos::View< PetscReal *, Kokkos::LayoutRight, Kokkos::HostSpace > grid_coords, Kokkos::View< PetscReal *, Kokkos::LayoutRight, Kokkos::HostSpace > grid_weights)
For description see source file.
Definition: pseudo_inverse.cpp:161
PetscErrorCode pseudo_inv_set_rhs(const DM dm, Mat MM, Kokkos::View< PetscReal *, Kokkos::LayoutRight, Kokkos::HostSpace > grid_weights, const PetscInt Ngrid, Vec rho, Vec rhs)
For description see source file.
Definition: pseudo_inverse.cpp:432
Kokkos::View< double ***, Kokkos::LayoutRight, Kokkos::HostSpace > ptl_coords_h
host version of ptl_coords
Definition: pseudo_inverse.hpp:57
DM swarm
PETSc particle swarm object.
Definition: pseudo_inverse.hpp:38
Kokkos::View< double ***, Kokkos::LayoutRight, Kokkos::HostSpace > grid_alpha_weights_h
host version of grid_alpha_weights
Definition: pseudo_inverse.hpp:78
Pseudo_inverse()
Definition: pseudo_inverse.hpp:90
Kokkos::View< bool **, Kokkos::LayoutRight, Kokkos::HostSpace > vgrid_filled_h
host version of vgrid_filled
Definition: pseudo_inverse.hpp:81
void send_grid_host_to_device()
Definition: pseudo_inverse.hpp:215
Definition: pseudo_inverse.hpp:37
PetscErrorCode pseudo_inv_grid_to_particles(const DM dm, DM sw, Vec rhs, Mat M_p, Mat PM_p, const PetscInt Ngrid, const PetscInt Nparticles, const PetscInt order, Kokkos::View< PetscReal *, Kokkos::LayoutRight, Kokkos::HostSpace > particle_weights, int &ksp_num_its, bool &ksp_converged)
Definition: pseudo_inverse.cpp:505
void zero_out_grid_weights()
Definition: pseudo_inverse.hpp:226
Kokkos::View< double ***, Kokkos::LayoutRight, Device > ptl_coords
Stores particle velocity coordinates: v_para(non-adiabatic species, node, 2*(particle index on node))...
Definition: pseudo_inverse.hpp:56
Definition: magnetic_field.F90:1
Definition: plasma.hpp:13
Vec rhs
PETSc right-hand side vector.
Definition: pseudo_inverse.hpp:40
int expected_n_grid_points
Expected number of grid points returned from PETSc.
Definition: pseudo_inverse.hpp:83
Kokkos::View< double ***, Kokkos::LayoutRight, Kokkos::HostSpace > grid_weights_h
host version of grid_weights
Definition: pseudo_inverse.hpp:75
KOKKOS_INLINE_FUNCTION bool pseudo_inv_do_node(const Pseudo_inverse< Device > &pseudo_inv, int isp_non_ad, int node)
Definition: pseudo_inverse.hpp:299
void set_vgrid_filled(Kokkos::View< bool **, Kokkos::LayoutRight, DeviceType > vgrid_filled_in)
Definition: pseudo_inverse.hpp:221
void send_device_to_host()
Definition: pseudo_inverse.hpp:182
bool pseudo_inv_do_node_h(Pseudo_inverse< DeviceType > &pseudo_inv, int isp_non_ad, int node)
Definition: pseudo_inverse.hpp:303
void send_ptl_device_to_host()
Definition: pseudo_inverse.hpp:206
Kokkos::View< PseudoInversePetscObjects *, Kokkos::LayoutRight, HostType > obj
Stores number of threads PseudoInversePetscObjects.
Definition: pseudo_inverse.hpp:52
PetscErrorCode pseudo_inv_create_Swarm(const DM dm, DM *sw)
For description see source file.
Definition: pseudo_inverse.cpp:301
Kokkos::View< double ***, Kokkos::LayoutRight, Kokkos::HostSpace > ptl_alpha_weights_h
host version of ptl_alpha_weights
Definition: pseudo_inverse.hpp:63
Kokkos::View< double ***, Kokkos::LayoutRight, Kokkos::HostSpace > grid_coords_h
host version of grid_coords
Definition: pseudo_inverse.hpp:72