[petsc-dev] DMPlexComputeCellGeometryFVM from fortran

Adrian Croucher a.croucher at auckland.ac.nz
Wed Mar 19 19:35:44 CDT 2014


On 19/03/14 02:14, Matthew Knepley wrote:
>
>     hi, I'm trying to use the DMPlexComputeCellGeometryFVM() function
>     from a fortran program, but when I compile and link it complains:
>
>     undefined reference to `dmplexcomputecellgeometryfvm_'
>
>     I am running the 'next' branch of PETSc (having previously tried
>     the master branch) and configured with
>     --with-fortran-interfaces=1. At the top of my code I have #include
>     <finclude/petsc.h90>.
>
>     I noticed the fortran interfaces to some other functions (e.g.
>     DMPlexCreateFromDAG) also seemed to be missing in older versions
>     (e.g. the 'maint' branch)- is this one perhaps still to be added too?
>
>
> This is my fault. I will push it to 'next' today.

Thanks for doing that. I can access DMPlexComputeCellGeometryFVM() from 
fortran now.

However at runtime it is giving me a segfault- perhaps I am not 
declaring the centroid and normal arrays properly? I also tried 
declaring them as ordinary fortran arrays instead of pointers, but that 
also segfaulted. It works OK when I call it from C. My fortran code is 
below- any clues much appreciated. The error it gives is:

[0]PETSC ERROR: Caught signal number 11 SEGV: Segmentation Violation, 
probably memory access out of range
[0]PETSC ERROR: Try option -start_in_debugger or -on_error_attach_debugger
[0]PETSC ERROR: or see 
http://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind[0]PETSC 
ERROR: or try http://valgrind.org on G
NU/linux and Apple Mac OS X to find memory corruption errors
[0]PETSC ERROR: likely location of problem given in stack below
[0]PETSC ERROR: ---------------------  Stack Frames 
------------------------------------
[0]PETSC ERROR: Note: The EXACT line numbers in the stack are not available,
[0]PETSC ERROR:       INSTEAD the line number of the start of the function
[0]PETSC ERROR:       is given.
[0]PETSC ERROR: [0] F90Array1dAccess line 70 
/home/acro018/software/PETSc/code/src/sys/f90-src/f90_cwrap.c
[0]PETSC ERROR: --------------------- Error Message 
--------------------------------------------------------------
[0]PETSC ERROR: Signal received
[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-4938-g89a3bfc  
GIT Date: 2014-03-19 17:00:11 +0000
[0]PETSC ERROR: testdag on a linux-gnu-c-opt named des108 by acro018 Thu 
Mar 20 13:22:15 2014
[0]PETSC ERROR: Configure options --download-netcdf --download-exodusii 
--with-hdf5-dir=/usr/lib --with-numpy-dir=
[0]PETSC ERROR: #1 User provided function() line 0 in  unknown file

- Adrian

--
program main

! setting up 3-D DMPlex using DMPlexCreateFromDAG(), DMPlexInterpolate() and
! DMPlexComputeCellGeometryFVM()

   implicit none

#include <finclude/petsc.h90>

   DM :: dm, dmi
   PetscInt :: dim = 3
   PetscInt :: depth = 1
   PetscErrorCode :: ierr
   PetscInt, allocatable :: numPoints(:)
   PetscInt, allocatable :: coneSize(:)
   PetscInt, allocatable :: cones(:)
   PetscInt, allocatable :: coneOrientations(:)
   PetscScalar, allocatable :: vertexCoords(:)
   PetscReal ::  vol
   PetscReal, pointer :: centroid(:), normal(:)
   PetscInt :: i

   numPoints = [12, 2]
   coneSize  = [8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0]
   cones = [2,5,4,3,6,7,8,9,  3,4,11,10,7,12,13,8]
   coneOrientations = [0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0]
   vertexCoords = [ &
        -0.5,0.0,0.0, 0.0,0.0,0.0, 0.0,1.0,0.0, -0.5,1.0,0.0, &
        -0.5,0.0,1.0, 0.0,0.0,1.0, 0.0,1.0,1.0, -0.5,1.0,1.0, &
         0.5,0.0,0.0, 0.5,1.0,0.0, 0.5,0.0,1.0,  0.5,1.0,1.0 ]

   call PetscInitialize(PETSC_NULL_CHARACTER,ierr); CHKERRQ(ierr)

   call DMPlexCreate(PETSC_COMM_WORLD, dm, ierr); CHKERRQ(ierr)
   call PetscObjectSetName(dm, "testplex", ierr); CHKERRQ(ierr)
   call DMPlexSetDimension(dm, dim, ierr); CHKERRQ(ierr)

   call DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, &
        coneOrientations, vertexCoords, ierr); CHKERRQ(ierr)

   call DMPlexInterpolate(dm, dmi, ierr); CHKERRQ(ierr)
   call DMPlexCopyCoordinates(dm, dmi, ierr); CHKERRQ(ierr)
   call DMDestroy(dm, ierr); CHKERRQ(ierr)
   dm = dmi

   call DMView(dm, PETSC_VIEWER_STDOUT_WORLD, ierr); CHKERRQ(ierr)

   do i = 0, 1
      call DMPlexComputeCellGeometryFVM(dm, i, vol, centroid, normal, ierr)
      CHKERRQ(ierr)
      write(*, '(a, i2, a, f8.4, a, 3(f8.4, 1x))'), &
           'cell:', i, 'volume:', vol, 'centroid:', centroid(1), 
centroid(2), centroid(3)
   end do

   call DMDestroy(dm, ierr); CHKERRQ(ierr)
   call PetscFinalize(ierr); CHKERRQ(ierr)
   deallocate(numPoints, coneSize, cones, coneOrientations, vertexCoords)

end program main

-- 
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/20140320/ac551507/attachment.html>


More information about the petsc-dev mailing list