Actual source code: veccomp0.h
slepc-3.7.3 2016-09-29
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: #include <petsc/private/vecimpl.h>
24: #if defined(__WITH_MPI__)
26: #else
28: #endif
36: PetscErrorCode __SUF__(VecDot_Comp)(Vec a,Vec b,PetscScalar *z)
37: {
38: PetscScalar sum = 0.0,work;
39: PetscInt i;
41: Vec_Comp *as = (Vec_Comp*)a->data,*bs = (Vec_Comp*)b->data;
44: SlepcValidVecComp(a);
45: SlepcValidVecComp(b);
46: if (as->x[0]->ops->dot_local) {
47: for (i=0,sum=0.0;i<as->n->n;i++) {
48: as->x[i]->ops->dot_local(as->x[i],bs->x[i],&work);
49: sum += work;
50: }
51: #if defined(__WITH_MPI__)
52: work = sum;
53: MPI_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)a));
54: #endif
55: } else {
56: for (i=0,sum=0.0;i<as->n->n;i++) {
57: VecDot(as->x[i],bs->x[i],&work);
58: sum += work;
59: }
60: }
61: *z = sum;
62: return(0);
63: }
67: PetscErrorCode __SUF__(VecMDot_Comp)(Vec a,PetscInt n,const Vec b[],PetscScalar *z)
68: {
69: PetscScalar *work,*work0,*r;
71: Vec_Comp *as = (Vec_Comp*)a->data;
72: Vec *bx;
73: PetscInt i,j;
76: SlepcValidVecComp(a);
77: for (i=0;i<n;i++) SlepcValidVecComp(b[i]);
79: if (as->n->n == 0) {
80: *z = 0;
81: return(0);
82: }
84: PetscMalloc2(n,&work0,n,&bx);
86: #if defined(__WITH_MPI__)
87: if (as->x[0]->ops->mdot_local) {
88: r = work0; work = z;
89: } else
90: #endif
91: {
92: r = z; work = work0;
93: }
95: /* z[i] <- a.x' * b[i].x */
96: for (i=0;i<n;i++) bx[i] = ((Vec_Comp*)b[i]->data)->x[0];
97: if (as->x[0]->ops->mdot_local) {
98: as->x[0]->ops->mdot_local(as->x[0],n,bx,r);
99: } else {
100: VecMDot(as->x[0],n,bx,r);
101: }
102: for (j=0;j<as->n->n;j++) {
103: for (i=0;i<n;i++) bx[i] = ((Vec_Comp*)b[i]->data)->x[j];
104: if (as->x[0]->ops->mdot_local) {
105: as->x[j]->ops->mdot_local(as->x[j],n,bx,work);
106: } else {
107: VecMDot(as->x[j],n,bx,work);
108: }
109: for (i=0;i<n;i++) r[i] += work[i];
110: }
112: /* If def(__WITH_MPI__) and exists mdot_local */
113: #if defined(__WITH_MPI__)
114: if (as->x[0]->ops->mdot_local) {
115: /* z[i] <- Allreduce(work[i]) */
116: MPI_Allreduce(r,z,n,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)a));
117: }
118: #endif
120: PetscFree2(work0,bx);
121: return(0);
122: }
126: PetscErrorCode __SUF__(VecTDot_Comp)(Vec a,Vec b,PetscScalar *z)
127: {
128: PetscScalar sum = 0.0,work;
129: PetscInt i;
131: Vec_Comp *as = (Vec_Comp*)a->data,*bs = (Vec_Comp*)b->data;
134: SlepcValidVecComp(a);
135: SlepcValidVecComp(b);
136: if (as->x[0]->ops->tdot_local) {
137: for (i=0,sum=0.0;i<as->n->n;i++) {
138: as->x[i]->ops->tdot_local(as->x[i],bs->x[i],&work);
139: sum += work;
140: }
141: #if defined(__WITH_MPI__)
142: work = sum;
143: MPI_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)a));
144: #endif
145: } else {
146: for (i=0,sum=0.0;i<as->n->n;i++) {
147: VecTDot(as->x[i],bs->x[i],&work);
148: sum += work;
149: }
150: }
151: *z = sum;
152: return(0);
153: }
157: PetscErrorCode __SUF__(VecMTDot_Comp)(Vec a,PetscInt n,const Vec b[],PetscScalar *z)
158: {
159: PetscScalar *work,*work0,*r;
161: Vec_Comp *as = (Vec_Comp*)a->data;
162: Vec *bx;
163: PetscInt i,j;
166: SlepcValidVecComp(a);
167: for (i=0;i<n;i++) SlepcValidVecComp(b[i]);
169: if (as->n->n == 0) {
170: *z = 0;
171: return(0);
172: }
174: PetscMalloc2(n,&work0,n,&bx);
176: #if defined(__WITH_MPI__)
177: if (as->x[0]->ops->mtdot_local) {
178: r = work0; work = z;
179: } else
180: #endif
181: {
182: r = z; work = work0;
183: }
185: /* z[i] <- a.x' * b[i].x */
186: for (i=0;i<n;i++) bx[i] = ((Vec_Comp*)b[i]->data)->x[0];
187: if (as->x[0]->ops->mtdot_local) {
188: as->x[0]->ops->mtdot_local(as->x[0],n,bx,r);
189: } else {
190: VecMTDot(as->x[0],n,bx,r);
191: }
192: for (j=0;j<as->n->n;j++) {
193: for (i=0;i<n;i++) bx[i] = ((Vec_Comp*)b[i]->data)->x[j];
194: if (as->x[0]->ops->mtdot_local) {
195: as->x[j]->ops->mtdot_local(as->x[j],n,bx,work);
196: } else {
197: VecMTDot(as->x[j],n,bx,work);
198: }
199: for (i=0;i<n;i++) r[i] += work[i];
200: }
202: /* If def(__WITH_MPI__) and exists mtdot_local */
203: #if defined(__WITH_MPI__)
204: if (as->x[0]->ops->mtdot_local) {
205: /* z[i] <- Allreduce(work[i]) */
206: MPI_Allreduce(r,z,n,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)a));
207: }
208: #endif
210: PetscFree2(work0,bx);
211: return(0);
212: }
216: PetscErrorCode __SUF__(VecNorm_Comp)(Vec a,NormType t,PetscReal *norm)
217: {
218: PetscReal work[3],s=0.0;
220: Vec_Comp *as = (Vec_Comp*)a->data;
221: PetscInt i;
224: SlepcValidVecComp(a);
225: /* Initialize norm */
226: switch (t) {
227: case NORM_1: case NORM_INFINITY: *norm = 0.0; break;
228: case NORM_2: case NORM_FROBENIUS: *norm = 1.0; s = 0.0; break;
229: case NORM_1_AND_2: norm[0] = 0.0; norm[1] = 1.0; s = 0.0; break;
230: }
231: for (i=0;i<as->n->n;i++) {
232: if (as->x[0]->ops->norm_local) {
233: as->x[0]->ops->norm_local(as->x[i],t,work);
234: } else {
235: VecNorm(as->x[i],t,work);
236: }
237: /* norm+= work */
238: switch (t) {
239: case NORM_1: *norm+= *work; break;
240: case NORM_2: case NORM_FROBENIUS: AddNorm2(norm,&s,*work); break;
241: case NORM_1_AND_2: norm[0]+= work[0]; AddNorm2(&norm[1],&s,work[1]); break;
242: case NORM_INFINITY: *norm = PetscMax(*norm,*work); break;
243: }
244: }
246: /* If def(__WITH_MPI__) and exists norm_local */
247: #if defined(__WITH_MPI__)
248: if (as->x[0]->ops->norm_local) {
249: PetscReal work0[3];
250: /* norm <- Allreduce(work) */
251: switch (t) {
252: case NORM_1:
253: work[0] = *norm;
254: MPI_Allreduce(work,norm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)a));
255: break;
256: case NORM_2: case NORM_FROBENIUS:
257: work[0] = *norm; work[1] = s;
258: MPI_Allreduce(work,work0,1,MPIU_NORM2,MPIU_NORM2_SUM,PetscObjectComm((PetscObject)a));
259: *norm = GetNorm2(work0[0],work0[1]);
260: break;
261: case NORM_1_AND_2:
262: work[0] = norm[0]; work[1] = norm[1]; work[2] = s;
263: MPI_Allreduce(work,work0,1,MPIU_NORM1_AND_2,MPIU_NORM2_SUM,PetscObjectComm((PetscObject)a));
264: norm[0] = work0[0];
265: norm[1] = GetNorm2(work0[1],work0[2]);
266: break;
267: case NORM_INFINITY:
268: work[0] = *norm;
269: MPI_Allreduce(work,norm,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)a));
270: break;
271: }
272: }
273: #else
274: /* Norm correction */
275: switch (t) {
276: case NORM_2: case NORM_FROBENIUS: *norm = GetNorm2(*norm,s); break;
277: case NORM_1_AND_2: norm[1] = GetNorm2(norm[1],s); break;
278: default: ;
279: }
280: #endif
281: return(0);
282: }
286: PetscErrorCode __SUF__(VecDotNorm2_Comp)(Vec v,Vec w,PetscScalar *dp,PetscScalar *nm)
287: {
288: PetscErrorCode ierr;
289: PetscScalar dp0,nm0,dp1,nm1;
290: const PetscScalar *vx,*wx;
291: Vec_Comp *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data;
292: PetscInt i,n;
293: PetscBool t0,t1;
294: #if defined(__WITH_MPI__)
295: PetscScalar work[4];
296: #endif
299: /* Compute recursively the local part */
300: dp0 = nm0 = 0.0;
301: PetscObjectTypeCompare((PetscObject)v,VECCOMP,&t0);
302: PetscObjectTypeCompare((PetscObject)w,VECCOMP,&t1);
303: if (t0 && t1) {
304: SlepcValidVecComp(v);
305: SlepcValidVecComp(w);
306: for (i=0;i<vs->n->n;i++) {
307: VecDotNorm2_Comp_Seq(vs->x[i],ws->x[i],&dp1,&nm1);
308: dp0 += dp1;
309: nm0 += nm1;
310: }
311: } else if (!t0 && !t1) {
312: VecGetLocalSize(v,&n);
313: VecGetArrayRead(v,&vx);
314: VecGetArrayRead(w,&wx);
315: for (i=0;i<n;i++) {
316: dp0 += vx[i]*PetscConj(wx[i]);
317: nm0 += wx[i]*PetscConj(wx[i]);
318: }
319: VecRestoreArrayRead(v,&vx);
320: VecRestoreArrayRead(w,&wx);
321: } else SETERRQ(PetscObjectComm((PetscObject)v),PETSC_ERR_ARG_INCOMP,"Incompatible vector types");
323: #if defined(__WITH_MPI__)
324: /* [dp, nm] <- Allreduce([dp0, nm0]) */
325: work[0] = dp0; work[1] = nm0;
326: MPI_Allreduce(work,&work[2],2,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)v));
327: *dp = work[2]; *nm = work[3];
328: #else
329: *dp = dp0; *nm = nm0;
330: #endif
331: return(0);
332: }