Actual source code: pvecimpl.h
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 first_assembly_done; /* Is the first time assembly done? */
30: PetscBool use_status; /* Use MPI_Status to determine number of items in each message */
31: PetscMPIInt nsendranks;
32: PetscMPIInt nrecvranks;
33: PetscMPIInt *sendranks;
34: PetscMPIInt *recvranks;
35: VecAssemblyHeader *sendhdr,*recvhdr;
36: VecAssemblyFrame *sendptrs; /* pointers to the main messages */
37: MPI_Request *sendreqs;
38: MPI_Request *recvreqs;
39: PetscSegBuffer segrecvint;
40: PetscSegBuffer segrecvscalar;
41: PetscSegBuffer segrecvframe;
42: #if defined(PETSC_HAVE_NVSHMEM)
43: PetscBool use_nvshmem; /* Try to use NVSHMEM in communication of, for example, VecNorm */
44: #endif
45: } Vec_MPI;
47: PETSC_INTERN PetscErrorCode VecDot_MPI(Vec,Vec,PetscScalar*);
48: PETSC_INTERN PetscErrorCode VecMDot_MPI(Vec,PetscInt,const Vec[],PetscScalar*);
49: PETSC_INTERN PetscErrorCode VecTDot_MPI(Vec,Vec,PetscScalar*);
50: PETSC_INTERN PetscErrorCode VecMTDot_MPI(Vec,PetscInt,const Vec[],PetscScalar*);
51: PETSC_INTERN PetscErrorCode VecNorm_MPI(Vec,NormType,PetscReal*);
52: PETSC_INTERN PetscErrorCode VecMax_MPI(Vec,PetscInt*,PetscReal*);
53: PETSC_INTERN PetscErrorCode VecMin_MPI(Vec,PetscInt*,PetscReal*);
54: PETSC_INTERN PetscErrorCode VecDestroy_MPI(Vec);
55: PETSC_INTERN PetscErrorCode VecView_MPI_Binary(Vec,PetscViewer);
56: PETSC_INTERN PetscErrorCode VecView_MPI_Draw_LG(Vec,PetscViewer);
57: PETSC_INTERN PetscErrorCode VecView_MPI_Socket(Vec,PetscViewer);
58: PETSC_INTERN PetscErrorCode VecView_MPI_HDF5(Vec,PetscViewer);
59: PETSC_INTERN PetscErrorCode VecView_MPI_ADIOS(Vec,PetscViewer);
60: PETSC_EXTERN PetscErrorCode VecView_MPI(Vec,PetscViewer);
61: PETSC_INTERN PetscErrorCode VecGetSize_MPI(Vec,PetscInt*);
62: PETSC_INTERN PetscErrorCode VecGetValues_MPI(Vec,PetscInt,const PetscInt [], PetscScalar []);
63: PETSC_INTERN PetscErrorCode VecSetValues_MPI(Vec,PetscInt,const PetscInt [],const PetscScalar[],InsertMode);
64: PETSC_INTERN PetscErrorCode VecSetValuesBlocked_MPI(Vec,PetscInt,const PetscInt [],const PetscScalar[],InsertMode);
65: PETSC_INTERN PetscErrorCode VecAssemblyBegin_MPI(Vec);
66: PETSC_INTERN PetscErrorCode VecAssemblyEnd_MPI(Vec);
67: PETSC_INTERN PetscErrorCode VecAssemblyReset_MPI(Vec);
68: PETSC_INTERN PetscErrorCode VecCreate_MPI_Private(Vec,PetscBool,PetscInt,const PetscScalar[]);
69: PETSC_EXTERN PetscErrorCode VecCreate_MPI(Vec);
70: PETSC_INTERN PetscErrorCode VecDuplicate_MPI(Vec,Vec*);
72: #endif