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