xref: /petsc/src/vec/vec/impls/mpi/pvecimpl.h (revision bbfbdd31d4e286473f4179abdcd08de394213aa7)
1 #pragma once
2 
3 #include <../src/vec/vec/impls/dvecimpl.h>
4 
5 typedef struct {
6   PetscInt insertmode;
7   PetscInt count;
8   PetscInt bcount;
9 } VecAssemblyHeader;
10 
11 typedef struct {
12   PetscInt    *ints;
13   PetscInt    *intb;
14   PetscScalar *scalars;
15   PetscScalar *scalarb;
16   char         pendings;
17   char         pendingb;
18 } VecAssemblyFrame;
19 
20 typedef struct {
21   VECHEADER
22   PetscInt   nghost;      /* number of ghost points on this process */
23   IS         ghost;       /* global indices of ghost values */
24   Vec        localrep;    /* local representation of vector */
25   VecScatter localupdate; /* scatter to update ghost values */
26 
27   PetscBool          assembly_subset;     /* Subsequent assemblies will set a subset (perhaps equal) of off-process entries set on first assembly */
28   PetscBool          first_assembly_done; /* Is the first time assembly done? */
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 #if defined(PETSC_HAVE_NVSHMEM)
42   PetscBool use_nvshmem; /* Try to use NVSHMEM in communication of, for example, VecNorm */
43 #endif
44 
45   /* COO fields, assuming m is the vector's local size */
46   PetscCount  coo_n;
47   PetscCount  tot1;  /* Total local entries in COO arrays */
48   PetscCount *jmap1; /* [m+1]: i-th entry of the vector has jmap1[i+1]-jmap1[i] repeats in COO arrays */
49   PetscCount *perm1; /* [tot1]: permutation array for local entries */
50 
51   PetscCount  nnz2;  /* Unique entries in recvbuf */
52   PetscCount *imap2; /* [nnz2]: i-th unique entry in recvbuf is imap2[i]-th entry in the vector */
53   PetscCount *jmap2; /* [nnz2+1] */
54   PetscCount *perm2; /* [recvlen] */
55 
56   PetscSF      coo_sf;
57   PetscCount   sendlen, recvlen;  /* Lengths (in unit of PetscScalar) of send/recvbuf */
58   PetscCount  *Cperm;             /* [sendlen]: permutation array to fill sendbuf[]. 'C' for communication */
59   PetscScalar *sendbuf, *recvbuf; /* Buffers for remote values in VecSetValuesCOO() */
60 } Vec_MPI;
61 
62 PETSC_INTERN PetscErrorCode VecMTDot_MPI(Vec, PetscInt, const Vec[], PetscScalar *);
63 PETSC_INTERN PetscErrorCode VecView_MPI_Binary(Vec, PetscViewer);
64 PETSC_INTERN PetscErrorCode VecView_MPI_Draw_LG(Vec, PetscViewer);
65 PETSC_INTERN PetscErrorCode VecView_MPI_Socket(Vec, PetscViewer);
66 PETSC_INTERN PetscErrorCode VecView_MPI_HDF5(Vec, PetscViewer);
67 PETSC_INTERN PetscErrorCode VecView_MPI_ADIOS(Vec, PetscViewer);
68 PETSC_INTERN PetscErrorCode VecGetSize_MPI(Vec, PetscInt *);
69 PETSC_INTERN PetscErrorCode VecGetValues_MPI(Vec, PetscInt, const PetscInt[], PetscScalar[]);
70 PETSC_INTERN PetscErrorCode VecSetValues_MPI(Vec, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
71 PETSC_INTERN PetscErrorCode VecSetValuesBlocked_MPI(Vec, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
72 PETSC_INTERN PetscErrorCode VecAssemblyBegin_MPI(Vec);
73 PETSC_INTERN PetscErrorCode VecAssemblyEnd_MPI(Vec);
74 PETSC_INTERN PetscErrorCode VecAssemblyReset_MPI(Vec);
75 PETSC_EXTERN PetscErrorCode VecCreate_MPI(Vec);
76 PETSC_INTERN PetscErrorCode VecMDot_MPI_GEMV(Vec, PetscInt, const Vec[], PetscScalar *);
77 PETSC_INTERN PetscErrorCode VecMTDot_MPI_GEMV(Vec, PetscInt, const Vec[], PetscScalar *);
78 
79 PETSC_INTERN PetscErrorCode VecDuplicate_MPI(Vec, Vec *);
80 PETSC_INTERN PetscErrorCode VecDuplicateWithArray_MPI(Vec, const PetscScalar *, Vec *);
81 PETSC_INTERN PetscErrorCode VecSetPreallocationCOO_MPI(Vec, PetscCount, const PetscInt[]);
82 PETSC_INTERN PetscErrorCode VecSetValuesCOO_MPI(Vec, const PetscScalar[], InsertMode);
83 
84 PETSC_INTERN PetscErrorCode VecDot_MPI(Vec, Vec, PetscScalar *);
85 PETSC_INTERN PetscErrorCode VecMDot_MPI(Vec, PetscInt, const Vec[], PetscScalar *);
86 PETSC_INTERN PetscErrorCode VecTDot_MPI(Vec, Vec, PetscScalar *);
87 PETSC_INTERN PetscErrorCode VecNorm_MPI(Vec, NormType, PetscReal *);
88 PETSC_INTERN PetscErrorCode VecMax_MPI(Vec, PetscInt *, PetscReal *);
89 PETSC_INTERN PetscErrorCode VecMin_MPI(Vec, PetscInt *, PetscReal *);
90 PETSC_INTERN PetscErrorCode VecMaxPointwiseDivide_MPI(Vec, Vec, PetscReal *);
91 PETSC_INTERN PetscErrorCode VecPlaceArray_MPI(Vec, const PetscScalar *);
92 PETSC_INTERN PetscErrorCode VecResetArray_MPI(Vec);
93 PETSC_INTERN PetscErrorCode VecCreate_MPI_Private(Vec, PetscBool, PetscInt, const PetscScalar[]);
94 
95 PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);
96 PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecDestroy_MPI(Vec);
97 
VecMXDot_MPI_Default(Vec xin,PetscInt nv,const Vec y[],PetscScalar * z,PetscErrorCode (* VecMXDot_SeqFn)(Vec,PetscInt,const Vec[],PetscScalar *))98 static inline PetscErrorCode VecMXDot_MPI_Default(Vec xin, PetscInt nv, const Vec y[], PetscScalar *z, PetscErrorCode (*VecMXDot_SeqFn)(Vec, PetscInt, const Vec[], PetscScalar *))
99 {
100   PetscFunctionBegin;
101   PetscCall(VecMXDot_SeqFn(xin, nv, y, z));
102   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, z, nv, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)xin)));
103   PetscFunctionReturn(PETSC_SUCCESS);
104 }
105 
VecXDot_MPI_Default(Vec xin,Vec yin,PetscScalar * z,PetscErrorCode (* VecXDot_SeqFn)(Vec,Vec,PetscScalar *))106 static inline PetscErrorCode VecXDot_MPI_Default(Vec xin, Vec yin, PetscScalar *z, PetscErrorCode (*VecXDot_SeqFn)(Vec, Vec, PetscScalar *))
107 {
108   PetscFunctionBegin;
109   PetscCall(VecXDot_SeqFn(xin, yin, z));
110   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, z, 1, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)xin)));
111   PetscFunctionReturn(PETSC_SUCCESS);
112 }
113 
VecMinMax_MPI_Default(Vec xin,PetscInt * idx,PetscReal * z,PetscErrorCode (* VecMinMax_SeqFn)(Vec,PetscInt *,PetscReal *),const MPI_Op ops[2])114 static inline PetscErrorCode VecMinMax_MPI_Default(Vec xin, PetscInt *idx, PetscReal *z, PetscErrorCode (*VecMinMax_SeqFn)(Vec, PetscInt *, PetscReal *), const MPI_Op ops[2])
115 {
116   PetscFunctionBegin;
117   /* Find the local max */
118   PetscCall(VecMinMax_SeqFn(xin, idx, z));
119   if (PetscDefined(HAVE_MPIUNI)) PetscFunctionReturn(PETSC_SUCCESS);
120   /* Find the global max */
121   if (idx) {
122     struct {
123       PetscReal v;
124       PetscInt  i;
125     } in = {*z, *idx + xin->map->rstart};
126 
127     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &in, 1, MPIU_REAL_INT, ops[0], PetscObjectComm((PetscObject)xin)));
128     *z   = in.v;
129     *idx = in.i;
130   } else {
131     /* User does not need idx */
132     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, z, 1, MPIU_REAL, ops[1], PetscObjectComm((PetscObject)xin)));
133   }
134   PetscFunctionReturn(PETSC_SUCCESS);
135 }
136 
VecDotNorm2_MPI_Default(Vec s,Vec t,PetscScalar * dp,PetscScalar * nm,PetscErrorCode (* VecDotNorm2_SeqFn)(Vec,Vec,PetscScalar *,PetscScalar *))137 static inline PetscErrorCode VecDotNorm2_MPI_Default(Vec s, Vec t, PetscScalar *dp, PetscScalar *nm, PetscErrorCode (*VecDotNorm2_SeqFn)(Vec, Vec, PetscScalar *, PetscScalar *))
138 {
139   PetscFunctionBegin;
140   PetscCall(VecDotNorm2_SeqFn(s, t, dp, nm));
141   {
142     PetscScalar sum[] = {*dp, *nm};
143 
144     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &sum, 2, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)s)));
145     *dp = sum[0];
146     *nm = sum[1];
147   }
148   PetscFunctionReturn(PETSC_SUCCESS);
149 }
150 
VecNorm_MPI_Default(Vec xin,NormType type,PetscReal * z,PetscErrorCode (* VecNorm_SeqFn)(Vec,NormType,PetscReal *))151 static inline PetscErrorCode VecNorm_MPI_Default(Vec xin, NormType type, PetscReal *z, PetscErrorCode (*VecNorm_SeqFn)(Vec, NormType, PetscReal *))
152 {
153   PetscMPIInt zn = 1;
154   MPI_Op      op = MPIU_SUM;
155 
156   PetscFunctionBegin;
157   PetscCall(VecNorm_SeqFn(xin, type, z));
158   switch (type) {
159   case NORM_1_AND_2:
160     // the 2 norm needs to be squared below before being summed. NORM_2 stores the norm in the
161     // first slot but while NORM_1_AND_2 stores it in the second
162     z[1] *= z[1];
163     zn = 2;
164     break;
165   case NORM_2:
166   case NORM_FROBENIUS:
167     z[0] *= z[0];
168   case NORM_1:
169     break;
170   case NORM_INFINITY:
171     op = MPIU_MAX;
172     break;
173   }
174   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, z, zn, MPIU_REAL, op, PetscObjectComm((PetscObject)xin)));
175   if (type == NORM_2 || type == NORM_FROBENIUS || type == NORM_1_AND_2) z[type == NORM_1_AND_2] = PetscSqrtReal(z[type == NORM_1_AND_2]);
176   PetscFunctionReturn(PETSC_SUCCESS);
177 }
178 
VecErrorWeightedNorms_MPI_Default(Vec U,Vec Y,Vec E,NormType wnormtype,PetscReal atol,Vec vatol,PetscReal rtol,Vec vrtol,PetscReal ignore_max,PetscReal * norm,PetscInt * norm_loc,PetscReal * norma,PetscInt * norma_loc,PetscReal * normr,PetscInt * normr_loc,PetscErrorCode (* SeqFn)(Vec,Vec,Vec,NormType,PetscReal,Vec,PetscReal,Vec,PetscReal,PetscReal *,PetscInt *,PetscReal *,PetscInt *,PetscReal *,PetscInt *))179 static inline PetscErrorCode VecErrorWeightedNorms_MPI_Default(Vec U, Vec Y, Vec E, NormType wnormtype, PetscReal atol, Vec vatol, PetscReal rtol, Vec vrtol, PetscReal ignore_max, PetscReal *norm, PetscInt *norm_loc, PetscReal *norma, PetscInt *norma_loc, PetscReal *normr, PetscInt *normr_loc, PetscErrorCode (*SeqFn)(Vec, Vec, Vec, NormType, PetscReal, Vec, PetscReal, Vec, PetscReal, PetscReal *, PetscInt *, PetscReal *, PetscInt *, PetscReal *, PetscInt *))
180 {
181   PetscReal loc[6];
182 
183   PetscFunctionBegin;
184   PetscCall(SeqFn(U, Y, E, wnormtype, atol, vatol, rtol, vrtol, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc));
185   if (wnormtype == NORM_2) {
186     loc[0] = PetscSqr(*norm);
187     loc[1] = PetscSqr(*norma);
188     loc[2] = PetscSqr(*normr);
189   } else {
190     loc[0] = *norm;
191     loc[1] = *norma;
192     loc[2] = *normr;
193   }
194   loc[3] = (PetscReal)*norm_loc;
195   loc[4] = (PetscReal)*norma_loc;
196   loc[5] = (PetscReal)*normr_loc;
197   if (wnormtype == NORM_2) {
198     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, loc, 6, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)U)));
199   } else {
200     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, loc, 3, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)U)));
201     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, loc + 3, 3, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)U)));
202   }
203   if (wnormtype == NORM_2) {
204     *norm  = PetscSqrtReal(loc[0]);
205     *norma = PetscSqrtReal(loc[1]);
206     *normr = PetscSqrtReal(loc[2]);
207   } else {
208     *norm  = loc[0];
209     *norma = loc[1];
210     *normr = loc[2];
211   }
212   *norm_loc  = loc[3];
213   *norma_loc = loc[4];
214   *normr_loc = loc[5];
215   PetscFunctionReturn(PETSC_SUCCESS);
216 }
217