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