[petsc-dev] Issues with DMComplexSetCone
Chris Eldred
chris.eldred at gmail.com
Mon Aug 27 10:58:57 CDT 2012
I need to add:
#include <finclude/petscdmcomplex.h>
in order for Blaise's code to work- but it otherwise works fine.
Thanks!
On Fri, Aug 24, 2012 at 2:57 PM, Matthew Knepley <knepley at gmail.com> wrote:
> On Fri, Aug 24, 2012 at 3:48 PM, Blaise Bourdin <bourdin at lsu.edu> wrote:
>>
>> I am attaching a debugged version.
>
>
> Checked in as DMComplex ex1f90.F
>
> Matt
>
>>
>> Blaise
>>
>>
>>
>> On Aug 24, 2012, at 3:30 PM, Matthew Knepley <knepley at gmail.com> wrote:
>>
>> On Fri, Aug 24, 2012 at 2:16 PM, Chris Eldred <chris.eldred at gmail.com>
>> wrote:
>>>
>>> The following code snippet does not appear to be working:
>>
>>
>> 1) I hate Fortran. A lot. It is not debuggable.
>>
>> 2) Here is my code that should work, but the do loop somehow does not
>> work:
>>
>> program main
>> implicit none
>> !
>> #include <finclude/petsc.h90>
>> #include <finclude/petscdmcomplex.h90>
>> !
>> DM dm
>> PetscInt, target, dimension(4) :: EC
>> PetscInt, pointer :: pEC(:)
>> PetscInt c, firstCell, numCells
>> PetscErrorCode ierr
>>
>> call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
>> call DMComplexCreate(PETSC_COMM_WORLD, dm, ierr)
>> firstCell = 0
>> numCells = 5
>> call DMComplexSetChart(dm, 0, numCells, ierr)
>> do c=firstCell,numCells
>> write(*,*) 'c',c
>> call DMComplexSetConeSize(dm, c, 4)
>> end do
>> call DMSetUp(dm, ierr)
>>
>> EC(1) = 1
>> EC(2) = 2
>> EC(3) = 3
>> EC(4) = 4
>> pEC => EC
>> write(*,*) 'cell',c,pEC
>> c = 0
>> call DMComplexSetCone(dm, c , pEC, ierr)
>> CHKERRQ(ierr)
>> call DMComplexGetCone(dm, c , pEC, ierr)
>> CHKERRQ(ierr)
>> write(*,*) 'cell',c,pEC
>>
>> call DMDestroy(dm,ierr)
>> call PetscFinalize(ierr)
>> end
>>
>>
>> Matt
>>
>>>
>>> PetscInt, dimension(4) :: EC
>>>
>>> EC(1) = 1
>>> EC(2) = 2
>>> EC(3) = 3
>>> EC(4) = 4
>>> write(*,*) 'cell',EC
>>> call DMComplexSetCone(model_mesh, lcell , EC, ierr)
>>> CHKERRQ(ierr)
>>> call DMComplexGetCone(model_mesh, lcell , EC, ierr)
>>> CHKERRQ(ierr)
>>> write(*,*) 'cell', EC
>>>
>>> I get the following as output:
>>>
>>> cell 1 2 3 4
>>> cell 18673792 0 -1 -1
>>>
>>> DMComplexGetConeSize returns 4 as expected. I am using the latest
>>> version of petsc-dev.
>>>
>>> Later on there is a segmentation fault violation, probably related to
>>> this issue.
>>>
>>> Any ideas?
>>>
>>> --
>>> Chris Eldred
>>> DOE Computational Science Graduate Fellow
>>> Graduate Student, Atmospheric Science, Colorado State University
>>> B.S. Applied Computational Physics, Carnegie Mellon University, 2009
>>> chris.eldred at gmail.com
>>
>>
>>
>>
>> --
>> 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
>>
>>
>> --
>> Department of Mathematics and Center for Computation & Technology
>> Louisiana State University, Baton Rouge, LA 70803, USA
>> Tel. +1 (225) 578 1612, Fax +1 (225) 578 4276
>> http://www.math.lsu.edu/~bourdin
>>
>>
>>
>>
>>
>>
>>
>>
>
>
>
> --
> 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
--
Chris Eldred
DOE Computational Science Graduate Fellow
Graduate Student, Atmospheric Science, Colorado State University
B.S. Applied Computational Physics, Carnegie Mellon University, 2009
chris.eldred at gmail.com
More information about the petsc-dev
mailing list