xref: /petsc/src/vec/vec/impls/mpi/pvec2.c (revision f05ece33fd483fcf0ae47f3ecf85fefd993a83b6)
1 
2 /*
3      Code for some of the parallel vector primatives.
4 */
5 #include <../src/vec/vec/impls/mpi/pvecimpl.h>
6 #include <petscblaslapack.h>
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "VecMDot_MPI"
10 PetscErrorCode VecMDot_MPI(Vec xin,PetscInt nv,const Vec y[],PetscScalar *z)
11 {
12   PetscScalar    awork[128],*work = awork;
13   PetscErrorCode ierr;
14 
15   PetscFunctionBegin;
16   if (nv > 128) {
17     ierr = PetscMalloc(nv*sizeof(PetscScalar),&work);CHKERRQ(ierr);
18   }
19   ierr = VecMDot_Seq(xin,nv,y,work);CHKERRQ(ierr);
20   ierr = MPI_Allreduce(work,z,nv,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr);
21   if (nv > 128) {
22     ierr = PetscFree(work);CHKERRQ(ierr);
23   }
24   PetscFunctionReturn(0);
25 }
26 
27 #undef __FUNCT__
28 #define __FUNCT__ "VecMTDot_MPI"
29 PetscErrorCode VecMTDot_MPI(Vec xin,PetscInt nv,const Vec y[],PetscScalar *z)
30 {
31   PetscScalar    awork[128],*work = awork;
32   PetscErrorCode ierr;
33 
34   PetscFunctionBegin;
35   if (nv > 128) {
36     ierr = PetscMalloc(nv*sizeof(PetscScalar),&work);CHKERRQ(ierr);
37   }
38   ierr = VecMTDot_Seq(xin,nv,y,work);CHKERRQ(ierr);
39   ierr = MPI_Allreduce(work,z,nv,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr);
40   if (nv > 128) {
41     ierr = PetscFree(work);CHKERRQ(ierr);
42   }
43   PetscFunctionReturn(0);
44 }
45 
46 #include <../src/vec/vec/impls/seq/ftn-kernels/fnorm.h>
47 #undef __FUNCT__
48 #define __FUNCT__ "VecNorm_MPI"
49 PetscErrorCode VecNorm_MPI(Vec xin,NormType type,PetscReal *z)
50 {
51   PetscReal         sum,work = 0.0;
52   const PetscScalar *xx;
53   PetscErrorCode    ierr;
54   PetscInt          n   = xin->map->n;
55   PetscBLASInt      one = 1,bn;
56 
57   PetscFunctionBegin;
58   ierr = PetscBLASIntCast(n,&bn);CHKERRQ(ierr);
59   if (type == NORM_2 || type == NORM_FROBENIUS) {
60     ierr = VecGetArrayRead(xin,&xx);CHKERRQ(ierr);
61     work = PetscRealPart(BLASdot_(&bn,xx,&one,xx,&one));
62     ierr = VecRestoreArrayRead(xin,&xx);CHKERRQ(ierr);
63     ierr = MPI_Allreduce(&work,&sum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr);
64     *z   = PetscSqrtReal(sum);
65     ierr = PetscLogFlops(2.0*xin->map->n);CHKERRQ(ierr);
66   } else if (type == NORM_1) {
67     /* Find the local part */
68     ierr = VecNorm_Seq(xin,NORM_1,&work);CHKERRQ(ierr);
69     /* Find the global max */
70     ierr = MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr);
71   } else if (type == NORM_INFINITY) {
72     /* Find the local max */
73     ierr = VecNorm_Seq(xin,NORM_INFINITY,&work);CHKERRQ(ierr);
74     /* Find the global max */
75     ierr = MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr);
76   } else if (type == NORM_1_AND_2) {
77     PetscReal temp[2];
78     ierr = VecNorm_Seq(xin,NORM_1,temp);CHKERRQ(ierr);
79     ierr = VecNorm_Seq(xin,NORM_2,temp+1);CHKERRQ(ierr);
80     temp[1] = temp[1]*temp[1];
81     ierr = MPI_Allreduce(temp,z,2,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr);
82     z[1] = PetscSqrtReal(z[1]);
83   }
84   PetscFunctionReturn(0);
85 }
86 
87 /*
88        These two functions are the MPI reduction operation used for max and min with index
89    A call to MPI_Op_create() converts the function Vec[Max,Min]_Local() to the MPI operator Vec[Max,Min]_Local_Op.
90    These are marked PETSC_EXTERN since the function pointers are passed to MPI.
91 */
92 MPI_Op VecMax_Local_Op = 0;
93 MPI_Op VecMin_Local_Op = 0;
94 
95 #undef __FUNCT__
96 #define __FUNCT__ "VecMax_Local"
97 PETSC_EXTERN void MPIAPI VecMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
98 {
99   PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
100 
101   PetscFunctionBegin;
102   if (*datatype != MPIU_REAL) {
103     (*PetscErrorPrintf)("Can only handle MPIU_REAL data types");
104     MPI_Abort(MPI_COMM_SELF,1);
105   }
106   if (xin[0] > xout[0]) {
107     xout[0] = xin[0];
108     xout[1] = xin[1];
109   } else if (xin[0] == xout[0]) {
110     xout[1] = PetscMin(xin[1],xout[1]);
111   }
112   PetscFunctionReturnVoid(); /* cannot return a value */
113 }
114 
115 #undef __FUNCT__
116 #define __FUNCT__ "VecMin_Local"
117 PETSC_EXTERN void MPIAPI VecMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
118 {
119   PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
120 
121   PetscFunctionBegin;
122   if (*datatype != MPIU_REAL) {
123     (*PetscErrorPrintf)("Can only handle MPIU_REAL data types");
124     MPI_Abort(MPI_COMM_SELF,1);
125   }
126   if (xin[0] < xout[0]) {
127     xout[0] = xin[0];
128     xout[1] = xin[1];
129   } else if (xin[0] == xout[0]) {
130     xout[1] = PetscMin(xin[1],xout[1]);
131   }
132   PetscFunctionReturnVoid();
133 }
134 
135 #undef __FUNCT__
136 #define __FUNCT__ "VecMax_MPI"
137 PetscErrorCode VecMax_MPI(Vec xin,PetscInt *idx,PetscReal *z)
138 {
139   PetscErrorCode ierr;
140   PetscReal      work;
141 
142   PetscFunctionBegin;
143   /* Find the local max */
144   ierr = VecMax_Seq(xin,idx,&work);CHKERRQ(ierr);
145 
146   /* Find the global max */
147   if (!idx) {
148     ierr = MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr);
149   } else {
150     PetscReal work2[2],z2[2];
151     PetscInt  rstart;
152     rstart   = xin->map->rstart;
153     work2[0] = work;
154     work2[1] = *idx + rstart;
155     ierr     = MPI_Allreduce(work2,z2,2,MPIU_REAL,VecMax_Local_Op,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr);
156     *z       = z2[0];
157     *idx     = (PetscInt)z2[1];
158   }
159   PetscFunctionReturn(0);
160 }
161 
162 #undef __FUNCT__
163 #define __FUNCT__ "VecMin_MPI"
164 PetscErrorCode VecMin_MPI(Vec xin,PetscInt *idx,PetscReal *z)
165 {
166   PetscErrorCode ierr;
167   PetscReal      work;
168 
169   PetscFunctionBegin;
170   /* Find the local Min */
171   ierr = VecMin_Seq(xin,idx,&work);CHKERRQ(ierr);
172 
173   /* Find the global Min */
174   if (!idx) {
175     ierr = MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr);
176   } else {
177     PetscReal work2[2],z2[2];
178     PetscInt  rstart;
179 
180     ierr = VecGetOwnershipRange(xin,&rstart,NULL);CHKERRQ(ierr);
181     work2[0] = work;
182     work2[1] = *idx + rstart;
183     ierr = MPI_Allreduce(work2,z2,2,MPIU_REAL,VecMin_Local_Op,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr);
184     *z   = z2[0];
185     *idx = (PetscInt)z2[1];
186   }
187   PetscFunctionReturn(0);
188 }
189 
190 
191 
192 
193 
194 
195 
196 
197