xref: /petsc/src/vec/vec/utils/vinv.c (revision 3e08d2bebf9958b2a7617e5d4703c5f9c1d28ec4)
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 #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   PetscValidRealPointer(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   PetscValidRealPointer(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] = PetscSqrtReal(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   PetscValidRealPointer(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   PetscValidRealPointer(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 Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,v->map->bs);
938   if (!v->ops->stridegather) SETERRQ(((PetscObject)s)->comm,PETSC_ERR_SUP,"Not implemented for this Vec class");
939   ierr = (*v->ops->stridegather)(v,start,s,addv);CHKERRQ(ierr);
940   PetscFunctionReturn(0);
941 }
942 
943 #undef __FUNCT__
944 #define __FUNCT__ "VecStrideScatter"
945 /*@
946    VecStrideScatter - Scatters a single component from a vector into a multi-component vector.
947 
948    Collective on Vec
949 
950    Input Parameter:
951 +  s - the single-component vector
952 .  start - starting point of the subvector (defined by a stride)
953 -  addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES
954 
955    Output Parameter:
956 .  v - the location where the subvector is scattered (the multi-component vector)
957 
958    Notes:
959    One must call VecSetBlockSize() on the multi-component vector before this
960    routine to set the stride  information, or use a vector created from a multicomponent DMDA.
961 
962    The parallel layout of the vector and the subvector must be the same;
963    i.e., nlocal of v = stride*(nlocal of s)
964 
965    Not optimized; could be easily
966 
967    Level: advanced
968 
969    Concepts: scatter^into strided vector
970 
971 .seealso: VecStrideNorm(), VecStrideGather(), VecStrideMin(), VecStrideMax(), VecStrideGatherAll(),
972           VecStrideScatterAll()
973 @*/
974 PetscErrorCode  VecStrideScatter(Vec s,PetscInt start,Vec v,InsertMode addv)
975 {
976   PetscErrorCode ierr;
977 
978   PetscFunctionBegin;
979   PetscValidHeaderSpecific(s,VEC_CLASSID,1);
980   PetscValidHeaderSpecific(v,VEC_CLASSID,3);
981   if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start);
982   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);
983   if (!v->ops->stridescatter) SETERRQ(((PetscObject)s)->comm,PETSC_ERR_SUP,"Not implemented for this Vec class");
984   ierr = (*v->ops->stridescatter)(s,start,v,addv);CHKERRQ(ierr);
985   PetscFunctionReturn(0);
986 }
987 
988 #undef __FUNCT__
989 #define __FUNCT__ "VecStrideGather_Default"
990 PetscErrorCode  VecStrideGather_Default(Vec v,PetscInt start,Vec s,InsertMode addv)
991 {
992   PetscErrorCode ierr;
993   PetscInt       i,n,bs,ns;
994   PetscScalar    *x,*y;
995 
996   PetscFunctionBegin;
997   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
998   ierr = VecGetLocalSize(s,&ns);CHKERRQ(ierr);
999   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1000   ierr = VecGetArray(s,&y);CHKERRQ(ierr);
1001 
1002   bs   = v->map->bs;
1003   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);
1004   x += start;
1005   n =  n/bs;
1006 
1007   if (addv == INSERT_VALUES) {
1008     for (i=0; i<n; i++) {
1009       y[i] = x[bs*i];
1010     }
1011   } else if (addv == ADD_VALUES) {
1012     for (i=0; i<n; i++) {
1013       y[i] += x[bs*i];
1014     }
1015 #if !defined(PETSC_USE_COMPLEX)
1016   } else if (addv == MAX_VALUES) {
1017     for (i=0; i<n; i++) {
1018       y[i] = PetscMax(y[i],x[bs*i]);
1019     }
1020 #endif
1021   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type");
1022 
1023   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1024   ierr = VecRestoreArray(s,&y);CHKERRQ(ierr);
1025   PetscFunctionReturn(0);
1026 }
1027 
1028 #undef __FUNCT__
1029 #define __FUNCT__ "VecStrideScatter_Default"
1030 PetscErrorCode  VecStrideScatter_Default(Vec s,PetscInt start,Vec v,InsertMode addv)
1031 {
1032   PetscErrorCode ierr;
1033   PetscInt       i,n,bs,ns;
1034   PetscScalar    *x,*y;
1035 
1036   PetscFunctionBegin;
1037   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1038   ierr = VecGetLocalSize(s,&ns);CHKERRQ(ierr);
1039   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1040   ierr = VecGetArray(s,&y);CHKERRQ(ierr);
1041 
1042   bs   = v->map->bs;
1043   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);
1044   x += start;
1045   n =  n/bs;
1046 
1047   if (addv == INSERT_VALUES) {
1048     for (i=0; i<n; i++) {
1049       x[bs*i] = y[i];
1050     }
1051   } else if (addv == ADD_VALUES) {
1052     for (i=0; i<n; i++) {
1053       x[bs*i] += y[i];
1054     }
1055 #if !defined(PETSC_USE_COMPLEX)
1056   } else if (addv == MAX_VALUES) {
1057     for (i=0; i<n; i++) {
1058       x[bs*i] = PetscMax(y[i],x[bs*i]);
1059     }
1060 #endif
1061   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type");
1062 
1063   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1064   ierr = VecRestoreArray(s,&y);CHKERRQ(ierr);
1065   PetscFunctionReturn(0);
1066 }
1067 
1068 #undef __FUNCT__
1069 #define __FUNCT__ "VecReciprocal_Default"
1070 PetscErrorCode VecReciprocal_Default(Vec v)
1071 {
1072   PetscErrorCode ierr;
1073   PetscInt       i,n;
1074   PetscScalar    *x;
1075 
1076   PetscFunctionBegin;
1077   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1078   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1079   for (i=0; i<n; i++) {
1080     if (x[i] != (PetscScalar)0.0) x[i] = (PetscScalar)1.0/x[i];
1081   }
1082   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1083   PetscFunctionReturn(0);
1084 }
1085 
1086 #undef __FUNCT__
1087 #define __FUNCT__ "VecExp"
1088 /*@
1089   VecExp - Replaces each component of a vector by e^x_i
1090 
1091   Not collective
1092 
1093   Input Parameter:
1094 . v - The vector
1095 
1096   Output Parameter:
1097 . v - The vector of exponents
1098 
1099   Level: beginner
1100 
1101 .seealso:  VecLog(), VecAbs(), VecSqrtAbs(), VecReciprocal()
1102 
1103 .keywords: vector, sqrt, square root
1104 @*/
1105 PetscErrorCode  VecExp(Vec v)
1106 {
1107   PetscScalar    *x;
1108   PetscInt       i, n;
1109   PetscErrorCode ierr;
1110 
1111   PetscFunctionBegin;
1112   PetscValidHeaderSpecific(v, VEC_CLASSID,1);
1113   if (v->ops->exp) {
1114     ierr = (*v->ops->exp)(v);CHKERRQ(ierr);
1115   } else {
1116     ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
1117     ierr = VecGetArray(v, &x);CHKERRQ(ierr);
1118     for (i = 0; i < n; i++) {
1119       x[i] = PetscExpScalar(x[i]);
1120     }
1121     ierr = VecRestoreArray(v, &x);CHKERRQ(ierr);
1122   }
1123   PetscFunctionReturn(0);
1124 }
1125 
1126 #undef __FUNCT__
1127 #define __FUNCT__ "VecLog"
1128 /*@
1129   VecLog - Replaces each component of a vector by log(x_i), the natural logarithm
1130 
1131   Not collective
1132 
1133   Input Parameter:
1134 . v - The vector
1135 
1136   Output Parameter:
1137 . v - The vector of logs
1138 
1139   Level: beginner
1140 
1141 .seealso:  VecExp(), VecAbs(), VecSqrtAbs(), VecReciprocal()
1142 
1143 .keywords: vector, sqrt, square root
1144 @*/
1145 PetscErrorCode  VecLog(Vec v)
1146 {
1147   PetscScalar    *x;
1148   PetscInt       i, n;
1149   PetscErrorCode ierr;
1150 
1151   PetscFunctionBegin;
1152   PetscValidHeaderSpecific(v, VEC_CLASSID,1);
1153   if (v->ops->log) {
1154     ierr = (*v->ops->log)(v);CHKERRQ(ierr);
1155   } else {
1156     ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
1157     ierr = VecGetArray(v, &x);CHKERRQ(ierr);
1158     for (i = 0; i < n; i++) {
1159       x[i] = PetscLogScalar(x[i]);
1160     }
1161     ierr = VecRestoreArray(v, &x);CHKERRQ(ierr);
1162   }
1163   PetscFunctionReturn(0);
1164 }
1165 
1166 #undef __FUNCT__
1167 #define __FUNCT__ "VecSqrtAbs"
1168 /*@
1169   VecSqrtAbs - Replaces each component of a vector by the square root of its magnitude.
1170 
1171   Not collective
1172 
1173   Input Parameter:
1174 . v - The vector
1175 
1176   Output Parameter:
1177 . v - The vector square root
1178 
1179   Level: beginner
1180 
1181   Note: The actual function is sqrt(|x_i|)
1182 
1183 .seealso: VecLog(), VecExp(), VecReciprocal(), VecAbs()
1184 
1185 .keywords: vector, sqrt, square root
1186 @*/
1187 PetscErrorCode  VecSqrtAbs(Vec v)
1188 {
1189   PetscScalar    *x;
1190   PetscInt       i, n;
1191   PetscErrorCode ierr;
1192 
1193   PetscFunctionBegin;
1194   PetscValidHeaderSpecific(v, VEC_CLASSID,1);
1195   if (v->ops->sqrt) {
1196     ierr = (*v->ops->sqrt)(v);CHKERRQ(ierr);
1197   } else {
1198     ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
1199     ierr = VecGetArray(v, &x);CHKERRQ(ierr);
1200     for (i = 0; i < n; i++) {
1201       x[i] = PetscSqrtReal(PetscAbsScalar(x[i]));
1202     }
1203     ierr = VecRestoreArray(v, &x);CHKERRQ(ierr);
1204   }
1205   PetscFunctionReturn(0);
1206 }
1207 
1208 #undef __FUNCT__
1209 #define __FUNCT__ "VecDotNorm2"
1210 /*@
1211   VecDotNorm2 - computes the inner product of two vectors and the 2-norm squared of the second vector
1212 
1213   Collective on Vec
1214 
1215   Input Parameter:
1216 + s - first vector
1217 - t - second vector
1218 
1219   Output Parameter:
1220 + dp - s't
1221 - nm - t't
1222 
1223   Level: advanced
1224 
1225   Developer Notes: Even though the second return argument is a norm and hence could be a PetscReal value it is returned as PetscScalar
1226 
1227 .seealso:   VecDot(), VecNorm(), VecDotBegin(), VecNormBegin(), VecDotEnd(), VecNormEnd()
1228 
1229 .keywords: vector, sqrt, square root
1230 @*/
1231 PetscErrorCode  VecDotNorm2(Vec s,Vec t,PetscScalar *dp, PetscScalar *nm)
1232 {
1233   PetscScalar    *sx, *tx, dpx = 0.0, nmx = 0.0,work[2],sum[2];
1234   PetscInt       i, n;
1235   PetscErrorCode ierr;
1236 
1237   PetscFunctionBegin;
1238   PetscValidHeaderSpecific(s, VEC_CLASSID,1);
1239   PetscValidHeaderSpecific(t, VEC_CLASSID,2);
1240   PetscValidScalarPointer(dp,3);
1241   PetscValidScalarPointer(nm,4);
1242   PetscValidType(s,1);
1243   PetscValidType(t,2);
1244   PetscCheckSameTypeAndComm(s,1,t,2);
1245   if (s->map->N != t->map->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
1246   if (s->map->n != t->map->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");
1247 
1248   ierr = PetscLogEventBarrierBegin(VEC_DotNormBarrier,s,t,0,0,((PetscObject)s)->comm);CHKERRQ(ierr);
1249   if (s->ops->dotnorm2) {
1250     ierr = (*s->ops->dotnorm2)(s,t,dp,nm);CHKERRQ(ierr);
1251   } else {
1252     PetscReal tmp = 0.0;
1253     ierr = VecGetLocalSize(s, &n);CHKERRQ(ierr);
1254     ierr = VecGetArray(s, &sx);CHKERRQ(ierr);
1255     ierr = VecGetArray(t, &tx);CHKERRQ(ierr);
1256 
1257     ierr = VecNorm_MPI(s,NORM_2, &tmp);
1258     ierr = VecNorm_MPI(t,NORM_2, &tmp);
1259 
1260     for (i = 0; i<n; i++) {
1261       dpx += sx[i]*PetscConj(tx[i]);
1262       nmx += tx[i]*PetscConj(tx[i]);
1263     }
1264     work[0] = dpx;
1265     work[1] = nmx;
1266     ierr = MPI_Allreduce(&work,&sum,2,MPIU_SCALAR,MPIU_SUM,((PetscObject)s)->comm);CHKERRQ(ierr);
1267     *dp  = sum[0];
1268     *nm  = sum[1];
1269 
1270     ierr = VecRestoreArray(t, &tx);CHKERRQ(ierr);
1271     ierr = VecRestoreArray(s, &sx);CHKERRQ(ierr);
1272     ierr = PetscLogFlops(4.0*n);CHKERRQ(ierr);
1273   }
1274   ierr = PetscLogEventBarrierEnd(VEC_DotNormBarrier,s,t,0,0,((PetscObject)s)->comm);CHKERRQ(ierr);
1275   PetscFunctionReturn(0);
1276 }
1277 
1278 #undef __FUNCT__
1279 #define __FUNCT__ "VecSum"
1280 /*@
1281    VecSum - Computes the sum of all the components of a vector.
1282 
1283    Collective on Vec
1284 
1285    Input Parameter:
1286 .  v - the vector
1287 
1288    Output Parameter:
1289 .  sum - the result
1290 
1291    Level: beginner
1292 
1293    Concepts: sum^of vector entries
1294 
1295 .seealso: VecNorm()
1296 @*/
1297 PetscErrorCode  VecSum(Vec v,PetscScalar *sum)
1298 {
1299   PetscErrorCode ierr;
1300   PetscInt       i,n;
1301   PetscScalar    *x,lsum = 0.0;
1302 
1303   PetscFunctionBegin;
1304   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
1305   PetscValidScalarPointer(sum,2);
1306   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1307   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1308   for (i=0; i<n; i++) {
1309     lsum += x[i];
1310   }
1311   ierr = MPI_Allreduce(&lsum,sum,1,MPIU_SCALAR,MPIU_SUM,((PetscObject)v)->comm);CHKERRQ(ierr);
1312   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1313   PetscFunctionReturn(0);
1314 }
1315 
1316 #undef __FUNCT__
1317 #define __FUNCT__ "VecShift"
1318 /*@
1319    VecShift - Shifts all of the components of a vector by computing
1320    x[i] = x[i] + shift.
1321 
1322    Logically Collective on Vec
1323 
1324    Input Parameters:
1325 +  v - the vector
1326 -  shift - the shift
1327 
1328    Output Parameter:
1329 .  v - the shifted vector
1330 
1331    Level: intermediate
1332 
1333    Concepts: vector^adding constant
1334 
1335 @*/
1336 PetscErrorCode  VecShift(Vec v,PetscScalar shift)
1337 {
1338   PetscErrorCode ierr;
1339   PetscInt       i,n;
1340   PetscScalar    *x;
1341 
1342   PetscFunctionBegin;
1343   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
1344   PetscValidLogicalCollectiveScalar(v,shift,2);
1345   if (v->ops->shift) {
1346     ierr = (*v->ops->shift)(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] += shift;
1352     }
1353     ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1354   }
1355   PetscFunctionReturn(0);
1356 }
1357 
1358 #undef __FUNCT__
1359 #define __FUNCT__ "VecAbs"
1360 /*@
1361    VecAbs - Replaces every element in a vector with its absolute value.
1362 
1363    Logically Collective on Vec
1364 
1365    Input Parameters:
1366 .  v - the vector
1367 
1368    Level: intermediate
1369 
1370    Concepts: vector^absolute value
1371 
1372 @*/
1373 PetscErrorCode  VecAbs(Vec v)
1374 {
1375   PetscErrorCode ierr;
1376   PetscInt       i,n;
1377   PetscScalar    *x;
1378 
1379   PetscFunctionBegin;
1380   PetscValidHeaderSpecific(v,VEC_CLASSID,1);
1381   if (v->ops->abs) {
1382     ierr = (*v->ops->abs)(v);CHKERRQ(ierr);
1383   } else {
1384     ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1385     ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1386     for (i=0; i<n; i++) {
1387       x[i] = PetscAbsScalar(x[i]);
1388     }
1389     ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1390   }
1391   PetscFunctionReturn(0);
1392 }
1393 
1394 #undef __FUNCT__
1395 #define __FUNCT__ "VecPermute"
1396 /*@
1397   VecPermute - Permutes a vector in place using the given ordering.
1398 
1399   Input Parameters:
1400 + vec   - The vector
1401 . order - The ordering
1402 - inv   - The flag for inverting the permutation
1403 
1404   Level: beginner
1405 
1406   Note: This function does not yet support parallel Index Sets
1407 
1408 .seealso: MatPermute()
1409 .keywords: vec, permute
1410 @*/
1411 PetscErrorCode  VecPermute(Vec x, IS row, PetscBool  inv)
1412 {
1413   PetscScalar    *array, *newArray;
1414   const PetscInt *idx;
1415   PetscInt       i;
1416   PetscErrorCode ierr;
1417 
1418   PetscFunctionBegin;
1419   ierr = ISGetIndices(row, &idx);CHKERRQ(ierr);
1420   ierr = VecGetArray(x, &array);CHKERRQ(ierr);
1421   ierr = PetscMalloc(x->map->n*sizeof(PetscScalar), &newArray);CHKERRQ(ierr);
1422 #ifdef PETSC_USE_DEBUG
1423   for (i = 0; i < x->map->n; i++) {
1424     if ((idx[i] < 0) || (idx[i] >= x->map->n)) {
1425       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT, "Permutation index %D is out of bounds: %D", i, idx[i]);
1426     }
1427   }
1428 #endif
1429   if (!inv) {
1430     for (i = 0; i < x->map->n; i++) newArray[i]      = array[idx[i]];
1431   } else {
1432     for (i = 0; i < x->map->n; i++) newArray[idx[i]] = array[i];
1433   }
1434   ierr = VecRestoreArray(x, &array);CHKERRQ(ierr);
1435   ierr = ISRestoreIndices(row, &idx);CHKERRQ(ierr);
1436   ierr = VecReplaceArray(x, newArray);CHKERRQ(ierr);
1437   PetscFunctionReturn(0);
1438 }
1439 
1440 #undef __FUNCT__
1441 #define __FUNCT__ "VecEqual"
1442 /*@
1443    VecEqual - Compares two vectors.
1444 
1445    Collective on Vec
1446 
1447    Input Parameters:
1448 +  vec1 - the first vector
1449 -  vec2 - the second vector
1450 
1451    Output Parameter:
1452 .  flg - PETSC_TRUE if the vectors are equal; PETSC_FALSE otherwise.
1453 
1454    Level: intermediate
1455 
1456    Concepts: equal^two vectors
1457    Concepts: vector^equality
1458 
1459 @*/
1460 PetscErrorCode  VecEqual(Vec vec1,Vec vec2,PetscBool  *flg)
1461 {
1462   PetscScalar    *v1,*v2;
1463   PetscErrorCode ierr;
1464   PetscInt       n1,n2,N1,N2;
1465   PetscInt       state1,state2;
1466   PetscBool      flg1;
1467 
1468   PetscFunctionBegin;
1469   PetscValidHeaderSpecific(vec1,VEC_CLASSID,1);
1470   PetscValidHeaderSpecific(vec2,VEC_CLASSID,2);
1471   PetscValidPointer(flg,3);
1472   if (vec1 == vec2) {
1473     *flg = PETSC_TRUE;
1474   } else {
1475     ierr = VecGetSize(vec1,&N1);CHKERRQ(ierr);
1476     ierr = VecGetSize(vec2,&N2);CHKERRQ(ierr);
1477     if (N1 != N2) {
1478       flg1 = PETSC_FALSE;
1479     } else {
1480       ierr = VecGetLocalSize(vec1,&n1);CHKERRQ(ierr);
1481       ierr = VecGetLocalSize(vec2,&n2);CHKERRQ(ierr);
1482       if (n1 != n2) {
1483         flg1 = PETSC_FALSE;
1484       } else {
1485         ierr = PetscObjectStateQuery((PetscObject) vec1,&state1);CHKERRQ(ierr);
1486         ierr = PetscObjectStateQuery((PetscObject) vec2,&state2);CHKERRQ(ierr);
1487         ierr = VecGetArray(vec1,&v1);CHKERRQ(ierr);
1488         ierr = VecGetArray(vec2,&v2);CHKERRQ(ierr);
1489 #if defined(PETSC_USE_COMPLEX)
1490         {
1491           PetscInt k;
1492           flg1 = PETSC_TRUE;
1493           for (k=0; k<n1; k++){
1494             if (PetscRealPart(v1[k]) != PetscRealPart(v2[k]) || PetscImaginaryPart(v1[k]) != PetscImaginaryPart(v2[k])){
1495               flg1 = PETSC_FALSE;
1496               break;
1497             }
1498           }
1499         }
1500 #else
1501         ierr = PetscMemcmp(v1,v2,n1*sizeof(PetscScalar),&flg1);CHKERRQ(ierr);
1502 #endif
1503         ierr = VecRestoreArray(vec1,&v1);CHKERRQ(ierr);
1504         ierr = VecRestoreArray(vec2,&v2);CHKERRQ(ierr);
1505         ierr = PetscObjectSetState((PetscObject) vec1,state1);CHKERRQ(ierr);
1506         ierr = PetscObjectSetState((PetscObject) vec2,state2);CHKERRQ(ierr);
1507       }
1508     }
1509     /* combine results from all processors */
1510     ierr = MPI_Allreduce(&flg1,flg,1,MPI_INT,MPI_MIN,((PetscObject)vec1)->comm);CHKERRQ(ierr);
1511   }
1512   PetscFunctionReturn(0);
1513 }
1514 
1515 
1516 
1517