[petsc-users] VecView doesn't work properly with DA global vectors
Barry Smith
bsmith at mcs.anl.gov
Mon Mar 21 20:28:33 CDT 2011
On Mar 21, 2011, at 8:11 PM, Алексей Рязанов wrote:
> Hi
>
> Is it ok, that VecView doesn't work properly with vectors, obtained with DACreateGlobalVector function, or im doing something wrong?
> In my program I create 5x5 da and global vector, corresponding to this da.
> Then each of 4 process writes its rank to its own local elements using DAVecGet(Restore)Array functions.
> Then i use VecView to visualize my global vector.
> And i see, that each process owns elements with different values.
> Where is the mistake?
>
VecView() for DM global vectors present the vector entries in the global natural ordering. Since PETSc divides the domain with DA into squareish blocks when you print the values on a process (by calling DAVecGetArray) it is not in the global natural ordering. Hence the two are different.
Barry
> Thank you very much!
>
> Alexey Ryazanov.
>
>
>
> CODE:
>
> #include "petscksp.h"
> #include "petscda.h"
> static char help[] = "VecVeiw doesn't work properly with DAGlobalVectors";
> int main(int argc,char **args)
> {
> Vec GlobalVec;
> DA da;
> PetscInt rank, dof, stencil_width, M, N, m, n, p, x, y, z, i, j;
> PetscScalar **array;
> PetscErrorCode ierr;
>
> PetscInitialize(&argc,&args,(char *)0,help);
> MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
>
> M=N=5; dof=1; stencil_width=1;
> ierr = DACreate2d(PETSC_COMM_WORLD, DA_NONPERIODIC, DA_STENCIL_BOX, M, N, PETSC_DECIDE, PETSC_DECIDE, dof, stencil_width, PETSC_NULL, PETSC_NULL, &da);
> CHKERRQ(ierr);
> ierr = DACreateGlobalVector(da,&GlobalVec); CHKERRQ(ierr);
>
> ierr = DAGetCorners(da,&x,&y,&z,&m,&n,&p); CHKERRQ(ierr);
>
> ierr = DAVecGetArray(da, GlobalVec, &array); CHKERRQ(ierr);
> for(i=y; i<y+n ; i++)
> for(j=x; j<x+m; j++)
> array[i][j]=rank;
> ierr = DAVecRestoreArray(da, GlobalVec, &array); CHKERRQ(ierr);
>
> ierr = VecView(GlobalVec, PETSC_VIEWER_STDOUT_WORLD); CHKERRQ(ierr);
>
> ierr = VecDestroy(GlobalVec); CHKERRQ(ierr);
> ierr = DADestroy(da); CHKERRQ(ierr);
> ierr = PetscFinalize(); CHKERRQ(ierr);
> return 0;
> }
>
>
> OUTPUT:
>
> Process [0]
> 0
> 0
> 0
> 1
> 1
> 0
> 0
> 0
> 1
> Process [1]
> 1
> 0
> 0
> 0
> 1
> 1
> Process [2]
> 2
> 2
> 2
> 3
> 3
> 2
> Process [3]
> 2
> 2
> 3
> 3
>
More information about the petsc-users
mailing list