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