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