Actual source code: aodatabasic.c
1: /*$Id: aodatabasic.c,v 1.64 2001/08/07 03:04:35 balay Exp $*/
3: /*
4: The most basic AOData routines. These store the entire database on each processor.
5: These routines are very simple; note that we do not even use a private data structure
6: for AOData, and the private datastructure for AODataSegment is just used as a simple array.
8: These are made slightly complicated by having to be able to handle logical variables
9: stored in bit arrays. Thus,
10: - Before mallocing to hold a bit array, we shrunk the array length by a factor
11: of 8 using PetscBTLength()
12: - We use PetscBitMemcpy() to allow us to copy at the individual bit level;
13: for regular datatypes this just does a regular memcpy().
14: */
16: #include src/dm/ao/aoimpl.h
17: #include petscsys.h
18: #include petscbt.h
22: int AODataDestroy_Basic(AOData ao)
23: {
24: int ierr;
25: AODataKey *key = ao->keys,*nextkey;
26: AODataSegment *seg,*nextseg;
29: /* if memory was published with AMS then destroy it */
30: PetscObjectDepublish(ao);
32: while (key) {
33: PetscFree(key->name);
34: if (key->ltog) {
35: ISLocalToGlobalMappingDestroy(key->ltog);
36: }
37: seg = key->segments;
38: while (seg) {
39: PetscFree(seg->data);
40: PetscFree(seg->name);
41: nextseg = seg->next;
42: PetscFree(seg);
43: seg = nextseg;
44: }
45: PetscFree(key->rowners);
46: nextkey = key->next;
47: PetscFree(key);
48: key = nextkey;
49: }
50:
51: PetscLogObjectDestroy(ao);
52: PetscHeaderDestroy(ao);
53: return(0);
54: }
58: int AODataView_Basic_Binary(AOData ao,PetscViewer viewer)
59: {
60: int ierr,N,fd;
61: AODataSegment *segment;
62: AODataKey *key = ao->keys;
63: char paddedname[256];
67: PetscViewerBinaryGetDescriptor(viewer,&fd);
69: /* write out number of keys */
70: PetscBinaryWrite(fd,&ao->nkeys,1,PETSC_INT,0);
72: while (key) {
73: N = key->N;
74: /*
75: Write out name of key - use a fixed length for the name in the binary
76: file to make seeking easier
77: */
78: PetscMemzero(paddedname,256*sizeof(char));
79: PetscStrncpy(paddedname,key->name,255);
80: PetscBinaryWrite(fd,paddedname,256,PETSC_CHAR,0);
81: /* write out the number of indices */
82: PetscBinaryWrite(fd,&key->N,1,PETSC_INT,0);
83: /* write out number of segments */
84: PetscBinaryWrite(fd,&key->nsegments,1,PETSC_INT,0);
85:
86: /* loop over segments writing them out */
87: segment = key->segments;
88: while (segment) {
89: /*
90: Write out name of segment - use a fixed length for the name in the binary
91: file to make seeking easier
92: */
93: PetscMemzero(paddedname,256*sizeof(char));
94: PetscStrncpy(paddedname,segment->name,255);
95: PetscBinaryWrite(fd,paddedname,256,PETSC_CHAR,0);
96: PetscBinaryWrite(fd,&segment->bs,1,PETSC_INT,0);
97: PetscBinaryWrite(fd,&segment->datatype,1,PETSC_INT,0);
98: /* write out the data */
99: PetscBinaryWrite(fd,segment->data,N*segment->bs,segment->datatype,0);
100: segment = segment->next;
101: }
102: key = key->next;
103: }
105: return(0);
106: }
108: /*
109: All processors have the same data so processor 0 prints it
110: */
113: int AODataView_Basic_ASCII(AOData ao,PetscViewer viewer)
114: {
115: int ierr,j,k,l,rank,size,nkeys,nsegs,i,N,bs,zero = 0;
116: char **keynames,**segnames,*segvalue;
117: const char *stype,*dt;
118: AODataSegment *segment;
119: AODataKey *key = ao->keys;
120: PetscDataType dtype;
121: PetscViewerFormat format;
124: MPI_Comm_rank(ao->comm,&rank);
125: if (rank) return(0);
126: MPI_Comm_size(ao->comm,&size);
128: PetscViewerGetFormat(viewer,&format);
129: if (format == PETSC_VIEWER_ASCII_INFO) {
130: AODataGetInfo(ao,&nkeys,&keynames);
131: for (i=0; i<nkeys; i++) {
132: AODataKeyGetInfo(ao,keynames[i],&N,0,&nsegs,&segnames);
133: PetscViewerASCIIPrintf(viewer," %s: (%d)\n",keynames[i],N);
134: for (j=0; j<nsegs; j++) {
135: AODataSegmentGetInfo(ao,keynames[i],segnames[j],&bs,&dtype);
136: PetscDataTypeGetName(dtype,&stype);
137: if (dtype == PETSC_CHAR) {
138: AODataSegmentGet(ao,keynames[i],segnames[j],1,&zero,(void **)&segvalue);
139: PetscViewerASCIIPrintf(viewer," %s: (%d) %s -> %s\n",segnames[j],bs,stype,segvalue);
140: AODataSegmentRestore(ao,keynames[i],segnames[j],1,&zero,(void **)&segvalue);
141: } else {
142: PetscViewerASCIIPrintf(viewer," %s: (%d) %s\n",segnames[j],bs,stype);
143: }
144: }
145: }
146: PetscFree(keynames);
147: } else {
148: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
149: while (key) {
150: PetscViewerASCIIPrintf(viewer,"AOData Key: %s Length %d Ownership: ",key->name,key->N);
151: for (j=0; j<size+1; j++) {PetscViewerASCIIPrintf(viewer,"%d ",key->rowners[j]);}
152: PetscViewerASCIIPrintf(viewer,"\n");
154: segment = key->segments;
155: while (segment) {
156: PetscDataTypeGetName(segment->datatype,&dt);
157: PetscViewerASCIIPrintf(viewer," AOData Segment: %s Blocksize %d datatype %s\n",segment->name,segment->bs,dt);
158: if (segment->datatype == PETSC_INT) {
159: int *mdata = (int*)segment->data;
160: for (k=0; k<key->N; k++) {
161: PetscViewerASCIIPrintf(viewer," %d: ",k);
162: for (l=0; l<segment->bs; l++) {
163: PetscViewerASCIIPrintf(viewer," %d ",mdata[k*segment->bs + l]);
164: }
165: PetscViewerASCIIPrintf(viewer,"\n");
166: }
167: } else if (segment->datatype == PETSC_DOUBLE) {
168: PetscReal *mdata = (PetscReal*)segment->data;
169: for (k=0; k<key->N; k++) {
170: PetscViewerASCIIPrintf(viewer," %d: ",k);
171: for (l=0; l<segment->bs; l++) {
172: PetscViewerASCIIPrintf(viewer," %18.16e ",mdata[k*segment->bs + l]);
173: }
174: PetscViewerASCIIPrintf(viewer,"\n");
175: }
176: } else if (segment->datatype == PETSC_SCALAR) {
177: PetscScalar *mdata = (PetscScalar*)segment->data;
178: for (k=0; k<key->N; k++) {
179: PetscViewerASCIIPrintf(viewer," %d: ",k);
180: for (l=0; l<segment->bs; l++) {
181: #if !defined(PETSC_USE_COMPLEX)
182: PetscViewerASCIIPrintf(viewer," %18.16e ",mdata[k*segment->bs + l]);
183: #else
184: PetscScalar x = mdata[k*segment->bs + l];
185: if (PetscImaginaryPart(x) > 0.0) {
186: PetscViewerASCIIPrintf(viewer," %18.16e + %18.16e i \n",PetscRealPart(x),PetscImaginaryPart(x));
187: } else if (PetscImaginaryPart(x) < 0.0) {
188: PetscViewerASCIIPrintf(viewer," %18.16e - %18.16e i \n",PetscRealPart(x),-PetscImaginaryPart(x));
189: } else {
190: PetscViewerASCIIPrintf(viewer," %18.16e \n",PetscRealPart(x));
191: }
192: #endif
193: }
194: }
195: PetscViewerASCIIPrintf(viewer,"\n");
196: } else if (segment->datatype == PETSC_LOGICAL) {
197: PetscBT mdata = (PetscBT) segment->data;
198: for (k=0; k<key->N; k++) {
199: PetscViewerASCIIPrintf(viewer," %d: ",k);
200: for (l=0; l<segment->bs; l++) {
201: PetscViewerASCIIPrintf(viewer," %d ",(int)PetscBTLookup(mdata,k*segment->bs + l));
202: }
203: PetscViewerASCIIPrintf(viewer,"\n");
204: }
205: } else if (segment->datatype == PETSC_CHAR) {
206: char * mdata = (char*)segment->data;
207: for (k=0; k<key->N; k++) {
208: PetscViewerASCIIPrintf(viewer," %s ",mdata + k*segment->bs);
209: }
210: PetscViewerASCIIPrintf(viewer,"\n");
211: } else {
212: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown PETSc data format");
213: }
214: segment = segment->next;
215: }
216: key = key->next;
217: }
218: }
219: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
220: return(0);
221: }
225: int AODataView_Basic(AOData ao,PetscViewer viewer)
226: {
227: int rank,ierr;
228: PetscTruth isascii,isbinary;
231: MPI_Comm_rank(ao->comm,&rank);
232: if (rank) return(0);
234: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
235: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
236: if (isascii) {
237: AODataView_Basic_ASCII(ao,viewer);
238: } else if (isbinary) {
239: AODataView_Basic_Binary(ao,viewer);
240: } else {
241: SETERRQ1(1,"Viewer type %s not supported for AOData basic",((PetscObject)viewer)->type_name);
242: }
244: return(0);
245: }
249: int AODataKeyRemove_Basic(AOData aodata,const char name[])
250: {
251: AODataSegment *segment,*iseg;
252: AODataKey *key,*ikey;
253: int ierr;
254: PetscTruth flag;
257: AODataKeyFind_Private(aodata,name,&flag,&key);
258: if (flag != 1) return(0);
259: aodata->nkeys--;
261: segment = key->segments;
262: while (segment) {
263: iseg = segment->next;
264: PetscFree(segment->name);
265: PetscFree(segment->data);
266: PetscFree(segment);
267: segment = iseg;
268: }
269: ikey = aodata->keys;
270: if (key == ikey) {
271: aodata->keys = key->next;
272: goto finishup1;
273: }
274: while (ikey->next != key) {
275: ikey = ikey->next;
276: }
277: ikey->next = key->next;
279: finishup1:
281: PetscFree(key->name);
282: PetscFree(key->rowners);
283: PetscFree(key);
285: return(0);
286: }
290: int AODataSegmentRemove_Basic(AOData aodata,const char name[],const char segname[])
291: {
292: AODataSegment *segment,*iseg;
293: AODataKey *key;
294: int ierr;
295: PetscTruth flag;
298: AODataSegmentFind_Private(aodata,name,segname,&flag,&key,&iseg);
299: if (!flag) return(0);
300: key->nsegments--;
302: segment = key->segments;
303: if (segment == iseg) {
304: key->segments = segment->next;
305: goto finishup2;
306: }
307: while (segment->next != iseg) {
308: segment = segment->next;
309: }
310: segment->next = iseg->next;
311: segment = iseg;
312:
313: finishup2:
315: PetscFree(segment->name);
316: PetscFree(segment->data);
317: PetscFree(segment);
318: return(0);
319: }
324: int AODataSegmentAdd_Basic(AOData aodata,const char name[],const char segname[],int bs,int n,int *keys,void *data,PetscDataType dtype)
325: {
326: AODataSegment *segment,*iseg;
327: AODataKey *key;
328: int N,size,ierr,*lens,i,*disp,*akeys,datasize,*fkeys,Nb,j;
329: MPI_Datatype mtype;
330: char *adata;
331: MPI_Comm comm = aodata->comm;
332: PetscTruth flag;
335: AODataKeyFind_Private(aodata,name,&flag,&key);
336: if (!flag) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Key %s doesn't exist",name);
337: AODataSegmentFind_Private(aodata,name,segname,&flag,&key,&iseg);
338: if (flag) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Segment %s in key %s already exists",name,segname);
340: PetscNew(AODataSegment,&segment);
341: if (iseg) {
342: iseg->next = segment;
343: } else {
344: key->segments = segment;
345: }
346: segment->next = 0;
347: segment->bs = bs;
348: segment->datatype = dtype;
350: PetscDataTypeGetSize(dtype,&datasize);
352: /*
353: If keys not given, assume each processor provides entire data
354: */
355: if (!keys && n == key->N) {
356: char *fdata1;
357: if (dtype == PETSC_LOGICAL) Nb = PetscBTLength(key->N); else Nb = key->N;
358: PetscMalloc((Nb*bs+1)*datasize,&fdata1);
359: PetscBitMemcpy(fdata1,0,data,0,key->N*bs,dtype);
360: segment->data = (void*)fdata1;
361: } else if (!keys) {
362: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Keys not given, but not all data given on each processor");
363: } else {
364: /* transmit all lengths to all processors */
365: MPI_Comm_size(comm,&size);
366: PetscMalloc(2*size*sizeof(int),&lens);
367: disp = lens + size;
368: MPI_Allgather(&n,1,MPI_INT,lens,1,MPI_INT,comm);
369: N = 0;
370: for (i=0; i<size; i++) {
371: disp[i] = N;
372: N += lens[i];
373: }
374: if (N != key->N) {
375: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Did not provide correct number of keys for keyname");
376: }
378: /*
379: Allocate space for all keys and all data
380: */
381: PetscMalloc((N+1)*sizeof(int),&akeys);
382: PetscMalloc((N*bs+1)*datasize,&adata);
384: MPI_Allgatherv(keys,n,MPI_INT,akeys,lens,disp,MPI_INT,comm);
385: for (i=0; i<size; i++) {
386: disp[i] *= bs;
387: lens[i] *= bs;
388: }
389:
390: if (dtype != PETSC_LOGICAL) {
391: char *fdata2;
393: PetscDataTypeToMPIDataType(dtype,&mtype);
394: MPI_Allgatherv(data,n*bs,mtype,adata,lens,disp,mtype,comm);
395: PetscFree(lens);
397: /*
398: Now we have all the keys and data we need to put it in order
399: */
400: PetscMalloc((key->N+1)*sizeof(int),&fkeys);
401: PetscMemzero(fkeys,(key->N+1)*sizeof(int));
402: PetscMalloc((key->N*bs+1)*datasize,&fdata2);
404: for (i=0; i<N; i++) {
405: if (fkeys[akeys[i]] != 0) {
406: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Duplicate key");
407: }
408: if (fkeys[akeys[i]] >= N) {
409: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Key out of range");
410: }
411: fkeys[akeys[i]] = 1;
412: PetscBitMemcpy(fdata2,akeys[i]*bs,adata,i*bs,bs,dtype);
413: }
414: for (i=0; i<N; i++) {
415: if (!fkeys[i]) {
416: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Missing key");
417: }
418: }
419: segment->data = (void*)fdata2;
420: } else {
421: /*
422: For logical input the length is given by the user in bits; we need to
423: convert to bytes to send with MPI
424: */
425: PetscBT fdata3;
426: PetscBT mvalues = (PetscBT) data;
427: char *values;
428: PetscMalloc((n+1)*bs*sizeof(char),&values);
429: for (i=0; i<n; i++) {
430: for (j=0; j<bs; j++) {
431: if (PetscBTLookup(mvalues,i*bs+j)) values[i*bs+j] = 1; else values[i*bs+j] = 0;
432: }
433: }
435: MPI_Allgatherv(values,n*bs,MPI_BYTE,adata,lens,disp,MPI_BYTE,comm);
436: PetscFree(lens);
437: PetscFree(values);
439: /*
440: Now we have all the keys and data we need to put it in order
441: */
442: PetscMalloc((key->N+1)*sizeof(int),&fkeys);
443: PetscMemzero(fkeys,(key->N+1)*sizeof(int));
444: PetscBTCreate(N*bs,fdata3);
446: for (i=0; i<N; i++) {
447: if (fkeys[akeys[i]] != 0) {
448: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Duplicate key");
449: }
450: if (fkeys[akeys[i]] >= N) {
451: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Key out of range");
452: }
453: fkeys[akeys[i]] = 1;
454: for (j=0; j<bs; j++) {
455: if (adata[i*bs+j]) { PetscBTSet(fdata3,i*bs+j); }
456: }
457: }
458: for (i=0; i<N; i++) {
459: if (!fkeys[i]) {
460: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Missing key");
461: }
462: }
463: segment->data = (void*)fdata3;
464: }
465: PetscFree(akeys);
466: PetscFree(adata);
467: PetscFree(fkeys);
468: }
470: key->nsegments++;
472: PetscStrallocpy(segname,&segment->name);
473: return(0);
474: }
478: int AODataSegmentGetExtrema_Basic(AOData ao,const char name[],const char segname[],void *xmax,void *xmin)
479: {
480: AODataSegment *segment;
481: AODataKey *key;
482: int ierr,i,bs,n,j;
483: PetscTruth flag;
486: /* find the correct segment */
487: AODataSegmentFind_Private(ao,name,segname,&flag,&key,&segment);
488: if (!flag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Cannot locate segment");
490: n = key->N;
491: bs = segment->bs;
493: if (segment->datatype == PETSC_INT) {
494: int *vmax = (int*)xmax,*vmin = (int*)xmin,*values = (int*)segment->data;
495: for (j=0; j<bs; j++) {
496: vmax[j] = vmin[j] = values[j];
497: }
498: for (i=1; i<n; i++) {
499: for (j=0; j<bs; j++) {
500: vmax[j] = PetscMax(vmax[j],values[bs*i+j]);
501: vmin[j] = PetscMin(vmin[j],values[bs*i+j]);
502: }
503: }
504: } else if (segment->datatype == PETSC_DOUBLE) {
505: PetscReal *vmax = (PetscReal*)xmax,*vmin = (PetscReal*)xmin,*values = (PetscReal*)segment->data;
506: for (j=0; j<bs; j++) {
507: vmax[j] = vmin[j] = values[j];
508: }
509: for (i=1; i<n; i++) {
510: for (j=0; j<bs; j++) {
511: vmax[j] = PetscMax(vmax[j],values[bs*i+j]);
512: vmin[j] = PetscMin(vmin[j],values[bs*i+j]);
513: }
514: }
515: } else SETERRQ(PETSC_ERR_SUP,"Cannot find extrema for this data type");
517: return(0);
518: }
522: int AODataSegmentGet_Basic(AOData ao,const char name[],const char segname[],int n,int *keys,void **data)
523: {
524: AODataSegment *segment;
525: AODataKey *key;
526: int ierr,dsize,i,bs,nb;
527: char *idata,*odata;
528: PetscTruth flag;
529:
531: /* find the correct segment */
532: AODataSegmentFind_Private(ao,name,segname,&flag,&key,&segment);
533: if (!flag) SETERRQ2(PETSC_ERR_ARG_WRONG,"Cannot locate segment %s in key %s",segname,name);
535: PetscDataTypeGetSize(segment->datatype,&dsize);
536: bs = segment->bs;
537: if (segment->datatype == PETSC_LOGICAL) nb = PetscBTLength(n); else nb = n;
538: PetscMalloc((nb+1)*bs*dsize,&odata);
539: idata = (char*)segment->data;
540: for (i=0; i<n; i++) {
541: PetscBitMemcpy(odata,i*bs,idata,keys[i]*bs,bs,segment->datatype);
542: }
543: *data = (void*)odata;
544: return(0);
545: }
549: int AODataSegmentRestore_Basic(AOData aodata,const char name[],const char segname[],int n,int *keys,void **data)
550: {
554: PetscFree(*data);
555: return(0);
556: }
560: int AODataSegmentGetLocal_Basic(AOData ao,const char name[],const char segname[],int n,int *keys,void **data)
561: {
562: int ierr,*globals,*locals,bs;
563: PetscDataType dtype;
564: AODataKey *key;
565: PetscTruth flag;
568: AODataKeyFind_Private(ao,segname,&flag,&key);
569: if (!flag) SETERRQ(PETSC_ERR_ARG_WRONG,"Segment does not have corresponding key");
570: if (!key->ltog) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No local to global mapping set for key");
571: AODataSegmentGetInfo(ao,name,segname,&bs,&dtype);
572: if (dtype != PETSC_INT) SETERRQ(PETSC_ERR_ARG_WRONG,"Datatype of segment must be PETSC_INT");
574: /* get the values in global indexing */
575: AODataSegmentGet_Basic(ao,name,segname,n,keys,(void **)&globals);
576:
577: /* allocate space to store them in local indexing */
578: PetscMalloc((n+1)*bs*sizeof(int),&locals);
580: ISGlobalToLocalMappingApply(key->ltog,IS_GTOLM_MASK,n*bs,globals,PETSC_NULL,locals);
582: AODataSegmentRestore_Basic(ao,name,segname,n,keys,(void **)&globals);
584: *data = (void*)locals;
585: return(0);
586: }
590: int AODataSegmentRestoreLocal_Basic(AOData aodata,const char name[],const char segname[],int n,int *keys,void **data)
591: {
595: PetscFree(*data);
596: return(0);
597: }
599: EXTERN int AOBasicGetIndices_Private(AO,int **,int **);
603: int AODataKeyRemap_Basic(AOData aodata,const char keyname[],AO ao)
604: {
605: int ierr,*inew,k,*ii,nk,dsize,bs,nkb;
606: char *data,*tmpdata;
607: AODataKey *key;
608: AODataSegment *seg;
609: PetscTruth flag,match;
613: /* remap all the values in the segments that match the key */
614: key = aodata->keys;
615: while (key) {
616: seg = key->segments;
617: while (seg) {
618: PetscStrcmp(seg->name,keyname,&match);
619: if (!match) {
620: seg = seg->next;
621: continue;
622: }
623: if (seg->datatype != PETSC_INT) {
624: SETERRQ(PETSC_ERR_ARG_WRONG,"Segment name same as key but not integer type");
625: }
626: nk = seg->bs*key->N;
627: ii = (int*)seg->data;
628: AOPetscToApplication(ao,nk,ii);
629: seg = seg->next;
630: }
631: key = key->next;
632: }
633:
634: AOBasicGetIndices_Private(ao,&inew,PETSC_NULL);
635: /* reorder in the arrays all the values for the key */
636: AODataKeyFind_Private(aodata,keyname,&flag,&key);
637: if (!flag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Could not find key");
638: nk = key->N;
639: seg = key->segments;
640: while (seg) {
641: PetscDataTypeGetSize(seg->datatype,&dsize);
642: bs = seg->bs;
643: data = (char*)seg->data;
644: if (seg->datatype == PETSC_LOGICAL) nkb = PetscBTLength(nk*bs); else nkb = nk*bs;
645: PetscMalloc((nkb+1)*dsize,&tmpdata);
647: for (k=0; k<nk; k++) {
648: PetscBitMemcpy(tmpdata,inew[k]*bs,data,k*bs,bs,seg->datatype);
649: }
650: PetscMemcpy(data,tmpdata,nkb*dsize);
651: PetscFree(tmpdata);
652: seg = seg->next;
653: }
655: return(0);
656: }
660: int AODataKeyGetAdjacency_Basic(AOData aodata,const char keyname[],Mat *adj)
661: {
662: int ierr,cnt,i,j,*jj,*ii,nlocal,n,*nb,bs,ls;
663: AODataKey *key;
664: AODataSegment *seg;
665: PetscTruth flag;
668: AODataSegmentFind_Private(aodata,keyname,keyname,&flag,&key,&seg);
669: if (!flag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Cannot locate key with neighbor segment");
671: /*
672: Get the beginning of the neighbor list for this processor
673: */
674: bs = seg->bs;
675: nb = (int*)seg->data;
676: nb += bs*key->rstart;
677: nlocal = key->rend - key->rstart;
678: n = bs*key->N;
680: /*
681: Assemble the adjacency graph: first we determine total number of entries
682: */
683: cnt = 0;
684: for (i=0; i<bs*nlocal; i++) {
685: if (nb[i] >= 0) cnt++;
686: }
687: PetscMalloc((nlocal + 1)*sizeof(int),&ii);
688: PetscMalloc((cnt+1)*sizeof(int),&jj);
689: ii[0] = 0;
690: cnt = 0;
691: for (i=0; i<nlocal; i++) {
692: ls = 0;
693: for (j=0; j<bs; j++) {
694: if (nb[bs*i+j] >= 0) {
695: jj[cnt++] = nb[bs*i+j];
696: ls++;
697: }
698: }
699: /* now sort the column indices for this row */
700: PetscSortInt(ls,jj+cnt-ls);
701: ii[i+1] = cnt;
702: }
704: MatCreate(aodata->comm,nlocal,n,PETSC_DETERMINE,n,adj);
705: MatSetType(*adj,MATMPIADJ);
706: MatMPIAdjSetPreallocation(*adj,ii,jj,PETSC_NULL);
707: return(0);
708: }
712: int AODataSegmentPartition_Basic(AOData aodata,const char keyname[],const char segname[])
713: {
714: int ierr,size,bs,i,j,*idx,nc,*isc;
715: AO ao;
716: AODataKey *key,*keyseg;
717: AODataSegment *segment;
718: PetscTruth flag;
722: AODataKeyFind_Private(aodata,segname,&flag,&keyseg);
723: if (!flag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Cannot locate segment as a key");
724: PetscMalloc(keyseg->N*sizeof(int),&isc);
725: PetscMemzero(isc,keyseg->N*sizeof(int));
727: AODataSegmentFind_Private(aodata,keyname,segname,&flag,&key,&segment);
728: if (flag != 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Cannot locate segment");
729: MPI_Comm_size(aodata->comm,&size);
731: bs = segment->bs;
733: idx = (int*)segment->data;
734: nc = 0;
735: for (i=0; i<size; i++) {
736: for (j=bs*key->rowners[i]; j<bs*key->rowners[i+1]; j++) {
737: if (!isc[idx[j]]) {
738: isc[idx[j]] = ++nc;
739: }
740: }
741: }
742: for (i=0; i<keyseg->N; i++) {
743: isc[i]--;
744: }
746: AOCreateBasic(aodata->comm,keyseg->nlocal,isc+keyseg->rstart,PETSC_NULL,&ao);
747: PetscFree(isc);
749: AODataKeyRemap(aodata,segname,ao);
750: AODestroy(ao);
751: return(0);
752: }
756: int AODataKeyGetActive_Basic(AOData aodata,const char name[],const char segname[],int n,int *keys,int wl,IS *is)
757: {
758: int ierr,i,cnt,*fnd,bs;
759: AODataKey *key;
760: AODataSegment *segment;
761: PetscBT bt;
762: PetscTruth flag;
765: AODataSegmentFind_Private(aodata,name,segname,&flag,&key,&segment);
766: if (!flag) SETERRQ(1,"Cannot locate segment");
768: bt = (PetscBT) segment->data;
769: bs = segment->bs;
771: if (wl >= bs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Bit field (wl) argument too large");
773: /* count the active ones */
774: cnt = 0;
775: for (i=0; i<n; i++) {
776: if (PetscBTLookup(bt,keys[i]*bs+wl)) {
777: cnt++;
778: }
779: }
781: PetscMalloc((cnt+1)*sizeof(int),&fnd);
782: cnt = 0;
783: for (i=0; i<n; i++) {
784: if (PetscBTLookup(bt,keys[i]*bs+wl)) {
785: fnd[cnt++] = keys[i];
786: }
787: }
788:
789: ISCreateGeneral(aodata->comm,cnt,fnd,is);
790: PetscFree(fnd);
791: return(0);
792: }
796: int AODataKeyGetActiveLocal_Basic(AOData aodata,const char name[],const char segname[],int n,int *keys,int wl,IS *is)
797: {
798: int ierr,i,cnt,*fnd,bs,*locals;
799: AODataKey *key;
800: AODataSegment *segment;
801: PetscBT bt;
802: PetscTruth flag;
805: AODataSegmentFind_Private(aodata,name,segname,&flag,&key,&segment);
806: if (!flag) SETERRQ(1,"Cannot locate segment");
808: bt = (PetscBT) segment->data;
809: bs = segment->bs;
811: if (wl >= bs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Bit field (wl) argument too large");
813: /* count the active ones */
814: cnt = 0;
815: for (i=0; i<n; i++) {
816: if (PetscBTLookup(bt,keys[i]*bs+wl)) {
817: cnt++;
818: }
819: }
821: PetscMalloc((cnt+1)*sizeof(int),&fnd);
822: cnt = 0;
823: for (i=0; i<n; i++) {
824: if (PetscBTLookup(bt,keys[i]*bs+wl)) {
825: fnd[cnt++] = keys[i];
826: }
827: }
828:
829: PetscMalloc((n+1)*sizeof(int),&locals);
830: ISGlobalToLocalMappingApply(key->ltog,IS_GTOLM_MASK,cnt,fnd,PETSC_NULL,locals);
831: PetscFree(fnd);
832: ISCreateGeneral(aodata->comm,cnt,locals,is);
833: PetscFree(locals);
834: return(0);
835: }
837: EXTERN int AODataSegmentGetReduced_Basic(AOData,const char[],const char[],int,int*,IS *);
838: EXTERN int AODataPublish_Petsc(PetscObject);
840: static struct _AODataOps myops = {AODataSegmentAdd_Basic,
841: AODataSegmentGet_Basic,
842: AODataSegmentRestore_Basic,
843: AODataSegmentGetLocal_Basic,
844: AODataSegmentRestoreLocal_Basic,
845: AODataSegmentGetReduced_Basic,
846: AODataSegmentGetExtrema_Basic,
847: AODataKeyRemap_Basic,
848: AODataKeyGetAdjacency_Basic,
849: AODataKeyGetActive_Basic,
850: AODataKeyGetActiveLocal_Basic,
851: AODataSegmentPartition_Basic,
852: AODataKeyRemove_Basic,
853: AODataSegmentRemove_Basic,
854: AODataDestroy_Basic,
855: AODataView_Basic};
859: /*@C
860: AODataCreateBasic - Creates an AO datastructure.
862: Collective on MPI_Comm
864: Input Parameters:
865: . comm - MPI communicator that is to share AO
867: Output Parameter:
868: . aoout - the new database
870: Options Database Keys:
871: + -ao_data_view - Prints entire database at the conclusion of AODataSegmentAdd()
872: - -ao_data_view_info - Prints info about database at the conclusion of AODataSegmentAdd()
874: Level: intermediate
876: .keywords: AOData, create
878: .seealso: AODataSegmentAdd(), AODataDestroy()
879: @*/
880: int AODataCreateBasic(MPI_Comm comm,AOData *aoout)
881: {
882: AOData ao;
883: int ierr;
887: *aoout = 0;
888: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
889: DMInitializePackage(PETSC_NULL);
890: #endif
892: PetscHeaderCreate(ao,_p_AOData,struct _AODataOps,AODATA_COOKIE,AODATA_BASIC,"AOData",comm,AODataDestroy,AODataView);
893: PetscLogObjectCreate(ao);
894: PetscLogObjectMemory(ao,sizeof(struct _p_AOData));
896: PetscMemcpy(ao->ops,&myops,sizeof(myops));
897: ao->bops->publish = AODataPublish_Petsc;
899: ao->nkeys = 0;
900: ao->keys = 0;
901: ao->datacomplete = 0;
903: PetscPublishAll(ao);
904: *aoout = ao;
905: return(0);
906: }
910: /*@C
911: AODataLoadBasic - Loads an AO database from a file.
913: Collective on PetscViewer
915: Input Parameters:
916: . viewer - the binary file containing the data
918: Output Parameter:
919: . aoout - the new database
921: Options Database Keys:
922: + -ao_data_view - Prints entire database at the conclusion of AODataLoadBasic()
923: - -ao_data_view_info - Prints info about database at the conclusion of AODataLoadBasic()
925: Level: intermediate
927: .keywords: AOData, create, load, basic
929: .seealso: AODataSegmentAdd(), AODataDestroy(), AODataCreateBasic(), AODataView()
930: @*/
931: int AODataLoadBasic(PetscViewer viewer,AOData *aoout)
932: {
933: AOData ao;
934: int fd,nkeys,i,ierr,dsize,j,size,rank,Nb;
935: char paddedname[256];
936: AODataSegment *seg = 0;
937: AODataKey *key = 0;
938: MPI_Comm comm;
939: PetscTruth isbinary,flg1;
942: *aoout = 0;
943: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
944: DMInitializePackage(PETSC_NULL);
945: #endif
946: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
947: if (!isbinary) {
948: SETERRQ(PETSC_ERR_ARG_WRONG,"Viewer must be obtained from PetscViewerBinaryOpen()");
949: }
951: PetscObjectGetComm((PetscObject)viewer,&comm);
952: MPI_Comm_size(comm,&size);
953: MPI_Comm_rank(comm,&rank);
955: PetscViewerBinaryGetDescriptor(viewer,&fd);
957: /* read in number of segments */
958: PetscBinaryRead(fd,&nkeys,1,PETSC_INT);
960: PetscHeaderCreate(ao,_p_AOData,struct _AODataOps,AODATA_COOKIE,AODATA_BASIC,"AOData",comm,AODataDestroy,AODataView);
961: PetscLogObjectCreate(ao);
962: PetscLogObjectMemory(ao,sizeof(struct _p_AOData) + nkeys*sizeof(void *));
964: PetscMemcpy(ao->ops,&myops,sizeof(myops));
965: ao->bops->publish = AODataPublish_Petsc;
967: ao->nkeys = nkeys;
969: for (i=0; i<nkeys; i++) {
970: if (!i) {
971: PetscNew(AODataKey,&key);
972: ao->keys = key;
973: } else {
974: PetscNew(AODataKey,&key->next);
975: key = key->next;
976: }
977: key->ltog = 0;
978: key->next = 0;
980: /* read in key name */
981: PetscBinaryRead(fd,paddedname,256,PETSC_CHAR);
982: PetscStrallocpy(paddedname,&key->name);
983: PetscBinaryRead(fd,&key->N,1,PETSC_INT);
985: /* determine Nlocal and rowners for key */
986: key->nlocal = key->N/size + ((key->N % size) > rank);
987: PetscMalloc((size+1)*sizeof(int),&key->rowners);
988: MPI_Allgather(&key->nlocal,1,MPI_INT,key->rowners+1,1,MPI_INT,comm);
989: key->rowners[0] = 0;
990: for (j=2; j<=size; j++) {
991: key->rowners[j] += key->rowners[j-1];
992: }
993: key->rstart = key->rowners[rank];
994: key->rend = key->rowners[rank+1];
996: /* loop key's segments, reading them in */
997: PetscBinaryRead(fd,&key->nsegments,1,PETSC_INT);
999: for (j=0; j<key->nsegments; j++) {
1000: if (!j) {
1001: PetscNew(AODataSegment,&seg);
1002: key->segments = seg;
1003: } else {
1004: PetscNew(AODataSegment,&seg->next);
1005: seg = seg->next;
1006: }
1008: /* read in segment name */
1009: PetscBinaryRead(fd,paddedname,256,PETSC_CHAR);
1010: PetscStrallocpy(paddedname,&seg->name);
1012: /* read in segment blocksize and datatype */
1013: PetscBinaryRead(fd,&seg->bs,1,PETSC_INT);
1014: PetscBinaryRead(fd,&seg->datatype,1,PETSC_INT);
1016: /* allocate the space for the data */
1017: PetscDataTypeGetSize(seg->datatype,&dsize);
1018: if (seg->datatype == PETSC_LOGICAL) Nb = PetscBTLength(key->N*seg->bs); else Nb = key->N*seg->bs;
1019: PetscMalloc(Nb*dsize,&seg->data);
1020: /* read in the data */
1021: PetscBinaryRead(fd,seg->data,key->N*seg->bs,seg->datatype);
1022: seg->next = 0;
1023: }
1024: }
1025: *aoout = ao;
1027: PetscOptionsHasName(PETSC_NULL,"-ao_data_view",&flg1);
1028: if (flg1) {
1029: AODataView(ao,PETSC_VIEWER_STDOUT_(comm));
1030: }
1031: PetscOptionsHasName(PETSC_NULL,"-ao_data_view_info",&flg1);
1032: if (flg1) {
1033: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(comm),PETSC_VIEWER_ASCII_INFO);
1034: AODataView(ao,PETSC_VIEWER_STDOUT_(comm));
1035: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(comm));
1036: }
1037: PetscPublishAll(ao);
1038: return(0);
1039: }