1 2 /* 3 Code for manipulating distributed regular arrays in parallel. 4 */ 5 6 #include "private/daimpl.h" /*I "petscdm.h" I*/ 7 8 /* 9 This allows the DMDA vectors to properly tell MATLAB their dimensions 10 */ 11 #if defined(PETSC_HAVE_MATLAB_ENGINE) 12 #include "engine.h" /* MATLAB include file */ 13 #include "mex.h" /* MATLAB include file */ 14 EXTERN_C_BEGIN 15 #undef __FUNCT__ 16 #define __FUNCT__ "VecMatlabEnginePut_DA2d" 17 PetscErrorCode VecMatlabEnginePut_DA2d(PetscObject obj,void *mengine) 18 { 19 PetscErrorCode ierr; 20 PetscInt n,m; 21 Vec vec = (Vec)obj; 22 PetscScalar *array; 23 mxArray *mat; 24 DM da; 25 26 PetscFunctionBegin; 27 ierr = PetscObjectQuery((PetscObject)vec,"DMDA",(PetscObject*)&da);CHKERRQ(ierr); 28 if (!da) SETERRQ(((PetscObject)vec)->comm,PETSC_ERR_ARG_WRONGSTATE,"Vector not associated with a DMDA"); 29 ierr = DMDAGetGhostCorners(da,0,0,0,&m,&n,0);CHKERRQ(ierr); 30 31 ierr = VecGetArray(vec,&array);CHKERRQ(ierr); 32 #if !defined(PETSC_USE_COMPLEX) 33 mat = mxCreateDoubleMatrix(m,n,mxREAL); 34 #else 35 mat = mxCreateDoubleMatrix(m,n,mxCOMPLEX); 36 #endif 37 ierr = PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));CHKERRQ(ierr); 38 ierr = PetscObjectName(obj);CHKERRQ(ierr); 39 engPutVariable((Engine *)mengine,obj->name,mat); 40 41 ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr); 42 PetscFunctionReturn(0); 43 } 44 EXTERN_C_END 45 #endif 46 47 48 #undef __FUNCT__ 49 #define __FUNCT__ "DMCreateLocalVector_DA" 50 PetscErrorCode DMCreateLocalVector_DA(DM da,Vec* g) 51 { 52 PetscErrorCode ierr; 53 DM_DA *dd = (DM_DA*)da->data; 54 55 PetscFunctionBegin; 56 PetscValidHeaderSpecific(da,DM_CLASSID,1); 57 PetscValidPointer(g,2); 58 ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr); 59 ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr); 60 ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr); 61 ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr); 62 ierr = PetscObjectCompose((PetscObject)*g,"DMDA",(PetscObject)da);CHKERRQ(ierr); 63 #if defined(PETSC_HAVE_MATLAB_ENGINE) 64 if (dd->w == 1 && dd->dim == 2) { 65 ierr = PetscObjectComposeFunctionDynamic((PetscObject)*g,"PetscMatlabEnginePut_C","VecMatlabEnginePut_DA2d",VecMatlabEnginePut_DA2d);CHKERRQ(ierr); 66 } 67 #endif 68 PetscFunctionReturn(0); 69 } 70 71 #undef __FUNCT__ 72 #define __FUNCT__ "DMGetLocalVector" 73 /*@ 74 DMGetLocalVector - Gets a Seq PETSc vector that 75 may be used with the DMXXX routines. This vector has spaces for the ghost values. 76 77 Not Collective 78 79 Input Parameter: 80 . dm - the distributed array 81 82 Output Parameter: 83 . g - the local vector 84 85 Level: beginner 86 87 Note: 88 The vector values are NOT initialized and may have garbage in them, so you may need 89 to zero them. 90 91 The output parameter, g, is a regular PETSc vector that should be returned with 92 DMRestoreLocalVector() DO NOT call VecDestroy() on it. 93 94 VecStride*() operations can be useful when using DM with dof > 1 95 96 .keywords: distributed array, create, local, vector 97 98 .seealso: DMCreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(), 99 DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMGlobalToLocalBegin(), 100 DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMCreateLocalVector(), DMRestoreLocalVector(), 101 VecStrideMax(), VecStrideMin(), VecStrideNorm() 102 @*/ 103 PetscErrorCode DMGetLocalVector(DM dm,Vec* g) 104 { 105 PetscErrorCode ierr,i; 106 107 PetscFunctionBegin; 108 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 109 PetscValidPointer(g,2); 110 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 111 if (dm->localin[i]) { 112 *g = dm->localin[i]; 113 dm->localin[i] = PETSC_NULL; 114 goto alldone; 115 } 116 } 117 ierr = DMCreateLocalVector(dm,g);CHKERRQ(ierr); 118 119 alldone: 120 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 121 if (!dm->localout[i]) { 122 dm->localout[i] = *g; 123 break; 124 } 125 } 126 PetscFunctionReturn(0); 127 } 128 129 #undef __FUNCT__ 130 #define __FUNCT__ "DMRestoreLocalVector" 131 /*@ 132 DMRestoreLocalVector - Returns a Seq PETSc vector that 133 obtained from DMGetLocalVector(). Do not use with vector obtained via 134 DMCreateLocalVector(). 135 136 Not Collective 137 138 Input Parameter: 139 + dm - the distributed array 140 - g - the local vector 141 142 Level: beginner 143 144 .keywords: distributed array, create, local, vector 145 146 .seealso: DMCreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(), 147 DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMGlobalToLocalBegin(), 148 DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMCreateLocalVector(), DMGetLocalVector() 149 @*/ 150 PetscErrorCode DMRestoreLocalVector(DM dm,Vec* g) 151 { 152 PetscErrorCode ierr; 153 PetscInt i,j; 154 155 PetscFunctionBegin; 156 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 157 PetscValidPointer(g,2); 158 for (j=0; j<DM_MAX_WORK_VECTORS; j++) { 159 if (*g == dm->localout[j]) { 160 dm->localout[j] = PETSC_NULL; 161 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 162 if (!dm->localin[i]) { 163 dm->localin[i] = *g; 164 goto alldone; 165 } 166 } 167 } 168 } 169 ierr = VecDestroy(*g);CHKERRQ(ierr); 170 alldone: 171 PetscFunctionReturn(0); 172 } 173 174 #undef __FUNCT__ 175 #define __FUNCT__ "DMGetGlobalVector" 176 /*@ 177 DMGetGlobalVector - Gets a MPI PETSc vector that 178 may be used with the DMXXX routines. 179 180 Collective on DM 181 182 Input Parameter: 183 . dm - the distributed array 184 185 Output Parameter: 186 . g - the global vector 187 188 Level: beginner 189 190 Note: 191 The vector values are NOT initialized and may have garbage in them, so you may need 192 to zero them. 193 194 The output parameter, g, is a regular PETSc vector that should be returned with 195 DMRestoreGlobalVector() DO NOT call VecDestroy() on it. 196 197 VecStride*() operations can be useful when using DM with dof > 1 198 199 .keywords: distributed array, create, Global, vector 200 201 .seealso: DMCreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(), 202 DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMGlobalToLocalBegin(), 203 DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMCreateLocalVector(), DMRestoreLocalVector() 204 VecStrideMax(), VecStrideMin(), VecStrideNorm() 205 206 @*/ 207 PetscErrorCode DMGetGlobalVector(DM dm,Vec* g) 208 { 209 PetscErrorCode ierr; 210 PetscInt i; 211 212 PetscFunctionBegin; 213 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 214 PetscValidPointer(g,2); 215 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 216 if (dm->globalin[i]) { 217 *g = dm->globalin[i]; 218 dm->globalin[i] = PETSC_NULL; 219 goto alldone; 220 } 221 } 222 ierr = DMCreateGlobalVector(dm,g);CHKERRQ(ierr); 223 224 alldone: 225 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 226 if (!dm->globalout[i]) { 227 dm->globalout[i] = *g; 228 break; 229 } 230 } 231 PetscFunctionReturn(0); 232 } 233 234 #undef __FUNCT__ 235 #define __FUNCT__ "DMRestoreGlobalVector" 236 /*@ 237 DMRestoreGlobalVector - Returns a Seq PETSc vector that 238 obtained from DMGetGlobalVector(). Do not use with vector obtained via 239 DMCreateGlobalVector(). 240 241 Not Collective 242 243 Input Parameter: 244 + dm - the distributed array 245 - g - the global vector 246 247 Level: beginner 248 249 .keywords: distributed array, create, global, vector 250 251 .seealso: DMCreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(), 252 DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMGlobalToGlobalBegin(), 253 DMGlobalToGlobalEnd(), DMGlobalToGlobal(), DMCreateLocalVector(), DMGetGlobalVector() 254 @*/ 255 PetscErrorCode DMRestoreGlobalVector(DM dm,Vec* g) 256 { 257 PetscErrorCode ierr; 258 PetscInt i,j; 259 260 PetscFunctionBegin; 261 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 262 PetscValidPointer(g,2); 263 for (j=0; j<DM_MAX_WORK_VECTORS; j++) { 264 if (*g == dm->globalout[j]) { 265 dm->globalout[j] = PETSC_NULL; 266 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 267 if (!dm->globalin[i]) { 268 dm->globalin[i] = *g; 269 goto alldone; 270 } 271 } 272 } 273 } 274 ierr = VecDestroy(*g);CHKERRQ(ierr); 275 alldone: 276 PetscFunctionReturn(0); 277 } 278 279 /* ------------------------------------------------------------------- */ 280 #if defined(PETSC_HAVE_ADIC) 281 282 EXTERN_C_BEGIN 283 #include "adic/ad_utils.h" 284 EXTERN_C_END 285 286 #undef __FUNCT__ 287 #define __FUNCT__ "DMDAGetAdicArray" 288 /*@C 289 DMDAGetAdicArray - Gets an array of derivative types for a DMDA 290 291 Input Parameter: 292 + da - information about my local patch 293 - ghosted - do you want arrays for the ghosted or nonghosted patch 294 295 Output Parameters: 296 + vptr - array data structured to be passed to ad_FormFunctionLocal() 297 . array_start - actual start of 1d array of all values that adiC can access directly (may be null) 298 - tdof - total number of degrees of freedom represented in array_start (may be null) 299 300 Notes: 301 The vector values are NOT initialized and may have garbage in them, so you may need 302 to zero them. 303 304 Returns the same type of object as the DMDAVecGetArray() except its elements are 305 derivative types instead of PetscScalars 306 307 Level: advanced 308 309 .seealso: DMDARestoreAdicArray() 310 311 @*/ 312 PetscErrorCode DMDAGetAdicArray(DM da,PetscBool ghosted,void *vptr,void *array_start,PetscInt *tdof) 313 { 314 PetscErrorCode ierr; 315 PetscInt j,i,deriv_type_size,xs,ys,xm,ym,zs,zm,itdof; 316 char *iarray_start; 317 void **iptr = (void**)vptr; 318 DM_DA *dd = (DM_DA*)da->data; 319 320 PetscFunctionBegin; 321 PetscValidHeaderSpecific(da,DM_CLASSID,1); 322 if (ghosted) { 323 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 324 if (dd->adarrayghostedin[i]) { 325 *iptr = dd->adarrayghostedin[i]; 326 iarray_start = (char*)dd->adstartghostedin[i]; 327 itdof = dd->ghostedtdof; 328 dd->adarrayghostedin[i] = PETSC_NULL; 329 dd->adstartghostedin[i] = PETSC_NULL; 330 331 goto done; 332 } 333 } 334 xs = dd->Xs; 335 ys = dd->Ys; 336 zs = dd->Zs; 337 xm = dd->Xe-dd->Xs; 338 ym = dd->Ye-dd->Ys; 339 zm = dd->Ze-dd->Zs; 340 } else { 341 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 342 if (dd->adarrayin[i]) { 343 *iptr = dd->adarrayin[i]; 344 iarray_start = (char*)dd->adstartin[i]; 345 itdof = dd->tdof; 346 dd->adarrayin[i] = PETSC_NULL; 347 dd->adstartin[i] = PETSC_NULL; 348 349 goto done; 350 } 351 } 352 xs = dd->xs; 353 ys = dd->ys; 354 zs = dd->zs; 355 xm = dd->xe-dd->xs; 356 ym = dd->ye-dd->ys; 357 zm = dd->ze-dd->zs; 358 } 359 deriv_type_size = PetscADGetDerivTypeSize(); 360 361 switch (dd->dim) { 362 case 1: { 363 void *ptr; 364 itdof = xm; 365 366 ierr = PetscMalloc(xm*deriv_type_size,&iarray_start);CHKERRQ(ierr); 367 368 ptr = (void*)(iarray_start - xs*deriv_type_size); 369 *iptr = (void*)ptr; 370 break;} 371 case 2: { 372 void **ptr; 373 itdof = xm*ym; 374 375 ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*deriv_type_size,&iarray_start);CHKERRQ(ierr); 376 377 ptr = (void**)(iarray_start + xm*ym*deriv_type_size - ys*sizeof(void*)); 378 for(j=ys;j<ys+ym;j++) { 379 ptr[j] = iarray_start + deriv_type_size*(xm*(j-ys) - xs); 380 } 381 *iptr = (void*)ptr; 382 break;} 383 case 3: { 384 void ***ptr,**bptr; 385 itdof = xm*ym*zm; 386 387 ierr = PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*deriv_type_size,&iarray_start);CHKERRQ(ierr); 388 389 ptr = (void***)(iarray_start + xm*ym*zm*deriv_type_size - zs*sizeof(void*)); 390 bptr = (void**)(iarray_start + xm*ym*zm*deriv_type_size + zm*sizeof(void**)); 391 392 for(i=zs;i<zs+zm;i++) { 393 ptr[i] = bptr + ((i-zs)*ym - ys); 394 } 395 for (i=zs; i<zs+zm; i++) { 396 for (j=ys; j<ys+ym; j++) { 397 ptr[i][j] = iarray_start + deriv_type_size*(xm*ym*(i-zs) + xm*(j-ys) - xs); 398 } 399 } 400 401 *iptr = (void*)ptr; 402 break;} 403 default: 404 SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim); 405 } 406 407 done: 408 /* add arrays to the checked out list */ 409 if (ghosted) { 410 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 411 if (!dd->adarrayghostedout[i]) { 412 dd->adarrayghostedout[i] = *iptr ; 413 dd->adstartghostedout[i] = iarray_start; 414 dd->ghostedtdof = itdof; 415 break; 416 } 417 } 418 } else { 419 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 420 if (!dd->adarrayout[i]) { 421 dd->adarrayout[i] = *iptr ; 422 dd->adstartout[i] = iarray_start; 423 dd->tdof = itdof; 424 break; 425 } 426 } 427 } 428 if (i == DMDA_MAX_AD_ARRAYS+1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Too many DMDA ADIC arrays obtained"); 429 if (tdof) *tdof = itdof; 430 if (array_start) *(void**)array_start = iarray_start; 431 PetscFunctionReturn(0); 432 } 433 434 #undef __FUNCT__ 435 #define __FUNCT__ "DMDARestoreAdicArray" 436 /*@C 437 DMDARestoreAdicArray - Restores an array of derivative types for a DMDA 438 439 Input Parameter: 440 + da - information about my local patch 441 - ghosted - do you want arrays for the ghosted or nonghosted patch 442 443 Output Parameters: 444 + ptr - array data structured to be passed to ad_FormFunctionLocal() 445 . array_start - actual start of 1d array of all values that adiC can access directly 446 - tdof - total number of degrees of freedom represented in array_start 447 448 Level: advanced 449 450 .seealso: DMDAGetAdicArray() 451 452 @*/ 453 PetscErrorCode DMDARestoreAdicArray(DM da,PetscBool ghosted,void *ptr,void *array_start,PetscInt *tdof) 454 { 455 PetscInt i; 456 void **iptr = (void**)ptr,iarray_start = 0; 457 458 PetscFunctionBegin; 459 PetscValidHeaderSpecific(da,DM_CLASSID,1); 460 if (ghosted) { 461 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 462 if (dd->adarrayghostedout[i] == *iptr) { 463 iarray_start = dd->adstartghostedout[i]; 464 dd->adarrayghostedout[i] = PETSC_NULL; 465 dd->adstartghostedout[i] = PETSC_NULL; 466 break; 467 } 468 } 469 if (!iarray_start) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Could not find array in checkout list"); 470 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 471 if (!dd->adarrayghostedin[i]){ 472 dd->adarrayghostedin[i] = *iptr; 473 dd->adstartghostedin[i] = iarray_start; 474 break; 475 } 476 } 477 } else { 478 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 479 if (dd->adarrayout[i] == *iptr) { 480 iarray_start = dd->adstartout[i]; 481 dd->adarrayout[i] = PETSC_NULL; 482 dd->adstartout[i] = PETSC_NULL; 483 break; 484 } 485 } 486 if (!iarray_start) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Could not find array in checkout list"); 487 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 488 if (!dd->adarrayin[i]){ 489 dd->adarrayin[i] = *iptr; 490 dd->adstartin[i] = iarray_start; 491 break; 492 } 493 } 494 } 495 PetscFunctionReturn(0); 496 } 497 498 #undef __FUNCT__ 499 #define __FUNCT__ "ad_DAGetArray" 500 PetscErrorCode ad_DAGetArray(DM da,PetscBool ghosted,void *iptr) 501 { 502 PetscErrorCode ierr; 503 PetscFunctionBegin; 504 ierr = DMDAGetAdicArray(da,ghosted,iptr,0,0);CHKERRQ(ierr); 505 PetscFunctionReturn(0); 506 } 507 508 #undef __FUNCT__ 509 #define __FUNCT__ "ad_DARestoreArray" 510 PetscErrorCode ad_DARestoreArray(DM da,PetscBool ghosted,void *iptr) 511 { 512 PetscErrorCode ierr; 513 PetscFunctionBegin; 514 ierr = DMDARestoreAdicArray(da,ghosted,iptr,0,0);CHKERRQ(ierr); 515 PetscFunctionReturn(0); 516 } 517 518 #endif 519 520 #undef __FUNCT__ 521 #define __FUNCT__ "DMDAGetArray" 522 /*@C 523 DMDAGetArray - Gets a work array for a DMDA 524 525 Input Parameter: 526 + da - information about my local patch 527 - ghosted - do you want arrays for the ghosted or nonghosted patch 528 529 Output Parameters: 530 . vptr - array data structured 531 532 Note: The vector values are NOT initialized and may have garbage in them, so you may need 533 to zero them. 534 535 Level: advanced 536 537 .seealso: DMDARestoreArray(), DMDAGetAdicArray() 538 539 @*/ 540 PetscErrorCode DMDAGetArray(DM da,PetscBool ghosted,void *vptr) 541 { 542 PetscErrorCode ierr; 543 PetscInt j,i,xs,ys,xm,ym,zs,zm; 544 char *iarray_start; 545 void **iptr = (void**)vptr; 546 DM_DA *dd = (DM_DA*)da->data; 547 548 PetscFunctionBegin; 549 PetscValidHeaderSpecific(da,DM_CLASSID,1); 550 if (ghosted) { 551 for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 552 if (dd->arrayghostedin[i]) { 553 *iptr = dd->arrayghostedin[i]; 554 iarray_start = (char*)dd->startghostedin[i]; 555 dd->arrayghostedin[i] = PETSC_NULL; 556 dd->startghostedin[i] = PETSC_NULL; 557 558 goto done; 559 } 560 } 561 xs = dd->Xs; 562 ys = dd->Ys; 563 zs = dd->Zs; 564 xm = dd->Xe-dd->Xs; 565 ym = dd->Ye-dd->Ys; 566 zm = dd->Ze-dd->Zs; 567 } else { 568 for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 569 if (dd->arrayin[i]) { 570 *iptr = dd->arrayin[i]; 571 iarray_start = (char*)dd->startin[i]; 572 dd->arrayin[i] = PETSC_NULL; 573 dd->startin[i] = PETSC_NULL; 574 575 goto done; 576 } 577 } 578 xs = dd->xs; 579 ys = dd->ys; 580 zs = dd->zs; 581 xm = dd->xe-dd->xs; 582 ym = dd->ye-dd->ys; 583 zm = dd->ze-dd->zs; 584 } 585 586 switch (dd->dim) { 587 case 1: { 588 void *ptr; 589 590 ierr = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 591 592 ptr = (void*)(iarray_start - xs*sizeof(PetscScalar)); 593 *iptr = (void*)ptr; 594 break;} 595 case 2: { 596 void **ptr; 597 598 ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 599 600 ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*)); 601 for(j=ys;j<ys+ym;j++) { 602 ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs); 603 } 604 *iptr = (void*)ptr; 605 break;} 606 case 3: { 607 void ***ptr,**bptr; 608 609 ierr = PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 610 611 ptr = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*)); 612 bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**)); 613 for(i=zs;i<zs+zm;i++) { 614 ptr[i] = bptr + ((i-zs)*ym - ys); 615 } 616 for (i=zs; i<zs+zm; i++) { 617 for (j=ys; j<ys+ym; j++) { 618 ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs); 619 } 620 } 621 622 *iptr = (void*)ptr; 623 break;} 624 default: 625 SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim); 626 } 627 628 done: 629 /* add arrays to the checked out list */ 630 if (ghosted) { 631 for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 632 if (!dd->arrayghostedout[i]) { 633 dd->arrayghostedout[i] = *iptr ; 634 dd->startghostedout[i] = iarray_start; 635 break; 636 } 637 } 638 } else { 639 for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 640 if (!dd->arrayout[i]) { 641 dd->arrayout[i] = *iptr ; 642 dd->startout[i] = iarray_start; 643 break; 644 } 645 } 646 } 647 PetscFunctionReturn(0); 648 } 649 650 #undef __FUNCT__ 651 #define __FUNCT__ "DMDARestoreArray" 652 /*@C 653 DMDARestoreArray - Restores an array of derivative types for a DMDA 654 655 Input Parameter: 656 + da - information about my local patch 657 . ghosted - do you want arrays for the ghosted or nonghosted patch 658 - vptr - array data structured to be passed to ad_FormFunctionLocal() 659 660 Level: advanced 661 662 .seealso: DMDAGetArray(), DMDAGetAdicArray() 663 664 @*/ 665 PetscErrorCode DMDARestoreArray(DM da,PetscBool ghosted,void *vptr) 666 { 667 PetscInt i; 668 void **iptr = (void**)vptr,*iarray_start = 0; 669 DM_DA *dd = (DM_DA*)da->data; 670 671 PetscFunctionBegin; 672 PetscValidHeaderSpecific(da,DM_CLASSID,1); 673 if (ghosted) { 674 for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 675 if (dd->arrayghostedout[i] == *iptr) { 676 iarray_start = dd->startghostedout[i]; 677 dd->arrayghostedout[i] = PETSC_NULL; 678 dd->startghostedout[i] = PETSC_NULL; 679 break; 680 } 681 } 682 for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 683 if (!dd->arrayghostedin[i]){ 684 dd->arrayghostedin[i] = *iptr; 685 dd->startghostedin[i] = iarray_start; 686 break; 687 } 688 } 689 } else { 690 for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 691 if (dd->arrayout[i] == *iptr) { 692 iarray_start = dd->startout[i]; 693 dd->arrayout[i] = PETSC_NULL; 694 dd->startout[i] = PETSC_NULL; 695 break; 696 } 697 } 698 for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 699 if (!dd->arrayin[i]){ 700 dd->arrayin[i] = *iptr; 701 dd->startin[i] = iarray_start; 702 break; 703 } 704 } 705 } 706 PetscFunctionReturn(0); 707 } 708 709 #undef __FUNCT__ 710 #define __FUNCT__ "DMDAGetAdicMFArray" 711 /*@C 712 DMDAGetAdicMFArray - Gets an array of derivative types for a DMDA for matrix-free ADIC. 713 714 Input Parameter: 715 + da - information about my local patch 716 - ghosted - do you want arrays for the ghosted or nonghosted patch? 717 718 Output Parameters: 719 + vptr - array data structured to be passed to ad_FormFunctionLocal() 720 . array_start - actual start of 1d array of all values that adiC can access directly (may be null) 721 - tdof - total number of degrees of freedom represented in array_start (may be null) 722 723 Notes: 724 The vector values are NOT initialized and may have garbage in them, so you may need 725 to zero them. 726 727 This routine returns the same type of object as the DMDAVecGetArray(), except its 728 elements are derivative types instead of PetscScalars. 729 730 Level: advanced 731 732 .seealso: DMDARestoreAdicMFArray(), DMDAGetArray(), DMDAGetAdicArray() 733 734 @*/ 735 PetscErrorCode DMDAGetAdicMFArray(DM da,PetscBool ghosted,void *vptr,void *array_start,PetscInt *tdof) 736 { 737 PetscErrorCode ierr; 738 PetscInt j,i,xs,ys,xm,ym,zs,zm,itdof = 0; 739 char *iarray_start; 740 void **iptr = (void**)vptr; 741 DM_DA *dd = (DM_DA*)da->data; 742 743 PetscFunctionBegin; 744 PetscValidHeaderSpecific(da,DM_CLASSID,1); 745 if (ghosted) { 746 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 747 if (dd->admfarrayghostedin[i]) { 748 *iptr = dd->admfarrayghostedin[i]; 749 iarray_start = (char*)dd->admfstartghostedin[i]; 750 itdof = dd->ghostedtdof; 751 dd->admfarrayghostedin[i] = PETSC_NULL; 752 dd->admfstartghostedin[i] = PETSC_NULL; 753 754 goto done; 755 } 756 } 757 xs = dd->Xs; 758 ys = dd->Ys; 759 zs = dd->Zs; 760 xm = dd->Xe-dd->Xs; 761 ym = dd->Ye-dd->Ys; 762 zm = dd->Ze-dd->Zs; 763 } else { 764 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 765 if (dd->admfarrayin[i]) { 766 *iptr = dd->admfarrayin[i]; 767 iarray_start = (char*)dd->admfstartin[i]; 768 itdof = dd->tdof; 769 dd->admfarrayin[i] = PETSC_NULL; 770 dd->admfstartin[i] = PETSC_NULL; 771 772 goto done; 773 } 774 } 775 xs = dd->xs; 776 ys = dd->ys; 777 zs = dd->zs; 778 xm = dd->xe-dd->xs; 779 ym = dd->ye-dd->ys; 780 zm = dd->ze-dd->zs; 781 } 782 783 switch (dd->dim) { 784 case 1: { 785 void *ptr; 786 itdof = xm; 787 788 ierr = PetscMalloc(xm*2*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 789 790 ptr = (void*)(iarray_start - xs*2*sizeof(PetscScalar)); 791 *iptr = (void*)ptr; 792 break;} 793 case 2: { 794 void **ptr; 795 itdof = xm*ym; 796 797 ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*2*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 798 799 ptr = (void**)(iarray_start + xm*ym*2*sizeof(PetscScalar) - ys*sizeof(void*)); 800 for(j=ys;j<ys+ym;j++) { 801 ptr[j] = iarray_start + 2*sizeof(PetscScalar)*(xm*(j-ys) - xs); 802 } 803 *iptr = (void*)ptr; 804 break;} 805 case 3: { 806 void ***ptr,**bptr; 807 itdof = xm*ym*zm; 808 809 ierr = PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*2*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 810 811 ptr = (void***)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) - zs*sizeof(void*)); 812 bptr = (void**)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) + zm*sizeof(void**)); 813 for(i=zs;i<zs+zm;i++) { 814 ptr[i] = bptr + ((i-zs)*ym* - ys)*sizeof(void*); 815 } 816 for (i=zs; i<zs+zm; i++) { 817 for (j=ys; j<ys+ym; j++) { 818 ptr[i][j] = iarray_start + 2*sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs); 819 } 820 } 821 822 *iptr = (void*)ptr; 823 break;} 824 default: 825 SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim); 826 } 827 828 done: 829 /* add arrays to the checked out list */ 830 if (ghosted) { 831 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 832 if (!dd->admfarrayghostedout[i]) { 833 dd->admfarrayghostedout[i] = *iptr ; 834 dd->admfstartghostedout[i] = iarray_start; 835 dd->ghostedtdof = itdof; 836 break; 837 } 838 } 839 } else { 840 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 841 if (!dd->admfarrayout[i]) { 842 dd->admfarrayout[i] = *iptr ; 843 dd->admfstartout[i] = iarray_start; 844 dd->tdof = itdof; 845 break; 846 } 847 } 848 } 849 if (i == DMDA_MAX_AD_ARRAYS+1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Too many DMDA ADIC arrays obtained"); 850 if (tdof) *tdof = itdof; 851 if (array_start) *(void**)array_start = iarray_start; 852 PetscFunctionReturn(0); 853 } 854 855 #undef __FUNCT__ 856 #define __FUNCT__ "DMDAGetAdicMFArray4" 857 PetscErrorCode DMDAGetAdicMFArray4(DM da,PetscBool ghosted,void *vptr,void *array_start,PetscInt *tdof) 858 { 859 PetscErrorCode ierr; 860 PetscInt j,i,xs,ys,xm,ym,zs,zm,itdof = 0; 861 char *iarray_start; 862 void **iptr = (void**)vptr; 863 DM_DA *dd = (DM_DA*)da->data; 864 865 PetscFunctionBegin; 866 PetscValidHeaderSpecific(da,DM_CLASSID,1); 867 if (ghosted) { 868 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 869 if (dd->admfarrayghostedin[i]) { 870 *iptr = dd->admfarrayghostedin[i]; 871 iarray_start = (char*)dd->admfstartghostedin[i]; 872 itdof = dd->ghostedtdof; 873 dd->admfarrayghostedin[i] = PETSC_NULL; 874 dd->admfstartghostedin[i] = PETSC_NULL; 875 876 goto done; 877 } 878 } 879 xs = dd->Xs; 880 ys = dd->Ys; 881 zs = dd->Zs; 882 xm = dd->Xe-dd->Xs; 883 ym = dd->Ye-dd->Ys; 884 zm = dd->Ze-dd->Zs; 885 } else { 886 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 887 if (dd->admfarrayin[i]) { 888 *iptr = dd->admfarrayin[i]; 889 iarray_start = (char*)dd->admfstartin[i]; 890 itdof = dd->tdof; 891 dd->admfarrayin[i] = PETSC_NULL; 892 dd->admfstartin[i] = PETSC_NULL; 893 894 goto done; 895 } 896 } 897 xs = dd->xs; 898 ys = dd->ys; 899 zs = dd->zs; 900 xm = dd->xe-dd->xs; 901 ym = dd->ye-dd->ys; 902 zm = dd->ze-dd->zs; 903 } 904 905 switch (dd->dim) { 906 case 2: { 907 void **ptr; 908 itdof = xm*ym; 909 910 ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*5*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 911 912 ptr = (void**)(iarray_start + xm*ym*5*sizeof(PetscScalar) - ys*sizeof(void*)); 913 for(j=ys;j<ys+ym;j++) { 914 ptr[j] = iarray_start + 5*sizeof(PetscScalar)*(xm*(j-ys) - xs); 915 } 916 *iptr = (void*)ptr; 917 break;} 918 default: 919 SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim); 920 } 921 922 done: 923 /* add arrays to the checked out list */ 924 if (ghosted) { 925 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 926 if (!dd->admfarrayghostedout[i]) { 927 dd->admfarrayghostedout[i] = *iptr ; 928 dd->admfstartghostedout[i] = iarray_start; 929 dd->ghostedtdof = itdof; 930 break; 931 } 932 } 933 } else { 934 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 935 if (!dd->admfarrayout[i]) { 936 dd->admfarrayout[i] = *iptr ; 937 dd->admfstartout[i] = iarray_start; 938 dd->tdof = itdof; 939 break; 940 } 941 } 942 } 943 if (i == DMDA_MAX_AD_ARRAYS+1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Too many DMDA ADIC arrays obtained"); 944 if (tdof) *tdof = itdof; 945 if (array_start) *(void**)array_start = iarray_start; 946 PetscFunctionReturn(0); 947 } 948 949 #undef __FUNCT__ 950 #define __FUNCT__ "DMDAGetAdicMFArray9" 951 PetscErrorCode DMDAGetAdicMFArray9(DM da,PetscBool ghosted,void *vptr,void *array_start,PetscInt *tdof) 952 { 953 PetscErrorCode ierr; 954 PetscInt j,i,xs,ys,xm,ym,zs,zm,itdof = 0; 955 char *iarray_start; 956 void **iptr = (void**)vptr; 957 DM_DA *dd = (DM_DA*)da->data; 958 959 PetscFunctionBegin; 960 PetscValidHeaderSpecific(da,DM_CLASSID,1); 961 if (ghosted) { 962 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 963 if (dd->admfarrayghostedin[i]) { 964 *iptr = dd->admfarrayghostedin[i]; 965 iarray_start = (char*)dd->admfstartghostedin[i]; 966 itdof = dd->ghostedtdof; 967 dd->admfarrayghostedin[i] = PETSC_NULL; 968 dd->admfstartghostedin[i] = PETSC_NULL; 969 970 goto done; 971 } 972 } 973 xs = dd->Xs; 974 ys = dd->Ys; 975 zs = dd->Zs; 976 xm = dd->Xe-dd->Xs; 977 ym = dd->Ye-dd->Ys; 978 zm = dd->Ze-dd->Zs; 979 } else { 980 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 981 if (dd->admfarrayin[i]) { 982 *iptr = dd->admfarrayin[i]; 983 iarray_start = (char*)dd->admfstartin[i]; 984 itdof = dd->tdof; 985 dd->admfarrayin[i] = PETSC_NULL; 986 dd->admfstartin[i] = PETSC_NULL; 987 988 goto done; 989 } 990 } 991 xs = dd->xs; 992 ys = dd->ys; 993 zs = dd->zs; 994 xm = dd->xe-dd->xs; 995 ym = dd->ye-dd->ys; 996 zm = dd->ze-dd->zs; 997 } 998 999 switch (dd->dim) { 1000 case 2: { 1001 void **ptr; 1002 itdof = xm*ym; 1003 1004 ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*10*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 1005 1006 ptr = (void**)(iarray_start + xm*ym*10*sizeof(PetscScalar) - ys*sizeof(void*)); 1007 for(j=ys;j<ys+ym;j++) { 1008 ptr[j] = iarray_start + 10*sizeof(PetscScalar)*(xm*(j-ys) - xs); 1009 } 1010 *iptr = (void*)ptr; 1011 break;} 1012 default: 1013 SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim); 1014 } 1015 1016 done: 1017 /* add arrays to the checked out list */ 1018 if (ghosted) { 1019 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 1020 if (!dd->admfarrayghostedout[i]) { 1021 dd->admfarrayghostedout[i] = *iptr ; 1022 dd->admfstartghostedout[i] = iarray_start; 1023 dd->ghostedtdof = itdof; 1024 break; 1025 } 1026 } 1027 } else { 1028 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 1029 if (!dd->admfarrayout[i]) { 1030 dd->admfarrayout[i] = *iptr ; 1031 dd->admfstartout[i] = iarray_start; 1032 dd->tdof = itdof; 1033 break; 1034 } 1035 } 1036 } 1037 if (i == DMDA_MAX_AD_ARRAYS+1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Too many DMDA ADIC arrays obtained"); 1038 if (tdof) *tdof = itdof; 1039 if (array_start) *(void**)array_start = iarray_start; 1040 PetscFunctionReturn(0); 1041 } 1042 1043 #undef __FUNCT__ 1044 #define __FUNCT__ "DMDAGetAdicMFArrayb" 1045 /*@C 1046 DMDAGetAdicMFArrayb - Gets an array of derivative types for a DMDA for matrix-free ADIC. 1047 1048 Input Parameter: 1049 + da - information about my local patch 1050 - ghosted - do you want arrays for the ghosted or nonghosted patch? 1051 1052 Output Parameters: 1053 + vptr - array data structured to be passed to ad_FormFunctionLocal() 1054 . array_start - actual start of 1d array of all values that adiC can access directly (may be null) 1055 - tdof - total number of degrees of freedom represented in array_start (may be null) 1056 1057 Notes: 1058 The vector values are NOT initialized and may have garbage in them, so you may need 1059 to zero them. 1060 1061 This routine returns the same type of object as the DMDAVecGetArray(), except its 1062 elements are derivative types instead of PetscScalars. 1063 1064 Level: advanced 1065 1066 .seealso: DMDARestoreAdicMFArray(), DMDAGetArray(), DMDAGetAdicArray() 1067 1068 @*/ 1069 PetscErrorCode DMDAGetAdicMFArrayb(DM da,PetscBool ghosted,void *vptr,void *array_start,PetscInt *tdof) 1070 { 1071 PetscErrorCode ierr; 1072 PetscInt j,i,xs,ys,xm,ym,zs,zm,itdof = 0; 1073 char *iarray_start; 1074 void **iptr = (void**)vptr; 1075 DM_DA *dd = (DM_DA*)da->data; 1076 PetscInt bs = dd->w,bs1 = bs+1; 1077 1078 PetscFunctionBegin; 1079 PetscValidHeaderSpecific(da,DM_CLASSID,1); 1080 if (ghosted) { 1081 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 1082 if (dd->admfarrayghostedin[i]) { 1083 *iptr = dd->admfarrayghostedin[i]; 1084 iarray_start = (char*)dd->admfstartghostedin[i]; 1085 itdof = dd->ghostedtdof; 1086 dd->admfarrayghostedin[i] = PETSC_NULL; 1087 dd->admfstartghostedin[i] = PETSC_NULL; 1088 1089 goto done; 1090 } 1091 } 1092 xs = dd->Xs; 1093 ys = dd->Ys; 1094 zs = dd->Zs; 1095 xm = dd->Xe-dd->Xs; 1096 ym = dd->Ye-dd->Ys; 1097 zm = dd->Ze-dd->Zs; 1098 } else { 1099 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 1100 if (dd->admfarrayin[i]) { 1101 *iptr = dd->admfarrayin[i]; 1102 iarray_start = (char*)dd->admfstartin[i]; 1103 itdof = dd->tdof; 1104 dd->admfarrayin[i] = PETSC_NULL; 1105 dd->admfstartin[i] = PETSC_NULL; 1106 1107 goto done; 1108 } 1109 } 1110 xs = dd->xs; 1111 ys = dd->ys; 1112 zs = dd->zs; 1113 xm = dd->xe-dd->xs; 1114 ym = dd->ye-dd->ys; 1115 zm = dd->ze-dd->zs; 1116 } 1117 1118 switch (dd->dim) { 1119 case 1: { 1120 void *ptr; 1121 itdof = xm; 1122 1123 ierr = PetscMalloc(xm*bs1*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 1124 1125 ptr = (void*)(iarray_start - xs*bs1*sizeof(PetscScalar)); 1126 *iptr = (void*)ptr; 1127 break;} 1128 case 2: { 1129 void **ptr; 1130 itdof = xm*ym; 1131 1132 ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*bs1*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 1133 1134 ptr = (void**)(iarray_start + xm*ym*bs1*sizeof(PetscScalar) - ys*sizeof(void*)); 1135 for(j=ys;j<ys+ym;j++) { 1136 ptr[j] = iarray_start + bs1*sizeof(PetscScalar)*(xm*(j-ys) - xs); 1137 } 1138 *iptr = (void*)ptr; 1139 break;} 1140 case 3: { 1141 void ***ptr,**bptr; 1142 itdof = xm*ym*zm; 1143 1144 ierr = PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*bs1*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 1145 1146 ptr = (void***)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) - zs*sizeof(void*)); 1147 bptr = (void**)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) + zm*sizeof(void**)); 1148 for(i=zs;i<zs+zm;i++) { 1149 ptr[i] = bptr + ((i-zs)*ym* - ys)*sizeof(void*); 1150 } 1151 for (i=zs; i<zs+zm; i++) { 1152 for (j=ys; j<ys+ym; j++) { 1153 ptr[i][j] = iarray_start + bs1*sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs); 1154 } 1155 } 1156 1157 *iptr = (void*)ptr; 1158 break;} 1159 default: 1160 SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim); 1161 } 1162 1163 done: 1164 /* add arrays to the checked out list */ 1165 if (ghosted) { 1166 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 1167 if (!dd->admfarrayghostedout[i]) { 1168 dd->admfarrayghostedout[i] = *iptr ; 1169 dd->admfstartghostedout[i] = iarray_start; 1170 dd->ghostedtdof = itdof; 1171 break; 1172 } 1173 } 1174 } else { 1175 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 1176 if (!dd->admfarrayout[i]) { 1177 dd->admfarrayout[i] = *iptr ; 1178 dd->admfstartout[i] = iarray_start; 1179 dd->tdof = itdof; 1180 break; 1181 } 1182 } 1183 } 1184 if (i == DMDA_MAX_AD_ARRAYS+1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Too many DMDA ADIC arrays obtained"); 1185 if (tdof) *tdof = itdof; 1186 if (array_start) *(void**)array_start = iarray_start; 1187 PetscFunctionReturn(0); 1188 } 1189 1190 #undef __FUNCT__ 1191 #define __FUNCT__ "DMDARestoreAdicMFArray" 1192 /*@C 1193 DMDARestoreAdicMFArray - Restores an array of derivative types for a DMDA. 1194 1195 Input Parameter: 1196 + da - information about my local patch 1197 - ghosted - do you want arrays for the ghosted or nonghosted patch? 1198 1199 Output Parameters: 1200 + ptr - array data structure to be passed to ad_FormFunctionLocal() 1201 . array_start - actual start of 1d array of all values that adiC can access directly 1202 - tdof - total number of degrees of freedom represented in array_start 1203 1204 Level: advanced 1205 1206 .seealso: DMDAGetAdicArray() 1207 1208 @*/ 1209 PetscErrorCode DMDARestoreAdicMFArray(DM da,PetscBool ghosted,void *vptr,void *array_start,PetscInt *tdof) 1210 { 1211 PetscInt i; 1212 void **iptr = (void**)vptr,*iarray_start = 0; 1213 DM_DA *dd = (DM_DA*)da->data; 1214 1215 PetscFunctionBegin; 1216 PetscValidHeaderSpecific(da,DM_CLASSID,1); 1217 if (ghosted) { 1218 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 1219 if (dd->admfarrayghostedout[i] == *iptr) { 1220 iarray_start = dd->admfstartghostedout[i]; 1221 dd->admfarrayghostedout[i] = PETSC_NULL; 1222 dd->admfstartghostedout[i] = PETSC_NULL; 1223 break; 1224 } 1225 } 1226 if (!iarray_start) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Could not find array in checkout list"); 1227 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 1228 if (!dd->admfarrayghostedin[i]){ 1229 dd->admfarrayghostedin[i] = *iptr; 1230 dd->admfstartghostedin[i] = iarray_start; 1231 break; 1232 } 1233 } 1234 } else { 1235 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 1236 if (dd->admfarrayout[i] == *iptr) { 1237 iarray_start = dd->admfstartout[i]; 1238 dd->admfarrayout[i] = PETSC_NULL; 1239 dd->admfstartout[i] = PETSC_NULL; 1240 break; 1241 } 1242 } 1243 if (!iarray_start) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Could not find array in checkout list"); 1244 for (i=0; i<DMDA_MAX_AD_ARRAYS; i++) { 1245 if (!dd->admfarrayin[i]){ 1246 dd->admfarrayin[i] = *iptr; 1247 dd->admfstartin[i] = iarray_start; 1248 break; 1249 } 1250 } 1251 } 1252 PetscFunctionReturn(0); 1253 } 1254 1255 #undef __FUNCT__ 1256 #define __FUNCT__ "admf_DAGetArray" 1257 PetscErrorCode admf_DAGetArray(DM da,PetscBool ghosted,void *iptr) 1258 { 1259 PetscErrorCode ierr; 1260 PetscFunctionBegin; 1261 ierr = DMDAGetAdicMFArray(da,ghosted,iptr,0,0);CHKERRQ(ierr); 1262 PetscFunctionReturn(0); 1263 } 1264 1265 #undef __FUNCT__ 1266 #define __FUNCT__ "admf_DARestoreArray" 1267 PetscErrorCode admf_DARestoreArray(DM da,PetscBool ghosted,void *iptr) 1268 { 1269 PetscErrorCode ierr; 1270 PetscFunctionBegin; 1271 ierr = DMDARestoreAdicMFArray(da,ghosted,iptr,0,0);CHKERRQ(ierr); 1272 PetscFunctionReturn(0); 1273 } 1274 1275