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