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