Actual source code: pvec2.c
petsc-3.3-p5 2012-12-01
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: const PetscScalar *xx;
53: PetscErrorCode ierr;
54: PetscInt n = xin->map->n;
55: PetscBLASInt one = 1,bn = PetscBLASIntCast(n);
58: if (type == NORM_2 || type == NORM_FROBENIUS) {
59: VecGetArrayRead(xin,&xx);
60: work = PetscRealPart(BLASdot_(&bn,xx,&one,xx,&one));
61: VecRestoreArrayRead(xin,&xx);
62: MPI_Allreduce(&work,&sum,1,MPIU_REAL,MPIU_SUM,((PetscObject)xin)->comm);
63: *z = PetscSqrtReal(sum);
64: PetscLogFlops(2.0*xin->map->n);
65: } else if (type == NORM_1) {
66: /* Find the local part */
67: VecNorm_Seq(xin,NORM_1,&work);
68: /* Find the global max */
69: MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_SUM,((PetscObject)xin)->comm);
70: } else if (type == NORM_INFINITY) {
71: /* Find the local max */
72: VecNorm_Seq(xin,NORM_INFINITY,&work);
73: /* Find the global max */
74: MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_MAX,((PetscObject)xin)->comm);
75: } else if (type == NORM_1_AND_2) {
76: PetscReal temp[2];
77: VecNorm_Seq(xin,NORM_1,temp);
78: VecNorm_Seq(xin,NORM_2,temp+1);
79: temp[1] = temp[1]*temp[1];
80: MPI_Allreduce(temp,z,2,MPIU_REAL,MPIU_SUM,((PetscObject)xin)->comm);
81: z[1] = PetscSqrtReal(z[1]);
82: }
83: return(0);
84: }
86: /*
87: These two functions are the MPI reduction operation used for max and min with index
88: The call below to MPI_Op_create() converts the function Vec[Max,Min]_Local() to the
89: MPI operator Vec[Max,Min]_Local_Op.
90: */
91: MPI_Op VecMax_Local_Op = 0;
92: MPI_Op VecMin_Local_Op = 0;
94: EXTERN_C_BEGIN
97: void MPIAPI VecMax_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
98: {
99: PetscReal *xin = (PetscReal *)in,*xout = (PetscReal*)out;
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: EXTERN_C_END
116: EXTERN_C_BEGIN
119: void MPIAPI VecMin_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
120: {
121: PetscReal *xin = (PetscReal *)in,*xout = (PetscReal*)out;
124: if (*datatype != MPIU_REAL) {
125: (*PetscErrorPrintf)("Can only handle MPIU_REAL data types");
126: MPI_Abort(MPI_COMM_SELF,1);
127: }
128: if (xin[0] < xout[0]) {
129: xout[0] = xin[0];
130: xout[1] = xin[1];
131: } else if (xin[0] == xout[0]) {
132: xout[1] = PetscMin(xin[1],xout[1]);
133: }
134: PetscFunctionReturnVoid();
135: }
136: EXTERN_C_END
140: PetscErrorCode VecMax_MPI(Vec xin,PetscInt *idx,PetscReal *z)
141: {
143: PetscReal work;
146: /* Find the local max */
147: VecMax_Seq(xin,idx,&work);
149: /* Find the global max */
150: if (!idx) {
151: MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_MAX,((PetscObject)xin)->comm);
152: } else {
153: PetscReal work2[2],z2[2];
154: PetscInt rstart;
155: rstart = xin->map->rstart;
156: work2[0] = work;
157: work2[1] = *idx + rstart;
158: MPI_Allreduce(work2,z2,2,MPIU_REAL,VecMax_Local_Op,((PetscObject)xin)->comm);
159: *z = z2[0];
160: *idx = (PetscInt)z2[1];
161: }
162: return(0);
163: }
167: PetscErrorCode VecMin_MPI(Vec xin,PetscInt *idx,PetscReal *z)
168: {
170: PetscReal work;
173: /* Find the local Min */
174: VecMin_Seq(xin,idx,&work);
176: /* Find the global Min */
177: if (!idx) {
178: MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_MIN,((PetscObject)xin)->comm);
179: } else {
180: PetscReal work2[2],z2[2];
181: PetscInt rstart;
183: VecGetOwnershipRange(xin,&rstart,PETSC_NULL);
184: work2[0] = work;
185: work2[1] = *idx + rstart;
186: MPI_Allreduce(work2,z2,2,MPIU_REAL,VecMin_Local_Op,((PetscObject)xin)->comm);
187: *z = z2[0];
188: *idx = (PetscInt)z2[1];
189: }
190: return(0);
191: }