xref: /petsc/src/dm/impls/plex/plexexodusii.c (revision 552f735842aa6127d09b62f740a903cdd0631614)
1*552f7358SJed Brown #define PETSCDM_DLL
2*552f7358SJed Brown #include <petsc-private/pleximpl.h>    /*I   "petscdmplex.h"   I*/
3*552f7358SJed Brown 
4*552f7358SJed Brown #ifdef PETSC_HAVE_EXODUSII
5*552f7358SJed Brown #include<netcdf.h>
6*552f7358SJed Brown #include<exodusII.h>
7*552f7358SJed Brown #endif
8*552f7358SJed Brown 
9*552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexInterpolate_2D(DM, DM *);
10*552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexInterpolate_3D(DM, DM *);
11*552f7358SJed Brown 
12*552f7358SJed Brown #undef __FUNCT__
13*552f7358SJed Brown #define __FUNCT__ "DMPlexCreateExodus"
14*552f7358SJed Brown /*@
15*552f7358SJed Brown   DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file.
16*552f7358SJed Brown 
17*552f7358SJed Brown   Collective on comm
18*552f7358SJed Brown 
19*552f7358SJed Brown   Input Parameters:
20*552f7358SJed Brown + comm  - The MPI communicator
21*552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open
22*552f7358SJed Brown - interpolate - Create faces and edges in the mesh
23*552f7358SJed Brown 
24*552f7358SJed Brown   Output Parameter:
25*552f7358SJed Brown . dm  - The DM object representing the mesh
26*552f7358SJed Brown 
27*552f7358SJed Brown   Level: beginner
28*552f7358SJed Brown 
29*552f7358SJed Brown .keywords: mesh,ExodusII
30*552f7358SJed Brown .seealso: MeshCreate(), MeshCreateExodus()
31*552f7358SJed Brown @*/
32*552f7358SJed Brown PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm)
33*552f7358SJed Brown {
34*552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
35*552f7358SJed Brown   PetscMPIInt    num_proc, rank;
36*552f7358SJed Brown   PetscSection   coordSection;
37*552f7358SJed Brown   Vec            coordinates;
38*552f7358SJed Brown   PetscScalar   *coords;
39*552f7358SJed Brown   PetscInt       coordSize, v;
40*552f7358SJed Brown   PetscErrorCode ierr;
41*552f7358SJed Brown   /* Read from ex_get_init() */
42*552f7358SJed Brown   char           title[PETSC_MAX_PATH_LEN+1];
43*552f7358SJed Brown   int            dim = 0, numVertices = 0, numCells = 0;
44*552f7358SJed Brown   int            num_cs = 0, num_vs = 0, num_fs = 0;
45*552f7358SJed Brown #endif
46*552f7358SJed Brown 
47*552f7358SJed Brown   PetscFunctionBegin;
48*552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
49*552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
50*552f7358SJed Brown   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
51*552f7358SJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
52*552f7358SJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
53*552f7358SJed Brown   /* Open EXODUS II file and read basic informations on rank 0, then broadcast to all processors */
54*552f7358SJed Brown   if (!rank) {
55*552f7358SJed Brown     ierr = ex_get_init(exoid, title, &dim, &numVertices, &numCells, &num_cs, &num_vs, &num_fs);
56*552f7358SJed Brown     if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ExodusII ex_get_init() failed with error code %D\n",ierr);
57*552f7358SJed Brown     if (!num_cs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set\n");
58*552f7358SJed Brown   }
59*552f7358SJed Brown   ierr = MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm);CHKERRQ(ierr);
60*552f7358SJed Brown   ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
61*552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, title);CHKERRQ(ierr);
62*552f7358SJed Brown   ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr);
63*552f7358SJed Brown   ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr);
64*552f7358SJed Brown 
65*552f7358SJed Brown   /* Read cell sets information */
66*552f7358SJed Brown   if (!rank) {
67*552f7358SJed Brown     PetscInt *cone;
68*552f7358SJed Brown     int  c, cs, c_loc, v, v_loc;
69*552f7358SJed Brown     /* Read from ex_get_elem_blk_ids() */
70*552f7358SJed Brown     int *cs_id;
71*552f7358SJed Brown     /* Read from ex_get_elem_block() */
72*552f7358SJed Brown     char buffer[PETSC_MAX_PATH_LEN+1];
73*552f7358SJed Brown     int  num_cell_in_set, num_vertex_per_cell, num_attr;
74*552f7358SJed Brown     /* Read from ex_get_elem_conn() */
75*552f7358SJed Brown     int *cs_connect;
76*552f7358SJed Brown 
77*552f7358SJed Brown     /* Get cell sets IDs */
78*552f7358SJed Brown     ierr = PetscMalloc(num_cs * sizeof(int), &cs_id);CHKERRQ(ierr);
79*552f7358SJed Brown     ierr = ex_get_elem_blk_ids(exoid, cs_id);CHKERRQ(ierr);
80*552f7358SJed Brown     /* Read the cell set connectivity table and build mesh topology
81*552f7358SJed Brown        EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */
82*552f7358SJed Brown     /* First set sizes */
83*552f7358SJed Brown     for (cs = 0, c = 0; cs < num_cs; cs++) {
84*552f7358SJed Brown       ierr = ex_get_elem_block(exoid, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, &num_attr);CHKERRQ(ierr);
85*552f7358SJed Brown       for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
86*552f7358SJed Brown         ierr = DMPlexSetConeSize(*dm, c, num_vertex_per_cell);CHKERRQ(ierr);
87*552f7358SJed Brown       }
88*552f7358SJed Brown     }
89*552f7358SJed Brown     ierr = DMSetUp(*dm);CHKERRQ(ierr);
90*552f7358SJed Brown     for (cs = 0, c = 0; cs < num_cs; cs++) {
91*552f7358SJed Brown       ierr = ex_get_elem_block(exoid, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, &num_attr);CHKERRQ(ierr);
92*552f7358SJed Brown       ierr = PetscMalloc2(num_vertex_per_cell*num_cell_in_set,int,&cs_connect,num_vertex_per_cell,PetscInt,&cone);CHKERRQ(ierr);
93*552f7358SJed Brown       ierr = ex_get_elem_conn(exoid, cs_id[cs], cs_connect);CHKERRQ(ierr);
94*552f7358SJed Brown       /* EXO uses Fortran-based indexing, sieve uses C-style and numbers cell first then vertices. */
95*552f7358SJed Brown       for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
96*552f7358SJed Brown         for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) {
97*552f7358SJed Brown           cone[v_loc] = cs_connect[v]+numCells-1;
98*552f7358SJed Brown         }
99*552f7358SJed Brown         ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr);
100*552f7358SJed Brown         ierr = DMPlexSetLabelValue(*dm, "Cell Sets", c, cs_id[cs]);CHKERRQ(ierr);
101*552f7358SJed Brown       }
102*552f7358SJed Brown       ierr = PetscFree2(cs_connect,cone);CHKERRQ(ierr);
103*552f7358SJed Brown     }
104*552f7358SJed Brown     ierr = PetscFree(cs_id);CHKERRQ(ierr);
105*552f7358SJed Brown   }
106*552f7358SJed Brown   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
107*552f7358SJed Brown   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
108*552f7358SJed Brown   if (interpolate) {
109*552f7358SJed Brown     DM idm;
110*552f7358SJed Brown 
111*552f7358SJed Brown     switch(dim) {
112*552f7358SJed Brown     case 2:
113*552f7358SJed Brown       ierr = DMPlexInterpolate_2D(*dm, &idm);CHKERRQ(ierr);break;
114*552f7358SJed Brown     case 3:
115*552f7358SJed Brown       ierr = DMPlexInterpolate_3D(*dm, &idm);CHKERRQ(ierr);break;
116*552f7358SJed Brown     default:
117*552f7358SJed Brown       SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No mesh interpolation support for dimension %D", dim);
118*552f7358SJed Brown     }
119*552f7358SJed Brown     /* Maintain Cell Sets label */
120*552f7358SJed Brown     {
121*552f7358SJed Brown       DMLabel label;
122*552f7358SJed Brown 
123*552f7358SJed Brown       ierr = DMPlexRemoveLabel(*dm, "Cell Sets", &label);CHKERRQ(ierr);
124*552f7358SJed Brown       if (label) {ierr = DMPlexAddLabel(idm, label);CHKERRQ(ierr);}
125*552f7358SJed Brown     }
126*552f7358SJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
127*552f7358SJed Brown     *dm  = idm;
128*552f7358SJed Brown   }
129*552f7358SJed Brown 
130*552f7358SJed Brown   /* Create vertex set label */
131*552f7358SJed Brown   if (!rank && (num_vs > 0)) {
132*552f7358SJed Brown     int  vs, v;
133*552f7358SJed Brown     /* Read from ex_get_node_set_ids() */
134*552f7358SJed Brown     int *vs_id;
135*552f7358SJed Brown     /* Read from ex_get_node_set_param() */
136*552f7358SJed Brown     int  num_vertex_in_set, num_attr;
137*552f7358SJed Brown     /* Read from ex_get_node_set() */
138*552f7358SJed Brown     int *vs_vertex_list;
139*552f7358SJed Brown 
140*552f7358SJed Brown     /* Get vertex set ids */
141*552f7358SJed Brown     ierr = PetscMalloc(num_vs * sizeof(int), &vs_id);CHKERRQ(ierr);
142*552f7358SJed Brown     ierr = ex_get_node_set_ids(exoid, vs_id);CHKERRQ(ierr);
143*552f7358SJed Brown     for(vs = 0; vs < num_vs; ++vs) {
144*552f7358SJed Brown       ierr = ex_get_node_set_param(exoid, vs_id[vs], &num_vertex_in_set, &num_attr);CHKERRQ(ierr);
145*552f7358SJed Brown       ierr = PetscMalloc(num_vertex_in_set * sizeof(int), &vs_vertex_list);CHKERRQ(ierr);
146*552f7358SJed Brown       ierr = ex_get_node_set(exoid, vs_id[vs], vs_vertex_list);CHKERRQ(ierr);
147*552f7358SJed Brown       for(v = 0; v < num_vertex_in_set; ++v) {
148*552f7358SJed Brown         ierr = DMPlexSetLabelValue(*dm, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]);CHKERRQ(ierr);
149*552f7358SJed Brown       }
150*552f7358SJed Brown       ierr = PetscFree(vs_vertex_list);CHKERRQ(ierr);
151*552f7358SJed Brown     }
152*552f7358SJed Brown     ierr = PetscFree(vs_id);CHKERRQ(ierr);
153*552f7358SJed Brown   }
154*552f7358SJed Brown   /* Read coordinates */
155*552f7358SJed Brown   ierr = DMPlexGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
156*552f7358SJed Brown   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
157*552f7358SJed Brown   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
158*552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
159*552f7358SJed Brown   for (v = numCells; v < numCells+numVertices; ++v) {
160*552f7358SJed Brown     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
161*552f7358SJed Brown     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
162*552f7358SJed Brown   }
163*552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
164*552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
165*552f7358SJed Brown   ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr);
166*552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
167*552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
168*552f7358SJed Brown   ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr);
169*552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
170*552f7358SJed Brown   if (rank == 0) {
171*552f7358SJed Brown     float *x, *y, *z;
172*552f7358SJed Brown 
173*552f7358SJed Brown     ierr = PetscMalloc3(numVertices,float,&x,numVertices,float,&y,numVertices,float,&z);CHKERRQ(ierr);
174*552f7358SJed Brown     ierr = ex_get_coord(exoid, x, y, z);CHKERRQ(ierr);
175*552f7358SJed Brown     if (dim > 0) {for (v = 0; v < numVertices; ++v) {coords[v*dim+0] = x[v];}}
176*552f7358SJed Brown     if (dim > 1) {for (v = 0; v < numVertices; ++v) {coords[v*dim+1] = y[v];}}
177*552f7358SJed Brown     if (dim > 2) {for (v = 0; v < numVertices; ++v) {coords[v*dim+2] = z[v];}}
178*552f7358SJed Brown     ierr = PetscFree3(x,y,z);CHKERRQ(ierr);
179*552f7358SJed Brown   }
180*552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
181*552f7358SJed Brown   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
182*552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
183*552f7358SJed Brown 
184*552f7358SJed Brown   /* Create side set label */
185*552f7358SJed Brown   if (!rank && interpolate && (num_fs > 0)) {
186*552f7358SJed Brown     int  fs, f, voff;
187*552f7358SJed Brown     /* Read from ex_get_side_set_ids() */
188*552f7358SJed Brown     int *fs_id;
189*552f7358SJed Brown     /* Read from ex_get_side_set_param() */
190*552f7358SJed Brown     int  num_side_in_set, num_dist_fact_in_set;
191*552f7358SJed Brown     /* Read from ex_get_side_set_node_list() */
192*552f7358SJed Brown     int *fs_vertex_count_list, *fs_vertex_list;
193*552f7358SJed Brown 
194*552f7358SJed Brown     /* Get side set ids */
195*552f7358SJed Brown     ierr = PetscMalloc(num_fs * sizeof(int), &fs_id);CHKERRQ(ierr);
196*552f7358SJed Brown     ierr = ex_get_side_set_ids(exoid, fs_id);CHKERRQ(ierr);
197*552f7358SJed Brown     for(fs = 0; fs < num_fs; ++fs) {
198*552f7358SJed Brown       ierr = ex_get_side_set_param(exoid, fs_id[fs], &num_side_in_set, &num_dist_fact_in_set);CHKERRQ(ierr);
199*552f7358SJed Brown       ierr = PetscMalloc2(num_side_in_set,int,&fs_vertex_count_list,num_side_in_set*4,int,&fs_vertex_list);CHKERRQ(ierr);
200*552f7358SJed Brown       ierr = ex_get_side_set_node_list(exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list);CHKERRQ(ierr);
201*552f7358SJed Brown       for(f = 0, voff = 0; f < num_side_in_set; ++f) {
202*552f7358SJed Brown         const PetscInt *faces    = PETSC_NULL;
203*552f7358SJed Brown         PetscInt        faceSize = fs_vertex_count_list[f], numFaces;
204*552f7358SJed Brown         PetscInt        faceVertices[4], v;
205*552f7358SJed Brown 
206*552f7358SJed Brown         if (faceSize > 4) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize);
207*552f7358SJed Brown         for(v = 0; v < faceSize; ++v, ++voff) {
208*552f7358SJed Brown           faceVertices[v] = fs_vertex_list[voff]+numCells-1;
209*552f7358SJed Brown         }
210*552f7358SJed Brown         ierr = DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr);
211*552f7358SJed Brown         if (numFaces != 1) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces);
212*552f7358SJed Brown         ierr = DMPlexSetLabelValue(*dm, "Face Sets", faces[0], fs_id[fs]);CHKERRQ(ierr);
213*552f7358SJed Brown         ierr = DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr);
214*552f7358SJed Brown       }
215*552f7358SJed Brown       ierr = PetscFree2(fs_vertex_count_list,fs_vertex_list);CHKERRQ(ierr);
216*552f7358SJed Brown     }
217*552f7358SJed Brown     ierr = PetscFree(fs_id);CHKERRQ(ierr);
218*552f7358SJed Brown   }
219*552f7358SJed Brown #else
220*552f7358SJed Brown   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
221*552f7358SJed Brown #endif
222*552f7358SJed Brown   PetscFunctionReturn(0);
223*552f7358SJed Brown }
224