[petsc-users] Question on DMPlexCreateSection for Fortran
Danyang Su
danyang.su at gmail.com
Thu Feb 22 17:10:15 CST 2018
Hi Matt,
Just to let you know that after updating to PETSc 3.8.3 version, the
DMPlexCreateSection in my code now works.
One more question, what is the PETSC_NULL_XXX for IS pointer, as shown
below, in C code, it just pass NULL, but in fortran, what is the name of
null object for pBcCompIS and pBcPointIS.
call DMPlexCreateSection(dmda_flow%da,dmda_flow%dim, &
numFields,pNumComp,pNumDof, &
numBC,pBcField, &
pBcCompIS,pBcPointIS, &
PETSC_NULL_IS,section,ierr)
CHKERRQ(ierr)
Thanks,
Danyang
On 18-02-21 09:22 AM, Danyang Su wrote:
>
> Hi Matt,
>
> To test the Segmentation Violation problem in my code, I modified the
> example ex1f90.F to reproduce the problem I have in my own code.
>
> If use DMPlexCreateBoxMesh to generate the mesh, the code works fine.
> However, if I use DMPlexCreateGmshFromFile, using the same mesh
> exported from "DMPlexCreateBoxMesh", it gives Segmentation Violation
> error.
>
> Did I miss something in the input mesh file? My first guess is the
> label "marker" used in the code, but I couldn't find any place to set
> this label.
>
> Would you please let me know how to solve this problem. My code is
> done in a similar way as ex1f90, it reads mesh from external file or
> creates from cell list, distributes the mesh (these already work), and
> then creates sections and sets ndof to the nodes.
>
> Thanks,
>
> Danyang
>
>
> On 18-02-20 10:07 AM, Danyang Su wrote:
>> On 18-02-20 09:52 AM, Matthew Knepley wrote:
>>> On Tue, Feb 20, 2018 at 12:30 PM, Danyang Su <danyang.su at gmail.com
>>> <mailto:danyang.su at gmail.com>> wrote:
>>>
>>> Hi All,
>>>
>>> I tried to compile the DMPlexCreateSection code but got error
>>> information as shown below.
>>>
>>> Error: Symbol 'petsc_null_is' at (1) has no IMPLICIT type
>>>
>>> I tried to use PETSC_NULL_OBJECT instead of PETSC_NULL_IS, then
>>> the code can be compiled but run into Segmentation Violation
>>> error in DMPlexCreateSection.
>>>
>>> From the webpage
>>>
>>> http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMPLEX/DMPlexCreateSection.html
>>>
>>>
>>> The F90 version is DMPlexCreateSectionF90. Doing this with F77
>>> arrays would have been too painful.
>> Hi Matt,
>>
>> Sorry, I still cannot compile the code if use DMPlexCreateSectionF90
>> instead of DMPlexCreateSection. Would you please tell me in more
>> details?
>>
>> undefined reference to `dmplexcreatesectionf90_'
>>
>> then I #include <petsc/finclude/petscdmplex.h90>, but this throws
>> more error during compilation.
>>
>>
>> Included at
>> /home/dsu/Soft/PETSc/petsc-3.7.5/include/petsc/finclude/petscdmplex.h90:6:
>> Included at ../../solver/solver_ddmethod.F90:62:
>>
>> PETSCSECTION_HIDE section
>> 1
>> Error: Unclassifiable statement at (1)
>> /home/dsu/Soft/PETSc/petsc-3.7.5/include/petsc/finclude/ftn-custom/petscdmplex.h90:167.10:
>> Included at
>> /home/dsu/Soft/PETSc/petsc-3.7.5/include/petsc/finclude/petscdmplex.h90:6:
>> Included at ../../solver/solver_ddmethod.F90:62:
>>
>> PETSCSECTION_HIDE section
>> 1
>> Error: Unclassifiable statement at (1)
>> /home/dsu/Soft/PETSc/petsc-3.7.5/include/petsc/finclude/ftn-custom/petscdmplex.h90:179.10:
>> Included at
>> /home/dsu/Soft/PETSc/petsc-3.7.5/include/petsc/finclude/petscdmplex.h90:6:
>> Included at ../../solver/solver_ddmethod.F90:62:
>>
>>>
>>> Thanks,
>>>
>>> Matt
>>>
>>> dmda_flow%da is distributed dm object that works fine.
>>>
>>> The fortran example I follow is
>>> http://www.mcs.anl.gov/petsc/petsc-dev/src/dm/impls/plex/examples/tutorials/ex1f90.F90
>>> <http://www.mcs.anl.gov/petsc/petsc-dev/src/dm/impls/plex/examples/tutorials/ex1f90.F90>.
>>>
>>>
>>> What parameters should I use if passing null to bcField,
>>> bcComps, bcPoints and perm.
>>>
>>> PetscErrorCode
>>> <https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/Vec/PetscErrorCode.html#PetscErrorCode> DMPlexCreateSection
>>> <https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/DMPLEX/DMPlexCreateSection.html#DMPlexCreateSection>(DM
>>> <https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/DM/DM.html#DM> dm,PetscInt
>>> <https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/Vec/PetscInt.html#PetscInt> dim,PetscInt
>>> <https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/Vec/PetscInt.html#PetscInt> numFields,constPetscInt
>>> <https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/Vec/PetscInt.html#PetscInt> numComp[],constPetscInt
>>> <https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/Vec/PetscInt.html#PetscInt> numDof[],PetscInt
>>> <https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/Vec/PetscInt.html#PetscInt> numBC,constPetscInt
>>> <https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/Vec/PetscInt.html#PetscInt> bcField[],
>>> constIS
>>> <https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/IS/IS.html#IS> bcComps[], constIS
>>> <https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/IS/IS.html#IS> bcPoints[],IS
>>> <https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/IS/IS.html#IS> perm,PetscSection
>>> <https://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/IS/PetscSection.html#PetscSection> *section)
>>>
>>> #include <petsc/finclude/petscis.h>
>>> #include <petsc/finclude/petscis.h90>
>>> #include <petsc/finclude/petscdmplex.h>
>>>
>>> ...
>>>
>>> #ifdef USG
>>> numFields = 1
>>> numComp(1) = 1
>>> pNumComp => numComp
>>>
>>> do i = 1, numFields*(dmda_flow%dim+1)
>>> numDof(i) = 0
>>> end do
>>> numDof(0*(dmda_flow%dim+1)+1) = dmda_flow%dof
>>> pNumDof => numDof
>>>
>>> numBC = 0
>>>
>>> call DMPlexCreateSection(dmda_flow%da,dmda_flow%dim, &
>>> numFields,pNumComp,pNumDof, &
>>> numBC,PETSC_NULL_INTEGER, &
>>> PETSC_NULL_IS,PETSC_NULL_IS, & !Error here
>>> PETSC_NULL_IS,section,ierr)
>>> CHKERRQ(ierr)
>>>
>>> call PetscSectionSetFieldName(section,0,'flow',ierr)
>>> CHKERRQ(ierr)
>>>
>>> call DMSetDefaultSection(dmda_flow%da,section,ierr)
>>> CHKERRQ(ierr)
>>>
>>> call PetscSectionDestroy(section,ierr)
>>> CHKERRQ(ierr)
>>> #endif
>>>
>>> Thanks,
>>>
>>> Danyang
>>>
>>>
>>>
>>>
>>> --
>>> 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
>>>
>>> https://www.cse.buffalo.edu/~knepley/
>>> <http://www.caam.rice.edu/%7Emk51/>
>>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20180222/0d1623b3/attachment.html>
More information about the petsc-users
mailing list