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