Actual source code: pvec2.c
2: /*
3: Code for some of the parallel vector primatives.
4: */
5: #include <../src/vec/vec/impls/mpi/pvecimpl.h>
6: #include <petscblaslapack.h>
10: PetscErrorCode VecMDot_MPI(Vec xin,PetscInt nv,const Vec y[],PetscScalar *z)
11: {
12: PetscScalar awork[128],*work = awork;
16: if (nv > 128) {
17: PetscMalloc(nv*sizeof(PetscScalar),&work);
18: }
19: VecMDot_Seq(xin,nv,y,work);
20: MPI_Allreduce(work,z,nv,MPIU_SCALAR,MPIU_SUM,((PetscObject)xin)->comm);
21: if (nv > 128) {
22: PetscFree(work);
23: }
24: return(0);
25: }
29: PetscErrorCode VecMTDot_MPI(Vec xin,PetscInt nv,const Vec y[],PetscScalar *z)
30: {
31: PetscScalar awork[128],*work = awork;
35: if (nv > 128) {
36: PetscMalloc(nv*sizeof(PetscScalar),&work);
37: }
38: VecMTDot_Seq(xin,nv,y,work);
39: MPI_Allreduce(work,z,nv,MPIU_SCALAR,MPIU_SUM,((PetscObject)xin)->comm);
40: if (nv > 128) {
41: PetscFree(work);
42: }
43: return(0);
44: }
46: #include <../src/vec/vec/impls/seq/ftn-kernels/fnorm.h>
49: PetscErrorCode VecNorm_MPI(Vec xin,NormType type,PetscReal *z)
50: {
51: PetscReal sum,work = 0.0;
52: PetscScalar *xx = *(PetscScalar**)xin->data;
54: PetscInt n = xin->map->n;
55: PetscBLASInt one = 1,bn = PetscBLASIntCast(n);
58: if (type == NORM_2 || type == NORM_FROBENIUS) {
59: work = BLASnrm2_(&bn,xx,&one);
60: work *= work;
61: MPI_Allreduce(&work,&sum,1,MPIU_REAL,MPIU_SUM,((PetscObject)xin)->comm);
62: *z = sqrt(sum);
63: PetscLogFlops(2.0*xin->map->n);
64: } else if (type == NORM_1) {
65: /* Find the local part */
66: VecNorm_Seq(xin,NORM_1,&work);
67: /* Find the global max */
68: MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_SUM,((PetscObject)xin)->comm);
69: } else if (type == NORM_INFINITY) {
70: /* Find the local max */
71: VecNorm_Seq(xin,NORM_INFINITY,&work);
72: /* Find the global max */
73: MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_MAX,((PetscObject)xin)->comm);
74: } else if (type == NORM_1_AND_2) {
75: PetscReal temp[2];
76: VecNorm_Seq(xin,NORM_1,temp);
77: VecNorm_Seq(xin,NORM_2,temp+1);
78: temp[1] = temp[1]*temp[1];
79: MPI_Allreduce(temp,z,2,MPIU_REAL,MPIU_SUM,((PetscObject)xin)->comm);
80: z[1] = sqrt(z[1]);
81: }
82: return(0);
83: }
85: /*
86: These two functions are the MPI reduction operation used for max and min with index
87: The call below to MPI_Op_create() converts the function Vec[Max,Min]_Local() to the
88: MPI operator Vec[Max,Min]_Local_Op.
89: */
90: MPI_Op VecMax_Local_Op = 0;
91: MPI_Op VecMin_Local_Op = 0;
96: void MPIAPI VecMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
97: {
98: PetscReal *xin = (PetscReal *)in,*xout = (PetscReal*)out;
101: if (*datatype != MPIU_REAL) {
102: (*PetscErrorPrintf)("Can only handle MPIU_REAL data types");
103: MPI_Abort(MPI_COMM_WORLD,1);
104: }
105: if (xin[0] > xout[0]) {
106: xout[0] = xin[0];
107: xout[1] = xin[1];
108: } else if (xin[0] == xout[0]) {
109: xout[1] = PetscMin(xin[1],xout[1]);
110: }
111: PetscFunctionReturnVoid(); /* cannot return a value */
112: }
118: void MPIAPI VecMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
119: {
120: PetscReal *xin = (PetscReal *)in,*xout = (PetscReal*)out;
123: if (*datatype != MPIU_REAL) {
124: (*PetscErrorPrintf)("Can only handle MPIU_REAL data types");
125: MPI_Abort(MPI_COMM_WORLD,1);
126: }
127: if (xin[0] < xout[0]) {
128: xout[0] = xin[0];
129: xout[1] = xin[1];
130: } else if (xin[0] == xout[0]) {
131: xout[1] = PetscMin(xin[1],xout[1]);
132: }
133: PetscFunctionReturnVoid();
134: }
139: PetscErrorCode VecMax_MPI(Vec xin,PetscInt *idx,PetscReal *z)
140: {
142: PetscReal work;
145: /* Find the local max */
146: VecMax_Seq(xin,idx,&work);
148: /* Find the global max */
149: if (!idx) {
150: MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_MAX,((PetscObject)xin)->comm);
151: } else {
152: PetscReal work2[2],z2[2];
153: PetscInt rstart;
154: rstart = xin->map->rstart;
155: work2[0] = work;
156: work2[1] = *idx + rstart;
157: MPI_Allreduce(work2,z2,2,MPIU_REAL,VecMax_Local_Op,((PetscObject)xin)->comm);
158: *z = z2[0];
159: *idx = (PetscInt)z2[1];
160: }
161: return(0);
162: }
166: PetscErrorCode VecMin_MPI(Vec xin,PetscInt *idx,PetscReal *z)
167: {
169: PetscReal work;
172: /* Find the local Min */
173: VecMin_Seq(xin,idx,&work);
175: /* Find the global Min */
176: if (!idx) {
177: MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_MIN,((PetscObject)xin)->comm);
178: } else {
179: PetscReal work2[2],z2[2];
180: PetscInt rstart;
182: VecGetOwnershipRange(xin,&rstart,PETSC_NULL);
183: work2[0] = work;
184: work2[1] = *idx + rstart;
185: MPI_Allreduce(work2,z2,2,MPIU_REAL,VecMin_Local_Op,((PetscObject)xin)->comm);
186: *z = z2[0];
187: *idx = (PetscInt)z2[1];
188: }
189: return(0);
190: }