Actual source code: vinv.c
1: /*$Id: vinv.c,v 1.71 2001/09/11 16:31:37 bsmith Exp $*/
2: /*
3: Some useful vector utility functions.
4: */
5: #include src/vec/vecimpl.h
9: /*@C
10: VecStrideScale - Scales a subvector of a vector defined
11: by a starting point and a stride.
13: Collective on Vec
15: Input Parameter:
16: + v - the vector
17: . start - starting point of the subvector (defined by a stride)
18: - scale - value to multiply each subvector entry by
20: Notes:
21: One must call VecSetBlockSize() before this routine to set the stride
22: information, or use a vector created from a multicomponent DA.
24: This will only work if the desire subvector is a stride subvector
26: Level: advanced
28: Concepts: scale^on stride of vector
29: Concepts: stride^scale
31: .seealso: VecNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideScale()
32: @*/
33: int VecStrideScale(Vec v,int start,PetscScalar *scale)
34: {
35: int i,n,ierr,bs;
36: PetscScalar *x,xscale = *scale;
40: VecGetLocalSize(v,&n);
41: VecGetArrayFast(v,&x);
43: bs = v->bs;
44: if (start >= bs) {
45: SETERRQ2(1,"Start of stride subvector (%d) is too large for stride\n\
46: Have you set the vector blocksize (%d) correctly with VecSetBlockSize()?",start,bs);
47: }
48: x += start;
50: for (i=0; i<n; i+=bs) {
51: x[i] *= xscale;
52: }
53: x -= start;
54: VecRestoreArrayFast(v,&x);
55: return(0);
56: }
60: /*@C
61: VecStrideNorm - Computes the norm of subvector of a vector defined
62: by a starting point and a stride.
64: Collective on Vec
66: Input Parameter:
67: + v - the vector
68: . start - starting point of the subvector (defined by a stride)
69: - ntype - type of norm, one of NORM_1, NORM_2, NORM_INFINITY
71: Output Parameter:
72: . norm - the norm
74: Notes:
75: One must call VecSetBlockSize() before this routine to set the stride
76: information, or use a vector created from a multicomponent DA.
78: If x is the array representing the vector x then this computes the norm
79: of the array (x[start],x[start+stride],x[start+2*stride], ....)
81: This is useful for computing, say the norm of the pressure variable when
82: the pressure is stored (interlaced) with other variables, say density etc.
84: This will only work if the desire subvector is a stride subvector
86: Level: advanced
88: Concepts: norm^on stride of vector
89: Concepts: stride^norm
91: .seealso: VecNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax()
92: @*/
93: int VecStrideNorm(Vec v,int start,NormType ntype,PetscReal *nrm)
94: {
95: int i,n,ierr,bs;
96: PetscScalar *x;
97: PetscReal tnorm;
98: MPI_Comm comm;
102: VecGetLocalSize(v,&n);
103: VecGetArrayFast(v,&x);
104: PetscObjectGetComm((PetscObject)v,&comm);
106: bs = v->bs;
107: if (start >= bs) {
108: SETERRQ2(1,"Start of stride subvector (%d) is too large for stride\n\
109: Have you set the vector blocksize (%d) correctly with VecSetBlockSize()?",start,bs);
110: }
111: x += start;
113: if (ntype == NORM_2) {
114: PetscScalar sum = 0.0;
115: for (i=0; i<n; i+=bs) {
116: sum += x[i]*(PetscConj(x[i]));
117: }
118: tnorm = PetscRealPart(sum);
119: MPI_Allreduce(&tnorm,nrm,1,MPIU_REAL,MPI_SUM,comm);
120: *nrm = sqrt(*nrm);
121: } else if (ntype == NORM_1) {
122: tnorm = 0.0;
123: for (i=0; i<n; i+=bs) {
124: tnorm += PetscAbsScalar(x[i]);
125: }
126: MPI_Allreduce(&tnorm,nrm,1,MPIU_REAL,MPI_SUM,comm);
127: } else if (ntype == NORM_INFINITY) {
128: PetscReal tmp;
129: tnorm = 0.0;
131: for (i=0; i<n; i+=bs) {
132: if ((tmp = PetscAbsScalar(x[i])) > tnorm) tnorm = tmp;
133: /* check special case of tmp == NaN */
134: if (tmp != tmp) {tnorm = tmp; break;}
135: }
136: MPI_Allreduce(&tnorm,nrm,1,MPIU_REAL,MPI_MAX,comm);
137: } else {
138: SETERRQ(1,"Unknown norm type");
139: }
141: VecRestoreArrayFast(v,&x);
142: return(0);
143: }
147: /*@C
148: VecStrideMax - Computes the maximum of subvector of a vector defined
149: by a starting point and a stride and optionally its location.
151: Collective on Vec
153: Input Parameter:
154: + v - the vector
155: - start - starting point of the subvector (defined by a stride)
157: Output Parameter:
158: + index - the location where the maximum occurred (not supported, pass PETSC_NULL,
159: if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
160: - nrm - the max
162: Notes:
163: One must call VecSetBlockSize() before this routine to set the stride
164: information, or use a vector created from a multicomponent DA.
166: If xa is the array representing the vector x, then this computes the max
167: of the array (xa[start],xa[start+stride],xa[start+2*stride], ....)
169: This is useful for computing, say the maximum of the pressure variable when
170: the pressure is stored (interlaced) with other variables, e.g., density, etc.
171: This will only work if the desire subvector is a stride subvector.
173: Level: advanced
175: Concepts: maximum^on stride of vector
176: Concepts: stride^maximum
178: .seealso: VecMax(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin()
179: @*/
180: int VecStrideMax(Vec v,int start,int *idex,PetscReal *nrm)
181: {
182: int i,n,ierr,bs;
183: PetscScalar *x;
184: PetscReal max,tmp;
185: MPI_Comm comm;
189: if (idex) {
190: SETERRQ(1,"No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
191: }
192: VecGetLocalSize(v,&n);
193: VecGetArrayFast(v,&x);
194: PetscObjectGetComm((PetscObject)v,&comm);
196: bs = v->bs;
197: if (start >= bs) {
198: SETERRQ2(1,"Start of stride subvector (%d) is too large for stride\n\
199: Have you set the vector blocksize (%d) correctly with VecSetBlockSize()?",start,bs);
200: }
201: x += start;
203: if (!n) {
204: max = PETSC_MIN;
205: } else {
206: #if defined(PETSC_USE_COMPLEX)
207: max = PetscRealPart(x[0]);
208: #else
209: max = x[0];
210: #endif
211: for (i=bs; i<n; i+=bs) {
212: #if defined(PETSC_USE_COMPLEX)
213: if ((tmp = PetscRealPart(x[i])) > max) { max = tmp;}
214: #else
215: if ((tmp = x[i]) > max) { max = tmp; }
216: #endif
217: }
218: }
219: MPI_Allreduce(&max,nrm,1,MPIU_REAL,MPI_MAX,comm);
221: VecRestoreArrayFast(v,&x);
222: return(0);
223: }
227: /*@C
228: VecStrideMin - Computes the minimum of subvector of a vector defined
229: by a starting point and a stride and optionally its location.
231: Collective on Vec
233: Input Parameter:
234: + v - the vector
235: - start - starting point of the subvector (defined by a stride)
237: Output Parameter:
238: + idex - the location where the minimum occurred (not supported, pass PETSC_NULL,
239: if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
240: - nrm - the min
242: Level: advanced
244: Notes:
245: One must call VecSetBlockSize() before this routine to set the stride
246: information, or use a vector created from a multicomponent DA.
248: If xa is the array representing the vector x, then this computes the min
249: of the array (xa[start],xa[start+stride],xa[start+2*stride], ....)
251: This is useful for computing, say the minimum of the pressure variable when
252: the pressure is stored (interlaced) with other variables, e.g., density, etc.
253: This will only work if the desire subvector is a stride subvector.
255: Concepts: minimum^on stride of vector
256: Concepts: stride^minimum
258: .seealso: VecMin(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMax()
259: @*/
260: int VecStrideMin(Vec v,int start,int *idex,PetscReal *nrm)
261: {
262: int i,n,ierr,bs;
263: PetscScalar *x;
264: PetscReal min,tmp;
265: MPI_Comm comm;
269: if (idex) {
270: SETERRQ(1,"No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
271: }
272: VecGetLocalSize(v,&n);
273: VecGetArrayFast(v,&x);
274: PetscObjectGetComm((PetscObject)v,&comm);
276: bs = v->bs;
277: if (start >= bs) {
278: SETERRQ2(1,"Start of stride subvector (%d) is too large for stride\n\
279: Have you set the vector blocksize (%d) correctly with VecSetBlockSize()?",start,bs);
280: }
281: x += start;
283: if (!n) {
284: min = PETSC_MAX;
285: } else {
286: #if defined(PETSC_USE_COMPLEX)
287: min = PetscRealPart(x[0]);
288: #else
289: min = x[0];
290: #endif
291: for (i=bs; i<n; i+=bs) {
292: #if defined(PETSC_USE_COMPLEX)
293: if ((tmp = PetscRealPart(x[i])) < min) { min = tmp;}
294: #else
295: if ((tmp = x[i]) < min) { min = tmp; }
296: #endif
297: }
298: }
299: MPI_Allreduce(&min,nrm,1,MPIU_REAL,MPI_MIN,comm);
301: VecRestoreArrayFast(v,&x);
302: return(0);
303: }
307: /*@C
308: VecStrideScaleAll - Scales the subvectors of a vector defined
309: by a starting point and a stride.
311: Collective on Vec
313: Input Parameter:
314: + v - the vector
315: - scales - values to multiply each subvector entry by
317: Notes:
318: One must call VecSetBlockSize() before this routine to set the stride
319: information, or use a vector created from a multicomponent DA.
322: Level: advanced
324: Concepts: scale^on stride of vector
325: Concepts: stride^scale
327: .seealso: VecNorm(), VecStrideScale(), VecScale(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax()
328: @*/
329: int VecStrideScaleAll(Vec v,PetscScalar *scales)
330: {
331: int i,j,n,ierr,bs;
332: PetscScalar *x;
336: VecGetLocalSize(v,&n);
337: VecGetArrayFast(v,&x);
339: bs = v->bs;
341: /* need to provide optimized code for each bs */
342: for (i=0; i<n; i+=bs) {
343: for (j=0; j<bs; j++) {
344: x[i+j] *= scales[j];
345: }
346: }
347: VecRestoreArrayFast(v,&x);
348: return(0);
349: }
353: /*@C
354: VecStrideNormAll - Computes the norms subvectors of a vector defined
355: by a starting point and a stride.
357: Collective on Vec
359: Input Parameter:
360: + v - the vector
361: - ntype - type of norm, one of NORM_1, NORM_2, NORM_INFINITY
363: Output Parameter:
364: . nrm - the norms
366: Notes:
367: One must call VecSetBlockSize() before this routine to set the stride
368: information, or use a vector created from a multicomponent DA.
370: If x is the array representing the vector x then this computes the norm
371: of the array (x[start],x[start+stride],x[start+2*stride], ....)
373: This is useful for computing, say the norm of the pressure variable when
374: the pressure is stored (interlaced) with other variables, say density etc.
376: This will only work if the desire subvector is a stride subvector
378: Level: advanced
380: Concepts: norm^on stride of vector
381: Concepts: stride^norm
383: .seealso: VecNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax()
384: @*/
385: int VecStrideNormAll(Vec v,NormType ntype,PetscReal *nrm)
386: {
387: int i,j,n,ierr,bs;
388: PetscScalar *x;
389: PetscReal tnorm[128];
390: MPI_Comm comm;
394: VecGetLocalSize(v,&n);
395: VecGetArrayFast(v,&x);
396: PetscObjectGetComm((PetscObject)v,&comm);
398: bs = v->bs;
399: if (bs > 128) SETERRQ(1,"Currently supports only blocksize up to 128");
401: if (ntype == NORM_2) {
402: PetscScalar sum[128];
403: for (j=0; j<bs; j++) sum[j] = 0.0;
404: for (i=0; i<n; i+=bs) {
405: for (j=0; j<bs; j++) {
406: sum[j] += x[i+j]*(PetscConj(x[i+j]));
407: }
408: }
409: for (j=0; j<bs; j++) {
410: tnorm[j] = PetscRealPart(sum[j]);
411: }
412: MPI_Allreduce(tnorm,nrm,bs,MPIU_REAL,MPI_SUM,comm);
413: for (j=0; j<bs; j++) {
414: nrm[j] = sqrt(nrm[j]);
415: }
416: } else if (ntype == NORM_1) {
417: for (j=0; j<bs; j++) {
418: tnorm[j] = 0.0;
419: }
420: for (i=0; i<n; i+=bs) {
421: for (j=0; j<bs; j++) {
422: tnorm[j] += PetscAbsScalar(x[i+j]);
423: }
424: }
425: MPI_Allreduce(tnorm,nrm,bs,MPIU_REAL,MPI_SUM,comm);
426: } else if (ntype == NORM_INFINITY) {
427: PetscReal tmp;
428: for (j=0; j<bs; j++) {
429: tnorm[j] = 0.0;
430: }
432: for (i=0; i<n; i+=bs) {
433: for (j=0; j<bs; j++) {
434: if ((tmp = PetscAbsScalar(x[i+j])) > tnorm[j]) tnorm[j] = tmp;
435: /* check special case of tmp == NaN */
436: if (tmp != tmp) {tnorm[j] = tmp; break;}
437: }
438: }
439: MPI_Allreduce(tnorm,nrm,bs,MPIU_REAL,MPI_MAX,comm);
440: } else {
441: SETERRQ(1,"Unknown norm type");
442: }
444: VecRestoreArrayFast(v,&x);
445: return(0);
446: }
450: /*@C
451: VecStrideMaxAll - Computes the maximums of subvectors of a vector defined
452: by a starting point and a stride and optionally its location.
454: Collective on Vec
456: Input Parameter:
457: . v - the vector
459: Output Parameter:
460: + index - the location where the maximum occurred (not supported, pass PETSC_NULL,
461: if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
462: - nrm - the maximums
464: Notes:
465: One must call VecSetBlockSize() before this routine to set the stride
466: information, or use a vector created from a multicomponent DA.
468: This is useful for computing, say the maximum of the pressure variable when
469: the pressure is stored (interlaced) with other variables, e.g., density, etc.
470: This will only work if the desire subvector is a stride subvector.
472: Level: advanced
474: Concepts: maximum^on stride of vector
475: Concepts: stride^maximum
477: .seealso: VecMax(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin()
478: @*/
479: int VecStrideMaxAll(Vec v,int *idex,PetscReal *nrm)
480: {
481: int i,j,n,ierr,bs;
482: PetscScalar *x;
483: PetscReal max[128],tmp;
484: MPI_Comm comm;
488: if (idex) {
489: SETERRQ(1,"No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
490: }
491: VecGetLocalSize(v,&n);
492: VecGetArrayFast(v,&x);
493: PetscObjectGetComm((PetscObject)v,&comm);
495: bs = v->bs;
496: if (bs > 128) SETERRQ(1,"Currently supports only blocksize up to 128");
498: if (!n) {
499: for (j=0; j<bs; j++) {
500: max[j] = PETSC_MIN;
501: }
502: } else {
503: for (j=0; j<bs; j++) {
504: #if defined(PETSC_USE_COMPLEX)
505: max[j] = PetscRealPart(x[j]);
506: #else
507: max[j] = x[j];
508: #endif
509: }
510: for (i=bs; i<n; i+=bs) {
511: for (j=0; j<bs; j++) {
512: #if defined(PETSC_USE_COMPLEX)
513: if ((tmp = PetscRealPart(x[i+j])) > max[j]) { max[j] = tmp;}
514: #else
515: if ((tmp = x[i+j]) > max[j]) { max[j] = tmp; }
516: #endif
517: }
518: }
519: }
520: MPI_Allreduce(max,nrm,bs,MPIU_REAL,MPI_MAX,comm);
522: VecRestoreArrayFast(v,&x);
523: return(0);
524: }
528: /*@C
529: VecStrideMinAll - Computes the minimum of subvector of a vector defined
530: by a starting point and a stride and optionally its location.
532: Collective on Vec
534: Input Parameter:
535: . v - the vector
537: Output Parameter:
538: + idex - the location where the minimum occurred (not supported, pass PETSC_NULL,
539: if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
540: - nrm - the minimums
542: Level: advanced
544: Notes:
545: One must call VecSetBlockSize() before this routine to set the stride
546: information, or use a vector created from a multicomponent DA.
548: This is useful for computing, say the minimum of the pressure variable when
549: the pressure is stored (interlaced) with other variables, e.g., density, etc.
550: This will only work if the desire subvector is a stride subvector.
552: Concepts: minimum^on stride of vector
553: Concepts: stride^minimum
555: .seealso: VecMin(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMax()
556: @*/
557: int VecStrideMinAll(Vec v,int *idex,PetscReal *nrm)
558: {
559: int i,n,ierr,bs,j;
560: PetscScalar *x;
561: PetscReal min[128],tmp;
562: MPI_Comm comm;
566: if (idex) {
567: SETERRQ(1,"No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
568: }
569: VecGetLocalSize(v,&n);
570: VecGetArrayFast(v,&x);
571: PetscObjectGetComm((PetscObject)v,&comm);
573: bs = v->bs;
574: if (bs > 128) SETERRQ(1,"Currently supports only blocksize up to 128");
576: if (!n) {
577: for (j=0; j<bs; j++) {
578: min[j] = PETSC_MAX;
579: }
580: } else {
581: for (j=0; j<bs; j++) {
582: #if defined(PETSC_USE_COMPLEX)
583: min[j] = PetscRealPart(x[j]);
584: #else
585: min[j] = x[j];
586: #endif
587: }
588: for (i=bs; i<n; i+=bs) {
589: for (j=0; j<bs; j++) {
590: #if defined(PETSC_USE_COMPLEX)
591: if ((tmp = PetscRealPart(x[i+j])) < min[j]) { min[j] = tmp;}
592: #else
593: if ((tmp = x[i+j]) < min[j]) { min[j] = tmp; }
594: #endif
595: }
596: }
597: }
598: MPI_Allreduce(min,nrm,bs,MPIU_REAL,MPI_MIN,comm);
600: VecRestoreArrayFast(v,&x);
601: return(0);
602: }
604: /*----------------------------------------------------------------------------------------------*/
607: /*@
608: VecStrideGatherAll - Gathers all the single components from a multi-component vector into
609: seperate vectors.
611: Collective on Vec
613: Input Parameter:
614: + v - the vector
615: - addv - one of ADD_VALUES,SET_VALUES,MAX_VALUES
617: Output Parameter:
618: . s - the location where the subvectors are stored
620: Notes:
621: One must call VecSetBlockSize() before this routine to set the stride
622: information, or use a vector created from a multicomponent DA.
624: If x is the array representing the vector x then this gathers
625: the arrays (x[start],x[start+stride],x[start+2*stride], ....)
626: for start=0,1,2,...bs-1
628: The parallel layout of the vector and the subvector must be the same;
629: i.e., nlocal of v = stride*(nlocal of s)
631: Not optimized; could be easily
633: Level: advanced
635: Concepts: gather^into strided vector
637: .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGather(),
638: VecStrideScatterAll()
639: @*/
640: int VecStrideGatherAll(Vec v,Vec *s,InsertMode addv)
641: {
642: int i,n,ierr,bs,j,k,*bss,nv,jj,nvc;
643: PetscScalar *x,**y;
648: VecGetLocalSize(v,&n);
649: VecGetArrayFast(v,&x);
650: bs = v->bs;
652: PetscMalloc(bs*sizeof(PetscReal*),&y);
653: PetscMalloc(bs*sizeof(int),&bss);
654: nv = 0;
655: nvc = 0;
656: for (i=0; i<bs; i++) {
657: VecGetBlockSize(s[i],&bss[i]);
658: VecGetArrayFast(s[i],&y[i]);
659: nvc += bss[i];
660: nv++;
661: if (nvc > bs) SETERRQ(1,"Number of subvectors in subvectors > number of vectors in main vector");
662: if (nvc == bs) break;
663: }
665: n = n/bs;
667: jj = 0;
668: if (addv == INSERT_VALUES) {
669: for (j=0; j<nv; j++) {
670: for (k=0; k<bss[j]; k++) {
671: for (i=0; i<n; i++) {
672: y[j][i*bss[j] + k] = x[bs*i+jj+k];
673: }
674: }
675: jj += bss[j];
676: }
677: } else if (addv == ADD_VALUES) {
678: for (j=0; j<nv; j++) {
679: for (k=0; k<bss[j]; k++) {
680: for (i=0; i<n; i++) {
681: y[j][i*bss[j] + k] += x[bs*i+jj+k];
682: }
683: }
684: jj += bss[j];
685: }
686: #if !defined(PETSC_USE_COMPLEX)
687: } else if (addv == MAX_VALUES) {
688: for (j=0; j<nv; j++) {
689: for (k=0; k<bss[j]; k++) {
690: for (i=0; i<n; i++) {
691: y[j][i*bss[j] + k] = PetscMax(y[j][i*bss[j] + k],x[bs*i+jj+k]);
692: }
693: }
694: jj += bss[j];
695: }
696: #endif
697: } else {
698: SETERRQ(1,"Unknown insert type");
699: }
701: VecRestoreArrayFast(v,&x);
702: for (i=0; i<nv; i++) {
703: VecRestoreArrayFast(s[i],&y[i]);
704: }
705: PetscFree(y);
706: PetscFree(bss);
707: return(0);
708: }
712: /*@
713: VecStrideScatterAll - Scatters all the single components from seperate vectors into
714: a multi-component vector.
716: Collective on Vec
718: Input Parameter:
719: + s - the location where the subvectors are stored
720: - addv - one of ADD_VALUES,SET_VALUES,MAX_VALUES
722: Output Parameter:
723: . v - the multicomponent vector
725: Notes:
726: One must call VecSetBlockSize() before this routine to set the stride
727: information, or use a vector created from a multicomponent DA.
729: The parallel layout of the vector and the subvector must be the same;
730: i.e., nlocal of v = stride*(nlocal of s)
732: Not optimized; could be easily
734: Level: advanced
736: Concepts: scatter^into strided vector
738: .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGather(),
739: VecStrideScatterAll()
740: @*/
741: int VecStrideScatterAll(Vec *s,Vec v,InsertMode addv)
742: {
743: int i,n,ierr,bs,j,jj,k,*bss,nv,nvc;
744: PetscScalar *x,**y;
749: VecGetLocalSize(v,&n);
750: VecGetArrayFast(v,&x);
751: bs = v->bs;
753: PetscMalloc(bs*sizeof(PetscReal*),&y);
754: PetscMalloc(bs*sizeof(int*),&bss);
755: nv = 0;
756: nvc = 0;
757: for (i=0; i<bs; i++) {
758: VecGetBlockSize(s[i],&bss[i]);
759: VecGetArrayFast(s[i],&y[i]);
760: nvc += bss[i];
761: nv++;
762: if (nvc > bs) SETERRQ(1,"Number of subvectors in subvectors > number of vectors in main vector");
763: if (nvc == bs) break;
764: }
766: n = n/bs;
768: jj = 0;
769: if (addv == INSERT_VALUES) {
770: for (j=0; j<nv; j++) {
771: for (k=0; k<bss[j]; k++) {
772: for (i=0; i<n; i++) {
773: x[bs*i+jj+k] = y[j][i*bss[j] + k];
774: }
775: }
776: jj += bss[j];
777: }
778: } else if (addv == ADD_VALUES) {
779: for (j=0; j<nv; j++) {
780: for (k=0; k<bss[j]; k++) {
781: for (i=0; i<n; i++) {
782: x[bs*i+jj+k] += y[j][i*bss[j] + k];
783: }
784: }
785: jj += bss[j];
786: }
787: #if !defined(PETSC_USE_COMPLEX)
788: } else if (addv == MAX_VALUES) {
789: for (j=0; j<nv; j++) {
790: for (k=0; k<bss[j]; k++) {
791: for (i=0; i<n; i++) {
792: x[bs*i+jj+k] = PetscMax(x[bs*i+jj+k],y[j][i*bss[j] + k]);
793: }
794: }
795: jj += bss[j];
796: }
797: #endif
798: } else {
799: SETERRQ(1,"Unknown insert type");
800: }
802: VecRestoreArrayFast(v,&x);
803: for (i=0; i<nv; i++) {
804: VecRestoreArrayFast(s[i],&y[i]);
805: }
806: PetscFree(y);
807: PetscFree(bss);
808: return(0);
809: }
813: /*@
814: VecStrideGather - Gathers a single component from a multi-component vector into
815: another vector.
817: Collective on Vec
819: Input Parameter:
820: + v - the vector
821: . start - starting point of the subvector (defined by a stride)
822: - addv - one of ADD_VALUES,SET_VALUES,MAX_VALUES
824: Output Parameter:
825: . s - the location where the subvector is stored
827: Notes:
828: One must call VecSetBlockSize() before this routine to set the stride
829: information, or use a vector created from a multicomponent DA.
831: If x is the array representing the vector x then this gathers
832: the array (x[start],x[start+stride],x[start+2*stride], ....)
834: The parallel layout of the vector and the subvector must be the same;
835: i.e., nlocal of v = stride*(nlocal of s)
837: Not optimized; could be easily
839: Level: advanced
841: Concepts: gather^into strided vector
843: .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGatherAll(),
844: VecStrideScatterAll()
845: @*/
846: int VecStrideGather(Vec v,int start,Vec s,InsertMode addv)
847: {
848: int i,n,ierr,bs,ns;
849: PetscScalar *x,*y;
854: VecGetLocalSize(v,&n);
855: VecGetLocalSize(s,&ns);
856: VecGetArrayFast(v,&x);
857: VecGetArrayFast(s,&y);
859: bs = v->bs;
860: if (start >= bs) {
861: SETERRQ2(1,"Start of stride subvector (%d) is too large for stride\n\
862: Have you set the vector blocksize (%d) correctly with VecSetBlockSize()?",start,bs);
863: }
864: if (n != ns*bs) {
865: SETERRQ2(1,"Subvector length * blocksize %d not correct for gather from original vector %d",ns*bs,n);
866: }
867: x += start;
868: n = n/bs;
870: if (addv == INSERT_VALUES) {
871: for (i=0; i<n; i++) {
872: y[i] = x[bs*i];
873: }
874: } else if (addv == ADD_VALUES) {
875: for (i=0; i<n; i++) {
876: y[i] += x[bs*i];
877: }
878: #if !defined(PETSC_USE_COMPLEX)
879: } else if (addv == MAX_VALUES) {
880: for (i=0; i<n; i++) {
881: y[i] = PetscMax(y[i],x[bs*i]);
882: }
883: #endif
884: } else {
885: SETERRQ(1,"Unknown insert type");
886: }
888: VecRestoreArrayFast(v,&x);
889: VecRestoreArrayFast(s,&y);
890: return(0);
891: }
895: /*@
896: VecStrideScatter - Scatters a single component from a vector into a multi-component vector.
898: Collective on Vec
900: Input Parameter:
901: + s - the single-component vector
902: . start - starting point of the subvector (defined by a stride)
903: - addv - one of ADD_VALUES,SET_VALUES,MAX_VALUES
905: Output Parameter:
906: . v - the location where the subvector is scattered (the multi-component vector)
908: Notes:
909: One must call VecSetBlockSize() on the multi-component vector before this
910: routine to set the stride information, or use a vector created from a multicomponent DA.
912: The parallel layout of the vector and the subvector must be the same;
913: i.e., nlocal of v = stride*(nlocal of s)
915: Not optimized; could be easily
917: Level: advanced
919: Concepts: scatter^into strided vector
921: .seealso: VecStrideNorm(), VecStrideGather(), VecStrideMin(), VecStrideMax(), VecStrideGatherAll(),
922: VecStrideScatterAll()
923: @*/
924: int VecStrideScatter(Vec s,int start,Vec v,InsertMode addv)
925: {
926: int i,n,ierr,bs,ns;
927: PetscScalar *x,*y;
932: VecGetLocalSize(v,&n);
933: VecGetLocalSize(s,&ns);
934: VecGetArrayFast(v,&x);
935: VecGetArrayFast(s,&y);
937: bs = v->bs;
938: if (start >= bs) {
939: SETERRQ2(1,"Start of stride subvector (%d) is too large for stride\n\
940: Have you set the vector blocksize (%d) correctly with VecSetBlockSize()?",start,bs);
941: }
942: if (n != ns*bs) {
943: SETERRQ2(1,"Subvector length * blocksize %d not correct for scatter to multicomponent vector %d",ns*bs,n);
944: }
945: x += start;
946: n = n/bs;
949: if (addv == INSERT_VALUES) {
950: for (i=0; i<n; i++) {
951: x[bs*i] = y[i];
952: }
953: } else if (addv == ADD_VALUES) {
954: for (i=0; i<n; i++) {
955: x[bs*i] += y[i];
956: }
957: #if !defined(PETSC_USE_COMPLEX)
958: } else if (addv == MAX_VALUES) {
959: for (i=0; i<n; i++) {
960: x[bs*i] = PetscMax(y[i],x[bs*i]);
961: }
962: #endif
963: } else {
964: SETERRQ(1,"Unknown insert type");
965: }
968: VecRestoreArrayFast(v,&x);
969: VecRestoreArrayFast(s,&y);
970: return(0);
971: }
975: int VecReciprocal_Default(Vec v)
976: {
977: int i,n,ierr;
978: PetscScalar *x;
981: VecGetLocalSize(v,&n);
982: VecGetArrayFast(v,&x);
983: for (i=0; i<n; i++) {
984: if (x[i] != 0.0) x[i] = 1.0/x[i];
985: }
986: VecRestoreArrayFast(v,&x);
987: return(0);
988: }
992: /*@
993: VecSqrt - Replaces each component of a vector by the square root of its magnitude.
995: Not collective
997: Input Parameter:
998: . v - The vector
1000: Output Parameter:
1001: . v - The vector square root
1003: Level: beginner
1005: Note: The actual function is sqrt(|x_i|)
1007: .keywords: vector, sqrt, square root
1008: @*/
1009: int VecSqrt(Vec v)
1010: {
1011: PetscScalar *x;
1012: int i, n;
1013: int ierr;
1017: VecGetLocalSize(v, &n);
1018: VecGetArrayFast(v, &x);
1019: for(i = 0; i < n; i++) {
1020: x[i] = sqrt(PetscAbsScalar(x[i]));
1021: }
1022: VecRestoreArrayFast(v, &x);
1023: return(0);
1024: }
1028: /*@
1029: VecSum - Computes the sum of all the components of a vector.
1031: Collective on Vec
1033: Input Parameter:
1034: . v - the vector
1036: Output Parameter:
1037: . sum - the result
1039: Level: beginner
1041: Concepts: sum^of vector entries
1043: .seealso: VecNorm()
1044: @*/
1045: int VecSum(Vec v,PetscScalar *sum)
1046: {
1047: int i,n,ierr;
1048: PetscScalar *x,lsum = 0.0;
1052: VecGetLocalSize(v,&n);
1053: VecGetArrayFast(v,&x);
1054: for (i=0; i<n; i++) {
1055: lsum += x[i];
1056: }
1057: MPI_Allreduce(&lsum,sum,1,MPIU_SCALAR,PetscSum_Op,v->comm);
1058: VecRestoreArrayFast(v,&x);
1059: return(0);
1060: }
1064: /*@
1065: VecShift - Shifts all of the components of a vector by computing
1066: x[i] = x[i] + shift.
1068: Collective on Vec
1070: Input Parameters:
1071: + v - the vector
1072: - shift - the shift
1074: Output Parameter:
1075: . v - the shifted vector
1077: Level: intermediate
1079: Concepts: vector^adding constant
1081: @*/
1082: int VecShift(const PetscScalar *shift,Vec v)
1083: {
1084: int i,n,ierr;
1085: PetscScalar *x,lsum = *shift;
1089: VecGetLocalSize(v,&n);
1090: VecGetArrayFast(v,&x);
1091: for (i=0; i<n; i++) {
1092: x[i] += lsum;
1093: }
1094: VecRestoreArrayFast(v,&x);
1095: return(0);
1096: }
1100: /*@
1101: VecAbs - Replaces every element in a vector with its absolute value.
1103: Collective on Vec
1105: Input Parameters:
1106: . v - the vector
1108: Level: intermediate
1110: Concepts: vector^absolute value
1112: @*/
1113: int VecAbs(Vec v)
1114: {
1115: int i,n,ierr;
1116: PetscScalar *x;
1120: VecGetLocalSize(v,&n);
1121: VecGetArrayFast(v,&x);
1122: for (i=0; i<n; i++) {
1123: x[i] = PetscAbsScalar(x[i]);
1124: }
1125: VecRestoreArrayFast(v,&x);
1126: return(0);
1127: }
1131: /*@
1132: VecPermute - Permutes a vector in place using the given ordering.
1134: Input Parameters:
1135: + vec - The vector
1136: . order - The ordering
1137: - inv - The flag for inverting the permutation
1139: Level: beginner
1141: Note: This function does not yet support parallel Index Sets
1143: .seealso: MatPermute()
1144: .keywords: vec, permute
1145: @*/
1146: int VecPermute(Vec x, IS row, PetscTruth inv)
1147: {
1148: PetscScalar *array, *newArray;
1149: int *idx;
1150: int i;
1151: int ierr;
1154: ISGetIndices(row, &idx);
1155: VecGetArrayFast(x, &array);
1156: PetscMalloc((x->n+1) * sizeof(PetscScalar), &newArray);
1157: #ifdef PETSC_USE_BOPT_g
1158: for(i = 0; i < x->n; i++) {
1159: if ((idx[i] < 0) || (idx[i] >= x->n)) {
1160: SETERRQ2(PETSC_ERR_ARG_CORRUPT, "Permutation index %d is out of bounds: %d", i, idx[i]);
1161: }
1162: }
1163: #endif
1164: if (inv == PETSC_FALSE) {
1165: for(i = 0; i < x->n; i++) newArray[i] = array[idx[i]];
1166: } else {
1167: for(i = 0; i < x->n; i++) newArray[idx[i]] = array[i];
1168: }
1169: VecRestoreArrayFast(x, &array);
1170: ISRestoreIndices(row, &idx);
1171: VecReplaceArray(x, newArray);
1172: return(0);
1173: }
1177: /*@
1178: VecEqual - Compares two vectors.
1180: Collective on Vec
1182: Input Parameters:
1183: + vec1 - the first matrix
1184: - vec2 - the second matrix
1186: Output Parameter:
1187: . flg - PETSC_TRUE if the vectors are equal; PETSC_FALSE otherwise.
1189: Level: intermediate
1191: Concepts: equal^two vectors
1192: Concepts: vector^equality
1194: @*/
1195: int VecEqual(Vec vec1,Vec vec2,PetscTruth *flg)
1196: {
1197: PetscScalar *v1,*v2;
1198: int n1,n2,ierr;
1199: PetscTruth flg1;
1202: VecGetSize(vec1,&n1);
1203: VecGetSize(vec2,&n2);
1204: if (vec1 == vec2) {
1205: flg1 = PETSC_TRUE;
1206: } else if (n1 != n2) {
1207: flg1 = PETSC_FALSE;
1208: } else {
1209: VecGetArrayFast(vec1,&v1);
1210: VecGetArrayFast(vec2,&v2);
1211: PetscMemcmp(v1,v2,n1*sizeof(PetscScalar),&flg1);
1212: VecRestoreArrayFast(vec1,&v1);
1213: VecRestoreArrayFast(vec2,&v2);
1214: }
1216: /* combine results from all processors */
1217: MPI_Allreduce(&flg1,flg,1,MPI_INT,MPI_MIN,vec1->comm);
1218:
1220: return(0);
1221: }