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