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