Actual source code: ex22f.F
1: !
2: ! Laplacian in 3D. Modeled by the partial differential equation
3: !
4: ! Laplacian u = 1,0 < x,y,z < 1,
5: !
6: ! with boundary conditions
7: !
8: ! u = 1 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1.
9: !
10: ! This uses multigrid to solve the linear system
12: program main
13: implicit none
15: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
16: ! Include files
17: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
18: !
19: ! petsc.h - base PETSc routines petscvec.h - vectors
20: ! petscsys.h - system routines petscmat.h - matrices
21: ! petscksp.h - Krylov subspace methods petscpc.h - preconditioners
23: #include include/finclude/petsc.h
24: #include include/finclude/petscvec.h
25: #include include/finclude/petscmat.h
26: #include include/finclude/petscpc.h
27: #include include/finclude/petscksp.h
28: #include include/finclude/petscda.h
29: DMMG dmmg
30: PetscErrorCode ierr
31: double precision norm
32: DA da
33: external ComputeRHS,ComputeJacobian
34: Vec x
35: PetscInt i1,i3
37: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
39: i3 = 3
40: i1 = 1
41: call DMMGCreate(PETSC_COMM_WORLD,i3,PETSC_NULL_INTEGER,dmmg,ierr)
42: call DACreate3d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR, &
43: & i3,i3,i3,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,i1,i1, &
44: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
45: & PETSC_NULL_INTEGER,da,ierr)
46: call DMMGSetDM(dmmg,da,ierr)
47: call DADestroy(da,ierr)
50: call DMMGSetKSP(dmmg,ComputeRHS,ComputeJacobian,ierr)
52: call DMMGSolve(dmmg,ierr)
54: call DMMGGetx(dmmg,x,ierr)
55: ! call VecView(x,PETSC_VIEWER_STDOUT_WORLD,ierr)
57: call DMMGDestroy(dmmg,ierr)
58: call PetscFinalize(ierr)
60: end
62: subroutine ComputeRHS(dmmg,b,ierr)
63: implicit none
65: #include include/finclude/petsc.h
66: #include include/finclude/petscvec.h
67: #include include/finclude/petscda.h
69: PetscErrorCode ierr
70: PetscInt mx,my,mz
71: PetscScalar h
72: Vec b
73: DMMG dmmg
74: DA da
76: call DMMGGetDA(dmmg,da,ierr)
77: call DAGetInfo(da,PETSC_NULL_INTEGER,mx,my,mz, &
78: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
79: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
80: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
81: & PETSC_NULL_INTEGER,ierr)
82: h = 1.d0/((mx-1)*(my-1)*(mz-1))
84: call VecSet(h,b,ierr)
85: return
86: end
88:
89: subroutine ComputeJacobian(dmmg,jac,ierr)
90: implicit none
92: #include include/finclude/petsc.h
93: #include include/finclude/petscvec.h
94: #include include/finclude/petscmat.h
95: #include include/finclude/petscda.h
97: DMMG dmmg
98: Mat jac
99: PetscErrorCode ierr
101: DA da
102: PetscInt i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs,i1,i7
103: PetscScalar v(7),Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy
104: MatStencil row(4),col(4,7)
106: i1 = 1
107: i7 = 7
108: call DMMGGetDA(dmmg,da,ierr)
109: call DAGetInfo(da,PETSC_NULL_INTEGER,mx,my,mz, &
110: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
111: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
112: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
113: & PETSC_NULL_INTEGER,ierr)
115: Hx = 1.d0 / (mx-1)
116: Hy = 1.d0 / (my-1)
117: Hz = 1.d0 / (mz-1)
118: HxHydHz = Hx*Hy/Hz
119: HxHzdHy = Hx*Hz/Hy
120: HyHzdHx = Hy*Hz/Hx
121: call DAGetCorners(da,xs,ys,zs,xm,ym,zm,ierr)
122:
123: do 10,k=zs,zs+zm-1
124: do 20,j=ys,ys+ym-1
125: do 30,i=xs,xs+xm-1
126: row(MatStencil_i) = i
127: row(MatStencil_j) = j
128: row(MatStencil_k) = k
129: if (i.eq.0 .or. j.eq.0 .or. k.eq.0 .or. i.eq.mx-1 .or. &
130: & j.eq.my-1 .or. k.eq.mz-1) then
131: v(1) = 2.d0*(HxHydHz + HxHzdHy + HyHzdHx)
132: call MatSetValuesStencil(jac,i1,row,i1,row,v,INSERT_VALUES, &
133: & ierr)
134: else
135: v(1) = -HxHydHz
136: col(MatStencil_i,1) = i
137: col(MatStencil_j,1) = j
138: col(MatStencil_k,1) = k-1
139: v(2) = -HxHzdHy
140: col(MatStencil_i,2) = i
141: col(MatStencil_j,2) = j-1
142: col(MatStencil_k,2) = k
143: v(3) = -HyHzdHx
144: col(MatStencil_i,3) = i-1
145: col(MatStencil_j,3) = j
146: col(MatStencil_k,3) = k
147: v(4) = 2.d0*(HxHydHz + HxHzdHy + HyHzdHx)
148: col(MatStencil_i,4) = i
149: col(MatStencil_j,4) = j
150: col(MatStencil_k,4) = k
151: v(5) = -HyHzdHx
152: col(MatStencil_i,5) = i+1
153: col(MatStencil_j,5) = j
154: col(MatStencil_k,5) = k
155: v(6) = -HxHzdHy
156: col(MatStencil_i,6) = i
157: col(MatStencil_j,6) = j+1
158: col(MatStencil_k,6) = k
159: v(7) = -HxHydHz
160: col(MatStencil_i,7) = i
161: col(MatStencil_j,7) = j
162: col(MatStencil_k,7) = k+1
163: call MatSetValuesStencil(jac,i1,row,i7,col,v,INSERT_VALUES, &
164: & ierr)
165: endif
166: 30 continue
167: 20 continue
168: 10 continue
170: call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
171: call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
172: return
173: end