Actual source code: section.c

petsc-3.4.5 2014-06-29
  1: #include <petsc-private/meshimpl.h>   /*I      "petscdmmesh.h"   I*/
  2: #include <petscdmmesh_viewers.hh>

  4: /* Logging support */
  5: PetscClassId  SECTIONREAL_CLASSID;
  6: PetscLogEvent SectionReal_View;
  7: PetscClassId  SECTIONINT_CLASSID;
  8: PetscLogEvent SectionInt_View;
  9: PetscClassId  SECTIONPAIR_CLASSID;
 10: PetscLogEvent SectionPair_View;

 14: PetscErrorCode SectionRealView_Sieve(SectionReal section, PetscViewer viewer)
 15: {
 16:   PetscBool      iascii, isbinary, isdraw;

 20:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
 21:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERBINARY, &isbinary);
 22:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW, &isdraw);

 24:   if (iascii) {
 25:     ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
 26:     ALE::Obj<PETSC_MESH_TYPE>                    b;
 27:     const char                                   *name;

 29:     SectionRealGetSection(section, s);
 30:     SectionRealGetBundle(section, b);
 31:     PetscObjectGetName((PetscObject) section, &name);
 32:     SectionView_Sieve_Ascii(b, s, name, viewer);
 33:   }
 34:   return(0);
 35: }

 39: /*@C
 40:    SectionRealView - Views a Section object.

 42:    Collective on Section

 44:    Input Parameters:
 45: +  section - the Section
 46: -  viewer - an optional visualization context

 48:    Notes:
 49:    The available visualization contexts include
 50: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
 51: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
 52:          output where only the first processor opens
 53:          the file.  All other processors send their
 54:          data to the first processor to print.

 56:    You can change the format the section is printed using the
 57:    option PetscViewerSetFormat().

 59:    The user can open alternative visualization contexts with
 60: +    PetscViewerASCIIOpen() - Outputs section to a specified file
 61: .    PetscViewerBinaryOpen() - Outputs section in binary to a
 62:          specified file; corresponding input uses SectionLoad()
 63: .    PetscViewerDrawOpen() - Outputs section to an X window display

 65:    The user can call PetscViewerSetFormat() to specify the output
 66:    format of ASCII printed objects (when using PETSC_VIEWER_STDOUT_SELF,
 67:    PETSC_VIEWER_STDOUT_WORLD and PetscViewerASCIIOpen).  Available formats include
 68: +    PETSC_VIEWER_DEFAULT - default, prints section information
 69: -    PETSC_VIEWER_ASCII_VTK - outputs a VTK file describing the section

 71:    Level: beginner

 73:    Concepts: section^printing
 74:    Concepts: section^saving to disk

 76: .seealso: VecView(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerBinaryOpen(), PetscViewerCreate()
 77: @*/
 78: PetscErrorCode SectionRealView(SectionReal section, PetscViewer viewer)
 79: {

 85:   if (!viewer) {
 86:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)section),&viewer);
 87:   }

 91:   PetscLogEventBegin(SectionReal_View,0,0,0,0);
 92:   (*section->ops->view)(section, viewer);
 93:   PetscLogEventEnd(SectionReal_View,0,0,0,0);
 94:   return(0);
 95: }

 99: /*@C
100:   SectionRealDuplicate - Create an equivalent Section object

102:   Not collective

104:   Input Parameter:
105: . section - the section object

107:   Output Parameter:
108: . newSection - the duplicate

110:   Level: advanced

112: .seealso SectionRealCreate(), SectionRealSetSection()
113: @*/
114: PetscErrorCode  SectionRealDuplicate(SectionReal section, SectionReal *newSection)
115: {

121:   const ALE::Obj<PETSC_MESH_TYPE::real_section_type>& s = section->s;
122:   ALE::Obj<PETSC_MESH_TYPE::real_section_type>        t = new PETSC_MESH_TYPE::real_section_type(s->comm(), s->debug());

124:   t->setAtlas(s->getAtlas());
125:   t->allocateStorage();
126:   t->copyBC(s);
127:   SectionRealCreate(s->comm(), newSection);
128:   SectionRealSetSection(*newSection, t);
129:   SectionRealSetBundle(*newSection, section->b);
130:   return(0);
131: }

135: /*@C
136:   SectionRealGetSection - Gets the internal section object

138:   Not collective

140:   Input Parameter:
141: . section - the section object

143:   Output Parameter:
144: . s - the internal section object

146:   Level: advanced

148: .seealso SectionRealCreate(), SectionRealSetSection()
149: @*/
150: PetscErrorCode  SectionRealGetSection(SectionReal section, ALE::Obj<PETSC_MESH_TYPE::real_section_type>& s)
151: {
154:   s = section->s;
155:   return(0);
156: }

160: /*@C
161:   SectionRealSetSection - Sets the internal section object

163:   Not collective

165:   Input Parameters:
166: + section - the section object
167: - s - the internal section object

169:   Level: advanced

171: .seealso SectionRealCreate(), SectionRealGetSection()
172: @*/
173: PetscErrorCode  SectionRealSetSection(SectionReal section, const ALE::Obj<PETSC_MESH_TYPE::real_section_type>& s)
174: {

179:   if (!s.isNull()) {PetscObjectSetName((PetscObject) section, s->getName().c_str());}
180:   section->s = s;
181:   return(0);
182: }

186: /*@C
187:   SectionRealGetBundle - Gets the section bundle

189:   Not collective

191:   Input Parameter:
192: . section - the section object

194:   Output Parameter:
195: . b - the section bundle

197:   Level: advanced

199: .seealso SectionRealCreate(), SectionRealGetSection(), SectionRealSetSection()
200: @*/
201: PetscErrorCode  SectionRealGetBundle(SectionReal section, ALE::Obj<PETSC_MESH_TYPE>& b)
202: {
205:   b = section->b;
206:   return(0);
207: }

211: /*@C
212:   SectionRealSetBundle - Sets the section bundle

214:   Not collective

216:   Input Parameters:
217: + section - the section object
218: - b - the section bundle

220:   Level: advanced

222: .seealso SectionRealCreate(), SectionRealGetSection(), SectionRealSetSection()
223: @*/
224: PetscErrorCode  SectionRealSetBundle(SectionReal section, const ALE::Obj<PETSC_MESH_TYPE>& b)
225: {
228:   section->b = b;
229:   return(0);
230: }

234: /*@C
235:   SectionRealCreate - Creates a Section object, used to manage data for an unstructured problem
236:   described by a Sieve.

238:   Collective on MPI_Comm

240:   Input Parameter:
241: . comm - the processors that will share the global section

243:   Output Parameters:
244: . section - the section object

246:   Level: advanced

248: .seealso SectionRealDestroy(), SectionRealView()
249: @*/
250: PetscErrorCode  SectionRealCreate(MPI_Comm comm, SectionReal *section)
251: {
253:   SectionReal    s;

257:   *section = NULL;

259:   PetscHeaderCreate(s,_p_SectionReal,struct _SectionRealOps,SECTIONREAL_CLASSID,"SectionReal","Section","DM",comm,SectionRealDestroy,0);

261:   s->ops->view            = SectionRealView_Sieve;
262:   s->ops->restrictClosure = SectionRealRestrict;
263:   s->ops->update          = SectionRealUpdate;

265:   PetscObjectChangeTypeName((PetscObject) s, "sieve");

267:   new(&s->s) ALE::Obj<PETSC_MESH_TYPE::real_section_type>(PETSC_MESH_TYPE::real_section_type(comm));
268:   new(&s->b) ALE::Obj<PETSC_MESH_TYPE>(NULL);
269:   *section = s;
270:   return(0);
271: }

275: /*@
276:   SectionRealDestroy - Destroys a section.

278:   Collective on Section

280:   Input Parameter:
281: . section - the section object

283:   Level: advanced

285: .seealso SectionRealCreate(), SectionRealView()
286: @*/
287: PetscErrorCode  SectionRealDestroy(SectionReal *section)
288: {

292:   if (!*section) return(0);
294:   if (--((PetscObject)(*section))->refct > 0) return(0);
295:   PetscHeaderDestroy(section);
296:   return(0);
297: }

301: /*@
302:   SectionRealDistribute - Distributes the sections.

304:   Not Collective

306:   Input Parameters:
307: + serialSection - The original Section object
308: - parallelMesh - The parallel DMMesh

310:   Output Parameter:
311: . parallelSection - The distributed Section object

313:   Level: intermediate

315: .keywords: mesh, section, distribute
316: .seealso: DMMeshCreate()
317: @*/
318: PetscErrorCode SectionRealDistribute(SectionReal serialSection, DM parallelMesh, SectionReal *parallelSection)
319: {
320:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> oldSection;
321:   ALE::Obj<PETSC_MESH_TYPE>                    m;
322:   PetscErrorCode                               ierr;

325:   SectionRealGetSection(serialSection, oldSection);
326:   DMMeshGetMesh(parallelMesh, m);
327:   SectionRealCreate(oldSection->comm(), parallelSection);
328: #if defined(PETSC_OPT_SIEVE)
329:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> newSection;

331:   // We assume all integer sections are complete sections
332:   newSection->setName(oldSection->getName());
333:   newSection->setChart(m->getSieve()->getChart());
334:   //distributeSection(oldSection, partition, m->getRenumbering(), m->getDistSendOverlap(), m->getDistRecvOverlap(), newSection);
335:   SectionRealSetSection(*parallelSection, newSection);
336:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Not working because the partition is unavailable");
337: #else
338:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> newSection = ALE::Distribution<PETSC_MESH_TYPE>::distributeSection(oldSection, m, m->getDistSendOverlap(), m->getDistRecvOverlap());
339:   SectionRealSetSection(*parallelSection, newSection);
340: #endif
341:   return(0);
342: }

346: /*@C
347:   SectionRealRestrict - Restricts the SectionReal to a subset of the topology, returning an array of values.

349:   Not collective

351:   Input Parameters:
352: + section - the section object
353: - point - the Sieve point

355:   Output Parameter:
356: . values - The values associated with the submesh

358:   Level: advanced

360: .seealso SectionUpdate(), SectionCreate(), SectionView()
361: @*/
362: PetscErrorCode  SectionRealRestrict(SectionReal section, PetscInt point, PetscScalar *values[])
363: {
367:   try {
368:     *values = (PetscScalar*) section->s->restrictPoint(point);
369:   } catch(ALE::Exception e) {
370:     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid argument: %s", e.message());
371:   }
372:   return(0);
373: }

377: /*@C
378:   SectionRealUpdate - Updates the array of values associated to a subset of the topology in this Section.

380:   Not collective

382:   Input Parameters:
383: + section - the section object
384: . point - the Sieve point
385: . values - The values associated with the submesh
386: - mode - The insertion mode

388:   Level: advanced

390: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
391: @*/
392: PetscErrorCode  SectionRealUpdate(SectionReal section, PetscInt point, const PetscScalar values[], InsertMode mode)
393: {
397: #if defined(PETSC_USE_COMPLEX)
398:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "SectionReal does not support complex updates");
399: #else
400:   try {
401:     if (mode == INSERT_VALUES) {
402:       section->b->update(section->s, point, values);
403:     } else if (mode == ADD_VALUES) {
404:       section->b->updateAdd(section->s, point, values);
405:     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Invalid insertion mode: %d", mode);
406:   } catch(ALE::Exception e) {
407:     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid argument: %s", e.message());
408:   }
409: #endif
410:   return(0);
411: }

415: /*@C
416:   SectionRealRestrictClosure - Returns an array with the values in a given closure

418:   Not Collective

420:   Input Parameters:
421: + section - The section
422: . mesh    - The DMMesh object
423: - point   - The sieve point

425:   Output Parameter:
426: . array - The array full of values in the closure

428:   Level: intermediate

430: .keywords: mesh, elements
431: .seealso: DMMeshCreate()
432: @*/
433: PetscErrorCode SectionRealRestrictClosure(SectionReal section, DM dm, PetscInt point, const PetscScalar *values[])
434: {
435:   ALE::Obj<PETSC_MESH_TYPE>                    m;
436:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
437:   PetscErrorCode                               ierr;

440:   DMMeshGetMesh(dm, m);
441:   SectionRealGetSection(section, s);
442: #if defined(PETSC_USE_COMPLEX)
443:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "SectionReal does not support complex restriction");
444: #else
445:   try {
446:     *values = m->restrictClosure(s, point);
447:   } catch(ALE::Exception e) {
448:     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid argument: %s", e.message());
449:   }
450: #endif
451:   return(0);
452: }

456: /*@C
457:   SectionRealRestrictClosureWithArray - Returns an array with the values in a given closure

459:   Not Collective

461:   Input Parameters:
462: + section - The section
463: . mesh    - The DMMesh object
464: . point   - The sieve point
465: . n       - The array size
466: - array   - The array to fill up

468:   Output Parameter:
469: . array - The array full of values in the closure

471:   Level: intermediate

473: .keywords: mesh, elements
474: .seealso: DMMeshCreate()
475: @*/
476: PetscErrorCode SectionRealRestrictClosureWithArray(SectionReal section, DM dm, PetscInt point, PetscInt n, PetscScalar values[])
477: {
478:   ALE::Obj<PETSC_MESH_TYPE>                    m;
479:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
480:   PetscErrorCode                               ierr;

483:   DMMeshGetMesh(dm, m);
484:   SectionRealGetSection(section, s);
485: #if defined(PETSC_USE_COMPLEX)
486:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "SectionReal does not support complex restriction");
487: #else
488:   try {
489:     m->restrictClosure(s, point, values, n);
490:   } catch(ALE::Exception e) {
491:     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid argument: %s", e.message());
492:   }
493: #endif
494:   return(0);
495: }

499: /*@C
500:   SectionRealUpdateClosure - Updates the values in a given closure from the array

502:   Not Collective

504:   Input Parameters:
505: + section - The section
506: . mesh    - The DMMesh object
507: . point   - The sieve point
508: . array   - The array to fill up
509: - mode    - The insertion mode

511:   Output Parameter:
512: . array - The array full of values in the closure

514:   Level: intermediate

516: .keywords: mesh, elements
517: .seealso: DMMeshCreate()
518: @*/
519: PetscErrorCode SectionRealUpdateClosure(SectionReal section, DM dm, PetscInt point, PetscScalar values[], InsertMode mode)
520: {
521:   ALE::Obj<PETSC_MESH_TYPE>                    m;
522:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
523:   PetscErrorCode                               ierr;

526:   DMMeshGetMesh(dm, m);
527:   SectionRealGetSection(section, s);
528: #if defined(PETSC_USE_COMPLEX)
529:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "SectionReal does not support complex update");
530: #else
531:   try {
532:     if (mode == INSERT_VALUES) {
533:       m->update(s, point, values);
534:     } else if (mode == ADD_VALUES) {
535:       m->updateAdd(s, point, values);
536:     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Invalid insertion mode: %d", mode);
537:   } catch(ALE::Exception e) {
538:     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid argument: %s", e.message());
539:   }
540: #endif
541:   return(0);
542: }

546: /*@
547:   SectionRealComplete - Exchanges data across the mesh overlap.

549:   Not collective

551:   Input Parameter:
552: . section - the section object

554:   Level: advanced

556: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
557: @*/
558: PetscErrorCode SectionRealComplete(SectionReal section)
559: {
560:   Obj<PETSC_MESH_TYPE::real_section_type> s;
561:   Obj<PETSC_MESH_TYPE>                    b;
562:   PetscErrorCode                          ierr;

565:   SectionRealGetSection(section, s);
566:   SectionRealGetBundle(section, b);
567: #if 0
568:   ALE::Distribution<PETSC_MESH_TYPE>::completeSection(b, s);
569: #else
570:   ALE::Completion::completeSectionAdd(b->getSendOverlap(), b->getRecvOverlap(), s, s);
571: #endif
572:   return(0);
573: }

577: /*@
578:   SectionRealZero - Zero out the entries

580:   Not collective

582:   Input Parameter:
583: . section - the section object

585:   Level: advanced

587: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
588: @*/
589: PetscErrorCode SectionRealZero(SectionReal section)
590: {
591:   Obj<PETSC_MESH_TYPE::real_section_type> s;
592:   PetscErrorCode                          ierr;

595:   SectionRealGetSection(section, s);
596:   s->zero();
597:   return(0);
598: }

602: /*@
603:   SectionRealGetFiberDimension - Get the size of the vector space attached to the point

605:   Not collective

607:   Input Parameters:
608: + section - the section object
609: - point - the Sieve point

611:   Output Parameters:
612: . size - The fiber dimension

614:   Level: advanced

616: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
617: @*/
618: PetscErrorCode  SectionRealGetFiberDimension(SectionReal section, PetscInt point, PetscInt *size)
619: {
622:   *size = section->s->getFiberDimension(point);
623:   return(0);
624: }

628: /*@
629:   SectionRealSetFiberDimension - Set the size of the vector space attached to the point

631:   Not collective

633:   Input Parameters:
634: + section - the section object
635: . point - the Sieve point
636: - size - The fiber dimension

638:   Level: advanced

640: .seealso SectionRealSetFiberDimensionField(), SectionRealRestrict(), SectionRealCreate(), SectionRealView()
641: @*/
642: PetscErrorCode  SectionRealSetFiberDimension(SectionReal section, PetscInt point, const PetscInt size)
643: {
646:   section->s->setFiberDimension(point, size);
647:   return(0);
648: }

652: /*@
653:   SectionRealSetFiberDimensionField - Set the size of the vector space attached to the point for a given field

655:   Not collective

657:   Input Parameters:
658: + section - the section object
659: . point - the Sieve point
660: . size - The fiber dimension
661: - field - The field number

663:   Level: advanced

665: .seealso SectionRealSetFiberDimension(), SectionRealRestrict(), SectionRealCreate(), SectionRealView()
666: @*/
667: PetscErrorCode  SectionRealSetFiberDimensionField(SectionReal section, PetscInt point, const PetscInt size, const PetscInt field)
668: {
671:   section->s->setFiberDimension(point, size, field);
672:   return(0);
673: }

677: /*@
678:   SectionRealGetSize - Gets the number of local dofs in this Section

680:   Not collective

682:   Input Parameter:
683: . section - the section object

685:   Output Parameter:
686: . size - the section size

688:   Level: advanced

690: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
691: @*/
692: PetscErrorCode SectionRealGetSize(SectionReal section, PetscInt *size)
693: {
694:   Obj<PETSC_MESH_TYPE::real_section_type> s;
695:   PetscErrorCode                          ierr;

699:   SectionRealGetSection(section, s);
700:   *size = s->size();
701:   return(0);
702: }

706: /*@
707:   SectionRealAllocate - Allocate storage for this section

709:   Not collective

711:   Input Parameter:
712: . section - the section object

714:   Level: advanced

716: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
717: @*/
718: PetscErrorCode  SectionRealAllocate(SectionReal section)
719: {
722:   section->b->allocate(section->s);
723:   return(0);
724: }

728: /*@
729:   SectionRealCreateLocalVector - Creates a vector with the local piece of the Section

731:   Collective on DMMesh

733:   Input Parameter:
734: . section - the Section

736:   Output Parameter:
737: . localVec - the local vector

739:   Level: advanced

741:   Notes: The vector can safely be destroyed using VecDestroy().
742: .seealso DMMeshDestroy(), DMMeshCreate()
743: @*/
744: PetscErrorCode  SectionRealCreateLocalVector(SectionReal section, Vec *localVec)
745: {
746:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
747:   PetscErrorCode                               ierr;

750: #if defined(PETSC_USE_COMPLEX)
751:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "SectionReal does not support complex Vec");
752: #else
753:   SectionRealGetSection(section, s);
754:   VecCreateSeqWithArray(PETSC_COMM_SELF,1, s->getStorageSize(), s->restrictSpace(), localVec);
755: #endif
756:   return(0);
757: }

761: /*@
762:   SectionRealAddSpace - Add another field to this section

764:   Collective on DMMesh

766:   Input Parameter:
767: . section - the Section

769:   Level: advanced

771: .seealso SectionRealCreate(), SectionRealGetFibration()
772: @*/
773: PetscErrorCode  SectionRealAddSpace(SectionReal section)
774: {
775:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
776:   PetscErrorCode                               ierr;

779:   SectionRealGetSection(section, s);
780:   s->addSpace();
781:   return(0);
782: }

786: /*@
787:   SectionRealGetFibration - Creates a section for only the data associated with the given field

789:   Collective on DMMesh

791:   Input Parameter:
792: + section - the Section
793: - field- The field number

795:   Output Parameter:
796: . subsection - the section of the given field

798:   Level: advanced

800: .seealso SectionRealCreate()
801: @*/
802: PetscErrorCode  SectionRealGetFibration(SectionReal section, const PetscInt field, SectionReal *subsection)
803: {
804:   ALE::Obj<PETSC_MESH_TYPE>                    b;
805:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
806:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> t;
807:   MPI_Comm                                     comm;
808:   PetscErrorCode                               ierr;

811:   PetscObjectGetComm((PetscObject) section, &comm);
812:   SectionRealGetBundle(section, b);
813:   SectionRealGetSection(section, s);
814:   SectionRealCreate(comm, subsection);
815:   SectionRealSetBundle(*subsection, b);
816:   t    = s->getFibration(field);
817:   SectionRealSetSection(*subsection, t);
818:   return(0);
819: }

823: /*@C
824:   SectionRealToVecDM - Maps the given section to a Vec

826:   Collective on Section

828:   Input Parameters:
829: + section - the real Section
830: - mesh - The DMMesh

832:   Output Parameter:
833: . vec - the Vec

835:   Level: intermediate

837: .seealso VecCreate(), SectionRealCreate()
838: @*/
839: PetscErrorCode  SectionRealToVecDM(SectionReal section, DM dm, ScatterMode mode, Vec vec)
840: {
841:   Vec            localVec;
842:   VecScatter     scatter;

847:   SectionRealCreateLocalVector(section, &localVec);
848:   DMMeshGetGlobalScatter(dm, &scatter);
849:   if (mode == SCATTER_FORWARD) {
850:     VecScatterBegin(scatter, localVec, vec, INSERT_VALUES, mode);
851:     VecScatterEnd(scatter, localVec, vec, INSERT_VALUES, mode);
852:   } else {
853:     VecScatterBegin(scatter, vec, localVec, INSERT_VALUES, mode);
854:     VecScatterEnd(scatter, vec, localVec, INSERT_VALUES, mode);
855:   }
856:   VecDestroy(&localVec);
857:   return(0);
858: }

862: /*@
863:   SectionRealToVec - Map between unassembled local Section storage and a globally assembled Vec

865:   Collective on VecScatter

867:   Input Parameters:
868: + section - the Section
869: . scatter - the scatter from the Section to the Vec
870: . mode - the mode, SCATTER_FORWARD (Section to Vec) or SCATTER_REVERSE (Vec to Section)
871: - vec - the Vec

873:   Level: advanced

875: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
876: @*/
877: PetscErrorCode  SectionRealToVec(SectionReal section, VecScatter scatter, ScatterMode mode, Vec vec)
878: {
879:   Vec            localVec;

884:   SectionRealCreateLocalVector(section, &localVec);
885:   if (mode == SCATTER_FORWARD) {
886:     VecScatterBegin(scatter, localVec, vec, INSERT_VALUES, mode);
887:     VecScatterEnd(scatter, localVec, vec, INSERT_VALUES, mode);
888:   } else {
889:     VecScatterBegin(scatter, vec, localVec, INSERT_VALUES, mode);
890:     VecScatterEnd(scatter, vec, localVec, INSERT_VALUES, mode);
891:   }
892:   VecDestroy(&localVec);
893:   return(0);
894: }

898: /*@
899:   SectionRealClear - Dellocate storage for this section

901:   Not collective

903:   Input Parameter:
904: . section - the section object

906:   Level: advanced

908: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
909: @*/
910: PetscErrorCode  SectionRealClear(SectionReal section)
911: {
914:   section->s->clear();
915:   return(0);
916: }

920: /*@
921:   SectionRealSet - Sets all the values to the given value

923:   Not collective

925:   Input Parameters:
926: + section - the real Section
927: - val - the value

929:   Level: intermediate

931: .seealso VecNorm(), SectionRealCreate()
932: @*/
933: PetscErrorCode  SectionRealSet(SectionReal section, PetscReal val)
934: {
935:   Obj<PETSC_MESH_TYPE::real_section_type> s;
936:   PetscErrorCode                          ierr;

939:   SectionRealGetSection(section, s);
940:   s->set(val);
941:   return(0);
942: }

946: /*@C
947:   SectionRealNorm - Computes the vector norm.

949:   Collective on Section

951:   Input Parameters:
952: +  section - the real Section
953: -  type - one of NORM_1, NORM_2, NORM_INFINITY.  Also available
954:           NORM_1_AND_2, which computes both norms and stores them
955:           in a two element array.

957:   Output Parameter:
958: . val - the norm

960:   Notes:
961: $     NORM_1 denotes sum_i |x_i|
962: $     NORM_2 denotes sqrt(sum_i (x_i)^2)
963: $     NORM_INFINITY denotes max_i |x_i|

965:   Level: intermediate

967: .seealso VecNorm(), SectionRealCreate()
968: @*/
969: PetscErrorCode  SectionRealNorm(SectionReal section, DM dm, NormType type, PetscReal *val)
970: {
971:   Obj<PETSC_MESH_TYPE>                    m;
972:   Obj<PETSC_MESH_TYPE::real_section_type> s;
973:   Vec                                     v;
974:   PetscErrorCode                          ierr;

978:   DMMeshGetMesh(dm, m);
979:   SectionRealGetSection(section, s);
980:   const ALE::Obj<PETSC_MESH_TYPE::order_type>& order = m->getFactory()->getGlobalOrder(m, s->getName(), s);
981:   VecCreate(m->comm(), &v);
982:   VecSetSizes(v, order->getLocalSize(), order->getGlobalSize());
983:   VecSetFromOptions(v);
984:   SectionRealToVecDM(section, dm, SCATTER_FORWARD, v);
985:   VecNorm(v, type, val);
986:   VecDestroy(&v);
987:   return(0);
988: }

992: /*@
993:   SectionRealAXPY -

995:   Collective on Section

997:   Input Parameters:
998: +  section - the real Section
999: .  alpha - a scalar
1000: -  X - the other real Section

1002:   Output Parameter:
1003: . section - the difference

1005:   Level: intermediate

1007: .seealso VecNorm(), SectionRealCreate()
1008: @*/
1009: PetscErrorCode  SectionRealAXPY(SectionReal section, DM dm, PetscScalar alpha, SectionReal X)
1010: {
1011:   Obj<PETSC_MESH_TYPE>                    m;
1012:   Obj<PETSC_MESH_TYPE::real_section_type> s;
1013:   Obj<PETSC_MESH_TYPE::real_section_type> sX;
1014:   Vec                                     v, x;
1015:   PetscErrorCode                          ierr;

1019:   DMMeshGetMesh(dm, m);
1020:   SectionRealGetSection(section, s);
1021:   SectionRealGetSection(X, sX);
1022:   const ALE::Obj<PETSC_MESH_TYPE::order_type>& order = m->getFactory()->getGlobalOrder(m, s->getName(), s);
1023:   VecCreate(m->comm(), &v);
1024:   VecSetSizes(v, order->getLocalSize(), order->getGlobalSize());
1025:   VecSetFromOptions(v);
1026:   VecDuplicate(v, &x);
1027:   SectionRealToVecDM(section, dm, SCATTER_FORWARD, v);
1028:   SectionRealToVecDM(X,       dm, SCATTER_FORWARD, x);
1029:   VecAXPY(v, alpha, x);
1030:   SectionRealToVecDM(section, dm, SCATTER_REVERSE, v);
1031:   VecDestroy(&v);
1032:   VecDestroy(&x);
1033:   return(0);
1034: }

1038: /*@C
1039:   DMMeshGetVertexSectionReal - Create a Section over the vertices with the specified fiber dimension

1041:   Collective on DMMesh

1043:   Input Parameters:
1044: + mesh - The DMMesh object
1045: . name - The name of the section
1046: - fiberDim - The number of degrees of freedom per vertex

1048:   Output Parameter:
1049: . section - The section

1051:   Level: intermediate

1053: .keywords: mesh, section, vertex
1054: .seealso: DMMeshCreate(), SectionRealCreate()
1055: @*/
1056: PetscErrorCode DMMeshGetVertexSectionReal(DM dm, const char name[], PetscInt fiberDim, SectionReal *section)
1057: {
1058:   ALE::Obj<PETSC_MESH_TYPE>                    m;
1059:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
1060:   PetscErrorCode                               ierr;

1063:   DMMeshGetMesh(dm, m);
1064:   SectionRealCreate(m->comm(), section);
1065:   PetscObjectSetName((PetscObject) *section, name);
1066:   SectionRealSetBundle(*section, m);
1067:   SectionRealGetSection(*section, s);
1068:   s->setChart(m->getSieve()->getChart());
1069:   s->setName(name);
1070:   s->setFiberDimension(m->depthStratum(0), fiberDim);
1071:   m->allocate(s);
1072:   return(0);
1073: }

1077: /*@C
1078:   DMMeshGetCellSectionReal - Create a Section over the cells with the specified fiber dimension

1080:   Collective on DMMesh

1082:   Input Parameters:
1083: + mesh - The DMMesh object
1084: . name - The name of the section
1085: - fiberDim - The number of degrees of freedom per cell

1087:   Output Parameter:
1088: . section - The section

1090:   Level: intermediate

1092: .keywords: mesh, section, cell
1093: .seealso: DMMeshCreate(), SectionRealCreate(), DMMeshGetVertexSectionReal(), DMMeshCreateSectionRealIS()
1094: @*/
1095: PetscErrorCode DMMeshGetCellSectionReal(DM dm, const char name[], PetscInt fiberDim, SectionReal *section)
1096: {
1097:   ALE::Obj<PETSC_MESH_TYPE>                    m;
1098:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
1099:   PetscErrorCode                               ierr;

1102:   DMMeshGetMesh(dm, m);
1103:   SectionRealCreate(m->comm(), section);
1104:   PetscObjectSetName((PetscObject) *section, name);
1105:   SectionRealSetBundle(*section, m);
1106:   SectionRealGetSection(*section, s);
1107:   s->setChart(m->getSieve()->getChart());
1108:   s->setName(name);
1109:   s->setFiberDimension(m->heightStratum(0), fiberDim);
1110:   m->allocate(s);
1111:   return(0);
1112: }

1116: /*@C
1117:   DMMeshCreateSectionRealIS - Create a Section over the points specified in an IS.

1119:   Collective on DMMesh

1121:   Input Parameters:
1122:   + dm - the DMMesh object
1123:   . is - The IS describing the points associated with the degrees of freedom
1124:   . name - The name of the section
1125:   - fiberDim - The number of degrees of freedom per point of the IS

1127:   Output Parameter:
1128:   . section - The section

1130:   Level: intermediate

1132: .keywords: mesh, section, cell
1133: .seealso: DMMeshCreate(), SectionRealCreate(), DMMeshGetVertexSectionReal(), DMMeshGetCellSectionReal()

1135: @*/
1136: PetscErrorCode DMMeshCreateSectionRealIS(DM dm,IS is,const char name[],PetscInt fiberDim,SectionReal *section)
1137: {
1138:   MPI_Comm       comm;
1139:   PetscSection   s;
1140:   PetscInt       pStart,pEnd;
1141:   const PetscInt *points;
1142:   PetscInt       numpoints,p;

1146:   PetscObjectGetComm((PetscObject) dm,&comm);

1148:   PetscSectionCreate(comm,&s);
1149:   DMMeshGetChart(dm,&pStart,&pEnd);
1150:   PetscSectionSetChart(s,pStart,pEnd);

1152:   ISGetLocalSize(is,&numpoints);
1153:   ISGetIndices(is,&points);
1154:   for (p =0; p < numpoints; p++) {
1155:     PetscSectionSetDof(s,points[p],fiberDim);
1156:   }
1157:   ISRestoreIndices(is,&points);
1158:   DMMeshSetSection(dm,name,s);
1159:   PetscSectionDestroy(&s);
1160:   DMMeshGetSectionReal(dm,name,section);
1161:   return(0);
1162: }


1167: /*@C
1168:   DMMeshCreateGlobalRealVector - Creates a vector of the correct size to be gathered into by the mesh.

1170:   Collective on DMMesh

1172:   Input Parameters:
1173: + mesh - the mesh object
1174: - section - The SectionReal

1176:   Output Parameters:
1177: . gvec - the global vector

1179:   Level: advanced

1181: .seealso DMMeshDestroy(), DMMeshCreate(), DMMeshCreateGlobalRealVector()
1182: @*/
1183: PetscErrorCode DMMeshCreateGlobalRealVector(DM dm, SectionReal section, Vec *gvec)
1184: {
1185:   ALE::Obj<PETSC_MESH_TYPE>                    m;
1186:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
1187:   const char                                   *name;
1188:   PetscErrorCode                               ierr;

1191:   DMMeshGetMesh(dm, m);
1192:   SectionRealGetSection(section, s);
1193:   PetscObjectGetName((PetscObject) section, &name);
1194:   const ALE::Obj<PETSC_MESH_TYPE::order_type>& order = m->getFactory()->getGlobalOrder(m, name, s);

1196:   VecCreate(m->comm(), gvec);
1197:   VecSetSizes(*gvec, order->getLocalSize(), order->getGlobalSize());
1198:   VecSetFromOptions(*gvec);
1199:   return(0);
1200: }

1204: PetscErrorCode SectionIntView_Sieve(SectionInt section, PetscViewer viewer)
1205: {
1206:   PetscBool      iascii;

1210:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);

1212:   if (iascii) {
1213:     ALE::Obj<PETSC_MESH_TYPE::int_section_type> s;
1214:     ALE::Obj<PETSC_MESH_TYPE>                   b;
1215:     const char                                  *name;

1217:     SectionIntGetSection(section, s);
1218:     SectionIntGetBundle(section, b);
1219:     PetscObjectGetName((PetscObject) section, &name);
1220:     SectionView_Sieve_Ascii(b, s, name, viewer);
1221:   }
1222:   return(0);
1223: }

1227: /*@C
1228:    SectionIntView - Views a Section object.

1230:    Collective on Section

1232:    Input Parameters:
1233: +  section - the Section
1234: -  viewer - an optional visualization context

1236:    Notes:
1237:    The available visualization contexts include
1238: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1239: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1240:          output where only the first processor opens
1241:          the file.  All other processors send their
1242:          data to the first processor to print.

1244:    You can change the format the section is printed using the
1245:    option PetscViewerSetFormat().

1247:    The user can open alternative visualization contexts with
1248: +    PetscViewerASCIIOpen() - Outputs section to a specified file
1249: .    PetscViewerBinaryOpen() - Outputs section in binary to a
1250:          specified file; corresponding input uses SectionLoad()
1251: .    PetscViewerDrawOpen() - Outputs section to an X window display

1253:    The user can call PetscViewerSetFormat() to specify the output
1254:    format of ASCII printed objects (when using PETSC_VIEWER_STDOUT_SELF,
1255:    PETSC_VIEWER_STDOUT_WORLD and PetscViewerASCIIOpen).  Available formats include
1256: +    PETSC_VIEWER_DEFAULT - default, prints section information
1257: -    PETSC_VIEWER_ASCII_VTK - outputs a VTK file describing the section

1259:    Level: beginner

1261:    Concepts: section^printing
1262:    Concepts: section^saving to disk

1264: .seealso: VecView(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerBinaryOpen(), PetscViewerCreate()
1265: @*/
1266: PetscErrorCode SectionIntView(SectionInt section, PetscViewer viewer)
1267: {

1273:   if (!viewer) {
1274:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)section),&viewer);
1275:   }

1279:   PetscLogEventBegin(SectionInt_View,0,0,0,0);
1280:   (*section->ops->view)(section, viewer);
1281:   PetscLogEventEnd(SectionInt_View,0,0,0,0);
1282:   return(0);
1283: }

1287: /*@C
1288:   SectionIntGetSection - Gets the internal section object

1290:   Not collective

1292:   Input Parameter:
1293: . section - the section object

1295:   Output Parameter:
1296: . s - the internal section object

1298:   Level: advanced

1300: .seealso SectionIntCreate(), SectionIntSetSection()
1301: @*/
1302: PetscErrorCode  SectionIntGetSection(SectionInt section, ALE::Obj<PETSC_MESH_TYPE::int_section_type>& s)
1303: {
1306:   s = section->s;
1307:   return(0);
1308: }

1312: /*@C
1313:   SectionIntSetSection - Sets the internal section object

1315:   Not collective

1317:   Input Parameters:
1318: + section - the section object
1319: - s - the internal section object

1321:   Level: advanced

1323: .seealso SectionIntCreate(), SectionIntGetSection()
1324: @*/
1325: PetscErrorCode  SectionIntSetSection(SectionInt section, const ALE::Obj<PETSC_MESH_TYPE::int_section_type>& s)
1326: {

1331:   if (!s.isNull()) {PetscObjectSetName((PetscObject) section, s->getName().c_str());}
1332:   section->s = s;
1333:   return(0);
1334: }

1338: /*@C
1339:   SectionIntGetBundle - Gets the section bundle

1341:   Not collective

1343:   Input Parameter:
1344: . section - the section object

1346:   Output Parameter:
1347: . b - the section bundle

1349:   Level: advanced

1351: .seealso SectionIntCreate(), SectionIntGetSection(), SectionIntSetSection()
1352: @*/
1353: PetscErrorCode  SectionIntGetBundle(SectionInt section, ALE::Obj<PETSC_MESH_TYPE>& b)
1354: {
1357:   b = section->b;
1358:   return(0);
1359: }

1363: /*@C
1364:   SectionIntSetBundle - Sets the section bundle

1366:   Not collective

1368:   Input Parameters:
1369: + section - the section object
1370: - b - the section bundle

1372:   Level: advanced

1374: .seealso SectionIntCreate(), SectionIntGetSection(), SectionIntSetSection()
1375: @*/
1376: PetscErrorCode  SectionIntSetBundle(SectionInt section, const ALE::Obj<PETSC_MESH_TYPE>& b)
1377: {
1380:   section->b = b;
1381:   return(0);
1382: }

1386: /*@C
1387:   SectionIntCreate - Creates a Section object, used to manage data for an unstructured problem
1388:   described by a Sieve.

1390:   Collective on MPI_Comm

1392:   Input Parameter:
1393: . comm - the processors that will share the global section

1395:   Output Parameters:
1396: . section - the section object

1398:   Level: advanced

1400: .seealso SectionIntDestroy(), SectionIntView()
1401: @*/
1402: PetscErrorCode  SectionIntCreate(MPI_Comm comm, SectionInt *section)
1403: {
1405:   SectionInt     s;

1409:   *section = NULL;

1411:   PetscHeaderCreate(s,_p_SectionInt,struct _SectionIntOps,SECTIONINT_CLASSID,"SectionInt","Section","DM",comm,SectionIntDestroy,0);

1413:   s->ops->view            = SectionIntView_Sieve;
1414:   s->ops->restrictClosure = SectionIntRestrict;
1415:   s->ops->update          = SectionIntUpdate;

1417:   PetscObjectChangeTypeName((PetscObject) s, "sieve");

1419:   new(&s->s) ALE::Obj<PETSC_MESH_TYPE::int_section_type>(PETSC_MESH_TYPE::int_section_type(comm));
1420:   new(&s->b) ALE::Obj<PETSC_MESH_TYPE>(NULL);
1421:   *section = s;
1422:   return(0);
1423: }

1427: /*@
1428:   SectionIntDestroy - Destroys a section.

1430:   Collective on Section

1432:   Input Parameter:
1433: . section - the section object

1435:   Level: advanced

1437: .seealso SectionIntCreate(), SectionIntView()
1438: @*/
1439: PetscErrorCode  SectionIntDestroy(SectionInt *section)
1440: {

1444:   if (!*section) return(0);
1446:   if (--((PetscObject)(*section))->refct > 0) return(0);
1447:   PetscHeaderDestroy(section);
1448:   return(0);
1449: }

1453: /*@
1454:   SectionIntDistribute - Distributes the sections.

1456:   Not Collective

1458:   Input Parameters:
1459: + serialSection - The original Section object
1460: - parallelMesh - The parallel DMMesh

1462:   Output Parameter:
1463: . parallelSection - The distributed Section object

1465:   Level: intermediate

1467: .keywords: mesh, section, distribute
1468: .seealso: DMMeshCreate()
1469: @*/
1470: PetscErrorCode SectionIntDistribute(SectionInt serialSection, DM parallelMesh, SectionInt *parallelSection)
1471: {
1472:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> oldSection;
1473:   ALE::Obj<PETSC_MESH_TYPE>                   m;
1474:   PetscErrorCode                              ierr;

1477:   SectionIntGetSection(serialSection, oldSection);
1478:   DMMeshGetMesh(parallelMesh, m);
1479:   SectionIntCreate(oldSection->comm(), parallelSection);
1480: #if defined(PETSC_OPT_SIEVE)
1481:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> newSection;

1483:   // We assume all integer sections are complete sections
1484:   newSection->setName(oldSection->getName());
1485:   newSection->setChart(m->getSieve()->getChart());
1486:   //distributeSection(oldSection, partition, m->getRenumbering(), m->getDistSendOverlap(), m->getDistRecvOverlap(), newSection);
1487:   SectionIntSetSection(*parallelSection, newSection);
1488:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Not working because the partition is unavailable");
1489: #else
1490:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> newSection = ALE::Distribution<PETSC_MESH_TYPE>::distributeSection(oldSection, m, m->getDistSendOverlap(), m->getDistRecvOverlap());
1491:   SectionIntSetSection(*parallelSection, newSection);
1492: #endif
1493:   return(0);
1494: }

1498: /*@C
1499:   SectionIntRestrict - Restricts the SectionInt to a subset of the topology, returning an array of values.

1501:   Not collective

1503:   Input Parameters:
1504: + section - the section object
1505: - point - the Sieve point

1507:   Output Parameter:
1508: . values - The values associated with the submesh

1510:   Level: advanced

1512: .seealso SectionIntUpdate(), SectionIntCreate(), SectionIntView()
1513: @*/
1514: PetscErrorCode  SectionIntRestrict(SectionInt section, PetscInt point, PetscInt *values[])
1515: {
1519:   *values = (PetscInt*) section->b->restrictClosure(section->s, point);
1520:   return(0);
1521: }

1525: /*@C
1526:   SectionIntUpdate - Updates the array of values associated to a subset of the topology in this Section.

1528:   Not collective

1530:   Input Parameters:
1531: + section - the section object
1532: . point - the Sieve point
1533: . values - The values associated with the submesh
1534: - mode - The insertion mode

1536:   Level: advanced

1538: .seealso SectionIntRestrict(), SectionIntCreate(), SectionIntView()
1539: @*/
1540: PetscErrorCode  SectionIntUpdate(SectionInt section, PetscInt point, const PetscInt values[], InsertMode mode)
1541: {
1545:   if (mode == INSERT_VALUES) {
1546:     section->b->update(section->s, point, values);
1547:   } else if (mode == ADD_VALUES) {
1548:     section->b->updateAdd(section->s, point, values);
1549:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Invalid insertion mode: %d", mode);
1550:   return(0);
1551: }

1555: /*@C
1556:   SectionIntRestrictClosureWithArray - Returns an array with the values in a given closure

1558:   Not Collective

1560:   Input Parameters:
1561: + section - The section
1562: . mesh    - The DMMesh object
1563: . point   - The sieve point
1564: . n       - The array size
1565: - array   - The array to fill up

1567:   Output Parameter:
1568: . array - The array full of values in the closure

1570:   Level: intermediate

1572: .keywords: mesh, elements
1573: .seealso: DMMeshCreate()
1574: @*/
1575: PetscErrorCode SectionIntRestrictClosureWithArray(SectionInt section, DM dm, PetscInt point, PetscInt n, PetscInt values[])
1576: {
1577:   ALE::Obj<PETSC_MESH_TYPE>                   m;
1578:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> s;
1579:   PetscErrorCode                              ierr;

1582:   DMMeshGetMesh(dm, m);
1583:   SectionIntGetSection(section, s);
1584:   m->restrictClosure(s, point, values, n);
1585:   return(0);
1586: }

1590: /*@C
1591:   SectionIntUpdateClosure - Updates the values in a given closure from the array

1593:   Not Collective

1595:   Input Parameters:
1596: + section - The section
1597: . mesh    - The DMMesh object
1598: . point   - The sieve point
1599: . array   - The array to fill up
1600: - mode    - The insertion mode

1602:   Output Parameter:
1603: . array - The array full of values in the closure

1605:   Level: intermediate

1607: .keywords: mesh, elements
1608: .seealso: DMMeshCreate()
1609: @*/
1610: PetscErrorCode SectionIntUpdateClosure(SectionInt section, DM dm, PetscInt point, PetscInt values[], InsertMode mode)
1611: {
1612:   ALE::Obj<PETSC_MESH_TYPE>                   m;
1613:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> s;
1614:   PetscErrorCode                              ierr;

1617:   DMMeshGetMesh(dm, m);
1618:   SectionIntGetSection(section, s);
1619:   if (mode == INSERT_VALUES) {
1620:     m->update(s, point, values);
1621:   } else if (mode == ADD_VALUES) {
1622:     m->updateAdd(s, point, values);
1623:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Invalid insertion mode: %d", mode);
1624:   return(0);
1625: }

1629: /*@
1630:   SectionIntComplete - Exchanges data across the mesh overlap.

1632:   Not collective

1634:   Input Parameter:
1635: . section - the section object

1637:   Level: advanced

1639: .seealso SectionIntRestrict(), SectionIntCreate(), SectionIntView()
1640: @*/
1641: PetscErrorCode SectionIntComplete(SectionInt section)
1642: {
1643:   Obj<PETSC_MESH_TYPE::int_section_type> s;
1644:   Obj<PETSC_MESH_TYPE>                   b;
1645:   PetscErrorCode                         ierr;

1648:   SectionIntGetSection(section, s);
1649:   SectionIntGetBundle(section, b);
1650: #if 0
1651:   ALE::Distribution<PETSC_MESH_TYPE>::completeSection(b, s);
1652: #else
1653:   ALE::Completion::completeSectionAdd(b->getSendOverlap(), b->getRecvOverlap(), s, s);
1654: #endif
1655:   return(0);
1656: }

1660: /*@
1661:   SectionIntZero - Zero out the entries

1663:   Not collective

1665:   Input Parameter:
1666: . section - the section object

1668:   Level: advanced

1670: .seealso SectionIntRestrict(), SectionIntCreate(), SectionIntView()
1671: @*/
1672: PetscErrorCode SectionIntZero(SectionInt section)
1673: {
1674:   Obj<PETSC_MESH_TYPE::int_section_type> s;
1675:   PetscErrorCode                         ierr;

1678:   SectionIntGetSection(section, s);
1679:   s->zero();
1680:   return(0);
1681: }

1685: /*@
1686:   SectionIntGetFiberDimension - Get the size of the vector space attached to the point

1688:   Not collective

1690:   Input Parameters:
1691: + section - the section object
1692: - point - the Sieve point

1694:   Output Parameters:
1695: . size - The fiber dimension

1697:   Level: advanced

1699: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
1700: @*/
1701: PetscErrorCode  SectionIntGetFiberDimension(SectionInt section, PetscInt point, PetscInt *size)
1702: {
1705:   *size = section->s->getFiberDimension(point);
1706:   return(0);
1707: }

1711: /*@
1712:   SectionIntSetFiberDimension - Set the size of the vector space attached to the point

1714:   Not collective

1716:   Input Parameters:
1717: + section - the section object
1718: . point - the Sieve point
1719: - size - The fiber dimension

1721:   Level: advanced

1723: .seealso SectionIntRestrict(), SectionIntCreate(), SectionIntView()
1724: @*/
1725: PetscErrorCode  SectionIntSetFiberDimension(SectionInt section, PetscInt point, const PetscInt size)
1726: {
1729:   section->s->setFiberDimension(point, size);
1730:   return(0);
1731: }

1735: /*@
1736:   SectionIntSetFiberDimensionField - Set the size of the vector space attached to the point for a given field

1738:   Not collective

1740:   Input Parameters:
1741: + section - the section object
1742: . point - the Sieve point
1743: . size - The fiber dimension
1744: - field - The field number

1746:   Level: advanced

1748: .seealso SectionIntSetFiberDimension(), SectionIntRestrict(), SectionIntCreate(), SectionIntView()
1749: @*/
1750: PetscErrorCode  SectionIntSetFiberDimensionField(SectionInt section, PetscInt point, const PetscInt size, const PetscInt field)
1751: {
1754:   section->s->setFiberDimension(point, size, field);
1755:   return(0);
1756: }

1760: /*@
1761:   SectionIntGetSize - Gets the number of local dofs in this Section

1763:   Not collective

1765:   Input Parameter:
1766: . section - the section object

1768:   Output Parameter:
1769: . size - the section size

1771:   Level: advanced

1773: .seealso SectionIntRestrict(), SectionIntCreate(), SectionIntView()
1774: @*/
1775: PetscErrorCode SectionIntGetSize(SectionInt section, PetscInt *size)
1776: {
1777:   Obj<PETSC_MESH_TYPE::int_section_type> s;
1778:   PetscErrorCode                         ierr;

1782:   SectionIntGetSection(section, s);
1783:   *size = s->size();
1784:   return(0);
1785: }

1789: /*@
1790:   SectionIntAllocate - Allocate storage for this section

1792:   Not collective

1794:   Input Parameter:
1795: . section - the section object

1797:   Level: advanced

1799: .seealso SectionIntRestrict(), SectionIntCreate(), SectionIntView()
1800: @*/
1801: PetscErrorCode  SectionIntAllocate(SectionInt section)
1802: {
1805:   section->b->allocate(section->s);
1806:   return(0);
1807: }

1811: /*@C
1812:   SectionIntClear - Dellocate storage for this section

1814:   Not collective

1816:   Input Parameter:
1817: . section - the section object

1819:   Level: advanced

1821: .seealso SectionIntRestrict(), SectionIntCreate(), SectionIntView()
1822: @*/
1823: PetscErrorCode  SectionIntClear(SectionInt section)
1824: {
1827:   section->s->clear();
1828:   return(0);
1829: }

1833: /*@
1834:   SectionIntSet - Sets all the values to the given value

1836:   Not collective

1838:   Input Parameters:
1839: + section - the real Section
1840: - val - the value

1842:   Level: intermediate

1844: .seealso VecNorm(), SectionIntCreate()
1845: @*/
1846: PetscErrorCode  SectionIntSet(SectionInt section, PetscInt val)
1847: {
1848:   Obj<PETSC_MESH_TYPE::int_section_type> s;
1849:   PetscErrorCode                         ierr;

1852:   SectionIntGetSection(section, s);
1853:   s->set(val);
1854:   return(0);
1855: }

1859: /*@
1860:   SectionIntAddSpace - Add another field to this section

1862:   Collective on DMMesh

1864:   Input Parameter:
1865: . section - the Section

1867:   Level: advanced

1869: .seealso SectionIntCreate(), SectionIntGetFibration()
1870: @*/
1871: PetscErrorCode  SectionIntAddSpace(SectionInt section)
1872: {
1873:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> s;
1874:   PetscErrorCode                              ierr;

1877:   SectionIntGetSection(section, s);
1878:   s->addSpace();
1879:   return(0);
1880: }

1884: /*@
1885:   SectionIntGetFibration - Creates a section for only the data associated with the given field

1887:   Collective on DMMesh

1889:   Input Parameter:
1890: + section - the Section
1891: - field- The field number

1893:   Output Parameter:
1894: . subsection - the section of the given field

1896:   Level: advanced

1898: .seealso SectionIntCreate(), SectionIntAddSpace()
1899: @*/
1900: PetscErrorCode  SectionIntGetFibration(SectionInt section, const PetscInt field, SectionInt *subsection)
1901: {
1902:   ALE::Obj<PETSC_MESH_TYPE>                   b;
1903:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> s;
1904:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> t;
1905:   MPI_Comm                                    comm;
1906:   PetscErrorCode                              ierr;

1909:   PetscObjectGetComm((PetscObject) section, &comm);
1910:   SectionIntGetBundle(section, b);
1911:   SectionIntGetSection(section, s);
1912:   SectionIntCreate(comm, subsection);
1913:   SectionIntSetBundle(*subsection, b);
1914:   t    = s->getFibration(field);
1915:   SectionIntSetSection(*subsection, t);
1916:   return(0);
1917: }

1921: /*@C
1922:   DMMeshGetVertexSectionInt - Create a Section over the vertices with the specified fiber dimension

1924:   Collective on DMMesh

1926:   Input Parameters:
1927: + mesh - The DMMesh object
1928: - fiberDim - The section name

1930:   Output Parameter:
1931: . section - The section

1933:   Level: intermediate

1935: .keywords: mesh, section, vertex
1936: .seealso: DMMeshCreate(), SectionIntCreate()
1937: @*/
1938: PetscErrorCode DMMeshGetVertexSectionInt(DM dm, const char name[], PetscInt fiberDim, SectionInt *section)
1939: {
1940:   ALE::Obj<PETSC_MESH_TYPE>                   m;
1941:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> s;
1942:   PetscErrorCode                              ierr;

1945:   DMMeshGetMesh(dm, m);
1946:   SectionIntCreate(m->comm(), section);
1947:   PetscObjectSetName((PetscObject) *section, name);
1948:   SectionIntSetBundle(*section, m);
1949:   SectionIntGetSection(*section, s);
1950: #if defined(PETSC_OPT_SIEVE)
1951:   s->setChart(m->getSieve()->getChart());
1952: #endif
1953:   s->setName(name);
1954:   s->setFiberDimension(m->depthStratum(0), fiberDim);
1955:   m->allocate(s);
1956:   return(0);
1957: }

1961: /*@C
1962:   DMMeshGetCellSectionInt - Create a Section over the cells with the specified fiber dimension

1964:   Collective on DMMesh

1966:   Input Parameters:
1967: + mesh - The DMMesh object
1968: - fiberDim - The section name

1970:   Output Parameter:
1971: . section - The section

1973:   Level: intermediate

1975: .keywords: mesh, section, cell
1976: .seealso: DMMeshCreate(), SectionIntCreate()
1977: @*/
1978: PetscErrorCode DMMeshGetCellSectionInt(DM dm, const char name[], PetscInt fiberDim, SectionInt *section)
1979: {
1980:   ALE::Obj<PETSC_MESH_TYPE>                   m;
1981:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> s;
1982:   PetscErrorCode                              ierr;

1985:   DMMeshGetMesh(dm, m);
1986:   SectionIntCreate(m->comm(), section);
1987:   PetscObjectSetName((PetscObject) *section, name);
1988:   SectionIntSetBundle(*section, m);
1989:   SectionIntGetSection(*section, s);
1990: #if defined(PETSC_OPT_SIEVE)
1991:   s->setChart(m->getSieve()->getChart());
1992: #endif
1993:   s->setName(name);
1994:   s->setFiberDimension(m->heightStratum(0), fiberDim);
1995:   m->allocate(s);
1996:   return(0);
1997: }

2001: /*@C
2002:   DMMeshSetupSection - Layout Section based upon discretization and boundary condition information in the Mesh

2004:   Collective on DMMesh

2006:   Input Parameters:
2007: . mesh - The DMMesh object

2009:   Output Parameter:
2010: . section - The section

2012:   Level: intermediate

2014: .keywords: mesh, section, cell
2015: .seealso: DMMeshCreate(), SectionRealCreate()
2016: @*/
2017: PetscErrorCode DMMeshSetupSection(DM dm, SectionReal section)
2018: {
2019:   ALE::Obj<PETSC_MESH_TYPE>                    m;
2020:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
2021:   PetscErrorCode                               ierr;

2024:   DMMeshGetMesh(dm, m);
2025:   SectionRealGetSection(section, s);
2026:   try {
2027:     m->setupField(s);
2028:   } catch(ALE::Exception e) {
2029:     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "%s", e.message());
2030:   }
2031:   return(0);
2032: }