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