Actual source code: pvecimpl.h

petsc-3.6.4 2016-04-12
Report Typos and Errors

  5: #include <../src/vec/vec/impls/dvecimpl.h>

  7: typedef struct {
  8:   VECHEADER
  9:   MPI_Request *send_waits,*recv_waits;  /* for communication during VecAssembly() */
 10:   PetscInt    nsends,nrecvs;
 11:   PetscScalar *svalues,*rvalues;
 12:   PetscInt    rmax;
 13:   PetscInt    nghost;                   /* length of local portion including ghost padding */
 14:   Vec         localrep;                 /* local representation of vector */
 15:   VecScatter  localupdate;              /* scatter to update ghost values */
 16: } Vec_MPI;

 18: PETSC_INTERN PetscErrorCode VecMDot_MPI(Vec,PetscInt,const Vec[],PetscScalar*);
 19: PETSC_INTERN PetscErrorCode VecMTDot_MPI(Vec,PetscInt,const Vec[],PetscScalar*);
 20: PETSC_INTERN PetscErrorCode VecNorm_MPI(Vec,NormType,PetscReal*);
 21: PETSC_INTERN PetscErrorCode VecMax_MPI(Vec,PetscInt*,PetscReal*);
 22: PETSC_INTERN PetscErrorCode VecMin_MPI(Vec,PetscInt*,PetscReal*);
 23: PETSC_INTERN PetscErrorCode VecDestroy_MPI(Vec);
 24: PETSC_INTERN PetscErrorCode VecView_MPI_Binary(Vec,PetscViewer);
 25: PETSC_INTERN PetscErrorCode VecView_MPI_Netcdf(Vec,PetscViewer);
 26: PETSC_INTERN PetscErrorCode VecView_MPI_Draw_LG(Vec,PetscViewer);
 27: PETSC_INTERN PetscErrorCode VecView_MPI_Socket(Vec,PetscViewer);
 28: PETSC_INTERN PetscErrorCode VecView_MPI_HDF5(Vec,PetscViewer);
 29: PETSC_EXTERN PetscErrorCode VecView_MPI(Vec,PetscViewer);
 30: PETSC_INTERN PetscErrorCode VecGetSize_MPI(Vec,PetscInt*);
 31: PETSC_INTERN PetscErrorCode VecGetValues_MPI(Vec,PetscInt,const PetscInt [], PetscScalar []);
 32: PETSC_INTERN PetscErrorCode VecSetValues_MPI(Vec,PetscInt,const PetscInt [],const PetscScalar[],InsertMode);
 33: PETSC_INTERN PetscErrorCode VecSetValuesBlocked_MPI(Vec,PetscInt,const PetscInt [],const PetscScalar[],InsertMode);
 34: PETSC_INTERN PetscErrorCode VecAssemblyBegin_MPI(Vec);
 35: PETSC_INTERN PetscErrorCode VecAssemblyEnd_MPI(Vec);
 36: PETSC_INTERN PetscErrorCode VecCreate_MPI_Private(Vec,PetscBool,PetscInt,const PetscScalar[]);



 40: #endif