<div dir="ltr">Hi, when looping over the nodal coordinates, is there a way to get the corresponding DOF in the index set for the DM Field(s) that is used in the residual vector? <br><div>I would like to update the jacobian/residual entries obtained from <span style="color:rgb(5,80,174);font-family:ui-monospace,SFMono-Regular,"SF Mono",Menlo,Consolas,"Liberation Mono",monospace;font-size:12px;white-space:pre">DMCreateFieldIS</span> that correspond to particular nodal coordinates that are not on the boundary.</div><div><br></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Fri, Jul 10, 2015 at 11:14 AM Matthew Knepley <<a href="mailto:knepley@gmail.com">knepley@gmail.com</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Fri, Jul 10, 2015 at 6:48 AM, Alejandro D Otero <span dir="ltr"><<a href="mailto:aotero@fi.uba.ar" target="_blank">aotero@fi.uba.ar</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">Hi Matt, thanks for your answer.<br><div><br>I got a vector from <span>getCoordinates(). How are their components indexed? is (p * dim + d) with p: node, dim: of the problem, x_d coordinate of the node, ok? <br>Which numbering for p? The local number of node, the number of point in the DAG of the dm, the original number of node? <br></span></div></div></blockquote><div><br></div><div>All data layouts for Plex are described by a PetscSection (which is how it should be in the rest of PETSc as well). The PetscSection is a map</div><div>from mesh points to (# dof, offset) pairs. Thus, the way I get coordinates is the following:</div><div><br></div><div>DMGetCoordinatesLocal(dm, &coordinates);</div><div>DMGetCoordinateDM(dm, &cdm);</div><div>DMGetDefaultSection(cdm, &cs);</div><div>PetscSectionGetChart(cs, &pStart, &pEnd);</div><div>VecGetArrayRead(coordinates, &coords);</div><div>for (p = pStart; p < pEnd; ++p) {</div><div> PetscInt dof, off;</div><div><br></div><div> PetscSectionGetDof(cs, p, &dof);</div><div> PetscSectionGetOffset(cs, p, &off);</div><div> for (d = 0; d < dof; ++d) <compute with coords[off+d]></div><div>}</div><div>VecRestoreArrayRead(coordinates, &coords);</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div><span></span></div><div><span>I am trying a simple square mesh with 16 4-node square elements parted into 2 process. Total of 25 nodes.<br>The distributed dm seems alright to me. Each process gets a dm with 8 elements an 15 nodes, which means that the 5 shared nodes are local to each process. But one of the process gives negative values for the shared nodes. How need them to be mapped to get the right number.<br></span></div></div></blockquote><div><br></div><div>Where do you get negative point numbers? I encode off-process point numbers as -(remote point + 1).</div><div> <br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div><span></span></div><div><span>It seems I'm using a wrong approach to this. May be I need to get the coordinates in a somehow different way. I'm working on a from-scratch implementation of a FEM code based on petsc4py. I want to code the problem matrices assembly from elemental matrices. I've already done this sequentially, but I got stuck when trying to compute elemental matrices in parallel because I don't understand the way of obtaining the coordinates of the nodes in for each element.<br></span></div></div></blockquote><div><br></div><div>The way I get coordinates for the element c is</div><div><br></div><div> PetscScalar *coords = NULL;</div><div> PetscInt csize;</div><div><br></div><div><div> ierr = DMPlexVecGetClosure(dm, cs, coordinates, c, &csize, &coords);CHKERRQ(ierr);</div></div><div><div> ierr = DMPlexVecRestoreClosure(dm, cs, coordinates, c, &csize, &coords);CHKERRQ(ierr);</div></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:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div><span></span></div><div><span>Again, thanks for your help,<br><br></span></div><div><span>A.<br></span></div><div><span><br><br></span></div><div><span><br><br><br><br><br><br></span></div></div><div class="gmail_extra"><br><div class="gmail_quote">On Thu, Jul 9, 2015 at 5:18 PM, Matthew Knepley <span dir="ltr"><<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div class="gmail_extra"><div class="gmail_quote"><span>On Thu, Jul 9, 2015 at 7:42 AM, Alejandro D Otero <span dir="ltr"><<a href="mailto:aotero@fi.uba.ar" target="_blank">aotero@fi.uba.ar</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div><div><div><div><div><div>Hi, sorry if this is an obvious question, but I cannot figure out how to recover finite element nodes coordinates once I have distributed a mesh stored as a dmplex. I am using petsc4py as interface to petsc rutines.<br><br></div><div>I first created a dmplex using:<br>dm.createFromCellList()<br><br></div><div>In a sequential run I got the coordinates with:<br>Coords = dm.getCoordinates()<br></div><div><br></div>which gave a sequential vector with the coordinates of the mesh nodes.<br><br></div>When I distribute the mesh with:<br></div>dm.distribute()<br><br></div>each mpi process has it own dm but the indexing of the vector resulting from getCoordinates() or getCoordinatesLocal() seems not consistent with the local numbering of the cells and nodes.<br></div></div></div></blockquote><div><br></div></span><div>When the mesh is distributed, the vertices are renumbered. Thus the coordinates you get out are</div><div>for reordered local vertices, but they are consistent with the local topology (cells still contain the</div><div>right vertices) and the overlap mapping (SF still connects the shared vertices).</div><div><br></div><div>What do you need it to do?</div><div><br></div><div> Thanks,</div><div><br></div><div> Matt</div><span><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div><div></div>Which is the correct way of doing this in PETSc philosophy? <br><br></div>Thanks in advance, <br><div><div><div>Alejandro<br></div></div></div></div>
</blockquote></span></div><span><font color="#888888"><br><br clear="all"><span><font color="#888888"><div><br></div>-- <br><div>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>
</font></span></font></span></div></div>
</blockquote></div><br></div>
</blockquote></div><br><br clear="all"><div><br></div>-- <br><div>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></div>
</blockquote></div>