[petsc-users] Declaring struct to represent field for dof > 1 for DM in Fortran

TAY wee-beng zonexo at gmail.com
Tue Jul 10 13:22:20 CDT 2012


On 10/7/2012 6:05 PM, Matthew Knepley wrote:
> On Tue, Jul 10, 2012 at 11:03 AM, TAY wee-beng <zonexo at gmail.com 
> <mailto:zonexo at gmail.com>> wrote:
>
>
>     Yours sincerely,
>
>     TAY wee-beng
>
>     On 10/7/2012 2:07 PM, Matthew Knepley wrote:
>>     On Tue, Jul 10, 2012 at 5:39 AM, TAY wee-beng <zonexo at gmail.com
>>     <mailto:zonexo at gmail.com>> wrote:
>>
>>         Hi,
>>
>>         I read in the manual in page 50 that it's recommended to
>>         declare struct to represent field for dof > 1 for DM.
>>
>>
>>     We mean C struct. C makes it easy (just use a pointer type cast).
>>     Fortran makes it hard unfortunately.
>>
>>        Matt
>     Ok, I'll try to use another mtd.
>
>     Btw, if I declare:
>
>     /PetscScalar,pointer :: array2(:,:,:)
>
>     with DMDACreate2d using dof = 2,
>
>     call DMDAVecGetArrayF90(da,x_local,array2,ierr)
>
>     access array2 ....
>
>     call DMDAVecRestoreArrayF90(da,x_local,array2,ierr)/
>
>     How is the memory for "array2" allocated ? Is it allocated all the
>     time, or only between the DMDAVecGetArrayF90 and
>     DMDAVecRestoreArrayF90?
>
>     Also, can I "reuse" array2? For e.g., now for y_local:
>
>     /call DMDAVecGetArrayF90(da,y_local,array2,ierr)
>
>     access array2 ..../ /
>
>     call DMDAVecRestoreArrayF90(da,y_local,array2,ierr)/
>
>
> The right thing to do here is to implement DMDAVecGetArratDOFF90().
>
>     Matt

Hi,

Do you mean DMDAVecGetArrayDOFF90 ? I tried to compile but it gives the 
error during linking:

1>dm_test2d.obj : error LNK2019: unresolved external symbol 
DMDAVECGETARRAYDOFF90 referenced in function MAIN__

Also from the manual of DMDAVecGetArray, it says:
/
Fortran Notes: From Fortran use DMDAVecGetArrayF90() and pass for the 
array type PetscScalar 
<http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sys/PetscScalar.html#PetscScalar>,pointer 
:: array(:,...,:) of the appropriate dimension. For a DMDA created with 
a dof of 1 use the dimension of the DMDA, for a DMDA created with a dof 
greater than 1 use one more than the dimension of the DMDA. The order of 
the indices is array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) (when dof is 1) 
otherwise array(1:dof,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values 
are obtained from DMDAGetCorners 
<http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMDAGetCorners.html#DMDAGetCorners>() 
for a global array or DMDAGetGhostCorners 
<http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMDAGetGhostCorners.html#DMDAGetGhostCorners>() 
for a local array. Include finclude/petscdmda.h90 to access this routine. /

I just tried with dof = 2 and there's no problem. However, the manual 
says that for dof > 1, the array is 
/array(1:dof,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)/.

Should it be /array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)/ instead? 
I had problems with the former, but the latter works fine.

Also, I'm still not sure how the memory is allocated. If I have:

/Vec x_local

PetscScalar,pointer :: array2(:,:,:)

with DMDACreate2d using dof = 2,

call DMDAVecGetArrayF90(da,x_local,array2,ierr)

access array2 ....

call DMDAVecRestoreArrayF90(da,x_local,array2,ierr)/

How is the memory for "array2" allocated ? Is it allocated all the time, 
or only between the DMDAVecGetArrayF90 and DMDAVecRestoreArrayF90?

Thanks!


>     Thank you!
>
>>
>>         I'm using Fortran and for testing, I use dof = 1 and write as:
>>
>>         /type field
>>
>>         //PetscScalar//u        (or real(8) :: u)
>>
>>         end type field
>>
>>         type(field), pointer :: field_u(:,:)/
>>
>>         When I tried to use :
>>
>>         /call DMDAVecGetArrayF90(da,x_local,field_u,ierr)/
>>
>>         I got the error : There is no matching specific subroutine
>>         for this generic subroutine call. [DMDAVECGETARRAYF90]
>>
>>         The da, x_local has been defined with the specific DM
>>         routines. It worked if I use :
>>
>>         /PetscScalar,pointer :: array(:,:) and
>>
>>         call DMDAVecGetArrayF90(da,x_local,array,ierr)/
>>
>>         May I know what did I do wrong?
>>
>>
>>         -- 
>>         Yours sincerely,
>>
>>         TAY wee-beng
>>
>>
>>
>>
>>     -- 
>>     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-users/attachments/20120710/d310e07c/attachment-0001.html>


More information about the petsc-users mailing list