Actual source code: vinv.c
1: /*$Id: vinv.c,v 1.67 2001/03/23 23:21:18 balay Exp $*/
2: /*
3: Some useful vector utility functions.
4: */
5: #include petscvec.h
6: #include src/vec/vecimpl.h
8: /*@C
9: VecStrideNorm - Computes the norm of subvector of a vector defined
10: by a starting point and a stride.
12: Collective on Vec
14: Input Parameter:
15: + v - the vector
16: . start - starting point of the subvector (defined by a stride)
17: - ntype - type of norm, one of NORM_1, NORM_2, NORM_INFINITY
19: Output Parameter:
20: . norm - the norm
22: Notes:
23: One must call VecSetBlockSize() before this routine to set the stride
24: information, or use a vector created from a multicomponent DA.
26: If x is the array representing the vector x then this computes the norm
27: of the array (x[start],x[start+stride],x[start+2*stride], ....)
29: This is useful for computing, say the norm of the pressure variable when
30: the pressure is stored (interlaced) with other variables, say density etc.
32: This will only work if the desire subvector is a stride subvector
34: Level: advanced
36: Concepts: norm^on stride of vector
37: Concepts: stride^norm
39: .seealso: VecNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax()
40: @*/
41: int VecStrideNorm(Vec v,int start,NormType ntype,PetscReal *norm)
42: {
43: int i,n,ierr,bs;
44: Scalar *x;
45: PetscReal tnorm;
46: MPI_Comm comm;
50: VecGetLocalSize(v,&n);
51: VecGetArray(v,&x);
52: PetscObjectGetComm((PetscObject)v,&comm);
54: bs = v->bs;
55: if (start >= bs) {
56: SETERRQ2(1,"Start of stride subvector (%d) is too large for striden
57: Have you set the vector blocksize (%d) correctly with VecSetBlockSize()?",start,bs);
58: }
59: x += start;
61: if (ntype == NORM_2) {
62: Scalar sum = 0.0;
63: for (i=0; i<n; i+=bs) {
64: sum += x[i]*(PetscConj(x[i]));
65: }
66: tnorm = PetscRealPart(sum);
67: ierr = MPI_Allreduce(&tnorm,norm,1,MPI_DOUBLE,MPI_SUM,comm);
68: *norm = sqrt(*norm);
69: } else if (ntype == NORM_1) {
70: tnorm = 0.0;
71: for (i=0; i<n; i+=bs) {
72: tnorm += PetscAbsScalar(x[i]);
73: }
74: ierr = MPI_Allreduce(&tnorm,norm,1,MPI_DOUBLE,MPI_SUM,comm);
75: } else if (ntype == NORM_INFINITY) {
76: PetscReal tmp;
77: tnorm = 0.0;
79: for (i=0; i<n; i+=bs) {
80: if ((tmp = PetscAbsScalar(x[i])) > tnorm) tnorm = tmp;
81: /* check special case of tmp == NaN */
82: if (tmp != tmp) {tnorm = tmp; break;}
83: }
84: ierr = MPI_Allreduce(&tnorm,norm,1,MPI_DOUBLE,MPI_MAX,comm);
85: } else {
86: SETERRQ(1,"Unknown norm type");
87: }
89: VecRestoreArray(v,&x);
90: return(0);
91: }
93: /*@C
94: VecStrideMax - Computes the maximum of subvector of a vector defined
95: by a starting point and a stride and optionally its location.
97: Collective on Vec
99: Input Parameter:
100: + v - the vector
101: - start - starting point of the subvector (defined by a stride)
103: Output Parameter:
104: + index - the location where the maximum occurred (not supported, pass PETSC_NULL,
105: if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
106: - norm - the max
108: Notes:
109: One must call VecSetBlockSize() before this routine to set the stride
110: information, or use a vector created from a multicomponent DA.
112: If xa is the array representing the vector x, then this computes the max
113: of the array (xa[start],xa[start+stride],xa[start+2*stride], ....)
115: This is useful for computing, say the maximum of the pressure variable when
116: the pressure is stored (interlaced) with other variables, e.g., density, etc.
117: This will only work if the desire subvector is a stride subvector.
119: Level: advanced
121: Concepts: maximum^on stride of vector
122: Concepts: stride^maximum
124: .seealso: VecMax(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin()
125: @*/
126: int VecStrideMax(Vec v,int start,int *index,PetscReal *norm)
127: {
128: int i,n,ierr,bs;
129: Scalar *x;
130: PetscReal max,tmp;
131: MPI_Comm comm;
135: if (index) {
136: SETERRQ(1,"No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
137: }
138: VecGetLocalSize(v,&n);
139: VecGetArray(v,&x);
140: PetscObjectGetComm((PetscObject)v,&comm);
142: bs = v->bs;
143: if (start >= bs) {
144: SETERRQ2(1,"Start of stride subvector (%d) is too large for striden
145: Have you set the vector blocksize (%d) correctly with VecSetBlockSize()?",start,bs);
146: }
147: x += start;
149: if (!n) {
150: max = PETSC_MIN;
151: } else {
152: #if defined(PETSC_USE_COMPLEX)
153: max = PetscRealPart(x[0]);
154: #else
155: max = x[0];
156: #endif
157: for (i=bs; i<n; i+=bs) {
158: #if defined(PETSC_USE_COMPLEX)
159: if ((tmp = PetscRealPart(x[i])) > max) { max = tmp;}
160: #else
161: if ((tmp = x[i]) > max) { max = tmp; }
162: #endif
163: }
164: }
165: ierr = MPI_Allreduce(&max,norm,1,MPI_DOUBLE,MPI_MAX,comm);
167: VecRestoreArray(v,&x);
168: return(0);
169: }
171: /*@C
172: VecStrideMin - Computes the minimum of subvector of a vector defined
173: by a starting point and a stride and optionally its location.
175: Collective on Vec
177: Input Parameter:
178: + v - the vector
179: - start - starting point of the subvector (defined by a stride)
181: Output Parameter:
182: + index - the location where the minimum occurred (not supported, pass PETSC_NULL,
183: if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
184: - norm - the min
186: Level: intermediate
188: Notes:
189: One must call VecSetBlockSize() before this routine to set the stride
190: information, or use a vector created from a multicomponent DA.
192: If xa is the array representing the vector x, then this computes the min
193: of the array (xa[start],xa[start+stride],xa[start+2*stride], ....)
195: This is useful for computing, say the minimum of the pressure variable when
196: the pressure is stored (interlaced) with other variables, e.g., density, etc.
197: This will only work if the desire subvector is a stride subvector.
199: Concepts: minimum^on stride of vector
200: Concepts: stride^minimum
202: .seealso: VecMin(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMax()
203: @*/
204: int VecStrideMin(Vec v,int start,int *index,PetscReal *norm)
205: {
206: int i,n,ierr,bs;
207: Scalar *x;
208: PetscReal min,tmp;
209: MPI_Comm comm;
213: if (index) {
214: SETERRQ(1,"No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
215: }
216: VecGetLocalSize(v,&n);
217: VecGetArray(v,&x);
218: PetscObjectGetComm((PetscObject)v,&comm);
220: bs = v->bs;
221: if (start >= bs) {
222: SETERRQ2(1,"Start of stride subvector (%d) is too large for striden
223: Have you set the vector blocksize (%d) correctly with VecSetBlockSize()?",start,bs);
224: }
225: x += start;
227: if (!n) {
228: min = PETSC_MAX;
229: } else {
230: #if defined(PETSC_USE_COMPLEX)
231: min = PetscRealPart(x[0]);
232: #else
233: min = x[0];
234: #endif
235: for (i=bs; i<n; i+=bs) {
236: #if defined(PETSC_USE_COMPLEX)
237: if ((tmp = PetscRealPart(x[i])) < min) { min = tmp;}
238: #else
239: if ((tmp = x[i]) < min) { min = tmp; }
240: #endif
241: }
242: }
243: ierr = MPI_Allreduce(&min,norm,1,MPI_DOUBLE,MPI_MIN,comm);
245: VecRestoreArray(v,&x);
246: return(0);
247: }
249: /*@
250: VecStrideGatherAll - Gathers all the single components from a multi-component vector into
251: seperate vectors.
253: Collective on Vec
255: Input Parameter:
256: + v - the vector
257: - addv - one of ADD_VALUES,SET_VALUES,MAX_VALUES
259: Output Parameter:
260: . s - the location where the subvectors are stored
262: Notes:
263: One must call VecSetBlockSize() before this routine to set the stride
264: information, or use a vector created from a multicomponent DA.
266: If x is the array representing the vector x then this gathers
267: the arrays (x[start],x[start+stride],x[start+2*stride], ....)
268: for start=0,1,2,...bs-1
270: The parallel layout of the vector and the subvector must be the same;
271: i.e., nlocal of v = stride*(nlocal of s)
273: Level: advanced
275: Concepts: gather^into strided vector
277: .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGather(),
278: VecStrideScatterAll()
279: @*/
280: int VecStrideGatherAll(Vec v,Vec *s,InsertMode addv)
281: {
282: int i,n,ierr,bs,ns,j;
283: Scalar *x,**y;
288: VecGetLocalSize(v,&n);
289: VecGetLocalSize(*s,&ns);
290: VecGetArray(v,&x);
291: bs = v->bs;
293: PetscMalloc(bs*sizeof(PetscReal*),&y);
294: for (i=0; i<bs; i++) {
295: VecGetArray(s[i],&y[i]);
296: }
298: if (n != ns*bs) {
299: SETERRQ2(1,"Subvector length * blocksize %d not correct for gather from original vector %d",ns*bs,n);
300: }
301: n = n/bs;
303: if (addv == INSERT_VALUES) {
304: for (j=0; j<bs; j++) {
305: for (i=0; i<n; i++) {
306: y[j][i] = x[bs*i+j];
307: }
308: }
309: } else if (addv == ADD_VALUES) {
310: for (j=0; j<bs; j++) {
311: for (i=0; i<n; i++) {
312: y[j][i] += x[bs*i+j];
313: }
314: }
315: #if !defined(PETSC_USE_COMPLEX)
316: } else if (addv == MAX_VALUES) {
317: for (j=0; j<bs; j++) {
318: for (i=0; i<n; i++) {
319: y[j][i] = PetscMax(y[j][i],x[bs*i+j]);
320: }
321: }
322: #endif
323: } else {
324: SETERRQ(1,"Unknown insert type");
325: }
327: VecRestoreArray(v,&x);
328: for (i=0; i<bs; i++) {
329: VecRestoreArray(s[i],&y[i]);
330: }
331: PetscFree(y);
332: return(0);
333: }
334: /*@
335: VecStrideScatterAll - Scatters all the single components from seperate vectors into
336: a multi-component vector.
338: Collective on Vec
340: Input Parameter:
341: + s - the location where the subvectors are stored
342: - addv - one of ADD_VALUES,SET_VALUES,MAX_VALUES
344: Output Parameter:
345: . v - the multicomponent vector
347: Notes:
348: One must call VecSetBlockSize() before this routine to set the stride
349: information, or use a vector created from a multicomponent DA.
351: The parallel layout of the vector and the subvector must be the same;
352: i.e., nlocal of v = stride*(nlocal of s)
354: Level: advanced
356: Concepts: scatter^into strided vector
358: .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGather(),
359: VecStrideScatterAll()
360: @*/
361: int VecStrideScatterAll(Vec *s,Vec v,InsertMode addv)
362: {
363: int i,n,ierr,bs,ns,j;
364: Scalar *x,**y;
369: VecGetLocalSize(v,&n);
370: VecGetLocalSize(*s,&ns);
371: VecGetArray(v,&x);
372: bs = v->bs;
374: PetscMalloc(bs*sizeof(PetscReal*),&y);
375: for (i=0; i<bs; i++) {
376: VecGetArray(s[i],&y[i]);
377: }
379: if (n != ns*bs) {
380: SETERRQ2(1,"Subvector length * blocksize %d not correct for gather from original vector %d",ns*bs,n);
381: }
382: n = n/bs;
384: if (addv == INSERT_VALUES) {
385: for (j=0; j<bs; j++) {
386: for (i=0; i<n; i++) {
387: x[bs*i+j] = y[j][i];
388: }
389: }
390: } else if (addv == ADD_VALUES) {
391: for (j=0; j<bs; j++) {
392: for (i=0; i<n; i++) {
393: x[bs*i+j] += y[j][i];
394: }
395: }
396: #if !defined(PETSC_USE_COMPLEX)
397: } else if (addv == MAX_VALUES) {
398: for (j=0; j<bs; j++) {
399: for (i=0; i<n; i++) {
400: x[bs*i+j] = PetscMax(y[j][i],x[bs*i+j]);
401: }
402: }
403: #endif
404: } else {
405: SETERRQ(1,"Unknown insert type");
406: }
408: VecRestoreArray(v,&x);
409: for (i=0; i<bs; i++) {
410: VecRestoreArray(s[i],&y[i]);
411: }
412: PetscFree(y);
413: return(0);
414: }
416: /*@
417: VecStrideGather - Gathers a single component from a multi-component vector into
418: another vector.
420: Collective on Vec
422: Input Parameter:
423: + v - the vector
424: . start - starting point of the subvector (defined by a stride)
425: - addv - one of ADD_VALUES,SET_VALUES,MAX_VALUES
427: Output Parameter:
428: . s - the location where the subvector is stored
430: Notes:
431: One must call VecSetBlockSize() before this routine to set the stride
432: information, or use a vector created from a multicomponent DA.
434: If x is the array representing the vector x then this gathers
435: the array (x[start],x[start+stride],x[start+2*stride], ....)
437: The parallel layout of the vector and the subvector must be the same;
438: i.e., nlocal of v = stride*(nlocal of s)
440: Level: advanced
442: Concepts: gather^into strided vector
444: .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGatherAll(),
445: VecStrideScatterAll()
446: @*/
447: int VecStrideGather(Vec v,int start,Vec s,InsertMode addv)
448: {
449: int i,n,ierr,bs,ns;
450: Scalar *x,*y;
455: VecGetLocalSize(v,&n);
456: VecGetLocalSize(s,&ns);
457: VecGetArray(v,&x);
458: VecGetArray(s,&y);
460: bs = v->bs;
461: if (start >= bs) {
462: SETERRQ2(1,"Start of stride subvector (%d) is too large for striden
463: Have you set the vector blocksize (%d) correctly with VecSetBlockSize()?",start,bs);
464: }
465: if (n != ns*bs) {
466: SETERRQ2(1,"Subvector length * blocksize %d not correct for gather from original vector %d",ns*bs,n);
467: }
468: x += start;
469: n = n/bs;
471: if (addv == INSERT_VALUES) {
472: for (i=0; i<n; i++) {
473: y[i] = x[bs*i];
474: }
475: } else if (addv == ADD_VALUES) {
476: for (i=0; i<n; i++) {
477: y[i] += x[bs*i];
478: }
479: #if !defined(PETSC_USE_COMPLEX)
480: } else if (addv == MAX_VALUES) {
481: for (i=0; i<n; i++) {
482: y[i] = PetscMax(y[i],x[bs*i]);
483: }
484: #endif
485: } else {
486: SETERRQ(1,"Unknown insert type");
487: }
489: VecRestoreArray(v,&x);
490: VecRestoreArray(s,&y);
491: return(0);
492: }
494: /*@
495: VecStrideScatter - Scatters a single component from a vector into a multi-component vector.
497: Collective on Vec
499: Input Parameter:
500: + s - the single-component vector
501: . start - starting point of the subvector (defined by a stride)
502: - addv - one of ADD_VALUES,SET_VALUES,MAX_VALUES
504: Output Parameter:
505: . v - the location where the subvector is scattered (the multi-component vector)
507: Notes:
508: One must call VecSetBlockSize() on the multi-component vector before this
509: routine to set the stride information, or use a vector created from a multicomponent DA.
511: The parallel layout of the vector and the subvector must be the same;
512: i.e., nlocal of v = stride*(nlocal of s)
514: Level: advanced
516: Concepts: scatter^into strided vector
518: .seealso: VecStrideNorm(), VecStrideGather(), VecStrideMin(), VecStrideMax(), VecStrideGatherAll(),
519: VecStrideScatterAll()
520: @*/
521: int VecStrideScatter(Vec s,int start,Vec v,InsertMode addv)
522: {
523: int i,n,ierr,bs,ns;
524: Scalar *x,*y;
529: VecGetLocalSize(v,&n);
530: VecGetLocalSize(s,&ns);
531: VecGetArray(v,&x);
532: VecGetArray(s,&y);
534: bs = v->bs;
535: if (start >= bs) {
536: SETERRQ2(1,"Start of stride subvector (%d) is too large for striden
537: Have you set the vector blocksize (%d) correctly with VecSetBlockSize()?",start,bs);
538: }
539: if (n != ns*bs) {
540: SETERRQ2(1,"Subvector length * blocksize %d not correct for scatter to multicomponent vector %d",ns*bs,n);
541: }
542: x += start;
543: n = n/bs;
546: if (addv == INSERT_VALUES) {
547: for (i=0; i<n; i++) {
548: x[bs*i] = y[i];
549: }
550: } else if (addv == ADD_VALUES) {
551: for (i=0; i<n; i++) {
552: x[bs*i] += y[i];
553: }
554: #if !defined(PETSC_USE_COMPLEX)
555: } else if (addv == MAX_VALUES) {
556: for (i=0; i<n; i++) {
557: x[bs*i] = PetscMax(y[i],x[bs*i]);
558: }
559: #endif
560: } else {
561: SETERRQ(1,"Unknown insert type");
562: }
565: VecRestoreArray(v,&x);
566: VecRestoreArray(s,&y);
567: return(0);
568: }
570: int VecReciprocal_Default(Vec v)
571: {
572: int i,n,ierr;
573: Scalar *x;
576: VecGetLocalSize(v,&n);
577: VecGetArray(v,&x);
578: for (i=0; i<n; i++) {
579: if (x[i] != 0.0) x[i] = 1.0/x[i];
580: }
581: VecRestoreArray(v,&x);
582: return(0);
583: }
585: /*@
586: VecSum - Computes the sum of all the components of a vector.
588: Collective on Vec
590: Input Parameter:
591: . v - the vector
593: Output Parameter:
594: . sum - the result
596: Level: beginner
598: Concepts: sum^of vector entries
600: .seealso: VecNorm()
601: @*/
602: int VecSum(Vec v,Scalar *sum)
603: {
604: int i,n,ierr;
605: Scalar *x,lsum = 0.0;
609: VecGetLocalSize(v,&n);
610: VecGetArray(v,&x);
611: for (i=0; i<n; i++) {
612: lsum += x[i];
613: }
614: MPI_Allreduce(&lsum,sum,1,MPIU_SCALAR,PetscSum_Op,v->comm);
615: VecRestoreArray(v,&x);
616: return(0);
617: }
619: /*@
620: VecShift - Shifts all of the components of a vector by computing
621: x[i] = x[i] + shift.
623: Collective on Vec
625: Input Parameters:
626: + v - the vector
627: - shift - the shift
629: Output Parameter:
630: . v - the shifted vector
632: Level: intermediate
634: Concepts: vector^adding constant
636: @*/
637: int VecShift(const Scalar *shift,Vec v)
638: {
639: int i,n,ierr;
640: Scalar *x,lsum = *shift;
644: VecGetLocalSize(v,&n);
645: VecGetArray(v,&x);
646: for (i=0; i<n; i++) {
647: x[i] += lsum;
648: }
649: VecRestoreArray(v,&x);
650: return(0);
651: }
653: /*@
654: VecAbs - Replaces every element in a vector with its absolute value.
656: Collective on Vec
658: Input Parameters:
659: . v - the vector
661: Level: intermediate
663: Concepts: vector^absolute value
665: @*/
666: int VecAbs(Vec v)
667: {
668: int i,n,ierr;
669: Scalar *x;
673: VecGetLocalSize(v,&n);
674: VecGetArray(v,&x);
675: for (i=0; i<n; i++) {
676: x[i] = PetscAbsScalar(x[i]);
677: }
678: VecRestoreArray(v,&x);
679: return(0);
680: }
683: /*@
684: VecEqual - Compares two vectors.
686: Collective on Vec
688: Input Parameters:
689: + vec1 - the first matrix
690: - vec2 - the second matrix
692: Output Parameter:
693: . flg - PETSC_TRUE if the vectors are equal; PETSC_FALSE otherwise.
695: Level: intermediate
697: Concepts: equal^two vectors
698: Concepts: vector^equality
700: @*/
701: int VecEqual(Vec vec1,Vec vec2,PetscTruth *flg)
702: {
703: Scalar * v1,*v2;
704: int n1,n2,ierr;
705: PetscTruth flg1;
708: VecGetSize(vec1,&n1);
709: VecGetSize(vec2,&n2);
710: if (n1 != n2) {
711: flg1 = PETSC_FALSE;
712: } else {
713: VecGetArray(vec1,&v1);
714: VecGetArray(vec2,&v2);
715: PetscMemcmp(v1,v2,n1*sizeof(Scalar),&flg1);
716: VecRestoreArray(vec1,&v1);
717: VecRestoreArray(vec2,&v2);
718: }
720: /* combine results from all processors */
721: MPI_Allreduce(&flg1,flg,1,MPI_INT,MPI_MIN,vec1->comm);
722:
724: return(0);
725: }