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