xref: /petsc/src/vec/vec/utils/vinv.c (revision fa8c639a052448572ef4c09d337391668dacca91)
1 
2 /*
3      Some useful vector utility functions.
4 */
5 #include <petsc-private/vecimpl.h>          /*I "petscvec.h" I*/
6 #include <../src/vec/vec/impls/mpi/pvecimpl.h>
7 extern MPI_Op VecMax_Local_Op;
8 extern MPI_Op VecMin_Local_Op;
9 
10 #undef __FUNCT__
11 #define __FUNCT__ "VecStrideSet"
12 /*@
13    VecStrideSet - Sets a subvector of a vector defined
14    by a starting point and a stride with a given value
15 
16    Logically Collective on Vec
17 
18    Input Parameter:
19 +  v - the vector
20 .  start - starting point of the subvector (defined by a stride)
21 -  s - value to multiply each subvector entry by
22 
23    Notes:
24    One must call VecSetBlockSize() before this routine to set the stride
25    information, or use a vector created from a multicomponent DMDA.
26 
27    This will only work if the desire subvector is a stride subvector
28 
29    Level: advanced
30 
31    Concepts: scale^on stride of vector
32    Concepts: stride^scale
33 
34 .seealso: VecNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideScale()
35 @*/
36 PetscErrorCode  VecStrideSet(Vec v,PetscInt start,PetscScalar s)
37 {
38   PetscErrorCode ierr;
39   PetscInt       i,n,bs;
40   PetscScalar    *x;
41 
42   PetscFunctionBegin;
43   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
44   PetscValidLogicalCollectiveInt(v,start,2);
45   PetscValidLogicalCollectiveScalar(v,s,3);
46 
47   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
48   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
49   bs   = v->map->bs;
50   if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start);
51   else if (start >= bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Start of stride subvector (%D) is too large for stride\n  Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
52   x += start;
53 
54   for (i=0; i<n; i+=bs) {
55     x[i] = s;
56   }
57   x -= start;
58 
59   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
60   PetscFunctionReturn(0);
61 }
62 
63 #undef __FUNCT__
64 #define __FUNCT__ "VecStrideScale"
65 /*@
66    VecStrideScale - Scales a subvector of a vector defined
67    by a starting point and a stride.
68 
69    Logically Collective on Vec
70 
71    Input Parameter:
72 +  v - the vector
73 .  start - starting point of the subvector (defined by a stride)
74 -  scale - value to multiply each subvector entry by
75 
76    Notes:
77    One must call VecSetBlockSize() before this routine to set the stride
78    information, or use a vector created from a multicomponent DMDA.
79 
80    This will only work if the desire subvector is a stride subvector
81 
82    Level: advanced
83 
84    Concepts: scale^on stride of vector
85    Concepts: stride^scale
86 
87 .seealso: VecNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideScale()
88 @*/
89 PetscErrorCode  VecStrideScale(Vec v,PetscInt start,PetscScalar scale)
90 {
91   PetscErrorCode ierr;
92   PetscInt       i,n,bs;
93   PetscScalar    *x;
94 
95   PetscFunctionBegin;
96   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
97   PetscValidLogicalCollectiveInt(v,start,2);
98   PetscValidLogicalCollectiveScalar(v,scale,3);
99 
100   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
101   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
102   bs   = v->map->bs;
103   if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start);
104   else if (start >= bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Start of stride subvector (%D) is too large for stride\n  Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
105   x += start;
106 
107   for (i=0; i<n; i+=bs) {
108     x[i] *= scale;
109   }
110   x -= start;
111 
112   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
113   PetscFunctionReturn(0);
114 }
115 
116 #undef __FUNCT__
117 #define __FUNCT__ "VecStrideNorm"
118 /*@
119    VecStrideNorm - Computes the norm of subvector of a vector defined
120    by a starting point and a stride.
121 
122    Collective on Vec
123 
124    Input Parameter:
125 +  v - the vector
126 .  start - starting point of the subvector (defined by a stride)
127 -  ntype - type of norm, one of NORM_1, NORM_2, NORM_INFINITY
128 
129    Output Parameter:
130 .  norm - the norm
131 
132    Notes:
133    One must call VecSetBlockSize() before this routine to set the stride
134    information, or use a vector created from a multicomponent DMDA.
135 
136    If x is the array representing the vector x then this computes the norm
137    of the array (x[start],x[start+stride],x[start+2*stride], ....)
138 
139    This is useful for computing, say the norm of the pressure variable when
140    the pressure is stored (interlaced) with other variables, say density etc.
141 
142    This will only work if the desire subvector is a stride subvector
143 
144    Level: advanced
145 
146    Concepts: norm^on stride of vector
147    Concepts: stride^norm
148 
149 .seealso: VecNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax()
150 @*/
151 PetscErrorCode  VecStrideNorm(Vec v,PetscInt start,NormType ntype,PetscReal *nrm)
152 {
153   PetscErrorCode ierr;
154   PetscInt       i,n,bs;
155   PetscScalar    *x;
156   PetscReal      tnorm;
157   MPI_Comm       comm;
158 
159   PetscFunctionBegin;
160   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
161   PetscValidRealPointer(nrm,3);
162   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
163   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
164   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
165 
166   bs = v->map->bs;
167   if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start);
168   else if (start >= bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Start of stride subvector (%D) is too large for stride\n Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
169   x += start;
170 
171   if (ntype == NORM_2) {
172     PetscScalar sum = 0.0;
173     for (i=0; i<n; i+=bs) {
174       sum += x[i]*(PetscConj(x[i]));
175     }
176     tnorm = PetscRealPart(sum);
177     ierr  = MPI_Allreduce(&tnorm,nrm,1,MPIU_REAL,MPIU_SUM,comm);CHKERRQ(ierr);
178     *nrm  = PetscSqrtReal(*nrm);
179   } else if (ntype == NORM_1) {
180     tnorm = 0.0;
181     for (i=0; i<n; i+=bs) {
182       tnorm += PetscAbsScalar(x[i]);
183     }
184     ierr = MPI_Allreduce(&tnorm,nrm,1,MPIU_REAL,MPIU_SUM,comm);CHKERRQ(ierr);
185   } else if (ntype == NORM_INFINITY) {
186     PetscReal tmp;
187     tnorm = 0.0;
188 
189     for (i=0; i<n; i+=bs) {
190       if ((tmp = PetscAbsScalar(x[i])) > tnorm) tnorm = tmp;
191       /* check special case of tmp == NaN */
192       if (tmp != tmp) {tnorm = tmp; break;}
193     }
194     ierr = MPI_Allreduce(&tnorm,nrm,1,MPIU_REAL,MPIU_MAX,comm);CHKERRQ(ierr);
195   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown norm type");
196   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
197   PetscFunctionReturn(0);
198 }
199 
200 #undef __FUNCT__
201 #define __FUNCT__ "VecStrideMax"
202 /*@
203    VecStrideMax - Computes the maximum of subvector of a vector defined
204    by a starting point and a stride and optionally its location.
205 
206    Collective on Vec
207 
208    Input Parameter:
209 +  v - the vector
210 -  start - starting point of the subvector (defined by a stride)
211 
212    Output Parameter:
213 +  index - the location where the maximum occurred  (pass PETSC_NULL if not required)
214 -  nrm - the max
215 
216    Notes:
217    One must call VecSetBlockSize() before this routine to set the stride
218    information, or use a vector created from a multicomponent DMDA.
219 
220    If xa is the array representing the vector x, then this computes the max
221    of the array (xa[start],xa[start+stride],xa[start+2*stride], ....)
222 
223    This is useful for computing, say the maximum of the pressure variable when
224    the pressure is stored (interlaced) with other variables, e.g., density, etc.
225    This will only work if the desire subvector is a stride subvector.
226 
227    Level: advanced
228 
229    Concepts: maximum^on stride of vector
230    Concepts: stride^maximum
231 
232 .seealso: VecMax(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin()
233 @*/
234 PetscErrorCode  VecStrideMax(Vec v,PetscInt start,PetscInt *idex,PetscReal *nrm)
235 {
236   PetscErrorCode ierr;
237   PetscInt       i,n,bs,id;
238   PetscScalar    *x;
239   PetscReal      max,tmp;
240   MPI_Comm       comm;
241 
242   PetscFunctionBegin;
243   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
244   PetscValidRealPointer(nrm,3);
245 
246   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
247   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
248   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
249 
250   bs = v->map->bs;
251   if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start);
252   else if (start >= bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Start of stride subvector (%D) is too large for stride\n Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
253   x += start;
254 
255   id = -1;
256   if (!n) {
257     max = PETSC_MIN_REAL;
258   } else {
259     id  = 0;
260     max = PetscRealPart(x[0]);
261     for (i=bs; i<n; i+=bs) {
262       if ((tmp = PetscRealPart(x[i])) > max) { max = tmp; id = i;}
263     }
264   }
265   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
266 
267   if (!idex) {
268     ierr = MPI_Allreduce(&max,nrm,1,MPIU_REAL,MPIU_MAX,comm);CHKERRQ(ierr);
269   } else {
270     PetscReal in[2],out[2];
271     PetscInt  rstart;
272 
273     ierr  = VecGetOwnershipRange(v,&rstart,PETSC_NULL);CHKERRQ(ierr);
274     in[0] = max;
275     in[1] = rstart+id+start;
276     ierr  = MPI_Allreduce(in,out,2,MPIU_REAL,VecMax_Local_Op,((PetscObject)v)->comm);CHKERRQ(ierr);
277     *nrm  = out[0];
278     *idex = (PetscInt)out[1];
279   }
280   PetscFunctionReturn(0);
281 }
282 
283 #undef __FUNCT__
284 #define __FUNCT__ "VecStrideMin"
285 /*@
286    VecStrideMin - Computes the minimum of subvector of a vector defined
287    by a starting point and a stride and optionally its location.
288 
289    Collective on Vec
290 
291    Input Parameter:
292 +  v - the vector
293 -  start - starting point of the subvector (defined by a stride)
294 
295    Output Parameter:
296 +  idex - the location where the minimum occurred. (pass PETSC_NULL if not required)
297 -  nrm - the min
298 
299    Level: advanced
300 
301    Notes:
302    One must call VecSetBlockSize() before this routine to set the stride
303    information, or use a vector created from a multicomponent DMDA.
304 
305    If xa is the array representing the vector x, then this computes the min
306    of the array (xa[start],xa[start+stride],xa[start+2*stride], ....)
307 
308    This is useful for computing, say the minimum of the pressure variable when
309    the pressure is stored (interlaced) with other variables, e.g., density, etc.
310    This will only work if the desire subvector is a stride subvector.
311 
312    Concepts: minimum^on stride of vector
313    Concepts: stride^minimum
314 
315 .seealso: VecMin(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMax()
316 @*/
317 PetscErrorCode  VecStrideMin(Vec v,PetscInt start,PetscInt *idex,PetscReal *nrm)
318 {
319   PetscErrorCode ierr;
320   PetscInt       i,n,bs,id;
321   PetscScalar    *x;
322   PetscReal      min,tmp;
323   MPI_Comm       comm;
324 
325   PetscFunctionBegin;
326   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
327   PetscValidRealPointer(nrm,4);
328 
329   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
330   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
331   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
332 
333   bs = v->map->bs;
334   if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start);
335   else if (start >= bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Start of stride subvector (%D) is too large for stride\nHave you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
336   x += start;
337 
338   id = -1;
339   if (!n) {
340     min = PETSC_MAX_REAL;
341   } else {
342     id = 0;
343     min = PetscRealPart(x[0]);
344     for (i=bs; i<n; i+=bs) {
345       if ((tmp = PetscRealPart(x[i])) < min) { min = tmp; id = i;}
346     }
347   }
348   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
349 
350   if (!idex) {
351     ierr = MPI_Allreduce(&min,nrm,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr);
352   } else {
353     PetscReal in[2],out[2];
354     PetscInt  rstart;
355 
356     ierr  = VecGetOwnershipRange(v,&rstart,PETSC_NULL);CHKERRQ(ierr);
357     in[0] = min;
358     in[1] = rstart+id;
359     ierr  = MPI_Allreduce(in,out,2,MPIU_REAL,VecMin_Local_Op,((PetscObject)v)->comm);CHKERRQ(ierr);
360     *nrm  = out[0];
361     *idex = (PetscInt)out[1];
362   }
363   PetscFunctionReturn(0);
364 }
365 
366 #undef __FUNCT__
367 #define __FUNCT__ "VecStrideScaleAll"
368 /*@
369    VecStrideScaleAll - Scales the subvectors of a vector defined
370    by a starting point and a stride.
371 
372    Logically Collective on Vec
373 
374    Input Parameter:
375 +  v - the vector
376 -  scales - values to multiply each subvector entry by
377 
378    Notes:
379    One must call VecSetBlockSize() before this routine to set the stride
380    information, or use a vector created from a multicomponent DMDA.
381 
382 
383    Level: advanced
384 
385    Concepts: scale^on stride of vector
386    Concepts: stride^scale
387 
388 .seealso: VecNorm(), VecStrideScale(), VecScale(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax()
389 @*/
390 PetscErrorCode  VecStrideScaleAll(Vec v,const PetscScalar *scales)
391 {
392   PetscErrorCode ierr;
393   PetscInt       i,j,n,bs;
394   PetscScalar    *x;
395 
396   PetscFunctionBegin;
397   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
398   PetscValidScalarPointer(scales,2);
399   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
400   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
401 
402   bs = v->map->bs;
403 
404   /* need to provide optimized code for each bs */
405   for (i=0; i<n; i+=bs) {
406     for (j=0; j<bs; j++) {
407       x[i+j] *= scales[j];
408     }
409   }
410   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
411   PetscFunctionReturn(0);
412 }
413 
414 #undef __FUNCT__
415 #define __FUNCT__ "VecStrideNormAll"
416 /*@
417    VecStrideNormAll - Computes the norms  subvectors of a vector defined
418    by a starting point and a stride.
419 
420    Collective on Vec
421 
422    Input Parameter:
423 +  v - the vector
424 -  ntype - type of norm, one of NORM_1, NORM_2, NORM_INFINITY
425 
426    Output Parameter:
427 .  nrm - the norms
428 
429    Notes:
430    One must call VecSetBlockSize() before this routine to set the stride
431    information, or use a vector created from a multicomponent DMDA.
432 
433    If x is the array representing the vector x then this computes the norm
434    of the array (x[start],x[start+stride],x[start+2*stride], ....)
435 
436    This is useful for computing, say the norm of the pressure variable when
437    the pressure is stored (interlaced) with other variables, say density etc.
438 
439    This will only work if the desire subvector is a stride subvector
440 
441    Level: advanced
442 
443    Concepts: norm^on stride of vector
444    Concepts: stride^norm
445 
446 .seealso: VecNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax()
447 @*/
448 PetscErrorCode  VecStrideNormAll(Vec v,NormType ntype,PetscReal nrm[])
449 {
450   PetscErrorCode ierr;
451   PetscInt       i,j,n,bs;
452   PetscScalar    *x;
453   PetscReal      tnorm[128];
454   MPI_Comm       comm;
455 
456   PetscFunctionBegin;
457   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
458   PetscValidRealPointer(nrm,3);
459   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
460   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
461   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
462 
463   bs = v->map->bs;
464   if (bs > 128) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently supports only blocksize up to 128");
465 
466   if (ntype == NORM_2) {
467     PetscScalar sum[128];
468     for (j=0; j<bs; j++) sum[j] = 0.0;
469     for (i=0; i<n; i+=bs) {
470       for (j=0; j<bs; j++) {
471         sum[j] += x[i+j]*(PetscConj(x[i+j]));
472       }
473     }
474     for (j=0; j<bs; j++) {
475       tnorm[j]  = PetscRealPart(sum[j]);
476     }
477     ierr   = MPI_Allreduce(tnorm,nrm,bs,MPIU_REAL,MPIU_SUM,comm);CHKERRQ(ierr);
478     for (j=0; j<bs; j++) {
479       nrm[j] = PetscSqrtReal(nrm[j]);
480     }
481   } else if (ntype == NORM_1) {
482     for (j=0; j<bs; j++) {
483       tnorm[j] = 0.0;
484     }
485     for (i=0; i<n; i+=bs) {
486       for (j=0; j<bs; j++) {
487         tnorm[j] += PetscAbsScalar(x[i+j]);
488       }
489     }
490     ierr   = MPI_Allreduce(tnorm,nrm,bs,MPIU_REAL,MPIU_SUM,comm);CHKERRQ(ierr);
491   } else if (ntype == NORM_INFINITY) {
492     PetscReal tmp;
493     for (j=0; j<bs; j++) {
494       tnorm[j] = 0.0;
495     }
496 
497     for (i=0; i<n; i+=bs) {
498       for (j=0; j<bs; j++) {
499         if ((tmp = PetscAbsScalar(x[i+j])) > tnorm[j]) tnorm[j] = tmp;
500         /* check special case of tmp == NaN */
501         if (tmp != tmp) {tnorm[j] = tmp; break;}
502       }
503     }
504     ierr   = MPI_Allreduce(tnorm,nrm,bs,MPIU_REAL,MPIU_MAX,comm);CHKERRQ(ierr);
505   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown norm type");
506   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
507   PetscFunctionReturn(0);
508 }
509 
510 #undef __FUNCT__
511 #define __FUNCT__ "VecStrideMaxAll"
512 /*@
513    VecStrideMaxAll - Computes the maximums of subvectors of a vector defined
514    by a starting point and a stride and optionally its location.
515 
516    Collective on Vec
517 
518    Input Parameter:
519 .  v - the vector
520 
521    Output Parameter:
522 +  index - the location where the maximum occurred (not supported, pass PETSC_NULL,
523            if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
524 -  nrm - the maximums
525 
526    Notes:
527    One must call VecSetBlockSize() before this routine to set the stride
528    information, or use a vector created from a multicomponent DMDA.
529 
530    This is useful for computing, say the maximum of the pressure variable when
531    the pressure is stored (interlaced) with other variables, e.g., density, etc.
532    This will only work if the desire subvector is a stride subvector.
533 
534    Level: advanced
535 
536    Concepts: maximum^on stride of vector
537    Concepts: stride^maximum
538 
539 .seealso: VecMax(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin()
540 @*/
541 PetscErrorCode  VecStrideMaxAll(Vec v,PetscInt idex[],PetscReal nrm[])
542 {
543   PetscErrorCode ierr;
544   PetscInt       i,j,n,bs;
545   PetscScalar    *x;
546   PetscReal      max[128],tmp;
547   MPI_Comm       comm;
548 
549   PetscFunctionBegin;
550   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
551   PetscValidRealPointer(nrm,3);
552   if (idex) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
553   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
554   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
555   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
556 
557   bs = v->map->bs;
558   if (bs > 128) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently supports only blocksize up to 128");
559 
560   if (!n) {
561     for (j=0; j<bs; j++) {
562       max[j] = PETSC_MIN_REAL;
563     }
564   } else {
565     for (j=0; j<bs; j++) {
566       max[j] = PetscRealPart(x[j]);
567     }
568     for (i=bs; i<n; i+=bs) {
569       for (j=0; j<bs; j++) {
570         if ((tmp = PetscRealPart(x[i+j])) > max[j]) { max[j] = tmp;}
571       }
572     }
573   }
574   ierr   = MPI_Allreduce(max,nrm,bs,MPIU_REAL,MPIU_MAX,comm);CHKERRQ(ierr);
575 
576   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
577   PetscFunctionReturn(0);
578 }
579 
580 #undef __FUNCT__
581 #define __FUNCT__ "VecStrideMinAll"
582 /*@
583    VecStrideMinAll - Computes the minimum of subvector of a vector defined
584    by a starting point and a stride and optionally its location.
585 
586    Collective on Vec
587 
588    Input Parameter:
589 .  v - the vector
590 
591    Output Parameter:
592 +  idex - the location where the minimum occurred (not supported, pass PETSC_NULL,
593            if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
594 -  nrm - the minimums
595 
596    Level: advanced
597 
598    Notes:
599    One must call VecSetBlockSize() before this routine to set the stride
600    information, or use a vector created from a multicomponent DMDA.
601 
602    This is useful for computing, say the minimum of the pressure variable when
603    the pressure is stored (interlaced) with other variables, e.g., density, etc.
604    This will only work if the desire subvector is a stride subvector.
605 
606    Concepts: minimum^on stride of vector
607    Concepts: stride^minimum
608 
609 .seealso: VecMin(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMax()
610 @*/
611 PetscErrorCode  VecStrideMinAll(Vec v,PetscInt idex[],PetscReal nrm[])
612 {
613   PetscErrorCode ierr;
614   PetscInt       i,n,bs,j;
615   PetscScalar    *x;
616   PetscReal      min[128],tmp;
617   MPI_Comm       comm;
618 
619   PetscFunctionBegin;
620   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
621   PetscValidRealPointer(nrm,3);
622   if (idex) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
623   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
624   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
625   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
626 
627   bs = v->map->bs;
628   if (bs > 128) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently supports only blocksize up to 128");
629 
630   if (!n) {
631     for (j=0; j<bs; j++) {
632       min[j] = PETSC_MAX_REAL;
633     }
634   } else {
635     for (j=0; j<bs; j++) {
636       min[j] = PetscRealPart(x[j]);
637     }
638     for (i=bs; i<n; i+=bs) {
639       for (j=0; j<bs; j++) {
640         if ((tmp = PetscRealPart(x[i+j])) < min[j]) { min[j] = tmp;}
641       }
642     }
643   }
644   ierr   = MPI_Allreduce(min,nrm,bs,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr);
645 
646   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
647   PetscFunctionReturn(0);
648 }
649 
650 /*----------------------------------------------------------------------------------------------*/
651 #undef __FUNCT__
652 #define __FUNCT__ "VecStrideGatherAll"
653 /*@
654    VecStrideGatherAll - Gathers all the single components from a multi-component vector into
655    separate vectors.
656 
657    Collective on Vec
658 
659    Input Parameter:
660 +  v - the vector
661 -  addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES
662 
663    Output Parameter:
664 .  s - the location where the subvectors are stored
665 
666    Notes:
667    One must call VecSetBlockSize() before this routine to set the stride
668    information, or use a vector created from a multicomponent DMDA.
669 
670    If x is the array representing the vector x then this gathers
671    the arrays (x[start],x[start+stride],x[start+2*stride], ....)
672    for start=0,1,2,...bs-1
673 
674    The parallel layout of the vector and the subvector must be the same;
675    i.e., nlocal of v = stride*(nlocal of s)
676 
677    Not optimized; could be easily
678 
679    Level: advanced
680 
681    Concepts: gather^into strided vector
682 
683 .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGather(),
684           VecStrideScatterAll()
685 @*/
686 PetscErrorCode  VecStrideGatherAll(Vec v,Vec s[],InsertMode addv)
687 {
688   PetscErrorCode ierr;
689   PetscInt       i,n,n2,bs,j,k,*bss = PETSC_NULL,nv,jj,nvc;
690   PetscScalar    *x,**y;
691 
692   PetscFunctionBegin;
693   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
694   PetscValidPointer(s,2);
695   PetscValidHeaderSpecific(*s,VEC_CLASSID,2);
696   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
697   ierr = VecGetLocalSize(s[0],&n2);CHKERRQ(ierr);
698   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
699   bs   = v->map->bs;
700   if (bs < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Input vector does not have a valid blocksize set");
701 
702   ierr = PetscMalloc2(bs,PetscReal*,&y,bs,PetscInt,&bss);CHKERRQ(ierr);
703   nv   = 0;
704   nvc  = 0;
705   for (i=0; i<bs; i++) {
706     ierr = VecGetBlockSize(s[i],&bss[i]);CHKERRQ(ierr);
707     if (bss[i] < 1) bss[i] = 1; /* if user never set it then assume 1  Re: [PETSC #8241] VecStrideGatherAll */
708     ierr = VecGetArray(s[i],&y[i]);CHKERRQ(ierr);
709     nvc += bss[i];
710     nv++;
711     if (nvc > bs)  SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subvectors in subvectors > number of vectors in main vector");
712     if (nvc == bs) break;
713   }
714 
715   n =  n/bs;
716 
717   jj = 0;
718   if (addv == INSERT_VALUES) {
719     for (j=0; j<nv; j++) {
720       for (k=0; k<bss[j]; k++) {
721         for (i=0; i<n; i++) {
722           y[j][i*bss[j] + k] = x[bs*i+jj+k];
723         }
724       }
725       jj += bss[j];
726     }
727   } else if (addv == ADD_VALUES) {
728     for (j=0; j<nv; j++) {
729       for (k=0; k<bss[j]; k++) {
730         for (i=0; i<n; i++) {
731           y[j][i*bss[j] + k] += x[bs*i+jj+k];
732         }
733       }
734       jj += bss[j];
735     }
736 #if !defined(PETSC_USE_COMPLEX)
737   } else if (addv == MAX_VALUES) {
738     for (j=0; j<nv; j++) {
739       for (k=0; k<bss[j]; k++) {
740         for (i=0; i<n; i++) {
741           y[j][i*bss[j] + k] = PetscMax(y[j][i*bss[j] + k],x[bs*i+jj+k]);
742         }
743       }
744       jj += bss[j];
745     }
746 #endif
747   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type");
748 
749   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
750   for (i=0; i<nv; i++) {
751     ierr = VecRestoreArray(s[i],&y[i]);CHKERRQ(ierr);
752   }
753 
754   ierr = PetscFree2(y,bss);CHKERRQ(ierr);
755   PetscFunctionReturn(0);
756 }
757 
758 #undef __FUNCT__
759 #define __FUNCT__ "VecStrideScatterAll"
760 /*@
761    VecStrideScatterAll - Scatters all the single components from separate vectors into
762      a multi-component vector.
763 
764    Collective on Vec
765 
766    Input Parameter:
767 +  s - the location where the subvectors are stored
768 -  addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES
769 
770    Output Parameter:
771 .  v - the multicomponent vector
772 
773    Notes:
774    One must call VecSetBlockSize() before this routine to set the stride
775    information, or use a vector created from a multicomponent DMDA.
776 
777    The parallel layout of the vector and the subvector must be the same;
778    i.e., nlocal of v = stride*(nlocal of s)
779 
780    Not optimized; could be easily
781 
782    Level: advanced
783 
784    Concepts:  scatter^into strided vector
785 
786 .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGather(),
787           VecStrideScatterAll()
788 @*/
789 PetscErrorCode  VecStrideScatterAll(Vec s[],Vec v,InsertMode addv)
790 {
791   PetscErrorCode ierr;
792   PetscInt       i,n,n2,bs,j,jj,k,*bss = PETSC_NULL,nv,nvc;
793   PetscScalar    *x,**y;
794 
795   PetscFunctionBegin;
796   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
797   PetscValidPointer(s,2);
798   PetscValidHeaderSpecific(*s,VEC_CLASSID,2);
799   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
800   ierr = VecGetLocalSize(s[0],&n2);CHKERRQ(ierr);
801   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
802   bs   = v->map->bs;
803   if (bs < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Input vector does not have a valid blocksize set");
804 
805   ierr = PetscMalloc2(bs,PetscScalar**,&y,bs,PetscInt,&bss);CHKERRQ(ierr);
806   nv   = 0;
807   nvc  = 0;
808   for (i=0; i<bs; i++) {
809     ierr = VecGetBlockSize(s[i],&bss[i]);CHKERRQ(ierr);
810     if (bss[i] < 1) bss[i] = 1; /* if user never set it then assume 1  Re: [PETSC #8241] VecStrideGatherAll */
811     ierr = VecGetArray(s[i],&y[i]);CHKERRQ(ierr);
812     nvc += bss[i];
813     nv++;
814     if (nvc > bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subvectors in subvectors > number of vectors in main vector");
815     if (nvc == bs) break;
816   }
817 
818   n =  n/bs;
819 
820   jj = 0;
821   if (addv == INSERT_VALUES) {
822     for (j=0; j<nv; j++) {
823       for (k=0; k<bss[j]; k++) {
824         for (i=0; i<n; i++) {
825           x[bs*i+jj+k] = y[j][i*bss[j] + k];
826         }
827       }
828       jj += bss[j];
829     }
830   } else if (addv == ADD_VALUES) {
831     for (j=0; j<nv; j++) {
832       for (k=0; k<bss[j]; k++) {
833         for (i=0; i<n; i++) {
834           x[bs*i+jj+k] += y[j][i*bss[j] + k];
835         }
836       }
837       jj += bss[j];
838     }
839 #if !defined(PETSC_USE_COMPLEX)
840   } else if (addv == MAX_VALUES) {
841     for (j=0; j<nv; j++) {
842       for (k=0; k<bss[j]; k++) {
843         for (i=0; i<n; i++) {
844           x[bs*i+jj+k] = PetscMax(x[bs*i+jj+k],y[j][i*bss[j] + k]);
845         }
846       }
847       jj += bss[j];
848     }
849 #endif
850   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type");
851 
852   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
853   for (i=0; i<nv; i++) {
854     ierr = VecRestoreArray(s[i],&y[i]);CHKERRQ(ierr);
855   }
856   ierr = PetscFree2(y,bss);CHKERRQ(ierr);
857   PetscFunctionReturn(0);
858 }
859 
860 #undef __FUNCT__
861 #define __FUNCT__ "VecStrideGather"
862 /*@
863    VecStrideGather - Gathers a single component from a multi-component vector into
864    another vector.
865 
866    Collective on Vec
867 
868    Input Parameter:
869 +  v - the vector
870 .  start - starting point of the subvector (defined by a stride)
871 -  addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES
872 
873    Output Parameter:
874 .  s - the location where the subvector is stored
875 
876    Notes:
877    One must call VecSetBlockSize() before this routine to set the stride
878    information, or use a vector created from a multicomponent DMDA.
879 
880    If x is the array representing the vector x then this gathers
881    the array (x[start],x[start+stride],x[start+2*stride], ....)
882 
883    The parallel layout of the vector and the subvector must be the same;
884    i.e., nlocal of v = stride*(nlocal of s)
885 
886    Not optimized; could be easily
887 
888    Level: advanced
889 
890    Concepts: gather^into strided vector
891 
892 .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGatherAll(),
893           VecStrideScatterAll()
894 @*/
895 PetscErrorCode  VecStrideGather(Vec v,PetscInt start,Vec s,InsertMode addv)
896 {
897   PetscErrorCode ierr;
898 
899   PetscFunctionBegin;
900   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
901   PetscValidHeaderSpecific(s,VEC_CLASSID,3);
902   if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start);
903   if (start >= v->map->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Start of stride subvector (%D) is too large for stride\n Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,v->map->bs);
904   if (!v->ops->stridegather) SETERRQ(((PetscObject)s)->comm,PETSC_ERR_SUP,"Not implemented for this Vec class");
905   ierr = (*v->ops->stridegather)(v,start,s,addv);CHKERRQ(ierr);
906   PetscFunctionReturn(0);
907 }
908 
909 #undef __FUNCT__
910 #define __FUNCT__ "VecStrideScatter"
911 /*@
912    VecStrideScatter - Scatters a single component from a vector into a multi-component vector.
913 
914    Collective on Vec
915 
916    Input Parameter:
917 +  s - the single-component vector
918 .  start - starting point of the subvector (defined by a stride)
919 -  addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES
920 
921    Output Parameter:
922 .  v - the location where the subvector is scattered (the multi-component vector)
923 
924    Notes:
925    One must call VecSetBlockSize() on the multi-component vector before this
926    routine to set the stride  information, or use a vector created from a multicomponent DMDA.
927 
928    The parallel layout of the vector and the subvector must be the same;
929    i.e., nlocal of v = stride*(nlocal of s)
930 
931    Not optimized; could be easily
932 
933    Level: advanced
934 
935    Concepts: scatter^into strided vector
936 
937 .seealso: VecStrideNorm(), VecStrideGather(), VecStrideMin(), VecStrideMax(), VecStrideGatherAll(),
938           VecStrideScatterAll()
939 @*/
940 PetscErrorCode  VecStrideScatter(Vec s,PetscInt start,Vec v,InsertMode addv)
941 {
942   PetscErrorCode ierr;
943 
944   PetscFunctionBegin;
945   PetscValidHeaderSpecific(s,VEC_CLASSID,1);
946   PetscValidHeaderSpecific(v,VEC_CLASSID,3);
947   if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start);
948   if (start >= v->map->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Start of stride subvector (%D) is too large for stride\n Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,v->map->bs);
949   if (!v->ops->stridescatter) SETERRQ(((PetscObject)s)->comm,PETSC_ERR_SUP,"Not implemented for this Vec class");
950   ierr = (*v->ops->stridescatter)(s,start,v,addv);CHKERRQ(ierr);
951   PetscFunctionReturn(0);
952 }
953 
954 #undef __FUNCT__
955 #define __FUNCT__ "VecStrideGather_Default"
956 PetscErrorCode  VecStrideGather_Default(Vec v,PetscInt start,Vec s,InsertMode addv)
957 {
958   PetscErrorCode ierr;
959   PetscInt       i,n,bs,ns;
960   PetscScalar    *x,*y;
961 
962   PetscFunctionBegin;
963   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
964   ierr = VecGetLocalSize(s,&ns);CHKERRQ(ierr);
965   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
966   ierr = VecGetArray(s,&y);CHKERRQ(ierr);
967 
968   bs = v->map->bs;
969   if (n != ns*bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Subvector length * blocksize %D not correct for gather from original vector %D",ns*bs,n);
970   x += start;
971   n  =  n/bs;
972 
973   if (addv == INSERT_VALUES) {
974     for (i=0; i<n; i++) {
975       y[i] = x[bs*i];
976     }
977   } else if (addv == ADD_VALUES) {
978     for (i=0; i<n; i++) {
979       y[i] += x[bs*i];
980     }
981 #if !defined(PETSC_USE_COMPLEX)
982   } else if (addv == MAX_VALUES) {
983     for (i=0; i<n; i++) {
984       y[i] = PetscMax(y[i],x[bs*i]);
985     }
986 #endif
987   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type");
988 
989   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
990   ierr = VecRestoreArray(s,&y);CHKERRQ(ierr);
991   PetscFunctionReturn(0);
992 }
993 
994 #undef __FUNCT__
995 #define __FUNCT__ "VecStrideScatter_Default"
996 PetscErrorCode  VecStrideScatter_Default(Vec s,PetscInt start,Vec v,InsertMode addv)
997 {
998   PetscErrorCode ierr;
999   PetscInt       i,n,bs,ns;
1000   PetscScalar    *x,*y;
1001 
1002   PetscFunctionBegin;
1003   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1004   ierr = VecGetLocalSize(s,&ns);CHKERRQ(ierr);
1005   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1006   ierr = VecGetArray(s,&y);CHKERRQ(ierr);
1007 
1008   bs = v->map->bs;
1009   if (n != ns*bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Subvector length * blocksize %D not correct for scatter to multicomponent vector %D",ns*bs,n);
1010   x += start;
1011   n  =  n/bs;
1012 
1013   if (addv == INSERT_VALUES) {
1014     for (i=0; i<n; i++) {
1015       x[bs*i] = y[i];
1016     }
1017   } else if (addv == ADD_VALUES) {
1018     for (i=0; i<n; i++) {
1019       x[bs*i] += y[i];
1020     }
1021 #if !defined(PETSC_USE_COMPLEX)
1022   } else if (addv == MAX_VALUES) {
1023     for (i=0; i<n; i++) {
1024       x[bs*i] = PetscMax(y[i],x[bs*i]);
1025     }
1026 #endif
1027   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type");
1028 
1029   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1030   ierr = VecRestoreArray(s,&y);CHKERRQ(ierr);
1031   PetscFunctionReturn(0);
1032 }
1033 
1034 #undef __FUNCT__
1035 #define __FUNCT__ "VecReciprocal_Default"
1036 PetscErrorCode VecReciprocal_Default(Vec v)
1037 {
1038   PetscErrorCode ierr;
1039   PetscInt       i,n;
1040   PetscScalar    *x;
1041 
1042   PetscFunctionBegin;
1043   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1044   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1045   for (i=0; i<n; i++) {
1046     if (x[i] != (PetscScalar)0.0) x[i] = (PetscScalar)1.0/x[i];
1047   }
1048   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1049   PetscFunctionReturn(0);
1050 }
1051 
1052 #undef __FUNCT__
1053 #define __FUNCT__ "VecExp"
1054 /*@
1055   VecExp - Replaces each component of a vector by e^x_i
1056 
1057   Not collective
1058 
1059   Input Parameter:
1060 . v - The vector
1061 
1062   Output Parameter:
1063 . v - The vector of exponents
1064 
1065   Level: beginner
1066 
1067 .seealso:  VecLog(), VecAbs(), VecSqrtAbs(), VecReciprocal()
1068 
1069 .keywords: vector, sqrt, square root
1070 @*/
1071 PetscErrorCode  VecExp(Vec v)
1072 {
1073   PetscScalar    *x;
1074   PetscInt       i, n;
1075   PetscErrorCode ierr;
1076 
1077   PetscFunctionBegin;
1078   PetscValidHeaderSpecific(v, VEC_CLASSID,1);
1079   if (v->ops->exp) {
1080     ierr = (*v->ops->exp)(v);CHKERRQ(ierr);
1081   } else {
1082     ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
1083     ierr = VecGetArray(v, &x);CHKERRQ(ierr);
1084     for (i = 0; i < n; i++) {
1085       x[i] = PetscExpScalar(x[i]);
1086     }
1087     ierr = VecRestoreArray(v, &x);CHKERRQ(ierr);
1088   }
1089   PetscFunctionReturn(0);
1090 }
1091 
1092 #undef __FUNCT__
1093 #define __FUNCT__ "VecLog"
1094 /*@
1095   VecLog - Replaces each component of a vector by log(x_i), the natural logarithm
1096 
1097   Not collective
1098 
1099   Input Parameter:
1100 . v - The vector
1101 
1102   Output Parameter:
1103 . v - The vector of logs
1104 
1105   Level: beginner
1106 
1107 .seealso:  VecExp(), VecAbs(), VecSqrtAbs(), VecReciprocal()
1108 
1109 .keywords: vector, sqrt, square root
1110 @*/
1111 PetscErrorCode  VecLog(Vec v)
1112 {
1113   PetscScalar    *x;
1114   PetscInt       i, n;
1115   PetscErrorCode ierr;
1116 
1117   PetscFunctionBegin;
1118   PetscValidHeaderSpecific(v, VEC_CLASSID,1);
1119   if (v->ops->log) {
1120     ierr = (*v->ops->log)(v);CHKERRQ(ierr);
1121   } else {
1122     ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
1123     ierr = VecGetArray(v, &x);CHKERRQ(ierr);
1124     for (i = 0; i < n; i++) {
1125       x[i] = PetscLogScalar(x[i]);
1126     }
1127     ierr = VecRestoreArray(v, &x);CHKERRQ(ierr);
1128   }
1129   PetscFunctionReturn(0);
1130 }
1131 
1132 #undef __FUNCT__
1133 #define __FUNCT__ "VecSqrtAbs"
1134 /*@
1135   VecSqrtAbs - Replaces each component of a vector by the square root of its magnitude.
1136 
1137   Not collective
1138 
1139   Input Parameter:
1140 . v - The vector
1141 
1142   Output Parameter:
1143 . v - The vector square root
1144 
1145   Level: beginner
1146 
1147   Note: The actual function is sqrt(|x_i|)
1148 
1149 .seealso: VecLog(), VecExp(), VecReciprocal(), VecAbs()
1150 
1151 .keywords: vector, sqrt, square root
1152 @*/
1153 PetscErrorCode  VecSqrtAbs(Vec v)
1154 {
1155   PetscScalar    *x;
1156   PetscInt       i, n;
1157   PetscErrorCode ierr;
1158 
1159   PetscFunctionBegin;
1160   PetscValidHeaderSpecific(v, VEC_CLASSID,1);
1161   if (v->ops->sqrt) {
1162     ierr = (*v->ops->sqrt)(v);CHKERRQ(ierr);
1163   } else {
1164     ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
1165     ierr = VecGetArray(v, &x);CHKERRQ(ierr);
1166     for (i = 0; i < n; i++) {
1167       x[i] = PetscSqrtReal(PetscAbsScalar(x[i]));
1168     }
1169     ierr = VecRestoreArray(v, &x);CHKERRQ(ierr);
1170   }
1171   PetscFunctionReturn(0);
1172 }
1173 
1174 #undef __FUNCT__
1175 #define __FUNCT__ "VecDotNorm2"
1176 /*@
1177   VecDotNorm2 - computes the inner product of two vectors and the 2-norm squared of the second vector
1178 
1179   Collective on Vec
1180 
1181   Input Parameter:
1182 + s - first vector
1183 - t - second vector
1184 
1185   Output Parameter:
1186 + dp - s'conj(t)
1187 - nm - t'conj(t)
1188 
1189   Level: advanced
1190 
1191   Notes: conj(x) is the complex conjugate of x when x is complex
1192 
1193 
1194 .seealso:   VecDot(), VecNorm(), VecDotBegin(), VecNormBegin(), VecDotEnd(), VecNormEnd()
1195 
1196 .keywords: vector, sqrt, square root
1197 @*/
1198 PetscErrorCode  VecDotNorm2(Vec s,Vec t,PetscScalar *dp, PetscReal *nm)
1199 {
1200   PetscScalar    *sx, *tx, dpx = 0.0, nmx = 0.0,work[2],sum[2];
1201   PetscInt       i, n;
1202   PetscErrorCode ierr;
1203 
1204   PetscFunctionBegin;
1205   PetscValidHeaderSpecific(s, VEC_CLASSID,1);
1206   PetscValidHeaderSpecific(t, VEC_CLASSID,2);
1207   PetscValidScalarPointer(dp,3);
1208   PetscValidScalarPointer(nm,4);
1209   PetscValidType(s,1);
1210   PetscValidType(t,2);
1211   PetscCheckSameTypeAndComm(s,1,t,2);
1212   if (s->map->N != t->map->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
1213   if (s->map->n != t->map->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");
1214 
1215   ierr = PetscLogEventBarrierBegin(VEC_DotNormBarrier,s,t,0,0,((PetscObject)s)->comm);CHKERRQ(ierr);
1216   if (s->ops->dotnorm2) {
1217     ierr = (*s->ops->dotnorm2)(s,t,dp,&dpx);CHKERRQ(ierr);
1218     *nm  = PetscRealPart(dpx);CHKERRQ(ierr);
1219   } else {
1220     ierr = VecGetLocalSize(s, &n);CHKERRQ(ierr);
1221     ierr = VecGetArray(s, &sx);CHKERRQ(ierr);
1222     ierr = VecGetArray(t, &tx);CHKERRQ(ierr);
1223 
1224     for (i = 0; i<n; i++) {
1225       dpx += sx[i]*PetscConj(tx[i]);
1226       nmx += tx[i]*PetscConj(tx[i]);
1227     }
1228     work[0] = dpx;
1229     work[1] = nmx;
1230     ierr = MPI_Allreduce(work,sum,2,MPIU_SCALAR,MPIU_SUM,((PetscObject)s)->comm);CHKERRQ(ierr);
1231     *dp  = sum[0];
1232     *nm  = PetscRealPart(sum[1]);
1233 
1234     ierr = VecRestoreArray(t, &tx);CHKERRQ(ierr);
1235     ierr = VecRestoreArray(s, &sx);CHKERRQ(ierr);
1236     ierr = PetscLogFlops(4.0*n);CHKERRQ(ierr);
1237   }
1238   ierr = PetscLogEventBarrierEnd(VEC_DotNormBarrier,s,t,0,0,((PetscObject)s)->comm);CHKERRQ(ierr);
1239   PetscFunctionReturn(0);
1240 }
1241 
1242 #undef __FUNCT__
1243 #define __FUNCT__ "VecSum"
1244 /*@
1245    VecSum - Computes the sum of all the components of a vector.
1246 
1247    Collective on Vec
1248 
1249    Input Parameter:
1250 .  v - the vector
1251 
1252    Output Parameter:
1253 .  sum - the result
1254 
1255    Level: beginner
1256 
1257    Concepts: sum^of vector entries
1258 
1259 .seealso: VecNorm()
1260 @*/
1261 PetscErrorCode  VecSum(Vec v,PetscScalar *sum)
1262 {
1263   PetscErrorCode ierr;
1264   PetscInt       i,n;
1265   PetscScalar    *x,lsum = 0.0;
1266 
1267   PetscFunctionBegin;
1268   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
1269   PetscValidScalarPointer(sum,2);
1270   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1271   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1272   for (i=0; i<n; i++) {
1273     lsum += x[i];
1274   }
1275   ierr = MPI_Allreduce(&lsum,sum,1,MPIU_SCALAR,MPIU_SUM,((PetscObject)v)->comm);CHKERRQ(ierr);
1276   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1277   PetscFunctionReturn(0);
1278 }
1279 
1280 #undef __FUNCT__
1281 #define __FUNCT__ "VecShift"
1282 /*@
1283    VecShift - Shifts all of the components of a vector by computing
1284    x[i] = x[i] + shift.
1285 
1286    Logically Collective on Vec
1287 
1288    Input Parameters:
1289 +  v - the vector
1290 -  shift - the shift
1291 
1292    Output Parameter:
1293 .  v - the shifted vector
1294 
1295    Level: intermediate
1296 
1297    Concepts: vector^adding constant
1298 
1299 @*/
1300 PetscErrorCode  VecShift(Vec v,PetscScalar shift)
1301 {
1302   PetscErrorCode ierr;
1303   PetscInt       i,n;
1304   PetscScalar    *x;
1305 
1306   PetscFunctionBegin;
1307   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
1308   PetscValidLogicalCollectiveScalar(v,shift,2);
1309   if (v->ops->shift) {
1310     ierr = (*v->ops->shift)(v);CHKERRQ(ierr);
1311   } else {
1312     ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1313     ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1314     for (i=0; i<n; i++) {
1315       x[i] += shift;
1316     }
1317     ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1318   }
1319   PetscFunctionReturn(0);
1320 }
1321 
1322 #undef __FUNCT__
1323 #define __FUNCT__ "VecAbs"
1324 /*@
1325    VecAbs - Replaces every element in a vector with its absolute value.
1326 
1327    Logically Collective on Vec
1328 
1329    Input Parameters:
1330 .  v - the vector
1331 
1332    Level: intermediate
1333 
1334    Concepts: vector^absolute value
1335 
1336 @*/
1337 PetscErrorCode  VecAbs(Vec v)
1338 {
1339   PetscErrorCode ierr;
1340   PetscInt       i,n;
1341   PetscScalar    *x;
1342 
1343   PetscFunctionBegin;
1344   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
1345   if (v->ops->abs) {
1346     ierr = (*v->ops->abs)(v);CHKERRQ(ierr);
1347   } else {
1348     ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1349     ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1350     for (i=0; i<n; i++) {
1351       x[i] = PetscAbsScalar(x[i]);
1352     }
1353     ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1354   }
1355   PetscFunctionReturn(0);
1356 }
1357 
1358 #undef __FUNCT__
1359 #define __FUNCT__ "VecPermute"
1360 /*@
1361   VecPermute - Permutes a vector in place using the given ordering.
1362 
1363   Input Parameters:
1364 + vec   - The vector
1365 . order - The ordering
1366 - inv   - The flag for inverting the permutation
1367 
1368   Level: beginner
1369 
1370   Note: This function does not yet support parallel Index Sets with non-local permutations
1371 
1372 .seealso: MatPermute()
1373 .keywords: vec, permute
1374 @*/
1375 PetscErrorCode  VecPermute(Vec x, IS row, PetscBool inv)
1376 {
1377   PetscScalar    *array, *newArray;
1378   const PetscInt *idx;
1379   PetscInt       i,rstart,rend;
1380   PetscErrorCode ierr;
1381 
1382   PetscFunctionBegin;
1383   ierr = VecGetOwnershipRange(x,&rstart,&rend);CHKERRQ(ierr);
1384   ierr = ISGetIndices(row, &idx);CHKERRQ(ierr);
1385   ierr = VecGetArray(x, &array);CHKERRQ(ierr);
1386   ierr = PetscMalloc(x->map->n*sizeof(PetscScalar), &newArray);CHKERRQ(ierr);
1387 #if defined(PETSC_USE_DEBUG)
1388   for (i = 0; i < x->map->n; i++) {
1389     if ((idx[i] < rstart) || (idx[i] >= rend)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT, "Permutation index %D is out of bounds: %D", i, idx[i]);
1390   }
1391 #endif
1392   if (!inv) {
1393     for (i = 0; i < x->map->n; i++) newArray[i] = array[idx[i]-rstart];
1394   } else {
1395     for (i = 0; i < x->map->n; i++) newArray[idx[i]-rstart] = array[i];
1396   }
1397   ierr = VecRestoreArray(x, &array);CHKERRQ(ierr);
1398   ierr = ISRestoreIndices(row, &idx);CHKERRQ(ierr);
1399   ierr = VecReplaceArray(x, newArray);CHKERRQ(ierr);
1400   PetscFunctionReturn(0);
1401 }
1402 
1403 #undef __FUNCT__
1404 #define __FUNCT__ "VecEqual"
1405 /*@
1406    VecEqual - Compares two vectors.
1407 
1408    Collective on Vec
1409 
1410    Input Parameters:
1411 +  vec1 - the first vector
1412 -  vec2 - the second vector
1413 
1414    Output Parameter:
1415 .  flg - PETSC_TRUE if the vectors are equal; PETSC_FALSE otherwise.
1416 
1417    Level: intermediate
1418 
1419    Concepts: equal^two vectors
1420    Concepts: vector^equality
1421 
1422 @*/
1423 PetscErrorCode  VecEqual(Vec vec1,Vec vec2,PetscBool  *flg)
1424 {
1425   PetscScalar    *v1,*v2;
1426   PetscErrorCode ierr;
1427   PetscInt       n1,n2,N1,N2;
1428   PetscInt       state1,state2;
1429   PetscBool      flg1;
1430 
1431   PetscFunctionBegin;
1432   PetscValidHeaderSpecific(vec1,VEC_CLASSID,1);
1433   PetscValidHeaderSpecific(vec2,VEC_CLASSID,2);
1434   PetscValidPointer(flg,3);
1435   if (vec1 == vec2) {
1436     *flg = PETSC_TRUE;
1437   } else {
1438     ierr = VecGetSize(vec1,&N1);CHKERRQ(ierr);
1439     ierr = VecGetSize(vec2,&N2);CHKERRQ(ierr);
1440     if (N1 != N2) {
1441       flg1 = PETSC_FALSE;
1442     } else {
1443       ierr = VecGetLocalSize(vec1,&n1);CHKERRQ(ierr);
1444       ierr = VecGetLocalSize(vec2,&n2);CHKERRQ(ierr);
1445       if (n1 != n2) {
1446         flg1 = PETSC_FALSE;
1447       } else {
1448         ierr = PetscObjectStateQuery((PetscObject) vec1,&state1);CHKERRQ(ierr);
1449         ierr = PetscObjectStateQuery((PetscObject) vec2,&state2);CHKERRQ(ierr);
1450         ierr = VecGetArray(vec1,&v1);CHKERRQ(ierr);
1451         ierr = VecGetArray(vec2,&v2);CHKERRQ(ierr);
1452 #if defined(PETSC_USE_COMPLEX)
1453         {
1454           PetscInt k;
1455           flg1 = PETSC_TRUE;
1456           for (k=0; k<n1; k++) {
1457             if (PetscRealPart(v1[k]) != PetscRealPart(v2[k]) || PetscImaginaryPart(v1[k]) != PetscImaginaryPart(v2[k])) {
1458               flg1 = PETSC_FALSE;
1459               break;
1460             }
1461           }
1462         }
1463 #else
1464         ierr = PetscMemcmp(v1,v2,n1*sizeof(PetscScalar),&flg1);CHKERRQ(ierr);
1465 #endif
1466         ierr = VecRestoreArray(vec1,&v1);CHKERRQ(ierr);
1467         ierr = VecRestoreArray(vec2,&v2);CHKERRQ(ierr);
1468         ierr = PetscObjectSetState((PetscObject) vec1,state1);CHKERRQ(ierr);
1469         ierr = PetscObjectSetState((PetscObject) vec2,state2);CHKERRQ(ierr);
1470       }
1471     }
1472     /* combine results from all processors */
1473     ierr = MPI_Allreduce(&flg1,flg,1,MPIU_BOOL,MPI_MIN,((PetscObject)vec1)->comm);CHKERRQ(ierr);
1474   }
1475   PetscFunctionReturn(0);
1476 }
1477 
1478 
1479 
1480