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