Actual source code: triconvert.c
1: #ifdef PETSC_RCS_HEADER
2: static char vcid[] = "$Id: triconvert.c,v 1.9 2001/08/07 21:31:52 bsmith Exp $";
3: #endif
4: /*
5: Converts triangulated grid data file.node and file.ele generated
6: by the Triangle code of Shewchuk to a PETSc AOData base.
8: */
10: #include petscao.h
11: #include petscbt.h
13: int main(int argc,char **args)
14: {
15: int ierr,nvertex,dim,nstuff,nbound,dummy,bound,i,ncell,*cell,*cell_edge,*edge_vertex;
16: int nedge,j,k,*edge_cell,*cell_cell,nc,ncp;
17: int i0,i1;
18: int *shift0,*shift1;
19: char filebase[256],filename[256];
20: FILE *file;
21: AOData ao;
22: PetscReal *vertex,ddummy;
23: PetscBT vertex_boundary;
24: PetscTruth flag;
26: PetscInitialize(&argc,&args,0,0);
27: PetscMemzero(filename,256*sizeof(char));
29: PetscOptionsGetString(0,"-f",filebase,246,&flag);
30: if (!flag) {
31: SETERRQ(1,"Must provide filebase name with -f");
32: }
34: /*
35: Create empty database
36: */
37: AODataCreateBasic(PETSC_COMM_SELF,&ao);
39: /* -----------------------------------------------------------------------------------
40: Read in vertex information
41: */
42: PetscStrcpy(filename,filebase);
43: PetscStrcat(filename,".node");
44: file = fopen(filename,"r");
45: if (!file) {
46: SETERRQ1(1,"Unable to open node file: %s",filename);
47: }
48: fscanf(file,"%d %d %d %dn",&nvertex,&dim,&nstuff,&nbound);
49: if (dim != 2) {
50: SETERRQ(1,"Triangulation is not in two dimensions");
51: }
52: PetscMalloc(2*nvertex*sizeof(PetscReal),&vertex);
53: ierr = PetscBTCreate(nvertex,vertex_boundary);
55: if (nstuff == 1) {
56: for (i=0; i<nvertex; i++) {
57: fscanf(file,"%d %le %le %le %dn",&dummy,vertex+2*i,vertex+2*i+1,&ddummy,&bound);
58: if (bound) PetscBTSet(vertex_boundary,i);
59: }
60: } else if (nstuff == 0) {
61: for (i=0; i<nvertex; i++) {
62: fscanf(file,"%d %le %le %dn",&dummy,vertex+2*i,vertex+2*i+1,&bound);
63: if (bound) PetscBTSet(vertex_boundary,i);
64: }
65: } else SETERRQ(1,"No support yet for that number of vertex quantities");
66: fclose(file);
68: /*
69: Put vertex into database
70: */
71: AODataKeyAdd(ao,"vertex",nvertex,nvertex);
72: AODataSegmentAdd(ao,"vertex","values",2,nvertex,0,vertex,PETSC_DOUBLE);
73: AODataSegmentAdd(ao,"vertex","boundary",1,nvertex,0,vertex_boundary,PETSC_LOGICAL);
74: PetscBTDestroy(vertex_boundary);
76: /* -----------------------------------------------------------------------------------
77: Read in triangle information
78: */
79: PetscStrcpy(filename,filebase);
80: PetscStrcat(filename,".ele");
81: file = fopen(filename,"r");
82: if (!file) {
83: SETERRQ(1,"Unable to open element file");
84: }
85: fscanf(file,"%d %d %dn",&ncell,&nc,&nstuff);ncp = nc;
87: PetscMalloc(nc*ncell*sizeof(int),&cell);
88: if (nstuff == 0) {
89: if (nc == 3) {
90: for (i=0; i<ncell; i++) {
91: fscanf(file,"%d %d %d %dn",&dummy,cell+3*i,cell+3*i+1,cell+3*i+2);
92: }
93: } else if (nc == 6) {
94: for (i=0; i<ncell; i++) {
95: fscanf(file,"%d %d %d %d %d %d %dn",&dummy,cell+6*i,cell+6*i+2,cell+6*i+4,
96: cell+6*i+3,cell+6*i+5,cell+6*i+1);
97: }
98: }
99: } else {
100: SETERRQ(1,"No support yet for that number of element quantities");
101: }
102: fclose(file);
103: for (i=0; i<nc*ncell; i++) {
104: cell[i]--; /* shift because triangle starts at 1, not 0 */
105: }
107: /*
108: Put cell into database
109: */
110: AODataKeyAdd(ao,"cell",ncell,ncell);
111: AODataSegmentAdd(ao,"cell","vertex",nc,ncell,0,cell,PETSC_INT);
113: PetscMalloc(nc*ncell*sizeof(int),&cell_edge);
114: PetscMalloc(2*nc*ncell*sizeof(int),&edge_cell);
115: PetscMalloc(2*nc*ncell*sizeof(int),&edge_vertex);
116: PetscMalloc(3*ncell*sizeof(int),&cell_cell);
118: /*
119: Determine edges
120: */
121: PetscMalloc(nc*sizeof(int),&shift0);
122: PetscMalloc(nc*sizeof(int),&shift1);
123: for (i=0; i<nc; i++) {
124: shift0[i] = i;
125: shift1[i] = (i + 1) % nc;
126: }
129: nedge = 0;
130: for (i=0; i<ncell; i++) {
131: for (k=0; k<nc; k++) {
132: i0 = cell[nc*i+shift0[k]];
133: i1 = cell[nc*i+shift1[k]];
134: for (j=0; j<nedge; j++) {
135: if ((i0 == edge_vertex[2*j+1] && i1 == edge_vertex[2*j]) ||
136: (i1 == edge_vertex[2*j+1] && i0 == edge_vertex[2*j])) {
137: cell_edge[nc*i+k] = j;
138: edge_cell[2*j+1] = i;
139: goto found;
140: }
141: }
142: /*
143: Add a new edge to the list
144: */
145: edge_cell[2*nedge] = i;
146: edge_cell[2*nedge+1] = -1;
147: edge_vertex[2*nedge] = i0;
148: edge_vertex[2*nedge+1] = i1;
149: cell_edge[nc*i+k] = nedge;
150: nedge++;
151:
152: found:;
153: }
154: }
156: /*
157: Determine cell neighbors
158: */
159: PetscStrcpy(filename,filebase);
160: PetscStrcat(filename,".neigh");
161: file = fopen(filename,"r");
162: if (file) {
163: fscanf(file,"%d %dn",&ncell,&nc);
164: if (nc != 3) SETERRQ(PETSC_ERR_SUP,"Can only handle three neighbors");
165: for (i=0; i<ncell; i++) {
166: fscanf(file,"%d %d %d %dn",&dummy,cell_cell+3*i+1,cell_cell+3*i+2,cell_cell+3*i);
167: }
168: for (i=0; i<3*ncell; i++) {
169: cell_cell[i]--; /* shift because triangle starts at 1, not 0 */
170: }
171: fclose(file);
173: } else { /* no neighbor list file given, generate manually only works for nc == 3 */
174: if (nc != 3) SETERRQ(PETSC_ERR_SUP,"No neighbor file given and cannot determine neighbors");
176: for (i=0; i<ncell; i++) {
177: for (k=0; k<3; k++) {
178: i0 = cell_edge[3*i+k];
179: if (edge_cell[2*i0] != i) cell_cell[3*i+k] = edge_cell[2*i0];
180: else cell_cell[3*i+k] = edge_cell[2*i0+1];
181: }
182: }
183: }
185: AODataSegmentAdd(ao,"cell","cell",3,ncell,0,cell_cell,PETSC_INT);
186: AODataSegmentAdd(ao,"cell","edge",nc,ncell,0,cell_edge,PETSC_INT);
188: AODataKeyAdd(ao,"edge",nedge,nedge);
189: AODataSegmentAdd(ao,"edge","vertex",2,nedge,0,edge_vertex,PETSC_INT);
191: PetscFree(vertex);
192: PetscFree(cell);
193: PetscFree(shift0);
194: PetscFree(shift1);
196: /*
197: Add information about cell shape and element type to database
198: */
199: AODataKeyAdd(ao,"info",1,1);
200: AODataSegmentAdd(ao,"info","shape",10,1,0,(void*)"triangular",PETSC_CHAR);
201: if (ncp == 3) {
202: AODataSegmentAdd(ao,"info","element",6,1,0,(void*)"linear",PETSC_CHAR);
203: } else if (ncp == 6) {
204: AODataSegmentAdd(ao,"info","element",13,1,0,(void*)"quadratic",PETSC_CHAR);
205: }
207: /* AODataView(ao,0); */
209: { PetscViewer binary;
210: PetscViewerBinaryOpen(PETSC_COMM_WORLD,filebase,PETSC_BINARY_CREATE,&binary);
211: AODataView(ao,binary);
212: PetscViewerDestroy(binary);
213: }
214:
215: return 0;
216: }