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