Actual source code: mesh.c
1: #define PETSCDM_DLL
2:
3: #include petscmesh.h
4: #include petscmat.h
7: typedef struct _MeshOps *MeshOps;
8: struct _MeshOps {
9: PetscErrorCode (*view)(Mesh,PetscViewer);
10: PetscErrorCode (*createglobalvector)(Mesh,Vec*);
11: PetscErrorCode (*getcoloring)(Mesh,ISColoringType,ISColoring*);
12: PetscErrorCode (*getmatrix)(Mesh,MatType,Mat*);
13: PetscErrorCode (*getinterpolation)(Mesh,Mesh,Mat*,Vec*);
14: PetscErrorCode (*refine)(Mesh,MPI_Comm,Mesh*);
15: };
17: struct _p_Mesh {
18: PETSCHEADER(struct _MeshOps);
19: ALE::Obj<ALE::Two::Mesh> m;
20: Vec globalvector;
21: PetscInt bs,n,N,Nghosts,*ghosts;
22: PetscInt d_nz,o_nz,*d_nnz,*o_nnz;
23: };
27: PetscErrorCode WriteVTKHeader(PetscViewer viewer)
28: {
32: PetscViewerASCIIPrintf(viewer,"# vtk DataFile Version 2.0\n");
33: PetscViewerASCIIPrintf(viewer,"Simplicial Mesh Example\n");
34: PetscViewerASCIIPrintf(viewer,"ASCII\n");
35: PetscViewerASCIIPrintf(viewer,"DATASET UNSTRUCTURED_GRID\n");
36: return(0);
37: }
41: PetscErrorCode WriteVTKVertices(ALE::Obj<ALE::def::Mesh> mesh, PetscViewer viewer)
42: {
43: ALE::Obj<ALE::def::Mesh::coordinate_type> coordinates = mesh->getCoordinates();
44: const double *array = coordinates->restrict(0);
45: int dim = mesh->getDimension();
46: int numVertices = mesh->getTopology()->depthStratum(0)->size();
50: PetscViewerASCIIPrintf(viewer,"POINTS %d double\n", numVertices);
51: for(int v = 0; v < numVertices; v++) {
52: for(int d = 0; d < dim; d++) {
53: if (d > 0) {
54: PetscViewerASCIIPrintf(viewer," ");
55: }
56: PetscViewerASCIIPrintf(viewer,"%G",array[v*dim+d]);
57: }
58: for(int d = dim; d < 3; d++) {
59: PetscViewerASCIIPrintf(viewer," 0.0");
60: }
61: PetscViewerASCIIPrintf(viewer,"\n");
62: }
63: return(0);
64: }
68: PetscErrorCode WriteVTKVertices_New(ALE::Obj<ALE::Two::Mesh> mesh, PetscViewer viewer)
69: {
70: const double *array = mesh->getCoordinates()->restrict(ALE::Two::Mesh::field_type::patch_type());
71: int dim = mesh->getDimension();
72: int numVertices = mesh->getTopology()->depthStratum(0)->size();
76: PetscViewerASCIIPrintf(viewer, "POINTS %d double\n", numVertices);
77: for(int v = 0; v < numVertices; v++) {
78: for(int d = 0; d < dim; d++) {
79: if (d > 0) {
80: PetscViewerASCIIPrintf(viewer, " ");
81: }
82: PetscViewerASCIIPrintf(viewer, "%G", array[v*dim+d]);
83: }
84: for(int d = dim; d < 3; d++) {
85: PetscViewerASCIIPrintf(viewer, " 0.0");
86: }
87: PetscViewerASCIIPrintf(viewer, "\n");
88: }
89: return(0);
90: }
94: PetscErrorCode WriteVTKElements(ALE::Obj<ALE::def::Mesh> mesh, PetscViewer viewer)
95: {
96: MPI_Comm comm = mesh->getComm();
97: ALE::Obj<ALE::def::Mesh::sieve_type> topology = mesh->getTopology();
98: ALE::Obj<ALE::def::Mesh::sieve_type::heightSequence> elements = topology->heightStratum(0);
99: // FIX: Needs to be global
100: int numElements = elements->size();
101: int corners = topology->nCone(*elements->begin(), topology->depth())->size();
102: ALE::Obj<ALE::def::Mesh::bundle_type> vertexBundle = mesh->getBundle(0);
103: PetscMPIInt rank, size;
104: PetscErrorCode ierr;
107: MPI_Comm_rank(comm, &rank);
108: MPI_Comm_size(comm, &size);
109: PetscViewerASCIIPrintf(viewer,"CELLS %d %d\n", numElements, numElements*(corners+1));
110: if (rank == 0) {
111: for(ALE::def::Mesh::sieve_type::heightSequence::iterator e_itor = elements->begin(); e_itor != elements->end(); ++e_itor) {
112: ALE::Obj<ALE::def::PointArray> cone = topology->nCone(*e_itor, topology->depth());
114: PetscViewerASCIIPrintf(viewer, "%d ", corners);
115: for(ALE::def::PointArray::iterator c_itor = cone->begin(); c_itor != cone->end(); ++c_itor) {
116: PetscViewerASCIIPrintf(viewer, " %d", vertexBundle->getIndex(0, *c_itor).prefix);
117: }
118: PetscViewerASCIIPrintf(viewer, "\n");
119: }
120: #ifdef PARALLEL
121: for(int p = 1; p < size; p++) {
122: MPI_Status status;
123: int *remoteVertices;
124: int numLocalElements;
126: MPI_Recv(&numLocalElements, 1, MPI_INT, p, 1, comm, &status);
127: PetscMalloc(numLocalElements*corners * sizeof(int), &remoteVertices);
128: MPI_Recv(remoteVertices, numLocalElements*corners, MPI_INT, p, 1, comm, &status);
129: for(int e = 0; e < numLocalElements; e++) {
130: PetscViewerASCIIPrintf(viewer, "%d ", corners);
131: for(int c = 0; c < corners; c++) {
132: PetscViewerASCIIPrintf(viewer, " %d", remoteVertices[e*corners+c]);
133: }
134: PetscViewerASCIIPrintf(viewer, "\n");
135: }
136: PetscFree(remoteVertices);
137: }
138: #endif
139: } else {
140: #ifdef PARALLEL
141: int numLocalElements = elements.size(), offset = 0;
142: int *array;
144: PetscMalloc(numLocalElements*corners * sizeof(int), &array);
145: for(ALE::Point_set::iterator e_itor = elements.begin(); e_itor != elements.end(); e_itor++) {
146: ALE::Point_set cone = topology->nCone(*e_itor, dim);
147: for(ALE::Point_set::iterator c_itor = cone.begin(); c_itor != cone.end(); c_itor++) {
148: ALE::Point index = vertexBundle->getGlobalFiberInterval(*c_itor);
150: array[offset++] = index.prefix;
151: }
152: }
153: if (offset != numLocalElements*corners) {
154: SETERRQ2(PETSC_ERR_PLIB, "Invalid number of vertices to send %d shuold be %d", offset, numLocalElements*corners);
155: }
156: MPI_Send(&numLocalElements, 1, MPI_INT, 0, 1, comm);
157: MPI_Send(array, numLocalElements*corners, MPI_INT, 0, 1, comm);
158: PetscFree(array);
159: #endif
160: }
161: PetscViewerASCIIPrintf(viewer, "CELL_TYPES %d\n", numElements);
162: if (corners == 2) {
163: // VTK_LINE
164: for(int e = 0; e < numElements; e++) {
165: PetscViewerASCIIPrintf(viewer, "3\n");
166: }
167: } else if (corners == 3) {
168: // VTK_TRIANGLE
169: for(int e = 0; e < numElements; e++) {
170: PetscViewerASCIIPrintf(viewer, "5\n");
171: }
172: } else if (corners == 4) {
173: // VTK_TETRA
174: for(int e = 0; e < numElements; e++) {
175: PetscViewerASCIIPrintf(viewer, "10\n");
176: }
177: }
178: return(0);
179: }
183: PetscErrorCode WriteVTKElements_New(ALE::Obj<ALE::Two::Mesh> mesh, PetscViewer viewer)
184: {
185: MPI_Comm comm = mesh->getComm();
186: ALE::Obj<ALE::Two::Mesh::sieve_type> topology = mesh->getTopology();
187: ALE::Obj<ALE::Two::Mesh::sieve_type::heightSequence> elements = topology->heightStratum(0);
188: // FIX: Needs to be global
189: int numElements = elements->size();
190: int corners = topology->nCone(*elements->begin(), topology->depth())->size();
191: ALE::Obj<ALE::Two::Mesh::bundle_type> vertexBundle = mesh->getBundle(0);
192: PetscMPIInt rank, size;
193: PetscErrorCode ierr;
196: MPI_Comm_rank(comm, &rank);
197: MPI_Comm_size(comm, &size);
198: PetscViewerASCIIPrintf(viewer,"CELLS %d %d\n", numElements, numElements*(corners+1));
199: if (rank == 0) {
200: for(ALE::Two::Mesh::sieve_type::heightSequence::iterator e_itor = elements->begin(); e_itor != elements->end(); ++e_itor) {
201: ALE::Obj<ALE::def::PointArray> cone = topology->nCone(*e_itor, topology->depth());
202: ALE::Two::Mesh::bundle_type::patch_type patch;
204: PetscViewerASCIIPrintf(viewer, "%d ", corners);
205: for(ALE::def::PointArray::iterator c_itor = cone->begin(); c_itor != cone->end(); ++c_itor) {
206: PetscViewerASCIIPrintf(viewer, " %d", vertexBundle->getIndex(patch, *c_itor).prefix);
207: }
208: PetscViewerASCIIPrintf(viewer, "\n");
209: }
210: #ifdef PARALLEL
211: for(int p = 1; p < size; p++) {
212: MPI_Status status;
213: int *remoteVertices;
214: int numLocalElements;
216: MPI_Recv(&numLocalElements, 1, MPI_INT, p, 1, comm, &status);
217: PetscMalloc(numLocalElements*corners * sizeof(int), &remoteVertices);
218: MPI_Recv(remoteVertices, numLocalElements*corners, MPI_INT, p, 1, comm, &status);
219: for(int e = 0; e < numLocalElements; e++) {
220: PetscViewerASCIIPrintf(viewer, "%d ", corners);
221: for(int c = 0; c < corners; c++) {
222: PetscViewerASCIIPrintf(viewer, " %d", remoteVertices[e*corners+c]);
223: }
224: PetscViewerASCIIPrintf(viewer, "\n");
225: }
226: PetscFree(remoteVertices);
227: }
228: #endif
229: } else {
230: #ifdef PARALLEL
231: int numLocalElements = elements.size(), offset = 0;
232: int *array;
234: PetscMalloc(numLocalElements*corners * sizeof(int), &array);
235: for(ALE::Point_set::iterator e_itor = elements.begin(); e_itor != elements.end(); e_itor++) {
236: ALE::Point_set cone = topology->nCone(*e_itor, dim);
237: for(ALE::Point_set::iterator c_itor = cone.begin(); c_itor != cone.end(); c_itor++) {
238: ALE::Point index = vertexBundle->getGlobalFiberInterval(*c_itor);
240: array[offset++] = index.prefix;
241: }
242: }
243: if (offset != numLocalElements*corners) {
244: SETERRQ2(PETSC_ERR_PLIB, "Invalid number of vertices to send %d shuold be %d", offset, numLocalElements*corners);
245: }
246: MPI_Send(&numLocalElements, 1, MPI_INT, 0, 1, comm);
247: MPI_Send(array, numLocalElements*corners, MPI_INT, 0, 1, comm);
248: PetscFree(array);
249: #endif
250: }
251: PetscViewerASCIIPrintf(viewer, "CELL_TYPES %d\n", numElements);
252: if (corners == 2) {
253: // VTK_LINE
254: for(int e = 0; e < numElements; e++) {
255: PetscViewerASCIIPrintf(viewer, "3\n");
256: }
257: } else if (corners == 3) {
258: // VTK_TRIANGLE
259: for(int e = 0; e < numElements; e++) {
260: PetscViewerASCIIPrintf(viewer, "5\n");
261: }
262: } else if (corners == 4) {
263: // VTK_TETRA
264: for(int e = 0; e < numElements; e++) {
265: PetscViewerASCIIPrintf(viewer, "10\n");
266: }
267: }
268: return(0);
269: }
273: PetscErrorCode WritePCICEVertices(ALE::Obj<ALE::def::Mesh> mesh, PetscViewer viewer)
274: {
275: ALE::Obj<ALE::def::Mesh::coordinate_type> coordinates = mesh->getCoordinates();
276: const double *array = coordinates->restrict(0);
277: int dim = mesh->getDimension();
278: int numVertices = mesh->getTopology()->depthStratum(0)->size();
282: PetscViewerASCIIPrintf(viewer,"%D\n", numVertices);
283: for(int v = 0; v < numVertices; v++) {
284: PetscViewerASCIIPrintf(viewer,"%7D ", v+1);
285: for(int d = 0; d < dim; d++) {
286: if (d > 0) {
287: PetscViewerASCIIPrintf(viewer," ");
288: }
289: PetscViewerASCIIPrintf(viewer,"% 12.5E",array[v*dim+d]);
290: }
291: PetscViewerASCIIPrintf(viewer,"\n");
292: }
293: return(0);
294: }
298: PetscErrorCode WritePCICEElements(ALE::Obj<ALE::def::Mesh> mesh, PetscViewer viewer)
299: {
300: MPI_Comm comm = mesh->getComm();
301: int dim = mesh->getDimension();
302: ALE::Obj<ALE::def::Mesh::sieve_type> topology = mesh->getTopology();
303: ALE::Obj<ALE::def::Mesh::sieve_type> orientation = mesh->getOrientation();
304: ALE::Obj<ALE::def::Mesh::sieve_type::heightSequence> elements = topology->heightStratum(0);
305: ALE::Obj<ALE::def::Mesh::bundle_type> vertexBundle = mesh->getBundle(0);
306: // FIX: Needs to be global
307: int numElements = elements->size();
308: int corners = topology->nCone(*elements->begin(), dim)->size();
309: PetscMPIInt rank, size;
310: int elementCount = 1;
311: PetscErrorCode ierr;
314: MPI_Comm_rank(comm, &rank);
315: MPI_Comm_size(comm, &size);
316: if (corners != dim+1) {
317: SETERRQ(PETSC_ERR_SUP, "PCICE only supports simplicies");
318: }
319: if (rank == 0) {
320: PetscViewerASCIIPrintf(viewer, "%d\n", numElements);
321: for(ALE::def::Mesh::sieve_type::heightSequence::iterator e_itor = elements->begin(); e_itor != elements->end(); ++e_itor) {
322: // FIX: Need to globalize
323: ALE::Obj<ALE::def::Mesh::bundle_type::IndexArray> intervals = vertexBundle->getOrderedIndices(0, orientation->cone(*e_itor));
325: PetscViewerASCIIPrintf(viewer, "%7d", elementCount++);
326: for(ALE::def::Mesh::bundle_type::IndexArray::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
327: PetscViewerASCIIPrintf(viewer, " %7d", i_itor->prefix+1);
328: }
329: PetscViewerASCIIPrintf(viewer, "\n");
330: }
331: #ifdef PARALLEL
332: for(int p = 1; p < size; p++) {
333: MPI_Status status;
334: int *remoteVertices;
335: int numLocalElements;
337: MPI_Recv(&numLocalElements, 1, MPI_INT, p, 1, comm, &status);
338: PetscMalloc(numLocalElements*corners * sizeof(int), &remoteVertices);
339: MPI_Recv(remoteVertices, numLocalElements*corners, MPI_INT, p, 1, comm, &status);
340: for(int e = 0; e < numLocalElements; e++) {
341: PetscViewerASCIIPrintf(viewer, "%7d", elementCount++);
342: for(int c = 0; c < corners; c++) {
343: PetscViewerASCIIPrintf(viewer, " %7d", remoteVertices[e*corners+c]);
344: }
345: PetscViewerASCIIPrintf(viewer, "\n");
346: }
347: PetscFree(remoteVertices);
348: }
349: #endif
350: } else {
351: #ifdef PARALLEL
352: int numLocalElements = elements.size(), offset = 0;
353: int *array;
355: PetscMalloc(numLocalElements*corners * sizeof(int), &array);
356: for(ALE::Point_set::iterator e_itor = elements.begin(); e_itor != elements.end(); e_itor++) {
357: ALE::Obj<ALE::Point_array> intervals = vertexBundle.getGlobalOrderedClosureIndices(orientation->cone(*e_itor));
359: for(ALE::Point_array::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
360: array[offset++] = (*i_itor).prefix+1;
361: }
362: }
363: if (offset != numLocalElements*corners) {
364: SETERRQ2(PETSC_ERR_PLIB, "Invalid number of vertices to send %d shuold be %d", offset, numLocalElements*corners);
365: }
366: MPI_Send(&numLocalElements, 1, MPI_INT, 0, 1, comm);
367: MPI_Send(array, numLocalElements*corners, MPI_INT, 0, 1, comm);
368: PetscFree(array);
369: #endif
370: }
371: return(0);
372: }
376: PetscErrorCode WritePyLithVertices(ALE::Obj<ALE::def::Mesh> mesh, PetscViewer viewer)
377: {
378: ALE::Obj<ALE::def::Mesh::coordinate_type> coordinates = mesh->getCoordinates();
379: const double *array = coordinates->restrict(0);
380: int dim = mesh->getDimension();
381: int numVertices = mesh->getTopology()->depthStratum(0)->size();
385: PetscViewerASCIIPrintf(viewer,"#\n");
386: PetscViewerASCIIPrintf(viewer,"coord_units = km\n");
387: PetscViewerASCIIPrintf(viewer,"#\n");
388: PetscViewerASCIIPrintf(viewer,"#\n");
389: PetscViewerASCIIPrintf(viewer,"# Node X-coord Y-coord Z-coord\n");
390: PetscViewerASCIIPrintf(viewer,"#\n");
391: for(int v = 0; v < numVertices; v++) {
392: PetscViewerASCIIPrintf(viewer,"%7D ", v+1);
393: for(int d = 0; d < dim; d++) {
394: if (d > 0) {
395: PetscViewerASCIIPrintf(viewer," ");
396: }
397: PetscViewerASCIIPrintf(viewer,"% 16.8E",array[v*dim+d]);
398: }
399: PetscViewerASCIIPrintf(viewer,"\n");
400: }
401: return(0);
402: }
406: PetscErrorCode WritePyLithElements(ALE::Obj<ALE::def::Mesh> mesh, PetscViewer viewer)
407: {
408: MPI_Comm comm = mesh->getComm();
409: int dim = mesh->getDimension();
410: ALE::Obj<ALE::def::Mesh::sieve_type> topology = mesh->getTopology();
411: ALE::Obj<ALE::def::Mesh::sieve_type> orientation = mesh->getOrientation();
412: ALE::Obj<ALE::def::Mesh::sieve_type::heightSequence> elements = topology->heightStratum(0);
413: ALE::Obj<ALE::def::Mesh::bundle_type> vertexBundle = mesh->getBundle(0);
414: // FIX: Needs to be global
415: int corners = topology->nCone(*elements->begin(), topology->depth())->size();
416: PetscMPIInt rank, size;
417: int elementCount = 1;
418: PetscErrorCode ierr;
421: MPI_Comm_rank(comm, &rank);
422: MPI_Comm_size(comm, &size);
423: if (dim != 3) {
424: SETERRQ(PETSC_ERR_SUP, "PyLith only supports 3D meshes.");
425: }
426: if (corners != 4) {
427: SETERRQ(PETSC_ERR_SUP, "We only support linear tetrahedra for PyLith.");
428: }
429: PetscViewerASCIIPrintf(viewer,"#\n");
430: PetscViewerASCIIPrintf(viewer,"# N ETP MAT INF N1 N2 N3 N4 N5 N6 N7 N8\n");
431: PetscViewerASCIIPrintf(viewer,"#\n");
432: if (rank == 0) {
433: for(ALE::def::Mesh::sieve_type::heightSequence::iterator e_itor = elements->begin(); e_itor != elements->end(); ++e_itor) {
434: // FIX: Need to globalize
435: ALE::Obj<ALE::def::Mesh::bundle_type::IndexArray> intervals = vertexBundle->getOrderedIndices(0, orientation->cone(*e_itor));
437: // Only linear tetrahedra, 1 material, no infinite elements
438: PetscViewerASCIIPrintf(viewer, "%7d %3d %3d %3d", elementCount++, 5, 1, 0);
439: for(ALE::def::Mesh::bundle_type::IndexArray::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
440: PetscViewerASCIIPrintf(viewer, " %6d", (*i_itor).prefix+1);
441: }
442: PetscViewerASCIIPrintf(viewer, "\n");
443: }
444: #ifdef PARALLEL
445: for(int p = 1; p < size; p++) {
446: MPI_Status status;
447: int *remoteVertices;
448: int numLocalElements;
450: MPI_Recv(&numLocalElements, 1, MPI_INT, p, 1, comm, &status);
451: PetscMalloc(numLocalElements*corners * sizeof(int), &remoteVertices);
452: MPI_Recv(remoteVertices, numLocalElements*corners, MPI_INT, p, 1, comm, &status);
453: for(int e = 0; e < numLocalElements; e++) {
454: // Only linear tetrahedra, 1 material, no infinite elements
455: PetscViewerASCIIPrintf(viewer, "%7d %3d %3d %3d", elementCount++, 5, 1, 0);
456: for(int c = 0; c < corners; c++) {
457: PetscViewerASCIIPrintf(viewer, " %6d", remoteVertices[e*corners+c]);
458: }
459: PetscViewerASCIIPrintf(viewer, "\n");
460: }
461: PetscFree(remoteVertices);
462: }
463: #endif
464: } else {
465: #ifdef PARALLEL
466: int numLocalElements = elements.size(), offset = 0;
467: int *array;
469: PetscMalloc(numLocalElements*corners * sizeof(int), &array);
470: for(ALE::Point_set::iterator e_itor = elements.begin(); e_itor != elements.end(); e_itor++) {
471: ALE::Obj<ALE::Point_array> intervals = vertexBundle.getGlobalOrderedClosureIndices(orientation->cone(*e_itor));
473: for(ALE::Point_array::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
474: array[offset++] = (*i_itor).prefix+1;
475: }
476: }
477: if (offset != numLocalElements*corners) {
478: SETERRQ2(PETSC_ERR_PLIB, "Invalid number of vertices to send %d shuold be %d", offset, numLocalElements*corners);
479: }
480: MPI_Send(&numLocalElements, 1, MPI_INT, 0, 1, comm);
481: MPI_Send(array, numLocalElements*corners, MPI_INT, 0, 1, comm);
482: PetscFree(array);
483: #endif
484: }
485: return(0);
486: }
490: PetscErrorCode WritePyLithVerticesLocal(ALE::Obj<ALE::def::Mesh> mesh, PetscViewer viewer)
491: {
492: ALE::Obj<ALE::def::Mesh::coordinate_type> coordinates = mesh->getCoordinates();
493: const double *array = coordinates->restrict(0);
494: int dim = mesh->getDimension();
495: int numVertices = mesh->getTopology()->depthStratum(0)->size();
499: PetscViewerASCIIPrintf(viewer,"#\n");
500: PetscViewerASCIIPrintf(viewer,"coord_units = km\n");
501: PetscViewerASCIIPrintf(viewer,"#\n");
502: PetscViewerASCIIPrintf(viewer,"#\n");
503: PetscViewerASCIIPrintf(viewer,"# Node X-coord Y-coord Z-coord\n");
504: PetscViewerASCIIPrintf(viewer,"#\n");
505: for(int v = 0; v < numVertices; v++) {
506: PetscViewerASCIIPrintf(viewer,"%7D ", v+1);
507: for(int d = 0; d < dim; d++) {
508: if (d > 0) {
509: PetscViewerASCIIPrintf(viewer," ");
510: }
511: PetscViewerASCIIPrintf(viewer,"% 16.8E",array[v*dim+d]);
512: }
513: PetscViewerASCIIPrintf(viewer,"\n");
514: }
515: return(0);
516: }
520: PetscErrorCode WritePyLithElementsLocal(ALE::Obj<ALE::def::Mesh> mesh, PetscViewer viewer)
521: {
522: int dim = mesh->getDimension();
523: ALE::Obj<ALE::def::Mesh::sieve_type> topology = mesh->getTopology();
524: ALE::Obj<ALE::def::Mesh::sieve_type> orientation = mesh->getOrientation();
525: ALE::Obj<ALE::def::Mesh::sieve_type::heightSequence> elements = topology->heightStratum(0);
526: ALE::Obj<ALE::def::Mesh::bundle_type> vertexBundle = mesh->getBundle(0);
527: int corners = topology->nCone(*elements->begin(), topology->depth())->size();
528: int elementCount = 1;
532: if (dim != 3) {
533: SETERRQ(PETSC_ERR_SUP, "PyLith only supports 3D meshes.");
534: }
535: if (corners != 4) {
536: SETERRQ(PETSC_ERR_SUP, "We only support linear tetrahedra for PyLith.");
537: }
538: PetscViewerASCIIPrintf(viewer,"#\n");
539: PetscViewerASCIIPrintf(viewer,"# N ETP MAT INF N1 N2 N3 N4 N5 N6 N7 N8\n");
540: PetscViewerASCIIPrintf(viewer,"#\n");
541: for(ALE::def::Mesh::sieve_type::heightSequence::iterator e_itor = elements->begin(); e_itor != elements->end(); ++e_itor) {
542: ALE::Obj<ALE::def::Mesh::bundle_type::IndexArray> intervals = vertexBundle->getOrderedIndices(0, orientation->cone(*e_itor));
544: // Only linear tetrahedra, 1 material, no infinite elements
545: PetscViewerASCIIPrintf(viewer, "%7d %3d %3d %3d", elementCount++, 5, 1, 0);
546: for(ALE::def::Mesh::bundle_type::IndexArray::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
547: PetscViewerASCIIPrintf(viewer, " %6d", (*i_itor).prefix+1);
548: }
549: PetscViewerASCIIPrintf(viewer, "\n");
550: }
551: return(0);
552: }
556: PetscErrorCode MeshView_Sieve_Ascii(ALE::Obj<ALE::def::Mesh> mesh, PetscViewer viewer)
557: {
558: PetscViewerFormat format;
559: PetscErrorCode ierr;
562: PetscViewerGetFormat(viewer, &format);
563: if (format == PETSC_VIEWER_ASCII_VTK) {
564: WriteVTKHeader(viewer);
565: WriteVTKVertices(mesh, viewer);
566: WriteVTKElements(mesh, viewer);
567: } else if (format == PETSC_VIEWER_ASCII_PCICE) {
568: char *filename;
569: char coordFilename[2048];
570: PetscTruth isConnect;
571: size_t len;
573: PetscViewerFileGetName(viewer, &filename);
574: PetscStrlen(filename, &len);
575: PetscStrcmp(&(filename[len-5]), ".lcon", &isConnect);
576: if (!isConnect) {
577: SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid element connectivity filename: %s", filename);
578: }
579: WritePCICEElements(mesh, viewer);
580: PetscStrncpy(coordFilename, filename, len-5);
581: coordFilename[len-5] = '\0';
582: PetscStrcat(coordFilename, ".nodes");
583: PetscViewerFileSetName(viewer, coordFilename);
584: WritePCICEVertices(mesh, viewer);
585: } else if (format == PETSC_VIEWER_ASCII_PYLITH) {
586: char *filename;
587: char connectFilename[2048];
588: char coordFilename[2048];
590: PetscViewerFileGetName(viewer, &filename);
591: PetscViewerFileSetMode(viewer, FILE_MODE_WRITE);
592: PetscStrcpy(connectFilename, filename);
593: PetscStrcat(connectFilename, ".connect");
594: PetscViewerFileSetName(viewer, connectFilename);
595: WritePyLithElements(mesh, viewer);
596: PetscStrcpy(coordFilename, filename);
597: PetscStrcat(coordFilename, ".coord");
598: PetscViewerFileSetName(viewer, coordFilename);
599: WritePyLithVertices(mesh, viewer);
600: PetscViewerFileSetMode(viewer, FILE_MODE_READ);
601: PetscExceptionTry1(PetscViewerFileSetName(viewer, filename), PETSC_ERR_FILE_OPEN);
602: if (PetscExceptionValue(ierr)) {
603: /* this means that a caller above me has also tryed this exception so I don't handle it here, pass it up */
604: } else if (PetscExceptionCaught(ierr, PETSC_ERR_FILE_OPEN)) {
605: 0;
606: }
607:
608: #if 0
609: } else if (format == PETSC_VIEWER_ASCII_PYLITH_LOCAL) {
610: PetscViewer connectViewer, coordViewer;
611: char *filename;
612: char localFilename[2048];
613: MPI_Comm comm;
614: PetscMPIInt rank;
616: PetscObjectGetComm((PetscObject) mesh, &comm);
617: MPI_Comm_rank(comm, &rank);
618: PetscViewerFileGetName(viewer, &filename);
620: sprintf(localFilename, "%s.%d.connect", filename, rank);
621: PetscViewerCreate(PETSC_COMM_SELF, &connectViewer);
622: PetscViewerSetType(connectViewer, PETSC_VIEWER_ASCII);
623: PetscViewerSetFormat(connectViewer, PETSC_VIEWER_ASCII_PYLITH);
624: PetscViewerFileSetName(connectViewer, localFilename);
625: WritePyLithElementsLocal(mesh, connectViewer);
626: PetscViewerDestroy(connectViewer);
628: sprintf(localFilename, "%s.%d.coord", filename, rank);
629: PetscViewerCreate(PETSC_COMM_SELF, &coordViewer);
630: PetscViewerSetType(coordViewer, PETSC_VIEWER_ASCII);
631: PetscViewerSetFormat(coordViewer, PETSC_VIEWER_ASCII_PYLITH);
632: PetscViewerFileSetName(coordViewer, localFilename);
633: WritePyLithVerticesLocal(mesh, coordViewer);
634: PetscViewerDestroy(coordViewer);
635: #endif
636: } else {
637: int dim = mesh->getDimension();
639: PetscViewerASCIIPrintf(viewer, "Mesh in %d dimensions:\n", dim);
640: for(int d = 0; d <= dim; d++) {
641: // FIX: Need to globalize
642: PetscViewerASCIIPrintf(viewer, " %d %d-cells\n", mesh->getTopology()->depthStratum(d)->size(), d);
643: }
644: }
645: PetscViewerFlush(viewer);
646: return(0);
647: }
651: PetscErrorCode MeshView_Sieve_Ascii(ALE::Obj<ALE::Two::Mesh> mesh, PetscViewer viewer)
652: {
653: PetscViewerFormat format;
654: PetscErrorCode ierr;
657: PetscViewerGetFormat(viewer, &format);
658: if (format == PETSC_VIEWER_ASCII_VTK) {
659: WriteVTKHeader(viewer);
660: WriteVTKVertices_New(mesh, viewer);
661: WriteVTKElements_New(mesh, viewer);
662: } else {
663: int dim = mesh->getDimension();
665: PetscViewerASCIIPrintf(viewer, "Mesh in %d dimensions:\n", dim);
666: for(int d = 0; d <= dim; d++) {
667: // FIX: Need to globalize
668: PetscViewerASCIIPrintf(viewer, " %d %d-cells\n", mesh->getTopology()->depthStratum(d)->size(), d);
669: }
670: }
671: PetscViewerFlush(viewer);
672: return(0);
673: }
677: PetscErrorCode MeshView_Sieve_Newer(ALE::Obj<ALE::Two::Mesh> mesh, PetscViewer viewer)
678: {
679: PetscTruth iascii, isbinary, isdraw;
683: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII, &iascii);
684: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_BINARY, &isbinary);
685: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_DRAW, &isdraw);
687: if (iascii){
688: MeshView_Sieve_Ascii(mesh, viewer);
689: } else if (isbinary) {
690: SETERRQ(PETSC_ERR_SUP, "Binary viewer not implemented for Mesh");
691: } else if (isdraw){
692: SETERRQ(PETSC_ERR_SUP, "Draw viewer not implemented for Mesh");
693: } else {
694: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by this mesh object", ((PetscObject)viewer)->type_name);
695: }
696: return(0);
697: }
701: PetscErrorCode MeshView_Sieve_New(ALE::Obj<ALE::def::Mesh> mesh, PetscViewer viewer)
702: {
703: PetscTruth iascii, isbinary, isdraw;
707: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII, &iascii);
708: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_BINARY, &isbinary);
709: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_DRAW, &isdraw);
711: if (iascii){
712: MeshView_Sieve_Ascii(mesh, viewer);
713: } else if (isbinary) {
714: SETERRQ(PETSC_ERR_SUP, "Binary viewer not implemented for Mesh");
715: } else if (isdraw){
716: SETERRQ(PETSC_ERR_SUP, "Draw viewer not implemented for Mesh");
717: } else {
718: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by this mesh object", ((PetscObject)viewer)->type_name);
719: }
720: return(0);
721: }
725: PetscErrorCode MeshView_Sieve(Mesh mesh, PetscViewer viewer)
726: {
727: PetscTruth iascii, isbinary, isdraw;
731: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII, &iascii);
732: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_BINARY, &isbinary);
733: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_DRAW, &isdraw);
735: if (iascii){
736: SETERRQ(PETSC_ERR_SUP, "Ascii viewer not implemented for Mesh");
737: } else if (isbinary) {
738: SETERRQ(PETSC_ERR_SUP, "Binary viewer not implemented for Mesh");
739: } else if (isdraw){
740: SETERRQ(PETSC_ERR_SUP, "Draw viewer not implemented for Mesh");
741: } else {
742: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by this mesh object", ((PetscObject)viewer)->type_name);
743: }
744: return(0);
745: }
749: /*@C
750: MeshView - Views a Mesh object.
752: Collective on Mesh
754: Input Parameters:
755: + mesh - the mesh
756: - viewer - an optional visualization context
758: Notes:
759: The available visualization contexts include
760: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
761: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
762: output where only the first processor opens
763: the file. All other processors send their
764: data to the first processor to print.
766: You can change the format the mesh is printed using the
767: option PetscViewerSetFormat().
769: The user can open alternative visualization contexts with
770: + PetscViewerASCIIOpen() - Outputs mesh to a specified file
771: . PetscViewerBinaryOpen() - Outputs mesh in binary to a
772: specified file; corresponding input uses MeshLoad()
773: . PetscViewerDrawOpen() - Outputs mesh to an X window display
775: The user can call PetscViewerSetFormat() to specify the output
776: format of ASCII printed objects (when using PETSC_VIEWER_STDOUT_SELF,
777: PETSC_VIEWER_STDOUT_WORLD and PetscViewerASCIIOpen). Available formats include
778: + PETSC_VIEWER_ASCII_DEFAULT - default, prints mesh information
779: - PETSC_VIEWER_ASCII_VTK - outputs a VTK file describing the mesh
781: Level: beginner
783: Concepts: mesh^printing
784: Concepts: mesh^saving to disk
786: .seealso: PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerBinaryOpen(),
787: MeshLoad(), PetscViewerCreate()
788: @*/
789: PetscErrorCode PETSCDM_DLLEXPORT MeshView(Mesh mesh, PetscViewer viewer)
790: {
796: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(mesh->comm);
800: (*mesh->ops->view)(mesh, viewer);
801: return(0);
802: }
806: /*@C
807: MeshLoad - Create a mesh topology from the saved data in a viewer.
809: Collective on Viewer
811: Input Parameter:
812: . viewer - The viewer containing the data
814: Output Parameters:
815: . mesh - the mesh object
817: Level: advanced
819: .seealso MeshView()
821: @*/
822: PetscErrorCode PETSCDM_DLLEXPORT MeshLoad(PetscViewer viewer, Mesh *mesh)
823: {
824: SETERRQ(PETSC_ERR_SUP, "");
825: }
829: /*@C
830: MeshGetMesh - Gets the internal mesh object
832: Not collective
834: Input Parameter:
835: . mesh - the mesh object
837: Output Parameter:
838: . m - the internal mesh object
839:
840: Level: advanced
842: .seealso MeshCreate(), MeshSetMesh()
844: @*/
845: PetscErrorCode PETSCDM_DLLEXPORT MeshGetMesh(Mesh mesh, ALE::Obj<ALE::Two::Mesh> *m)
846: {
849: if (m) {
851: *m = mesh->m;
852: }
853: return(0);
854: }
858: /*@C
859: MeshSetMesh - Sets the internal mesh object
861: Not collective
863: Input Parameters:
864: + mesh - the mesh object
865: - boundary - the internal mesh object
866:
867: Level: advanced
869: .seealso MeshCreate(), MeshGetMesh()
871: @*/
872: PetscErrorCode PETSCDM_DLLEXPORT MeshSetMesh(Mesh mesh, ALE::Obj<ALE::Two::Mesh> m)
873: {
876: mesh->m = m;
877: return(0);
878: }
882: /*@C
883: MeshGetMatrix - Creates a matrix with the correct parallel layout required for
884: computing the Jacobian on a function defined using the informatin in Mesh.
886: Collective on Mesh
888: Input Parameter:
889: + mesh - the mesh object
890: - mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, MATSEQSBAIJ, MATMPISBAIJ,
891: or any type which inherits from one of these (such as MATAIJ, MATLUSOL, etc.).
893: Output Parameters:
894: . J - matrix with the correct nonzero preallocation
895: (obviously without the correct Jacobian values)
897: Level: advanced
899: Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
900: do not need to do it yourself.
902: .seealso ISColoringView(), ISColoringGetIS(), MatFDColoringCreate(), DASetBlockFills()
904: @*/
905: PetscErrorCode PETSCDM_DLLEXPORT MeshGetMatrix(Mesh mesh, MatType mtype,Mat *J)
906: {
907: ALE::Obj<ALE::Two::Mesh> m;
908: #if 0
909: ISLocalToGlobalMapping lmap;
910: PetscInt *globals,rstart,i;
911: #endif
912: PetscInt localSize, globalSize;
913: PetscErrorCode ierr;
916: MeshGetMesh(mesh, &m);
917: localSize = m->getField("u")->getSize(ALE::Two::Mesh::field_type::patch_type());
918: globalSize = localSize;
920: MatCreate(mesh->comm,J);
921: MatSetSizes(*J,localSize,localSize,globalSize,globalSize);
922: MatSetType(*J,mtype);
923: MatSetBlockSize(*J,1);
924: MatSeqAIJSetPreallocation(*J,mesh->d_nz,mesh->d_nnz);
925: MatMPIAIJSetPreallocation(*J,mesh->d_nz,mesh->d_nnz,mesh->o_nz,mesh->o_nnz);
926: MatSeqBAIJSetPreallocation(*J,mesh->bs,mesh->d_nz,mesh->d_nnz);
927: MatMPIBAIJSetPreallocation(*J,mesh->bs,mesh->d_nz,mesh->d_nnz,mesh->o_nz,mesh->o_nnz);
929: #if 0
930: PetscMalloc((mesh->n+mesh->Nghosts+1)*sizeof(PetscInt),&globals);
931: MatGetOwnershipRange(*J,&rstart,PETSC_NULL);
932: for (i=0; i<mesh->n; i++) {
933: globals[i] = rstart + i;
934: }
935: PetscMemcpy(globals+mesh->n,mesh->ghosts,mesh->Nghosts*sizeof(PetscInt));
936: ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,mesh->n+mesh->Nghosts,globals,&lmap);
937: PetscFree(globals);
938: MatSetLocalToGlobalMapping(*J,lmap);
939: ISLocalToGlobalMappingDestroy(lmap);
940: #endif
941: return(0);
942: }
946: /*@C
947: MeshSetGhosts - Sets the global indices of other processes elements that will
948: be ghosts on this process
950: Not Collective
952: Input Parameters:
953: + mesh - the Mesh object
954: . bs - block size
955: . nlocal - number of local (non-ghost) entries
956: . Nghosts - number of ghosts on this process
957: - ghosts - indices of all the ghost points
959: Level: advanced
961: .seealso MeshDestroy(), MeshCreateGlobalVector(), MeshGetGlobalIndices()
963: @*/
964: PetscErrorCode PETSCDM_DLLEXPORT MeshSetGhosts(Mesh mesh,PetscInt bs,PetscInt nlocal,PetscInt Nghosts,const PetscInt ghosts[])
965: {
970: PetscFree(mesh->ghosts);
971: PetscMalloc((1+Nghosts)*sizeof(PetscInt),&mesh->ghosts);
972: PetscMemcpy(mesh->ghosts,ghosts,Nghosts*sizeof(PetscInt));
973: mesh->bs = bs;
974: mesh->n = nlocal;
975: mesh->Nghosts = Nghosts;
976: return(0);
977: }
981: /*@C
982: MeshSetPreallocation - sets the matrix memory preallocation for matrices computed by Mesh
984: Not Collective
986: Input Parameters:
987: + mesh - the Mesh object
988: . d_nz - maximum number of nonzeros in any row of diagonal block
989: . d_nnz - number of nonzeros in each row of diagonal block
990: . o_nz - maximum number of nonzeros in any row of off-diagonal block
991: . o_nnz - number of nonzeros in each row of off-diagonal block
994: Level: advanced
996: .seealso MeshDestroy(), MeshCreateGlobalVector(), MeshGetGlobalIndices(), MatMPIAIJSetPreallocation(),
997: MatMPIBAIJSetPreallocation()
999: @*/
1000: PetscErrorCode PETSCDM_DLLEXPORT MeshSetPreallocation(Mesh mesh,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
1001: {
1004: mesh->d_nz = d_nz;
1005: mesh->d_nnz = (PetscInt*)d_nnz;
1006: mesh->o_nz = o_nz;
1007: mesh->o_nnz = (PetscInt*)o_nnz;
1008: return(0);
1009: }
1013: /*@C
1014: MeshCreate - Creates a DM object, used to manage data for an unstructured problem
1015: described by a Sieve.
1017: Collective on MPI_Comm
1019: Input Parameter:
1020: . comm - the processors that will share the global vector
1022: Output Parameters:
1023: . mesh - the mesh object
1025: Level: advanced
1027: .seealso MeshDestroy(), MeshCreateGlobalVector(), MeshGetGlobalIndices()
1029: @*/
1030: PetscErrorCode PETSCDM_DLLEXPORT MeshCreate(MPI_Comm comm,Mesh *mesh)
1031: {
1033: Mesh p;
1037: *mesh = PETSC_NULL;
1038: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
1039: DMInitializePackage(PETSC_NULL);
1040: #endif
1042: PetscHeaderCreate(p,_p_Mesh,struct _MeshOps,DA_COOKIE,0,"Mesh",comm,MeshDestroy,0);
1043: p->ops->view = MeshView_Sieve;
1044: p->ops->createglobalvector = MeshCreateGlobalVector;
1045: p->ops->getmatrix = MeshGetMatrix;
1047: PetscObjectChangeTypeName((PetscObject) p, "sieve");
1049: p->globalvector = PETSC_NULL;
1050: *mesh = p;
1051: return(0);
1052: }
1056: /*@C
1057: MeshDestroy - Destroys a mesh.
1059: Collective on Mesh
1061: Input Parameter:
1062: . mesh - the mesh object
1064: Level: advanced
1066: .seealso MeshCreate(), MeshCreateGlobalVector(), MeshGetGlobalIndices()
1068: @*/
1069: PetscErrorCode PETSCDM_DLLEXPORT MeshDestroy(Mesh mesh)
1070: {
1071: PetscErrorCode ierr;
1074: if (--mesh->refct > 0) return(0);
1075: if (mesh->globalvector) {VecDestroy(mesh->globalvector);}
1076: PetscHeaderDestroy(mesh);
1077: return(0);
1078: }
1082: inline void ExpandInterval(ALE::Point interval, PetscInt indices[], PetscInt *indx)
1083: {
1084: for(int i = 0; i < interval.index; i++) {
1085: indices[(*indx)++] = interval.prefix + i;
1086: }
1087: }
1091: inline void ExpandInterval_New(ALE::def::Point interval, PetscInt indices[], PetscInt *indx)
1092: {
1093: for(int i = 0; i < interval.index; i++) {
1094: indices[(*indx)++] = interval.prefix + i;
1095: }
1096: }
1100: PetscErrorCode ExpandIntervals(ALE::Obj<ALE::def::Mesh::bundle_type::IndexArray> intervals, PetscInt *indices)
1101: {
1102: int k = 0;
1105: for(ALE::def::Mesh::bundle_type::IndexArray::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
1106: ExpandInterval_New(*i_itor, indices, &k);
1107: }
1108: return(0);
1109: }
1113: /*
1114: Creates a ghosted vector based upon the global ordering in the bundle.
1115: */
1116: PetscErrorCode MeshCreateVector(Mesh mesh, ALE::Obj<ALE::Two::Mesh> m, Vec *v)
1117: {
1118: ALE::Obj<ALE::Two::Mesh::field_type> field = m->getField("u");
1119: ALE::Two::Mesh::field_type::patch_type patch;
1120: // FIX: Must not include ghosts
1121: PetscInt localSize = field->getSize(ALE::Two::Mesh::field_type::patch_type());
1122: MPI_Comm comm = m->getComm();
1123: PetscMPIInt rank = m->getRank();
1124: PetscInt *ghostIndices = NULL;
1125: PetscInt ghostSize = 0;
1129: #ifdef PARALLEL
1130: ALE::Obj<ALE::PreSieve> globalIndices = bundle->getGlobalIndices();
1131: ALE::Obj<ALE::PreSieve> pointTypes = bundle->getPointTypes();
1132: ALE::Obj<ALE::Point_set> rentedPoints = pointTypes->cone(ALE::Point(rank, ALE::rentedPoint));
1134: for(ALE::Point_set::iterator e_itor = rentedPoints->begin(); e_itor != rentedPoints->end(); e_itor++) {
1135: ALE::Obj<ALE::Point_set> cone = globalIndices->cone(*e_itor);
1137: if (cone->size()) {
1138: ALE::Point interval = *cone->begin();
1140: ghostSize += interval.index;
1141: }
1142: }
1143: #endif
1144: if (ghostSize) {
1145: PetscMalloc(ghostSize * sizeof(PetscInt), &ghostIndices);
1146: }
1147: #ifdef PARALLEL
1148: PetscInt ghostIdx = 0;
1150: for(ALE::Point_set::iterator e_itor = rentedPoints->begin(); e_itor != rentedPoints->end(); e_itor++) {
1151: ALE::Obj<ALE::Point_set> cone = globalIndices->cone(*e_itor);
1153: if (cone->size()) {
1154: ALE::Point interval = *cone->begin();
1156: // Must insert into ghostIndices at the index given by localIndices
1157: // However, I think right now its correct because rentedPoints iterates in the same way in both methods
1158: ExpandInterval(interval, ghostIndices, &ghostIdx);
1159: }
1160: }
1161: #endif
1162: VecCreateGhost(comm, localSize, PETSC_DETERMINE, ghostSize, ghostIndices, v);
1163: if (m->debug) {
1164: PetscInt globalSize, g;
1166: VecGetSize(*v, &globalSize);
1167: PetscPrintf(comm, "Making an ordering over the vertices\n===============================\n");
1168: PetscSynchronizedPrintf(comm, "[%d] global size: %d localSize: %d ghostSize: %d\n", rank, globalSize, localSize, ghostSize);
1169: PetscSynchronizedPrintf(comm, "[%d] ghostIndices:", rank);
1170: for(g = 0; g < ghostSize; g++) {
1171: PetscSynchronizedPrintf(comm, "[%d] %d\n", rank, ghostIndices[g]);
1172: }
1173: PetscSynchronizedPrintf(comm, "\n");
1174: PetscSynchronizedFlush(comm);
1175: }
1176: PetscFree(ghostIndices);
1177: return(0);
1178: }
1182: /*@C
1183: MeshCreateGlobalVector - Creates a vector of the correct size to be gathered into
1184: by the mesh.
1186: Collective on Mesh
1188: Input Parameter:
1189: . mesh - the mesh object
1191: Output Parameters:
1192: . gvec - the global vector
1194: Level: advanced
1196: Notes: Once this has been created you cannot add additional arrays or vectors to be packed.
1198: .seealso MeshDestroy(), MeshCreate(), MeshGetGlobalIndices()
1200: @*/
1201: PetscErrorCode PETSCDM_DLLEXPORT MeshCreateGlobalVector(Mesh mesh,Vec *gvec)
1202: {
1206: /* Turned off caching for this method so that bundle can be reset to make different vectors */
1207: #if 0
1208: if (mesh->globalvector) {
1209: VecDuplicate(mesh->globalvector, gvec);
1210: return(0);
1211: }
1212: #endif
1213: #ifdef __cplusplus
1214: ALE::Obj<ALE::Two::Mesh> m;
1216: MeshGetMesh(mesh, &m);
1217: MeshCreateVector(mesh, m, gvec);
1218: #endif
1219: #if 0
1220: mesh->globalvector = *gvec;
1221: PetscObjectReference((PetscObject) mesh->globalvector);
1222: #endif
1223: return(0);
1224: }
1228: /*@C
1229: MeshGetGlobalIndices - Gets the global indices for all the local entries
1231: Collective on Mesh
1233: Input Parameter:
1234: . mesh - the mesh object
1236: Output Parameters:
1237: . idx - the individual indices for each packed vector/array
1238:
1239: Level: advanced
1241: Notes:
1242: The idx parameters should be freed by the calling routine with PetscFree()
1244: .seealso MeshDestroy(), MeshCreateGlobalVector(), MeshCreate()
1246: @*/
1247: PetscErrorCode PETSCDM_DLLEXPORT MeshGetGlobalIndices(Mesh mesh,PetscInt *idx[])
1248: {
1249: SETERRQ(PETSC_ERR_SUP, "");
1250: }
1252: EXTERN PetscErrorCode assembleFullField(VecScatter, Vec, Vec, InsertMode);
1256: /*@
1257: restrictVector - Insert values from a global vector into a local ghosted vector
1259: Collective on g
1261: Input Parameters:
1262: + g - The global vector
1263: . l - The local vector
1264: - mode - either ADD_VALUES or INSERT_VALUES, where
1265: ADD_VALUES adds values to any existing entries, and
1266: INSERT_VALUES replaces existing entries with new values
1268: Level: beginner
1270: .seealso: MatSetOption()
1271: @*/
1272: PetscErrorCode restrictVector(Vec g, Vec l, InsertMode mode)
1273: {
1274: //VecScatter injection;
1278: //PetscObjectQuery((PetscObject) g, "injection", (PetscObject *) &injection);
1279: //VecScatterBegin(g, l, mode, SCATTER_REVERSE, injection);
1280: //VecScatterEnd(g, l, mode, SCATTER_REVERSE, injection);
1281: if (mode == INSERT_VALUES) {
1282: VecCopy(g, l);
1283: } else {
1284: VecAXPY(l, 1.0, g);
1285: }
1286: return(0);
1287: }
1291: /*@
1292: assembleVectorComplete - Insert values from a local ghosted vector into a global vector
1294: Collective on g
1296: Input Parameters:
1297: + g - The global vector
1298: . l - The local vector
1299: - mode - either ADD_VALUES or INSERT_VALUES, where
1300: ADD_VALUES adds values to any existing entries, and
1301: INSERT_VALUES replaces existing entries with new values
1303: Level: beginner
1305: .seealso: MatSetOption()
1306: @*/
1307: PetscErrorCode assembleVectorComplete(Vec g, Vec l, InsertMode mode)
1308: {
1309: //VecScatter injection;
1313: //PetscObjectQuery((PetscObject) g, "injection", (PetscObject *) &injection);
1314: //VecScatterBegin(l, g, mode, SCATTER_FORWARD, injection);
1315: //VecScatterEnd(l, g, mode, SCATTER_FORWARD, injection);
1316: if (mode == INSERT_VALUES) {
1317: VecCopy(l, g);
1318: } else {
1319: VecAXPY(g, 1.0, l);
1320: }
1321: return(0);
1322: }
1326: /*@
1327: assembleVector - Insert values into a vector
1329: Collective on A
1331: Input Parameters:
1332: + b - the vector
1333: . e - The element number
1334: . v - The values
1335: - mode - either ADD_VALUES or INSERT_VALUES, where
1336: ADD_VALUES adds values to any existing entries, and
1337: INSERT_VALUES replaces existing entries with new values
1339: Level: beginner
1341: .seealso: VecSetOption()
1342: @*/
1343: PetscErrorCode assembleVector(Vec b, PetscInt e, PetscScalar v[], InsertMode mode)
1344: {
1345: Mesh mesh;
1346: ALE::Obj<ALE::Two::Mesh> m;
1347: ALE::Two::Mesh::field_type::patch_type patch;
1348: PetscInt firstElement;
1349: PetscErrorCode ierr;
1352: PetscObjectQuery((PetscObject) b, "mesh", (PetscObject *) &mesh);
1353: MeshGetMesh(mesh, &m);
1354: //firstElement = elementBundle->getLocalSizes()[bundle->getCommRank()];
1355: firstElement = 0;
1356: // Must relate b to field
1357: if (mode == INSERT_VALUES) {
1358: m->getField(std::string("x"))->update(patch, ALE::Two::Mesh::point_type(0, e + firstElement), v);
1359: } else {
1360: m->getField(std::string("x"))->updateAdd(patch, ALE::Two::Mesh::point_type(0, e + firstElement), v);
1361: }
1362: return(0);
1363: }
1367: PetscErrorCode updateOperator(Mat A, ALE::Obj<ALE::Two::Mesh::field_type> field, const ALE::Two::Mesh::point_type& e, PetscScalar array[], InsertMode mode)
1368: {
1369: ALE::Obj<ALE::Two::Mesh::field_type::IndexArray> intervals = field->getIndices("element", e);
1370: static PetscInt indicesSize = 0;
1371: static PetscInt *indices = NULL;
1372: PetscInt numIndices = 0;
1373: PetscErrorCode ierr;
1376: for(ALE::def::Mesh::bundle_type::IndexArray::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
1377: numIndices += (*i_itor).index;
1378: if (0) {
1379: //printf("[%d]interval (%d, %d)\n", mesh->getCommRank(), (*i_itor).prefix, (*i_itor).index);
1380: printf("[%d]interval (%d, %d)\n", 0, (*i_itor).prefix, (*i_itor).index);
1381: }
1382: }
1383: if (indicesSize && (indicesSize != numIndices)) {
1384: PetscFree(indices);
1385: indices = NULL;
1386: }
1387: if (!indices) {
1388: indicesSize = numIndices;
1389: PetscMalloc(indicesSize * sizeof(PetscInt), &indices);
1390: }
1391: ExpandIntervals(intervals, indices);
1392: if (0) {
1393: for(int i = 0; i < numIndices; i++) {
1394: //printf("[%d]indices[%d] = %d\n", mesh->getCommRank(), i, indices[i]);
1395: printf("[%d]indices[%d] = %d\n", 0, i, indices[i]);
1396: }
1397: }
1398: MatSetValues(A, numIndices, indices, numIndices, indices, array, mode);
1399: return(0);
1400: }
1402: PetscErrorCode assembleOperator_New(Mat A, ALE::Obj<ALE::def::Mesh::coordinate_type> field, ALE::Obj<ALE::def::Mesh::sieve_type> orientation, ALE::def::Mesh::sieve_type::point_type e, PetscScalar array[], InsertMode mode)
1403: {
1404: ALE::Obj<ALE::def::Mesh::bundle_type::IndexArray> intervals = field->getOrderedIndices(0, orientation->cone(e));
1405: static PetscInt indicesSize = 0;
1406: static PetscInt *indices = NULL;
1407: PetscInt numIndices = 0;
1408: PetscErrorCode ierr;
1411: for(ALE::def::Mesh::bundle_type::IndexArray::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
1412: numIndices += (*i_itor).index;
1413: if (0) {
1414: //printf("[%d]interval (%d, %d)\n", mesh->getCommRank(), (*i_itor).prefix, (*i_itor).index);
1415: printf("[%d]interval (%d, %d)\n", 0, (*i_itor).prefix, (*i_itor).index);
1416: }
1417: }
1418: if (indicesSize && (indicesSize != numIndices)) {
1419: PetscFree(indices);
1420: indices = NULL;
1421: }
1422: if (!indices) {
1423: indicesSize = numIndices;
1424: PetscMalloc(indicesSize * sizeof(PetscInt), &indices);
1425: }
1426: ExpandIntervals(intervals, indices);
1427: if (0) {
1428: for(int i = 0; i < numIndices; i++) {
1429: //printf("[%d]indices[%d] = %d\n", mesh->getCommRank(), i, indices[i]);
1430: printf("[%d]indices[%d] = %d\n", 0, i, indices[i]);
1431: }
1432: }
1433: MatSetValues(A, numIndices, indices, numIndices, indices, array, mode);
1434: return(0);
1435: }
1439: /*@
1440: assembleMatrix - Insert values into a matrix
1442: Collective on A
1444: Input Parameters:
1445: + A - the matrix
1446: . e - The element number
1447: . v - The values
1448: - mode - either ADD_VALUES or INSERT_VALUES, where
1449: ADD_VALUES adds values to any existing entries, and
1450: INSERT_VALUES replaces existing entries with new values
1452: Level: beginner
1454: .seealso: MatSetOption()
1455: @*/
1456: PetscErrorCode assembleMatrix(Mat A, PetscInt e, PetscScalar v[], InsertMode mode)
1457: {
1458: PetscObjectContainer c;
1459: ALE::def::Mesh *mesh;
1460: int firstElement = 0;
1461: PetscErrorCode ierr;
1464: PetscObjectQuery((PetscObject) A, "mesh", (PetscObject *) &c);
1465: PetscObjectContainerGetPointer(c, (void **) &mesh);
1466: // Need to globalize
1467: // firstElement = elementBundle->getLocalSizes()[bundle->getCommRank()];
1468: try {
1469: assembleOperator_New(A, mesh->getField(), mesh->getOrientation(), ALE::def::Mesh::sieve_type::point_type(0, e + firstElement), v, mode);
1470: } catch (ALE::Exception e) {
1471: std::cout << e << std::endl;
1472: }
1473: return(0);
1474: }