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