Actual source code: ctable.c
petsc-3.6.4 2016-04-12
2: /* Contributed by - Mark Adams */
4: #include <petscsys.h>
5: #include <petscctable.h>
9: /*
10: PetscTableCreate Creates a PETSc look up table
12: Input Parameters:
13: + n - expected number of keys
14: - maxkey- largest possible key
16: Notes: keys are between 1 and maxkey inclusive
18: */
19: PetscErrorCode PetscTableCreate(const PetscInt n,PetscInt maxkey,PetscTable *rta)
20: {
21: PetscTable ta;
25: if (n < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"n < 0");
26: PetscNew(&ta);
27: ta->tablesize = (3*n)/2 + 17;
28: if (ta->tablesize < n) ta->tablesize = PETSC_MAX_INT/4; /* overflow */
29: PetscCalloc1(ta->tablesize,&ta->keytable);
30: PetscMalloc1(ta->tablesize,&ta->table);
31: ta->head = 0;
32: ta->count = 0;
33: ta->maxkey = maxkey;
34: *rta = ta;
35: return(0);
36: }
40: /* PetscTableCreate() ********************************************
41: *
42: * hash table for non-zero data and keys
43: *
44: */
45: PetscErrorCode PetscTableCreateCopy(const PetscTable intable,PetscTable *rta)
46: {
48: PetscInt i;
49: PetscTable ta;
52: PetscNew(&ta);
53: ta->tablesize = intable->tablesize;
54: PetscMalloc1(ta->tablesize,&ta->keytable);
55: PetscMalloc1(ta->tablesize,&ta->table);
56: for (i = 0; i < ta->tablesize; i++) {
57: ta->keytable[i] = intable->keytable[i];
58: ta->table[i] = intable->table[i];
59: #if defined(PETSC_USE_DEBUG)
60: if (ta->keytable[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_COR,"ta->keytable[i] < 0");
61: #endif
62: }
63: ta->head = 0;
64: ta->count = intable->count;
65: ta->maxkey = intable->maxkey;
66: *rta = ta;
67: return(0);
68: }
72: /* PetscTableDestroy() ********************************************
73: *
74: *
75: */
76: PetscErrorCode PetscTableDestroy(PetscTable *ta)
77: {
81: if (!*ta) return(0);
82: PetscFree((*ta)->keytable);
83: PetscFree((*ta)->table);
84: PetscFree(*ta);
85: return(0);
86: }
90: /* PetscTableGetCount() ********************************************
91: */
92: PetscErrorCode PetscTableGetCount(const PetscTable ta,PetscInt *count)
93: {
95: *count = ta->count;
96: return(0);
97: }
101: /* PetscTableIsEmpty() ********************************************
102: */
103: PetscErrorCode PetscTableIsEmpty(const PetscTable ta,PetscInt *flag)
104: {
106: *flag = !(ta->count);
107: return(0);
108: }
112: /*
113: PetscTableAddExpand - called PetscTableAdd() if more space needed
115: */
116: PetscErrorCode PetscTableAddExpand(PetscTable ta,PetscInt key,PetscInt data,InsertMode imode)
117: {
119: PetscInt ii = 0;
120: const PetscInt tsize = ta->tablesize,tcount = ta->count;
121: PetscInt *oldtab = ta->table,*oldkt = ta->keytable,newk,ndata;
124: if (ta->tablesize == PETSC_MAX_INT/4) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_COR,"ta->tablesize < 0");
125: ta->tablesize = 2*tsize;
126: if (ta->tablesize <= tsize) ta->tablesize = PETSC_MAX_INT/4;
128: PetscMalloc1(ta->tablesize,&ta->table);
129: PetscCalloc1(ta->tablesize,&ta->keytable);
131: ta->count = 0;
132: ta->head = 0;
134: PetscTableAdd(ta,key,data,INSERT_VALUES);
135: /* rehash */
136: for (ii = 0; ii < tsize; ii++) {
137: newk = oldkt[ii];
138: if (newk) {
139: ndata = oldtab[ii];
140: PetscTableAdd(ta,newk,ndata,imode);
141: }
142: }
143: if (ta->count != tcount + 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_COR,"corrupted ta->count");
145: PetscFree(oldtab);
146: PetscFree(oldkt);
147: return(0);
148: }
153: /* PetscTableRemoveAll() ********************************************
154: *
155: *
156: */
157: PetscErrorCode PetscTableRemoveAll(PetscTable ta)
158: {
162: ta->head = 0;
163: if (ta->count) {
164: ta->count = 0;
166: PetscMemzero(ta->keytable,ta->tablesize*sizeof(PetscInt));
167: }
168: return(0);
169: }
175: /* PetscTableGetHeadPosition() ********************************************
176: *
177: */
178: PetscErrorCode PetscTableGetHeadPosition(PetscTable ta,PetscTablePosition *ppos)
179: {
180: PetscInt i = 0;
183: *ppos = NULL;
184: if (!ta->count) return(0);
186: /* find first valid place */
187: do {
188: if (ta->keytable[i]) {
189: *ppos = (PetscTablePosition)&ta->table[i];
190: break;
191: }
192: } while (i++ < ta->tablesize);
193: if (!*ppos) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_COR,"No head");
194: return(0);
195: }
199: /* PetscTableGetNext() ********************************************
200: *
201: * - iteration - PetscTablePosition is always valid (points to a data)
202: *
203: */
204: PetscErrorCode PetscTableGetNext(PetscTable ta,PetscTablePosition *rPosition,PetscInt *pkey,PetscInt *data)
205: {
206: PetscInt idex;
207: PetscTablePosition pos;
210: pos = *rPosition;
211: if (!pos) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Null position");
212: *data = *pos;
213: if (!*data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Null data");
214: idex = pos - ta->table;
215: *pkey = ta->keytable[idex];
216: if (!*pkey) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Null key");
218: /* get next */
219: do {
220: pos++; idex++;
221: if (idex >= ta->tablesize) {
222: pos = 0; /* end of list */
223: break;
224: } else if (ta->keytable[idex]) {
225: pos = ta->table + idex;
226: break;
227: }
228: } while (idex < ta->tablesize);
229: *rPosition = pos;
230: return(0);
231: }
236: PetscErrorCode PetscTableAddCountExpand(PetscTable ta,PetscInt key)
237: {
239: PetscInt ii = 0,hash = PetscHash(ta,key);
240: const PetscInt tsize = ta->tablesize,tcount = ta->count;
241: PetscInt *oldtab = ta->table,*oldkt = ta->keytable,newk,ndata;
244: /* before making the table larger check if key is already in table */
245: while (ii++ < ta->tablesize) {
246: if (ta->keytable[hash] == key) return(0);
247: hash = (hash == (ta->tablesize-1)) ? 0 : hash+1;
248: }
250: /* alloc new (bigger) table */
251: if (ta->tablesize == PETSC_MAX_INT/4) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_COR,"ta->tablesize < 0");
252: ta->tablesize = 2*tsize;
253: if (ta->tablesize <= tsize) ta->tablesize = PETSC_MAX_INT/4;
255: PetscMalloc1(ta->tablesize,&ta->table);
256: PetscCalloc1(ta->tablesize,&ta->keytable);
258: ta->count = 0;
259: ta->head = 0;
261: /* Build a new copy of the data */
262: for (ii = 0; ii < tsize; ii++) {
263: newk = oldkt[ii];
264: if (newk) {
265: ndata = oldtab[ii];
266: PetscTableAdd(ta,newk,ndata,INSERT_VALUES);
267: }
268: }
269: PetscTableAddCount(ta,key);
270: if (ta->count != tcount + 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_COR,"corrupted ta->count");
272: PetscFree(oldtab);
273: PetscFree(oldkt);
274: return(0);
275: }