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

Barry Smith bsmith at mcs.anl.gov
Mon Jul 2 19:10:40 CDT 2012


   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?

   Thanks

     Barry

On Jul 2, 2012, at 10:10 AM, Blaise Bourdin wrote:

> Hi,
> 
> There appears to be a bug in DMDAVecRestoreArrayF90. It is probably only triggered when the intel compilers. gfortran and intel seem to have very different internal implementations of fortran90 allocatable arrays. 
> 
> Developers, can you check if the attached patch makes sense? It will not fix the case of a 3d da with dof>1 since F90Array4dAccess is not implemented. Other than that, it seems to fix ex11f90 under linux and mac OS
> 
> <DMDAVecRestoreArrayF90.patch>
> 
> Blaise
> 
> 
> 
> On Jul 2, 2012, at 5:41 PM, TAY wee-beng wrote:
> 
>> On 2/7/2012 2:49 PM, Matthew Knepley wrote:
>>> On Mon, Jul 2, 2012 at 3:58 AM, TAY wee-beng <zonexo at gmail.com> wrote:
>>> Hi,
>>> 
>>> I have used DMDACreate2d for my code and then use:
>>> 
>>> call DMLocalToGlobalBegin(da,b_rhs_semi_local,INSERT_VALUES,b_rhs_semi_global,ierr)
>>> 
>>> call DMLocalToGlobalEnd(da,b_rhs_semi_local,INSERT_VALUES,b_rhs_semi_global,ierr)
>>> 
>>> to construct the global DM vector b_rhs_semi_global
>>> 
>>> Now I want to get the values with ghost values in a 2d array locally which is declared as:
>>> 
>>> real(8), allocatable :: array2d(:,:)
>>> 
>>> I guess I should use DMDAGetGhostCorners to get the corressponding indices and allocate it.  But what should I do next? How can I use something like VecGetArrayF90 to get to the pointer to access the local vector?
>>> 
>>> I can't use DMDAVecGetArrayF90/DMDAVecRestoreArrayF90 since I'm using intel fortran and they can't work. I can't use gfortran at the moment since I've problems with HYPRE with           gfortran in 3D.
>>> 
>>> Are you certain of this? That used to be true, but the current version should work for any F90.
>>> 
>>>    Matt
>> 
>> I just tested 3.3-p1 and it still doesn't work (example ex11f90 in dm). Is there a chance petsc-dev can work?
>>>  
>>> Thanks
>>> 
>>> -- 
>>> 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
>> 
>> 
> 
> -- 
> 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
> 
> 
> 
> 
> 
> 
> 



More information about the petsc-users mailing list