Actual source code: ex15f.F
petsc-3.7.7 2017-09-25
1: !
2: ! Solves a linear system in parallel with KSP. Also indicates
3: ! use of a user-provided preconditioner. Input parameters include:
4: ! -user_defined_pc : Activate a user-defined preconditioner
5: !
6: !
7: !/*T
8: ! Concepts: KSP^basic parallel example
9: ! Concepts: PC^setting a user-defined shell preconditioner
10: ! Processors: n
11: !T*/
12: !
13: ! -------------------------------------------------------------------------
15: program main
16: implicit none
18: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
19: ! Include files
20: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
21: !
22: ! petscsys.h - base PETSc routines petscvec.h - vectors
23: ! petscmat.h - matrices
24: ! petscksp.h - Krylov subspace methods petscpc.h - preconditioners
26: #include <petsc/finclude/petscsys.h>
27: #include <petsc/finclude/petscvec.h>
28: #include <petsc/finclude/petscmat.h>
29: #include <petsc/finclude/petscpc.h>
30: #include <petsc/finclude/petscksp.h>
32: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
33: ! Variable declarations
34: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35: !
36: ! Variables:
37: ! ksp - linear solver context
38: ! ksp - Krylov subspace method context
39: ! pc - preconditioner context
40: ! x, b, u - approx solution, right-hand-side, exact solution vectors
41: ! A - matrix that defines linear system
42: ! its - iterations for convergence
43: ! norm - norm of solution error
45: Vec x,b,u
46: Mat A
47: PC pc
48: KSP ksp
49: PetscScalar v,one,neg_one
50: PetscReal norm,tol
51: PetscErrorCode ierr
52: PetscInt i,j,II,JJ,Istart
53: PetscInt Iend,m,n,i1,its,five
54: PetscMPIInt rank
55: PetscBool user_defined_pc,flg
57: ! Note: Any user-defined Fortran routines MUST be declared as external.
59: external SampleShellPCSetUp, SampleShellPCApply
60: external SampleShellPCDestroy
62: ! Common block to store data for user-provided preconditioner
63: common /myshellpc/ diag
64: Vec diag
66: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67: ! Beginning of program
68: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
70: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
71: one = 1.0
72: neg_one = -1.0
73: i1 = 1
74: m = 8
75: n = 7
76: five = 5
77: call PetscOptionsGetInt(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER, &
78: & '-m',m,flg,ierr)
79: call PetscOptionsGetInt(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER, &
80: & '-n',n,flg,ierr)
81: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
83: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: ! Compute the matrix and right-hand-side vector that define
85: ! the linear system, Ax = b.
86: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88: ! Create parallel matrix, specifying only its global dimensions.
89: ! When using MatCreate(), the matrix format can be specified at
90: ! runtime. Also, the parallel partitioning of the matrix is
91: ! determined by PETSc at runtime.
93: call MatCreate(PETSC_COMM_WORLD,A,ierr)
94: call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,ierr)
95: call MatSetType(A, MATAIJ,ierr)
96: call MatSetFromOptions(A,ierr)
97: call MatMPIAIJSetPreallocation(A,five,PETSC_NULL_INTEGER,five, &
98: & PETSC_NULL_INTEGER,ierr)
99: call MatSeqAIJSetPreallocation(A,five,PETSC_NULL_INTEGER,ierr)
101: ! Currently, all PETSc parallel matrix formats are partitioned by
102: ! contiguous chunks of rows across the processors. Determine which
103: ! rows of the matrix are locally owned.
105: call MatGetOwnershipRange(A,Istart,Iend,ierr)
107: ! Set matrix elements for the 2-D, five-point stencil in parallel.
108: ! - Each processor needs to insert only elements that it owns
109: ! locally (but any non-local elements will be sent to the
110: ! appropriate processor during matrix assembly).
111: ! - Always specify global row and columns of matrix entries.
112: ! - Note that MatSetValues() uses 0-based row and column numbers
113: ! in Fortran as well as in C.
115: do 10, II=Istart,Iend-1
116: v = -1.0
117: i = II/n
118: j = II - i*n
119: if (i.gt.0) then
120: JJ = II - n
121: call MatSetValues(A,i1,II,i1,JJ,v,ADD_VALUES,ierr)
122: endif
123: if (i.lt.m-1) then
124: JJ = II + n
125: call MatSetValues(A,i1,II,i1,JJ,v,ADD_VALUES,ierr)
126: endif
127: if (j.gt.0) then
128: JJ = II - 1
129: call MatSetValues(A,i1,II,i1,JJ,v,ADD_VALUES,ierr)
130: endif
131: if (j.lt.n-1) then
132: JJ = II + 1
133: call MatSetValues(A,i1,II,i1,JJ,v,ADD_VALUES,ierr)
134: endif
135: v = 4.0
136: call MatSetValues(A,i1,II,i1,II,v,ADD_VALUES,ierr)
137: 10 continue
139: ! Assemble matrix, using the 2-step process:
140: ! MatAssemblyBegin(), MatAssemblyEnd()
141: ! Computations can be done while messages are in transition,
142: ! by placing code between these two statements.
144: call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
145: call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
147: ! Create parallel vectors.
148: ! - Here, the parallel partitioning of the vector is determined by
149: ! PETSc at runtime. We could also specify the local dimensions
150: ! if desired -- or use the more general routine VecCreate().
151: ! - When solving a linear system, the vectors and matrices MUST
152: ! be partitioned accordingly. PETSc automatically generates
153: ! appropriately partitioned matrices and vectors when MatCreate()
154: ! and VecCreate() are used with the same communicator.
155: ! - Note: We form 1 vector from scratch and then duplicate as needed.
157: call VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,m*n,u,ierr)
158: call VecDuplicate(u,b,ierr)
159: call VecDuplicate(b,x,ierr)
161: ! Set exact solution; then compute right-hand-side vector.
163: call VecSet(u,one,ierr)
164: call MatMult(A,u,b,ierr)
166: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167: ! Create the linear solver and set various options
168: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
170: ! Create linear solver context
172: call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)
174: ! Set operators. Here the matrix that defines the linear system
175: ! also serves as the preconditioning matrix.
177: call KSPSetOperators(ksp,A,A,ierr)
179: ! Set linear solver defaults for this problem (optional).
180: ! - By extracting the KSP and PC contexts from the KSP context,
181: ! we can then directly directly call any KSP and PC routines
182: ! to set various options.
184: call KSPGetPC(ksp,pc,ierr)
185: tol = 1.e-7
186: call KSPSetTolerances(ksp,tol,PETSC_DEFAULT_REAL, &
187: & PETSC_DEFAULT_REAL,PETSC_DEFAULT_INTEGER,ierr)
189: !
190: ! Set a user-defined shell preconditioner if desired
191: !
192: call PetscOptionsHasName(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER, &
193: & '-user_defined_pc',user_defined_pc,ierr)
195: if (user_defined_pc) then
197: ! (Required) Indicate to PETSc that we are using a shell preconditioner
198: call PCSetType(pc,PCSHELL,ierr)
200: ! (Required) Set the user-defined routine for applying the preconditioner
201: call PCShellSetApply(pc,SampleShellPCApply,ierr)
203: ! (Optional) Do any setup required for the preconditioner
204: call PCShellSetSetUp(pc,SampleShellPCSetUp,ierr)
206: ! (Optional) Frees any objects we created for the preconditioner
207: call PCShellSetDestroy(pc,SampleShellPCDestroy,ierr)
209: else
210: call PCSetType(pc,PCJACOBI,ierr)
211: endif
213: ! Set runtime options, e.g.,
214: ! -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
215: ! These options will override those specified above as long as
216: ! KSPSetFromOptions() is called _after_ any other customization
217: ! routines.
219: call KSPSetFromOptions(ksp,ierr)
221: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222: ! Solve the linear system
223: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225: call KSPSolve(ksp,b,x,ierr)
227: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
228: ! Check solution and clean up
229: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
231: ! Check the error
233: call VecAXPY(x,neg_one,u,ierr)
234: call VecNorm(x,NORM_2,norm,ierr)
235: call KSPGetIterationNumber(ksp,its,ierr)
237: if (rank .eq. 0) then
238: if (norm .gt. 1.e-12) then
239: write(6,100) norm,its
240: else
241: write(6,110) its
242: endif
243: endif
244: 100 format('Norm of error ',1pe11.4,' iterations ',i5)
245: 110 format('Norm of error < 1.e-12,iterations ',i5)
247: ! Free work space. All PETSc objects should be destroyed when they
248: ! are no longer needed.
250: call KSPDestroy(ksp,ierr)
251: call VecDestroy(u,ierr)
252: call VecDestroy(x,ierr)
253: call VecDestroy(b,ierr)
254: call MatDestroy(A,ierr)
256: ! Always call PetscFinalize() before exiting a program.
258: call PetscFinalize(ierr)
259: end
261: !/***********************************************************************/
262: !/* Routines for a user-defined shell preconditioner */
263: !/***********************************************************************/
265: !
266: ! SampleShellPCSetUp - This routine sets up a user-defined
267: ! preconditioner context.
268: !
269: ! Input Parameters:
270: ! pc - preconditioner object
271: !
272: ! Output Parameter:
273: ! ierr - error code (nonzero if error has been detected)
274: !
275: ! Notes:
276: ! In this example, we define the shell preconditioner to be Jacobi
277: ! method. Thus, here we create a work vector for storing the reciprocal
278: ! of the diagonal of the preconditioner matrix; this vector is then
279: ! used within the routine SampleShellPCApply().
280: !
281: subroutine SampleShellPCSetUp(pc,ierr)
283: implicit none
285: #include <petsc/finclude/petscsys.h>
286: #include <petsc/finclude/petscvec.h>
287: #include <petsc/finclude/petscmat.h>
288: PC pc
290: Mat pmat
291: integer ierr
293: ! Common block to store data for user-provided preconditioner
294: ! Normally we would recommend storing all the work data (like diag) in
295: ! the context set with PCShellSetContext()
297: common /myshellpc/ diag
298: Vec diag
300: call PCGetOperators(pc,PETSC_NULL_OBJECT,pmat,ierr)
301: call MatCreateVecs(pmat,diag,PETSC_NULL_OBJECT,ierr)
302: call MatGetDiagonal(pmat,diag,ierr)
303: call VecReciprocal(diag,ierr)
305: end
307: ! -------------------------------------------------------------------
308: !
309: ! SampleShellPCApply - This routine demonstrates the use of a
310: ! user-provided preconditioner.
311: !
312: ! Input Parameters:
313: ! pc - preconditioner object
314: ! x - input vector
315: !
316: ! Output Parameters:
317: ! y - preconditioned vector
318: ! ierr - error code (nonzero if error has been detected)
319: !
320: ! Notes:
321: ! This code implements the Jacobi preconditioner, merely as an
322: ! example of working with a PCSHELL. Note that the Jacobi method
323: ! is already provided within PETSc.
324: !
325: subroutine SampleShellPCApply(pc,x,y,ierr)
327: implicit none
329: #include <petsc/finclude/petscsys.h>
330: #include <petsc/finclude/petscvec.h>
332: PC pc
333: Vec x,y
334: integer ierr
336: ! Common block to store data for user-provided preconditioner
337: common /myshellpc/ diag
338: Vec diag
340: call VecPointwiseMult(y,x,diag,ierr)
342: end
344: !/***********************************************************************/
345: !/* Routines for a user-defined shell preconditioner */
346: !/***********************************************************************/
348: !
349: ! SampleShellPCDestroy - This routine destroys (frees the memory of) any
350: ! objects we made for the preconditioner
351: !
352: ! Input Parameters:
353: ! pc - for this example we use the actual PC as our shell context
354: !
355: ! Output Parameter:
356: ! ierr - error code (nonzero if error has been detected)
357: !
359: subroutine SampleShellPCDestroy(pc,ierr)
361: implicit none
363: #include <petsc/finclude/petscsys.h>
364: #include <petsc/finclude/petscvec.h>
365: #include <petsc/finclude/petscmat.h>
366: PC pc
367: PetscErrorCode ierr
369: ! Common block to store data for user-provided preconditioner
370: ! Normally we would recommend storing all the work data (like diag) in
371: ! the context set with PCShellSetContext()
373: common /myshellpc/ diag
374: Vec diag
376: call VecDestroy(diag,ierr)
378: end