<html>
<head>
<meta http-equiv="Content-Type" content="text/html; charset=UTF-8">
</head>
<body text="#000000" bgcolor="#FFFFFF">
<p><br>
</p>
<div class="moz-cite-prefix">On 2018-12-03 12:03 p.m., Matthew
Knepley wrote:<br>
</div>
<blockquote type="cite"
cite="mid:CAMYG4GmRXGMm7fP+aCfEpQiD_tX-fHpy2RLBeukFBk3O8zruGQ@mail.gmail.com">
<meta http-equiv="content-type" content="text/html; charset=UTF-8">
<div dir="ltr">
<div class="gmail_quote">
<div dir="ltr">On Mon, Dec 3, 2018 at 2:27 PM Danyang Su <<a
href="mailto:danyang.su@gmail.com" moz-do-not-send="true">danyang.su@gmail.com</a>>
wrote:<br>
</div>
<blockquote class="gmail_quote" style="margin:0 0 0
.8ex;border-left:1px #ccc solid;padding-left:1ex">
<div text="#000000" bgcolor="#FFFFFF">
<p>Hi Matt,</p>
<p>Thanks. <br>
</p>
<p>BTW: DmPlexGetVertexNumbering now can work using the
latest develop version. But the index is not in natural
ordering when DMSetUseNatural is called. That's why I
want to use PetscSFDistributeSection to check if I miss
anything in the code.</p>
</div>
</blockquote>
<div>Can you explain that a little more? Maybe you can just
push forward what you want using the migrationSF.</div>
</div>
</div>
</blockquote>
<p>Hi Matt,</p>
<p>Since I cannot figure what is wrong or missing in my code, I
followed an old ex26.c example in src/dm/impls/plex/examples/tests
to create similar code as shown below to test global to natural
ordering. The code may be ugly with unnecessary functions in it.
Using DmPlexGetVertexNumbering, I can get the value but it is not
in natural order, instead, it is still in default PETSc order
without calling DMSetUseNatural(dm,PETSC_TRUE,ierr).</p>
if (rank == 0) then<br>
call
DMPlexCreateFromCellList(Petsc_Comm_World,ndim,num_cells,
num_nodes,num_nodes_per_cell, &<br>
Petsc_False,dmplex_cells,ndim, dmplex_verts,dm,ierr)<br>
CHKERRQ(ierr)<br>
else<br>
call DMPlexCreateFromCellList(Petsc_Comm_World,ndim,0,
0,num_nodes_per_cell, &<br>
Petsc_False,dmplex_cells,ndim,dmplex_verts,dm,ierr)<br>
CHKERRQ(ierr)<br>
end if<br>
<br>
if (nprocs > 1) then<br>
call DMSetUseNatural(dm,PETSC_TRUE,ierr)<br>
CHKERRQ(ierr)<br>
end if<br>
<br>
call DMPlexDistribute(dm,stencil_width,
&<br>
migrationsf,distributedMesh,ierr)<br>
CHKERRQ(ierr)<br>
<br>
if (distributedMesh /= PETSC_NULL_DM) then<br>
call
PetscSFCreateInverseSF(migrationsf,migrationsf_inv,ierr)<br>
CHKERRQ(ierr)<br>
<br>
call
DMCreateGlobalToNatural(distributedMesh,migrationsf,migrationsf_inv,ierr)<br>
CHKERRQ(ierr)<br>
<br>
call DMGetSection(distributedMesh,section,ierr)<br>
CHKERRQ(ierr)<br>
<br>
call
PetscSectionCreate(Petsc_Comm_World,section_seq,ierr)<br>
CHKERRQ(ierr)<br>
<br>
call
PetscSFDistributeSection(migrationsf_inv,section, &<br>
PETSC_NULL_INTEGER,section_seq,ierr)<br>
CHKERRQ(ierr)<br>
<br>
call
DMPlexCreateGlobalToNaturalSF(distributedMesh, &<br>
section_seq,migrationsf,sf_natural,ierr)<br>
CHKERRQ(ierr)<br>
<br>
call DMSetUseNatural(distributedMesh,PETSC_TRUE,ierr)<br>
CHKERRQ(ierr)<br>
<br>
call PetscSFDestroy(migrationsf,ierr)<br>
CHKERRQ(ierr)<br>
<br>
call PetscSFDestroy(migrationsf_inv,ierr)<br>
CHKERRQ(ierr)<br>
<p> end if</p>
<p>Thanks,</p>
<p>Danyang<br>
</p>
<blockquote type="cite"
cite="mid:CAMYG4GmRXGMm7fP+aCfEpQiD_tX-fHpy2RLBeukFBk3O8zruGQ@mail.gmail.com">
<div dir="ltr">
<div class="gmail_quote">
<div><br>
</div>
<div> Thanks,</div>
<div><br>
</div>
<div> Matt <br>
</div>
<blockquote class="gmail_quote" style="margin:0 0 0
.8ex;border-left:1px #ccc solid;padding-left:1ex">
<div text="#000000" bgcolor="#FFFFFF">
<p>Regards,</p>
<p>Danyang<br>
</p>
<div class="m_-4663990655342689142moz-cite-prefix">On
2018-12-03 5:22 a.m., Matthew Knepley wrote:<br>
</div>
<blockquote type="cite">
<div dir="ltr">I need to write a custom Fortran stub for
this one. I will get it done as soon as possible.
<div><br>
</div>
<div> Thanks,</div>
<div><br>
</div>
<div> Matt</div>
</div>
<br>
<div class="gmail_quote">
<div dir="ltr">On Sat, Dec 1, 2018 at 7:16 PM Danyang
Su via petsc-users <<a
href="mailto:petsc-users@mcs.anl.gov"
target="_blank" moz-do-not-send="true">petsc-users@mcs.anl.gov</a>>
wrote:<br>
</div>
<blockquote class="gmail_quote" style="margin:0 0 0
.8ex;border-left:1px #ccc solid;padding-left:1ex">Hi
All,<br>
<br>
I got a simple compilation error when use
PetscSFDistributeSection in <br>
Fortran. It looks like the required head files are
included and the <br>
parameters are correctly defined. However, when
compile the code, I got <br>
error undefined reference to
`petscsfdistributesection_'. The code is <br>
shown below. Did I miss anything here?<br>
<br>
#include <petsc/finclude/petscsys.h><br>
#include <petsc/finclude/petscvec.h><br>
#include <petsc/finclude/petscdm.h><br>
#include <petsc/finclude/petscdmplex.h><br>
use petscsys<br>
use petscvec<br>
use petscdm<br>
use petscdmplex<br>
<br>
implicit none<br>
<br>
PetscSection :: section, section_seq<br>
PetscSF :: migrationsf_inv, sf_natural<br>
Vec :: vec_global, vec_natural<br>
PetscErrorCode :: ierr<br>
<br>
...<br>
<br>
call
PetscSFDistributeSection(migrationsf_inv,section,
&<br>
PETSC_NULL_INTEGER,section_seq,ierr)<br>
CHKERRQ(ierr)<br>
<br>
<br>
call
PetscSFDistributeSection(migrationsf_inv,section,
&<br>
PETSC_NULL_INTEGER,section_seq,ierr)<br>
CHKERRQ(ierr)<br>
<br>
Thanks,<br>
<br>
Danyang<br>
<br>
</blockquote>
</div>
<br clear="all">
<div><br>
</div>
-- <br>
<div dir="ltr"
class="m_-4663990655342689142gmail_signature"
data-smartmail="gmail_signature">
<div dir="ltr">
<div>
<div dir="ltr">
<div>
<div dir="ltr">
<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><br>
</div>
<div><a
href="http://www.cse.buffalo.edu/~knepley/"
target="_blank" moz-do-not-send="true">https://www.cse.buffalo.edu/~knepley/</a><br>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</blockquote>
</div>
</blockquote>
</div>
<br clear="all">
<div><br>
</div>
-- <br>
<div dir="ltr" class="gmail_signature"
data-smartmail="gmail_signature">
<div dir="ltr">
<div>
<div dir="ltr">
<div>
<div dir="ltr">
<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><br>
</div>
<div><a href="http://www.cse.buffalo.edu/~knepley/"
target="_blank" moz-do-not-send="true">https://www.cse.buffalo.edu/~knepley/</a><br>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</blockquote>
</body>
</html>