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