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_CLASSID,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_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start); 51 } else if (start >= bs) { 52 SETERRQ2(PETSC_COMM_SELF,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_CLASSID,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_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start); 119 } else if (start >= bs) { 120 SETERRQ2(PETSC_COMM_SELF,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_COMM_SELF,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_CLASSID,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_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start); 210 } else if (start >= bs) { 211 SETERRQ2(PETSC_COMM_SELF,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_CLASSID,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_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start); 306 } else if (start >= bs) { 307 SETERRQ2(PETSC_COMM_SELF,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_CLASSID,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_CLASSID,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_COMM_SELF,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_COMM_SELF,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_CLASSID,1); 537 PetscValidDoublePointer(nrm,3); 538 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"); 539 ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 540 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 541 ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr); 542 543 bs = v->map->bs; 544 if (bs > 128) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently supports only blocksize up to 128"); 545 546 if (!n) { 547 for (j=0; j<bs; j++) { 548 max[j] = PETSC_MIN; 549 } 550 } else { 551 for (j=0; j<bs; j++) { 552 #if defined(PETSC_USE_COMPLEX) 553 max[j] = PetscRealPart(x[j]); 554 #else 555 max[j] = x[j]; 556 #endif 557 } 558 for (i=bs; i<n; i+=bs) { 559 for (j=0; j<bs; j++) { 560 #if defined(PETSC_USE_COMPLEX) 561 if ((tmp = PetscRealPart(x[i+j])) > max[j]) { max[j] = tmp;} 562 #else 563 if ((tmp = x[i+j]) > max[j]) { max[j] = tmp; } 564 #endif 565 } 566 } 567 } 568 ierr = MPI_Allreduce(max,nrm,bs,MPIU_REAL,MPI_MAX,comm);CHKERRQ(ierr); 569 570 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 571 PetscFunctionReturn(0); 572 } 573 574 #undef __FUNCT__ 575 #define __FUNCT__ "VecStrideMinAll" 576 /*@ 577 VecStrideMinAll - Computes the minimum of subvector of a vector defined 578 by a starting point and a stride and optionally its location. 579 580 Collective on Vec 581 582 Input Parameter: 583 . v - the vector 584 585 Output Parameter: 586 + idex - the location where the minimum occurred (not supported, pass PETSC_NULL, 587 if you need this, send mail to petsc-maint@mcs.anl.gov to request it) 588 - nrm - the minimums 589 590 Level: advanced 591 592 Notes: 593 One must call VecSetBlockSize() before this routine to set the stride 594 information, or use a vector created from a multicomponent DA. 595 596 This is useful for computing, say the minimum of the pressure variable when 597 the pressure is stored (interlaced) with other variables, e.g., density, etc. 598 This will only work if the desire subvector is a stride subvector. 599 600 Concepts: minimum^on stride of vector 601 Concepts: stride^minimum 602 603 .seealso: VecMin(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMax() 604 @*/ 605 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideMinAll(Vec v,PetscInt idex[],PetscReal nrm[]) 606 { 607 PetscErrorCode ierr; 608 PetscInt i,n,bs,j; 609 PetscScalar *x; 610 PetscReal min[128],tmp; 611 MPI_Comm comm; 612 613 PetscFunctionBegin; 614 PetscValidHeaderSpecific(v,VEC_CLASSID,1); 615 PetscValidDoublePointer(nrm,3); 616 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"); 617 ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 618 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 619 ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr); 620 621 bs = v->map->bs; 622 if (bs > 128) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently supports only blocksize up to 128"); 623 624 if (!n) { 625 for (j=0; j<bs; j++) { 626 min[j] = PETSC_MAX; 627 } 628 } else { 629 for (j=0; j<bs; j++) { 630 #if defined(PETSC_USE_COMPLEX) 631 min[j] = PetscRealPart(x[j]); 632 #else 633 min[j] = x[j]; 634 #endif 635 } 636 for (i=bs; i<n; i+=bs) { 637 for (j=0; j<bs; j++) { 638 #if defined(PETSC_USE_COMPLEX) 639 if ((tmp = PetscRealPart(x[i+j])) < min[j]) { min[j] = tmp;} 640 #else 641 if ((tmp = x[i+j]) < min[j]) { min[j] = tmp; } 642 #endif 643 } 644 } 645 } 646 ierr = MPI_Allreduce(min,nrm,bs,MPIU_REAL,MPI_MIN,comm);CHKERRQ(ierr); 647 648 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 649 PetscFunctionReturn(0); 650 } 651 652 /*----------------------------------------------------------------------------------------------*/ 653 #undef __FUNCT__ 654 #define __FUNCT__ "VecStrideGatherAll" 655 /*@ 656 VecStrideGatherAll - Gathers all the single components from a multi-component vector into 657 separate vectors. 658 659 Collective on Vec 660 661 Input Parameter: 662 + v - the vector 663 - addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES 664 665 Output Parameter: 666 . s - the location where the subvectors are stored 667 668 Notes: 669 One must call VecSetBlockSize() before this routine to set the stride 670 information, or use a vector created from a multicomponent DA. 671 672 If x is the array representing the vector x then this gathers 673 the arrays (x[start],x[start+stride],x[start+2*stride], ....) 674 for start=0,1,2,...bs-1 675 676 The parallel layout of the vector and the subvector must be the same; 677 i.e., nlocal of v = stride*(nlocal of s) 678 679 Not optimized; could be easily 680 681 Level: advanced 682 683 Concepts: gather^into strided vector 684 685 .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGather(), 686 VecStrideScatterAll() 687 @*/ 688 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideGatherAll(Vec v,Vec s[],InsertMode addv) 689 { 690 PetscErrorCode ierr; 691 PetscInt i,n,n2,bs,j,k,*bss = PETSC_NULL,nv,jj,nvc; 692 PetscScalar *x,**y; 693 694 PetscFunctionBegin; 695 PetscValidHeaderSpecific(v,VEC_CLASSID,1); 696 PetscValidPointer(s,2); 697 PetscValidHeaderSpecific(*s,VEC_CLASSID,2); 698 ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 699 ierr = VecGetLocalSize(s[0],&n2);CHKERRQ(ierr); 700 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 701 bs = v->map->bs; 702 if (bs < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Input vector does not have a valid blocksize set"); 703 704 ierr = PetscMalloc2(bs,PetscReal*,&y,bs,PetscInt,&bss);CHKERRQ(ierr); 705 nv = 0; 706 nvc = 0; 707 for (i=0; i<bs; i++) { 708 ierr = VecGetBlockSize(s[i],&bss[i]);CHKERRQ(ierr); 709 if (bss[i] < 1) bss[i] = 1; /* if user never set it then assume 1 Re: [PETSC #8241] VecStrideGatherAll */ 710 ierr = VecGetArray(s[i],&y[i]);CHKERRQ(ierr); 711 nvc += bss[i]; 712 nv++; 713 if (nvc > bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subvectors in subvectors > number of vectors in main vector"); 714 if (nvc == bs) break; 715 } 716 717 n = n/bs; 718 719 jj = 0; 720 if (addv == INSERT_VALUES) { 721 for (j=0; j<nv; j++) { 722 for (k=0; k<bss[j]; k++) { 723 for (i=0; i<n; i++) { 724 y[j][i*bss[j] + k] = x[bs*i+jj+k]; 725 } 726 } 727 jj += bss[j]; 728 } 729 } else if (addv == ADD_VALUES) { 730 for (j=0; j<nv; j++) { 731 for (k=0; k<bss[j]; k++) { 732 for (i=0; i<n; i++) { 733 y[j][i*bss[j] + k] += x[bs*i+jj+k]; 734 } 735 } 736 jj += bss[j]; 737 } 738 #if !defined(PETSC_USE_COMPLEX) 739 } else if (addv == MAX_VALUES) { 740 for (j=0; j<nv; j++) { 741 for (k=0; k<bss[j]; k++) { 742 for (i=0; i<n; i++) { 743 y[j][i*bss[j] + k] = PetscMax(y[j][i*bss[j] + k],x[bs*i+jj+k]); 744 } 745 } 746 jj += bss[j]; 747 } 748 #endif 749 } else { 750 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type"); 751 } 752 753 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 754 for (i=0; i<nv; i++) { 755 ierr = VecRestoreArray(s[i],&y[i]);CHKERRQ(ierr); 756 } 757 758 ierr = PetscFree2(y,bss);CHKERRQ(ierr); 759 PetscFunctionReturn(0); 760 } 761 762 #undef __FUNCT__ 763 #define __FUNCT__ "VecStrideScatterAll" 764 /*@ 765 VecStrideScatterAll - Scatters all the single components from separate vectors into 766 a multi-component vector. 767 768 Collective on Vec 769 770 Input Parameter: 771 + s - the location where the subvectors are stored 772 - addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES 773 774 Output Parameter: 775 . v - the multicomponent vector 776 777 Notes: 778 One must call VecSetBlockSize() before this routine to set the stride 779 information, or use a vector created from a multicomponent DA. 780 781 The parallel layout of the vector and the subvector must be the same; 782 i.e., nlocal of v = stride*(nlocal of s) 783 784 Not optimized; could be easily 785 786 Level: advanced 787 788 Concepts: scatter^into strided vector 789 790 .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGather(), 791 VecStrideScatterAll() 792 @*/ 793 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideScatterAll(Vec s[],Vec v,InsertMode addv) 794 { 795 PetscErrorCode ierr; 796 PetscInt i,n,n2,bs,j,jj,k,*bss = PETSC_NULL,nv,nvc; 797 PetscScalar *x,**y; 798 799 PetscFunctionBegin; 800 PetscValidHeaderSpecific(v,VEC_CLASSID,1); 801 PetscValidPointer(s,2); 802 PetscValidHeaderSpecific(*s,VEC_CLASSID,2); 803 ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 804 ierr = VecGetLocalSize(s[0],&n2);CHKERRQ(ierr); 805 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 806 bs = v->map->bs; 807 if (bs < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Input vector does not have a valid blocksize set"); 808 809 ierr = PetscMalloc2(bs,PetscScalar**,&y,bs,PetscInt,&bss);CHKERRQ(ierr); 810 nv = 0; 811 nvc = 0; 812 for (i=0; i<bs; i++) { 813 ierr = VecGetBlockSize(s[i],&bss[i]);CHKERRQ(ierr); 814 if (bss[i] < 1) bss[i] = 1; /* if user never set it then assume 1 Re: [PETSC #8241] VecStrideGatherAll */ 815 ierr = VecGetArray(s[i],&y[i]);CHKERRQ(ierr); 816 nvc += bss[i]; 817 nv++; 818 if (nvc > bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subvectors in subvectors > number of vectors in main vector"); 819 if (nvc == bs) break; 820 } 821 822 n = n/bs; 823 824 jj = 0; 825 if (addv == INSERT_VALUES) { 826 for (j=0; j<nv; j++) { 827 for (k=0; k<bss[j]; k++) { 828 for (i=0; i<n; i++) { 829 x[bs*i+jj+k] = y[j][i*bss[j] + k]; 830 } 831 } 832 jj += bss[j]; 833 } 834 } else if (addv == ADD_VALUES) { 835 for (j=0; j<nv; j++) { 836 for (k=0; k<bss[j]; k++) { 837 for (i=0; i<n; i++) { 838 x[bs*i+jj+k] += y[j][i*bss[j] + k]; 839 } 840 } 841 jj += bss[j]; 842 } 843 #if !defined(PETSC_USE_COMPLEX) 844 } else if (addv == MAX_VALUES) { 845 for (j=0; j<nv; j++) { 846 for (k=0; k<bss[j]; k++) { 847 for (i=0; i<n; i++) { 848 x[bs*i+jj+k] = PetscMax(x[bs*i+jj+k],y[j][i*bss[j] + k]); 849 } 850 } 851 jj += bss[j]; 852 } 853 #endif 854 } else { 855 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type"); 856 } 857 858 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 859 for (i=0; i<nv; i++) { 860 ierr = VecRestoreArray(s[i],&y[i]);CHKERRQ(ierr); 861 } 862 ierr = PetscFree2(y,bss);CHKERRQ(ierr); 863 PetscFunctionReturn(0); 864 } 865 866 #undef __FUNCT__ 867 #define __FUNCT__ "VecStrideGather" 868 /*@ 869 VecStrideGather - Gathers a single component from a multi-component vector into 870 another vector. 871 872 Collective on Vec 873 874 Input Parameter: 875 + v - the vector 876 . start - starting point of the subvector (defined by a stride) 877 - addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES 878 879 Output Parameter: 880 . s - the location where the subvector is stored 881 882 Notes: 883 One must call VecSetBlockSize() before this routine to set the stride 884 information, or use a vector created from a multicomponent DA. 885 886 If x is the array representing the vector x then this gathers 887 the array (x[start],x[start+stride],x[start+2*stride], ....) 888 889 The parallel layout of the vector and the subvector must be the same; 890 i.e., nlocal of v = stride*(nlocal of s) 891 892 Not optimized; could be easily 893 894 Level: advanced 895 896 Concepts: gather^into strided vector 897 898 .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGatherAll(), 899 VecStrideScatterAll() 900 @*/ 901 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideGather(Vec v,PetscInt start,Vec s,InsertMode addv) 902 { 903 PetscErrorCode ierr; 904 905 PetscFunctionBegin; 906 PetscValidHeaderSpecific(v,VEC_CLASSID,1); 907 PetscValidHeaderSpecific(s,VEC_CLASSID,3); 908 if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start); 909 if (start >= v->map->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Start of stride subvector (%D) is too large for stride\n\ 910 Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,v->map->bs); 911 if (!v->ops->stridegather) SETERRQ(((PetscObject)s)->comm,PETSC_ERR_SUP,"Not implemented for this Vec class"); 912 ierr = (*v->ops->stridegather)(v,start,s,addv);CHKERRQ(ierr); 913 PetscFunctionReturn(0); 914 } 915 916 #undef __FUNCT__ 917 #define __FUNCT__ "VecStrideScatter" 918 /*@ 919 VecStrideScatter - Scatters a single component from a vector into a multi-component vector. 920 921 Collective on Vec 922 923 Input Parameter: 924 + s - the single-component vector 925 . start - starting point of the subvector (defined by a stride) 926 - addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES 927 928 Output Parameter: 929 . v - the location where the subvector is scattered (the multi-component vector) 930 931 Notes: 932 One must call VecSetBlockSize() on the multi-component vector before this 933 routine to set the stride information, or use a vector created from a multicomponent DA. 934 935 The parallel layout of the vector and the subvector must be the same; 936 i.e., nlocal of v = stride*(nlocal of s) 937 938 Not optimized; could be easily 939 940 Level: advanced 941 942 Concepts: scatter^into strided vector 943 944 .seealso: VecStrideNorm(), VecStrideGather(), VecStrideMin(), VecStrideMax(), VecStrideGatherAll(), 945 VecStrideScatterAll() 946 @*/ 947 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideScatter(Vec s,PetscInt start,Vec v,InsertMode addv) 948 { 949 PetscErrorCode ierr; 950 951 PetscFunctionBegin; 952 PetscValidHeaderSpecific(v,VEC_CLASSID,1); 953 PetscValidHeaderSpecific(s,VEC_CLASSID,3); 954 if (start < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start); 955 if (start >= v->map->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Start of stride subvector (%D) is too large for stride\n\ 956 Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,v->map->bs); 957 if (!v->ops->stridescatter) SETERRQ(((PetscObject)s)->comm,PETSC_ERR_SUP,"Not implemented for this Vec class"); 958 ierr = (*v->ops->stridescatter)(s,start,v,addv);CHKERRQ(ierr); 959 PetscFunctionReturn(0); 960 } 961 962 #undef __FUNCT__ 963 #define __FUNCT__ "VecStrideGather_Default" 964 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideGather_Default(Vec v,PetscInt start,Vec s,InsertMode addv) 965 { 966 PetscErrorCode ierr; 967 PetscInt i,n,bs,ns; 968 PetscScalar *x,*y; 969 970 PetscFunctionBegin; 971 ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 972 ierr = VecGetLocalSize(s,&ns);CHKERRQ(ierr); 973 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 974 ierr = VecGetArray(s,&y);CHKERRQ(ierr); 975 976 bs = v->map->bs; 977 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); 978 x += start; 979 n = n/bs; 980 981 if (addv == INSERT_VALUES) { 982 for (i=0; i<n; i++) { 983 y[i] = x[bs*i]; 984 } 985 } else if (addv == ADD_VALUES) { 986 for (i=0; i<n; i++) { 987 y[i] += x[bs*i]; 988 } 989 #if !defined(PETSC_USE_COMPLEX) 990 } else if (addv == MAX_VALUES) { 991 for (i=0; i<n; i++) { 992 y[i] = PetscMax(y[i],x[bs*i]); 993 } 994 #endif 995 } else { 996 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type"); 997 } 998 999 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1000 ierr = VecRestoreArray(s,&y);CHKERRQ(ierr); 1001 PetscFunctionReturn(0); 1002 } 1003 1004 #undef __FUNCT__ 1005 #define __FUNCT__ "VecStrideScatter_Default" 1006 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideScatter_Default(Vec s,PetscInt start,Vec v,InsertMode addv) 1007 { 1008 PetscErrorCode ierr; 1009 PetscInt i,n,bs,ns; 1010 PetscScalar *x,*y; 1011 1012 PetscFunctionBegin; 1013 ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 1014 ierr = VecGetLocalSize(s,&ns);CHKERRQ(ierr); 1015 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1016 ierr = VecGetArray(s,&y);CHKERRQ(ierr); 1017 1018 bs = v->map->bs; 1019 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); 1020 x += start; 1021 n = n/bs; 1022 1023 if (addv == INSERT_VALUES) { 1024 for (i=0; i<n; i++) { 1025 x[bs*i] = y[i]; 1026 } 1027 } else if (addv == ADD_VALUES) { 1028 for (i=0; i<n; i++) { 1029 x[bs*i] += y[i]; 1030 } 1031 #if !defined(PETSC_USE_COMPLEX) 1032 } else if (addv == MAX_VALUES) { 1033 for (i=0; i<n; i++) { 1034 x[bs*i] = PetscMax(y[i],x[bs*i]); 1035 } 1036 #endif 1037 } else { 1038 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type"); 1039 } 1040 1041 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1042 ierr = VecRestoreArray(s,&y);CHKERRQ(ierr); 1043 PetscFunctionReturn(0); 1044 } 1045 1046 #undef __FUNCT__ 1047 #define __FUNCT__ "VecReciprocal_Default" 1048 PetscErrorCode VecReciprocal_Default(Vec v) 1049 { 1050 PetscErrorCode ierr; 1051 PetscInt i,n; 1052 PetscScalar *x; 1053 1054 PetscFunctionBegin; 1055 ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 1056 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1057 for (i=0; i<n; i++) { 1058 if (x[i] != 0.0) x[i] = 1.0/x[i]; 1059 } 1060 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1061 PetscFunctionReturn(0); 1062 } 1063 1064 #undef __FUNCT__ 1065 #define __FUNCT__ "VecExp" 1066 /*@ 1067 VecExp - Replaces each component of a vector by e^x_i 1068 1069 Not collective 1070 1071 Input Parameter: 1072 . v - The vector 1073 1074 Output Parameter: 1075 . v - The vector of exponents 1076 1077 Level: beginner 1078 1079 .seealso: VecLog(), VecAbs(), VecSqrtAbs(), VecReciprocal() 1080 1081 .keywords: vector, sqrt, square root 1082 @*/ 1083 PetscErrorCode PETSCVEC_DLLEXPORT VecExp(Vec v) 1084 { 1085 PetscScalar *x; 1086 PetscInt i, n; 1087 PetscErrorCode ierr; 1088 1089 PetscFunctionBegin; 1090 PetscValidHeaderSpecific(v, VEC_CLASSID,1); 1091 if (v->ops->exp) { 1092 ierr = (*v->ops->exp)(v);CHKERRQ(ierr); 1093 } else { 1094 ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr); 1095 ierr = VecGetArray(v, &x);CHKERRQ(ierr); 1096 for(i = 0; i < n; i++) { 1097 x[i] = PetscExpScalar(x[i]); 1098 } 1099 ierr = VecRestoreArray(v, &x);CHKERRQ(ierr); 1100 } 1101 PetscFunctionReturn(0); 1102 } 1103 1104 #undef __FUNCT__ 1105 #define __FUNCT__ "VecLog" 1106 /*@ 1107 VecLog - Replaces each component of a vector by log(x_i), the natural logarithm 1108 1109 Not collective 1110 1111 Input Parameter: 1112 . v - The vector 1113 1114 Output Parameter: 1115 . v - The vector of logs 1116 1117 Level: beginner 1118 1119 .seealso: VecExp(), VecAbs(), VecSqrtAbs(), VecReciprocal() 1120 1121 .keywords: vector, sqrt, square root 1122 @*/ 1123 PetscErrorCode PETSCVEC_DLLEXPORT VecLog(Vec v) 1124 { 1125 PetscScalar *x; 1126 PetscInt i, n; 1127 PetscErrorCode ierr; 1128 1129 PetscFunctionBegin; 1130 PetscValidHeaderSpecific(v, VEC_CLASSID,1); 1131 if (v->ops->log) { 1132 ierr = (*v->ops->log)(v);CHKERRQ(ierr); 1133 } else { 1134 ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr); 1135 ierr = VecGetArray(v, &x);CHKERRQ(ierr); 1136 for(i = 0; i < n; i++) { 1137 x[i] = PetscLogScalar(x[i]); 1138 } 1139 ierr = VecRestoreArray(v, &x);CHKERRQ(ierr); 1140 } 1141 PetscFunctionReturn(0); 1142 } 1143 1144 #undef __FUNCT__ 1145 #define __FUNCT__ "VecSqrtAbs" 1146 /*@ 1147 VecSqrtAbs - Replaces each component of a vector by the square root of its magnitude. 1148 1149 Not collective 1150 1151 Input Parameter: 1152 . v - The vector 1153 1154 Output Parameter: 1155 . v - The vector square root 1156 1157 Level: beginner 1158 1159 Note: The actual function is sqrt(|x_i|) 1160 1161 .seealso: VecLog(), VecExp(), VecReciprocal(), VecAbs() 1162 1163 .keywords: vector, sqrt, square root 1164 @*/ 1165 PetscErrorCode PETSCVEC_DLLEXPORT VecSqrtAbs(Vec v) 1166 { 1167 PetscScalar *x; 1168 PetscInt i, n; 1169 PetscErrorCode ierr; 1170 1171 PetscFunctionBegin; 1172 PetscValidHeaderSpecific(v, VEC_CLASSID,1); 1173 if (v->ops->sqrt) { 1174 ierr = (*v->ops->sqrt)(v);CHKERRQ(ierr); 1175 } else { 1176 ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr); 1177 ierr = VecGetArray(v, &x);CHKERRQ(ierr); 1178 for(i = 0; i < n; i++) { 1179 x[i] = sqrt(PetscAbsScalar(x[i])); 1180 } 1181 ierr = VecRestoreArray(v, &x);CHKERRQ(ierr); 1182 } 1183 PetscFunctionReturn(0); 1184 } 1185 1186 #undef __FUNCT__ 1187 #define __FUNCT__ "VecDotNorm2" 1188 /*@ 1189 VecDotNorm2 - computes the inner product of two vectors and the 2-norm squared of the second vector 1190 1191 Collective on Vec 1192 1193 Input Parameter: 1194 + s - first vector 1195 - t - second vector 1196 1197 Output Parameter: 1198 + dp - s't 1199 - nm - t't 1200 1201 Level: advanced 1202 1203 .seealso: VecDot(), VecNorm(), VecDotBegin(), VecNormBegin(), VecDotEnd(), VecNormEnd() 1204 1205 .keywords: vector, sqrt, square root 1206 @*/ 1207 PetscErrorCode PETSCVEC_DLLEXPORT VecDotNorm2(Vec s,Vec t,PetscScalar *dp, PetscScalar *nm) 1208 { 1209 PetscScalar *sx, *tx, dpx = 0.0, nmx = 0.0,work[2],sum[2]; 1210 PetscInt i, n; 1211 PetscErrorCode ierr; 1212 1213 PetscFunctionBegin; 1214 PetscValidHeaderSpecific(s, VEC_CLASSID,1); 1215 PetscValidHeaderSpecific(t, VEC_CLASSID,2); 1216 1217 ierr = PetscLogEventBarrierBegin(VEC_DotNormBarrier,s,t,0,0,((PetscObject)s)->comm);CHKERRQ(ierr); 1218 ierr = VecGetLocalSize(s, &n);CHKERRQ(ierr); 1219 ierr = VecGetArray(s, &sx);CHKERRQ(ierr); 1220 ierr = VecGetArray(t, &tx);CHKERRQ(ierr); 1221 1222 for (i = 0; i<n; i++) { 1223 dpx += sx[i]*PetscConj(tx[i]); 1224 nmx += tx[i]*PetscConj(tx[i]); 1225 } 1226 work[0] = dpx; 1227 work[1] = nmx; 1228 ierr = MPI_Allreduce(&work,&sum,2,MPIU_SCALAR,MPIU_SUM,((PetscObject)s)->comm);CHKERRQ(ierr); 1229 *dp = sum[0]; 1230 *nm = sum[1]; 1231 1232 ierr = VecRestoreArray(t, &tx);CHKERRQ(ierr); 1233 ierr = VecRestoreArray(s, &sx);CHKERRQ(ierr); 1234 ierr = PetscLogFlops(4.0*n);CHKERRQ(ierr); 1235 ierr = PetscLogEventBarrierEnd(VEC_DotNormBarrier,s,t,0,0,((PetscObject)s)->comm);CHKERRQ(ierr); 1236 PetscFunctionReturn(0); 1237 } 1238 1239 #undef __FUNCT__ 1240 #define __FUNCT__ "VecSum" 1241 /*@ 1242 VecSum - Computes the sum of all the components of a vector. 1243 1244 Collective on Vec 1245 1246 Input Parameter: 1247 . v - the vector 1248 1249 Output Parameter: 1250 . sum - the result 1251 1252 Level: beginner 1253 1254 Concepts: sum^of vector entries 1255 1256 .seealso: VecNorm() 1257 @*/ 1258 PetscErrorCode PETSCVEC_DLLEXPORT VecSum(Vec v,PetscScalar *sum) 1259 { 1260 PetscErrorCode ierr; 1261 PetscInt i,n; 1262 PetscScalar *x,lsum = 0.0; 1263 1264 PetscFunctionBegin; 1265 PetscValidHeaderSpecific(v,VEC_CLASSID,1); 1266 PetscValidScalarPointer(sum,2); 1267 ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 1268 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1269 for (i=0; i<n; i++) { 1270 lsum += x[i]; 1271 } 1272 ierr = MPI_Allreduce(&lsum,sum,1,MPIU_SCALAR,MPIU_SUM,((PetscObject)v)->comm);CHKERRQ(ierr); 1273 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1274 PetscFunctionReturn(0); 1275 } 1276 1277 #undef __FUNCT__ 1278 #define __FUNCT__ "VecShift" 1279 /*@ 1280 VecShift - Shifts all of the components of a vector by computing 1281 x[i] = x[i] + shift. 1282 1283 Collective on Vec 1284 1285 Input Parameters: 1286 + v - the vector 1287 - shift - the shift 1288 1289 Output Parameter: 1290 . v - the shifted vector 1291 1292 Level: intermediate 1293 1294 Concepts: vector^adding constant 1295 1296 @*/ 1297 PetscErrorCode PETSCVEC_DLLEXPORT VecShift(Vec v,PetscScalar shift) 1298 { 1299 PetscErrorCode ierr; 1300 PetscInt i,n; 1301 PetscScalar *x; 1302 1303 PetscFunctionBegin; 1304 PetscValidHeaderSpecific(v,VEC_CLASSID,1); 1305 if (v->ops->shift) { 1306 ierr = (*v->ops->shift)(v);CHKERRQ(ierr); 1307 } else { 1308 ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 1309 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1310 for (i=0; i<n; i++) { 1311 x[i] += shift; 1312 } 1313 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1314 } 1315 PetscFunctionReturn(0); 1316 } 1317 1318 #undef __FUNCT__ 1319 #define __FUNCT__ "VecAbs" 1320 /*@ 1321 VecAbs - Replaces every element in a vector with its absolute value. 1322 1323 Collective on Vec 1324 1325 Input Parameters: 1326 . v - the vector 1327 1328 Level: intermediate 1329 1330 Concepts: vector^absolute value 1331 1332 @*/ 1333 PetscErrorCode PETSCVEC_DLLEXPORT VecAbs(Vec v) 1334 { 1335 PetscErrorCode ierr; 1336 PetscInt i,n; 1337 PetscScalar *x; 1338 1339 PetscFunctionBegin; 1340 PetscValidHeaderSpecific(v,VEC_CLASSID,1); 1341 if (v->ops->abs) { 1342 ierr = (*v->ops->abs)(v);CHKERRQ(ierr); 1343 } else { 1344 ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 1345 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1346 for (i=0; i<n; i++) { 1347 x[i] = PetscAbsScalar(x[i]); 1348 } 1349 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1350 } 1351 PetscFunctionReturn(0); 1352 } 1353 1354 #undef __FUNCT__ 1355 #define __FUNCT__ "VecPermute" 1356 /*@ 1357 VecPermute - Permutes a vector in place using the given ordering. 1358 1359 Input Parameters: 1360 + vec - The vector 1361 . order - The ordering 1362 - inv - The flag for inverting the permutation 1363 1364 Level: beginner 1365 1366 Note: This function does not yet support parallel Index Sets 1367 1368 .seealso: MatPermute() 1369 .keywords: vec, permute 1370 @*/ 1371 PetscErrorCode PETSCVEC_DLLEXPORT VecPermute(Vec x, IS row, PetscTruth inv) 1372 { 1373 PetscScalar *array, *newArray; 1374 const PetscInt *idx; 1375 PetscInt i; 1376 PetscErrorCode ierr; 1377 1378 PetscFunctionBegin; 1379 ierr = ISGetIndices(row, &idx);CHKERRQ(ierr); 1380 ierr = VecGetArray(x, &array);CHKERRQ(ierr); 1381 ierr = PetscMalloc(x->map->n*sizeof(PetscScalar), &newArray);CHKERRQ(ierr); 1382 #ifdef PETSC_USE_DEBUG 1383 for(i = 0; i < x->map->n; i++) { 1384 if ((idx[i] < 0) || (idx[i] >= x->map->n)) { 1385 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT, "Permutation index %D is out of bounds: %D", i, idx[i]); 1386 } 1387 } 1388 #endif 1389 if (!inv) { 1390 for(i = 0; i < x->map->n; i++) newArray[i] = array[idx[i]]; 1391 } else { 1392 for(i = 0; i < x->map->n; i++) newArray[idx[i]] = array[i]; 1393 } 1394 ierr = VecRestoreArray(x, &array);CHKERRQ(ierr); 1395 ierr = ISRestoreIndices(row, &idx);CHKERRQ(ierr); 1396 ierr = VecReplaceArray(x, newArray);CHKERRQ(ierr); 1397 PetscFunctionReturn(0); 1398 } 1399 1400 #undef __FUNCT__ 1401 #define __FUNCT__ "VecEqual" 1402 /*@ 1403 VecEqual - Compares two vectors. 1404 1405 Collective on Vec 1406 1407 Input Parameters: 1408 + vec1 - the first vector 1409 - vec2 - the second vector 1410 1411 Output Parameter: 1412 . flg - PETSC_TRUE if the vectors are equal; PETSC_FALSE otherwise. 1413 1414 Level: intermediate 1415 1416 Concepts: equal^two vectors 1417 Concepts: vector^equality 1418 1419 @*/ 1420 PetscErrorCode PETSCVEC_DLLEXPORT VecEqual(Vec vec1,Vec vec2,PetscTruth *flg) 1421 { 1422 PetscScalar *v1,*v2; 1423 PetscErrorCode ierr; 1424 PetscInt n1,n2,N1,N2; 1425 PetscInt state1,state2; 1426 PetscTruth flg1; 1427 1428 PetscFunctionBegin; 1429 PetscValidHeaderSpecific(vec1,VEC_CLASSID,1); 1430 PetscValidHeaderSpecific(vec2,VEC_CLASSID,2); 1431 PetscValidPointer(flg,3); 1432 if (vec1 == vec2) { 1433 *flg = PETSC_TRUE; 1434 } else { 1435 ierr = VecGetSize(vec1,&N1);CHKERRQ(ierr); 1436 ierr = VecGetSize(vec2,&N2);CHKERRQ(ierr); 1437 if (N1 != N2) { 1438 flg1 = PETSC_FALSE; 1439 } else { 1440 ierr = VecGetLocalSize(vec1,&n1);CHKERRQ(ierr); 1441 ierr = VecGetLocalSize(vec2,&n2);CHKERRQ(ierr); 1442 if (n1 != n2) { 1443 flg1 = PETSC_FALSE; 1444 } else { 1445 ierr = PetscObjectStateQuery((PetscObject) vec1,&state1);CHKERRQ(ierr); 1446 ierr = PetscObjectStateQuery((PetscObject) vec2,&state2);CHKERRQ(ierr); 1447 ierr = VecGetArray(vec1,&v1);CHKERRQ(ierr); 1448 ierr = VecGetArray(vec2,&v2);CHKERRQ(ierr); 1449 #if defined(PETSC_USE_COMPLEX) 1450 PetscInt k; 1451 for (k=0; k<n1; k++){ 1452 if (PetscRealPart(v1[k]) != PetscRealPart(v2[k]) || PetscImaginaryPart(v1[k]) != PetscImaginaryPart(v2[k])){ 1453 flg1 = PETSC_FALSE; 1454 break; 1455 } 1456 } 1457 #else 1458 ierr = PetscMemcmp(v1,v2,n1*sizeof(PetscScalar),&flg1);CHKERRQ(ierr); 1459 #endif 1460 ierr = VecRestoreArray(vec1,&v1);CHKERRQ(ierr); 1461 ierr = VecRestoreArray(vec2,&v2);CHKERRQ(ierr); 1462 ierr = PetscObjectSetState((PetscObject) vec1,state1);CHKERRQ(ierr); 1463 ierr = PetscObjectSetState((PetscObject) vec2,state2);CHKERRQ(ierr); 1464 } 1465 } 1466 /* combine results from all processors */ 1467 ierr = MPI_Allreduce(&flg1,flg,1,MPI_INT,MPI_MIN,((PetscObject)vec1)->comm);CHKERRQ(ierr); 1468 } 1469 PetscFunctionReturn(0); 1470 } 1471 1472 1473 1474