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