Actual source code: da2.c
1: /*$Id: da2.c,v 1.180 2001/09/07 20:12:17 bsmith Exp $*/
2:
3: #include src/dm/da/daimpl.h
5: int DAGetOwnershipRange(DA da,int **lx,int **ly,int **lz)
6: {
9: if (lx) *lx = da->lx;
10: if (ly) *ly = da->ly;
11: if (lz) *lz = da->lz;
12: return(0);
13: }
15: int DAView_2d(DA da,PetscViewer viewer)
16: {
17: int rank,ierr;
18: PetscTruth isascii,isdraw,isbinary;
21: MPI_Comm_rank(da->comm,&rank);
23: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
24: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
25: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
26: if (isascii) {
27: PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %d N %d m %d n %d w %d s %dn",rank,da->M,
28: da->N,da->m,da->n,da->w,da->s);
29: PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %d %d, Y range of indices: %d %dn",da->xs,da->xe,da->ys,da->ye);
30: PetscViewerFlush(viewer);
31: } else if (isdraw) {
32: PetscDraw draw;
33: double ymin = -1*da->s-1,ymax = da->N+da->s;
34: double xmin = -1*da->s-1,xmax = da->M+da->s;
35: double x,y;
36: int base,*idx;
37: char node[10];
38: PetscTruth isnull;
39:
40: PetscViewerDrawGetDraw(viewer,0,&draw);
41: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
42: PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);
43: PetscDrawSynchronizedClear(draw);
45: /* first processor draw all node lines */
46: if (!rank) {
47: ymin = 0.0; ymax = da->N - 1;
48: for (xmin=0; xmin<da->M; xmin++) {
49: PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);
50: }
51: xmin = 0.0; xmax = da->M - 1;
52: for (ymin=0; ymin<da->N; ymin++) {
53: PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);
54: }
55: }
56: PetscDrawSynchronizedFlush(draw);
57: PetscDrawPause(draw);
59: /* draw my box */
60: ymin = da->ys; ymax = da->ye - 1; xmin = da->xs/da->w;
61: xmax =(da->xe-1)/da->w;
62: PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);
63: PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);
64: PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);
65: PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);
67: /* put in numbers */
68: base = (da->base)/da->w;
69: for (y=ymin; y<=ymax; y++) {
70: for (x=xmin; x<=xmax; x++) {
71: sprintf(node,"%d",base++);
72: PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);
73: }
74: }
76: PetscDrawSynchronizedFlush(draw);
77: PetscDrawPause(draw);
78: /* overlay ghost numbers, useful for error checking */
79: /* put in numbers */
81: base = 0; idx = da->idx;
82: ymin = da->Ys; ymax = da->Ye; xmin = da->Xs; xmax = da->Xe;
83: for (y=ymin; y<ymax; y++) {
84: for (x=xmin; x<xmax; x++) {
85: if ((base % da->w) == 0) {
86: sprintf(node,"%d",idx[base]/da->w);
87: PetscDrawString(draw,x/da->w,y,PETSC_DRAW_BLUE,node);
88: }
89: base++;
90: }
91: }
92: PetscDrawSynchronizedFlush(draw);
93: PetscDrawPause(draw);
94: } else if (isbinary) {
95: DAView_Binary(da,viewer);
96: } else {
97: SETERRQ1(1,"Viewer type %s not supported for DA2d",((PetscObject)viewer)->type_name);
98: }
99: return(0);
100: }
102: #if defined(PETSC_HAVE_AMS)
103: /*
104: This function tells the AMS the layout of the vectors, it is called
105: in the VecPublish_xx routines.
106: */
107: EXTERN_C_BEGIN
108: int AMSSetFieldBlock_DA(AMS_Memory amem,char *name,Vec vec)
109: {
110: int ierr,dof,dim,ends[4],shift = 0,starts[] = {0,0,0,0};
111: DA da = 0;
112: PetscTruth isseq,ismpi;
115: if (((PetscObject)vec)->amem < 0) return(0); /* return if not published */
117: PetscObjectQuery((PetscObject)vec,"DA",(PetscObject*)&da);
118: if (!da) return(0);
119: DAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0);
120: if (dof > 1) {dim++; shift = 1; ends[0] = dof;}
122: PetscTypeCompare((PetscObject)vec,VECSEQ,&isseq);
123: PetscTypeCompare((PetscObject)vec,VECMPI,&ismpi);
124: if (isseq) {
125: DAGetGhostCorners(da,0,0,0,ends+shift,ends+shift+1,ends+shift+2);
126: ends[shift] += starts[shift]-1;
127: ends[shift+1] += starts[shift+1]-1;
128: ends[shift+2] += starts[shift+2]-1;
129: AMS_Memory_set_field_block(amem,name,dim,starts,ends);
130: if (ierr) {
131: char *message;
132: AMS_Explain_error(ierr,&message);
133: SETERRQ(ierr,message);
134: }
135: } else if (ismpi) {
136: DAGetCorners(da,starts+shift,starts+shift+1,starts+shift+2,
137: ends+shift,ends+shift+1,ends+shift+2);
138: ends[shift] += starts[shift]-1;
139: ends[shift+1] += starts[shift+1]-1;
140: ends[shift+2] += starts[shift+2]-1;
141: AMS_Memory_set_field_block(amem,name,dim,starts,ends);
142: if (ierr) {
143: char *message;
144: AMS_Explain_error(ierr,&message);
145: SETERRQ(ierr,message);
146: }
147: } else {
148: SETERRQ1(1,"Wrong vector type %s for this call",((PetscObject)vec)->type_name);
149: }
151: return(0);
152: }
153: EXTERN_C_END
154: #endif
156: int DAPublish_Petsc(PetscObject obj)
157: {
158: #if defined(PETSC_HAVE_AMS)
159: DA v = (DA) obj;
160: int ierr;
161: #endif
165: #if defined(PETSC_HAVE_AMS)
166: /* if it is already published then return */
167: if (v->amem >=0) return(0);
169: PetscObjectPublishBaseBegin(obj);
170: PetscObjectPublishBaseEnd(obj);
171: #endif
173: return(0);
174: }
176: /*
177: This allows the DA vectors to properly tell Matlab their dimensions
178: */
179: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
180: #include "engine.h" /* Matlab include file */
181: #include "mex.h" /* Matlab include file */
182: EXTERN_C_BEGIN
183: int VecMatlabEnginePut_DA2d(PetscObject obj,void *engine)
184: {
185: int ierr,n,m;
186: Vec vec = (Vec)obj;
187: PetscScalar *array;
188: mxArray *mat;
189: DA da;
192: PetscObjectQuery((PetscObject)vec,"DA",(PetscObject*)&da);
193: if (!da) SETERRQ(1,"Vector not associated with a DA");
194: DAGetGhostCorners(da,0,0,0,&m,&n,0);
196: VecGetArray(vec,&array);
197: #if !defined(PETSC_USE_COMPLEX)
198: mat = mxCreateDoubleMatrix(m,n,mxREAL);
199: #else
200: mat = mxCreateDoubleMatrix(m,n,mxCOMPLEX);
201: #endif
202: PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));
203: PetscObjectName(obj);
204: mxSetName(mat,obj->name);
205: engPutArray((Engine *)engine,mat);
206:
207: VecRestoreArray(vec,&array);
208: return(0);
209: }
210: EXTERN_C_END
211: #endif
213: /*@C
214: DACreate2d - Creates an object that will manage the communication of two-dimensional
215: regular array data that is distributed across some processors.
217: Collective on MPI_Comm
219: Input Parameters:
220: + comm - MPI communicator
221: . wrap - type of periodicity should the array have.
222: Use one of DA_NONPERIODIC, DA_XPERIODIC, DA_YPERIODIC, or DA_XYPERIODIC.
223: . stencil_type - stencil type. Use either DA_STENCIL_BOX or DA_STENCIL_STAR.
224: . M,N - global dimension in each direction of the array
225: . m,n - corresponding number of processors in each dimension
226: (or PETSC_DECIDE to have calculated)
227: . dof - number of degrees of freedom per node
228: . s - stencil width
229: - lx, ly - arrays containing the number of nodes in each cell along
230: the x and y coordinates, or PETSC_NULL. If non-null, these
231: must be of length as m and n, and the corresponding
232: m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
233: must be M, and the sum of the ly[] entries must be N.
235: Output Parameter:
236: . inra - the resulting distributed array object
238: Options Database Key:
239: + -da_view - Calls DAView() at the conclusion of DACreate2d()
240: . -da_grid_x <nx> - number of grid points in x direction, if M < 0
241: . -da_grid_y <ny> - number of grid points in y direction, if N < 0
242: . -da_processors_x <nx> - number of processors in x direction
243: . -da_processors_y <ny> - number of processors in y direction
244: - -da_noao - do not compute natural to PETSc ordering object
246: Level: beginner
248: Notes:
249: If you are having problems with running out of memory than run with the option -da_noao
251: The stencil type DA_STENCIL_STAR with width 1 corresponds to the
252: standard 5-pt stencil, while DA_STENCIL_BOX with width 1 denotes
253: the standard 9-pt stencil.
255: The array data itself is NOT stored in the DA, it is stored in Vec objects;
256: The appropriate vector objects can be obtained with calls to DACreateGlobalVector()
257: and DACreateLocalVector() and calls to VecDuplicate() if more are needed.
259: .keywords: distributed array, create, two-dimensional
261: .seealso: DADestroy(), DAView(), DACreate1d(), DACreate3d(), DAGlobalToLocalBegin(),
262: DAGlobalToLocalEnd(), DALocalToGlobal(), DALocalToLocalBegin(), DALocalToLocalEnd(),
263: DAGetInfo(), DACreateGlobalVector(), DACreateLocalVector(), DACreateNaturalVector(), DALoad(), DAView()
265: @*/
266: int DACreate2d(MPI_Comm comm,DAPeriodicType wrap,DAStencilType stencil_type,
267: int M,int N,int m,int n,int dof,int s,int *lx,int *ly,DA *inra)
268: {
269: int rank,size,xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,ierr,start,end;
270: int up,down,left,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
271: int xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
272: int s_x,s_y; /* s proportionalized to w */
273: int *flx = 0,*fly = 0;
274: int sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0,refine_x = 2, refine_y = 2,tM = M,tN = N;
275: PetscTruth flg1,flg2;
276: DA da;
277: Vec local,global;
278: VecScatter ltog,gtol;
279: IS to,from;
283: *inra = 0;
284: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
285: DMInitializePackage(PETSC_NULL);
286: #endif
288: if (dof < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %d",dof);
289: if (s < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %d",s);
291: PetscOptionsBegin(comm,PETSC_NULL,"2d DA Options","DA");
292: if (M < 0){
293: tM = -M;
294: PetscOptionsInt("-da_grid_x","Number of grid points in x direction","DACreate2d",tM,&tM,PETSC_NULL);
295: }
296: if (N < 0){
297: tN = -N;
298: PetscOptionsInt("-da_grid_y","Number of grid points in y direction","DACreate2d",tN,&tN,PETSC_NULL);
299: }
300: PetscOptionsInt("-da_processors_x","Number of processors in x direction","DACreate2d",m,&m,PETSC_NULL);
301: PetscOptionsInt("-da_processors_y","Number of processors in y direction","DACreate2d",n,&n,PETSC_NULL);
302: PetscOptionsInt("-da_refine_x","Refinement ratio in x direction","DACreate2d",refine_x,&refine_x,PETSC_NULL);
303: PetscOptionsInt("-da_refine_y","Refinement ratio in y direction","DACreate2d",refine_y,&refine_y,PETSC_NULL);
304: PetscOptionsEnd();
305: M = tM; N = tN;
307: PetscHeaderCreate(da,_p_DA,struct _DAOps,DA_COOKIE,0,"DA",comm,DADestroy,DAView);
308: PetscLogObjectCreate(da);
309: da->bops->publish = DAPublish_Petsc;
310: da->ops->createglobalvector = DACreateGlobalVector;
311: da->ops->getinterpolation = DAGetInterpolation;
312: da->ops->getcoloring = DAGetColoring;
313: da->ops->getmatrix = DAGetMatrix;
314: da->ops->refine = DARefine;
315: PetscLogObjectMemory(da,sizeof(struct _p_DA));
316: da->dim = 2;
317: da->interptype = DA_Q1;
318: da->refine_x = refine_x;
319: da->refine_y = refine_y;
320: PetscMalloc(dof*sizeof(char*),&da->fieldname);
321: PetscMemzero(da->fieldname,dof*sizeof(char*));
323: MPI_Comm_size(comm,&size);
324: MPI_Comm_rank(comm,&rank);
326: if (m != PETSC_DECIDE) {
327: if (m < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %d",m);}
328: else if (m > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %d %d",m,size);}
329: }
330: if (n != PETSC_DECIDE) {
331: if (n < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %d",n);}
332: else if (n > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %d %d",n,size);}
333: }
335: if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
336: /* try for squarish distribution */
337: /* This should use MPI_Dims_create instead */
338: m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)N)));
339: if (!m) m = 1;
340: while (m > 0) {
341: n = size/m;
342: if (m*n == size) break;
343: m--;
344: }
345: if (M > N && m < n) {int _m = m; m = n; n = _m;}
346: if (m*n != size) SETERRQ(PETSC_ERR_PLIB,"Internally Created Bad Partition");
347: } else if (m*n != size) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
349: if (M < m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %d %d",M,m);
350: if (N < n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %d %d",N,n);
352: /*
353: We should create an MPI Cartesian topology here, with reorder
354: set to true. That would create a NEW communicator that we would
355: need to use for operations on this distributed array
356: */
357: PetscOptionsHasName(PETSC_NULL,"-da_partition_nodes_at_end",&flg2);
359: /*
360: Determine locally owned region
361: xs is the first local node number, x is the number of local nodes
362: */
363: if (lx) { /* user sets distribution */
364: x = lx[rank % m];
365: xs = 0;
366: for (i=0; i<(rank % m); i++) {
367: xs += lx[i];
368: }
369: left = xs;
370: for (i=(rank % m); i<m; i++) {
371: left += lx[i];
372: }
373: if (left != M) {
374: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %d %d",left,M);
375: }
376: } else if (flg2) {
377: x = (M + rank%m)/m;
378: if (x < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column width is too thin for stencil! %d %d",x,s);
379: if (M/m == x) { xs = (rank % m)*x; }
380: else { xs = (rank % m)*(x-1) + (M+(rank % m))%(x*m); }
381: SETERRQ(PETSC_ERR_SUP,"-da_partition_nodes_at_end not supported");
382: } else { /* Normal PETSc distribution */
383: x = M/m + ((M % m) > (rank % m));
384: if (x < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column width is too thin for stencil! %d %d",x,s);
385: if ((M % m) > (rank % m)) { xs = (rank % m)*x; }
386: else { xs = (M % m)*(x+1) + ((rank % m)-(M % m))*x; }
387: PetscMalloc(m*sizeof(int),&lx);
388: flx = lx;
389: for (i=0; i<m; i++) {
390: lx[i] = M/m + ((M % m) > i);
391: }
392: }
394: /*
395: Determine locally owned region
396: ys is the first local node number, y is the number of local nodes
397: */
398: if (ly) { /* user sets distribution */
399: y = ly[rank/m];
400: ys = 0;
401: for (i=0; i<(rank/m); i++) {
402: ys += ly[i];
403: }
404: left = ys;
405: for (i=(rank/m); i<n; i++) {
406: left += ly[i];
407: }
408: if (left != N) {
409: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %d %d",left,N);
410: }
411: } else if (flg2) {
412: y = (N + rank/m)/n;
413: if (y < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row width is too thin for stencil! %d %d",y,s);
414: if (N/n == y) { ys = (rank/m)*y; }
415: else { ys = (rank/m)*(y-1) + (N+(rank/m))%(y*n); }
416: SETERRQ(PETSC_ERR_SUP,"-da_partition_nodes_at_end not supported");
417: } else { /* Normal PETSc distribution */
418: y = N/n + ((N % n) > (rank/m));
419: if (y < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row width is too thin for stencil! %d %d",y,s);
420: if ((N % n) > (rank/m)) { ys = (rank/m)*y; }
421: else { ys = (N % n)*(y+1) + ((rank/m)-(N % n))*y; }
422: PetscMalloc(n*sizeof(int),&ly);
423: fly = ly;
424: for (i=0; i<n; i++) {
425: ly[i] = N/n + ((N % n) > i);
426: }
427: }
429: xe = xs + x;
430: ye = ys + y;
432: /* determine ghost region */
433: /* Assume No Periodicity */
434: if (xs-s > 0) Xs = xs - s; else Xs = 0;
435: if (ys-s > 0) Ys = ys - s; else Ys = 0;
436: if (xe+s <= M) Xe = xe + s; else Xe = M;
437: if (ye+s <= N) Ye = ye + s; else Ye = N;
439: /* X Periodic */
440: if (DAXPeriodic(wrap)){
441: Xs = xs - s;
442: Xe = xe + s;
443: }
445: /* Y Periodic */
446: if (DAYPeriodic(wrap)){
447: Ys = ys - s;
448: Ye = ye + s;
449: }
451: /* Resize all X parameters to reflect w */
452: x *= dof;
453: xs *= dof;
454: xe *= dof;
455: Xs *= dof;
456: Xe *= dof;
457: s_x = s*dof;
458: s_y = s;
460: /* determine starting point of each processor */
461: nn = x*y;
462: PetscMalloc((2*size+1)*sizeof(int),&bases);
463: ldims = (int*)(bases+size+1);
464: MPI_Allgather(&nn,1,MPI_INT,ldims,1,MPI_INT,comm);
465: bases[0] = 0;
466: for (i=1; i<=size; i++) {
467: bases[i] = ldims[i-1];
468: }
469: for (i=1; i<=size; i++) {
470: bases[i] += bases[i-1];
471: }
473: /* allocate the base parallel and sequential vectors */
474: VecCreateMPI(comm,x*y,PETSC_DECIDE,&global);
475: VecSetBlockSize(global,dof);
476: VecCreateSeq(PETSC_COMM_SELF,(Xe-Xs)*(Ye-Ys),&local);
477: VecSetBlockSize(local,dof);
480: /* generate appropriate vector scatters */
481: /* local to global inserts non-ghost point region into global */
482: VecGetOwnershipRange(global,&start,&end);
483: ISCreateStride(comm,x*y,start,1,&to);
485: left = xs - Xs; down = ys - Ys; up = down + y;
486: PetscMalloc(x*(up - down)*sizeof(int),&idx);
487: count = 0;
488: for (i=down; i<up; i++) {
489: for (j=0; j<x; j++) {
490: idx[count++] = left + i*(Xe-Xs) + j;
491: }
492: }
493: ISCreateGeneral(comm,count,idx,&from);
494: PetscFree(idx);
496: VecScatterCreate(local,from,global,to,<og);
497: PetscLogObjectParent(da,to);
498: PetscLogObjectParent(da,from);
499: PetscLogObjectParent(da,ltog);
500: ISDestroy(from);
501: ISDestroy(to);
503: /* global to local must include ghost points */
504: if (stencil_type == DA_STENCIL_BOX) {
505: ISCreateStride(comm,(Xe-Xs)*(Ye-Ys),0,1,&to);
506: } else {
507: /* must drop into cross shape region */
508: /* ---------|
509: | top |
510: |--- ---|
511: | middle |
512: | |
513: ---- ----
514: | bottom |
515: -----------
516: Xs xs xe Xe */
517: /* bottom */
518: left = xs - Xs; down = ys - Ys; up = down + y;
519: count = down*(xe-xs) + (up-down)*(Xe-Xs) + (Ye-Ys-up)*(xe-xs);
520: ierr = PetscMalloc(count*sizeof(int),&idx);
521: count = 0;
522: for (i=0; i<down; i++) {
523: for (j=0; j<xe-xs; j++) {
524: idx[count++] = left + i*(Xe-Xs) + j;
525: }
526: }
527: /* middle */
528: for (i=down; i<up; i++) {
529: for (j=0; j<Xe-Xs; j++) {
530: idx[count++] = i*(Xe-Xs) + j;
531: }
532: }
533: /* top */
534: for (i=up; i<Ye-Ys; i++) {
535: for (j=0; j<xe-xs; j++) {
536: idx[count++] = left + i*(Xe-Xs) + j;
537: }
538: }
539: ISCreateGeneral(comm,count,idx,&to);
540: PetscFree(idx);
541: }
544: /* determine who lies on each side of us stored in n6 n7 n8
545: n3 n5
546: n0 n1 n2
547: */
549: /* Assume the Non-Periodic Case */
550: n1 = rank - m;
551: if (rank % m) {
552: n0 = n1 - 1;
553: } else {
554: n0 = -1;
555: }
556: if ((rank+1) % m) {
557: n2 = n1 + 1;
558: n5 = rank + 1;
559: n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
560: } else {
561: n2 = -1; n5 = -1; n8 = -1;
562: }
563: if (rank % m) {
564: n3 = rank - 1;
565: n6 = n3 + m; if (n6 >= m*n) n6 = -1;
566: } else {
567: n3 = -1; n6 = -1;
568: }
569: n7 = rank + m; if (n7 >= m*n) n7 = -1;
572: /* Modify for Periodic Cases */
573: if (wrap == DA_YPERIODIC) { /* Handle Top and Bottom Sides */
574: if (n1 < 0) n1 = rank + m * (n-1);
575: if (n7 < 0) n7 = rank - m * (n-1);
576: if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
577: if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
578: if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
579: if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
580: } else if (wrap == DA_XPERIODIC) { /* Handle Left and Right Sides */
581: if (n3 < 0) n3 = rank + (m-1);
582: if (n5 < 0) n5 = rank - (m-1);
583: if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
584: if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
585: if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
586: if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
587: } else if (wrap == DA_XYPERIODIC) {
589: /* Handle all four corners */
590: if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
591: if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
592: if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
593: if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
595: /* Handle Top and Bottom Sides */
596: if (n1 < 0) n1 = rank + m * (n-1);
597: if (n7 < 0) n7 = rank - m * (n-1);
598: if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
599: if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
600: if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
601: if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
603: /* Handle Left and Right Sides */
604: if (n3 < 0) n3 = rank + (m-1);
605: if (n5 < 0) n5 = rank - (m-1);
606: if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
607: if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
608: if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
609: if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
610: }
612: if (stencil_type == DA_STENCIL_STAR) {
613: /* save corner processor numbers */
614: sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
615: n0 = n2 = n6 = n8 = -1;
616: }
618: PetscMalloc((x+2*s_x)*(y+2*s_y)*sizeof(int),&idx);
619: PetscLogObjectMemory(da,(x+2*s_x)*(y+2*s_y)*sizeof(int));
620: nn = 0;
622: xbase = bases[rank];
623: for (i=1; i<=s_y; i++) {
624: if (n0 >= 0) { /* left below */
625: x_t = lx[n0 % m]*dof;
626: y_t = ly[(n0/m)];
627: s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
628: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
629: }
630: if (n1 >= 0) { /* directly below */
631: x_t = x;
632: y_t = ly[(n1/m)];
633: s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
634: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
635: }
636: if (n2 >= 0) { /* right below */
637: x_t = lx[n2 % m]*dof;
638: y_t = ly[(n2/m)];
639: s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
640: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
641: }
642: }
644: for (i=0; i<y; i++) {
645: if (n3 >= 0) { /* directly left */
646: x_t = lx[n3 % m]*dof;
647: /* y_t = y; */
648: s_t = bases[n3] + (i+1)*x_t - s_x;
649: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
650: }
652: for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
654: if (n5 >= 0) { /* directly right */
655: x_t = lx[n5 % m]*dof;
656: /* y_t = y; */
657: s_t = bases[n5] + (i)*x_t;
658: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
659: }
660: }
662: for (i=1; i<=s_y; i++) {
663: if (n6 >= 0) { /* left above */
664: x_t = lx[n6 % m]*dof;
665: /* y_t = ly[(n6/m)]; */
666: s_t = bases[n6] + (i)*x_t - s_x;
667: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
668: }
669: if (n7 >= 0) { /* directly above */
670: x_t = x;
671: /* y_t = ly[(n7/m)]; */
672: s_t = bases[n7] + (i-1)*x_t;
673: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
674: }
675: if (n8 >= 0) { /* right above */
676: x_t = lx[n8 % m]*dof;
677: /* y_t = ly[(n8/m)]; */
678: s_t = bases[n8] + (i-1)*x_t;
679: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
680: }
681: }
683: base = bases[rank];
684: ISCreateGeneral(comm,nn,idx,&from);
685: VecScatterCreate(global,from,local,to,>ol);
686: PetscLogObjectParent(da,to);
687: PetscLogObjectParent(da,from);
688: PetscLogObjectParent(da,gtol);
689: ISDestroy(to);
690: ISDestroy(from);
692: if (stencil_type == DA_STENCIL_STAR) {
693: /*
694: Recompute the local to global mappings, this time keeping the
695: information about the cross corner processor numbers.
696: */
697: n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
698: nn = 0;
699: xbase = bases[rank];
700: for (i=1; i<=s_y; i++) {
701: if (n0 >= 0) { /* left below */
702: x_t = lx[n0 % m]*dof;
703: y_t = ly[(n0/m)];
704: s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
705: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
706: }
707: if (n1 >= 0) { /* directly below */
708: x_t = x;
709: y_t = ly[(n1/m)];
710: s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
711: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
712: }
713: if (n2 >= 0) { /* right below */
714: x_t = lx[n2 % m]*dof;
715: y_t = ly[(n2/m)];
716: s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
717: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
718: }
719: }
721: for (i=0; i<y; i++) {
722: if (n3 >= 0) { /* directly left */
723: x_t = lx[n3 % m]*dof;
724: /* y_t = y; */
725: s_t = bases[n3] + (i+1)*x_t - s_x;
726: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
727: }
729: for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
731: if (n5 >= 0) { /* directly right */
732: x_t = lx[n5 % m]*dof;
733: /* y_t = y; */
734: s_t = bases[n5] + (i)*x_t;
735: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
736: }
737: }
739: for (i=1; i<=s_y; i++) {
740: if (n6 >= 0) { /* left above */
741: x_t = lx[n6 % m]*dof;
742: /* y_t = ly[(n6/m)]; */
743: s_t = bases[n6] + (i)*x_t - s_x;
744: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
745: }
746: if (n7 >= 0) { /* directly above */
747: x_t = x;
748: /* y_t = ly[(n7/m)]; */
749: s_t = bases[n7] + (i-1)*x_t;
750: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
751: }
752: if (n8 >= 0) { /* right above */
753: x_t = lx[n8 % m]*dof;
754: /* y_t = ly[(n8/m)]; */
755: s_t = bases[n8] + (i-1)*x_t;
756: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
757: }
758: }
759: }
761: da->M = M; da->N = N; da->m = m; da->n = n; da->w = dof; da->s = s;
762: da->xs = xs; da->xe = xe; da->ys = ys; da->ye = ye; da->zs = 0; da->ze = 1;
763: da->Xs = Xs; da->Xe = Xe; da->Ys = Ys; da->Ye = Ye; da->Zs = 0; da->Ze = 1;
764: da->P = 1; da->p = 1;
766: PetscLogObjectParent(da,global);
767: PetscLogObjectParent(da,local);
769: da->global = global;
770: da->local = local;
771: da->gtol = gtol;
772: da->ltog = ltog;
773: da->idx = idx;
774: da->Nl = nn;
775: da->base = base;
776: da->wrap = wrap;
777: da->ops->view = DAView_2d;
778: da->stencil_type = stencil_type;
780: /*
781: Set the local to global ordering in the global vector, this allows use
782: of VecSetValuesLocal().
783: */
784: ISLocalToGlobalMappingCreate(comm,nn,idx,&da->ltogmap);
785: VecSetLocalToGlobalMapping(da->global,da->ltogmap);
786: ISLocalToGlobalMappingBlock(da->ltogmap,da->w,&da->ltogmapb);
787: VecSetLocalToGlobalMappingBlock(da->global,da->ltogmapb);
788: PetscLogObjectParent(da,da->ltogmap);
790: *inra = da;
792: /* recalculate the idx including missed ghost points */
793: /* Assume the Non-Periodic Case */
794: n1 = rank - m;
795: if (rank % m) {
796: n0 = n1 - 1;
797: } else {
798: n0 = -1;
799: }
800: if ((rank+1) % m) {
801: n2 = n1 + 1;
802: n5 = rank + 1;
803: n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
804: } else {
805: n2 = -1; n5 = -1; n8 = -1;
806: }
807: if (rank % m) {
808: n3 = rank - 1;
809: n6 = n3 + m; if (n6 >= m*n) n6 = -1;
810: } else {
811: n3 = -1; n6 = -1;
812: }
813: n7 = rank + m; if (n7 >= m*n) n7 = -1;
816: /* Modify for Periodic Cases */
817: if (wrap == DA_YPERIODIC) { /* Handle Top and Bottom Sides */
818: if (n1 < 0) n1 = rank + m * (n-1);
819: if (n7 < 0) n7 = rank - m * (n-1);
820: if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
821: if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
822: if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
823: if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
824: } else if (wrap == DA_XPERIODIC) { /* Handle Left and Right Sides */
825: if (n3 < 0) n3 = rank + (m-1);
826: if (n5 < 0) n5 = rank - (m-1);
827: if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
828: if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
829: if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
830: if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
831: } else if (wrap == DA_XYPERIODIC) {
833: /* Handle all four corners */
834: if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
835: if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
836: if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
837: if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
839: /* Handle Top and Bottom Sides */
840: if (n1 < 0) n1 = rank + m * (n-1);
841: if (n7 < 0) n7 = rank - m * (n-1);
842: if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
843: if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
844: if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
845: if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
847: /* Handle Left and Right Sides */
848: if (n3 < 0) n3 = rank + (m-1);
849: if (n5 < 0) n5 = rank - (m-1);
850: if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
851: if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
852: if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
853: if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
854: }
856: nn = 0;
858: xbase = bases[rank];
859: for (i=1; i<=s_y; i++) {
860: if (n0 >= 0) { /* left below */
861: x_t = lx[n0 % m]*dof;
862: y_t = ly[(n0/m)];
863: s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
864: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
865: }
866: if (n1 >= 0) { /* directly below */
867: x_t = x;
868: y_t = ly[(n1/m)];
869: s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
870: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
871: }
872: if (n2 >= 0) { /* right below */
873: x_t = lx[n2 % m]*dof;
874: y_t = ly[(n2/m)];
875: s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
876: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
877: }
878: }
880: for (i=0; i<y; i++) {
881: if (n3 >= 0) { /* directly left */
882: x_t = lx[n3 % m]*dof;
883: /* y_t = y; */
884: s_t = bases[n3] + (i+1)*x_t - s_x;
885: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
886: }
888: for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
890: if (n5 >= 0) { /* directly right */
891: x_t = lx[n5 % m]*dof;
892: /* y_t = y; */
893: s_t = bases[n5] + (i)*x_t;
894: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
895: }
896: }
898: for (i=1; i<=s_y; i++) {
899: if (n6 >= 0) { /* left above */
900: x_t = lx[n6 % m]*dof;
901: /* y_t = ly[(n6/m)]; */
902: s_t = bases[n6] + (i)*x_t - s_x;
903: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
904: }
905: if (n7 >= 0) { /* directly above */
906: x_t = x;
907: /* y_t = ly[(n7/m)]; */
908: s_t = bases[n7] + (i-1)*x_t;
909: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
910: }
911: if (n8 >= 0) { /* right above */
912: x_t = lx[n8 % m]*dof;
913: /* y_t = ly[(n8/m)]; */
914: s_t = bases[n8] + (i-1)*x_t;
915: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
916: }
917: }
918: PetscFree(bases);
920: /* construct the local to local scatter context */
921: /*
922: We simply remap the values in the from part of
923: global to local to read from an array with the ghost values
924: rather then from the plan array.
925: */
926: VecScatterCopy(gtol,&da->ltol);
927: PetscLogObjectParent(da,da->ltol);
928: left = xs - Xs; down = ys - Ys; up = down + y;
929: PetscMalloc(x*(up - down)*sizeof(int),&idx);
930: count = 0;
931: for (i=down; i<up; i++) {
932: for (j=0; j<x; j++) {
933: idx[count++] = left + i*(Xe-Xs) + j;
934: }
935: }
936: VecScatterRemap(da->ltol,idx,PETSC_NULL);
937: PetscFree(idx);
939: /*
940: Build the natural ordering to PETSc ordering mappings.
941: */
942: PetscOptionsHasName(PETSC_NULL,"-da_noao",&flg1);
943: if (!flg1) {
944: IS ispetsc,isnatural;
945: int *lidx,lict = 0,Nlocal = (da->xe-da->xs)*(da->ye-da->ys);
947: ISCreateStride(comm,Nlocal,da->base,1,&ispetsc);
949: PetscMalloc(Nlocal*sizeof(int),&lidx);
950: for (j=ys; j<ye; j++) {
951: for (i=xs; i<xe; i++) {
952: /* global number in natural ordering */
953: lidx[lict++] = i + j*M*dof;
954: }
955: }
956: ISCreateGeneral(comm,Nlocal,lidx,&isnatural);
957: PetscFree(lidx);
960: AOCreateBasicIS(isnatural,ispetsc,&da->ao);
961: PetscLogObjectParent(da,da->ao);
962: ISDestroy(ispetsc);
963: ISDestroy(isnatural);
964: } else {
965: da->ao = PETSC_NULL;
966: }
968: if (!flx) {
969: PetscMalloc(m*sizeof(int),&flx);
970: PetscMemcpy(flx,lx,m*sizeof(int));
971: }
972: if (!fly) {
973: PetscMalloc(n*sizeof(int),&fly);
974: PetscMemcpy(fly,ly,n*sizeof(int));
975: }
976: da->lx = flx;
977: da->ly = fly;
979: PetscOptionsHasName(PETSC_NULL,"-da_view",&flg1);
980: if (flg1) {DAView(da,PETSC_VIEWER_STDOUT_(da->comm));}
981: PetscOptionsHasName(PETSC_NULL,"-da_view_draw",&flg1);
982: if (flg1) {DAView(da,PETSC_VIEWER_DRAW_(da->comm));}
983: PetscOptionsHasName(PETSC_NULL,"-help",&flg1);
984: if (flg1) {DAPrintHelp(da);}
986: PetscPublishAll(da);
987: #if defined(PETSC_HAVE_AMS)
988: PetscObjectComposeFunctionDynamic((PetscObject)global,"AMSSetFieldBlock_C",
989: "AMSSetFieldBlock_DA",AMSSetFieldBlock_DA);
990: PetscObjectComposeFunctionDynamic((PetscObject)local,"AMSSetFieldBlock_C",
991: "AMSSetFieldBlock_DA",AMSSetFieldBlock_DA);
992: if (((PetscObject)global)->amem > -1) {
993: AMSSetFieldBlock_DA(((PetscObject)global)->amem,"values",global);
994: }
995: #endif
996: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
997: if (dof == 1) {
998: PetscObjectComposeFunctionDynamic((PetscObject)local,"PetscMatlabEnginePut_C","VecMatlabEnginePut_DA2d",VecMatlabEnginePut_DA2d);
999: }
1000: #endif
1001: VecSetOperation(global,VECOP_VIEW,(void(*)(void))VecView_MPI_DA);
1002: VecSetOperation(global,VECOP_LOADINTOVECTOR,(void(*)(void))VecLoadIntoVector_Binary_DA);
1003: return(0);
1004: }
1006: /*@
1007: DAPrintHelp - Prints command line options for DA.
1009: Collective on DA
1011: Input Parameters:
1012: . da - the distributed array
1014: Level: intermediate
1016: .seealso: DACreate1d(), DACreate2d(), DACreate3d()
1018: .keywords: DA, help
1020: @*/
1021: int DAPrintHelp(DA da)
1022: {
1023: static PetscTruth called = PETSC_FALSE;
1024: MPI_Comm comm;
1025: int ierr;
1030: comm = da->comm;
1031: if (!called) {
1032: (*PetscHelpPrintf)(comm,"General Distributed Array (DA) options:n");
1033: (*PetscHelpPrintf)(comm," -da_view: print DA distribution to screenn");
1034: (*PetscHelpPrintf)(comm," -da_view_draw: display DA in windown");
1035: called = PETSC_TRUE;
1036: }
1037: return(0);
1038: }
1040: /*@C
1041: DARefine - Creates a new distributed array that is a refinement of a given
1042: distributed array.
1044: Collective on DA
1046: Input Parameter:
1047: + da - initial distributed array
1048: - comm - communicator to contain refined DA, must be either same as the da communicator or include the
1049: da communicator and be 2, 4, or 8 times larger. Currently ignored
1051: Output Parameter:
1052: . daref - refined distributed array
1054: Level: advanced
1056: Note:
1057: Currently, refinement consists of just doubling the number of grid spaces
1058: in each dimension of the DA.
1060: .keywords: distributed array, refine
1062: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy()
1063: @*/
1064: int DARefine(DA da,MPI_Comm comm,DA *daref)
1065: {
1066: int M,N,P,ierr;
1067: DA da2;
1072: if (DAXPeriodic(da->wrap) || da->interptype == DA_Q0){
1073: M = da->refine_x*da->M;
1074: } else {
1075: M = 1 + da->refine_x*(da->M - 1);
1076: }
1077: if (DAYPeriodic(da->wrap) || da->interptype == DA_Q0){
1078: N = da->refine_y*da->N;
1079: } else {
1080: N = 1 + da->refine_y*(da->N - 1);
1081: }
1082: if (DAZPeriodic(da->wrap) || da->interptype == DA_Q0){
1083: P = da->refine_z*da->P;
1084: } else {
1085: P = 1 + da->refine_z*(da->P - 1);
1086: }
1087: if (da->dim == 1) {
1088: DACreate1d(da->comm,da->wrap,M,da->w,da->s,PETSC_NULL,&da2);
1089: } else if (da->dim == 2) {
1090: DACreate2d(da->comm,da->wrap,da->stencil_type,M,N,da->m,da->n,da->w,da->s,PETSC_NULL,PETSC_NULL,&da2);
1091: } else if (da->dim == 3) {
1092: DACreate3d(da->comm,da->wrap,da->stencil_type,M,N,P,da->m,da->n,da->p,da->w,da->s,0,0,0,&da2);
1093: }
1094: *daref = da2;
1095: return(0);
1096: }
1098: /*
1099: M is number of grid points
1100: m is number of processors
1102: */
1103: int DASplitComm2d(MPI_Comm comm,int M,int N,int sw,MPI_Comm *outcomm)
1104: {
1105: int ierr,m,n = 0,csize,size,rank,x = 0,y = 0;
1108: MPI_Comm_size(comm,&size);
1109: MPI_Comm_rank(comm,&rank);
1111: csize = 4*size;
1112: do {
1113: if (csize % 4) SETERRQ4(1,"Cannot split communicator of size %d tried %d %d %d",size,csize,x,y);
1114: csize = csize/4;
1115:
1116: m = (int)(0.5 + sqrt(((double)M)*((double)csize)/((double)N)));
1117: if (!m) m = 1;
1118: while (m > 0) {
1119: n = csize/m;
1120: if (m*n == csize) break;
1121: m--;
1122: }
1123: if (M > N && m < n) {int _m = m; m = n; n = _m;}
1125: x = M/m + ((M % m) > ((csize-1) % m));
1126: y = (N + (csize-1)/m)/n;
1127: } while ((x < 4 || y < 4) && csize > 1);
1128: if (size != csize) {
1129: MPI_Group entire_group,sub_group;
1130: int i,*groupies;
1132: ierr = MPI_Comm_group(comm,&entire_group);
1133: PetscMalloc(csize*sizeof(int),&groupies);
1134: for (i=0; i<csize; i++) {
1135: groupies[i] = (rank/csize)*csize + i;
1136: }
1137: ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);
1138: ierr = PetscFree(groupies);
1139: ierr = MPI_Comm_create(comm,sub_group,outcomm);
1140: ierr = MPI_Group_free(&entire_group);
1141: ierr = MPI_Group_free(&sub_group);
1142: PetscLogInfo(0,"Creating redundant coarse problems of size %dn",csize);
1143: } else {
1144: *outcomm = comm;
1145: }
1146: return(0);
1147: }
1149: /*@C
1150: DASetLocalFunction - Caches in a DA a local function
1152: Collective on DA
1154: Input Parameter:
1155: + da - initial distributed array
1156: - lf - the local function
1158: Level: intermediate
1160: .keywords: distributed array, refine
1162: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunctioni()
1163: @*/
1164: int DASetLocalFunction(DA da,DALocalFunction1 lf)
1165: {
1168: da->lf = lf;
1169: return(0);
1170: }
1172: /*@C
1173: DASetLocalFunctioni - Caches in a DA a local function that evaluates a single component
1175: Collective on DA
1177: Input Parameter:
1178: + da - initial distributed array
1179: - lfi - the local function
1181: Level: intermediate
1183: .keywords: distributed array, refine
1185: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction()
1186: @*/
1187: int DASetLocalFunctioni(DA da,int (*lfi)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
1188: {
1191: da->lfi = lfi;
1192: return(0);
1193: }
1195: /*MC
1196: DASetLocalAdicFunction - Caches in a DA a local function computed by ADIC/ADIFOR
1198: Collective on DA
1200: Synopsis:
1201: int int DASetLocalAdicFunction(DA da,DALocalFunction1 ad_lf)
1202:
1203: Input Parameter:
1204: + da - initial distributed array
1205: - ad_lf - the local function as computed by ADIC/ADIFOR
1207: Level: intermediate
1209: .keywords: distributed array, refine
1211: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
1212: DASetLocalJacobian()
1213: M*/
1215: int DASetLocalAdicFunction_Private(DA da,DALocalFunction1 ad_lf)
1216: {
1219: da->adic_lf = ad_lf;
1220: return(0);
1221: }
1223: /*MC
1224: DASetLocalAdicFunctioni - Caches in a DA a local functioni computed by ADIC/ADIFOR
1226: Collective on DA
1228: Synopsis:
1229: int int DASetLocalAdicFunctioni(DA da,int (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*)
1230:
1231: Input Parameter:
1232: + da - initial distributed array
1233: - ad_lfi - the local function as computed by ADIC/ADIFOR
1235: Level: intermediate
1237: .keywords: distributed array, refine
1239: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
1240: DASetLocalJacobian(), DASetLocalFunctioni()
1241: M*/
1243: int DASetLocalAdicFunctioni_Private(DA da,int (*ad_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*))
1244: {
1247: da->adic_lfi = ad_lfi;
1248: return(0);
1249: }
1251: /*MC
1252: DASetLocalAdicMFFunctioni - Caches in a DA a local functioni computed by ADIC/ADIFOR
1254: Collective on DA
1256: Synopsis:
1257: int int DASetLocalAdicFunctioni(DA da,int (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*)
1258:
1259: Input Parameter:
1260: + da - initial distributed array
1261: - admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
1263: Level: intermediate
1265: .keywords: distributed array, refine
1267: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
1268: DASetLocalJacobian(), DASetLocalFunctioni()
1269: M*/
1271: int DASetLocalAdicMFFunctioni_Private(DA da,int (*admf_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*))
1272: {
1275: da->adicmf_lfi = admf_lfi;
1276: return(0);
1277: }
1279: /*MC
1280: DASetLocalAdicMFFunction - Caches in a DA a local function computed by ADIC/ADIFOR
1282: Collective on DA
1284: Synopsis:
1285: int int DASetLocalAdicMFFunction(DA da,DALocalFunction1 ad_lf)
1286:
1287: Input Parameter:
1288: + da - initial distributed array
1289: - ad_lf - the local function as computed by ADIC/ADIFOR
1291: Level: intermediate
1293: .keywords: distributed array, refine
1295: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
1296: DASetLocalJacobian()
1297: M*/
1299: int DASetLocalAdicMFFunction_Private(DA da,DALocalFunction1 ad_lf)
1300: {
1303: da->adicmf_lf = ad_lf;
1304: return(0);
1305: }
1307: /*@C
1308: DASetLocalJacobian - Caches in a DA a local Jacobian
1310: Collective on DA
1312:
1313: Input Parameter:
1314: + da - initial distributed array
1315: - lj - the local Jacobian
1317: Level: intermediate
1319: .keywords: distributed array, refine
1321: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction()
1322: @*/
1323: int DASetLocalJacobian(DA da,DALocalFunction1 lj)
1324: {
1327: da->lj = lj;
1328: return(0);
1329: }
1331: /*@C
1332: DAGetLocalFunction - Gets from a DA a local function and its ADIC/ADIFOR Jacobian
1334: Collective on DA
1336: Input Parameter:
1337: . da - initial distributed array
1339: Output Parameters:
1340: . lf - the local function
1342: Level: intermediate
1344: .keywords: distributed array, refine
1346: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DASetLocalFunction()
1347: @*/
1348: int DAGetLocalFunction(DA da,DALocalFunction1 *lf)
1349: {
1352: if (lf) *lf = da->lf;
1353: return(0);
1354: }
1356: /*@
1357: DAFormFunction1 - Evaluates a user provided function on each processor that
1358: share a DA
1360: Input Parameters:
1361: + da - the DA that defines the grid
1362: . vu - input vector
1363: . vfu - output vector
1364: - w - any user data
1366: Notes: Does NOT do ghost updates on vu upon entry
1368: Level: advanced
1370: .seealso: DAComputeJacobian1WithAdic()
1372: @*/
1373: int DAFormFunction1(DA da,Vec vu,Vec vfu,void *w)
1374: {
1375: int ierr;
1376: void *u,*fu;
1377: DALocalInfo info;
1378:
1381: DAGetLocalInfo(da,&info);
1382: DAVecGetArray(da,vu,(void**)&u);
1383: DAVecGetArray(da,vfu,(void**)&fu);
1385: (*da->lf)(&info,u,fu,w);
1387: DAVecRestoreArray(da,vu,(void**)&u);
1388: DAVecRestoreArray(da,vfu,(void**)&fu);
1389: return(0);
1390: }
1392: int DAFormFunctioniTest1(DA da,void *w)
1393: {
1394: Vec vu,fu,fui;
1395: int ierr,i,n;
1396: PetscScalar *ui,mone = -1.0;
1397: PetscRandom rnd;
1398: PetscReal norm;
1401: DAGetLocalVector(da,&vu);
1402: PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&rnd);
1403: VecSetRandom(rnd,vu);
1404: PetscRandomDestroy(rnd);
1406: DAGetGlobalVector(da,&fu);
1407: DAGetGlobalVector(da,&fui);
1408:
1409: DAFormFunction1(da,vu,fu,w);
1411: VecGetArray(fui,&ui);
1412: VecGetLocalSize(fui,&n);
1413: for (i=0; i<n; i++) {
1414: DAFormFunctioni1(da,i,vu,ui+i,w);
1415: }
1416: VecRestoreArray(fui,&ui);
1418: VecAXPY(&mone,fu,fui);
1419: VecNorm(fui,NORM_2,&norm);
1420: PetscPrintf(da->comm,"Norm of difference in vectors %gn",norm);
1421: VecView(fu,0);
1422: VecView(fui,0);
1424: DARestoreLocalVector(da,&vu);
1425: DARestoreGlobalVector(da,&fu);
1426: DARestoreGlobalVector(da,&fui);
1427: return(0);
1428: }
1430: /*@
1431: DAFormFunctioni1 - Evaluates a user provided function
1433: Input Parameters:
1434: + da - the DA that defines the grid
1435: . i - the component of the function we wish to compute (must be local)
1436: . vu - input vector
1437: . vfu - output value
1438: - w - any user data
1440: Notes: Does NOT do ghost updates on vu upon entry
1442: Level: advanced
1444: .seealso: DAComputeJacobian1WithAdic()
1446: @*/
1447: int DAFormFunctioni1(DA da,int i,Vec vu,PetscScalar *vfu,void *w)
1448: {
1449: int ierr;
1450: void *u;
1451: DALocalInfo info;
1452: MatStencil stencil;
1453:
1456: DAGetLocalInfo(da,&info);
1457: DAVecGetArray(da,vu,(void**)&u);
1459: /* figure out stencil value from i */
1460: stencil.c = i % info.dof;
1461: stencil.i = (i % (info.xm*info.dof))/info.dof;
1462: stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
1463: stencil.k = i/(info.xm*info.ym*info.dof);
1465: (*da->lfi)(&info,&stencil,u,vfu,w);
1467: DAVecRestoreArray(da,vu,(void**)&u);
1468: return(0);
1469: }
1471: #if defined(new)
1472: /*
1473: DAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
1474: function lives on a DA
1476: y ~= (F(u + ha) - F(u))/h,
1477: where F = nonlinear function, as set by SNESSetFunction()
1478: u = current iterate
1479: h = difference interval
1480: */
1481: int DAGetDiagonal_MFFD(DA da,Vec U,Vec a)
1482: {
1483: PetscScalar h,*aa,*ww,v;
1484: PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
1485: int ierr,gI,nI;
1486: MatStencil stencil;
1487: DALocalInfo info;
1488:
1490: (*ctx->func)(0,U,a,ctx->funcctx);
1491: (*ctx->funcisetbase)(U,ctx->funcctx);
1493: VecGetArray(U,&ww);
1494: VecGetArray(a,&aa);
1495:
1496: nI = 0;
1497: h = ww[gI];
1498: if (h == 0.0) h = 1.0;
1499: #if !defined(PETSC_USE_COMPLEX)
1500: if (h < umin && h >= 0.0) h = umin;
1501: else if (h < 0.0 && h > -umin) h = -umin;
1502: #else
1503: if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin;
1504: else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
1505: #endif
1506: h *= epsilon;
1507:
1508: ww[gI += h;
1509: ierr = (*ctx->funci)(i,w,&v,ctx->funcctx);
1510: aa[nI] = (v - aa[nI])/h;
1511: ww[gI] -= h;
1512: nI++;
1513: }
1514: VecRestoreArray(U,&ww);
1515: VecRestoreArray(a,&aa);
1516: return(0);
1517: }
1518: #endif
1520: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
1521: EXTERN_C_BEGIN
1522: #include "adic/ad_utils.h"
1523: EXTERN_C_END
1525: /*@
1526: DAComputeJacobian1WithAdic - Evaluates a adiC provided Jacobian function on each processor that
1527: share a DA
1529: Input Parameters:
1530: + da - the DA that defines the grid
1531: . vu - input vector (ghosted)
1532: . J - output matrix
1533: - w - any user data
1535: Level: advanced
1537: Notes: Does NOT do ghost updates on vu upon entry
1539: .seealso: DAFormFunction1()
1541: @*/
1542: int DAComputeJacobian1WithAdic(DA da,Vec vu,Mat J,void *w)
1543: {
1544: int ierr,gtdof,tdof;
1545: PetscScalar *ustart;
1546: DALocalInfo info;
1547: void *ad_u,*ad_f,*ad_ustart,*ad_fstart;
1548: ISColoring iscoloring;
1551: DAGetLocalInfo(da,&info);
1553: /* get space for derivative objects. */
1554: DAGetAdicArray(da,PETSC_TRUE,(void **)&ad_u,&ad_ustart,>dof);
1555: DAGetAdicArray(da,PETSC_FALSE,(void **)&ad_f,&ad_fstart,&tdof);
1556: VecGetArray(vu,&ustart);
1557: PetscADSetValArray(((DERIV_TYPE*)ad_ustart),gtdof,ustart);
1558: VecRestoreArray(vu,&ustart);
1560: PetscADResetIndep();
1561: DAGetColoring(da,IS_COLORING_GHOSTED,&iscoloring);
1562: PetscADSetIndepArrayColored(ad_ustart,gtdof,iscoloring->colors);
1563: PetscADIncrementTotalGradSize(iscoloring->n);
1564: ISColoringDestroy(iscoloring);
1565: PetscADSetIndepDone();
1567: (*da->adic_lf)(&info,ad_u,ad_f,w);
1569: /* stick the values into the matrix */
1570: MatSetValuesAdic(J,(PetscScalar**)ad_fstart);
1572: /* return space for derivative objects. */
1573: DARestoreAdicArray(da,PETSC_TRUE,(void **)&ad_u,&ad_ustart,>dof);
1574: DARestoreAdicArray(da,PETSC_FALSE,(void **)&ad_f,&ad_fstart,&tdof);
1575: return(0);
1576: }
1578: /*@C
1579: DAMultiplyByJacobian1WithAdic - Applies an ADIC-provided Jacobian function to a vector on
1580: each processor that shares a DA.
1582: Input Parameters:
1583: + da - the DA that defines the grid
1584: . vu - Jacobian is computed at this point (ghosted)
1585: . v - product is done on this vector (ghosted)
1586: . fu - output vector = J(vu)*v (not ghosted)
1587: - w - any user data
1589: Notes:
1590: This routine does NOT do ghost updates on vu upon entry.
1592: Level: advanced
1594: .seealso: DAFormFunction1()
1596: @*/
1597: int DAMultiplyByJacobian1WithAdic(DA da,Vec vu,Vec v,Vec f,void *w)
1598: {
1599: int ierr,i,gtdof,tdof;
1600: PetscScalar *avu,*av,*af,*ad_vustart,*ad_fstart;
1601: DALocalInfo info;
1602: void *ad_vu,*ad_f;
1605: DAGetLocalInfo(da,&info);
1607: /* get space for derivative objects. */
1608: DAGetAdicMFArray(da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,>dof);
1609: DAGetAdicMFArray(da,PETSC_FALSE,(void **)&ad_f,(void**)&ad_fstart,&tdof);
1611: /* copy input vector into derivative object */
1612: VecGetArray(vu,&avu);
1613: VecGetArray(v,&av);
1614: for (i=0; i<gtdof; i++) {
1615: ad_vustart[2*i] = avu[i];
1616: ad_vustart[2*i+1] = av[i];
1617: }
1618: VecRestoreArray(vu,&avu);
1619: VecRestoreArray(v,&av);
1621: PetscADResetIndep();
1622: PetscADIncrementTotalGradSize(1);
1623: PetscADSetIndepDone();
1625: (*da->adicmf_lf)(&info,ad_vu,ad_f,w);
1627: /* stick the values into the vector */
1628: VecGetArray(f,&af);
1629: for (i=0; i<tdof; i++) {
1630: af[i] = ad_fstart[2*i+1];
1631: }
1632: VecRestoreArray(f,&af);
1634: /* return space for derivative objects. */
1635: DARestoreAdicMFArray(da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,>dof);
1636: DARestoreAdicMFArray(da,PETSC_FALSE,(void **)&ad_f,(void**)&ad_fstart,&tdof);
1637: return(0);
1638: }
1641: #else
1643: int DAComputeJacobian1WithAdic(DA da,Vec vu,Mat J,void *w)
1644: {
1646: SETERRQ(1,"Must compile with bmake/PETSC_ARCH/packages flag PETSC_HAVE_ADIC for this routine");
1647: }
1649: int DAMultiplyByJacobian1WithAdic(DA da,Vec vu,Vec v,Vec f,void *w)
1650: {
1652: SETERRQ(1,"Must compile with bmake/PETSC_ARCH/packages flag PETSC_HAVE_ADIC for this routine");
1653: }
1655: #endif
1657: /*@
1658: DAComputeJacobian1 - Evaluates a local Jacobian function on each processor that
1659: share a DA
1661: Input Parameters:
1662: + da - the DA that defines the grid
1663: . vu - input vector (ghosted)
1664: . J - output matrix
1665: - w - any user data
1667: Notes: Does NOT do ghost updates on vu upon entry
1669: Level: advanced
1671: .seealso: DAFormFunction1()
1673: @*/
1674: int DAComputeJacobian1(DA da,Vec vu,Mat J,void *w)
1675: {
1676: int ierr;
1677: void *u;
1678: DALocalInfo info;
1681: DAGetLocalInfo(da,&info);
1682: DAVecGetArray(da,vu,&u);
1683: (*da->lj)(&info,u,J,w);
1684: DAVecRestoreArray(da,vu,&u);
1685: return(0);
1686: }
1689: /*
1690: DAComputeJacobian1WithAdifor - Evaluates a ADIFOR provided Jacobian local function on each processor that
1691: share a DA
1693: Input Parameters:
1694: + da - the DA that defines the grid
1695: . vu - input vector (ghosted)
1696: . J - output matrix
1697: - w - any user data
1699: Notes: Does NOT do ghost updates on vu upon entry
1701: .seealso: DAFormFunction1()
1703: */
1704: int DAComputeJacobian1WithAdifor(DA da,Vec vu,Mat J,void *w)
1705: {
1706: int i,ierr,Nc,*color,N;
1707: DALocalInfo info;
1708: PetscScalar *u,*g_u,*g_f,*f,*p_u;
1709: ISColoring iscoloring;
1710: void (*lf)(int *,DALocalInfo*,PetscScalar*,PetscScalar*,int*,PetscScalar*,PetscScalar*,int*,void*,int*) =
1711: (void (*)(int *,DALocalInfo*,PetscScalar*,PetscScalar*,int*,PetscScalar*,PetscScalar*,int*,void*,int*))*da->adifor_lf;
1714: DAGetColoring(da,IS_COLORING_GHOSTED,&iscoloring);
1715: Nc = iscoloring->n;
1716: DAGetLocalInfo(da,&info);
1717: N = info.gxm*info.gym*info.gzm*info.dof;
1719: /* get space for derivative objects. */
1720: ierr = PetscMalloc(Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar),&g_u);
1721: ierr = PetscMemzero(g_u,Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar));
1722: p_u = g_u;
1723: color = iscoloring->colors;
1724: for (i=0; i<N; i++) {
1725: p_u[*color++] = 1.0;
1726: p_u += Nc;
1727: }
1728: ISColoringDestroy(iscoloring);
1729: PetscMalloc(Nc*info.xm*info.ym*info.zm*info.dof*sizeof(PetscScalar),&g_f);
1730: PetscMalloc(info.xm*info.ym*info.zm*info.dof*sizeof(PetscScalar),&f);
1732: /* Seed the input array g_u with coloring information */
1733:
1734: VecGetArray(vu,&u);
1735: (lf)(&Nc,&info,u,g_u,&Nc,f,g_f,&Nc,w,&ierr);
1736: VecRestoreArray(vu,&u);
1738: /* stick the values into the matrix */
1739: /* PetscScalarView(Nc*info.xm*info.ym,g_f,0); */
1740: MatSetValuesAdifor(J,Nc,g_f);
1742: /* return space for derivative objects. */
1743: PetscFree(g_u);
1744: PetscFree(g_f);
1745: PetscFree(f);
1746: return(0);
1747: }
1749: /*@C
1750: DAMultiplyByJacobian1WithAD - Applies a Jacobian function supplied by ADIFOR or ADIC
1751: to a vector on each processor that shares a DA.
1753: Input Parameters:
1754: + da - the DA that defines the grid
1755: . vu - Jacobian is computed at this point (ghosted)
1756: . v - product is done on this vector (ghosted)
1757: . fu - output vector = J(vu)*v (not ghosted)
1758: - w - any user data
1760: Notes:
1761: This routine does NOT do ghost updates on vu and v upon entry.
1762:
1763: Automatically calls DAMultiplyByJacobian1WithAdifor() or DAMultiplyByJacobian1WithAdic()
1764: depending on whether DASetLocalAdicMFFunction() or DASetLocalAdiforMFFunction() was called.
1766: Level: advanced
1768: .seealso: DAFormFunction1(), DAMultiplyByJacobian1WithAdifor(), DAMultiplyByJacobian1WithAdic()
1770: @*/
1771: int DAMultiplyByJacobian1WithAD(DA da,Vec u,Vec v,Vec f,void *w)
1772: {
1773: int ierr;
1776: if (da->adicmf_lf) {
1777: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
1778: DAMultiplyByJacobian1WithAdic(da,u,v,f,w);
1779: #else
1780: SETERRQ(1,"Requires ADIC to be installed and cannot use complex numbers");
1781: #endif
1782: } else if (da->adiformf_lf) {
1783: DAMultiplyByJacobian1WithAdifor(da,u,v,f,w);
1784: } else {
1785: SETERRQ(1,"Must call DASetLocalAdiforMFFunction() or DASetLocalAdicMFFunction() before using");
1786: }
1787: return(0);
1788: }
1791: /*@C
1792: DAMultiplyByJacobian1WithAdifor - Applies a ADIFOR provided Jacobian function on each processor that
1793: share a DA to a vector
1795: Input Parameters:
1796: + da - the DA that defines the grid
1797: . vu - Jacobian is computed at this point (ghosted)
1798: . v - product is done on this vector (ghosted)
1799: . fu - output vector = J(vu)*v (not ghosted)
1800: - w - any user data
1802: Notes: Does NOT do ghost updates on vu and v upon entry
1804: Level: advanced
1806: .seealso: DAFormFunction1()
1808: @*/
1809: int DAMultiplyByJacobian1WithAdifor(DA da,Vec u,Vec v,Vec f,void *w)
1810: {
1811: int ierr;
1812: PetscScalar *au,*av,*af,*awork;
1813: Vec work;
1814: DALocalInfo info;
1815: void (*lf)(DALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,int*) =
1816: (void (*)(DALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,int*))*da->adiformf_lf;
1819: DAGetLocalInfo(da,&info);
1821: DAGetGlobalVector(da,&work);
1822: VecGetArray(u,&au);
1823: VecGetArray(v,&av);
1824: VecGetArray(f,&af);
1825: VecGetArray(work,&awork);
1826: (lf)(&info,au,av,awork,af,w,&ierr);
1827: VecRestoreArray(u,&au);
1828: VecRestoreArray(v,&av);
1829: VecRestoreArray(f,&af);
1830: VecRestoreArray(work,&awork);
1831: DARestoreGlobalVector(da,&work);
1833: return(0);
1834: }
1836: /*@C
1837: DASetInterpolationType - Sets the type of interpolation that will be
1838: returned by DAGetInterpolation()
1840: Collective on DA
1842: Input Parameter:
1843: + da - initial distributed array
1844: . ctype - DA_Q1 and DA_Q0 are currently the only supported forms
1846: Level: intermediate
1848: .keywords: distributed array, interpolation
1850: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DA, DAInterpolationType
1851: @*/
1852: int DASetInterpolationType(DA da,DAInterpolationType ctype)
1853: {
1856: da->interptype = ctype;
1857: return(0);
1858: }