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