[petsc-dev] DMViewFromOptions() in fortran

Matthew Knepley knepley at gmail.com
Fri Apr 4 10:54:59 CDT 2014


On Tue, Mar 25, 2014 at 10:11 PM, Matthew Knepley <knepley at gmail.com> wrote:

> On Tue, Mar 25, 2014 at 8:52 PM, Adrian Croucher <
> a.croucher at auckland.ac.nz> wrote:
>
>>
>>
>> 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.
>>
>
> Below, you should not need an if, and can just use DMView(). If that does
> not work in parallel, something is wrong.
> Likely on the interpretation of the default viewer. Let me try it.
>

I have written the Fortran interface for PetscObjectViewFromOptions() and
pushed to next, so
you should be able to use that now.

  Thanks,

     Matt


>    Matt
>
>
>> 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
>>
>>
>
>
> --
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> -- Norbert Wiener
>



-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20140404/556b78e3/attachment.html>


More information about the petsc-dev mailing list