<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Tue, Jan 28, 2014 at 7:03 PM, Geoffrey Irving <span dir="ltr"><<a href="mailto:irving@naml.us" target="_blank">irving@naml.us</a>></span> wrote:<br>
<blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div class="im">On Tue, Jan 28, 2014 at 4:18 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com">knepley@gmail.com</a>> wrote:<br>

> On Tue, Jan 28, 2014 at 6:04 PM, Geoffrey Irving <<a href="mailto:irving@naml.us">irving@naml.us</a>> wrote:<br>
>><br>
>> What's the cleanest way to save data from a petsc finite element<br>
>> simulation including boundary conditions?  First, is it correct that<br>
>> the local sections contains dofs for boundary conditions, while the<br>
>> global section does not?  If so, is the following the right plan?<br>
>><br>
>> 1. At start, save the DMPlex<br>
>> 2. At start, save the local section<br>
>> 3. Every frame, sync the global vector a local vector and write out in<br>
>> petsc binary format.<br>
>><br>
>> Are there standard binary formats for DMPlex and sections?  All I see<br>
>> is the following sad commented out call to DMPlexView_Binary.<br>
><br>
> We are having a big discussion about this on petsc-dev. I will write a longer<br>
> thing, but now the only thing that works all the way is VTK. You should just<br>
> be able to create the Viewer, do VecView for your fields, and then destroy<br>
> it.<br>
<br>
</div>Thanks, and apologies for missing that thread.<br>
<br>
Is there a working example of that kind of viewer use?<br>
DMPlexVTKWriteAll_ASCII seems to traverse a linked list of vectors<br>
attached to the viewer and write them all out, which is a bit<br>
confusing.</blockquote><div><br></div><div>In ex62:</div><div><br></div><div>  if (user.runType == RUN_FULL) {</div><div>    PetscViewer viewer;</div><div>    Vec         uLocal;</div><div>    const char *name;</div><div>
<br></div><div>    ierr = PetscViewerCreate(PETSC_COMM_WORLD, &viewer);CHKERRQ(ierr);</div><div>    ierr = PetscViewerSetType(viewer, PETSCVIEWERVTK);CHKERRQ(ierr);</div><div>    ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr);</div>
<div>    ierr = PetscViewerFileSetName(viewer, "ex62_sol.vtk");CHKERRQ(ierr);</div><div><br></div><div>    ierr = DMGetLocalVector(dm, &uLocal);CHKERRQ(ierr);</div><div>    ierr = PetscObjectGetName((PetscObject) u, &name);CHKERRQ(ierr);</div>
<div>    ierr = PetscObjectSetName((PetscObject) uLocal, name);CHKERRQ(ierr);</div><div>    ierr = DMGlobalToLocalBegin(dm, u, INSERT_VALUES, uLocal);CHKERRQ(ierr);</div><div>    ierr = DMGlobalToLocalEnd(dm, u, INSERT_VALUES, uLocal);CHKERRQ(ierr);</div>
<div>    ierr = VecView(uLocal, viewer);CHKERRQ(ierr);</div><div>    ierr = DMRestoreLocalVector(dm, &uLocal);CHKERRQ(ierr);</div><div><br></div><div>    ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);</div><div>
  }</div><div>  </div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div class="im">
> I create a global section with all the boundary values in it.<br>
<br>
</div>Do you mean this is what I'll need to do, or is that taken care of<br>
automatically in the VTK path?  It seems like boundary condition<br>
handling can't be automatic, but DMPlexVTKWriteAll_ASCII reads a bunch<br>
of section information out of the DM so I'm sure how one would<br>
substitute a global section + boundaries.</blockquote><div><br></div><div>You are right, this is not what I do anymore. I decided it was easier for all the Viewers</div><div>to just take local vectors. I handle all the "globalization" automatically. See the code snippet above.</div>
<div><br></div><div>  Thanks,</div><div><br></div><div>    Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">
<div class="im">
> I plan on having HDF5+Xdmf and Exodus.<br>
<br>
</div>Thanks,<br>
Geoffrey<br>
</blockquote></div><br><br clear="all"><div><br></div>-- <br>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>
-- Norbert Wiener
</div></div>