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