Actual source code: mesh.c

  1: 
 2:  #include src/dm/mesh/meshimpl.h
  3: #include <Distribution.hh>
 4:  #include src/dm/mesh/meshvtk.h
 5:  #include src/dm/mesh/meshpcice.h
 6:  #include src/dm/mesh/meshpylith.h

  8: /* Logging support */
  9: PetscCookie  MESH_COOKIE = 0;
 10: PetscEvent  Mesh_View = 0, Mesh_GetGlobalScatter = 0, Mesh_restrictVector = 0, Mesh_assembleVector = 0,
 11:             Mesh_assembleVectorComplete = 0, Mesh_assembleMatrix = 0, Mesh_updateOperator = 0;

 15: PetscErrorCode MeshView_Sieve_Ascii(const ALE::Obj<ALE::Mesh>& mesh, PetscViewer viewer)
 16: {
 17:   PetscViewerFormat format;
 18:   PetscErrorCode    ierr;

 21:   PetscViewerGetFormat(viewer, &format);
 22:   if (format == PETSC_VIEWER_ASCII_VTK) {
 23:     VTKViewer::writeHeader(viewer);
 24:     VTKViewer::writeVertices(mesh, viewer);
 25:     VTKViewer::writeElements(mesh, viewer);
 26:   } else if (format == PETSC_VIEWER_ASCII_PYLITH) {
 27:     char *filename;
 28:     char  connectFilename[2048];
 29:     char  coordFilename[2048];

 31:     PetscViewerFileGetName(viewer, &filename);
 32:     PetscViewerFileSetMode(viewer, FILE_MODE_WRITE);
 33:     PetscStrcpy(connectFilename, filename);
 34:     PetscStrcat(connectFilename, ".connect");
 35:     PetscViewerFileSetName(viewer, connectFilename);
 36:     ALE::PyLith::Viewer::writeElements(mesh, viewer);
 37:     PetscStrcpy(coordFilename, filename);
 38:     PetscStrcat(coordFilename, ".coord");
 39:     PetscViewerFileSetName(viewer, coordFilename);
 40:     ALE::PyLith::Viewer::writeVertices(mesh, viewer);
 41:     PetscViewerFileSetMode(viewer, FILE_MODE_READ);
 42:     PetscExceptionTry1(PetscViewerFileSetName(viewer, filename), PETSC_ERR_FILE_OPEN);
 43:     if (PetscExceptionValue(ierr)) {
 44:       /* this means that a caller above me has also tryed this exception so I don't handle it here, pass it up */
 45:     } else if (PetscExceptionCaught(ierr, PETSC_ERR_FILE_OPEN)) {
 46:       0;
 47:     }
 48: 
 49:   } else if (format == PETSC_VIEWER_ASCII_PYLITH_LOCAL) {
 50:     PetscViewer connectViewer, coordViewer, splitViewer;
 51:     char       *filename;
 52:     char        localFilename[2048];
 53:     int         rank = mesh->commRank();

 55:     PetscViewerFileGetName(viewer, &filename);

 57:     sprintf(localFilename, "%s.%d.connect", filename, rank);
 58:     PetscViewerCreate(PETSC_COMM_SELF, &connectViewer);
 59:     PetscViewerSetType(connectViewer, PETSC_VIEWER_ASCII);
 60:     PetscViewerSetFormat(connectViewer, PETSC_VIEWER_ASCII_PYLITH);
 61:     PetscViewerFileSetName(connectViewer, localFilename);
 62:     ALE::PyLith::Viewer::writeElementsLocal(mesh, connectViewer);
 63:     PetscViewerDestroy(connectViewer);

 65:     sprintf(localFilename, "%s.%d.coord", filename, rank);
 66:     PetscViewerCreate(PETSC_COMM_SELF, &coordViewer);
 67:     PetscViewerSetType(coordViewer, PETSC_VIEWER_ASCII);
 68:     PetscViewerSetFormat(coordViewer, PETSC_VIEWER_ASCII_PYLITH);
 69:     PetscViewerFileSetName(coordViewer, localFilename);
 70:     ALE::PyLith::Viewer::writeVerticesLocal(mesh, coordViewer);
 71:     PetscViewerDestroy(coordViewer);

 73:     if (!mesh->getSplitSection().isNull()) {
 74:       sprintf(localFilename, "%s.%d.split", filename, rank);
 75:       PetscViewerCreate(PETSC_COMM_SELF, &splitViewer);
 76:       PetscViewerSetType(splitViewer, PETSC_VIEWER_ASCII);
 77:       PetscViewerSetFormat(splitViewer, PETSC_VIEWER_ASCII_PYLITH);
 78:       PetscViewerFileSetName(splitViewer, localFilename);
 79:       ALE::PyLith::Viewer::writeSplitLocal(mesh, mesh->getSplitSection(), splitViewer);
 80:       PetscViewerDestroy(splitViewer);
 81:     }
 82:   } else if (format == PETSC_VIEWER_ASCII_PCICE) {
 83:     char      *filename;
 84:     char       coordFilename[2048];
 85:     PetscTruth isConnect;
 86:     size_t     len;

 88:     PetscViewerFileGetName(viewer, &filename);
 89:     PetscStrlen(filename, &len);
 90:     PetscStrcmp(&(filename[len-5]), ".lcon", &isConnect);
 91:     if (!isConnect) {
 92:       SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid element connectivity filename: %s", filename);
 93:     }
 94:     ALE::PCICE::Viewer::writeElements(mesh, viewer);
 95:     PetscStrncpy(coordFilename, filename, len-5);
 96:     coordFilename[len-5] = '\0';
 97:     PetscStrcat(coordFilename, ".nodes");
 98:     PetscViewerFileSetName(viewer, coordFilename);
 99:     ALE::PCICE::Viewer::writeVertices(mesh, viewer);
100:   } else {
101:     int dim = mesh->getDimension();

103:     PetscViewerASCIIPrintf(viewer, "Mesh in %d dimensions:\n", dim);
104:     for(int d = 0; d <= dim; d++) {
105:       // FIX: Need to globalize
106:       PetscViewerASCIIPrintf(viewer, "  %d %d-cells\n", mesh->getTopology()->depthStratum(d)->size(), d);
107:     }
108:   }
109:   PetscViewerFlush(viewer);
110:   return(0);
111: }

115: PetscErrorCode MeshView_Sieve(const ALE::Obj<ALE::Mesh>& mesh, PetscViewer viewer)
116: {
117:   PetscTruth     iascii, isbinary, isdraw;

121:   PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII, &iascii);
122:   PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_BINARY, &isbinary);
123:   PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_DRAW, &isdraw);

125:   if (iascii){
126:     MeshView_Sieve_Ascii(mesh, viewer);
127:   } else if (isbinary) {
128:     SETERRQ(PETSC_ERR_SUP, "Binary viewer not implemented for Mesh");
129:   } else if (isdraw){
130:     SETERRQ(PETSC_ERR_SUP, "Draw viewer not implemented for Mesh");
131:   } else {
132:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by this mesh object", ((PetscObject)viewer)->type_name);
133:   }
134:   return(0);
135: }

139: PetscErrorCode FieldView_Sieve_Ascii(const ALE::Obj<ALE::Mesh>& mesh, const std::string& name, PetscViewer viewer)
140: {
141:   // state 0: No header has been output
142:   // state 1: Only POINT_DATA has been output
143:   // state 2: Only CELL_DATA has been output
144:   // state 3: Output both, POINT_DATA last
145:   // state 4: Output both, CELL_DATA last
146:   PetscViewerFormat format;
147:   PetscErrorCode    ierr;

150:   PetscViewerGetFormat(viewer, &format);
151:   if (format == PETSC_VIEWER_ASCII_VTK || format == PETSC_VIEWER_ASCII_VTK_CELL) {
152:     static PetscInt   stateId     = -1;
153:     PetscInt          doOutput    = 0;
154:     PetscInt          outputState = 0;
155:     PetscInt          fiberDim    = 0;
156:     PetscTruth        hasState;

158:     if (stateId < 0) {
159:       PetscObjectComposedDataRegister(&stateId);
160:       PetscObjectComposedDataSetInt((PetscObject) viewer, stateId, 0);
161:     }
162:     PetscObjectComposedDataGetInt((PetscObject) viewer, stateId, outputState, hasState);
163:     if (format == PETSC_VIEWER_ASCII_VTK) {
164:       if (outputState == 0) {
165:         outputState = 1;
166:         doOutput = 1;
167:       } else if (outputState == 1) {
168:         doOutput = 0;
169:       } else if (outputState == 2) {
170:         outputState = 3;
171:         doOutput = 1;
172:       } else if (outputState == 3) {
173:         doOutput = 0;
174:       } else if (outputState == 4) {
175:         SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Tried to output POINT_DATA again after intervening CELL_DATA");
176:       }
177:       typedef ALE::New::Numbering<ALE::Mesh::topology_type> numbering_type;
178:       ALE::Obj<ALE::Mesh::section_type>   field     = mesh->getSection(name);
179:       ALE::Obj<numbering_type>            numbering = new numbering_type(mesh->getTopologyNew(), "depth", 0);

181:       numbering->construct();
182:       if (doOutput) {
183:         ALE::Mesh::section_type::patch_type patch = mesh->getTopologyNew()->getPatches().begin()->first;

185:         fiberDim = field->getAtlas()->size(patch, *mesh->getTopologyNew()->depthStratum(patch, 0)->begin());
186:         PetscViewerASCIIPrintf(viewer, "POINT_DATA %d\n", numbering->getGlobalSize());
187:       }
188:       VTKViewer::writeField(mesh, field, name, fiberDim, numbering, viewer);
189:     } else {
190:       if (outputState == 0) {
191:         outputState = 2;
192:         doOutput = 1;
193:       } else if (outputState == 1) {
194:         outputState = 4;
195:         doOutput = 1;
196:       } else if (outputState == 2) {
197:         doOutput = 0;
198:       } else if (outputState == 3) {
199:         SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Tried to output CELL_DATA again after intervening POINT_DATA");
200:       } else if (outputState == 4) {
201:         doOutput = 0;
202:       }
203:       typedef ALE::New::Numbering<ALE::Mesh::topology_type> numbering_type;
204:       ALE::Obj<ALE::Mesh::section_type>   field     = mesh->getSection(name);
205:       ALE::Obj<numbering_type>            numbering = new numbering_type(mesh->getTopologyNew(), "height", 0);

207:       numbering->construct();
208:       if (doOutput) {
209:         ALE::Mesh::section_type::patch_type patch = mesh->getTopologyNew()->getPatches().begin()->first;

211:         fiberDim = field->getAtlas()->size(patch, *mesh->getTopologyNew()->heightStratum(patch, 0)->begin());
212:         PetscViewerASCIIPrintf(viewer, "CELL_DATA %d\n", numbering->getGlobalSize());
213:       }
214:       VTKViewer::writeField(mesh, field, name, fiberDim, numbering, viewer);
215:     }
216:     PetscObjectComposedDataSetInt((PetscObject) viewer, stateId, outputState);
217:   } else {
218:   }
219:   return(0);
220: }

224: PetscErrorCode FieldView_Sieve(const ALE::Obj<ALE::Mesh>& mesh, const std::string& name, PetscViewer viewer)
225: {
226:   PetscTruth     iascii, isbinary, isdraw;

230:   PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII, &iascii);
231:   PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_BINARY, &isbinary);
232:   PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_DRAW, &isdraw);

234:   if (iascii){
235:     FieldView_Sieve_Ascii(mesh, name, viewer);
236:   } else if (isbinary) {
237:     SETERRQ(PETSC_ERR_SUP, "Binary viewer not implemented for Field");
238:   } else if (isdraw){
239:     SETERRQ(PETSC_ERR_SUP, "Draw viewer not implemented for Field");
240:   } else {
241:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by this field object", ((PetscObject)viewer)->type_name);
242:   }
243:   return(0);
244: }

248: /*@C
249:    MeshView - Views a Mesh object. 

251:    Collective on Mesh

253:    Input Parameters:
254: +  mesh - the mesh
255: -  viewer - an optional visualization context

257:    Notes:
258:    The available visualization contexts include
259: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
260: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
261:          output where only the first processor opens
262:          the file.  All other processors send their 
263:          data to the first processor to print. 

265:    You can change the format the mesh is printed using the 
266:    option PetscViewerSetFormat().

268:    The user can open alternative visualization contexts with
269: +    PetscViewerASCIIOpen() - Outputs mesh to a specified file
270: .    PetscViewerBinaryOpen() - Outputs mesh in binary to a
271:          specified file; corresponding input uses MeshLoad()
272: .    PetscViewerDrawOpen() - Outputs mesh to an X window display

274:    The user can call PetscViewerSetFormat() to specify the output
275:    format of ASCII printed objects (when using PETSC_VIEWER_STDOUT_SELF,
276:    PETSC_VIEWER_STDOUT_WORLD and PetscViewerASCIIOpen).  Available formats include
277: +    PETSC_VIEWER_ASCII_DEFAULT - default, prints mesh information
278: -    PETSC_VIEWER_ASCII_VTK - outputs a VTK file describing the mesh

280:    Level: beginner

282:    Concepts: mesh^printing
283:    Concepts: mesh^saving to disk

285: .seealso: PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerBinaryOpen(),
286:           MeshLoad(), PetscViewerCreate()
287: @*/
288: PetscErrorCode  MeshView(Mesh mesh, PetscViewer viewer)
289: {

295:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(mesh->comm);

300:   (*mesh->ops->view)(mesh->m, viewer);
302:   return(0);
303: }

307: /*@C
308:     MeshLoad - Create a mesh topology from the saved data in a viewer.

310:     Collective on Viewer

312:     Input Parameter:
313: .   viewer - The viewer containing the data

315:     Output Parameters:
316: .   mesh - the mesh object

318:     Level: advanced

320: .seealso MeshView()

322: @*/
323: PetscErrorCode  MeshLoad(PetscViewer viewer, Mesh *mesh)
324: {
325:   SETERRQ(PETSC_ERR_SUP, "");
326: }

330: /*@C
331:     MeshGetMesh - Gets the internal mesh object

333:     Not collective

335:     Input Parameter:
336: .    mesh - the mesh object

338:     Output Parameter:
339: .    m - the internal mesh object
340:  
341:     Level: advanced

343: .seealso MeshCreate(), MeshSetMesh()

345: @*/
346: PetscErrorCode  MeshGetMesh(Mesh mesh, ALE::Obj<ALE::Mesh> *m)
347: {
350:   if (m) {
352:     *m = mesh->m;
353:   }
354:   return(0);
355: }

359: /*@C
360:     MeshSetMesh - Sets the internal mesh object

362:     Not collective

364:     Input Parameters:
365: +    mesh - the mesh object
366: -    boundary - the internal mesh object
367:  
368:     Level: advanced

370: .seealso MeshCreate(), MeshGetMesh()

372: @*/
373: PetscErrorCode  MeshSetMesh(Mesh mesh, ALE::Obj<ALE::Mesh> m)
374: {
377:   mesh->m = m;
378:   return(0);
379: }

383: /*@C
384:     MeshGetMatrix - Creates a matrix with the correct parallel layout required for 
385:       computing the Jacobian on a function defined using the informatin in Mesh.

387:     Collective on Mesh

389:     Input Parameter:
390: +   mesh - the mesh object
391: -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, MATSEQSBAIJ, MATMPISBAIJ,
392:             or any type which inherits from one of these (such as MATAIJ, MATLUSOL, etc.).

394:     Output Parameters:
395: .   J  - matrix with the correct nonzero preallocation
396:         (obviously without the correct Jacobian values)

398:     Level: advanced

400:     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you 
401:        do not need to do it yourself.

403: .seealso ISColoringView(), ISColoringGetIS(), MatFDColoringCreate(), DASetBlockFills()

405: @*/
406: PetscErrorCode  MeshGetMatrix(Mesh mesh, MatType mtype,Mat *J)
407: {
408:   ALE::Obj<ALE::Mesh> m;
409: #if 0
410:   ISLocalToGlobalMapping lmap;
411:   PetscInt              *globals,rstart,i;
412: #endif
413:   PetscInt               localSize;
414:   PetscErrorCode         ierr;

417:   MeshGetMesh(mesh, &m);
418:   //localSize = m->getSection("u")->getGlobalOrder()->getSize(ALE::Mesh::field_type::patch_type());

420:   MatCreate(mesh->comm,J);
421:   MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);
422:   MatSetType(*J,mtype);
423:   MatSetBlockSize(*J,1);
424:   MatSeqAIJSetPreallocation(*J,mesh->d_nz,mesh->d_nnz);
425:   MatMPIAIJSetPreallocation(*J,mesh->d_nz,mesh->d_nnz,mesh->o_nz,mesh->o_nnz);
426:   MatSeqBAIJSetPreallocation(*J,mesh->bs,mesh->d_nz,mesh->d_nnz);
427:   MatMPIBAIJSetPreallocation(*J,mesh->bs,mesh->d_nz,mesh->d_nnz,mesh->o_nz,mesh->o_nnz);

429: #if 0
430:   PetscMalloc((mesh->n+mesh->Nghosts+1)*sizeof(PetscInt),&globals);
431:   MatGetOwnershipRange(*J,&rstart,PETSC_NULL);
432:   for (i=0; i<mesh->n; i++) {
433:     globals[i] = rstart + i;
434:   }
435:   PetscMemcpy(globals+mesh->n,mesh->ghosts,mesh->Nghosts*sizeof(PetscInt));
436:   ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,mesh->n+mesh->Nghosts,globals,&lmap);
437:   PetscFree(globals);
438:   MatSetLocalToGlobalMapping(*J,lmap);
439:   ISLocalToGlobalMappingDestroy(lmap);
440: #endif
441:   return(0);
442: }

446: /*@C
447:     MeshSetGhosts - Sets the global indices of other processes elements that will
448:       be ghosts on this process

450:     Not Collective

452:     Input Parameters:
453: +    mesh - the Mesh object
454: .    bs - block size
455: .    nlocal - number of local (non-ghost) entries
456: .    Nghosts - number of ghosts on this process
457: -    ghosts - indices of all the ghost points

459:     Level: advanced

461: .seealso MeshDestroy(), MeshCreateGlobalVector(), MeshGetGlobalIndices()

463: @*/
464: PetscErrorCode  MeshSetGhosts(Mesh mesh,PetscInt bs,PetscInt nlocal,PetscInt Nghosts,const PetscInt ghosts[])
465: {

470:   PetscFree(mesh->ghosts);
471:   PetscMalloc((1+Nghosts)*sizeof(PetscInt),&mesh->ghosts);
472:   PetscMemcpy(mesh->ghosts,ghosts,Nghosts*sizeof(PetscInt));
473:   mesh->bs      = bs;
474:   mesh->n       = nlocal;
475:   mesh->Nghosts = Nghosts;
476:   return(0);
477: }

481: /*@C
482:     MeshSetPreallocation - sets the matrix memory preallocation for matrices computed by Mesh

484:     Not Collective

486:     Input Parameters:
487: +    mesh - the Mesh object
488: .    d_nz - maximum number of nonzeros in any row of diagonal block
489: .    d_nnz - number of nonzeros in each row of diagonal block
490: .    o_nz - maximum number of nonzeros in any row of off-diagonal block
491: .    o_nnz - number of nonzeros in each row of off-diagonal block


494:     Level: advanced

496: .seealso MeshDestroy(), MeshCreateGlobalVector(), MeshGetGlobalIndices(), MatMPIAIJSetPreallocation(),
497:          MatMPIBAIJSetPreallocation()

499: @*/
500: PetscErrorCode  MeshSetPreallocation(Mesh mesh,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
501: {
504:   mesh->d_nz  = d_nz;
505:   mesh->d_nnz = (PetscInt*)d_nnz;
506:   mesh->o_nz  = o_nz;
507:   mesh->o_nnz = (PetscInt*)o_nnz;
508:   return(0);
509: }

513: /*@C
514:     MeshCreate - Creates a DM object, used to manage data for an unstructured problem
515:     described by a Sieve.

517:     Collective on MPI_Comm

519:     Input Parameter:
520: .   comm - the processors that will share the global vector

522:     Output Parameters:
523: .   mesh - the mesh object

525:     Level: advanced

527: .seealso MeshDestroy(), MeshCreateGlobalVector(), MeshGetGlobalIndices()

529: @*/
530: PetscErrorCode  MeshCreate(MPI_Comm comm,Mesh *mesh)
531: {
533:   Mesh         p;

537:   *mesh = PETSC_NULL;
538: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
539:   DMInitializePackage(PETSC_NULL);
540: #endif

542:   PetscHeaderCreate(p,_p_Mesh,struct _MeshOps,MESH_COOKIE,0,"Mesh",comm,MeshDestroy,0);
543:   p->ops->view               = MeshView_Sieve;
544:   p->ops->createglobalvector = MeshCreateGlobalVector;
545:   p->ops->getmatrix          = MeshGetMatrix;

547:   PetscObjectChangeTypeName((PetscObject) p, "sieve");

549:   p->m            = PETSC_NULL;
550:   p->globalvector = PETSC_NULL;
551:   *mesh = p;
552:   return(0);
553: }

557: /*@C
558:     MeshDestroy - Destroys a mesh.

560:     Collective on Mesh

562:     Input Parameter:
563: .   mesh - the mesh object

565:     Level: advanced

567: .seealso MeshCreate(), MeshCreateGlobalVector(), MeshGetGlobalIndices()

569: @*/
570: PetscErrorCode  MeshDestroy(Mesh mesh)
571: {
572:   PetscErrorCode     ierr;

575:   if (--mesh->refct > 0) return(0);
576:   if (mesh->globalvector) {VecDestroy(mesh->globalvector);}
577:   PetscHeaderDestroy(mesh);
578:   return(0);
579: }

583: inline void ExpandInterval(const ALE::Point& interval, int indices[], int& indx)
584: {
585:   const int end = interval.prefix + interval.index;
586:   for(int i = interval.prefix; i < end; i++) {
587:     indices[indx++] = i;
588:   }
589: }

593: inline void ExpandInterval_New(ALE::Point interval, PetscInt indices[], PetscInt *indx)
594: {
595:   for(int i = 0; i < interval.index; i++) {
596:     indices[(*indx)++] = interval.prefix + i;
597:   }
598: }

602: PetscErrorCode ExpandIntervals(ALE::Obj<ALE::Mesh::atlas_type::IndexArray> intervals, PetscInt *indices)
603: {
604:   int k = 0;

607:   for(ALE::Mesh::atlas_type::IndexArray::iterator i_itor = intervals->begin(); i_itor != intervals->end(); i_itor++) {
608:     ExpandInterval_New(*i_itor, indices, &k);
609:   }
610:   return(0);
611: }

615: /*
616:   Creates a ghosted vector based upon the global ordering in the bundle.
617: */
618: PetscErrorCode MeshCreateVector(ALE::Obj<ALE::Mesh> m, Vec *v)
619: {
620:   // FIX: Must not include ghosts
621:   PetscInt       localSize = 0;
622:   MPI_Comm       comm = m->comm();
623:   PetscMPIInt    rank = m->commRank();
624:   PetscInt      *ghostIndices = NULL;
625:   PetscInt       ghostSize = 0;

629: #ifdef PARALLEL
630:   ALE::Obj<ALE::PreSieve> globalIndices = bundle->getGlobalIndices();
631:   ALE::Obj<ALE::PreSieve> pointTypes = bundle->getPointTypes();
632:   ALE::Obj<ALE::Point_set> rentedPoints = pointTypes->cone(ALE::Point(rank, ALE::rentedPoint));

634:   for(ALE::Point_set::iterator e_itor = rentedPoints->begin(); e_itor != rentedPoints->end(); e_itor++) {
635:     ALE::Obj<ALE::Point_set> cone = globalIndices->cone(*e_itor);

637:     if (cone->size()) {
638:       ALE::Point interval = *cone->begin();

640:       ghostSize += interval.index;
641:     }
642:   }
643: #endif
644:   if (ghostSize) {
645:     PetscMalloc(ghostSize * sizeof(PetscInt), &ghostIndices);
646:   }
647: #ifdef PARALLEL
648:   PetscInt ghostIdx = 0;

650:   for(ALE::Point_set::iterator e_itor = rentedPoints->begin(); e_itor != rentedPoints->end(); e_itor++) {
651:     ALE::Obj<ALE::Point_set> cone = globalIndices->cone(*e_itor);

653:     if (cone->size()) {
654:       ALE::Point interval = *cone->begin();

656:       // Must insert into ghostIndices at the index given by localIndices
657:       //   However, I think right now its correct because rentedPoints iterates in the same way in both methods
658:       ExpandInterval(interval, ghostIndices, &ghostIdx);
659:     }
660:   }
661: #endif
662:   VecCreateGhost(comm, localSize, PETSC_DETERMINE, ghostSize, ghostIndices, v);
663:   if (m->debug) {
664:     PetscInt globalSize, g;

666:     VecGetSize(*v, &globalSize);
667:     PetscPrintf(comm, "Making an ordering over the vertices\n===============================\n");
668:     PetscSynchronizedPrintf(comm, "[%d]  global size: %d localSize: %d ghostSize: %d\n", rank, globalSize, localSize, ghostSize);
669:     PetscSynchronizedPrintf(comm, "[%d]  ghostIndices:", rank);
670:     for(g = 0; g < ghostSize; g++) {
671:       PetscSynchronizedPrintf(comm, "[%d] %d\n", rank, ghostIndices[g]);
672:     }
673:     PetscSynchronizedPrintf(comm, "\n");
674:     PetscSynchronizedFlush(comm);
675:   }
676:   PetscFree(ghostIndices);
677:   return(0);
678: }

682: /*@C
683:     MeshCreateGlobalVector - Creates a vector of the correct size to be gathered into 
684:         by the mesh.

686:     Collective on Mesh

688:     Input Parameter:
689: .    mesh - the mesh object

691:     Output Parameters:
692: .   gvec - the global vector

694:     Level: advanced

696:     Notes: Once this has been created you cannot add additional arrays or vectors to be packed.

698: .seealso MeshDestroy(), MeshCreate(), MeshGetGlobalIndices()

700: @*/
701: PetscErrorCode  MeshCreateGlobalVector(Mesh mesh,Vec *gvec)
702: {

706:   /* Turned off caching for this method so that bundle can be reset to make different vectors */
707: #if 0
708:   if (mesh->globalvector) {
709:     VecDuplicate(mesh->globalvector, gvec);
710:     return(0);
711:   }
712: #endif
713: #ifdef __cplusplus
714:   ALE::Obj<ALE::Mesh> m;

716:   MeshGetMesh(mesh, &m);
717:   MeshCreateVector(m, gvec);
718: #endif
719: #if 0
720:   mesh->globalvector = *gvec;
721:   PetscObjectReference((PetscObject) mesh->globalvector);
722: #endif
723:   return(0);
724: }

728: /*@C
729:     MeshGetGlobalIndices - Gets the global indices for all the local entries

731:     Collective on Mesh

733:     Input Parameter:
734: .    mesh - the mesh object

736:     Output Parameters:
737: .    idx - the individual indices for each packed vector/array
738:  
739:     Level: advanced

741:     Notes:
742:        The idx parameters should be freed by the calling routine with PetscFree()

744: .seealso MeshDestroy(), MeshCreateGlobalVector(), MeshCreate()

746: @*/
747: PetscErrorCode  MeshGetGlobalIndices(Mesh mesh,PetscInt *idx[])
748: {
749:   SETERRQ(PETSC_ERR_SUP, "");
750: }

754: PetscErrorCode  MeshGetGlobalScatter(ALE::Mesh *mesh, const char fieldName[], Vec g, VecScatter *scatter)
755: {
756:   typedef ALE::Mesh::atlas_type::index_type index_type;

761:   const ALE::Obj<ALE::Mesh::section_type>&   field       = mesh->getSection(std::string(fieldName));
762:   const ALE::Obj<ALE::Mesh::numbering_type>& globalOrder = mesh->getGlobalOrder(fieldName);
763:   const ALE::Mesh::section_type::patch_type  patch       = 0;
764:   const ALE::Mesh::atlas_type::chart_type&   chart       = field->getAtlas()->getChart(patch);
765:   int                                        localSize   = field->getAtlas()->size(patch);
766:   int *localIndices, *globalIndices;
767:   int  localIndx = 0, globalIndx = 0;
768:   Vec  localVec;
769:   IS   localIS, globalIS;

771:   // Loop over all local points
772:   PetscMalloc(localSize*sizeof(int), &localIndices);
773:   PetscMalloc(localSize*sizeof(int), &globalIndices);
774:   for(ALE::Mesh::atlas_type::chart_type::const_iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
775:     // Map local indices to global indices
776:     ExpandInterval(c_iter->second, localIndices, localIndx);
777:     ExpandInterval(index_type(globalOrder->getIndex(c_iter->first), c_iter->second.index), globalIndices, globalIndx);
778:   }
779:   if (localIndx  != localSize) SETERRQ2(PETSC_ERR_ARG_SIZ, "Invalid number of local indices %d, should be %d", localIndx, localSize);
780:   if (globalIndx != localSize) SETERRQ2(PETSC_ERR_ARG_SIZ, "Invalid number of global indices %d, should be %d", globalIndx, localSize);
781:   ISCreateGeneral(PETSC_COMM_SELF, localSize, localIndices,  &localIS);
782:   ISCreateGeneral(PETSC_COMM_SELF, localSize, globalIndices, &globalIS);
783:   PetscFree(localIndices);
784:   PetscFree(globalIndices);
785:   VecCreateSeqWithArray(PETSC_COMM_SELF, localSize, field->restrict(patch), &localVec);
786:   VecScatterCreate(localVec, localIS, g, globalIS, scatter);
787:   ISDestroy(globalIS);
788:   ISDestroy(localIS);
790:   return(0);
791: }

793: EXTERN PetscErrorCode assembleFullField(VecScatter, Vec, Vec, InsertMode);

797: /*@C
798:   restrictVector - Insert values from a global vector into a local ghosted vector

800:   Collective on g

802:   Input Parameters:
803: + g - The global vector
804: . l - The local vector
805: - mode - either ADD_VALUES or INSERT_VALUES, where
806:    ADD_VALUES adds values to any existing entries, and
807:    INSERT_VALUES replaces existing entries with new values

809:    Level: beginner

811: .seealso: MatSetOption()
812: @*/
813: PetscErrorCode restrictVector(Vec g, Vec l, InsertMode mode)
814: {
815:   VecScatter     injection;

820:   PetscObjectQuery((PetscObject) g, "injection", (PetscObject *) &injection);
821:   if (injection) {
822:     VecScatterBegin(g, l, mode, SCATTER_REVERSE, injection);
823:     VecScatterEnd(g, l, mode, SCATTER_REVERSE, injection);
824:   } else {
825:     if (mode == INSERT_VALUES) {
826:       VecCopy(g, l);
827:     } else {
828:       VecAXPY(l, 1.0, g);
829:     }
830:   }
832:   return(0);
833: }

837: /*@C
838:   assembleVectorComplete - Insert values from a local ghosted vector into a global vector

840:   Collective on g

842:   Input Parameters:
843: + g - The global vector
844: . l - The local vector
845: - mode - either ADD_VALUES or INSERT_VALUES, where
846:    ADD_VALUES adds values to any existing entries, and
847:    INSERT_VALUES replaces existing entries with new values

849:    Level: beginner

851: .seealso: MatSetOption()
852: @*/
853: PetscErrorCode assembleVectorComplete(Vec g, Vec l, InsertMode mode)
854: {
855:   VecScatter     injection;

860:   PetscObjectQuery((PetscObject) g, "injection", (PetscObject *) &injection);
861:   if (injection) {
862:     VecScatterBegin(l, g, mode, SCATTER_FORWARD, injection);
863:     VecScatterEnd(l, g, mode, SCATTER_FORWARD, injection);
864:   } else {
865:     if (mode == INSERT_VALUES) {
866:       VecCopy(l, g);
867:     } else {
868:       VecAXPY(g, 1.0, l);
869:     }
870:   }
872:   return(0);
873: }

877: /*@C
878:   assembleVector - Insert values into a vector

880:   Collective on A

882:   Input Parameters:
883: + b - the vector
884: . e - The element number
885: . v - The values
886: - mode - either ADD_VALUES or INSERT_VALUES, where
887:    ADD_VALUES adds values to any existing entries, and
888:    INSERT_VALUES replaces existing entries with new values

890:    Level: beginner

892: .seealso: VecSetOption()
893: @*/
894: PetscErrorCode assembleVector(Vec b, PetscInt e, PetscScalar v[], InsertMode mode)
895: {
896:   Mesh                     mesh;
897:   ALE::Obj<ALE::Mesh> m;
898:   ALE::Mesh::section_type::patch_type patch;
899:   PetscInt                 firstElement;
900:   PetscErrorCode           ierr;

904:   PetscObjectQuery((PetscObject) b, "mesh", (PetscObject *) &mesh);
905:   MeshGetMesh(mesh, &m);
906:   //firstElement = elementBundle->getLocalSizes()[bundle->getCommRank()];
907:   firstElement = 0;
908:   // Must relate b to field
909:   if (mode == INSERT_VALUES) {
910:     m->getSection(std::string("x"))->update(patch, ALE::Mesh::point_type(e + firstElement), v);
911:   } else {
912:     m->getSection(std::string("x"))->updateAdd(patch, ALE::Mesh::point_type(e + firstElement), v);
913:   }
915:   return(0);
916: }

920: PetscErrorCode updateOperator(Mat A, const ALE::Obj<ALE::Mesh::atlas_type>& atlas, const ALE::Obj<ALE::Mesh::numbering_type>& globalOrder, const ALE::Mesh::point_type& e, PetscScalar array[], InsertMode mode)
921: {
922:   ALE::Mesh::section_type::patch_type patch = 0;
923:   static PetscInt  indicesSize = 0;
924:   static PetscInt *indices = NULL;
925:   PetscInt         numIndices = 0;
926:   PetscErrorCode   ierr;

929:   const ALE::Obj<ALE::Mesh::atlas_type::IndexArray> intervals = atlas->getIndices(patch, e, globalOrder);

932:   if (atlas->debug()) {printf("[%d]mat for element %d\n", atlas->commRank(), e);}
933:   for(ALE::Mesh::atlas_type::IndexArray::iterator i_iter = intervals->begin(); i_iter != intervals->end(); ++i_iter) {
934:     numIndices += i_iter->index;
935:     if (atlas->debug()) {
936:       printf("[%d]mat interval (%d, %d)\n", atlas->commRank(), i_iter->prefix, i_iter->index);
937:     }
938:   }
939:   if (indicesSize && (indicesSize != numIndices)) {
940:     PetscFree(indices);
941:     indices = NULL;
942:   }
943:   if (!indices) {
944:     indicesSize = numIndices;
945:     PetscMalloc(indicesSize * sizeof(PetscInt), &indices);
946:   }
947:   ExpandIntervals(intervals, indices);
948:   if (atlas->debug()) {
949:     for(int i = 0; i < numIndices; i++) {
950:       printf("[%d]mat indices[%d] = %d\n", atlas->commRank(), i, indices[i]);
951:     }
952:     for(int i = 0; i < numIndices; i++) {
953:       printf("[%d]", atlas->commRank());
954:       for(int j = 0; j < numIndices; j++) {
955:         printf(" %g", array[i*numIndices+j]);
956:       }
957:       printf("\n");
958:     }
959:   }
960:   MatSetValues(A, numIndices, indices, numIndices, indices, array, mode);
962:   return(0);
963: }

967: /*@C
968:   assembleMatrix - Insert values into a matrix

970:   Collective on A

972:   Input Parameters:
973: + A - the matrix
974: . e - The element number
975: . v - The values
976: - mode - either ADD_VALUES or INSERT_VALUES, where
977:    ADD_VALUES adds values to any existing entries, and
978:    INSERT_VALUES replaces existing entries with new values

980:    Level: beginner

982: .seealso: MatSetOption()
983: @*/
984: PetscErrorCode assembleMatrix(Mat A, PetscInt e, PetscScalar v[], InsertMode mode)
985: {
986:   PetscObjectContainer c;
987:   ALE::Mesh           *mesh;
988:   PetscErrorCode       ierr;

992:   PetscObjectQuery((PetscObject) A, "mesh", (PetscObject *) &c);
993:   PetscObjectContainerGetPointer(c, (void **) &mesh);
994:   try {
995:     updateOperator(A, mesh->getSection("displacement")->getAtlas(), mesh->getGlobalOrder("displacement"), mesh->getLocalNumbering(mesh->getTopologyNew()->depth())->getPoint(e), v, mode);
996:   } catch (ALE::Exception e) {
997:     std::cout << e.msg() << std::endl;
998:   }
1000:   return(0);
1001: }

1005: PetscErrorCode preallocateMatrix(ALE::Mesh *mesh, const ALE::Obj<ALE::Mesh::section_type>& field, const ALE::Obj<ALE::Mesh::numbering_type>& globalOrder, Mat A)
1006: {
1007:   const ALE::Obj<ALE::Mesh::sieve_type>     adjGraph    = new ALE::Mesh::sieve_type(mesh->comm(), mesh->debug);
1008:   const ALE::Obj<ALE::Mesh::topology_type>  adjTopology = new ALE::Mesh::topology_type(mesh->comm(), mesh->debug);
1009:   const ALE::Mesh::section_type::patch_type patch       = 0;
1010:   const ALE::Obj<ALE::Mesh::sieve_type>&    sieve       = mesh->getTopologyNew()->getPatch(patch);
1011:   PetscInt       numLocalRows, firstRow, lastRow;
1012:   PetscInt      *dnz, *onz;

1016:   adjTopology->setPatch(patch, adjGraph);
1017:   MatGetLocalSize(A, &numLocalRows, PETSC_NULL);
1018:   MatGetOwnershipRange(A, &firstRow, &lastRow);
1019:   PetscMalloc2(numLocalRows, PetscInt, &dnz, numLocalRows, PetscInt, &onz);
1020:   /* Create local adjacency graph */
1021:   const ALE::Mesh::atlas_type::chart_type& chart = field->getAtlas()->getChart(patch);

1023:   for(ALE::Mesh::atlas_type::chart_type::const_iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
1024:     const ALE::Mesh::atlas_type::point_type& point = c_iter->first;

1026:     adjGraph->addCone(sieve->cone(sieve->support(point)), point);
1027:   }
1028:   /* Distribute adjacency graph */
1029: #if 1
1030:   const ALE::Obj<ALE::Mesh::send_overlap_type>& vertexSendOverlap = mesh->getVertexSendOverlap();
1031:   const ALE::Obj<ALE::Mesh::recv_overlap_type>& vertexRecvOverlap = mesh->getVertexRecvOverlap();

1033:   //ALE::New::Distribution<ALE::Mesh::topology_type>::coneCompletion(vertexSendOverlap, vertexRecvOverlap, adjTopology);
1034: #else
1035:   const ALE::Obj<ALE::Mesh::send_overlap_type>&             vertexSendOverlap = mesh->getVertexSendOverlap();
1036:   const ALE::Obj<ALE::Mesh::recv_overlap_type>&             vertexRecvOverlap = mesh->getVertexRecvOverlap();
1037:   const ALE::Obj<ALE::Mesh::numbering_type>&                vNumbering  = mesh->getLocalNumbering(0);
1038:   const ALE::Obj<ALE::Mesh::topology_type::label_sequence>& vertices    = field->getAtlas()->getTopology()->depthStratum(patch, 0);
1039:   int                                                       numVertices = vertices->size();
1040:   short                                                    *assignment  = new short[numVertices];
1041:   int                                                       rank        = mesh->commRank();

1043:   for(ALE::Mesh::topology_type::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
1044:     // FIX if (vertexOverlap->capContains(*v_iter)) {
1045:     if ((rank && vertexRecvOverlap->baseContains(*v_iter)) || (!rank && vertexSendOverlap->capContains(*v_iter))) {
1046:       int minRank = field->commSize();

1048:       if (!rank) {
1049:         const Obj<ALE::Mesh::send_overlap_type::traits::supportSequence>& sendPatches = vertexSendOverlap->support(*v_iter);

1051:         for(ALE::Mesh::send_overlap_type::traits::supportSequence::iterator p_iter = sendPatches->begin(); p_iter != sendPatches->end(); ++p_iter) {
1052:           if (*p_iter < minRank) minRank = *p_iter;
1053:         }
1054:       } else {
1055:         const Obj<ALE::Mesh::recv_overlap_type::traits::coneSequence>& recvPatches = vertexRecvOverlap->cone(*v_iter);

1057:         for(ALE::Mesh::recv_overlap_type::traits::coneSequence::iterator p_iter = recvPatches->begin(); p_iter != recvPatches->end(); ++p_iter) {
1058:           if (*p_iter < minRank) minRank = *p_iter;
1059:         }
1060:       }
1061:       if (minRank < rank) {
1062:         assignment[vNumbering->getIndex(*v_iter)] = minRank;
1063:       } else {
1064:         assignment[vNumbering->getIndex(*v_iter)] = rank;
1065:       }
1066:     } else {
1067:       assignment[vNumbering->getIndex(*v_iter)] = rank;
1068:     }
1069:   }
1070:   ALE::New::Distribution<ALE::Mesh::topology_type>::scatterTopology(serialTopology, numVertices, assignment, parallelTopology);
1071: #endif
1072:   /* Read out adjacency graph */
1073:   const ALE::Obj<ALE::Mesh::sieve_type> graph = adjTopology->getPatch(patch);

1075:   PetscMemzero(dnz, numLocalRows * sizeof(PetscInt));
1076:   PetscMemzero(onz, numLocalRows * sizeof(PetscInt));
1077:   for(ALE::Mesh::atlas_type::chart_type::const_iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
1078:     const ALE::Mesh::atlas_type::point_type& point = c_iter->first;

1080:     if (globalOrder->isLocal(point)) {
1081:       const Obj<ALE::Mesh::sieve_type::traits::coneSequence>& adj   = graph->cone(point);
1082:       const int                                               row   = globalOrder->getIndex(point);
1083:       const int                                               rSize = c_iter->second.index;

1085:       for(ALE::Mesh::sieve_type::traits::coneSequence::iterator v_iter = adj->begin(); v_iter != adj->end(); ++v_iter) {
1086:         const ALE::Mesh::atlas_type::point_type& neighbor = *v_iter;
1087:         const int                                col      = globalOrder->getIndex(neighbor);
1088:         const int                                cSize    = field->getAtlas()->getFiberDimension(patch, neighbor);
1089: 
1090:         if (col >= firstRow && col < lastRow) {
1091:           for(int r = 0; r < rSize; ++r) {dnz[row - firstRow + r] += cSize;}
1092:         } else {
1093:           for(int r = 0; r < rSize; ++r) {onz[row - firstRow + r] += cSize;}
1094:         }
1095:       }
1096:     }
1097:   }
1098:   int rank = mesh->commRank();
1099:   for(int r = 0; r < numLocalRows; r++) {
1100:     std::cout << "["<<rank<<"]: dnz["<<r<<"]: " << dnz[r] << " onz["<<r<<"]: " << onz[r] << std::endl;
1101:   }
1102:   MatSeqAIJSetPreallocation(A, 0, dnz);
1103:   MatMPIAIJSetPreallocation(A, 0, dnz, 0, onz);
1104:   PetscFree2(dnz, onz);
1105:   MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR);
1106:   return(0);
1107: }

1109: /******************************** C Wrappers **********************************/

1113: PetscErrorCode WriteVTKHeader(PetscViewer viewer)
1114: {
1115:   return VTKViewer::writeHeader(viewer);
1116: }

1120: PetscErrorCode WriteVTKVertices(Mesh mesh, PetscViewer viewer)
1121: {
1122:   ALE::Obj<ALE::Mesh> m;

1125:   MeshGetMesh(mesh, &m);
1126:   return VTKViewer::writeVertices(m, viewer);
1127: }

1131: PetscErrorCode WriteVTKElements(Mesh mesh, PetscViewer viewer)
1132: {
1133:   ALE::Obj<ALE::Mesh> m;

1136:   MeshGetMesh(mesh, &m);
1137:   return VTKViewer::writeElements(m, viewer);
1138: }

1142: PetscErrorCode WritePCICEVertices(Mesh mesh, PetscViewer viewer)
1143: {
1144:   ALE::Obj<ALE::Mesh> m;

1147:   MeshGetMesh(mesh, &m);
1148:   return ALE::PCICE::Viewer::writeVertices(m, viewer);
1149: }

1153: PetscErrorCode WritePCICEElements(Mesh mesh, PetscViewer viewer)
1154: {
1155:   ALE::Obj<ALE::Mesh> m;

1158:   MeshGetMesh(mesh, &m);
1159:   return ALE::PCICE::Viewer::writeElements(m, viewer);
1160: }

1164: PetscErrorCode WritePyLithVertices(Mesh mesh, PetscViewer viewer)
1165: {
1166:   ALE::Obj<ALE::Mesh> m;

1169:   MeshGetMesh(mesh, &m);
1170:   return ALE::PyLith::Viewer::writeVertices(m, viewer);
1171: }

1175: PetscErrorCode WritePyLithElements(Mesh mesh, PetscViewer viewer)
1176: {
1177:   ALE::Obj<ALE::Mesh> m;

1180:   MeshGetMesh(mesh, &m);
1181:   return ALE::PyLith::Viewer::writeElements(m, viewer);
1182: }

1186: PetscErrorCode WritePyLithVerticesLocal(Mesh mesh, PetscViewer viewer)
1187: {
1188:   ALE::Obj<ALE::Mesh> m;

1191:   MeshGetMesh(mesh, &m);
1192:   return ALE::PyLith::Viewer::writeVerticesLocal(m, viewer);
1193: }

1197: PetscErrorCode WritePyLithElementsLocal(Mesh mesh, PetscViewer viewer)
1198: {
1199:   ALE::Obj<ALE::Mesh> m;

1202:   MeshGetMesh(mesh, &m);
1203:   return ALE::PyLith::Viewer::writeElementsLocal(m, viewer);
1204: }

1208: PetscErrorCode WritePyLithSplitLocal(Mesh mesh, PetscViewer viewer)
1209: {
1210:   ALE::Obj<ALE::Mesh> m;

1213:   MeshGetMesh(mesh, &m);
1214:   return ALE::PyLith::Viewer::writeSplitLocal(m, m->getSplitSection(), viewer);
1215: }