xref: /petsc/src/vec/vec/impls/mpi/pvec2.c (revision a6abbb56490c78acb90244568e53007c89095626)
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,((PetscObject)xin)->comm);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,((PetscObject)xin)->comm);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   PetscScalar    *xx = *(PetscScalar**)xin->data;
53   PetscErrorCode ierr;
54   PetscInt       n = xin->map->n;
55   PetscBLASInt   one = 1,bn = PetscBLASIntCast(n);
56 
57   PetscFunctionBegin;
58   if (type == NORM_2 || type == NORM_FROBENIUS) {
59     work  = BLASdot_(&bn,xx,&one,xx,&one);
60     ierr = MPI_Allreduce(&work,&sum,1,MPIU_REAL,MPIU_SUM,((PetscObject)xin)->comm);CHKERRQ(ierr);
61     *z = PetscSqrtReal(sum);
62     ierr = PetscLogFlops(2.0*xin->map->n);CHKERRQ(ierr);
63   } else if (type == NORM_1) {
64     /* Find the local part */
65     ierr = VecNorm_Seq(xin,NORM_1,&work);CHKERRQ(ierr);
66     /* Find the global max */
67     ierr = MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_SUM,((PetscObject)xin)->comm);CHKERRQ(ierr);
68   } else if (type == NORM_INFINITY) {
69     /* Find the local max */
70     ierr = VecNorm_Seq(xin,NORM_INFINITY,&work);CHKERRQ(ierr);
71     /* Find the global max */
72     ierr = MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_MAX,((PetscObject)xin)->comm);CHKERRQ(ierr);
73   } else if (type == NORM_1_AND_2) {
74     PetscReal temp[2];
75     ierr = VecNorm_Seq(xin,NORM_1,temp);CHKERRQ(ierr);
76     ierr = VecNorm_Seq(xin,NORM_2,temp+1);CHKERRQ(ierr);
77     temp[1] = temp[1]*temp[1];
78     ierr = MPI_Allreduce(temp,z,2,MPIU_REAL,MPIU_SUM,((PetscObject)xin)->comm);CHKERRQ(ierr);
79     z[1] = PetscSqrtReal(z[1]);
80   }
81   PetscFunctionReturn(0);
82 }
83 
84 /*
85        These two functions are the MPI reduction operation used for max and min with index
86    The call below to MPI_Op_create() converts the function Vec[Max,Min]_Local() to the
87    MPI operator Vec[Max,Min]_Local_Op.
88 */
89 MPI_Op VecMax_Local_Op = 0;
90 MPI_Op VecMin_Local_Op = 0;
91 
92 EXTERN_C_BEGIN
93 #undef __FUNCT__
94 #define __FUNCT__ "VecMax_Local"
95 void  MPIAPI VecMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
96 {
97   PetscReal *xin = (PetscReal *)in,*xout = (PetscReal*)out;
98 
99   PetscFunctionBegin;
100   if (*datatype != MPIU_REAL) {
101     (*PetscErrorPrintf)("Can only handle MPIU_REAL data types");
102     MPI_Abort(MPI_COMM_WORLD,1);
103   }
104   if (xin[0] > xout[0]) {
105     xout[0] = xin[0];
106     xout[1] = xin[1];
107   } else if (xin[0] == xout[0]) {
108     xout[1] = PetscMin(xin[1],xout[1]);
109   }
110   PetscFunctionReturnVoid(); /* cannot return a value */
111 }
112 EXTERN_C_END
113 
114 EXTERN_C_BEGIN
115 #undef __FUNCT__
116 #define __FUNCT__ "VecMin_Local"
117 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_WORLD,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 EXTERN_C_END
135 
136 #undef __FUNCT__
137 #define __FUNCT__ "VecMax_MPI"
138 PetscErrorCode VecMax_MPI(Vec xin,PetscInt *idx,PetscReal *z)
139 {
140   PetscErrorCode ierr;
141   PetscReal      work;
142 
143   PetscFunctionBegin;
144   /* Find the local max */
145   ierr = VecMax_Seq(xin,idx,&work);CHKERRQ(ierr);
146 
147   /* Find the global max */
148   if (!idx) {
149     ierr = MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_MAX,((PetscObject)xin)->comm);CHKERRQ(ierr);
150   } else {
151     PetscReal work2[2],z2[2];
152     PetscInt  rstart;
153     rstart = xin->map->rstart;
154     work2[0] = work;
155     work2[1] = *idx + rstart;
156     ierr = MPI_Allreduce(work2,z2,2,MPIU_REAL,VecMax_Local_Op,((PetscObject)xin)->comm);CHKERRQ(ierr);
157     *z   = z2[0];
158     *idx = (PetscInt)z2[1];
159   }
160   PetscFunctionReturn(0);
161 }
162 
163 #undef __FUNCT__
164 #define __FUNCT__ "VecMin_MPI"
165 PetscErrorCode VecMin_MPI(Vec xin,PetscInt *idx,PetscReal *z)
166 {
167   PetscErrorCode ierr;
168   PetscReal      work;
169 
170   PetscFunctionBegin;
171   /* Find the local Min */
172   ierr = VecMin_Seq(xin,idx,&work);CHKERRQ(ierr);
173 
174   /* Find the global Min */
175   if (!idx) {
176     ierr = MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_MIN,((PetscObject)xin)->comm);CHKERRQ(ierr);
177   } else {
178     PetscReal work2[2],z2[2];
179     PetscInt  rstart;
180 
181     ierr = VecGetOwnershipRange(xin,&rstart,PETSC_NULL);CHKERRQ(ierr);
182     work2[0] = work;
183     work2[1] = *idx + rstart;
184     ierr = MPI_Allreduce(work2,z2,2,MPIU_REAL,VecMin_Local_Op,((PetscObject)xin)->comm);CHKERRQ(ierr);
185     *z   = z2[0];
186     *idx = (PetscInt)z2[1];
187   }
188   PetscFunctionReturn(0);
189 }
190 
191 
192 
193 
194 
195 
196 
197 
198