[petsc-dev] calling DMPlexPointLocalRef() from fortran
Adrian Croucher
a.croucher at auckland.ac.nz
Tue Apr 8 00:18:57 CDT 2014
On 08/04/14 06:52, Matthew Knepley wrote:
>
> I got your code running. There were a few code problems (fixed in my
> attached version), but the
> problem you are seeing above is a bug in our system for generating
> Fortran wrappers. The fix
> will be pushed out soon, but I am attaching dmf.c which you can put in
> src/dm/interface/ftn-auto
> and rebuild and things should work for you.
Thanks very much. That did fix the error I was getting yesterday (once I
also removed the call to DMSetPointSF()).
Next I am trying to visualize the partition vector via VTK, using the
subroutine below (again modelled on what is done in TS ex 33):
-------
subroutine VTK_output(filename, v)
! Writes the vector v to a VTK file
implicit none
#include <finclude/petsc.h90>
character(*), intent(in) :: filename
Vec, intent(in) :: v
! Locals:
PetscErrorCode :: ierr
PetscViewer :: viewer
call PetscViewerCreate(PETSC_COMM_WORLD, viewer, ierr); CHKERRQ(ierr)
call PetscViewerSetType(viewer, PETSCVIEWERVTK, ierr); CHKERRQ(ierr)
call PetscViewerFileSetName(viewer, filename, ierr); CHKERRQ(ierr)
call VecView(v, viewer, ierr); CHKERRQ(ierr)
call PetscViewerDestroy(viewer, ierr); CHKERRQ(ierr)
end subroutine VTK_output
------------
However it crashes with this error message:
[0]PETSC ERROR: --------------------- Error Message
--------------------------------------------------------------
[0]PETSC ERROR: Null argument, when expecting valid pointer
[0]PETSC ERROR: Null Object: Parameter # 1
[0]PETSC ERROR: See
http://http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble
shooting.
[0]PETSC ERROR: Petsc Development GIT revision: v3.4.4-5339-g8a3122f
GIT Date: 2014-04-07 17:35:29 -0500
[0]PETSC ERROR: testplex on a linux-gnu-c-opt named des108 by acro018
Tue Apr 8 17:12:35 2014
[0]PETSC ERROR: Configure options --download-netcdf --download-exodusii
--with-hdf5-dir=/usr/lib --with-numpy-dir=/usr/lib/pymodul
es/python2.7/numpy
--with-scientificpython-dir=/usr/lib/python2.7/dist-packages/scipy
--with-petsc4py-dir=/usr/local/lib/python2.7
/dist-packages/petsc4py --download-triangle --with-fortran-interfaces=1
--download-ptscotch --download-chaco --download-ctetgen
[0]PETSC ERROR: #1 VecGetDM() line 200 in
/home/acro018/software/PETSc/code/src/dm/interface/dm.c
[0]PETSC ERROR: #2 DMPlexInsertBoundaryValuesFVM() line 403 in
/home/acro018/software/PETSc/code/src/dm/impls/plex/plexfem.c
[0]PETSC ERROR: #3 VecView_Plex_Local() line 235 in
/home/acro018/software/PETSc/code/src/dm/impls/plex/plex.c
[0]PETSC ERROR: #4 VecView() line 666 in
/home/acro018/software/PETSc/code/src/vec/vec/interface/vector.c
Would this be expected to work in fortran at the moment? Or have I just
messed it up?
I also came across another function that seems to be lacking a fortran
interface: DMPlexConstructGhostCells(). Not totally sure if I need that,
but again was just trying to follow what is done in TS ex 33.
Cheers, Adrian
--
Dr Adrian Croucher
Department of Engineering Science
University of Auckland
New Zealand
tel 64-9-373-7599 ext 84611
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20140408/d3a05e0b/attachment.html>
More information about the petsc-dev
mailing list