[petsc-dev] DMGetDS in fortran
Adrian Croucher
a.croucher at auckland.ac.nz
Thu Nov 27 21:19:44 CST 2014
hi
I'm just getting back to using DMPlex after a few months busy on other
stuff, and looking at the TS ex11 problem, a few things have changed in
the meantime.
TS ex11 used to set up a PetscSection manually to define the data
layout, for constructing vectors on it. Now it seems this is handled by
creating a PetscFV object and adding it to the PetscDS in the DMPlex.
I set up a simplified code to test this in isolation, and it works in C
but not in Fortran (code below). It gets this error- looks like maybe
something wrong with the fortran interface:
[0]PETSC ERROR: --------------------- Error Message
--------------------------------------------------------------
[0]PETSC ERROR: Null argument, when expecting valid pointer
[0]PETSC ERROR: Null Pointer: Parameter # 2
[0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html
for trouble shooting.
[0]PETSC ERROR: Petsc Development GIT revision: v3.5.2-1240-gd4bc1ae
GIT Date: 2014-11-25 11:37:47 -0600
[0]PETSC ERROR: ./testplexf on a linux-gnu-c-opt named des108 by acro018
Fri Nov 28 16:02:22 2014
[0]PETSC ERROR: Configure options --download-netcdf --download-exodusii
--with-hdf5-dir=/usr/lib
--with-numpy-dir=/usr/lib/pymodules/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-triangl
e --with-fortran-interfaces=1 --download-ptscotch --download-chaco
--download-ctetgen
[0]PETSC ERROR: #1 DMGetDS() line 3407 in
/home/acro018/software/PETSc/code/src/dm/interface/dm.c
Also, any chance of a fortran interface to DMPlexCreateFromFile()? It
appears to be missing at the moment.
Regards, Adrian
--
program testplex
implicit none
#include <finclude/petscsys.h>
#include <finclude/petscvec.h>
#include <finclude/petscdm.h>
PetscInt, parameter :: dim = 3, dof = 1
MPI_Comm :: comm
DM :: dm, dist_dm, ghost_dm
PetscErrorCode :: ierr
character(40) :: filename = "column.exo"
PetscInt :: localsize
PetscFV :: fvm
PetscDS :: ds
Vec :: X
call PetscInitialize(PETSC_NULL_CHARACTER, ierr); CHKERRQ(ierr)
comm = PETSC_COMM_WORLD
call DMPlexCreateExodusFromFile(comm, filename, PETSC_TRUE, dm,
ierr); CHKERRQ(ierr)
call DMPlexSetAdjacencyUseCone(dm, PETSC_TRUE, ierr); CHKERRQ(ierr)
call DMPlexSetAdjacencyUseClosure(dm, PETSC_FALSE, ierr); CHKERRQ(ierr)
call DMPlexDistribute(dm, 0, PETSC_NULL_OBJECT, dist_dm, ierr)
CHKERRQ(ierr)
if (dist_dm /= 0) then
call DMDestroy(dm, ierr); CHKERRQ(ierr)
dm = dist_dm
end if
call DMPlexConstructGhostCells(dm, PETSC_NULL_CHARACTER,
PETSC_NULL_INTEGER, &
ghost_dm, ierr); CHKERRQ(ierr)
if (ghost_dm /= 0) then
call DMDestroy(dm, ierr); CHKERRQ(ierr);
dm = ghost_dm
end if
call DMView(dm, PETSC_VIEWER_STDOUT_WORLD, ierr); CHKERRQ(ierr)
call PetscFVCreate(comm, fvm, ierr); CHKERRQ(ierr)
call PetscFVSetNumComponents(fvm, dof, ierr); CHKERRQ(ierr)
call PetscFVSetSpatialDimension(fvm, dim, ierr); CHKERRQ(ierr)
call DMGetDS(dm, ds, ierr); CHKERRQ(ierr)
call PetscDSAddDiscretization(ds, fvm, ierr); CHKERRQ(ierr)
call DMCreateGlobalVector(dm, X, ierr); CHKERRQ(ierr)
call VecGetLocalSize(X, localsize, ierr); CHKERRQ(ierr);
print *, "Local size:", localsize
call DMDestroy(dm, ierr); CHKERRQ(ierr)
call PetscFinalize(ierr); CHKERRQ(ierr)
end program testplex
--
Dr Adrian Croucher
Senior Research Fellow
Department of Engineering Science
University of Auckland, New Zealand
email: a.croucher at auckland.ac.nz
tel: +64 (0)9 923 84611
More information about the petsc-dev
mailing list