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