1 #define PETSCDM_DLL 2 3 #include "private/daimpl.h" /*I "petscda.h" I*/ 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "DAView_2d" 7 PetscErrorCode DAView_2d(DA da,PetscViewer viewer) 8 { 9 PetscErrorCode ierr; 10 PetscMPIInt rank; 11 PetscBool iascii,isdraw; 12 DM_DA *dd = (DM_DA*)da->data; 13 14 PetscFunctionBegin; 15 ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr); 16 17 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 18 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 19 if (iascii) { 20 PetscViewerFormat format; 21 22 ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 23 if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) { 24 DALocalInfo info; 25 ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr); 26 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D N %D m %D n %D w %D s %D\n",rank,dd->M,dd->N,dd->m,dd->n,dd->w,dd->s);CHKERRQ(ierr); 27 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D\n",info.xs,info.xs+info.xm,info.ys,info.ys+info.ym);CHKERRQ(ierr); 28 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 29 } 30 } else if (isdraw) { 31 PetscDraw draw; 32 double ymin = -1*dd->s-1,ymax = dd->N+dd->s; 33 double xmin = -1*dd->s-1,xmax = dd->M+dd->s; 34 double x,y; 35 PetscInt base,*idx; 36 char node[10]; 37 PetscBool isnull; 38 39 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 40 ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 41 if (!dd->coordinates) { 42 ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr); 43 } 44 ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr); 45 46 /* first processor draw all node lines */ 47 if (!rank) { 48 ymin = 0.0; ymax = dd->N - 1; 49 for (xmin=0; xmin<dd->M; xmin++) { 50 ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr); 51 } 52 xmin = 0.0; xmax = dd->M - 1; 53 for (ymin=0; ymin<dd->N; ymin++) { 54 ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr); 55 } 56 } 57 ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 58 ierr = PetscDrawPause(draw);CHKERRQ(ierr); 59 60 /* draw my box */ 61 ymin = dd->ys; ymax = dd->ye - 1; xmin = dd->xs/dd->w; 62 xmax =(dd->xe-1)/dd->w; 63 ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr); 64 ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 65 ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 66 ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr); 67 68 /* put in numbers */ 69 base = (dd->base)/dd->w; 70 for (y=ymin; y<=ymax; y++) { 71 for (x=xmin; x<=xmax; x++) { 72 sprintf(node,"%d",(int)base++); 73 ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr); 74 } 75 } 76 77 ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 78 ierr = PetscDrawPause(draw);CHKERRQ(ierr); 79 /* overlay ghost numbers, useful for error checking */ 80 /* put in numbers */ 81 82 base = 0; idx = dd->idx; 83 ymin = dd->Ys; ymax = dd->Ye; xmin = dd->Xs; xmax = dd->Xe; 84 for (y=ymin; y<ymax; y++) { 85 for (x=xmin; x<xmax; x++) { 86 if ((base % dd->w) == 0) { 87 sprintf(node,"%d",(int)(idx[base]/dd->w)); 88 ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr); 89 } 90 base++; 91 } 92 } 93 ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 94 ierr = PetscDrawPause(draw);CHKERRQ(ierr); 95 } else { 96 SETERRQ1(((PetscObject)viewer)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DA2d",((PetscObject)viewer)->type_name); 97 } 98 PetscFunctionReturn(0); 99 } 100 101 /* 102 M is number of grid points 103 m is number of processors 104 105 */ 106 #undef __FUNCT__ 107 #define __FUNCT__ "DASplitComm2d" 108 PetscErrorCode PETSCDM_DLLEXPORT DASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm) 109 { 110 PetscErrorCode ierr; 111 PetscInt m,n = 0,x = 0,y = 0; 112 PetscMPIInt size,csize,rank; 113 114 PetscFunctionBegin; 115 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 116 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 117 118 csize = 4*size; 119 do { 120 if (csize % 4) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot split communicator of size %d tried %d %D %D",size,csize,x,y); 121 csize = csize/4; 122 123 m = (PetscInt)(0.5 + sqrt(((double)M)*((double)csize)/((double)N))); 124 if (!m) m = 1; 125 while (m > 0) { 126 n = csize/m; 127 if (m*n == csize) break; 128 m--; 129 } 130 if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 131 132 x = M/m + ((M % m) > ((csize-1) % m)); 133 y = (N + (csize-1)/m)/n; 134 } while ((x < 4 || y < 4) && csize > 1); 135 if (size != csize) { 136 MPI_Group entire_group,sub_group; 137 PetscMPIInt i,*groupies; 138 139 ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr); 140 ierr = PetscMalloc(csize*sizeof(PetscInt),&groupies);CHKERRQ(ierr); 141 for (i=0; i<csize; i++) { 142 groupies[i] = (rank/csize)*csize + i; 143 } 144 ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr); 145 ierr = PetscFree(groupies);CHKERRQ(ierr); 146 ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr); 147 ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr); 148 ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr); 149 ierr = PetscInfo1(0,"DASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr); 150 } else { 151 *outcomm = comm; 152 } 153 PetscFunctionReturn(0); 154 } 155 156 #undef __FUNCT__ 157 #define __FUNCT__ "DMGetElements_DA_2d_P1" 158 PetscErrorCode DMGetElements_DA_2d_P1(DA da,PetscInt *n,const PetscInt *e[]) 159 { 160 PetscErrorCode ierr; 161 DM_DA *dd = (DM_DA*)da->data; 162 PetscInt i,j,cnt,xs,xe = dd->xe,ys,ye = dd->ye,Xs = dd->Xs, Xe = dd->Xe, Ys = dd->Ys; 163 164 PetscFunctionBegin; 165 if (!dd->e) { 166 if (dd->xs == Xs) xs = dd->xs; else xs = dd->xs - 1; 167 if (dd->ys == Ys) ys = dd->ys; else ys = dd->ys - 1; 168 dd->ne = 2*(xe - xs - 1)*(ye - ys - 1); 169 ierr = PetscMalloc((1 + 3*dd->ne)*sizeof(PetscInt),&dd->e);CHKERRQ(ierr); 170 cnt = 0; 171 for (j=ys; j<ye-1; j++) { 172 for (i=xs; i<xe-1; i++) { 173 dd->e[cnt] = i - Xs + (j - Ys)*(Xe - Xs); 174 dd->e[cnt+1] = i - Xs + 1 + (j - Ys)*(Xe - Xs); 175 dd->e[cnt+2] = i - Xs + (j - Ys + 1)*(Xe - Xs); 176 177 dd->e[cnt+3] = i - Xs + 1 + (j - Ys + 1)*(Xe - Xs); 178 dd->e[cnt+4] = i - Xs + (j - Ys + 1)*(Xe - Xs); 179 dd->e[cnt+5] = i - Xs + 1 + (j - Ys)*(Xe - Xs); 180 cnt += 6; 181 } 182 } 183 } 184 *n = dd->ne; 185 *e = dd->e; 186 PetscFunctionReturn(0); 187 } 188 189 #undef __FUNCT__ 190 #define __FUNCT__ "DASetLocalFunction" 191 /*@C 192 DASetLocalFunction - Caches in a DA a local function. 193 194 Logically Collective on DA 195 196 Input Parameter: 197 + da - initial distributed array 198 - lf - the local function 199 200 Level: intermediate 201 202 Notes: The routine SNESDAFormFunction() uses this the cached function to evaluate the user provided function. 203 204 .keywords: distributed array, refine 205 206 .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunctioni() 207 @*/ 208 PetscErrorCode PETSCDM_DLLEXPORT DASetLocalFunction(DA da,DALocalFunction1 lf) 209 { 210 DM_DA *dd = (DM_DA*)da->data; 211 PetscFunctionBegin; 212 PetscValidHeaderSpecific(da,DM_CLASSID,1); 213 dd->lf = lf; 214 PetscFunctionReturn(0); 215 } 216 217 #undef __FUNCT__ 218 #define __FUNCT__ "DASetLocalFunctioni" 219 /*@C 220 DASetLocalFunctioni - Caches in a DA a local function that evaluates a single component 221 222 Logically Collective on DA 223 224 Input Parameter: 225 + da - initial distributed array 226 - lfi - the local function 227 228 Level: intermediate 229 230 .keywords: distributed array, refine 231 232 .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction() 233 @*/ 234 PetscErrorCode PETSCDM_DLLEXPORT DASetLocalFunctioni(DA da,PetscErrorCode (*lfi)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*)) 235 { 236 DM_DA *dd = (DM_DA*)da->data; 237 PetscFunctionBegin; 238 PetscValidHeaderSpecific(da,DM_CLASSID,1); 239 dd->lfi = lfi; 240 PetscFunctionReturn(0); 241 } 242 243 #undef __FUNCT__ 244 #define __FUNCT__ "DASetLocalFunctionib" 245 /*@C 246 DASetLocalFunctionib - Caches in a DA a block local function that evaluates a single component 247 248 Logically Collective on DA 249 250 Input Parameter: 251 + da - initial distributed array 252 - lfi - the local function 253 254 Level: intermediate 255 256 .keywords: distributed array, refine 257 258 .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction() 259 @*/ 260 PetscErrorCode PETSCDM_DLLEXPORT DASetLocalFunctionib(DA da,PetscErrorCode (*lfi)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*)) 261 { 262 DM_DA *dd = (DM_DA*)da->data; 263 PetscFunctionBegin; 264 PetscValidHeaderSpecific(da,DM_CLASSID,1); 265 dd->lfib = lfi; 266 PetscFunctionReturn(0); 267 } 268 269 #undef __FUNCT__ 270 #define __FUNCT__ "DASetLocalAdicFunction_Private" 271 PetscErrorCode DASetLocalAdicFunction_Private(DA da,DALocalFunction1 ad_lf) 272 { 273 DM_DA *dd = (DM_DA*)da->data; 274 PetscFunctionBegin; 275 PetscValidHeaderSpecific(da,DM_CLASSID,1); 276 dd->adic_lf = ad_lf; 277 PetscFunctionReturn(0); 278 } 279 280 /*MC 281 DASetLocalAdicFunctioni - Caches in a DA a local functioni computed by ADIC/ADIFOR 282 283 Synopsis: 284 PetscErrorCode DASetLocalAdicFunctioni(DA da,PetscInt (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*) 285 286 Logically Collective on DA 287 288 Input Parameter: 289 + da - initial distributed array 290 - ad_lfi - the local function as computed by ADIC/ADIFOR 291 292 Level: intermediate 293 294 .keywords: distributed array, refine 295 296 .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(), 297 DASetLocalJacobian(), DASetLocalFunctioni() 298 M*/ 299 300 #undef __FUNCT__ 301 #define __FUNCT__ "DASetLocalAdicFunctioni_Private" 302 PetscErrorCode DASetLocalAdicFunctioni_Private(DA da,PetscErrorCode (*ad_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*)) 303 { 304 DM_DA *dd = (DM_DA*)da->data; 305 PetscFunctionBegin; 306 PetscValidHeaderSpecific(da,DM_CLASSID,1); 307 dd->adic_lfi = ad_lfi; 308 PetscFunctionReturn(0); 309 } 310 311 /*MC 312 DASetLocalAdicMFFunctioni - Caches in a DA a local functioni computed by ADIC/ADIFOR 313 314 Synopsis: 315 PetscErrorCode DASetLocalAdicFunctioni(DA da,int (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*) 316 317 Logically Collective on DA 318 319 Input Parameter: 320 + da - initial distributed array 321 - admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR 322 323 Level: intermediate 324 325 .keywords: distributed array, refine 326 327 .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(), 328 DASetLocalJacobian(), DASetLocalFunctioni() 329 M*/ 330 331 #undef __FUNCT__ 332 #define __FUNCT__ "DASetLocalAdicMFFunctioni_Private" 333 PetscErrorCode DASetLocalAdicMFFunctioni_Private(DA da,PetscErrorCode (*admf_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*)) 334 { 335 DM_DA *dd = (DM_DA*)da->data; 336 PetscFunctionBegin; 337 PetscValidHeaderSpecific(da,DM_CLASSID,1); 338 dd->adicmf_lfi = admf_lfi; 339 PetscFunctionReturn(0); 340 } 341 342 /*MC 343 DASetLocalAdicFunctionib - Caches in a DA a block local functioni computed by ADIC/ADIFOR 344 345 Synopsis: 346 PetscErrorCode DASetLocalAdicFunctionib(DA da,PetscInt (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*) 347 348 Logically Collective on DA 349 350 Input Parameter: 351 + da - initial distributed array 352 - ad_lfi - the local function as computed by ADIC/ADIFOR 353 354 Level: intermediate 355 356 .keywords: distributed array, refine 357 358 .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(), 359 DASetLocalJacobian(), DASetLocalFunctionib() 360 M*/ 361 362 #undef __FUNCT__ 363 #define __FUNCT__ "DASetLocalAdicFunctionib_Private" 364 PetscErrorCode DASetLocalAdicFunctionib_Private(DA da,PetscErrorCode (*ad_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*)) 365 { 366 DM_DA *dd = (DM_DA*)da->data; 367 PetscFunctionBegin; 368 PetscValidHeaderSpecific(da,DM_CLASSID,1); 369 dd->adic_lfib = ad_lfi; 370 PetscFunctionReturn(0); 371 } 372 373 /*MC 374 DASetLocalAdicMFFunctionib - Caches in a DA a block local functioni computed by ADIC/ADIFOR 375 376 Synopsis: 377 PetscErrorCode DASetLocalAdicFunctionib(DA da,int (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*) 378 379 Logically Collective on DA 380 381 Input Parameter: 382 + da - initial distributed array 383 - admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR 384 385 Level: intermediate 386 387 .keywords: distributed array, refine 388 389 .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(), 390 DASetLocalJacobian(), DASetLocalFunctionib() 391 M*/ 392 393 #undef __FUNCT__ 394 #define __FUNCT__ "DASetLocalAdicMFFunctionib_Private" 395 PetscErrorCode DASetLocalAdicMFFunctionib_Private(DA da,PetscErrorCode (*admf_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*)) 396 { 397 DM_DA *dd = (DM_DA*)da->data; 398 PetscFunctionBegin; 399 PetscValidHeaderSpecific(da,DM_CLASSID,1); 400 dd->adicmf_lfib = admf_lfi; 401 PetscFunctionReturn(0); 402 } 403 404 /*MC 405 DASetLocalAdicMFFunction - Caches in a DA a local function computed by ADIC/ADIFOR 406 407 Synopsis: 408 PetscErrorCode DASetLocalAdicMFFunction(DA da,DALocalFunction1 ad_lf) 409 410 Logically Collective on DA 411 412 Input Parameter: 413 + da - initial distributed array 414 - ad_lf - the local function as computed by ADIC/ADIFOR 415 416 Level: intermediate 417 418 .keywords: distributed array, refine 419 420 .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(), 421 DASetLocalJacobian() 422 M*/ 423 424 #undef __FUNCT__ 425 #define __FUNCT__ "DASetLocalAdicMFFunction_Private" 426 PetscErrorCode DASetLocalAdicMFFunction_Private(DA da,DALocalFunction1 ad_lf) 427 { 428 DM_DA *dd = (DM_DA*)da->data; 429 PetscFunctionBegin; 430 PetscValidHeaderSpecific(da,DM_CLASSID,1); 431 dd->adicmf_lf = ad_lf; 432 PetscFunctionReturn(0); 433 } 434 435 /*@C 436 DASetLocalJacobian - Caches in a DA a local Jacobian computation function 437 438 Logically Collective on DA 439 440 441 Input Parameter: 442 + da - initial distributed array 443 - lj - the local Jacobian 444 445 Level: intermediate 446 447 Notes: The routine SNESDAFormFunction() uses this the cached function to evaluate the user provided function. 448 449 .keywords: distributed array, refine 450 451 .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction() 452 @*/ 453 #undef __FUNCT__ 454 #define __FUNCT__ "DASetLocalJacobian" 455 PetscErrorCode PETSCDM_DLLEXPORT DASetLocalJacobian(DA da,DALocalFunction1 lj) 456 { 457 DM_DA *dd = (DM_DA*)da->data; 458 PetscFunctionBegin; 459 PetscValidHeaderSpecific(da,DM_CLASSID,1); 460 dd->lj = lj; 461 PetscFunctionReturn(0); 462 } 463 464 #undef __FUNCT__ 465 #define __FUNCT__ "DAGetLocalFunction" 466 /*@C 467 DAGetLocalFunction - Gets from a DA a local function and its ADIC/ADIFOR Jacobian 468 469 Note Collective 470 471 Input Parameter: 472 . da - initial distributed array 473 474 Output Parameter: 475 . lf - the local function 476 477 Level: intermediate 478 479 .keywords: distributed array, refine 480 481 .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalJacobian(), DASetLocalFunction() 482 @*/ 483 PetscErrorCode PETSCDM_DLLEXPORT DAGetLocalFunction(DA da,DALocalFunction1 *lf) 484 { 485 DM_DA *dd = (DM_DA*)da->data; 486 PetscFunctionBegin; 487 PetscValidHeaderSpecific(da,DM_CLASSID,1); 488 if (lf) *lf = dd->lf; 489 PetscFunctionReturn(0); 490 } 491 492 #undef __FUNCT__ 493 #define __FUNCT__ "DAGetLocalJacobian" 494 /*@C 495 DAGetLocalJacobian - Gets from a DA a local jacobian 496 497 Not Collective 498 499 Input Parameter: 500 . da - initial distributed array 501 502 Output Parameter: 503 . lj - the local jacobian 504 505 Level: intermediate 506 507 .keywords: distributed array, refine 508 509 .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalJacobian() 510 @*/ 511 PetscErrorCode PETSCDM_DLLEXPORT DAGetLocalJacobian(DA da,DALocalFunction1 *lj) 512 { 513 DM_DA *dd = (DM_DA*)da->data; 514 PetscFunctionBegin; 515 PetscValidHeaderSpecific(da,DM_CLASSID,1); 516 if (lj) *lj = dd->lj; 517 PetscFunctionReturn(0); 518 } 519 520 #undef __FUNCT__ 521 #define __FUNCT__ "DAFormFunction" 522 /*@ 523 DAFormFunction - Evaluates a user provided function on each processor that 524 share a DA 525 526 Input Parameters: 527 + da - the DA that defines the grid 528 . vu - input vector 529 . vfu - output vector 530 - w - any user data 531 532 Notes: Does NOT do ghost updates on vu upon entry 533 534 This should eventually replace DAFormFunction1 535 536 Level: advanced 537 538 .seealso: DAComputeJacobian1WithAdic() 539 540 @*/ 541 PetscErrorCode PETSCDM_DLLEXPORT DAFormFunction(DA da,PetscErrorCode (*lf)(void),Vec vu,Vec vfu,void *w) 542 { 543 PetscErrorCode ierr; 544 void *u,*fu; 545 DALocalInfo info; 546 PetscErrorCode (*f)(DALocalInfo*,void*,void*,void*) = (PetscErrorCode (*)(DALocalInfo*,void*,void*,void*))lf; 547 548 PetscFunctionBegin; 549 ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr); 550 ierr = DAVecGetArray(da,vu,&u);CHKERRQ(ierr); 551 ierr = DAVecGetArray(da,vfu,&fu);CHKERRQ(ierr); 552 553 ierr = (*f)(&info,u,fu,w); 554 if (PetscExceptionValue(ierr)) { 555 PetscErrorCode pierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(pierr); 556 pierr = DAVecRestoreArray(da,vfu,&fu);CHKERRQ(pierr); 557 } 558 CHKERRQ(ierr); 559 560 ierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(ierr); 561 ierr = DAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr); 562 PetscFunctionReturn(0); 563 } 564 565 #undef __FUNCT__ 566 #define __FUNCT__ "DAFormFunctionLocal" 567 /*@C 568 DAFormFunctionLocal - This is a universal function evaluation routine for 569 a local DA function. 570 571 Collective on DA 572 573 Input Parameters: 574 + da - the DA context 575 . func - The local function 576 . X - input vector 577 . F - function vector 578 - ctx - A user context 579 580 Level: intermediate 581 582 .seealso: DASetLocalFunction(), DASetLocalJacobian(), DASetLocalAdicFunction(), DASetLocalAdicMFFunction(), 583 SNESSetFunction(), SNESSetJacobian() 584 585 @*/ 586 PetscErrorCode PETSCDM_DLLEXPORT DAFormFunctionLocal(DA da, DALocalFunction1 func, Vec X, Vec F, void *ctx) 587 { 588 Vec localX; 589 DALocalInfo info; 590 void *u; 591 void *fu; 592 PetscErrorCode ierr; 593 594 PetscFunctionBegin; 595 ierr = DAGetLocalVector(da,&localX);CHKERRQ(ierr); 596 /* 597 Scatter ghost points to local vector, using the 2-step process 598 DAGlobalToLocalBegin(), DAGlobalToLocalEnd(). 599 */ 600 ierr = DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 601 ierr = DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 602 ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr); 603 ierr = DAVecGetArray(da,localX,&u);CHKERRQ(ierr); 604 ierr = DAVecGetArray(da,F,&fu);CHKERRQ(ierr); 605 ierr = (*func)(&info,u,fu,ctx); 606 if (PetscExceptionValue(ierr)) { 607 PetscErrorCode pierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(pierr); 608 pierr = DAVecRestoreArray(da,F,&fu);CHKERRQ(pierr); 609 } 610 CHKERRQ(ierr); 611 ierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(ierr); 612 ierr = DAVecRestoreArray(da,F,&fu);CHKERRQ(ierr); 613 if (PetscExceptionValue(ierr)) { 614 PetscErrorCode pierr = DARestoreLocalVector(da,&localX);CHKERRQ(pierr); 615 } 616 CHKERRQ(ierr); 617 ierr = DARestoreLocalVector(da,&localX);CHKERRQ(ierr); 618 PetscFunctionReturn(0); 619 } 620 621 #undef __FUNCT__ 622 #define __FUNCT__ "DAFormFunctionLocalGhost" 623 /*@C 624 DAFormFunctionLocalGhost - This is a universal function evaluation routine for 625 a local DA function, but the ghost values of the output are communicated and added. 626 627 Collective on DA 628 629 Input Parameters: 630 + da - the DA context 631 . func - The local function 632 . X - input vector 633 . F - function vector 634 - ctx - A user context 635 636 Level: intermediate 637 638 .seealso: DASetLocalFunction(), DASetLocalJacobian(), DASetLocalAdicFunction(), DASetLocalAdicMFFunction(), 639 SNESSetFunction(), SNESSetJacobian() 640 641 @*/ 642 PetscErrorCode PETSCDM_DLLEXPORT DAFormFunctionLocalGhost(DA da, DALocalFunction1 func, Vec X, Vec F, void *ctx) 643 { 644 Vec localX, localF; 645 DALocalInfo info; 646 void *u; 647 void *fu; 648 PetscErrorCode ierr; 649 650 PetscFunctionBegin; 651 ierr = DAGetLocalVector(da,&localX);CHKERRQ(ierr); 652 ierr = DAGetLocalVector(da,&localF);CHKERRQ(ierr); 653 /* 654 Scatter ghost points to local vector, using the 2-step process 655 DAGlobalToLocalBegin(), DAGlobalToLocalEnd(). 656 */ 657 ierr = DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 658 ierr = DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 659 ierr = VecSet(F, 0.0);CHKERRQ(ierr); 660 ierr = VecSet(localF, 0.0);CHKERRQ(ierr); 661 ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr); 662 ierr = DAVecGetArray(da,localX,&u);CHKERRQ(ierr); 663 ierr = DAVecGetArray(da,localF,&fu);CHKERRQ(ierr); 664 ierr = (*func)(&info,u,fu,ctx); 665 if (PetscExceptionValue(ierr)) { 666 PetscErrorCode pierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(pierr); 667 pierr = DAVecRestoreArray(da,localF,&fu);CHKERRQ(pierr); 668 } 669 CHKERRQ(ierr); 670 ierr = DALocalToGlobalBegin(da,localF,F);CHKERRQ(ierr); 671 ierr = DALocalToGlobalEnd(da,localF,F);CHKERRQ(ierr); 672 ierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(ierr); 673 ierr = DAVecRestoreArray(da,localF,&fu);CHKERRQ(ierr); 674 if (PetscExceptionValue(ierr)) { 675 PetscErrorCode pierr = DARestoreLocalVector(da,&localX);CHKERRQ(pierr); 676 ierr = DARestoreLocalVector(da,&localF);CHKERRQ(ierr); 677 } 678 CHKERRQ(ierr); 679 ierr = DARestoreLocalVector(da,&localX);CHKERRQ(ierr); 680 ierr = DARestoreLocalVector(da,&localF);CHKERRQ(ierr); 681 PetscFunctionReturn(0); 682 } 683 684 #undef __FUNCT__ 685 #define __FUNCT__ "DAFormFunction1" 686 /*@ 687 DAFormFunction1 - Evaluates a user provided function on each processor that 688 share a DA 689 690 Input Parameters: 691 + da - the DA that defines the grid 692 . vu - input vector 693 . vfu - output vector 694 - w - any user data 695 696 Notes: Does NOT do ghost updates on vu upon entry 697 698 Level: advanced 699 700 .seealso: DAComputeJacobian1WithAdic() 701 702 @*/ 703 PetscErrorCode PETSCDM_DLLEXPORT DAFormFunction1(DA da,Vec vu,Vec vfu,void *w) 704 { 705 PetscErrorCode ierr; 706 void *u,*fu; 707 DALocalInfo info; 708 DM_DA *dd = (DM_DA*)da->data; 709 710 PetscFunctionBegin; 711 712 ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr); 713 ierr = DAVecGetArray(da,vu,&u);CHKERRQ(ierr); 714 ierr = DAVecGetArray(da,vfu,&fu);CHKERRQ(ierr); 715 716 CHKMEMQ; 717 ierr = (*dd->lf)(&info,u,fu,w); 718 if (PetscExceptionValue(ierr)) { 719 PetscErrorCode pierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(pierr); 720 pierr = DAVecRestoreArray(da,vfu,&fu);CHKERRQ(pierr); 721 } 722 CHKERRQ(ierr); 723 CHKMEMQ; 724 725 ierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(ierr); 726 ierr = DAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr); 727 PetscFunctionReturn(0); 728 } 729 730 #undef __FUNCT__ 731 #define __FUNCT__ "DAFormFunctioniTest1" 732 PetscErrorCode PETSCDM_DLLEXPORT DAFormFunctioniTest1(DA da,void *w) 733 { 734 Vec vu,fu,fui; 735 PetscErrorCode ierr; 736 PetscInt i,n; 737 PetscScalar *ui; 738 PetscRandom rnd; 739 PetscReal norm; 740 741 PetscFunctionBegin; 742 ierr = DAGetLocalVector(da,&vu);CHKERRQ(ierr); 743 ierr = PetscRandomCreate(PETSC_COMM_SELF,&rnd);CHKERRQ(ierr); 744 ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr); 745 ierr = VecSetRandom(vu,rnd);CHKERRQ(ierr); 746 ierr = PetscRandomDestroy(rnd);CHKERRQ(ierr); 747 748 ierr = DAGetGlobalVector(da,&fu);CHKERRQ(ierr); 749 ierr = DAGetGlobalVector(da,&fui);CHKERRQ(ierr); 750 751 ierr = DAFormFunction1(da,vu,fu,w);CHKERRQ(ierr); 752 753 ierr = VecGetArray(fui,&ui);CHKERRQ(ierr); 754 ierr = VecGetLocalSize(fui,&n);CHKERRQ(ierr); 755 for (i=0; i<n; i++) { 756 ierr = DAFormFunctioni1(da,i,vu,ui+i,w);CHKERRQ(ierr); 757 } 758 ierr = VecRestoreArray(fui,&ui);CHKERRQ(ierr); 759 760 ierr = VecAXPY(fui,-1.0,fu);CHKERRQ(ierr); 761 ierr = VecNorm(fui,NORM_2,&norm);CHKERRQ(ierr); 762 ierr = PetscPrintf(((PetscObject)da)->comm,"Norm of difference in vectors %G\n",norm);CHKERRQ(ierr); 763 ierr = VecView(fu,0);CHKERRQ(ierr); 764 ierr = VecView(fui,0);CHKERRQ(ierr); 765 766 ierr = DARestoreLocalVector(da,&vu);CHKERRQ(ierr); 767 ierr = DARestoreGlobalVector(da,&fu);CHKERRQ(ierr); 768 ierr = DARestoreGlobalVector(da,&fui);CHKERRQ(ierr); 769 PetscFunctionReturn(0); 770 } 771 772 #undef __FUNCT__ 773 #define __FUNCT__ "DAFormFunctioni1" 774 /*@ 775 DAFormFunctioni1 - Evaluates a user provided point-wise function 776 777 Input Parameters: 778 + da - the DA that defines the grid 779 . i - the component of the function we wish to compute (must be local) 780 . vu - input vector 781 . vfu - output value 782 - w - any user data 783 784 Notes: Does NOT do ghost updates on vu upon entry 785 786 Level: advanced 787 788 .seealso: DAComputeJacobian1WithAdic() 789 790 @*/ 791 PetscErrorCode PETSCDM_DLLEXPORT DAFormFunctioni1(DA da,PetscInt i,Vec vu,PetscScalar *vfu,void *w) 792 { 793 PetscErrorCode ierr; 794 void *u; 795 DALocalInfo info; 796 MatStencil stencil; 797 DM_DA *dd = (DM_DA*)da->data; 798 799 PetscFunctionBegin; 800 801 ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr); 802 ierr = DAVecGetArray(da,vu,&u);CHKERRQ(ierr); 803 804 /* figure out stencil value from i */ 805 stencil.c = i % info.dof; 806 stencil.i = (i % (info.xm*info.dof))/info.dof; 807 stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof); 808 stencil.k = i/(info.xm*info.ym*info.dof); 809 810 ierr = (*dd->lfi)(&info,&stencil,u,vfu,w);CHKERRQ(ierr); 811 812 ierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(ierr); 813 PetscFunctionReturn(0); 814 } 815 816 #undef __FUNCT__ 817 #define __FUNCT__ "DAFormFunctionib1" 818 /*@ 819 DAFormFunctionib1 - Evaluates a user provided point-block function 820 821 Input Parameters: 822 + da - the DA that defines the grid 823 . i - the component of the function we wish to compute (must be local) 824 . vu - input vector 825 . vfu - output value 826 - w - any user data 827 828 Notes: Does NOT do ghost updates on vu upon entry 829 830 Level: advanced 831 832 .seealso: DAComputeJacobian1WithAdic() 833 834 @*/ 835 PetscErrorCode PETSCDM_DLLEXPORT DAFormFunctionib1(DA da,PetscInt i,Vec vu,PetscScalar *vfu,void *w) 836 { 837 PetscErrorCode ierr; 838 void *u; 839 DALocalInfo info; 840 MatStencil stencil; 841 DM_DA *dd = (DM_DA*)da->data; 842 843 PetscFunctionBegin; 844 ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr); 845 ierr = DAVecGetArray(da,vu,&u);CHKERRQ(ierr); 846 847 /* figure out stencil value from i */ 848 stencil.c = i % info.dof; 849 if (stencil.c) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Point-block functions can only be called for the entire block"); 850 stencil.i = (i % (info.xm*info.dof))/info.dof; 851 stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof); 852 stencil.k = i/(info.xm*info.ym*info.dof); 853 854 ierr = (*dd->lfib)(&info,&stencil,u,vfu,w);CHKERRQ(ierr); 855 856 ierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(ierr); 857 PetscFunctionReturn(0); 858 } 859 860 #if defined(new) 861 #undef __FUNCT__ 862 #define __FUNCT__ "DAGetDiagonal_MFFD" 863 /* 864 DAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local 865 function lives on a DA 866 867 y ~= (F(u + ha) - F(u))/h, 868 where F = nonlinear function, as set by SNESSetFunction() 869 u = current iterate 870 h = difference interval 871 */ 872 PetscErrorCode DAGetDiagonal_MFFD(DA da,Vec U,Vec a) 873 { 874 PetscScalar h,*aa,*ww,v; 875 PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 876 PetscErrorCode ierr; 877 PetscInt gI,nI; 878 MatStencil stencil; 879 DALocalInfo info; 880 881 PetscFunctionBegin; 882 ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr); 883 ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr); 884 885 ierr = VecGetArray(U,&ww);CHKERRQ(ierr); 886 ierr = VecGetArray(a,&aa);CHKERRQ(ierr); 887 888 nI = 0; 889 h = ww[gI]; 890 if (h == 0.0) h = 1.0; 891 #if !defined(PETSC_USE_COMPLEX) 892 if (h < umin && h >= 0.0) h = umin; 893 else if (h < 0.0 && h > -umin) h = -umin; 894 #else 895 if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 896 else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 897 #endif 898 h *= epsilon; 899 900 ww[gI] += h; 901 ierr = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr); 902 aa[nI] = (v - aa[nI])/h; 903 ww[gI] -= h; 904 nI++; 905 } 906 ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr); 907 ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr); 908 PetscFunctionReturn(0); 909 } 910 #endif 911 912 #if defined(PETSC_HAVE_ADIC) 913 EXTERN_C_BEGIN 914 #include "adic/ad_utils.h" 915 EXTERN_C_END 916 917 #undef __FUNCT__ 918 #define __FUNCT__ "DAComputeJacobian1WithAdic" 919 /*@C 920 DAComputeJacobian1WithAdic - Evaluates a adiC provided Jacobian function on each processor that 921 share a DA 922 923 Input Parameters: 924 + da - the DA that defines the grid 925 . vu - input vector (ghosted) 926 . J - output matrix 927 - w - any user data 928 929 Level: advanced 930 931 Notes: Does NOT do ghost updates on vu upon entry 932 933 .seealso: DAFormFunction1() 934 935 @*/ 936 PetscErrorCode PETSCDM_DLLEXPORT DAComputeJacobian1WithAdic(DA da,Vec vu,Mat J,void *w) 937 { 938 PetscErrorCode ierr; 939 PetscInt gtdof,tdof; 940 PetscScalar *ustart; 941 DALocalInfo info; 942 void *ad_u,*ad_f,*ad_ustart,*ad_fstart; 943 ISColoring iscoloring; 944 945 PetscFunctionBegin; 946 ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr); 947 948 PetscADResetIndep(); 949 950 /* get space for derivative objects. */ 951 ierr = DAGetAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,>dof);CHKERRQ(ierr); 952 ierr = DAGetAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr); 953 ierr = VecGetArray(vu,&ustart);CHKERRQ(ierr); 954 ierr = DAGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr); 955 956 PetscADSetValueAndColor(ad_ustart,gtdof,iscoloring->colors,ustart); 957 958 ierr = VecRestoreArray(vu,&ustart);CHKERRQ(ierr); 959 ierr = ISColoringDestroy(iscoloring);CHKERRQ(ierr); 960 ierr = PetscADIncrementTotalGradSize(iscoloring->n);CHKERRQ(ierr); 961 PetscADSetIndepDone(); 962 963 ierr = PetscLogEventBegin(DA_LocalADFunction,0,0,0,0);CHKERRQ(ierr); 964 ierr = (*dd->adic_lf)(&info,ad_u,ad_f,w);CHKERRQ(ierr); 965 ierr = PetscLogEventEnd(DA_LocalADFunction,0,0,0,0);CHKERRQ(ierr); 966 967 /* stick the values into the matrix */ 968 ierr = MatSetValuesAdic(J,(PetscScalar**)ad_fstart);CHKERRQ(ierr); 969 970 /* return space for derivative objects. */ 971 ierr = DARestoreAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,>dof);CHKERRQ(ierr); 972 ierr = DARestoreAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr); 973 PetscFunctionReturn(0); 974 } 975 976 #undef __FUNCT__ 977 #define __FUNCT__ "DAMultiplyByJacobian1WithAdic" 978 /*@C 979 DAMultiplyByJacobian1WithAdic - Applies an ADIC-provided Jacobian function to a vector on 980 each processor that shares a DA. 981 982 Input Parameters: 983 + da - the DA that defines the grid 984 . vu - Jacobian is computed at this point (ghosted) 985 . v - product is done on this vector (ghosted) 986 . fu - output vector = J(vu)*v (not ghosted) 987 - w - any user data 988 989 Notes: 990 This routine does NOT do ghost updates on vu upon entry. 991 992 Level: advanced 993 994 .seealso: DAFormFunction1() 995 996 @*/ 997 PetscErrorCode PETSCDM_DLLEXPORT DAMultiplyByJacobian1WithAdic(DA da,Vec vu,Vec v,Vec f,void *w) 998 { 999 PetscErrorCode ierr; 1000 PetscInt i,gtdof,tdof; 1001 PetscScalar *avu,*av,*af,*ad_vustart,*ad_fstart; 1002 DALocalInfo info; 1003 void *ad_vu,*ad_f; 1004 1005 PetscFunctionBegin; 1006 ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr); 1007 1008 /* get space for derivative objects. */ 1009 ierr = DAGetAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,>dof);CHKERRQ(ierr); 1010 ierr = DAGetAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr); 1011 1012 /* copy input vector into derivative object */ 1013 ierr = VecGetArray(vu,&avu);CHKERRQ(ierr); 1014 ierr = VecGetArray(v,&av);CHKERRQ(ierr); 1015 for (i=0; i<gtdof; i++) { 1016 ad_vustart[2*i] = avu[i]; 1017 ad_vustart[2*i+1] = av[i]; 1018 } 1019 ierr = VecRestoreArray(vu,&avu);CHKERRQ(ierr); 1020 ierr = VecRestoreArray(v,&av);CHKERRQ(ierr); 1021 1022 PetscADResetIndep(); 1023 ierr = PetscADIncrementTotalGradSize(1);CHKERRQ(ierr); 1024 PetscADSetIndepDone(); 1025 1026 ierr = (*dd->adicmf_lf)(&info,ad_vu,ad_f,w);CHKERRQ(ierr); 1027 1028 /* stick the values into the vector */ 1029 ierr = VecGetArray(f,&af);CHKERRQ(ierr); 1030 for (i=0; i<tdof; i++) { 1031 af[i] = ad_fstart[2*i+1]; 1032 } 1033 ierr = VecRestoreArray(f,&af);CHKERRQ(ierr); 1034 1035 /* return space for derivative objects. */ 1036 ierr = DARestoreAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,>dof);CHKERRQ(ierr); 1037 ierr = DARestoreAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr); 1038 PetscFunctionReturn(0); 1039 } 1040 #endif 1041 1042 #undef __FUNCT__ 1043 #define __FUNCT__ "DAComputeJacobian1" 1044 /*@ 1045 DAComputeJacobian1 - Evaluates a local Jacobian function on each processor that 1046 share a DA 1047 1048 Input Parameters: 1049 + da - the DA that defines the grid 1050 . vu - input vector (ghosted) 1051 . J - output matrix 1052 - w - any user data 1053 1054 Notes: Does NOT do ghost updates on vu upon entry 1055 1056 Level: advanced 1057 1058 .seealso: DAFormFunction1() 1059 1060 @*/ 1061 PetscErrorCode PETSCDM_DLLEXPORT DAComputeJacobian1(DA da,Vec vu,Mat J,void *w) 1062 { 1063 PetscErrorCode ierr; 1064 void *u; 1065 DALocalInfo info; 1066 DM_DA *dd = (DM_DA*)da->data; 1067 1068 PetscFunctionBegin; 1069 ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr); 1070 ierr = DAVecGetArray(da,vu,&u);CHKERRQ(ierr); 1071 ierr = (*dd->lj)(&info,u,J,w);CHKERRQ(ierr); 1072 ierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(ierr); 1073 PetscFunctionReturn(0); 1074 } 1075 1076 1077 #undef __FUNCT__ 1078 #define __FUNCT__ "DAComputeJacobian1WithAdifor" 1079 /* 1080 DAComputeJacobian1WithAdifor - Evaluates a ADIFOR provided Jacobian local function on each processor that 1081 share a DA 1082 1083 Input Parameters: 1084 + da - the DA that defines the grid 1085 . vu - input vector (ghosted) 1086 . J - output matrix 1087 - w - any user data 1088 1089 Notes: Does NOT do ghost updates on vu upon entry 1090 1091 .seealso: DAFormFunction1() 1092 1093 */ 1094 PetscErrorCode PETSCDM_DLLEXPORT DAComputeJacobian1WithAdifor(DA da,Vec vu,Mat J,void *w) 1095 { 1096 PetscErrorCode ierr; 1097 PetscInt i,Nc,N; 1098 ISColoringValue *color; 1099 DALocalInfo info; 1100 PetscScalar *u,*g_u,*g_f,*f = 0,*p_u; 1101 ISColoring iscoloring; 1102 DM_DA *dd = (DM_DA*)da->data; 1103 void (*lf)(PetscInt*,DALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*) = 1104 (void (*)(PetscInt*,DALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*))*dd->adifor_lf; 1105 1106 PetscFunctionBegin; 1107 ierr = DAGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr); 1108 Nc = iscoloring->n; 1109 ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr); 1110 N = info.gxm*info.gym*info.gzm*info.dof; 1111 1112 /* get space for derivative objects. */ 1113 ierr = PetscMalloc(Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar),&g_u);CHKERRQ(ierr); 1114 ierr = PetscMemzero(g_u,Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar));CHKERRQ(ierr); 1115 p_u = g_u; 1116 color = iscoloring->colors; 1117 for (i=0; i<N; i++) { 1118 p_u[*color++] = 1.0; 1119 p_u += Nc; 1120 } 1121 ierr = ISColoringDestroy(iscoloring);CHKERRQ(ierr); 1122 ierr = PetscMalloc2(Nc*info.xm*info.ym*info.zm*info.dof,PetscScalar,&g_f,info.xm*info.ym*info.zm*info.dof,PetscScalar,&f);CHKERRQ(ierr); 1123 1124 /* Seed the input array g_u with coloring information */ 1125 1126 ierr = VecGetArray(vu,&u);CHKERRQ(ierr); 1127 (lf)(&Nc,&info,u,g_u,&Nc,f,g_f,&Nc,w,&ierr);CHKERRQ(ierr); 1128 ierr = VecRestoreArray(vu,&u);CHKERRQ(ierr); 1129 1130 /* stick the values into the matrix */ 1131 /* PetscScalarView(Nc*info.xm*info.ym,g_f,0); */ 1132 ierr = MatSetValuesAdifor(J,Nc,g_f);CHKERRQ(ierr); 1133 1134 /* return space for derivative objects. */ 1135 ierr = PetscFree(g_u);CHKERRQ(ierr); 1136 ierr = PetscFree2(g_f,f);CHKERRQ(ierr); 1137 PetscFunctionReturn(0); 1138 } 1139 1140 #undef __FUNCT__ 1141 #define __FUNCT__ "DAFormJacobianLocal" 1142 /*@C 1143 DAFormjacobianLocal - This is a universal Jacobian evaluation routine for 1144 a local DA function. 1145 1146 Collective on DA 1147 1148 Input Parameters: 1149 + da - the DA context 1150 . func - The local function 1151 . X - input vector 1152 . J - Jacobian matrix 1153 - ctx - A user context 1154 1155 Level: intermediate 1156 1157 .seealso: DASetLocalFunction(), DASetLocalJacobian(), DASetLocalAdicFunction(), DASetLocalAdicMFFunction(), 1158 SNESSetFunction(), SNESSetJacobian() 1159 1160 @*/ 1161 PetscErrorCode PETSCDM_DLLEXPORT DAFormJacobianLocal(DA da, DALocalFunction1 func, Vec X, Mat J, void *ctx) 1162 { 1163 Vec localX; 1164 DALocalInfo info; 1165 void *u; 1166 PetscErrorCode ierr; 1167 1168 PetscFunctionBegin; 1169 ierr = DAGetLocalVector(da,&localX);CHKERRQ(ierr); 1170 /* 1171 Scatter ghost points to local vector, using the 2-step process 1172 DAGlobalToLocalBegin(), DAGlobalToLocalEnd(). 1173 */ 1174 ierr = DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 1175 ierr = DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 1176 ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr); 1177 ierr = DAVecGetArray(da,localX,&u);CHKERRQ(ierr); 1178 ierr = (*func)(&info,u,J,ctx); 1179 if (PetscExceptionValue(ierr)) { 1180 PetscErrorCode pierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(pierr); 1181 } 1182 CHKERRQ(ierr); 1183 ierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(ierr); 1184 if (PetscExceptionValue(ierr)) { 1185 PetscErrorCode pierr = DARestoreLocalVector(da,&localX);CHKERRQ(pierr); 1186 } 1187 CHKERRQ(ierr); 1188 ierr = DARestoreLocalVector(da,&localX);CHKERRQ(ierr); 1189 PetscFunctionReturn(0); 1190 } 1191 1192 #undef __FUNCT__ 1193 #define __FUNCT__ "DAMultiplyByJacobian1WithAD" 1194 /*@C 1195 DAMultiplyByJacobian1WithAD - Applies a Jacobian function supplied by ADIFOR or ADIC 1196 to a vector on each processor that shares a DA. 1197 1198 Input Parameters: 1199 + da - the DA that defines the grid 1200 . vu - Jacobian is computed at this point (ghosted) 1201 . v - product is done on this vector (ghosted) 1202 . fu - output vector = J(vu)*v (not ghosted) 1203 - w - any user data 1204 1205 Notes: 1206 This routine does NOT do ghost updates on vu and v upon entry. 1207 1208 Automatically calls DAMultiplyByJacobian1WithAdifor() or DAMultiplyByJacobian1WithAdic() 1209 depending on whether DASetLocalAdicMFFunction() or DASetLocalAdiforMFFunction() was called. 1210 1211 Level: advanced 1212 1213 .seealso: DAFormFunction1(), DAMultiplyByJacobian1WithAdifor(), DAMultiplyByJacobian1WithAdic() 1214 1215 @*/ 1216 PetscErrorCode PETSCDM_DLLEXPORT DAMultiplyByJacobian1WithAD(DA da,Vec u,Vec v,Vec f,void *w) 1217 { 1218 PetscErrorCode ierr; 1219 DM_DA *dd = (DM_DA*)da->data; 1220 1221 PetscFunctionBegin; 1222 if (dd->adicmf_lf) { 1223 #if defined(PETSC_HAVE_ADIC) 1224 ierr = DAMultiplyByJacobian1WithAdic(da,u,v,f,w);CHKERRQ(ierr); 1225 #else 1226 SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP_SYS,"Requires ADIC to be installed and cannot use complex numbers"); 1227 #endif 1228 } else if (dd->adiformf_lf) { 1229 ierr = DAMultiplyByJacobian1WithAdifor(da,u,v,f,w);CHKERRQ(ierr); 1230 } else { 1231 SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ORDER,"Must call DASetLocalAdiforMFFunction() or DASetLocalAdicMFFunction() before using"); 1232 } 1233 PetscFunctionReturn(0); 1234 } 1235 1236 1237 #undef __FUNCT__ 1238 #define __FUNCT__ "DAMultiplyByJacobian1WithAdifor" 1239 /*@C 1240 DAMultiplyByJacobian1WithAdifor - Applies a ADIFOR provided Jacobian function on each processor that 1241 share a DA to a vector 1242 1243 Input Parameters: 1244 + da - the DA that defines the grid 1245 . vu - Jacobian is computed at this point (ghosted) 1246 . v - product is done on this vector (ghosted) 1247 . fu - output vector = J(vu)*v (not ghosted) 1248 - w - any user data 1249 1250 Notes: Does NOT do ghost updates on vu and v upon entry 1251 1252 Level: advanced 1253 1254 .seealso: DAFormFunction1() 1255 1256 @*/ 1257 PetscErrorCode PETSCDM_DLLEXPORT DAMultiplyByJacobian1WithAdifor(DA da,Vec u,Vec v,Vec f,void *w) 1258 { 1259 PetscErrorCode ierr; 1260 PetscScalar *au,*av,*af,*awork; 1261 Vec work; 1262 DALocalInfo info; 1263 DM_DA *dd = (DM_DA*)da->data; 1264 void (*lf)(DALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*) = 1265 (void (*)(DALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*))*dd->adiformf_lf; 1266 1267 PetscFunctionBegin; 1268 ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr); 1269 1270 ierr = DAGetGlobalVector(da,&work);CHKERRQ(ierr); 1271 ierr = VecGetArray(u,&au);CHKERRQ(ierr); 1272 ierr = VecGetArray(v,&av);CHKERRQ(ierr); 1273 ierr = VecGetArray(f,&af);CHKERRQ(ierr); 1274 ierr = VecGetArray(work,&awork);CHKERRQ(ierr); 1275 (lf)(&info,au,av,awork,af,w,&ierr);CHKERRQ(ierr); 1276 ierr = VecRestoreArray(u,&au);CHKERRQ(ierr); 1277 ierr = VecRestoreArray(v,&av);CHKERRQ(ierr); 1278 ierr = VecRestoreArray(f,&af);CHKERRQ(ierr); 1279 ierr = VecRestoreArray(work,&awork);CHKERRQ(ierr); 1280 ierr = DARestoreGlobalVector(da,&work);CHKERRQ(ierr); 1281 1282 PetscFunctionReturn(0); 1283 } 1284 1285 EXTERN_C_BEGIN 1286 #undef __FUNCT__ 1287 #define __FUNCT__ "DASetUp_2D" 1288 PetscErrorCode PETSCDM_DLLEXPORT DASetUp_2D(DA da) 1289 { 1290 DM_DA *dd = (DM_DA*)da->data; 1291 const PetscInt M = dd->M; 1292 const PetscInt N = dd->N; 1293 PetscInt m = dd->m; 1294 PetscInt n = dd->n; 1295 const PetscInt dof = dd->w; 1296 const PetscInt s = dd->s; 1297 const DAPeriodicType wrap = dd->wrap; 1298 const DAStencilType stencil_type = dd->stencil_type; 1299 PetscInt *lx = dd->lx; 1300 PetscInt *ly = dd->ly; 1301 MPI_Comm comm; 1302 PetscMPIInt rank,size; 1303 PetscInt xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,start,end; 1304 PetscInt up,down,left,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn; 1305 PetscInt xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count; 1306 PetscInt s_x,s_y; /* s proportionalized to w */ 1307 PetscInt sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0; 1308 Vec local,global; 1309 VecScatter ltog,gtol; 1310 IS to,from; 1311 PetscErrorCode ierr; 1312 1313 PetscFunctionBegin; 1314 ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 1315 1316 if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof); 1317 if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s); 1318 1319 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1320 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1321 1322 da->ops->getelements = DMGetElements_DA_2d_P1; 1323 1324 dd->dim = 2; 1325 dd->elementtype = DA_ELEMENT_P1; 1326 ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr); 1327 ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr); 1328 1329 if (m != PETSC_DECIDE) { 1330 if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m); 1331 else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size); 1332 } 1333 if (n != PETSC_DECIDE) { 1334 if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n); 1335 else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size); 1336 } 1337 1338 if (m == PETSC_DECIDE || n == PETSC_DECIDE) { 1339 if (n != PETSC_DECIDE) { 1340 m = size/n; 1341 } else if (m != PETSC_DECIDE) { 1342 n = size/m; 1343 } else { 1344 /* try for squarish distribution */ 1345 m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)N))); 1346 if (!m) m = 1; 1347 while (m > 0) { 1348 n = size/m; 1349 if (m*n == size) break; 1350 m--; 1351 } 1352 if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;} 1353 } 1354 if (m*n != size) SETERRQ(comm,PETSC_ERR_PLIB,"Unable to create partition, check the size of the communicator and input m and n "); 1355 } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition"); 1356 1357 if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m); 1358 if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n); 1359 1360 /* 1361 Determine locally owned region 1362 xs is the first local node number, x is the number of local nodes 1363 */ 1364 if (!lx) { 1365 ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr); 1366 lx = dd->lx; 1367 for (i=0; i<m; i++) { 1368 lx[i] = M/m + ((M % m) > i); 1369 } 1370 } 1371 x = lx[rank % m]; 1372 xs = 0; 1373 for (i=0; i<(rank % m); i++) { 1374 xs += lx[i]; 1375 } 1376 #if defined(PETSC_USE_DEBUG) 1377 left = xs; 1378 for (i=(rank % m); i<m; i++) { 1379 left += lx[i]; 1380 } 1381 if (left != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %D %D",left,M); 1382 #endif 1383 1384 /* 1385 Determine locally owned region 1386 ys is the first local node number, y is the number of local nodes 1387 */ 1388 if (!ly) { 1389 ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr); 1390 ly = dd->ly; 1391 for (i=0; i<n; i++) { 1392 ly[i] = N/n + ((N % n) > i); 1393 } 1394 } 1395 y = ly[rank/m]; 1396 ys = 0; 1397 for (i=0; i<(rank/m); i++) { 1398 ys += ly[i]; 1399 } 1400 #if defined(PETSC_USE_DEBUG) 1401 left = ys; 1402 for (i=(rank/m); i<n; i++) { 1403 left += ly[i]; 1404 } 1405 if (left != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %D %D",left,N); 1406 #endif 1407 1408 if (x < s) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s); 1409 if (y < s) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s); 1410 xe = xs + x; 1411 ye = ys + y; 1412 1413 /* determine ghost region */ 1414 /* Assume No Periodicity */ 1415 if (xs-s > 0) Xs = xs - s; else Xs = 0; 1416 if (ys-s > 0) Ys = ys - s; else Ys = 0; 1417 if (xe+s <= M) Xe = xe + s; else Xe = M; 1418 if (ye+s <= N) Ye = ye + s; else Ye = N; 1419 1420 /* X Periodic */ 1421 if (DAXPeriodic(wrap)){ 1422 Xs = xs - s; 1423 Xe = xe + s; 1424 } 1425 1426 /* Y Periodic */ 1427 if (DAYPeriodic(wrap)){ 1428 Ys = ys - s; 1429 Ye = ye + s; 1430 } 1431 1432 /* Resize all X parameters to reflect w */ 1433 x *= dof; 1434 xs *= dof; 1435 xe *= dof; 1436 Xs *= dof; 1437 Xe *= dof; 1438 s_x = s*dof; 1439 s_y = s; 1440 1441 /* determine starting point of each processor */ 1442 nn = x*y; 1443 ierr = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr); 1444 ierr = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr); 1445 bases[0] = 0; 1446 for (i=1; i<=size; i++) { 1447 bases[i] = ldims[i-1]; 1448 } 1449 for (i=1; i<=size; i++) { 1450 bases[i] += bases[i-1]; 1451 } 1452 1453 /* allocate the base parallel and sequential vectors */ 1454 dd->Nlocal = x*y; 1455 ierr = VecCreateMPIWithArray(comm,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr); 1456 ierr = VecSetBlockSize(global,dof);CHKERRQ(ierr); 1457 dd->nlocal = (Xe-Xs)*(Ye-Ys); 1458 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dd->nlocal,0,&local);CHKERRQ(ierr); 1459 ierr = VecSetBlockSize(local,dof);CHKERRQ(ierr); 1460 1461 1462 /* generate appropriate vector scatters */ 1463 /* local to global inserts non-ghost point region into global */ 1464 ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr); 1465 ierr = ISCreateStride(comm,x*y,start,1,&to);CHKERRQ(ierr); 1466 1467 left = xs - Xs; down = ys - Ys; up = down + y; 1468 ierr = PetscMalloc(x*(up - down)*sizeof(PetscInt),&idx);CHKERRQ(ierr); 1469 count = 0; 1470 for (i=down; i<up; i++) { 1471 for (j=0; j<x/dof; j++) { 1472 idx[count++] = (left + i*(Xe-Xs) + j*dof)/dof; 1473 } 1474 } 1475 ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr); 1476 1477 ierr = VecScatterCreate(local,from,global,to,<og);CHKERRQ(ierr); 1478 ierr = PetscLogObjectParent(dd,ltog);CHKERRQ(ierr); 1479 ierr = ISDestroy(from);CHKERRQ(ierr); 1480 ierr = ISDestroy(to);CHKERRQ(ierr); 1481 1482 /* global to local must include ghost points */ 1483 if (stencil_type == DA_STENCIL_BOX) { 1484 ierr = ISCreateStride(comm,(Xe-Xs)*(Ye-Ys),0,1,&to);CHKERRQ(ierr); 1485 } else { 1486 /* must drop into cross shape region */ 1487 /* ---------| 1488 | top | 1489 |--- ---| 1490 | middle | 1491 | | 1492 ---- ---- 1493 | bottom | 1494 ----------- 1495 Xs xs xe Xe */ 1496 /* bottom */ 1497 left = xs - Xs; down = ys - Ys; up = down + y; 1498 count = down*(xe-xs) + (up-down)*(Xe-Xs) + (Ye-Ys-up)*(xe-xs); 1499 ierr = PetscMalloc(count*sizeof(PetscInt)/dof,&idx);CHKERRQ(ierr); 1500 count = 0; 1501 for (i=0; i<down; i++) { 1502 for (j=0; j<xe-xs; j += dof) { 1503 idx[count++] = left + i*(Xe-Xs) + j; 1504 } 1505 } 1506 /* middle */ 1507 for (i=down; i<up; i++) { 1508 for (j=0; j<Xe-Xs; j += dof) { 1509 idx[count++] = i*(Xe-Xs) + j; 1510 } 1511 } 1512 /* top */ 1513 for (i=up; i<Ye-Ys; i++) { 1514 for (j=0; j<xe-xs; j += dof) { 1515 idx[count++] = (left + i*(Xe-Xs) + j); 1516 } 1517 } 1518 for (i=0; i<count; i++) idx[i] = idx[i]/dof; 1519 ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 1520 } 1521 1522 1523 /* determine who lies on each side of us stored in n6 n7 n8 1524 n3 n5 1525 n0 n1 n2 1526 */ 1527 1528 /* Assume the Non-Periodic Case */ 1529 n1 = rank - m; 1530 if (rank % m) { 1531 n0 = n1 - 1; 1532 } else { 1533 n0 = -1; 1534 } 1535 if ((rank+1) % m) { 1536 n2 = n1 + 1; 1537 n5 = rank + 1; 1538 n8 = rank + m + 1; if (n8 >= m*n) n8 = -1; 1539 } else { 1540 n2 = -1; n5 = -1; n8 = -1; 1541 } 1542 if (rank % m) { 1543 n3 = rank - 1; 1544 n6 = n3 + m; if (n6 >= m*n) n6 = -1; 1545 } else { 1546 n3 = -1; n6 = -1; 1547 } 1548 n7 = rank + m; if (n7 >= m*n) n7 = -1; 1549 1550 1551 /* Modify for Periodic Cases */ 1552 if (wrap == DA_YPERIODIC) { /* Handle Top and Bottom Sides */ 1553 if (n1 < 0) n1 = rank + m * (n-1); 1554 if (n7 < 0) n7 = rank - m * (n-1); 1555 if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 1556 if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1; 1557 if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 1558 if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1; 1559 } else if (wrap == DA_XPERIODIC) { /* Handle Left and Right Sides */ 1560 if (n3 < 0) n3 = rank + (m-1); 1561 if (n5 < 0) n5 = rank - (m-1); 1562 if ((n1 >= 0) && (n0 < 0)) n0 = rank-1; 1563 if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1; 1564 if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1; 1565 if ((n7 >= 0) && (n8 < 0)) n8 = rank+1; 1566 } else if (wrap == DA_XYPERIODIC) { 1567 1568 /* Handle all four corners */ 1569 if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1; 1570 if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0; 1571 if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m; 1572 if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1; 1573 1574 /* Handle Top and Bottom Sides */ 1575 if (n1 < 0) n1 = rank + m * (n-1); 1576 if (n7 < 0) n7 = rank - m * (n-1); 1577 if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 1578 if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1; 1579 if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 1580 if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1; 1581 1582 /* Handle Left and Right Sides */ 1583 if (n3 < 0) n3 = rank + (m-1); 1584 if (n5 < 0) n5 = rank - (m-1); 1585 if ((n1 >= 0) && (n0 < 0)) n0 = rank-1; 1586 if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1; 1587 if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1; 1588 if ((n7 >= 0) && (n8 < 0)) n8 = rank+1; 1589 } 1590 ierr = PetscMalloc(9*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr); 1591 dd->neighbors[0] = n0; 1592 dd->neighbors[1] = n1; 1593 dd->neighbors[2] = n2; 1594 dd->neighbors[3] = n3; 1595 dd->neighbors[4] = rank; 1596 dd->neighbors[5] = n5; 1597 dd->neighbors[6] = n6; 1598 dd->neighbors[7] = n7; 1599 dd->neighbors[8] = n8; 1600 1601 if (stencil_type == DA_STENCIL_STAR) { 1602 /* save corner processor numbers */ 1603 sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8; 1604 n0 = n2 = n6 = n8 = -1; 1605 } 1606 1607 ierr = PetscMalloc((x+2*s_x)*(y+2*s_y)*sizeof(PetscInt),&idx);CHKERRQ(ierr); 1608 ierr = PetscLogObjectMemory(da,(x+2*s_x)*(y+2*s_y)*sizeof(PetscInt));CHKERRQ(ierr); 1609 nn = 0; 1610 1611 xbase = bases[rank]; 1612 for (i=1; i<=s_y; i++) { 1613 if (n0 >= 0) { /* left below */ 1614 x_t = lx[n0 % m]*dof; 1615 y_t = ly[(n0/m)]; 1616 s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x; 1617 for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1618 } 1619 if (n1 >= 0) { /* directly below */ 1620 x_t = x; 1621 y_t = ly[(n1/m)]; 1622 s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t; 1623 for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1624 } 1625 if (n2 >= 0) { /* right below */ 1626 x_t = lx[n2 % m]*dof; 1627 y_t = ly[(n2/m)]; 1628 s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t; 1629 for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1630 } 1631 } 1632 1633 for (i=0; i<y; i++) { 1634 if (n3 >= 0) { /* directly left */ 1635 x_t = lx[n3 % m]*dof; 1636 /* y_t = y; */ 1637 s_t = bases[n3] + (i+1)*x_t - s_x; 1638 for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1639 } 1640 1641 for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */ 1642 1643 if (n5 >= 0) { /* directly right */ 1644 x_t = lx[n5 % m]*dof; 1645 /* y_t = y; */ 1646 s_t = bases[n5] + (i)*x_t; 1647 for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1648 } 1649 } 1650 1651 for (i=1; i<=s_y; i++) { 1652 if (n6 >= 0) { /* left above */ 1653 x_t = lx[n6 % m]*dof; 1654 /* y_t = ly[(n6/m)]; */ 1655 s_t = bases[n6] + (i)*x_t - s_x; 1656 for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1657 } 1658 if (n7 >= 0) { /* directly above */ 1659 x_t = x; 1660 /* y_t = ly[(n7/m)]; */ 1661 s_t = bases[n7] + (i-1)*x_t; 1662 for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1663 } 1664 if (n8 >= 0) { /* right above */ 1665 x_t = lx[n8 % m]*dof; 1666 /* y_t = ly[(n8/m)]; */ 1667 s_t = bases[n8] + (i-1)*x_t; 1668 for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1669 } 1670 } 1671 1672 base = bases[rank]; 1673 { 1674 PetscInt nnn = nn/dof,*iidx; 1675 ierr = PetscMalloc(nnn*sizeof(PetscInt),&iidx);CHKERRQ(ierr); 1676 for (i=0; i<nnn; i++) { 1677 iidx[i] = idx[dof*i]/dof; 1678 } 1679 ierr = ISCreateBlock(comm,dof,nnn,iidx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr); 1680 } 1681 ierr = VecScatterCreate(global,from,local,to,>ol);CHKERRQ(ierr); 1682 ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr); 1683 ierr = ISDestroy(to);CHKERRQ(ierr); 1684 ierr = ISDestroy(from);CHKERRQ(ierr); 1685 1686 if (stencil_type == DA_STENCIL_STAR) { 1687 /* 1688 Recompute the local to global mappings, this time keeping the 1689 information about the cross corner processor numbers. 1690 */ 1691 n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8; 1692 nn = 0; 1693 xbase = bases[rank]; 1694 for (i=1; i<=s_y; i++) { 1695 if (n0 >= 0) { /* left below */ 1696 x_t = lx[n0 % m]*dof; 1697 y_t = ly[(n0/m)]; 1698 s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x; 1699 for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1700 } 1701 if (n1 >= 0) { /* directly below */ 1702 x_t = x; 1703 y_t = ly[(n1/m)]; 1704 s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t; 1705 for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1706 } 1707 if (n2 >= 0) { /* right below */ 1708 x_t = lx[n2 % m]*dof; 1709 y_t = ly[(n2/m)]; 1710 s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t; 1711 for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1712 } 1713 } 1714 1715 for (i=0; i<y; i++) { 1716 if (n3 >= 0) { /* directly left */ 1717 x_t = lx[n3 % m]*dof; 1718 /* y_t = y; */ 1719 s_t = bases[n3] + (i+1)*x_t - s_x; 1720 for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1721 } 1722 1723 for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */ 1724 1725 if (n5 >= 0) { /* directly right */ 1726 x_t = lx[n5 % m]*dof; 1727 /* y_t = y; */ 1728 s_t = bases[n5] + (i)*x_t; 1729 for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1730 } 1731 } 1732 1733 for (i=1; i<=s_y; i++) { 1734 if (n6 >= 0) { /* left above */ 1735 x_t = lx[n6 % m]*dof; 1736 /* y_t = ly[(n6/m)]; */ 1737 s_t = bases[n6] + (i)*x_t - s_x; 1738 for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1739 } 1740 if (n7 >= 0) { /* directly above */ 1741 x_t = x; 1742 /* y_t = ly[(n7/m)]; */ 1743 s_t = bases[n7] + (i-1)*x_t; 1744 for (j=0; j<x_t; j++) { idx[nn++] = s_t++;} 1745 } 1746 if (n8 >= 0) { /* right above */ 1747 x_t = lx[n8 % m]*dof; 1748 /* y_t = ly[(n8/m)]; */ 1749 s_t = bases[n8] + (i-1)*x_t; 1750 for (j=0; j<s_x; j++) { idx[nn++] = s_t++;} 1751 } 1752 } 1753 } 1754 ierr = PetscFree2(bases,ldims);CHKERRQ(ierr); 1755 1756 dd->m = m; dd->n = n; 1757 dd->xs = xs; dd->xe = xe; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1; 1758 dd->Xs = Xs; dd->Xe = Xe; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1; 1759 1760 ierr = VecDestroy(local);CHKERRQ(ierr); 1761 ierr = VecDestroy(global);CHKERRQ(ierr); 1762 1763 dd->gtol = gtol; 1764 dd->ltog = ltog; 1765 dd->idx = idx; 1766 dd->Nl = nn; 1767 dd->base = base; 1768 da->ops->view = DAView_2d; 1769 1770 /* 1771 Set the local to global ordering in the global vector, this allows use 1772 of VecSetValuesLocal(). 1773 */ 1774 ierr = ISLocalToGlobalMappingCreate(comm,nn,idx,PETSC_OWN_POINTER,&dd->ltogmap);CHKERRQ(ierr); 1775 ierr = ISLocalToGlobalMappingBlock(dd->ltogmap,dd->w,&dd->ltogmapb);CHKERRQ(ierr); 1776 ierr = PetscLogObjectParent(da,dd->ltogmap);CHKERRQ(ierr); 1777 1778 dd->ltol = PETSC_NULL; 1779 dd->ao = PETSC_NULL; 1780 1781 PetscFunctionReturn(0); 1782 } 1783 EXTERN_C_END 1784 1785 #undef __FUNCT__ 1786 #define __FUNCT__ "DACreate2d" 1787 /*@C 1788 DACreate2d - Creates an object that will manage the communication of two-dimensional 1789 regular array data that is distributed across some processors. 1790 1791 Collective on MPI_Comm 1792 1793 Input Parameters: 1794 + comm - MPI communicator 1795 . wrap - type of periodicity should the array have. 1796 Use one of DA_NONPERIODIC, DA_XPERIODIC, DA_YPERIODIC, or DA_XYPERIODIC. 1797 . stencil_type - stencil type. Use either DA_STENCIL_BOX or DA_STENCIL_STAR. 1798 . M,N - global dimension in each direction of the array (use -M and or -N to indicate that it may be set to a different value 1799 from the command line with -da_grid_x <M> -da_grid_y <N>) 1800 . m,n - corresponding number of processors in each dimension 1801 (or PETSC_DECIDE to have calculated) 1802 . dof - number of degrees of freedom per node 1803 . s - stencil width 1804 - lx, ly - arrays containing the number of nodes in each cell along 1805 the x and y coordinates, or PETSC_NULL. If non-null, these 1806 must be of length as m and n, and the corresponding 1807 m and n cannot be PETSC_DECIDE. The sum of the lx[] entries 1808 must be M, and the sum of the ly[] entries must be N. 1809 1810 Output Parameter: 1811 . da - the resulting distributed array object 1812 1813 Options Database Key: 1814 + -da_view - Calls DAView() at the conclusion of DACreate2d() 1815 . -da_grid_x <nx> - number of grid points in x direction, if M < 0 1816 . -da_grid_y <ny> - number of grid points in y direction, if N < 0 1817 . -da_processors_x <nx> - number of processors in x direction 1818 . -da_processors_y <ny> - number of processors in y direction 1819 . -da_refine_x - refinement ratio in x direction 1820 - -da_refine_y - refinement ratio in y direction 1821 1822 Level: beginner 1823 1824 Notes: 1825 The stencil type DA_STENCIL_STAR with width 1 corresponds to the 1826 standard 5-pt stencil, while DA_STENCIL_BOX with width 1 denotes 1827 the standard 9-pt stencil. 1828 1829 The array data itself is NOT stored in the DA, it is stored in Vec objects; 1830 The appropriate vector objects can be obtained with calls to DACreateGlobalVector() 1831 and DACreateLocalVector() and calls to VecDuplicate() if more are needed. 1832 1833 .keywords: distributed array, create, two-dimensional 1834 1835 .seealso: DADestroy(), DAView(), DACreate1d(), DACreate3d(), DAGlobalToLocalBegin(), DAGetRefinementFactor(), 1836 DAGlobalToLocalEnd(), DALocalToGlobal(), DALocalToLocalBegin(), DALocalToLocalEnd(), DASetRefinementFactor(), 1837 DAGetInfo(), DACreateGlobalVector(), DACreateLocalVector(), DACreateNaturalVector(), DALoad(), DAView(), DAGetOwnershipRanges() 1838 1839 @*/ 1840 PetscErrorCode PETSCDM_DLLEXPORT DACreate2d(MPI_Comm comm,DAPeriodicType wrap,DAStencilType stencil_type, 1841 PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DA *da) 1842 { 1843 PetscErrorCode ierr; 1844 1845 PetscFunctionBegin; 1846 ierr = DACreate(comm, da);CHKERRQ(ierr); 1847 ierr = DASetDim(*da, 2);CHKERRQ(ierr); 1848 ierr = DASetSizes(*da, M, N, 1);CHKERRQ(ierr); 1849 ierr = DASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr); 1850 ierr = DASetPeriodicity(*da, wrap);CHKERRQ(ierr); 1851 ierr = DASetDof(*da, dof);CHKERRQ(ierr); 1852 ierr = DASetStencilType(*da, stencil_type);CHKERRQ(ierr); 1853 ierr = DASetStencilWidth(*da, s);CHKERRQ(ierr); 1854 ierr = DASetOwnershipRanges(*da, lx, ly, PETSC_NULL);CHKERRQ(ierr); 1855 /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */ 1856 ierr = DASetFromOptions(*da);CHKERRQ(ierr); 1857 ierr = DASetUp(*da);CHKERRQ(ierr); 1858 PetscFunctionReturn(0); 1859 } 1860