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