Actual source code: dalocal.c

  1: /*$Id: dalocal.c,v 1.35 2001/08/07 03:04:39 balay Exp $*/
  2: 
  3: /*
  4:   Code for manipulating distributed regular arrays in parallel.
  5: */

 7:  #include src/dm/da/daimpl.h

  9: /*
 10:    This allows the DA vectors to properly tell Matlab their dimensions
 11: */
 12: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
 13: #include "engine.h"   /* Matlab include file */
 14: #include "mex.h"      /* Matlab include file */
 15: EXTERN_C_BEGIN
 18: int VecMatlabEnginePut_DA2d(PetscObject obj,void *mengine)
 19: {
 20:   int          ierr,n,m;
 21:   Vec          vec = (Vec)obj;
 22:   PetscScalar  *array;
 23:   mxArray      *mat;
 24:   DA           da;

 27:   PetscObjectQuery((PetscObject)vec,"DA",(PetscObject*)&da);
 28:   if (!da) SETERRQ(1,"Vector not associated with a DA");
 29:   DAGetGhostCorners(da,0,0,0,&m,&n,0);

 31:   VecGetArray(vec,&array);
 32: #if !defined(PETSC_USE_COMPLEX)
 33:   mat  = mxCreateDoubleMatrix(m,n,mxREAL);
 34: #else
 35:   mat  = mxCreateDoubleMatrix(m,n,mxCOMPLEX);
 36: #endif
 37:   PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));
 38:   PetscObjectName(obj);
 39:   mxSetName(mat,obj->name);
 40:   engPutArray((Engine *)mengine,mat);
 41: 
 42:   VecRestoreArray(vec,&array);
 43:   return(0);
 44: }
 45: EXTERN_C_END
 46: #endif


 51: /*@C
 52:    DACreateLocalVector - Creates a Seq PETSc vector that
 53:    may be used with the DAXXX routines.

 55:    Not Collective

 57:    Input Parameter:
 58: .  da - the distributed array

 60:    Output Parameter:
 61: .  g - the local vector

 63:    Level: beginner

 65:    Note:
 66:    The output parameter, g, is a regular PETSc vector that should be destroyed
 67:    with a call to VecDestroy() when usage is finished.

 69: .keywords: distributed array, create, local, vector

 71: .seealso: DACreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
 72:           DACreate1d(), DACreate2d(), DACreate3d(), DAGlobalToLocalBegin(),
 73:           DAGlobalToLocalEnd(), DALocalToGlobal(), DAGetLocalVector(), DARestoreLocalVector()
 74: @*/
 75: int DACreateLocalVector(DA da,Vec* g)
 76: {

 81:   VecCreateSeq(PETSC_COMM_SELF,da->nlocal,g);
 82:   VecSetBlockSize(*g,da->w);
 83:   PetscObjectCompose((PetscObject)*g,"DA",(PetscObject)da);
 84: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
 85:   if (da->w == 1  && da->dim == 2) {
 86:     PetscObjectComposeFunctionDynamic((PetscObject)*g,"PetscMatlabEnginePut_C","VecMatlabEnginePut_DA2d",VecMatlabEnginePut_DA2d);
 87:   }
 88: #endif
 89:   return(0);
 90: }

 94: /*@C
 95:    DAGetLocalVector - Gets a Seq PETSc vector that
 96:    may be used with the DAXXX routines.

 98:    Not Collective

100:    Input Parameter:
101: .  da - the distributed array

103:    Output Parameter:
104: .  g - the local vector

106:    Level: beginner

108:    Note:
109:    The vector values are NOT initialized and may have garbage in them, so you may need
110:    to zero them.

112:    The output parameter, g, is a regular PETSc vector that should be returned with 
113:    DARestoreLocalVector() DO NOT call VecDestroy() on it.

115:    VecStride*() operations can be useful when using DA with dof > 1

117: .keywords: distributed array, create, local, vector

119: .seealso: DACreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
120:           DACreate1d(), DACreate2d(), DACreate3d(), DAGlobalToLocalBegin(),
121:           DAGlobalToLocalEnd(), DALocalToGlobal(), DACreateLocalVector(), DARestoreLocalVector(),
122:           VecStrideMax(), VecStrideMin(), VecStrideNorm()
123: @*/
124: int DAGetLocalVector(DA da,Vec* g)
125: {
126:   int ierr,i;

130:   for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
131:     if (da->localin[i]) {
132:       *g             = da->localin[i];
133:       da->localin[i] = PETSC_NULL;
134:       goto alldone;
135:     }
136:   }
137:   DACreateLocalVector(da,g);

139:   alldone:
140:   for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
141:     if (!da->localout[i]) {
142:       da->localout[i] = *g;
143:       break;
144:     }
145:   }
146:   return(0);
147: }

151: /*@C
152:    DARestoreLocalVector - Returns a Seq PETSc vector that
153:      obtained from DAGetLocalVector(). Do not use with vector obtained via
154:      DACreateLocalVector().

156:    Not Collective

158:    Input Parameter:
159: +  da - the distributed array
160: -  g - the local vector

162:    Level: beginner

164: .keywords: distributed array, create, local, vector

166: .seealso: DACreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
167:           DACreate1d(), DACreate2d(), DACreate3d(), DAGlobalToLocalBegin(),
168:           DAGlobalToLocalEnd(), DALocalToGlobal(), DACreateLocalVector(), DAGetLocalVector()
169: @*/
170: int DARestoreLocalVector(DA da,Vec* g)
171: {
172:   int ierr,i,j;

176:   for (j=0; j<DA_MAX_WORK_VECTORS; j++) {
177:     if (*g == da->localout[j]) {
178:       da->localout[j] = PETSC_NULL;
179:       for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
180:         if (!da->localin[i]) {
181:           da->localin[i] = *g;
182:           goto alldone;
183:         }
184:       }
185:     }
186:   }
187:   VecDestroy(*g);
188:   alldone:
189:   return(0);
190: }

194: /*@C
195:    DAGetGlobalVector - Gets a MPI PETSc vector that
196:    may be used with the DAXXX routines.

198:    Collective on DA

200:    Input Parameter:
201: .  da - the distributed array

203:    Output Parameter:
204: .  g - the global vector

206:    Level: beginner

208:    Note:
209:    The vector values are NOT initialized and may have garbage in them, so you may need
210:    to zero them.

212:    The output parameter, g, is a regular PETSc vector that should be returned with 
213:    DARestoreGlobalVector() DO NOT call VecDestroy() on it.

215:    VecStride*() operations can be useful when using DA with dof > 1

217: .keywords: distributed array, create, Global, vector

219: .seealso: DACreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
220:           DACreate1d(), DACreate2d(), DACreate3d(), DAGlobalToLocalBegin(),
221:           DAGlobalToLocalEnd(), DALocalToGlobal(), DACreateLocalVector(), DARestoreLocalVector()
222:           VecStrideMax(), VecStrideMin(), VecStrideNorm()

224: @*/
225: int DAGetGlobalVector(DA da,Vec* g)
226: {
227:   int ierr,i;

231:   for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
232:     if (da->globalin[i]) {
233:       *g             = da->globalin[i];
234:       da->globalin[i] = PETSC_NULL;
235:       goto alldone;
236:     }
237:   }
238:   DACreateGlobalVector(da,g);

240:   alldone:
241:   for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
242:     if (!da->globalout[i]) {
243:       da->globalout[i] = *g;
244:       break;
245:     }
246:   }
247:   return(0);
248: }

252: /*@C
253:    DARestoreGlobalVector - Returns a Seq PETSc vector that
254:      obtained from DAGetGlobalVector(). Do not use with vector obtained via
255:      DACreateGlobalVector().

257:    Not Collective

259:    Input Parameter:
260: +  da - the distributed array
261: -  g - the global vector

263:    Level: beginner

265: .keywords: distributed array, create, global, vector

267: .seealso: DACreateGlobalVector(), VecDuplicate(), VecDuplicateVecs(),
268:           DACreate1d(), DACreate2d(), DACreate3d(), DAGlobalToGlobalBegin(),
269:           DAGlobalToGlobalEnd(), DAGlobalToGlobal(), DACreateLocalVector(), DAGetGlobalVector()
270: @*/
271: int DARestoreGlobalVector(DA da,Vec* g)
272: {
273:   int ierr,i,j;

277:   for (j=0; j<DA_MAX_WORK_VECTORS; j++) {
278:     if (*g == da->globalout[j]) {
279:       da->globalout[j] = PETSC_NULL;
280:       for (i=0; i<DA_MAX_WORK_VECTORS; i++) {
281:         if (!da->globalin[i]) {
282:           da->globalin[i] = *g;
283:           goto alldone;
284:         }
285:       }
286:     }
287:   }
288:   VecDestroy(*g);
289:   alldone:
290:   return(0);
291: }

293: /* ------------------------------------------------------------------- */
294: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX)

296: EXTERN_C_BEGIN
297: #include "adic/ad_utils.h"
298: EXTERN_C_END

302: /*@C
303:      DAGetAdicArray - Gets an array of derivative types for a DA
304:           
305:     Input Parameter:
306: +    da - information about my local patch
307: -    ghosted - do you want arrays for the ghosted or nonghosted patch

309:     Output Parameters:
310: +    ptr - array data structured to be passed to ad_FormFunctionLocal()
311: .    array_start - actual start of 1d array of all values that adiC can access directly (may be null)
312: -    tdof - total number of degrees of freedom represented in array_start (may be null)

314:      Notes:
315:        The vector values are NOT initialized and may have garbage in them, so you may need
316:        to zero them.

318:        Returns the same type of object as the DAVecGetArray() except its elements are 
319:            derivative types instead of PetscScalars

321:      Level: advanced

323: .seealso: DARestoreAdicArray()

325: @*/
326: int DAGetAdicArray(DA da,PetscTruth ghosted,void **iptr,void **array_start,int *tdof)
327: {
328:   int  ierr,j,i,deriv_type_size,xs,ys,xm,ym,zs,zm,itdof;
329:   char *iarray_start;

333:   if (ghosted) {
334:     for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
335:       if (da->adarrayghostedin[i]) {
336:         *iptr                   = da->adarrayghostedin[i];
337:         iarray_start            = (char*)da->adstartghostedin[i];
338:         itdof                   = da->ghostedtdof;
339:         da->adarrayghostedin[i] = PETSC_NULL;
340:         da->adstartghostedin[i] = PETSC_NULL;
341: 
342:         goto done;
343:       }
344:     }
345:     xs = da->Xs;
346:     ys = da->Ys;
347:     zs = da->Zs;
348:     xm = da->Xe-da->Xs;
349:     ym = da->Ye-da->Ys;
350:     zm = da->Ze-da->Zs;
351:   } else {
352:     for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
353:       if (da->adarrayin[i]) {
354:         *iptr            = da->adarrayin[i];
355:         iarray_start     = (char*)da->adstartin[i];
356:         itdof            = da->tdof;
357:         da->adarrayin[i] = PETSC_NULL;
358:         da->adstartin[i] = PETSC_NULL;
359: 
360:         goto done;
361:       }
362:     }
363:     xs = da->xs;
364:     ys = da->ys;
365:     zs = da->zs;
366:     xm = da->xe-da->xs;
367:     ym = da->ye-da->ys;
368:     zm = da->ze-da->zs;
369:   }
370:   deriv_type_size = PetscADGetDerivTypeSize();

372:   switch (da->dim) {
373:     case 1: {
374:       void *ptr;
375:       itdof = xm;

377:       PetscMalloc(xm*deriv_type_size,&iarray_start);

379:       ptr   = (void*)(iarray_start - xs*deriv_type_size);
380:       *iptr = (void*)ptr;
381:       break;}
382:     case 2: {
383:       void **ptr;
384:       itdof = xm*ym;

386:       PetscMalloc((ym+1)*sizeof(void *)+xm*ym*deriv_type_size,&iarray_start);

388:       ptr  = (void**)(iarray_start + xm*ym*deriv_type_size - ys*sizeof(void*));
389:       for(j=ys;j<ys+ym;j++) {
390:         ptr[j] = iarray_start + deriv_type_size*(xm*(j-ys) - xs);
391:       }
392:       *iptr = (void*)ptr;
393:       break;}
394:     case 3: {
395:       void ***ptr,**bptr;
396:       itdof = xm*ym*zm;

398:       PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*deriv_type_size,&iarray_start);

400:       ptr  = (void***)(iarray_start + xm*ym*zm*deriv_type_size - zs*sizeof(void*));
401:       bptr = (void**)(iarray_start + xm*ym*zm*deriv_type_size + zm*sizeof(void**));

403:       for(i=zs;i<zs+zm;i++) {
404:         ptr[i] = bptr + ((i-zs)*ym - ys);
405:       }
406:       for (i=zs; i<zs+zm; i++) {
407:         for (j=ys; j<ys+ym; j++) {
408:           ptr[i][j] = iarray_start + deriv_type_size*(xm*ym*(i-zs) + xm*(j-ys) - xs);
409:         }
410:       }

412:       *iptr = (void*)ptr;
413:       break;}
414:     default:
415:       SETERRQ1(1,"Dimension %d not supported",da->dim);
416:   }

418:   done:
419:   /* add arrays to the checked out list */
420:   if (ghosted) {
421:     for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
422:       if (!da->adarrayghostedout[i]) {
423:         da->adarrayghostedout[i] = *iptr ;
424:         da->adstartghostedout[i] = iarray_start;
425:         da->ghostedtdof          = itdof;
426:         break;
427:       }
428:     }
429:   } else {
430:     for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
431:       if (!da->adarrayout[i]) {
432:         da->adarrayout[i] = *iptr ;
433:         da->adstartout[i] = iarray_start;
434:         da->tdof          = itdof;
435:         break;
436:       }
437:     }
438:   }
439:   if (i == DA_MAX_AD_ARRAYS+1) SETERRQ(1,"Too many DA ADIC arrays obtained");
440:   if (tdof)        *tdof = itdof;
441:   if (array_start) *array_start = iarray_start;
442:   return(0);
443: }

447: /*@C
448:      DARestoreAdicArray - Restores an array of derivative types for a DA
449:           
450:     Input Parameter:
451: +    da - information about my local patch
452: -    ghosted - do you want arrays for the ghosted or nonghosted patch

454:     Output Parameters:
455: +    ptr - array data structured to be passed to ad_FormFunctionLocal()
456: .    array_start - actual start of 1d array of all values that adiC can access directly
457: -    tdof - total number of degrees of freedom represented in array_start

459:      Level: advanced

461: .seealso: DAGetAdicArray()

463: @*/
464: int DARestoreAdicArray(DA da,PetscTruth ghosted,void **iptr,void **array_start,int *tdof)
465: {
466:   int  i;
467:   void *iarray_start = 0;
468: 
471:   if (ghosted) {
472:     for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
473:       if (da->adarrayghostedout[i] == *iptr) {
474:         iarray_start             = da->adstartghostedout[i];
475:         da->adarrayghostedout[i] = PETSC_NULL;
476:         da->adstartghostedout[i] = PETSC_NULL;
477:         break;
478:       }
479:     }
480:     if (!iarray_start) SETERRQ(1,"Could not find array in checkout list");
481:     for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
482:       if (!da->adarrayghostedin[i]){
483:         da->adarrayghostedin[i] = *iptr;
484:         da->adstartghostedin[i] = iarray_start;
485:         break;
486:       }
487:     }
488:   } else {
489:     for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
490:       if (da->adarrayout[i] == *iptr) {
491:         iarray_start      = da->adstartout[i];
492:         da->adarrayout[i] = PETSC_NULL;
493:         da->adstartout[i] = PETSC_NULL;
494:         break;
495:       }
496:     }
497:     if (!iarray_start) SETERRQ(1,"Could not find array in checkout list");
498:     for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
499:       if (!da->adarrayin[i]){
500:         da->adarrayin[i]   = *iptr;
501:         da->adstartin[i]   = iarray_start;
502:         break;
503:       }
504:     }
505:   }
506:   return(0);
507: }

511: int ad_DAGetArray(DA da,PetscTruth ghosted,void **iptr)
512: {
515:   DAGetAdicArray(da,ghosted,iptr,0,0);
516:   return(0);
517: }

521: int ad_DARestoreArray(DA da,PetscTruth ghosted,void **iptr)
522: {
525:   DARestoreAdicArray(da,ghosted,iptr,0,0);
526:   return(0);
527: }

529: #endif

533: /*@C
534:      DAGetArray - Gets a work array for a DA
535:           
536:     Input Parameter:
537: +    da - information about my local patch
538: -    ghosted - do you want arrays for the ghosted or nonghosted patch

540:     Output Parameters:
541: .    ptr - array data structured

543:     Note:  The vector values are NOT initialized and may have garbage in them, so you may need
544:            to zero them.

546:   Level: advanced

548: .seealso: DARestoreArray(), DAGetAdicArray()

550: @*/
551: int DAGetArray(DA da,PetscTruth ghosted,void **iptr)
552: {
553:   int  ierr,j,i,xs,ys,xm,ym,zs,zm;
554:   char *iarray_start;

558:   if (ghosted) {
559:     for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
560:       if (da->arrayghostedin[i]) {
561:         *iptr                 = da->arrayghostedin[i];
562:         iarray_start          = (char*)da->startghostedin[i];
563:         da->arrayghostedin[i] = PETSC_NULL;
564:         da->startghostedin[i] = PETSC_NULL;
565: 
566:         goto done;
567:       }
568:     }
569:     xs = da->Xs;
570:     ys = da->Ys;
571:     zs = da->Zs;
572:     xm = da->Xe-da->Xs;
573:     ym = da->Ye-da->Ys;
574:     zm = da->Ze-da->Zs;
575:   } else {
576:     for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
577:       if (da->arrayin[i]) {
578:         *iptr          = da->arrayin[i];
579:         iarray_start   = (char*)da->startin[i];
580:         da->arrayin[i] = PETSC_NULL;
581:         da->startin[i] = PETSC_NULL;
582: 
583:         goto done;
584:       }
585:     }
586:     xs = da->xs;
587:     ys = da->ys;
588:     zs = da->zs;
589:     xm = da->xe-da->xs;
590:     ym = da->ye-da->ys;
591:     zm = da->ze-da->zs;
592:   }

594:   switch (da->dim) {
595:     case 1: {
596:       void *ptr;

598:       PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);

600:       ptr   = (void*)(iarray_start - xs*sizeof(PetscScalar));
601:       *iptr = (void*)ptr;
602:       break;}
603:     case 2: {
604:       void **ptr;

606:       PetscMalloc((ym+1)*sizeof(void *)+xm*ym*sizeof(PetscScalar),&iarray_start);

608:       ptr  = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*));
609:       for(j=ys;j<ys+ym;j++) {
610:         ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs);
611:       }
612:       *iptr = (void*)ptr;
613:       break;}
614:     case 3: {
615:       void ***ptr,**bptr;

617:       PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);

619:       ptr  = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*));
620:       bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**));
621:       for(i=zs;i<zs+zm;i++) {
622:         ptr[i] = bptr + ((i-zs)*ym - ys);
623:       }
624:       for (i=zs; i<zs+zm; i++) {
625:         for (j=ys; j<ys+ym; j++) {
626:           ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
627:         }
628:       }

630:       *iptr = (void*)ptr;
631:       break;}
632:     default:
633:       SETERRQ1(1,"Dimension %d not supported",da->dim);
634:   }

636:   done:
637:   /* add arrays to the checked out list */
638:   if (ghosted) {
639:     for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
640:       if (!da->arrayghostedout[i]) {
641:         da->arrayghostedout[i] = *iptr ;
642:         da->startghostedout[i] = iarray_start;
643:         break;
644:       }
645:     }
646:   } else {
647:     for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
648:       if (!da->arrayout[i]) {
649:         da->arrayout[i] = *iptr ;
650:         da->startout[i] = iarray_start;
651:         break;
652:       }
653:     }
654:   }
655:   return(0);
656: }

660: /*@C
661:      DARestoreArray - Restores an array of derivative types for a DA
662:           
663:     Input Parameter:
664: +    da - information about my local patch
665: .    ghosted - do you want arrays for the ghosted or nonghosted patch
666: -    ptr - array data structured to be passed to ad_FormFunctionLocal()

668:      Level: advanced

670: .seealso: DAGetArray(), DAGetAdicArray()

672: @*/
673: int DARestoreArray(DA da,PetscTruth ghosted,void **iptr)
674: {
675:   int  i;
676:   void *iarray_start = 0;
677: 
680:   if (ghosted) {
681:     for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
682:       if (da->arrayghostedout[i] == *iptr) {
683:         iarray_start           = da->startghostedout[i];
684:         da->arrayghostedout[i] = PETSC_NULL;
685:         da->startghostedout[i] = PETSC_NULL;
686:         break;
687:       }
688:     }
689:     for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
690:       if (!da->arrayghostedin[i]){
691:         da->arrayghostedin[i] = *iptr;
692:         da->startghostedin[i] = iarray_start;
693:         break;
694:       }
695:     }
696:   } else {
697:     for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
698:       if (da->arrayout[i] == *iptr) {
699:         iarray_start    = da->startout[i];
700:         da->arrayout[i] = PETSC_NULL;
701:         da->startout[i] = PETSC_NULL;
702:         break;
703:       }
704:     }
705:     for (i=0; i<DA_MAX_WORK_ARRAYS; i++) {
706:       if (!da->arrayin[i]){
707:         da->arrayin[i]  = *iptr;
708:         da->startin[i]  = iarray_start;
709:         break;
710:       }
711:     }
712:   }
713:   return(0);
714: }

718: /*@C
719:      DAGetAdicMFArray - Gets an array of derivative types for a DA for matrix-free ADIC.
720:           
721:      Input Parameter:
722: +    da - information about my local patch
723: -    ghosted - do you want arrays for the ghosted or nonghosted patch?

725:      Output Parameters:
726: +    iptr - array data structured to be passed to ad_FormFunctionLocal()
727: .    array_start - actual start of 1d array of all values that adiC can access directly (may be null)
728: -    tdof - total number of degrees of freedom represented in array_start (may be null)

730:      Notes: 
731:      The vector values are NOT initialized and may have garbage in them, so you may need
732:      to zero them.

734:      This routine returns the same type of object as the DAVecGetArray(), except its
735:      elements are derivative types instead of PetscScalars.

737:      Level: advanced

739: .seealso: DARestoreAdicMFArray(), DAGetArray(), DAGetAdicArray()

741: @*/
742: int DAGetAdicMFArray(DA da,PetscTruth ghosted,void **iptr,void **array_start,int *tdof)
743: {
744:   int  ierr,j,i,xs,ys,xm,ym,zs,zm,itdof;
745:   char *iarray_start;

749:   if (ghosted) {
750:     for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
751:       if (da->admfarrayghostedin[i]) {
752:         *iptr                     = da->admfarrayghostedin[i];
753:         iarray_start              = (char*)da->admfstartghostedin[i];
754:         itdof                     = da->ghostedtdof;
755:         da->admfarrayghostedin[i] = PETSC_NULL;
756:         da->admfstartghostedin[i] = PETSC_NULL;
757: 
758:         goto done;
759:       }
760:     }
761:     xs = da->Xs;
762:     ys = da->Ys;
763:     zs = da->Zs;
764:     xm = da->Xe-da->Xs;
765:     ym = da->Ye-da->Ys;
766:     zm = da->Ze-da->Zs;
767:   } else {
768:     for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
769:       if (da->admfarrayin[i]) {
770:         *iptr              = da->admfarrayin[i];
771:         iarray_start       = (char*)da->admfstartin[i];
772:         itdof              = da->tdof;
773:         da->admfarrayin[i] = PETSC_NULL;
774:         da->admfstartin[i] = PETSC_NULL;
775: 
776:         goto done;
777:       }
778:     }
779:     xs = da->xs;
780:     ys = da->ys;
781:     zs = da->zs;
782:     xm = da->xe-da->xs;
783:     ym = da->ye-da->ys;
784:     zm = da->ze-da->zs;
785:   }

787:   switch (da->dim) {
788:     case 1: {
789:       void *ptr;
790:       itdof = xm;

792:       PetscMalloc(xm*2*sizeof(PetscScalar),&iarray_start);

794:       ptr   = (void*)(iarray_start - xs*2*sizeof(PetscScalar));
795:       *iptr = (void*)ptr;
796:       break;}
797:     case 2: {
798:       void **ptr;
799:       itdof = xm*ym;

801:       PetscMalloc((ym+1)*sizeof(void *)+xm*ym*2*sizeof(PetscScalar),&iarray_start);

803:       ptr  = (void**)(iarray_start + xm*ym*2*sizeof(PetscScalar) - ys*sizeof(void*));
804:       for(j=ys;j<ys+ym;j++) {
805:         ptr[j] = iarray_start + 2*sizeof(PetscScalar)*(xm*(j-ys) - xs);
806:       }
807:       *iptr = (void*)ptr;
808:       break;}
809:     case 3: {
810:       void ***ptr,**bptr;
811:       itdof = xm*ym*zm;

813:       PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*2*sizeof(PetscScalar),&iarray_start);

815:       ptr  = (void***)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) - zs*sizeof(void*));
816:       bptr = (void**)(iarray_start + xm*ym*zm*2*sizeof(PetscScalar) + zm*sizeof(void**));
817:       for(i=zs;i<zs+zm;i++) {
818:         ptr[i] = bptr + ((i-zs)*ym* - ys)*sizeof(void*);
819:       }
820:       for (i=zs; i<zs+zm; i++) {
821:         for (j=ys; j<ys+ym; j++) {
822:           ptr[i][j] = iarray_start + 2*sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
823:         }
824:       }

826:       *iptr = (void*)ptr;
827:       break;}
828:     default:
829:       SETERRQ1(1,"Dimension %d not supported",da->dim);
830:   }

832:   done:
833:   /* add arrays to the checked out list */
834:   if (ghosted) {
835:     for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
836:       if (!da->admfarrayghostedout[i]) {
837:         da->admfarrayghostedout[i] = *iptr ;
838:         da->admfstartghostedout[i] = iarray_start;
839:         da->ghostedtdof            = itdof;
840:         break;
841:       }
842:     }
843:   } else {
844:     for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
845:       if (!da->admfarrayout[i]) {
846:         da->admfarrayout[i] = *iptr ;
847:         da->admfstartout[i] = iarray_start;
848:         da->tdof            = itdof;
849:         break;
850:       }
851:     }
852:   }
853:   if (i == DA_MAX_AD_ARRAYS+1) SETERRQ(1,"Too many DA ADIC arrays obtained");
854:   if (tdof)        *tdof = itdof;
855:   if (array_start) *array_start = iarray_start;
856:   return(0);
857: }

861: /*@C
862:      DARestoreAdicMFArray - Restores an array of derivative types for a DA.
863:           
864:      Input Parameter:
865: +    da - information about my local patch
866: -    ghosted - do you want arrays for the ghosted or nonghosted patch?

868:      Output Parameters:
869: +    ptr - array data structure to be passed to ad_FormFunctionLocal()
870: .    array_start - actual start of 1d array of all values that adiC can access directly
871: -    tdof - total number of degrees of freedom represented in array_start

873:      Level: advanced

875: .seealso: DAGetAdicArray()

877: @*/
878: int DARestoreAdicMFArray(DA da,PetscTruth ghosted,void **iptr,void **array_start,int *tdof)
879: {
880:   int  i;
881:   void *iarray_start = 0;
882: 
885:   if (ghosted) {
886:     for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
887:       if (da->admfarrayghostedout[i] == *iptr) {
888:         iarray_start               = da->admfstartghostedout[i];
889:         da->admfarrayghostedout[i] = PETSC_NULL;
890:         da->admfstartghostedout[i] = PETSC_NULL;
891:         break;
892:       }
893:     }
894:     if (!iarray_start) SETERRQ(1,"Could not find array in checkout list");
895:     for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
896:       if (!da->admfarrayghostedin[i]){
897:         da->admfarrayghostedin[i] = *iptr;
898:         da->admfstartghostedin[i] = iarray_start;
899:         break;
900:       }
901:     }
902:   } else {
903:     for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
904:       if (da->admfarrayout[i] == *iptr) {
905:         iarray_start        = da->admfstartout[i];
906:         da->admfarrayout[i] = PETSC_NULL;
907:         da->admfstartout[i] = PETSC_NULL;
908:         break;
909:       }
910:     }
911:     if (!iarray_start) SETERRQ(1,"Could not find array in checkout list");
912:     for (i=0; i<DA_MAX_AD_ARRAYS; i++) {
913:       if (!da->admfarrayin[i]){
914:         da->admfarrayin[i] = *iptr;
915:         da->admfstartin[i] = iarray_start;
916:         break;
917:       }
918:     }
919:   }
920:   return(0);
921: }

925: int admf_DAGetArray(DA da,PetscTruth ghosted,void **iptr)
926: {
929:   DAGetAdicMFArray(da,ghosted,iptr,0,0);
930:   return(0);
931: }

935: int admf_DARestoreArray(DA da,PetscTruth ghosted,void **iptr)
936: {
939:   DARestoreAdicMFArray(da,ghosted,iptr,0,0);
940:   return(0);
941: }