[petsc-users] [petsc-dev] Getting 2d array with updated ghost values from DM global vector

TAY wee-beng zonexo at gmail.com
Sat Jul 7 10:36:08 CDT 2012


Thanks to everyone which helped to fix the bug.

Let me restate my problem.

I have used DMDACreate2d (dof = 2) for my code and solve the equation. Using

/call KSPSolve(ksp_semi,b_rhs_semi_global,velocity_global,ierr)

call 
DMGlobalToLocalBegin(da,velocity_global,INSERT_VALUES,velocity_local,ierr)
call 
DMGlobalToLocalEnd(da,velocity_global,INSERT_VALUES,velocity_local,ierr)

call VecGetArrayF90(velocity_local,p_velocity,ierr)

do j = start_ij(2),end_ij(2)

     do i = start_ij(1),end_ij(1)

         du(i,j)=p_velocity(ij)

         dv(i,j)=p_velocity(ij+1)

         ij=ij+2

     end do

end do/

/call VecRestoreArrayF90(velocity_local,p_velocity,ierr)/

The above du,dv correspond to the  region without the ghost cell. I 
would like to have the ghost cell values as well. I tried using:

/PetscScalar,pointer :: velocity_array(:,:,:)

call DMDAVecGetArrayF90(da,velocity_local,velocity_array,ierr)
call DMDAVecRestoreArrayF90(da,velocity_local,velocity_array,ierr)
/
I checked by looking at the values at :

velocity_array(0,34-1,47-1:48-1)     ->  0 is for 1st dof, -1 becos 
array starts from 0, 47-48 is where the grid divides

the But the values obtained is wrong. So how should I do it?

Yours sincerely,

TAY wee-beng

On 6/7/2012 1:24 PM, Blaise Bourdin wrote:
> On Jul 6, 2012, at 6:02 PM, TAY wee-beng wrote:
>
>> On 3/7/2012 1:23 PM, Satish Balay wrote:
>>> On Tue, 3 Jul 2012, Barry Smith wrote:
>>>
>>>> On Jul 3, 2012, at 3:08 AM, Blaise Bourdin wrote:
>>>>
>>>>> On Jul 3, 2012, at 4:10 AM, Barry Smith wrote:
>>>>>
>>>>>> Blaise,
>>>>>>
>>>>>>   I don't understand why the patch does anything:
>>>>>>
>>>>>> -  *ierr = VecRestoreArray(*v,0);if (*ierr) return;
>>>>>> +  PetscScalar *fa;
>>>>>> +  *ierr = F90Array1dAccess(a,PETSC_SCALAR,(void**)&fa PETSC_F90_2PTR_PARAM(ptrd));
>>>>>> +  *ierr = VecRestoreArray(*v,&fa);if (*ierr) return;
>>>>>> *ierr = F90Array1dDestroy(&a,PETSC_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
>>>>>>
>>>>>> All that passing &fa into VecRestoreArray() does is cause fa to be zeroed. Why would that have any affect on anything?
>>>>> Not sure either, I quite don't understand this code, but I noticed that the logic of VecRestoreArrayF90 was different from that of DMDAVecRestoreArrayF90
>>>>>
>>>>> src/vec/vec/interface/f90-custom/zvectorf90.c:33
>>>>> PetscScalar *fa;
>>>>> *__ierr = F90Array1dAccess(ptr,PETSC_SCALAR,(void**)&fa PETSC_F90_2PTR_PARAM(ptrd));if (*__ierr) return;
>>>>> *__ierr = F90Array1dDestroy(ptr,PETSC_SCALAR PETSC_F90_2PTR_PARAM(ptrd));if (*__ierr) return;
>>>>       It could be the above line is important; but the Accesser and restore array are not.
>>>>
>>>>       I'll have Satish apply the patch.
>>> pushed to petsc-3.3 [petsc-dev will get this update]
>>>
>>> Satish
>> Hi,
>>
>> I just tested with the latest petsc-dev but it doesn't work in intel linux for ex11f90. Has the patch been applied?
> Try to clone petsc-3.3  from the mercurial repositoryhttp://petsc.cs.iit.edu/petsc/releases/petsc-3.3/  It looks like the first patch has not made its way to the tarball or petsc-dev yet.
> You can also apply the patch manually:
> cd $PETSC_DIR
> patch -p1 < DMDAVecGetArrayF90.patch
>
>> Also, is there any chance of it working under 3d with multiple dof since that's what I'm using and I have other problems with gfortran. Lastly, if the patch is applied, it works with 3d da with 1 dof? Is that right?
> I am sending another patch that should take care of the 3d case with >1 dof to the developers list.
>
> Blaise
>


-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20120707/05d32df5/attachment.html>


More information about the petsc-users mailing list