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: #undef  __FUNCT__
 10: /* MeshDistribute_Private
 11:   This function distributes the mesh identically among processors.

 13:   Collective on Mesh

 15:   Output Parameters:
 16: + mesh - The mesh
 17: - out  - Structure for communicating

 19:   Level: developer

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

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

 51: #undef  __FUNCT__
 53: /* MeshTriangular1D_Create_Simple
 54:   This function creates a simple 1D unstructured mesh.

 56:   Collective on Mesh

 58:   Input Parameter:
 59: . numCorners  - The number of nodes in an element

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

 67:   Output Parameter:
 68: . mesh        - The new mesh

 70:   Level: developer

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

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

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

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

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

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

135:   MeshDistribute_Private(mesh);
136:   return(0);
137: }

139: static int getNewEdge(Mesh mesh, double leftPoint, double rightPoint, PointFunction area, double *newRightPoint) {
140:   PetscScalar maxLen;
141:   double      len, midPoint;
142:   int         ierr;

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

158: #undef  __FUNCT__
160: /* MeshTriangular1D_Refine_Simple
161:   This function refines a simple 1D unstructured mesh until all segments are less than a given length.

163:   Collective on Mesh

165:   Input Parameters:
166: + oldMesh - The previous mesh
167: - area    - The PointFunction giving the area limit in each segment

169:   Output Parameter:
170: . mesh    - The new mesh

172:   Level: developer

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

192:   PetscMalloc(maxNodes * sizeof(double), &nodes);
193:   mesh->numEdges    = 0;
194:   mesh->numCorners  = numCorners;
195:   mesh->numNodes    = 1;
196:   mesh->numVertices = 1;
197:   nodes[mesh->numNodes] = oldNodes[oldEdges[0]];
198:   for(oldEdge = 0; oldEdge < numOldEdges; oldEdge++) {
199:     leftPoint     = oldNodes[oldEdges[oldEdge*numCorners]];
200:     rightPoint    = oldNodes[oldEdges[oldEdge*numCorners+1]];
201:     newRightPoint = leftPoint;

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

206:       if (mesh->numNodes >= maxNodes) {
207:         PetscMalloc(maxNodes*2 * sizeof(double), &tri->nodes);
208:         PetscMemcpy(tri->nodes, nodes, maxNodes*2 * sizeof(double));
209:         PetscFree(nodes);
210:         maxNodes *= 2;
211:         nodes     = tri->nodes;
212:       }
213:       mesh->numVertices++;
214:       for(corner = 1; corner < mesh->numCorners-1; corner++, mesh->numNodes++) {
215:         nodes[mesh->numNodes] = nodes[mesh->numNodes-1] + (newRightPoint - leftPoint)/(mesh->numCorners-1);
216:       }
217:       nodes[mesh->numNodes++] = newRightPoint;
218:       mesh->numEdges++;

220:       leftPoint = newRightPoint;
221:     }
222:   }

224:   PetscMalloc(mesh->numNodes * sizeof(double), &tri->nodes);
225:   PetscMalloc(mesh->numNodes * sizeof(int),    &tri->markers);
226:   PetscLogObjectMemory(mesh, (mesh->numNodes) * sizeof(double));
227:   PetscLogObjectMemory(mesh, (mesh->numNodes) * sizeof(int));
228:   PetscMemcpy(tri->nodes, nodes, (mesh->numNodes) * sizeof(double));
229:   PetscMemzero(tri->markers, mesh->numNodes * sizeof(int));
230:   tri->markers[0]                = oldMarkers[0];
231:   tri->markers[mesh->numNodes-1] = oldMarkers[oldMesh->numNodes-1];
232:   PetscFree(nodes);

234:   PetscMalloc(mesh->numEdges*numCorners * sizeof(int), &tri->edges);
235:   PetscLogObjectMemory(mesh, (mesh->numEdges*mesh->numCorners) * sizeof(int));
236:   for(edge = 0, node = 0; edge < mesh->numEdges; edge++) {
237:     tri->edges[edge*mesh->numCorners+0] = node;
238:     tri->edges[edge*mesh->numCorners+1] = node + mesh->numCorners-1;
239:     for(corner = 2, node++; corner < mesh->numCorners; corner++, node++) {
240:       tri->edges[edge*mesh->numCorners+corner] = node;
241:     }
242:   }

244:   PetscMalloc(mesh->numEdges*2 * sizeof(int), &tri->neighbors);
245:   PetscLogObjectMemory(mesh, (mesh->numEdges*2) * sizeof(int));
246:   for(edge = 0; edge < mesh->numEdges; edge++) {
247:     tri->neighbors[edge*2+0] = edge-1;
248:     tri->neighbors[edge*2+1] = edge+1;
249:   }
250:   tri->neighbors[mesh->numEdges*2-1] = -1;
251:   return(0);
252: }

254: int MeshTriangular1D_Create_Periodic(MeshBoundary2D *bdCtx, int numCorners, Mesh mesh) {
255:   SETERRQ(PETSC_ERR_SUP, " ")
256: }

258: int MeshTriangular1D_Refine_Periodic(Mesh oldMesh, PointFunction area, Mesh mesh) {
259:   SETERRQ(PETSC_ERR_SUP, " ")
260: }