Actual source code: section.c

petsc-3.3-p7 2013-05-11
  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:   } else if (isbinary) {
 34:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Binary viewer not implemented for Section");
 35:   } else if (isdraw){
 36:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Draw viewer not implemented for Section");
 37:   } else {
 38:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by this section object", ((PetscObject)viewer)->type_name);
 39:   }
 40:   return(0);
 41: }

 45: /*@C
 46:    SectionRealView - Views a Section object. 

 48:    Collective on Section

 50:    Input Parameters:
 51: +  section - the Section
 52: -  viewer - an optional visualization context

 54:    Notes:
 55:    The available visualization contexts include
 56: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
 57: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
 58:          output where only the first processor opens
 59:          the file.  All other processors send their 
 60:          data to the first processor to print. 

 62:    You can change the format the section is printed using the 
 63:    option PetscViewerSetFormat().

 65:    The user can open alternative visualization contexts with
 66: +    PetscViewerASCIIOpen() - Outputs section to a specified file
 67: .    PetscViewerBinaryOpen() - Outputs section in binary to a
 68:          specified file; corresponding input uses SectionLoad()
 69: .    PetscViewerDrawOpen() - Outputs section to an X window display

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

 77:    Level: beginner

 79:    Concepts: section^printing
 80:    Concepts: section^saving to disk

 82: .seealso: VecView(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerBinaryOpen(), PetscViewerCreate()
 83: @*/
 84: PetscErrorCode SectionRealView(SectionReal section, PetscViewer viewer)
 85: {

 91:   if (!viewer) {
 92:     PetscViewerASCIIGetStdout(((PetscObject)section)->comm,&viewer);
 93:   }

 97:   PetscLogEventBegin(SectionReal_View,0,0,0,0);
 98:   (*section->ops->view)(section, viewer);
 99:   PetscLogEventEnd(SectionReal_View,0,0,0,0);
100:   return(0);
101: }

105: /*@C
106:   SectionRealDuplicate - Create an equivalent Section object

108:   Not collective

110:   Input Parameter:
111: . section - the section object

113:   Output Parameter:
114: . newSection - the duplicate
115:  
116:   Level: advanced

118: .seealso SectionRealCreate(), SectionRealSetSection()
119: @*/
120: PetscErrorCode  SectionRealDuplicate(SectionReal section, SectionReal *newSection)
121: {

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

130:   t->setAtlas(s->getAtlas());
131:   t->allocateStorage();
132:   t->copyBC(s);
133:   SectionRealCreate(s->comm(), newSection);
134:   SectionRealSetSection(*newSection, t);
135:   SectionRealSetBundle(*newSection, section->b);
136:   return(0);
137: }

141: /*@C
142:   SectionRealGetSection - Gets the internal section object

144:   Not collective

146:   Input Parameter:
147: . section - the section object

149:   Output Parameter:
150: . s - the internal section object
151:  
152:   Level: advanced

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

166: /*@C
167:   SectionRealSetSection - Sets the internal section object

169:   Not collective

171:   Input Parameters:
172: + section - the section object
173: - s - the internal section object
174:  
175:   Level: advanced

177: .seealso SectionRealCreate(), SectionRealGetSection()
178: @*/
179: PetscErrorCode  SectionRealSetSection(SectionReal section, const ALE::Obj<PETSC_MESH_TYPE::real_section_type>& s)
180: {

185:   if (!s.isNull()) {PetscObjectSetName((PetscObject) section, s->getName().c_str());}
186:   section->s = s;
187:   return(0);
188: }

192: /*@C
193:   SectionRealGetBundle - Gets the section bundle

195:   Not collective

197:   Input Parameter:
198: . section - the section object

200:   Output Parameter:
201: . b - the section bundle
202:  
203:   Level: advanced

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

217: /*@C
218:   SectionRealSetBundle - Sets the section bundle

220:   Not collective

222:   Input Parameters:
223: + section - the section object
224: - b - the section bundle
225:  
226:   Level: advanced

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

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

244:   Collective on MPI_Comm

246:   Input Parameter:
247: . comm - the processors that will share the global section

249:   Output Parameters:
250: . section - the section object

252:   Level: advanced

254: .seealso SectionRealDestroy(), SectionRealView()
255: @*/
256: PetscErrorCode  SectionRealCreate(MPI_Comm comm, SectionReal *section)
257: {
259:   SectionReal    s;

263:   *section = PETSC_NULL;

265:   PetscHeaderCreate(s,_p_SectionReal,struct _SectionRealOps,SECTIONREAL_CLASSID,0,"SectionReal","Section","DM",comm,SectionRealDestroy,0);
266:   s->ops->view     = SectionRealView_Sieve;
267:   s->ops->restrictClosure = SectionRealRestrict;
268:   s->ops->update   = SectionRealUpdate;

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

272:   new(&s->s) ALE::Obj<PETSC_MESH_TYPE::real_section_type>(PETSC_MESH_TYPE::real_section_type(comm));
273:   new(&s->b) ALE::Obj<PETSC_MESH_TYPE>(PETSC_NULL);
274:   *section = s;
275:   return(0);
276: }

280: /*@
281:   SectionRealDestroy - Destroys a section.

283:   Collective on Section

285:   Input Parameter:
286: . section - the section object

288:   Level: advanced

290: .seealso SectionRealCreate(), SectionRealView()
291: @*/
292: PetscErrorCode  SectionRealDestroy(SectionReal *section)
293: {

297:   if (!*section) return(0);
299:   if (--((PetscObject)(*section))->refct > 0) return(0);
300:   PetscHeaderDestroy(section);
301:   return(0);
302: }

306: /*@
307:   SectionRealDistribute - Distributes the sections.

309:   Not Collective

311:   Input Parameters:
312: + serialSection - The original Section object
313: - parallelMesh - The parallel DMMesh

315:   Output Parameter:
316: . parallelSection - The distributed Section object

318:   Level: intermediate

320: .keywords: mesh, section, distribute
321: .seealso: DMMeshCreate()
322: @*/
323: PetscErrorCode SectionRealDistribute(SectionReal serialSection, DM parallelMesh, SectionReal *parallelSection)
324: {
325:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> oldSection;
326:   ALE::Obj<PETSC_MESH_TYPE>               m;

330:   SectionRealGetSection(serialSection, oldSection);
331:   DMMeshGetMesh(parallelMesh, m);
332:   SectionRealCreate(oldSection->comm(), parallelSection);
333: #ifdef PETSC_OPT_SIEVE
334:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> newSection;

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

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

354:   Not collective

356:   Input Parameters:
357: + section - the section object
358: - point - the Sieve point

360:   Output Parameter:
361: . values - The values associated with the submesh

363:   Level: advanced

365: .seealso SectionUpdate(), SectionCreate(), SectionView()
366: @*/
367: PetscErrorCode  SectionRealRestrict(SectionReal section, PetscInt point, PetscScalar *values[])
368: {
372:   try {
373:     *values = (PetscScalar *) section->s->restrictPoint(point);
374:   } catch(ALE::Exception e) {
375:     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid argument: %s", e.message());
376:   }
377:   return(0);
378: }

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

385:   Not collective

387:   Input Parameters:
388: + section - the section object
389: . point - the Sieve point
390: . values - The values associated with the submesh
391: - mode - The insertion mode

393:   Level: advanced

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

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

425:   Not Collective

427:   Input Parameters:
428: + section - The section
429: . mesh    - The DMMesh object
430: - point   - The sieve point

432:   Output Parameter:
433: . array - The array full of values in the closure

435:   Level: intermediate

437: .keywords: mesh, elements
438: .seealso: DMMeshCreate()
439: @*/
440: PetscErrorCode SectionRealRestrictClosure(SectionReal section, DM dm, PetscInt point, const PetscScalar *values[])
441: {
442:   ALE::Obj<PETSC_MESH_TYPE> m;
443:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
444:   PetscErrorCode            ierr;

447:   DMMeshGetMesh(dm, m);
448:   SectionRealGetSection(section, s);
449: #ifdef PETSC_USE_COMPLEX
450:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "SectionReal does not support complex restriction");
451: #else
452:   try {
453:     *values = m->restrictClosure(s, point);
454:   } catch(ALE::Exception e) {
455:     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid argument: %s", e.message());
456:   }
457: #endif
458:   return(0);
459: }

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

466:   Not Collective

468:   Input Parameters:
469: + section - The section
470: . mesh    - The DMMesh object
471: . point   - The sieve point
472: . n       - The array size
473: - array   - The array to fill up

475:   Output Parameter:
476: . array - The array full of values in the closure

478:   Level: intermediate

480: .keywords: mesh, elements
481: .seealso: DMMeshCreate()
482: @*/
483: PetscErrorCode SectionRealRestrictClosure(SectionReal section, DM dm, PetscInt point, PetscInt n, PetscScalar values[])
484: {
485:   ALE::Obj<PETSC_MESH_TYPE> m;
486:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
487:   PetscErrorCode            ierr;

490:   DMMeshGetMesh(dm, m);
491:   SectionRealGetSection(section, s);
492: #ifdef PETSC_USE_COMPLEX
493:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "SectionReal does not support complex restriction");
494: #else
495:   try {
496:     m->restrictClosure(s, point, values, n);
497:   } catch(ALE::Exception e) {
498:     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid argument: %s", e.message());
499:   }
500: #endif
501:   return(0);
502: }

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

509:   Not Collective

511:   Input Parameters:
512: + section - The section
513: . mesh    - The DMMesh object
514: . point   - The sieve point
515: . array   - The array to fill up
516: - mode    - The insertion mode

518:   Output Parameter:
519: . array - The array full of values in the closure

521:   Level: intermediate

523: .keywords: mesh, elements
524: .seealso: DMMeshCreate()
525: @*/
526: PetscErrorCode SectionRealUpdateClosure(SectionReal section, DM dm, PetscInt point, PetscScalar values[], InsertMode mode)
527: {
528:   ALE::Obj<PETSC_MESH_TYPE> m;
529:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
530:   PetscErrorCode            ierr;

533:   DMMeshGetMesh(dm, m);
534:   SectionRealGetSection(section, s);
535: #ifdef PETSC_USE_COMPLEX
536:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "SectionReal does not support complex update");
537: #else
538:   try {
539:     if (mode == INSERT_VALUES) {
540:       m->update(s, point, values);
541:     } else if (mode == ADD_VALUES) {
542:       m->updateAdd(s, point, values);
543:     } else {
544:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Invalid insertion mode: %d", mode);
545:     }
546:   } catch(ALE::Exception e) {
547:     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid argument: %s", e.message());
548:   }
549: #endif
550:   return(0);
551: }

555: /*@
556:   SectionRealComplete - Exchanges data across the mesh overlap.

558:   Not collective

560:   Input Parameter:
561: . section - the section object

563:   Level: advanced

565: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
566: @*/
567: PetscErrorCode SectionRealComplete(SectionReal section)
568: {
569:   Obj<PETSC_MESH_TYPE::real_section_type> s;
570:   Obj<PETSC_MESH_TYPE>                    b;

574:   SectionRealGetSection(section, s);
575:   SectionRealGetBundle(section, b);
576: #if 0
577:   ALE::Distribution<PETSC_MESH_TYPE>::completeSection(b, s);
578: #else
579:   ALE::Completion::completeSectionAdd(b->getSendOverlap(), b->getRecvOverlap(), s, s);
580: #endif
581:   return(0);
582: }

586: /*@
587:   SectionRealZero - Zero out the entries

589:   Not collective

591:   Input Parameter:
592: . section - the section object

594:   Level: advanced

596: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
597: @*/
598: PetscErrorCode SectionRealZero(SectionReal section)
599: {
600:   Obj<PETSC_MESH_TYPE::real_section_type> s;

604:   SectionRealGetSection(section, s);
605:   s->zero();
606:   return(0);
607: }

611: /*@
612:   SectionRealGetFiberDimension - Get the size of the vector space attached to the point

614:   Not collective

616:   Input Parameters:
617: + section - the section object
618: - point - the Sieve point

620:   Output Parameters:
621: . size - The fiber dimension

623:   Level: advanced

625: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
626: @*/
627: PetscErrorCode  SectionRealGetFiberDimension(SectionReal section, PetscInt point, PetscInt *size)
628: {
631:   *size = section->s->getFiberDimension(point);
632:   return(0);
633: }

637: /*@
638:   SectionRealSetFiberDimension - Set the size of the vector space attached to the point

640:   Not collective

642:   Input Parameters:
643: + section - the section object
644: . point - the Sieve point
645: - size - The fiber dimension

647:   Level: advanced

649: .seealso SectionRealSetFiberDimensionField(), SectionRealRestrict(), SectionRealCreate(), SectionRealView()
650: @*/
651: PetscErrorCode  SectionRealSetFiberDimension(SectionReal section, PetscInt point, const PetscInt size)
652: {
655:   section->s->setFiberDimension(point, size);
656:   return(0);
657: }

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

664:   Not collective

666:   Input Parameters:
667: + section - the section object
668: . point - the Sieve point
669: . size - The fiber dimension
670: - field - The field number

672:   Level: advanced

674: .seealso SectionRealSetFiberDimension(), SectionRealRestrict(), SectionRealCreate(), SectionRealView()
675: @*/
676: PetscErrorCode  SectionRealSetFiberDimensionField(SectionReal section, PetscInt point, const PetscInt size, const PetscInt field)
677: {
680:   section->s->setFiberDimension(point, size, field);
681:   return(0);
682: }

686: /*@
687:   SectionRealGetSize - Gets the number of local dofs in this Section

689:   Not collective

691:   Input Parameter:
692: . section - the section object

694:   Output Parameter:
695: . size - the section size

697:   Level: advanced

699: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
700: @*/
701: PetscErrorCode SectionRealGetSize(SectionReal section, PetscInt *size)
702: {
703:   Obj<PETSC_MESH_TYPE::real_section_type> s;

708:   SectionRealGetSection(section, s);
709:   *size = s->size();
710:   return(0);
711: }

715: /*@
716:   SectionRealAllocate - Allocate storage for this section

718:   Not collective

720:   Input Parameter:
721: . section - the section object

723:   Level: advanced

725: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
726: @*/
727: PetscErrorCode  SectionRealAllocate(SectionReal section)
728: {
731:   section->b->allocate(section->s);
732:   return(0);
733: }

737: /*@
738:   SectionRealCreateLocalVector - Creates a vector with the local piece of the Section

740:   Collective on DMMesh

742:   Input Parameter:
743: . section - the Section  

745:   Output Parameter:
746: . localVec - the local vector

748:   Level: advanced

750:   Notes: The vector can safely be destroyed using VecDestroy().
751: .seealso DMMeshDestroy(), DMMeshCreate()
752: @*/
753: PetscErrorCode  SectionRealCreateLocalVector(SectionReal section, Vec *localVec)
754: {
755:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;

759: #ifdef PETSC_USE_COMPLEX
760:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "SectionReal does not support complex Vec");
761: #else
762:   SectionRealGetSection(section, s);
763:   VecCreateSeqWithArray(PETSC_COMM_SELF,1, s->getStorageSize(), s->restrictSpace(), localVec);
764: #endif
765:   return(0);
766: }

770: /*@
771:   SectionRealAddSpace - Add another field to this section

773:   Collective on DMMesh

775:   Input Parameter:
776: . section - the Section

778:   Level: advanced

780: .seealso SectionRealCreate(), SectionRealGetFibration()
781: @*/
782: PetscErrorCode  SectionRealAddSpace(SectionReal section)
783: {
784:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;

788:   SectionRealGetSection(section, s);
789:   s->addSpace();
790:   return(0);
791: }

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

798:   Collective on DMMesh

800:   Input Parameter:
801: + section - the Section
802: - field- The field number

804:   Output Parameter:
805: . subsection - the section of the given field

807:   Level: advanced

809: .seealso SectionRealCreate()
810: @*/
811: PetscErrorCode  SectionRealGetFibration(SectionReal section, const PetscInt field, SectionReal *subsection)
812: {
813:   ALE::Obj<PETSC_MESH_TYPE>                    b;
814:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
815:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> t;
816:   MPI_Comm       comm;

820:   PetscObjectGetComm((PetscObject) section, &comm);
821:   SectionRealGetBundle(section, b);
822:   SectionRealGetSection(section, s);
823:   SectionRealCreate(comm, subsection);
824:   SectionRealSetBundle(*subsection, b);
825:   t    = s->getFibration(field);
826:   SectionRealSetSection(*subsection, t);
827:   return(0);
828: }

832: /*@C
833:   SectionRealToVec - Maps the given section to a Vec

835:   Collective on Section

837:   Input Parameters:
838: + section - the real Section
839: - mesh - The DMMesh

841:   Output Parameter:
842: . vec - the Vec

844:   Level: intermediate

846: .seealso VecCreate(), SectionRealCreate()
847: @*/
848: PetscErrorCode  SectionRealToVec(SectionReal section, DM dm, ScatterMode mode, Vec vec)
849: {
850:   Vec            localVec;
851:   VecScatter     scatter;

856:   SectionRealCreateLocalVector(section, &localVec);
857:   DMMeshGetGlobalScatter(dm, &scatter);
858:   if (mode == SCATTER_FORWARD) {
859:     VecScatterBegin(scatter, localVec, vec, INSERT_VALUES, mode);
860:     VecScatterEnd(scatter, localVec, vec, INSERT_VALUES, mode);
861:   } else {
862:     VecScatterBegin(scatter, vec, localVec, INSERT_VALUES, mode);
863:     VecScatterEnd(scatter, vec, localVec, INSERT_VALUES, mode);
864:   }
865:   VecDestroy(&localVec);
866:   return(0);
867: }

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

874:   Collective on VecScatter

876:   Input Parameters:
877: + section - the Section
878: . scatter - the scatter from the Section to the Vec
879: . mode - the mode, SCATTER_FORWARD (Section to Vec) or SCATTER_REVERSE (Vec to Section)
880: - vec - the Vec

882:   Level: advanced

884: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
885: @*/
886: PetscErrorCode  SectionRealToVec(SectionReal section, VecScatter scatter, ScatterMode mode, Vec vec)
887: {
888:   Vec            localVec;

893:   SectionRealCreateLocalVector(section, &localVec);
894:   if (mode == SCATTER_FORWARD) {
895:     VecScatterBegin(scatter, localVec, vec, INSERT_VALUES, mode);
896:     VecScatterEnd(scatter, localVec, vec, INSERT_VALUES, mode);
897:   } else {
898:     VecScatterBegin(scatter, vec, localVec, INSERT_VALUES, mode);
899:     VecScatterEnd(scatter, vec, localVec, INSERT_VALUES, mode);
900:   }
901:   VecDestroy(&localVec);
902:   return(0);
903: }

907: /*@
908:   SectionRealClear - Dellocate storage for this section

910:   Not collective

912:   Input Parameter:
913: . section - the section object

915:   Level: advanced

917: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
918: @*/
919: PetscErrorCode  SectionRealClear(SectionReal section)
920: {
923:   section->s->clear();
924:   return(0);
925: }

929: /*@
930:   SectionRealSet - Sets all the values to the given value

932:   Not collective

934:   Input Parameters:
935: + section - the real Section
936: - val - the value

938:   Level: intermediate

940: .seealso VecNorm(), SectionRealCreate()
941: @*/
942: PetscErrorCode  SectionRealSet(SectionReal section, PetscReal val)
943: {
944:   Obj<PETSC_MESH_TYPE::real_section_type> s;

948:   SectionRealGetSection(section, s);
949:   s->set(val);
950:   return(0);
951: }

955: /*@C
956:   SectionRealNorm - Computes the vector norm.

958:   Collective on Section

960:   Input Parameters:
961: +  section - the real Section
962: -  type - one of NORM_1, NORM_2, NORM_INFINITY.  Also available
963:           NORM_1_AND_2, which computes both norms and stores them
964:           in a two element array.

966:   Output Parameter:
967: . val - the norm

969:   Notes:
970: $     NORM_1 denotes sum_i |x_i|
971: $     NORM_2 denotes sqrt(sum_i (x_i)^2)
972: $     NORM_INFINITY denotes max_i |x_i|

974:   Level: intermediate

976: .seealso VecNorm(), SectionRealCreate()
977: @*/
978: PetscErrorCode  SectionRealNorm(SectionReal section, DM dm, NormType type, PetscReal *val)
979: {
980:   Obj<PETSC_MESH_TYPE> m;
981:   Obj<PETSC_MESH_TYPE::real_section_type> s;
982:   Vec            v;

987:   DMMeshGetMesh(dm, m);
988:   SectionRealGetSection(section, s);
989:   const ALE::Obj<PETSC_MESH_TYPE::order_type>& order = m->getFactory()->getGlobalOrder(m, s->getName(), s);
990:   VecCreate(m->comm(), &v);
991:   VecSetSizes(v, order->getLocalSize(), order->getGlobalSize());
992:   VecSetFromOptions(v);
993:   SectionRealToVec(section, dm, SCATTER_FORWARD, v);
994:   VecNorm(v, type, val);
995:   VecDestroy(&v);
996:   return(0);
997: }

1001: /*@
1002:   SectionRealAXPY - 

1004:   Collective on Section

1006:   Input Parameters:
1007: +  section - the real Section
1008: .  alpha - a scalar
1009: -  X - the other real Section

1011:   Output Parameter:
1012: . section - the difference

1014:   Level: intermediate

1016: .seealso VecNorm(), SectionRealCreate()
1017: @*/
1018: PetscErrorCode  SectionRealAXPY(SectionReal section, DM dm, PetscScalar alpha, SectionReal X)
1019: {
1020:   Obj<PETSC_MESH_TYPE> m;
1021:   Obj<PETSC_MESH_TYPE::real_section_type> s;
1022:   Obj<PETSC_MESH_TYPE::real_section_type> sX;
1023:   Vec            v, x;

1028:   DMMeshGetMesh(dm, m);
1029:   SectionRealGetSection(section, s);
1030:   SectionRealGetSection(X, sX);
1031:   const ALE::Obj<PETSC_MESH_TYPE::order_type>& order = m->getFactory()->getGlobalOrder(m, s->getName(), s);
1032:   VecCreate(m->comm(), &v);
1033:   VecSetSizes(v, order->getLocalSize(), order->getGlobalSize());
1034:   VecSetFromOptions(v);
1035:   VecDuplicate(v, &x);
1036:   SectionRealToVec(section, dm, SCATTER_FORWARD, v);
1037:   SectionRealToVec(X,       dm, SCATTER_FORWARD, x);
1038:   VecAXPY(v, alpha, x);
1039:   SectionRealToVec(section, dm, SCATTER_REVERSE, v);
1040:   VecDestroy(&v);
1041:   VecDestroy(&x);
1042:   return(0);
1043: }

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

1050:   Collective on DMMesh

1052:   Input Parameters:
1053: + mesh - The DMMesh object
1054: . name - The name of the section
1055: - fiberDim - The number of degrees of freedom per vertex

1057:   Output Parameter:
1058: . section - The section

1060:   Level: intermediate

1062: .keywords: mesh, section, vertex
1063: .seealso: DMMeshCreate(), SectionRealCreate()
1064: @*/
1065: PetscErrorCode DMMeshGetVertexSectionReal(DM dm, const char name[], PetscInt fiberDim, SectionReal *section)
1066: {
1067:   ALE::Obj<PETSC_MESH_TYPE> m;
1068:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
1069:   PetscErrorCode      ierr;

1072:   DMMeshGetMesh(dm, m);
1073:   SectionRealCreate(m->comm(), section);
1074:   PetscObjectSetName((PetscObject) *section, name);
1075:   SectionRealSetBundle(*section, m);
1076:   SectionRealGetSection(*section, s);
1077:   s->setChart(m->getSieve()->getChart());
1078:   s->setName(name);
1079:   s->setFiberDimension(m->depthStratum(0), fiberDim);
1080:   m->allocate(s);
1081:   return(0);
1082: }

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

1089:   Collective on DMMesh

1091:   Input Parameters:
1092: + mesh - The DMMesh object
1093: . name - The name of the section
1094: - fiberDim - The number of degrees of freedom per cell

1096:   Output Parameter:
1097: . section - The section

1099:   Level: intermediate

1101: .keywords: mesh, section, cell
1102: .seealso: DMMeshCreate(), SectionRealCreate(), DMMeshGetVertexSectionReal(), DMMeshCreateSectionRealIS()
1103: @*/
1104: PetscErrorCode DMMeshGetCellSectionReal(DM dm, const char name[], PetscInt fiberDim, SectionReal *section)
1105: {
1106:   ALE::Obj<PETSC_MESH_TYPE> m;
1107:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
1108:   PetscErrorCode      ierr;

1111:   DMMeshGetMesh(dm, m);
1112:   SectionRealCreate(m->comm(), section);
1113:   PetscObjectSetName((PetscObject) *section, name);
1114:   SectionRealSetBundle(*section, m);
1115:   SectionRealGetSection(*section, s);
1116:   s->setChart(m->getSieve()->getChart());
1117:   s->setName(name);
1118:   s->setFiberDimension(m->heightStratum(0), fiberDim);
1119:   m->allocate(s);
1120:   return(0);
1121: }

1125: /*@C
1126:   DMMeshCreateSectionRealIS - Create a Section over the points specified in an IS.
1127:   
1128:   Collective on DMMesh
1129:   
1130:   Input Parameters:
1131:   + dm - the DMMesh object
1132:   . is - The IS describing the points associated with the degrees of freedom
1133:   . name - The name of the section
1134:   - fiberDim - The number of degrees of freedom per point of the IS
1135:   
1136:   Output Parameter:
1137:   . section - The section
1138:   
1139:   Level: intermediate

1141: .keywords: mesh, section, cell
1142: .seealso: DMMeshCreate(), SectionRealCreate(), DMMeshGetVertexSectionReal(), DMMeshGetCellSectionReal()

1144: @*/
1145: PetscErrorCode DMMeshCreateSectionRealIS(DM dm,IS is,const char name[],PetscInt fiberDim,SectionReal *section)
1146: {
1147:   MPI_Comm        comm;
1148:   PetscSection    s;
1149:   PetscInt        pStart,pEnd;
1150:   const PetscInt *points;
1151:   PetscInt        numpoints,p;
1152:   PetscErrorCode  ierr;
1153: 
1155:   PetscObjectGetComm((PetscObject) dm,&comm);
1156: 
1157:   PetscSectionCreate(comm,&s);
1158:   DMMeshGetChart(dm,&pStart,&pEnd);
1159:   PetscSectionSetChart(s,pStart,pEnd);
1160: 
1161:   ISGetLocalSize(is,&numpoints);
1162:   ISGetIndices(is,&points);
1163:   for (p =0; p < numpoints; p++) {
1164:     PetscSectionSetDof(s,points[p],fiberDim);
1165:   }
1166:   ISRestoreIndices(is,&points);
1167:   DMMeshSetSection(dm,name,s);
1168:   PetscSectionDestroy(&s);
1169:   DMMeshGetSectionReal(dm,name,section);
1170:   return(0);
1171: }


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

1179:   Collective on DMMesh

1181:   Input Parameters:
1182: + mesh - the mesh object
1183: - section - The SectionReal

1185:   Output Parameters:
1186: . gvec - the global vector

1188:   Level: advanced

1190: .seealso DMMeshDestroy(), DMMeshCreate(), DMMeshCreateGlobalRealVector()
1191: @*/
1192: PetscErrorCode DMMeshCreateGlobalRealVector(DM dm, SectionReal section, Vec *gvec)
1193: {
1194:   ALE::Obj<PETSC_MESH_TYPE> m;
1195:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
1196:   const char    *name;

1200:   DMMeshGetMesh(dm, m);
1201:   SectionRealGetSection(section, s);
1202:   PetscObjectGetName((PetscObject) section, &name);
1203:   const ALE::Obj<PETSC_MESH_TYPE::order_type>& order = m->getFactory()->getGlobalOrder(m, name, s);

1205:   VecCreate(m->comm(), gvec);
1206:   VecSetSizes(*gvec, order->getLocalSize(), order->getGlobalSize());
1207:   VecSetFromOptions(*gvec);
1208:   return(0);
1209: }

1213: PetscErrorCode SectionIntView_Sieve(SectionInt section, PetscViewer viewer)
1214: {
1215:   PetscBool      iascii, isbinary, isdraw;

1219:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
1220:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERBINARY, &isbinary);
1221:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW, &isdraw);

1223:   if (iascii){
1224:     ALE::Obj<PETSC_MESH_TYPE::int_section_type> s;
1225:     ALE::Obj<PETSC_MESH_TYPE>                   b;
1226:     const char                                  *name;

1228:     SectionIntGetSection(section, s);
1229:     SectionIntGetBundle(section, b);
1230:     PetscObjectGetName((PetscObject) section, &name);
1231:     SectionView_Sieve_Ascii(b, s, name, viewer);
1232:   } else if (isbinary) {
1233:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Binary viewer not implemented for Section");
1234:   } else if (isdraw){
1235:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Draw viewer not implemented for Section");
1236:   } else {
1237:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by this section object", ((PetscObject)viewer)->type_name);
1238:   }
1239:   return(0);
1240: }

1244: /*@C
1245:    SectionIntView - Views a Section object. 

1247:    Collective on Section

1249:    Input Parameters:
1250: +  section - the Section
1251: -  viewer - an optional visualization context

1253:    Notes:
1254:    The available visualization contexts include
1255: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1256: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1257:          output where only the first processor opens
1258:          the file.  All other processors send their 
1259:          data to the first processor to print. 

1261:    You can change the format the section is printed using the 
1262:    option PetscViewerSetFormat().

1264:    The user can open alternative visualization contexts with
1265: +    PetscViewerASCIIOpen() - Outputs section to a specified file
1266: .    PetscViewerBinaryOpen() - Outputs section in binary to a
1267:          specified file; corresponding input uses SectionLoad()
1268: .    PetscViewerDrawOpen() - Outputs section to an X window display

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

1276:    Level: beginner

1278:    Concepts: section^printing
1279:    Concepts: section^saving to disk

1281: .seealso: VecView(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerBinaryOpen(), PetscViewerCreate()
1282: @*/
1283: PetscErrorCode SectionIntView(SectionInt section, PetscViewer viewer)
1284: {

1290:   if (!viewer) {
1291:     PetscViewerASCIIGetStdout(((PetscObject)section)->comm,&viewer);
1292:   }

1296:   PetscLogEventBegin(SectionInt_View,0,0,0,0);
1297:   (*section->ops->view)(section, viewer);
1298:   PetscLogEventEnd(SectionInt_View,0,0,0,0);
1299:   return(0);
1300: }

1304: /*@C
1305:   SectionIntGetSection - Gets the internal section object

1307:   Not collective

1309:   Input Parameter:
1310: . section - the section object

1312:   Output Parameter:
1313: . s - the internal section object
1314:  
1315:   Level: advanced

1317: .seealso SectionIntCreate(), SectionIntSetSection()
1318: @*/
1319: PetscErrorCode  SectionIntGetSection(SectionInt section, ALE::Obj<PETSC_MESH_TYPE::int_section_type>& s)
1320: {
1323:   s = section->s;
1324:   return(0);
1325: }

1329: /*@C
1330:   SectionIntSetSection - Sets the internal section object

1332:   Not collective

1334:   Input Parameters:
1335: + section - the section object
1336: - s - the internal section object
1337:  
1338:   Level: advanced

1340: .seealso SectionIntCreate(), SectionIntGetSection()
1341: @*/
1342: PetscErrorCode  SectionIntSetSection(SectionInt section, const ALE::Obj<PETSC_MESH_TYPE::int_section_type>& s)
1343: {

1348:   if (!s.isNull()) {PetscObjectSetName((PetscObject) section, s->getName().c_str());}
1349:   section->s = s;
1350:   return(0);
1351: }

1355: /*@C
1356:   SectionIntGetBundle - Gets the section bundle

1358:   Not collective

1360:   Input Parameter:
1361: . section - the section object

1363:   Output Parameter:
1364: . b - the section bundle
1365:  
1366:   Level: advanced

1368: .seealso SectionIntCreate(), SectionIntGetSection(), SectionIntSetSection()
1369: @*/
1370: PetscErrorCode  SectionIntGetBundle(SectionInt section, ALE::Obj<PETSC_MESH_TYPE>& b)
1371: {
1374:   b = section->b;
1375:   return(0);
1376: }

1380: /*@C
1381:   SectionIntSetBundle - Sets the section bundle

1383:   Not collective

1385:   Input Parameters:
1386: + section - the section object
1387: - b - the section bundle
1388:  
1389:   Level: advanced

1391: .seealso SectionIntCreate(), SectionIntGetSection(), SectionIntSetSection()
1392: @*/
1393: PetscErrorCode  SectionIntSetBundle(SectionInt section, const ALE::Obj<PETSC_MESH_TYPE>& b)
1394: {
1397:   section->b = b;
1398:   return(0);
1399: }

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

1407:   Collective on MPI_Comm

1409:   Input Parameter:
1410: . comm - the processors that will share the global section

1412:   Output Parameters:
1413: . section - the section object

1415:   Level: advanced

1417: .seealso SectionIntDestroy(), SectionIntView()
1418: @*/
1419: PetscErrorCode  SectionIntCreate(MPI_Comm comm, SectionInt *section)
1420: {
1422:   SectionInt    s;

1426:   *section = PETSC_NULL;

1428:   PetscHeaderCreate(s,_p_SectionInt,struct _SectionIntOps,SECTIONINT_CLASSID,0,"SectionInt","Section","DM",comm,SectionIntDestroy,0);
1429:   s->ops->view     = SectionIntView_Sieve;
1430:   s->ops->restrictClosure = SectionIntRestrict;
1431:   s->ops->update   = SectionIntUpdate;

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

1435:   new(&s->s) ALE::Obj<PETSC_MESH_TYPE::int_section_type>(PETSC_MESH_TYPE::int_section_type(comm));
1436:   new(&s->b) ALE::Obj<PETSC_MESH_TYPE>(PETSC_NULL);
1437:   *section = s;
1438:   return(0);
1439: }

1443: /*@
1444:   SectionIntDestroy - Destroys a section.

1446:   Collective on Section

1448:   Input Parameter:
1449: . section - the section object

1451:   Level: advanced

1453: .seealso SectionIntCreate(), SectionIntView()
1454: @*/
1455: PetscErrorCode  SectionIntDestroy(SectionInt *section)
1456: {

1460:   if (!*section) return(0);
1462:   if (--((PetscObject)(*section))->refct > 0) return(0);
1463:   PetscHeaderDestroy(section);
1464:   return(0);
1465: }

1469: /*@
1470:   SectionIntDistribute - Distributes the sections.

1472:   Not Collective

1474:   Input Parameters:
1475: + serialSection - The original Section object
1476: - parallelMesh - The parallel DMMesh

1478:   Output Parameter:
1479: . parallelSection - The distributed Section object

1481:   Level: intermediate

1483: .keywords: mesh, section, distribute
1484: .seealso: DMMeshCreate()
1485: @*/
1486: PetscErrorCode SectionIntDistribute(SectionInt serialSection, DM parallelMesh, SectionInt *parallelSection)
1487: {
1488:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> oldSection;
1489:   ALE::Obj<PETSC_MESH_TYPE> m;
1490:   PetscErrorCode      ierr;

1493:   SectionIntGetSection(serialSection, oldSection);
1494:   DMMeshGetMesh(parallelMesh, m);
1495:   SectionIntCreate(oldSection->comm(), parallelSection);
1496: #ifdef PETSC_OPT_SIEVE
1497:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> newSection;

1499:   // We assume all integer sections are complete sections
1500:   newSection->setName(oldSection->getName());
1501:   newSection->setChart(m->getSieve()->getChart());
1502:   //distributeSection(oldSection, partition, m->getRenumbering(), m->getDistSendOverlap(), m->getDistRecvOverlap(), newSection);
1503:   SectionIntSetSection(*parallelSection, newSection);
1504:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Not working because the partition is unavailable");
1505: #else
1506:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> newSection = ALE::Distribution<PETSC_MESH_TYPE>::distributeSection(oldSection, m, m->getDistSendOverlap(), m->getDistRecvOverlap());
1507:   SectionIntSetSection(*parallelSection, newSection);
1508: #endif
1509:   return(0);
1510: }

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

1517:   Not collective

1519:   Input Parameters:
1520: + section - the section object
1521: - point - the Sieve point

1523:   Output Parameter:
1524: . values - The values associated with the submesh

1526:   Level: advanced

1528: .seealso SectionIntUpdate(), SectionIntCreate(), SectionIntView()
1529: @*/
1530: PetscErrorCode  SectionIntRestrict(SectionInt section, PetscInt point, PetscInt *values[])
1531: {
1535:   *values = (PetscInt *) section->b->restrictClosure(section->s, point);
1536:   return(0);
1537: }

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

1544:   Not collective

1546:   Input Parameters:
1547: + section - the section object
1548: . point - the Sieve point
1549: . values - The values associated with the submesh
1550: - mode - The insertion mode

1552:   Level: advanced

1554: .seealso SectionIntRestrict(), SectionIntCreate(), SectionIntView()
1555: @*/
1556: PetscErrorCode  SectionIntUpdate(SectionInt section, PetscInt point, const PetscInt values[], InsertMode mode)
1557: {
1561:   if (mode == INSERT_VALUES) {
1562:     section->b->update(section->s, point, values);
1563:   } else if (mode == ADD_VALUES) {
1564:     section->b->updateAdd(section->s, point, values);
1565:   } else {
1566:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Invalid insertion mode: %d", mode);
1567:   }
1568:   return(0);
1569: }

1573: /*@C
1574:   SectionIntRestrictClosure - Returns an array with the values in a given closure

1576:   Not Collective

1578:   Input Parameters:
1579: + section - The section
1580: . mesh    - The DMMesh object
1581: . point   - The sieve point
1582: . n       - The array size
1583: - array   - The array to fill up

1585:   Output Parameter:
1586: . array - The array full of values in the closure

1588:   Level: intermediate

1590: .keywords: mesh, elements
1591: .seealso: DMMeshCreate()
1592: @*/
1593: PetscErrorCode SectionIntRestrictClosure(SectionInt section, DM dm, PetscInt point, PetscInt n, PetscInt values[])
1594: {
1595:   ALE::Obj<PETSC_MESH_TYPE> m;
1596:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> s;
1597:   PetscErrorCode            ierr;

1600:   DMMeshGetMesh(dm, m);
1601:   SectionIntGetSection(section, s);
1602:   m->restrictClosure(s, point, values, n);
1603:   return(0);
1604: }

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

1611:   Not Collective

1613:   Input Parameters:
1614: + section - The section
1615: . mesh    - The DMMesh object
1616: . point   - The sieve point
1617: . array   - The array to fill up
1618: - mode    - The insertion mode

1620:   Output Parameter:
1621: . array - The array full of values in the closure

1623:   Level: intermediate

1625: .keywords: mesh, elements
1626: .seealso: DMMeshCreate()
1627: @*/
1628: PetscErrorCode SectionIntUpdateClosure(SectionInt section, DM dm, PetscInt point, PetscInt values[], InsertMode mode)
1629: {
1630:   ALE::Obj<PETSC_MESH_TYPE> m;
1631:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> s;
1632:   PetscErrorCode            ierr;

1635:   DMMeshGetMesh(dm, m);
1636:   SectionIntGetSection(section, s);
1637:   if (mode == INSERT_VALUES) {
1638:     m->update(s, point, values);
1639:   } else if (mode == ADD_VALUES) {
1640:     m->updateAdd(s, point, values);
1641:   } else {
1642:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Invalid insertion mode: %d", mode);
1643:   }
1644:   return(0);
1645: }

1649: /*@
1650:   SectionIntComplete - Exchanges data across the mesh overlap.

1652:   Not collective

1654:   Input Parameter:
1655: . section - the section object

1657:   Level: advanced

1659: .seealso SectionIntRestrict(), SectionIntCreate(), SectionIntView()
1660: @*/
1661: PetscErrorCode SectionIntComplete(SectionInt section)
1662: {
1663:   Obj<PETSC_MESH_TYPE::int_section_type> s;
1664:   Obj<PETSC_MESH_TYPE>                   b;

1668:   SectionIntGetSection(section, s);
1669:   SectionIntGetBundle(section, b);
1670: #if 0
1671:   ALE::Distribution<PETSC_MESH_TYPE>::completeSection(b, s);
1672: #else
1673:   ALE::Completion::completeSectionAdd(b->getSendOverlap(), b->getRecvOverlap(), s, s);
1674: #endif
1675:   return(0);
1676: }

1680: /*@
1681:   SectionIntZero - Zero out the entries

1683:   Not collective

1685:   Input Parameter:
1686: . section - the section object

1688:   Level: advanced

1690: .seealso SectionIntRestrict(), SectionIntCreate(), SectionIntView()
1691: @*/
1692: PetscErrorCode SectionIntZero(SectionInt section)
1693: {
1694:   Obj<PETSC_MESH_TYPE::int_section_type> s;

1698:   SectionIntGetSection(section, s);
1699:   s->zero();
1700:   return(0);
1701: }

1705: /*@
1706:   SectionIntGetFiberDimension - Get the size of the vector space attached to the point

1708:   Not collective

1710:   Input Parameters:
1711: + section - the section object
1712: - point - the Sieve point

1714:   Output Parameters:
1715: . size - The fiber dimension

1717:   Level: advanced

1719: .seealso SectionRealRestrict(), SectionRealCreate(), SectionRealView()
1720: @*/
1721: PetscErrorCode  SectionIntGetFiberDimension(SectionInt section, PetscInt point, PetscInt *size)
1722: {
1725:   *size = section->s->getFiberDimension(point);
1726:   return(0);
1727: }

1731: /*@
1732:   SectionIntSetFiberDimension - Set the size of the vector space attached to the point

1734:   Not collective

1736:   Input Parameters:
1737: + section - the section object
1738: . point - the Sieve point
1739: - size - The fiber dimension

1741:   Level: advanced

1743: .seealso SectionIntRestrict(), SectionIntCreate(), SectionIntView()
1744: @*/
1745: PetscErrorCode  SectionIntSetFiberDimension(SectionInt section, PetscInt point, const PetscInt size)
1746: {
1749:   section->s->setFiberDimension(point, size);
1750:   return(0);
1751: }

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

1758:   Not collective

1760:   Input Parameters:
1761: + section - the section object
1762: . point - the Sieve point
1763: . size - The fiber dimension
1764: - field - The field number

1766:   Level: advanced

1768: .seealso SectionIntSetFiberDimension(), SectionIntRestrict(), SectionIntCreate(), SectionIntView()
1769: @*/
1770: PetscErrorCode  SectionIntSetFiberDimensionField(SectionInt section, PetscInt point, const PetscInt size, const PetscInt field)
1771: {
1774:   section->s->setFiberDimension(point, size, field);
1775:   return(0);
1776: }

1780: /*@
1781:   SectionIntGetSize - Gets the number of local dofs in this Section

1783:   Not collective

1785:   Input Parameter:
1786: . section - the section object

1788:   Output Parameter:
1789: . size - the section size

1791:   Level: advanced

1793: .seealso SectionIntRestrict(), SectionIntCreate(), SectionIntView()
1794: @*/
1795: PetscErrorCode SectionIntGetSize(SectionInt section, PetscInt *size)
1796: {
1797:   Obj<PETSC_MESH_TYPE::int_section_type> s;

1802:   SectionIntGetSection(section, s);
1803:   *size = s->size();
1804:   return(0);
1805: }

1809: /*@
1810:   SectionIntAllocate - Allocate storage for this section

1812:   Not collective

1814:   Input Parameter:
1815: . section - the section object

1817:   Level: advanced

1819: .seealso SectionIntRestrict(), SectionIntCreate(), SectionIntView()
1820: @*/
1821: PetscErrorCode  SectionIntAllocate(SectionInt section)
1822: {
1825:   section->b->allocate(section->s);
1826:   return(0);
1827: }

1831: /*@C
1832:   SectionIntClear - Dellocate storage for this section

1834:   Not collective

1836:   Input Parameter:
1837: . section - the section object

1839:   Level: advanced

1841: .seealso SectionIntRestrict(), SectionIntCreate(), SectionIntView()
1842: @*/
1843: PetscErrorCode  SectionIntClear(SectionInt section)
1844: {
1847:   section->s->clear();
1848:   return(0);
1849: }

1853: /*@
1854:   SectionIntSet - Sets all the values to the given value

1856:   Not collective

1858:   Input Parameters:
1859: + section - the real Section
1860: - val - the value

1862:   Level: intermediate

1864: .seealso VecNorm(), SectionIntCreate()
1865: @*/
1866: PetscErrorCode  SectionIntSet(SectionInt section, PetscInt val)
1867: {
1868:   Obj<PETSC_MESH_TYPE::int_section_type> s;

1872:   SectionIntGetSection(section, s);
1873:   s->set(val);
1874:   return(0);
1875: }

1879: /*@
1880:   SectionIntAddSpace - Add another field to this section

1882:   Collective on DMMesh

1884:   Input Parameter:
1885: . section - the Section

1887:   Level: advanced

1889: .seealso SectionIntCreate(), SectionIntGetFibration()
1890: @*/
1891: PetscErrorCode  SectionIntAddSpace(SectionInt section)
1892: {
1893:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> s;

1897:   SectionIntGetSection(section, s);
1898:   s->addSpace();
1899:   return(0);
1900: }

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

1907:   Collective on DMMesh

1909:   Input Parameter:
1910: + section - the Section
1911: - field- The field number

1913:   Output Parameter:
1914: . subsection - the section of the given field

1916:   Level: advanced

1918: .seealso SectionIntCreate(), SectionIntAddSpace()
1919: @*/
1920: PetscErrorCode  SectionIntGetFibration(SectionInt section, const PetscInt field, SectionInt *subsection)
1921: {
1922:   ALE::Obj<PETSC_MESH_TYPE>                   b;
1923:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> s;
1924:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> t;
1925:   MPI_Comm       comm;

1929:   PetscObjectGetComm((PetscObject) section, &comm);
1930:   SectionIntGetBundle(section, b);
1931:   SectionIntGetSection(section, s);
1932:   SectionIntCreate(comm, subsection);
1933:   SectionIntSetBundle(*subsection, b);
1934:   t    = s->getFibration(field);
1935:   SectionIntSetSection(*subsection, t);
1936:   return(0);
1937: }

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

1944:   Collective on DMMesh

1946:   Input Parameters:
1947: + mesh - The DMMesh object
1948: - fiberDim - The section name

1950:   Output Parameter:
1951: . section - The section

1953:   Level: intermediate

1955: .keywords: mesh, section, vertex
1956: .seealso: DMMeshCreate(), SectionIntCreate()
1957: @*/
1958: PetscErrorCode DMMeshGetVertexSectionInt(DM dm, const char name[], PetscInt fiberDim, SectionInt *section)
1959: {
1960:   ALE::Obj<PETSC_MESH_TYPE> m;
1961:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> s;
1962:   PetscErrorCode      ierr;

1965:   DMMeshGetMesh(dm, m);
1966:   SectionIntCreate(m->comm(), section);
1967:   PetscObjectSetName((PetscObject) *section, name);
1968:   SectionIntSetBundle(*section, m);
1969:   SectionIntGetSection(*section, s);
1970: #ifdef PETSC_OPT_SIEVE
1971:   s->setChart(m->getSieve()->getChart());
1972: #endif
1973:   s->setName(name);
1974:   s->setFiberDimension(m->depthStratum(0), fiberDim);
1975:   m->allocate(s);
1976:   return(0);
1977: }

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

1984:   Collective on DMMesh

1986:   Input Parameters:
1987: + mesh - The DMMesh object
1988: - fiberDim - The section name

1990:   Output Parameter:
1991: . section - The section

1993:   Level: intermediate

1995: .keywords: mesh, section, cell
1996: .seealso: DMMeshCreate(), SectionIntCreate()
1997: @*/
1998: PetscErrorCode DMMeshGetCellSectionInt(DM dm, const char name[], PetscInt fiberDim, SectionInt *section)
1999: {
2000:   ALE::Obj<PETSC_MESH_TYPE> m;
2001:   ALE::Obj<PETSC_MESH_TYPE::int_section_type> s;
2002:   PetscErrorCode      ierr;

2005:   DMMeshGetMesh(dm, m);
2006:   SectionIntCreate(m->comm(), section);
2007:   PetscObjectSetName((PetscObject) *section, name);
2008:   SectionIntSetBundle(*section, m);
2009:   SectionIntGetSection(*section, s);
2010: #ifdef PETSC_OPT_SIEVE
2011:   s->setChart(m->getSieve()->getChart());
2012: #endif
2013:   s->setName(name);
2014:   s->setFiberDimension(m->heightStratum(0), fiberDim);
2015:   m->allocate(s);
2016:   return(0);
2017: }

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

2024:   Collective on DMMesh

2026:   Input Parameters:
2027: . mesh - The DMMesh object

2029:   Output Parameter:
2030: . section - The section

2032:   Level: intermediate

2034: .keywords: mesh, section, cell
2035: .seealso: DMMeshCreate(), SectionRealCreate()
2036: @*/
2037: PetscErrorCode DMMeshSetupSection(DM dm, SectionReal section)
2038: {
2039:   ALE::Obj<PETSC_MESH_TYPE> m;
2040:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
2041:   PetscErrorCode      ierr;

2044:   DMMeshGetMesh(dm, m);
2045:   SectionRealGetSection(section, s);
2046:   try {
2047:     m->setupField(s);
2048:   } catch(ALE::Exception e) {
2049:     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "%s", e.message());
2050:   }
2051:   return(0);
2052: }