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