Actual source code: tri1d_Simple.c

  1: #ifdef PETSC_RCS_HEADER
  2: static char vcid[] = "$Id: tri1d_Simple.c,v 1.6 2000/10/17 13:48:57 knepley Exp $";
  3: #endif

  5: #include "src/mesh/impls/triangular/1d/1dimpl.h"         /*I "mesh.h" I*/
  6: #include "tri1d_Simple.h"

  8: /* MeshDistribute_Private
  9:   This function distributes the mesh identically among processors.

 11:   Collective on Mesh

 13:   Output Parameters:
 14: + mesh - The mesh
 15: - out  - Structure for communicating

 17:   Level: developer

 19: .keywords mesh
 20: .seealso MeshTriangular1D_Create_Simple()
 21: */
 22: static int MeshDistribute_Private(Mesh mesh) {
 23:   Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
 24:   int              rank;
 25:   int              ierr;

 28:   MPI_Comm_rank(mesh->comm, &rank);
 29:   MPI_Bcast(&mesh->numBd,       1, MPI_INT, 0, mesh->comm);
 30:   MPI_Bcast(&mesh->numVertices, 1, MPI_INT, 0, mesh->comm);
 31:   MPI_Bcast(&mesh->numNodes,    1, MPI_INT, 0, mesh->comm);
 32:   MPI_Bcast(&mesh->numEdges,    1, MPI_INT, 0, mesh->comm);
 33:   MPI_Bcast(&mesh->numCorners,  1, MPI_INT, 0, mesh->comm);
 34:   if (rank != 0) {
 35:     PetscMalloc(mesh->numNodes*mesh->dim        * sizeof(double), &tri->nodes);
 36:     PetscMalloc(mesh->numNodes                  * sizeof(int),    &tri->markers);
 37:     PetscMalloc(mesh->numEdges*mesh->numCorners * sizeof(int),    &tri->edges);
 38:     PetscMalloc(mesh->numEdges*2                * sizeof(int),    &tri->neighbors);
 39:     PetscLogObjectMemory(mesh, (mesh->numNodes*mesh->dim) * sizeof(double));
 40:     PetscLogObjectMemory(mesh, (mesh->numNodes + mesh->numEdges*mesh->numCorners + mesh->numEdges*2) * sizeof(int));
 41:   }
 42:   MPI_Bcast(tri->nodes,     mesh->numNodes*mesh->dim,        MPI_DOUBLE, 0, mesh->comm);
 43:   MPI_Bcast(tri->markers,   mesh->numNodes,                  MPI_INT,    0, mesh->comm);
 44:   MPI_Bcast(tri->edges,     mesh->numEdges*mesh->numCorners, MPI_INT,    0, mesh->comm);
 45:   MPI_Bcast(tri->neighbors, mesh->numEdges*2,                MPI_INT,    0, mesh->comm);
 46:   return(0);
 47: }

 49: /* MeshTriangular1D_Create_Simple
 50:   This function creates a simple 1D unstructured mesh.

 52:   Collective on Mesh

 54:   Input Parameter:
 55: . numCorners  - The number of nodes in an element

 57:   Input Parameters from bdCtx:
 58: + numBD       - The number of closed boundaries in the geometry, or different markers
 59: . numVertices - The number of boundary points
 60: . vertices    - The (x,y) coordinates of the boundary points
 61: . markers     - The boundary markers for nodes, 0 indicates an interior point, each boundary must have a different marker

 63:   Output Parameter:
 64: . mesh        - The new mesh

 66:   Level: developer

 68: .keywords mesh, generation
 69: .seealso MeshTriangular1D_Refine_Simple(), MeshTriangular1D_Create_Periodic(), MeshTriangular1D_Refine_Periodic()
 70: */
 71: int MeshTriangular1D_Create_Simple(MeshBoundary2D *bdCtx, int numCorners, Mesh mesh) {
 72:   Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
 73:   int              rank;
 74:   int              ierr;

 77:   MPI_Comm_rank(mesh->comm, &rank);
 78:   if (rank == 0) {
 79:     double min = PetscMin(bdCtx->vertices[0], bdCtx->vertices[1]);
 80:     double max = PetscMax(bdCtx->vertices[0], bdCtx->vertices[1]);
 81:     double len = max - min;
 82:     double inc;
 83:     int    node, edge, corner;

 85:     /* I will interlace vertices and midnodes for easy refinement */
 86:     mesh->numBd       = bdCtx->numBd;
 87:     mesh->numCorners  = numCorners;
 88:     mesh->numVertices = 5;
 89:     if (mesh->numCorners == 2) {
 90:       mesh->numNodes  = 5;
 91:     } else if (mesh->numCorners == 3) {
 92:       mesh->numNodes  = 9;
 93:     } else {
 94:       SETERRQ1(PETSC_ERR_SUP, "Number of corners %d not supported", mesh->numCorners);
 95:     }

 97:     PetscMalloc(mesh->numNodes*mesh->dim * sizeof(double), &tri->nodes);
 98:     PetscMalloc(mesh->numNodes           * sizeof(int),    &tri->markers);
 99:     PetscLogObjectMemory(mesh, (mesh->numNodes*mesh->dim) * sizeof(double));
100:     PetscLogObjectMemory(mesh, (mesh->numNodes) * sizeof(int));
101:     tri->nodes[0]   = min;
102:     tri->markers[0] = bdCtx->markers[0];
103:     inc             = len/(mesh->numNodes-1);
104:     for(node = 1; node < mesh->numNodes-1; node++) {
105:       tri->nodes[node]   = tri->nodes[node-1] + inc;
106:       tri->markers[node] = 0;
107:     }
108:     tri->nodes[mesh->numNodes-1]   = max;
109:     tri->markers[mesh->numNodes-1] = bdCtx->markers[1];

111:     mesh->numEdges = mesh->numVertices-1;
112:     PetscMalloc(mesh->numEdges*mesh->numCorners * sizeof(int), &tri->edges);
113:     PetscLogObjectMemory(mesh, (mesh->numEdges*mesh->numCorners) * sizeof(int));
114:     for(edge = 0, node = 0; edge < mesh->numEdges; edge++) {
115:       tri->edges[edge*mesh->numCorners+0] = node;
116:       tri->edges[edge*mesh->numCorners+1] = node + mesh->numCorners-1;
117:       for(corner = 2, node++; corner < mesh->numCorners; corner++, node++) {
118:         tri->edges[edge*mesh->numCorners+corner] = node;
119:       }
120:     }

122:     PetscMalloc(mesh->numEdges*2 * sizeof(int), &tri->neighbors);
123:     PetscLogObjectMemory(mesh, (mesh->numEdges*2) * sizeof(int));
124:     for(edge = 0; edge < mesh->numEdges; edge++) {
125:       tri->neighbors[edge*2+0] = edge-1;
126:       tri->neighbors[edge*2+1] = edge+1;
127:     }
128:     tri->neighbors[mesh->numEdges*2-1] = -1;
129:   }

131:   MeshDistribute_Private(mesh);
132:   return(0);
133: }

135: static int getNewEdge(Mesh mesh, double leftPoint, double rightPoint, PointFunction area, double *newRightPoint) {
136:   PetscScalar maxLen;
137:   double      len, midPoint;
138:   int         ierr;

141:   len      = rightPoint - leftPoint;
142:   midPoint = leftPoint + len*0.5;
143:   (*area)(1, 1, &midPoint, PETSC_NULL, PETSC_NULL, &maxLen, mesh->usr);
144:   while (PetscAbsScalar(maxLen) < len) {
145:     rightPoint -= len*0.5;
146:     len         = rightPoint - leftPoint;
147:     midPoint    = leftPoint + len*0.5;
148:     (*area)(1, 1, &midPoint, PETSC_NULL, PETSC_NULL, &maxLen, mesh->usr);
149:   }
150:   *newRightPoint = rightPoint;
151:   return(0);
152: }

154: /* MeshTriangular1D_Refine_Simple
155:   This function refines a simple 1D unstructured mesh until all segments are less than a given length.

157:   Collective on Mesh

159:   Input Parameters:
160: + oldMesh - The previous mesh
161: - area    - The PointFunction giving the area limit in each segment

163:   Output Parameter:
164: . mesh    - The new mesh

166:   Level: developer

168: .keywords mesh, refinement
169: .seealso MeshTriangular1D_Create_Simple(), MeshTriangular1D_Create_Periodic(), MeshTriangular1D_Refine_Periodic()
170: */
171: int MeshTriangular1D_Refine_Simple(Mesh oldMesh, PointFunction area, Mesh mesh) {
172:   Mesh_Triangular *oldTri      = (Mesh_Triangular *) oldMesh->data;
173:   Mesh_Triangular *tri         = (Mesh_Triangular *) mesh->data;
174:   int              numCorners  = oldMesh->numCorners;
175:   int              numOldEdges = oldMesh->numEdges;
176:   double          *oldNodes    = oldTri->nodes;
177:   int             *oldMarkers  = oldTri->markers;
178:   int             *oldEdges    = oldTri->edges;
179:   int              maxNodes    = oldMesh->numNodes;
180:   double          *nodes;
181:   int              oldEdge, edge, corner, node;
182:   double           leftPoint, rightPoint, newRightPoint;
183:   int              ierr;

186:   PetscMalloc(maxNodes * sizeof(double), &nodes);
187:   mesh->numEdges    = 0;
188:   mesh->numCorners  = numCorners;
189:   mesh->numNodes    = 1;
190:   mesh->numVertices = 1;
191:   nodes[mesh->numNodes] = oldNodes[oldEdges[0]];
192:   for(oldEdge = 0; oldEdge < numOldEdges; oldEdge++) {
193:     leftPoint     = oldNodes[oldEdges[oldEdge*numCorners]];
194:     rightPoint    = oldNodes[oldEdges[oldEdge*numCorners+1]];
195:     newRightPoint = leftPoint;

197:     while(rightPoint != newRightPoint) {
198:       getNewEdge(oldMesh, leftPoint, rightPoint, area, &newRightPoint);

200:       if (mesh->numNodes >= maxNodes) {
201:         PetscMalloc(maxNodes*2 * sizeof(double), &tri->nodes);
202:         PetscMemcpy(tri->nodes, nodes, maxNodes*2 * sizeof(double));
203:         PetscFree(nodes);
204:         maxNodes *= 2;
205:         nodes     = tri->nodes;
206:       }
207:       mesh->numVertices++;
208:       for(corner = 1; corner < mesh->numCorners-1; corner++, mesh->numNodes++) {
209:         nodes[mesh->numNodes] = nodes[mesh->numNodes-1] + (newRightPoint - leftPoint)/(mesh->numCorners-1);
210:       }
211:       nodes[mesh->numNodes++] = newRightPoint;
212:       mesh->numEdges++;

214:       leftPoint = newRightPoint;
215:     }
216:   }

218:   PetscMalloc(mesh->numNodes * sizeof(double), &tri->nodes);
219:   PetscMalloc(mesh->numNodes * sizeof(int),    &tri->markers);
220:   PetscLogObjectMemory(mesh, (mesh->numNodes) * sizeof(double));
221:   PetscLogObjectMemory(mesh, (mesh->numNodes) * sizeof(int));
222:   PetscMemcpy(tri->nodes, nodes, (mesh->numNodes) * sizeof(double));
223:   PetscMemzero(tri->markers, mesh->numNodes * sizeof(int));
224:   tri->markers[0]                = oldMarkers[0];
225:   tri->markers[mesh->numNodes-1] = oldMarkers[oldMesh->numNodes-1];
226:   PetscFree(nodes);

228:   PetscMalloc(mesh->numEdges*numCorners * sizeof(int), &tri->edges);
229:   PetscLogObjectMemory(mesh, (mesh->numEdges*mesh->numCorners) * sizeof(int));
230:   for(edge = 0, node = 0; edge < mesh->numEdges; edge++) {
231:     tri->edges[edge*mesh->numCorners+0] = node;
232:     tri->edges[edge*mesh->numCorners+1] = node + mesh->numCorners-1;
233:     for(corner = 2, node++; corner < mesh->numCorners; corner++, node++) {
234:       tri->edges[edge*mesh->numCorners+corner] = node;
235:     }
236:   }

238:   PetscMalloc(mesh->numEdges*2 * sizeof(int), &tri->neighbors);
239:   PetscLogObjectMemory(mesh, (mesh->numEdges*2) * sizeof(int));
240:   for(edge = 0; edge < mesh->numEdges; edge++) {
241:     tri->neighbors[edge*2+0] = edge-1;
242:     tri->neighbors[edge*2+1] = edge+1;
243:   }
244:   tri->neighbors[mesh->numEdges*2-1] = -1;
245:   return(0);
246: }

248: int MeshTriangular1D_Create_Periodic(MeshBoundary2D *bdCtx, int numCorners, Mesh mesh) {
249:   SETERRQ(PETSC_ERR_SUP, " ")
250: }

252: int MeshTriangular1D_Refine_Periodic(Mesh oldMesh, PointFunction area, Mesh mesh) {
253:   SETERRQ(PETSC_ERR_SUP, " ")
254: }