[petsc-users] updating values in a DA Global array

Barry Smith bsmith at mcs.anl.gov
Tue Jun 22 14:23:40 CDT 2010


On Jun 22, 2010, at 7:59 AM, Mark Cheeseman wrote:

> Hi,
> 
> I am trying to write a PETSc program in FORTRAN90 where I need to update a single value in a global distributed array.  I know the global coordinates of the position that needs to be updated in the global array but I cannot get the mapping from the local vector correct.  In this case, I am working on a domain with global dimensions [arraysize(1),arraysize(2),arraysize(3)] and I want to alter a single point in the global distributed array, uGLOBAL, at the global position [arraysize(1)/2-1,arraysize(2)-1,3].  I cannot seem to be able to do this... what am I doing wrong?
> 
> ...
> DA da
> Vec uGLOBAL, uLOCAL, tmp
> PetscErrorCode ierr
> PetscScalar, pointer :: xx
> PetscInt rank, source_rank, i,j,k, row
> 
> ....
> 
> call MPI_Comm_rank( PETSC_COMM_WORLD, rank, ierr )
> call DACreate3d( PETSC_COMM_WORLD, DA_NONPERIODIC, DA_STENCIL_BOX,       &
>                       arraysize(1), arraysize(2), arraysize(3), PETSC_DECIDE, &
>                       PETSC_DECIDE, PETSC_DECIDE, 1, 5, PETSC_NULL_INTEGER,   &
>                       PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, da, ierr)
> call DACreateGlobalVector( da, pNOW, ierr )
> call DAGetCorners( da, xs, ys, zs, xl, yl, zl, ierr )
> 
>     do i = xs,xs+xl-1
>         if ( i.eq.arraysize(1)/2-1 ) then
>            src_loc(1) = i
>            do j = ys,ys+yl-1
>               if ( j.eq.arraysize(2)/2-1 ) then
>                  src_loc(2) = j
>                  do k = zs,zs+zl-1
>                     if ( k.eq.3 ) then
>                        src_loc(2) = j
>                        source_rank = rank
>                     endif
>                  enddo
>               endif
>            enddo
>         endif
>      enddo
> 
> call DAGetLocalVector( da, uLOCAL, ierr )
> call VecGetArrayF90( uLOCAL, xx, ierr )
> 
   Use VecGetArrayF90() directly on the global vector and set the value in there.  

   Or instead of this you can determine the global location of the entry in the "natural ordering" and then use DAGetAO() followed by AOApplicationToPetsc()  with the global location in the natural ordering to convert to the PETSc ordering and then call VecSetValues() with the global vector followed by VecAssemblyBegin() and VecAssemblyEnd().


> if  ( rank.eq.source_rank )  then
>      row = 15
>      xx(row) = pressure
> endif
> 
> call VecRestoreArrayF90( uLOCAL, xx, ierr )
> call DALocalToGlobal( da, uLOCAL, ADD_VALUES, uGLOBAL, ierr )
> call DARestoreLocalVector( da, uLOCAL, ierr )
> 

> 
> Any help would be greatly appreciated.  Is there a way to directly access the global distributed array (uGLOBAL) instead of working through the intermediate local array (uLOCAL)?  I have tried the approach
> 
> call DAVecGetArray( da, uGLOBAL, xx, ierr )
> xx = ....
> call DAVecRestoreArray( da, uGLOBAL, xx, ierr )
> 
> Unfortunately, the compile always fails with the message that the DAVecRestoreArray function cannot be found.

    These are not supported in Fortran.

  Barry


> 
> Thank you,
> Mark
> 
> -- 
> Mark Patrick Cheeseman
> 
> Computational Scientist
> KSL (KAUST Supercomputing Laboratory)
> Building 1, Office #126
> King Abdullah University of Science & Technology
> Thuwal 23955-6900
> Kingdom of Saudi Arabia
>         
> EMAIL   : mark.cheeseman at kaust.edu.sa
> PHONE : +966   (2) 808 0221 (office)
>               +966 (54) 470 1082 (mobile)
> SKYPE : mark.patrick.cheeseman

-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20100622/84fb6e4c/attachment.htm>


More information about the petsc-users mailing list