On Tue, Aug 28, 2012 at 4:39 PM, Chris Eldred <span dir="ltr"><<a href="mailto:chris.eldred@gmail.com" target="_blank">chris.eldred@gmail.com</a>></span> wrote:<br><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
DMComplexGetTransitiveClosure does not appear to have the calling<br>
sequence (DM dm, PetscInt p, PetscBool useCone, PetscInt *numPoints,<br>
PetscInt *points[]) indicated in the documentation- instead it has DM<br>
dm, PetscInt p, PetscBool useCone, PetscInt *points[]. In addition, it<br>
does not appear to be returning the full closure- instead it just<br>
gives cone + the point.<br></blockquote><div><br></div><div>1) I will turn this code into a running example</div><div><br></div><div>2) I will fix the documentation. The interface change for Fortran comes because the size</div>
<div>    is embedded in the array pointer. If you think I should still return the size in Fortran,</div><div>    let me know.</div><div><br></div><div>3) I will also check to make sure the documentation is clear that this thing returns each</div>
<div>    point AND its orientation.</div><div><br></div><div>Will write again when I have your thing running. It could be an error with the Fortran array</div><div>processing I put in.</div><div><br></div><div>  Thanks,</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">
Sample code is below:<br>
<br>
program main<br>
implicit none<br>
!<br>
! #echo '#include <finclude/petscsys.h>\n'<br>
#echo '#include <finclude/petsc.h90>'<br>
#echo '#include <finclude/petscdmcomplex.h>\n'<br>
#echo '#include <finclude/petscdmcomplex.h90>'<br>
! #echo '#include <finclude/petscdm.h>\n'<br>
! #echo '#include <finclude/petscdm.h90>\n'<br>
<br>
!<br>
DM dm<br>
PetscInt, target, dimension(3) :: EC<br>
PetscInt, target, dimension(2) :: VE<br>
PetscInt, pointer :: pEC(:), pVE(:)<br>
PetscInt, pointer :: nClosure(:)<br>
PetscInt, pointer :: pClosure(:), pCone(:), pCone2(:)<br>
PetscErrorCode ierr<br>
<br>
call PetscInitialize(PETSC_NULL_CHARACTER,ierr)<br>
CHKERRQ(ierr)<br>
call DMComplexCreate(PETSC_COMM_WORLD, dm, ierr)<br>
CHKERRQ(ierr)<br>
<br>
!Make Doublet Mesh from Fig 2 of Flexible Representation of<br>
Computational Meshes (except indexing is from 0 instead of 1)<br>
call DMComplexSetChart(dm, 0, 11, ierr)<br>
CHKERRQ(ierr)<br>
<br>
call DMComplexSetConeSize(dm, 9, 3,ierr)<br>
CHKERRQ(ierr)<br>
call DMComplexSetConeSize(dm, 10, 3,ierr)<br>
CHKERRQ(ierr)<br>
call DMComplexSetConeSize(dm, 4, 2,ierr)<br>
CHKERRQ(ierr)<br>
call DMComplexSetConeSize(dm, 5, 2,ierr)<br>
CHKERRQ(ierr)<br>
call DMComplexSetConeSize(dm, 6, 2,ierr)<br>
CHKERRQ(ierr)<br>
call DMComplexSetConeSize(dm, 7, 2,ierr)<br>
CHKERRQ(ierr)<br>
call DMComplexSetConeSize(dm, 8, 2,ierr)<br>
CHKERRQ(ierr)<br>
<br>
call DMSetUp(dm, ierr)<br>
CHKERRQ(ierr)<br>
<br>
EC(1) = 4<br>
EC(2) = 5<br>
EC(3) = 6<br>
pEC => EC<br>
call DMComplexSetCone(dm, 9 , pEC, ierr)<br>
CHKERRQ(ierr)<br>
EC(1) = 6<br>
EC(2) = 7<br>
EC(3) = 8<br>
pEC => EC<br>
call DMComplexSetCone(dm, 10 , pEC, ierr)<br>
CHKERRQ(ierr)<br>
<br>
VE(1) = 0<br>
VE(2) = 1<br>
pVE => VE<br>
call DMComplexSetCone(dm, 4 , pVE, ierr)<br>
CHKERRQ(ierr)<br>
VE(1) = 0<br>
VE(2) = 2<br>
pVE => VE<br>
call DMComplexSetCone(dm, 5 , pVE, ierr)<br>
CHKERRQ(ierr)<br>
VE(1) = 1<br>
VE(2) = 2<br>
pVE => VE<br>
call DMComplexSetCone(dm, 6 , pVE, ierr)<br>
CHKERRQ(ierr)<br>
VE(1) = 2<br>
VE(2) = 3<br>
pVE => VE<br>
call DMComplexSetCone(dm, 7 , pVE, ierr)<br>
CHKERRQ(ierr)<br>
VE(1) = 1<br>
VE(2) = 3<br>
pVE => VE<br>
call DMComplexSetCone(dm, 8 , pVE, ierr)<br>
CHKERRQ(ierr)<br>
<br>
call DMComplexSymmetrize(dm,ierr)<br>
CHKERRQ(ierr)<br>
call DMComplexStratify(dm,ierr)<br>
CHKERRQ(ierr)<br>
<br>
!Test Stuff<br>
<br>
call DMComplexGetTransitiveClosure(dm,10,PETSC_TRUE,nClosure,ierr)<br>
!pClosure,ierr)<br>
CHKERRQ(ierr)<br>
write(*,*) nClosure<br>
write(*,*) pClosure<br>
<br>
call DMDestroy(dm,ierr)<br>
CHKERRQ(ierr)<br>
call PetscFinalize(ierr)<br>
CHKERRQ(ierr)<br>
end<br>
<br>
<br>
The output is:<br>
          10           0           6           0           7<br>
0           8<br>
<br>
and then there is an error in DMDestroy<br>
<br>
[0]PETSC ERROR: --------------------- Error Message<br>
------------------------------------<br>
[0]PETSC ERROR: Object is in wrong state!<br>
[0]PETSC ERROR: Work array still checked out!<br>
[0]PETSC ERROR:<br>
------------------------------------------------------------------------<br>
[0]PETSC ERROR: Petsc Development HG revision:<br>
b6fe82991deee4a0d3f9a20654bc7750a6b1fe0f  HG Date: Mon Aug 27 13:17:10<br>
2012 -0500<br>
[0]PETSC ERROR: See docs/changes/index.html for recent updates.<br>
[0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.<br>
[0]PETSC ERROR: See docs/index.html for manual pages.<br>
[0]PETSC ERROR:<br>
------------------------------------------------------------------------<br>
[0]PETSC ERROR: ./test on a arch-linu named Puget-101334 by user Tue<br>
Aug 28 15:31:06 2012<br>
[0]PETSC ERROR: Libraries linked from<br>
/home/user/Desktop/LIBRARIES/petsc-dev/arch-linux2-cxx-debug/lib<br>
[0]PETSC ERROR: Configure run at Mon Aug 27 12:55:14 2012<br>
[0]PETSC ERROR: Configure options --download-boost --download-chaco<br>
--download-ctetgen --download-f-blas-lapack --download-fiat<br>
--download-generator --download-metis --download-ml --download-mpich<br>
--download-parmetis --download-scientificpython --download-triangle<br>
--with-clanguage=cxx --with-dynamic-loading --with-shared-libraries<br>
--with-sieve PETSC_ARCH=arch-linux2-cxx-debug<br>
[0]PETSC ERROR:<br>
------------------------------------------------------------------------<br>
[0]PETSC ERROR: DMDestroy() line 234 in<br>
/home/user/Desktop/LIBRARIES/petsc-dev/src/dm/interface/dm.c<br>
application called MPI_Abort(MPI_COMM_WORLD, 73) - process 0<br>
[unset]: aborting job:<br>
application called MPI_Abort(MPI_COMM_WORLD, 73) - process 0<br>
<br>
Ideas?<br>
<span class="HOEnZb"><font color="#888888"><br>
--<br>
Chris Eldred<br>
DOE Computational Science Graduate Fellow<br>
Graduate Student, Atmospheric Science, Colorado State University<br>
B.S. Applied Computational Physics, Carnegie Mellon University, 2009<br>
<a href="mailto:chris.eldred@gmail.com">chris.eldred@gmail.com</a><br>
</font></span></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<br>