1*c4762a1bSJed Brown program main 2*c4762a1bSJed Brown!----------------------------------------------------------------------- 3*c4762a1bSJed Brown! 4*c4762a1bSJed Brown! Tests DMDAGetVecGetArray() 5*c4762a1bSJed Brown!----------------------------------------------------------------------- 6*c4762a1bSJed Brown! 7*c4762a1bSJed Brown 8*c4762a1bSJed Brown#include <petsc/finclude/petscdm.h> 9*c4762a1bSJed Brown use petsc 10*c4762a1bSJed Brown implicit none 11*c4762a1bSJed Brown 12*c4762a1bSJed Brown Type(tVec) g 13*c4762a1bSJed Brown Type(tDM) ada 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown PetscScalar,pointer :: x1(:),x2(:,:) 16*c4762a1bSJed Brown PetscScalar,pointer :: x3(:,:,:),x4(:,:,:,:) 17*c4762a1bSJed Brown PetscErrorCode ierr 18*c4762a1bSJed Brown PetscInt m,n,p,dof,s,i,j,k,xs,xl 19*c4762a1bSJed Brown PetscInt ys,yl 20*c4762a1bSJed Brown PetscInt zs,zl,sw 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown m = 5 23*c4762a1bSJed Brown n = 6 24*c4762a1bSJed Brown p = 4; 25*c4762a1bSJed Brown s = 1 26*c4762a1bSJed Brown dof = 1 27*c4762a1bSJed Brown sw = 1 28*c4762a1bSJed Brown CALL PetscInitialize(PETSC_NULL_CHARACTER,ierr) 29*c4762a1bSJed Brown if (ierr .ne. 0) then 30*c4762a1bSJed Brown print*,'Unable to initialize PETSc' 31*c4762a1bSJed Brown stop 32*c4762a1bSJed Brown endif 33*c4762a1bSJed Brown call DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,m,dof,sw,PETSC_NULL_INTEGER,ada,ierr);CHKERRA(ierr) 34*c4762a1bSJed Brown call DMSetUp(ada,ierr);CHKERRA(ierr) 35*c4762a1bSJed Brown call DMGetGlobalVector(ada,g,ierr);CHKERRA(ierr) 36*c4762a1bSJed Brown call DMDAGetCorners(ada,xs,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,xl,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr);CHKERRA(ierr) 37*c4762a1bSJed Brown call DMDAVecGetArrayF90(ada,g,x1,ierr);CHKERRA(ierr) 38*c4762a1bSJed Brown do i=xs,xs+xl-1 39*c4762a1bSJed Brown! CHKMEMQ 40*c4762a1bSJed Brown x1(i) = i 41*c4762a1bSJed Brown! CHKMEMQ 42*c4762a1bSJed Brown enddo 43*c4762a1bSJed Brown call DMDAVecRestoreArrayF90(ada,g,x1,ierr);CHKERRA(ierr) 44*c4762a1bSJed Brown call VecView(g,PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr) 45*c4762a1bSJed Brown call DMRestoreGlobalVector(ada,g,ierr);CHKERRA(ierr) 46*c4762a1bSJed Brown call DMDestroy(ada,ierr);CHKERRA(ierr) 47*c4762a1bSJed Brown 48*c4762a1bSJed Brown call DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,m,n,PETSC_DECIDE,PETSC_DECIDE,dof,s, & 49*c4762a1bSJed Brown & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ada,ierr);CHKERRA(ierr) 50*c4762a1bSJed Brown call DMSetUp(ada,ierr);CHKERRA(ierr) 51*c4762a1bSJed Brown call DMGetGlobalVector(ada,g,ierr);CHKERRA(ierr) 52*c4762a1bSJed Brown call DMDAGetCorners(ada,xs,ys,PETSC_NULL_INTEGER,xl,yl,PETSC_NULL_INTEGER,ierr);CHKERRA(ierr) 53*c4762a1bSJed Brown call DMDAVecGetArrayF90(ada,g,x2,ierr);CHKERRA(ierr) 54*c4762a1bSJed Brown do i=xs,xs+xl-1 55*c4762a1bSJed Brown do j=ys,ys+yl-1 56*c4762a1bSJed Brown! CHKMEMQ 57*c4762a1bSJed Brown x2(i,j) = i + j 58*c4762a1bSJed Brown! CHKMEMQ 59*c4762a1bSJed Brown enddo 60*c4762a1bSJed Brown enddo 61*c4762a1bSJed Brown call DMDAVecRestoreArrayF90(ada,g,x2,ierr);CHKERRA(ierr) 62*c4762a1bSJed Brown call VecView(g,PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr) 63*c4762a1bSJed Brown call DMRestoreGlobalVector(ada,g,ierr);CHKERRA(ierr) 64*c4762a1bSJed Brown call DMDestroy(ada,ierr);CHKERRA(ierr) 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown call DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX, m,n,p,PETSC_DECIDE,PETSC_DECIDE, & 67*c4762a1bSJed Brown & PETSC_DECIDE,dof,s,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ada,ierr);CHKERRA(ierr) 68*c4762a1bSJed Brown call DMSetUp(ada,ierr);CHKERRA(ierr) 69*c4762a1bSJed Brown call DMGetGlobalVector(ada,g,ierr);CHKERRA(ierr) 70*c4762a1bSJed Brown call DMDAGetCorners(ada,xs,ys,zs,xl,yl,zl,ierr);CHKERRA(ierr) 71*c4762a1bSJed Brown call DMDAVecGetArrayF90(ada,g,x3,ierr);CHKERRA(ierr) 72*c4762a1bSJed Brown do i=xs,xs+xl-1 73*c4762a1bSJed Brown do j=ys,ys+yl-1 74*c4762a1bSJed Brown do k=zs,zs+zl-1 75*c4762a1bSJed Brown! CHKMEMQ 76*c4762a1bSJed Brown x3(i,j,k) = i + j + k 77*c4762a1bSJed Brown! CHKMEMQ 78*c4762a1bSJed Brown enddo 79*c4762a1bSJed Brown enddo 80*c4762a1bSJed Brown enddo 81*c4762a1bSJed Brown call DMDAVecRestoreArrayF90(ada,g,x3,ierr);CHKERRA(ierr) 82*c4762a1bSJed Brown call VecView(g,PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr) 83*c4762a1bSJed Brown call DMRestoreGlobalVector(ada,g,ierr);CHKERRA(ierr) 84*c4762a1bSJed Brown call DMDestroy(ada,ierr);CHKERRA(ierr) 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown! 87*c4762a1bSJed Brown! Same tests but now with DOF > 1, so dimensions of array are one higher 88*c4762a1bSJed Brown! 89*c4762a1bSJed Brown dof = 2 90*c4762a1bSJed Brown call DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,m,dof,sw,PETSC_NULL_INTEGER,ada,ierr);CHKERRA(ierr) 91*c4762a1bSJed Brown call DMSetUp(ada,ierr);CHKERRA(ierr) 92*c4762a1bSJed Brown call DMGetGlobalVector(ada,g,ierr);CHKERRA(ierr) 93*c4762a1bSJed Brown call DMDAGetCorners(ada,xs,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,xl,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr);CHKERRA(ierr) 94*c4762a1bSJed Brown call DMDAVecGetArrayF90(ada,g,x2,ierr);CHKERRA(ierr) 95*c4762a1bSJed Brown do i=xs,xs+xl-1 96*c4762a1bSJed Brown! CHKMEMQ 97*c4762a1bSJed Brown x2(0,i) = i 98*c4762a1bSJed Brown x2(1,i) = -i 99*c4762a1bSJed Brown! CHKMEMQ 100*c4762a1bSJed Brown enddo 101*c4762a1bSJed Brown call DMDAVecRestoreArrayF90(ada,g,x1,ierr);CHKERRA(ierr) 102*c4762a1bSJed Brown call VecView(g,PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr) 103*c4762a1bSJed Brown call DMRestoreGlobalVector(ada,g,ierr);CHKERRA(ierr) 104*c4762a1bSJed Brown call DMDestroy(ada,ierr);CHKERRA(ierr) 105*c4762a1bSJed Brown 106*c4762a1bSJed Brown dof = 2 107*c4762a1bSJed Brown call DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,m,n,PETSC_DECIDE,PETSC_DECIDE,dof,s, & 108*c4762a1bSJed Brown & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ada,ierr);CHKERRA(ierr) 109*c4762a1bSJed Brown call DMSetUp(ada,ierr);CHKERRA(ierr) 110*c4762a1bSJed Brown call DMGetGlobalVector(ada,g,ierr);CHKERRA(ierr) 111*c4762a1bSJed Brown call DMDAGetCorners(ada,xs,ys,PETSC_NULL_INTEGER,xl,yl,PETSC_NULL_INTEGER,ierr);CHKERRA(ierr) 112*c4762a1bSJed Brown call DMDAVecGetArrayF90(ada,g,x3,ierr);CHKERRA(ierr) 113*c4762a1bSJed Brown do i=xs,xs+xl-1 114*c4762a1bSJed Brown do j=ys,ys+yl-1 115*c4762a1bSJed Brown! CHKMEMQ 116*c4762a1bSJed Brown x3(0,i,j) = i + j 117*c4762a1bSJed Brown x3(1,i,j) = -(i + j) 118*c4762a1bSJed Brown! CHKMEMQ 119*c4762a1bSJed Brown enddo 120*c4762a1bSJed Brown enddo 121*c4762a1bSJed Brown call DMDAVecRestoreArrayF90(ada,g,x3,ierr);CHKERRA(ierr) 122*c4762a1bSJed Brown call VecView(g,PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr) 123*c4762a1bSJed Brown call DMRestoreGlobalVector(ada,g,ierr);CHKERRA(ierr) 124*c4762a1bSJed Brown call DMDestroy(ada,ierr);CHKERRA(ierr) 125*c4762a1bSJed Brown 126*c4762a1bSJed Brown dof = 3 127*c4762a1bSJed Brown call DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,m,n,p,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,dof,s, & 128*c4762a1bSJed Brown & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ada,ierr);CHKERRA(ierr) 129*c4762a1bSJed Brown call DMSetUp(ada,ierr);CHKERRA(ierr) 130*c4762a1bSJed Brown call DMGetGlobalVector(ada,g,ierr);CHKERRA(ierr) 131*c4762a1bSJed Brown call DMDAGetCorners(ada,xs,ys,zs,xl,yl,zl,ierr);CHKERRA(ierr) 132*c4762a1bSJed Brown call DMDAVecGetArrayF90(ada,g,x4,ierr);CHKERRA(ierr) 133*c4762a1bSJed Brown do i=xs,xs+xl-1 134*c4762a1bSJed Brown do j=ys,ys+yl-1 135*c4762a1bSJed Brown do k=zs,zs+zl-1 136*c4762a1bSJed Brown! CHKMEMQ 137*c4762a1bSJed Brown x4(0,i,j,k) = i + j + k 138*c4762a1bSJed Brown x4(1,i,j,k) = -(i + j + k) 139*c4762a1bSJed Brown x4(2,i,j,k) = i + j + k 140*c4762a1bSJed Brown! CHKMEMQ 141*c4762a1bSJed Brown enddo 142*c4762a1bSJed Brown enddo 143*c4762a1bSJed Brown enddo 144*c4762a1bSJed Brown call DMDAVecRestoreArrayF90(ada,g,x4,ierr);CHKERRA(ierr) 145*c4762a1bSJed Brown call VecView(g,PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr) 146*c4762a1bSJed Brown call DMRestoreGlobalVector(ada,g,ierr);CHKERRA(ierr) 147*c4762a1bSJed Brown call DMDestroy(ada,ierr);CHKERRA(ierr) 148*c4762a1bSJed Brown 149*c4762a1bSJed Brown CALL PetscFinalize(ierr) 150*c4762a1bSJed Brown END PROGRAM 151*c4762a1bSJed Brown 152*c4762a1bSJed Brown! 153*c4762a1bSJed Brown!/*TEST 154*c4762a1bSJed Brown! 155*c4762a1bSJed Brown! build: 156*c4762a1bSJed Brown! requires: !complex 157*c4762a1bSJed Brown! 158*c4762a1bSJed Brown! test: 159*c4762a1bSJed Brown! filter: Error: grep -v "Vec Object" | grep -v "Warning: ieee_inexact is signaling" 160*c4762a1bSJed Brown! 161*c4762a1bSJed Brown!TEST*/ 162