Actual source code: pvecimpl.h

petsc-3.10.5 2019-03-28
Report Typos and Errors


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

  7: typedef struct {
  8:   PetscInt insertmode;
  9:   PetscInt count;
 10:   PetscInt bcount;
 11: } VecAssemblyHeader;

 13: typedef struct {
 14:   PetscInt *ints;
 15:   PetscInt *intb;
 16:   PetscScalar *scalars;
 17:   PetscScalar *scalarb;
 18:   char        pendings;
 19:   char        pendingb;
 20: } VecAssemblyFrame;

 22: typedef struct {
 23:   VECHEADER
 24:   PetscInt    nghost;                   /* number of ghost points on this process */
 25:   Vec         localrep;                 /* local representation of vector */
 26:   VecScatter  localupdate;              /* scatter to update ghost values */

 28:   PetscBool   assembly_subset;          /* Subsequent assemblies will set a subset (perhaps equal) of off-process entries set on first assembly */
 29:   PetscBool   use_status;               /* Use MPI_Status to determine number of items in each message */
 30:   PetscMPIInt nsendranks;
 31:   PetscMPIInt nrecvranks;
 32:   PetscMPIInt *sendranks;
 33:   PetscMPIInt *recvranks;
 34:   VecAssemblyHeader *sendhdr,*recvhdr;
 35:   VecAssemblyFrame *sendptrs;   /* pointers to the main messages */
 36:   MPI_Request    *sendreqs;
 37:   MPI_Request    *recvreqs;
 38:   PetscSegBuffer segrecvint;
 39:   PetscSegBuffer segrecvscalar;
 40:   PetscSegBuffer segrecvframe;
 41: } Vec_MPI;

 43: PETSC_INTERN PetscErrorCode VecDot_MPI(Vec,Vec,PetscScalar*);
 44: PETSC_INTERN PetscErrorCode VecMDot_MPI(Vec,PetscInt,const Vec[],PetscScalar*);
 45: PETSC_INTERN PetscErrorCode VecTDot_MPI(Vec,Vec,PetscScalar*);
 46: PETSC_INTERN PetscErrorCode VecMTDot_MPI(Vec,PetscInt,const Vec[],PetscScalar*);
 47: PETSC_INTERN PetscErrorCode VecNorm_MPI(Vec,NormType,PetscReal*);
 48: PETSC_INTERN PetscErrorCode VecMax_MPI(Vec,PetscInt*,PetscReal*);
 49: PETSC_INTERN PetscErrorCode VecMin_MPI(Vec,PetscInt*,PetscReal*);
 50: PETSC_INTERN PetscErrorCode VecDestroy_MPI(Vec);
 51: PETSC_INTERN PetscErrorCode VecView_MPI_Binary(Vec,PetscViewer);
 52: PETSC_INTERN PetscErrorCode VecView_MPI_Draw_LG(Vec,PetscViewer);
 53: PETSC_INTERN PetscErrorCode VecView_MPI_Socket(Vec,PetscViewer);
 54: PETSC_INTERN PetscErrorCode VecView_MPI_HDF5(Vec,PetscViewer);
 55: PETSC_INTERN PetscErrorCode VecView_MPI_ADIOS(Vec,PetscViewer);
 56: PETSC_INTERN PetscErrorCode VecView_MPI_ADIOS2(Vec,PetscViewer);
 57: PETSC_EXTERN PetscErrorCode VecView_MPI(Vec,PetscViewer);
 58: PETSC_INTERN PetscErrorCode VecGetSize_MPI(Vec,PetscInt*);
 59: PETSC_INTERN PetscErrorCode VecGetValues_MPI(Vec,PetscInt,const PetscInt [], PetscScalar []);
 60: PETSC_INTERN PetscErrorCode VecSetValues_MPI(Vec,PetscInt,const PetscInt [],const PetscScalar[],InsertMode);
 61: PETSC_INTERN PetscErrorCode VecSetValuesBlocked_MPI(Vec,PetscInt,const PetscInt [],const PetscScalar[],InsertMode);
 62: PETSC_INTERN PetscErrorCode VecAssemblyBegin_MPI(Vec);
 63: PETSC_INTERN PetscErrorCode VecAssemblyEnd_MPI(Vec);
 64: PETSC_INTERN PetscErrorCode VecAssemblyReset_MPI(Vec);
 65: PETSC_INTERN PetscErrorCode VecCreate_MPI_Private(Vec,PetscBool,PetscInt,const PetscScalar[]);

 67: #endif