<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Sat, Dec 7, 2013 at 1:26 PM, Mark Lohry <span dir="ltr"><<a href="mailto:mlohry@gmail.com" target="_blank">mlohry@gmail.com</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
<div dir="ltr"><div>1) Is there somewhere on the doc site with a single list of all examples and what they do? As a new user I hugely appreciate all the great examples and documentation, but it's a little unwieldy to try and find relevant examples by navigating the function documents.</div>
</div></blockquote><div><br></div><div>No, we have not yet setup a nice system for this.</div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr">
2) I'm a little unclear as to the best-practice way to handle solution/grid related data that isn't actual solution DOFS or mesh points, and doesn't require global communication.<br><br>Say I want to work on an (ni X nj ) cell-centered finite volume layout with one DOF per cell. I start with<div>
<br></div><div><div>DMDACreate2d(PETSC_COMM_WORLD,</div><div><span style="white-space:pre-wrap"> </span> DMDA_BOUNDARY_PERIODIC, DMDA_BOUNDARY_NONE, /* periodic in i, no conn in j */</div><div><span style="white-space:pre-wrap"> </span> DMDA_STENCIL_BOX, /* 5-point stencil */</div>
<div><span style="white-space:pre-wrap"> </span> ni,nj, /* number of cells */</div><div><span style="white-space:pre-wrap"> </span> PETSC_DECIDE,PETSC_DECIDE, /* let petsc determine number of procs */</div>
<div><span style="white-space:pre-wrap"> </span> 1, /* one degree of freedom per cell */</div><div><span style="white-space:pre-wrap"> </span> 1, /* stencil width */</div><div><span style="white-space:pre-wrap"> </span> NULL,NULL,&da);</div>
</div><div><br></div><div>Now I can construct the solution vector like so:</div><div><br></div><div>struct Field { PetscScalar u; };</div><div>Field **soln;</div><div>Vec global_vec;</div><div>DMDACreateGlobalVector(da,&global_vec);</div>
<div>DMDAVecGetArray(da,global_vec,&soln); // converts global vector into Field** </div><div><br></div><div>and access the local elements through soln[j][i]. But say I want to store some grid information locally, say the FaceNormal for a cell face:</div>
<div>struct FaceNormal{ PetscScalar nx,ny; };</div><div>FaceNormal **IFace;</div><div><br></div><div>Now it seems that the CreateGlobalVector/CreateLocalVector are just going to allocate for 1 DOF because of how the DMDA was defined, so there wouldn't be sufficient storage for the 2 variables in FaceNormal. What's an appropriate way to do this? Clearly I could just malloc that myself, but it'd be nice to keep the same indexing as the petsc vectors.</div>
</div></blockquote><div><br></div><div>There are at least two answers:</div><div><br></div><div> 1) If you do not care about redundant storage, you could create another 2d DMDA using dof = 8, and store all the normals for</div>
<div> each cell. DMDAs take almost no memory. Your storage would be the Vec you use to store the normals.</div><div><br></div><div> 2) If you care about storage, then create another DMDA, but then create a PetscSection to assign 2 dof to every face, which you can do with</div>
<div><br></div><div> DMDACreateSection()</div><div><br></div><div> and use DMSetDefaultSection(). This will make DMCreateGlobalVector() create a Vec you can use to store all the face normals.</div><div><br></div>
<div> Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div>Thanks in advance,</div></div>
</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>