Actual source code: petscviennacl.h

  1: #pragma once

  3: #include <petscvec.h>

  5: #if defined(PETSC_HAVE_CUDA)
  6:   #define VIENNACL_WITH_CUDA
  7: #endif

  9: #if defined(PETSC_HAVE_OPENCL)
 10:   #define VIENNACL_WITH_OPENCL
 11: #endif

 13: #if defined(PETSC_HAVE_OPENMP)
 14:   #define VIENNACL_WITH_OPENMP
 15: #endif

 17: #include <viennacl/forwards.h>
 18: #include <viennacl/vector_proxy.hpp>
 19: #include <viennacl/vector.hpp>

 21: PETSC_EXTERN PetscErrorCode VecViennaCLGetArray(Vec v, viennacl::vector<PetscScalar> **a);
 22: PETSC_EXTERN PetscErrorCode VecViennaCLRestoreArray(Vec v, viennacl::vector<PetscScalar> **a);

 24: PETSC_EXTERN PetscErrorCode VecViennaCLGetArrayRead(Vec v, const viennacl::vector<PetscScalar> **a);
 25: PETSC_EXTERN PetscErrorCode VecViennaCLRestoreArrayRead(Vec v, const viennacl::vector<PetscScalar> **a);

 27: PETSC_EXTERN PetscErrorCode VecViennaCLGetArrayWrite(Vec v, viennacl::vector<PetscScalar> **a);
 28: PETSC_EXTERN PetscErrorCode VecViennaCLRestoreArrayWrite(Vec v, viennacl::vector<PetscScalar> **a);

 30: PETSC_EXTERN PetscErrorCode VecCreateSeqViennaCLWithArray(MPI_Comm, PetscInt, PetscInt, const viennacl::vector<PetscScalar> *, Vec *);
 31: PETSC_EXTERN PetscErrorCode VecCreateSeqViennaCLWithArrays(MPI_Comm, PetscInt, PetscInt, const PetscScalar *, const viennacl::vector<PetscScalar> *, Vec *);
 32: PETSC_EXTERN PetscErrorCode VecCreateMPIViennaCLWithArray(MPI_Comm, PetscInt, PetscInt, PetscInt, const viennacl::vector<PetscScalar> *, Vec *);
 33: PETSC_EXTERN PetscErrorCode VecCreateMPIViennaCLWithArrays(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscScalar *, const viennacl::vector<PetscScalar> *, Vec *);

 35: PETSC_EXTERN PetscErrorCode VecViennaCLPlaceArray(Vec, const viennacl::vector<PetscScalar> *);
 36: PETSC_EXTERN PetscErrorCode VecViennaCLResetArray(Vec);