[petsc-dev] DMViewFromOptions() in fortran
Adrian Croucher
a.croucher at auckland.ac.nz
Tue Mar 25 20:52:54 CDT 2014
>
> PETSC_STATIC_INLINE PetscErrorCode DMViewFromOptions(DM A,const char
> prefix[],const char name[]) {return
> PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
>
> so I think it might be better to just call the PetscObject version
> from Fortran.
I tried that but got a similar linker error from
PetscObjectViewFromOptions().
> Can you run an example, like SNES ex12, and give -dm_view in parallel?
Yes, that works. The C version of my test code also works- that's why I
figured it might be just a missing fortran interface. My fortran test
code is below.
Cheers, Adrian
--
program main
! distributing a DMPlex
implicit none
#include <finclude/petsc.h90>
DM :: dm, interp_dm, dist_dm
PetscInt, parameter :: dim = 3
PetscInt :: depth = 1
PetscErrorCode :: ierr
PetscInt :: rank, size
PetscInt, allocatable :: numPoints(:)
PetscInt, allocatable :: coneSize(:)
PetscInt, allocatable :: cones(:)
PetscInt, allocatable :: coneOrientations(:)
PetscScalar, allocatable :: vertexCoords(:)
character(5) :: partitioner = "chaco"
call PetscInitialize(PETSC_NULL_CHARACTER,ierr); CHKERRQ(ierr)
call MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr); CHKERRQ(ierr);
call MPI_Comm_size(PETSC_COMM_WORLD, size, 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)
if (rank == 0) then
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 DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, &
coneOrientations, vertexCoords, ierr)
else
numPoints = [0, 0, 0]
coneSize = [0]
cones = [0]
coneOrientations = [0]
vertexCoords = [0]
call DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, &
cones, coneOrientations, vertexCoords, ierr);
end if
CHKERRQ(ierr)
call DMPlexInterpolate(dm, interp_dm, ierr); CHKERRQ(ierr)
call DMPlexCopyCoordinates(dm, interp_dm, ierr); CHKERRQ(ierr)
call DMDestroy(dm, ierr); CHKERRQ(ierr)
dm = interp_dm
call DMPlexDistribute(dm, partitioner, 0, PETSC_NULL_OBJECT, dist_dm,
ierr)
CHKERRQ(ierr)
if (size > 1) then
call DMViewFromOptions(dist_dm, PETSC_NULL_OBJECT, "-dm_view", ierr)
CHKERRQ(ierr)
else ! one processor
call DMView(dm, PETSC_VIEWER_STDOUT_WORLD, ierr); CHKERRQ(ierr)
end if
call DMDestroy(dm, ierr); CHKERRQ(ierr)
dm = dist_dm
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/20140326/5e14e15d/attachment.html>
More information about the petsc-dev
mailing list