Actual source code: ex76f.F90
petsc-3.13.1 2020-05-02
1: !
2: ! Description: Solves a linear systems using PCHPDDM.
3: !
5: program main
6: #include <petsc/finclude/petscksp.h>
7: use petscksp
8: use petscisdef
9: implicit none
10: Vec x,b
11: Mat A,aux
12: KSP ksp
13: PC pc
14: IS is,sizes
15: PetscScalar one
16: PetscInt, pointer :: idx(:)
17: PetscMPIInt rank,size
18: PetscBool flg
19: character*(128) dir
20: character*(128) name
21: character*(8) fmt
22: character(1) crank,csize
23: PetscViewer viewer
24: PetscErrorCode ierr
26: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
27: if (ierr .ne. 0) then
28: print *,'Unable to initialize PETSc'
29: stop
30: endif
31: call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
32: if (size .ne. 4) then
33: print *,'This example requires 4 processes'
34: stop
35: endif
36: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
37: call MatCreate(PETSC_COMM_WORLD,A,ierr);CHKERRA(ierr)
38: call MatCreate(PETSC_COMM_SELF,aux,ierr);CHKERRA(ierr)
39: call ISCreate(PETSC_COMM_SELF,is,ierr);CHKERRA(ierr)
40: dir = '.'
41: call PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-load_dir',dir,flg,ierr);CHKERRA(ierr)
42: fmt = '(I1)'
43: write (crank,fmt) rank
44: write (csize,fmt) size
45: write (name,'(a)')trim(dir)//'/sizes_'//crank//'_'//csize//'.dat'
46: call PetscViewerBinaryOpen(PETSC_COMM_SELF,name,FILE_MODE_READ, viewer,ierr);CHKERRA(ierr)
47: call ISCreate(PETSC_COMM_SELF,sizes,ierr);CHKERRA(ierr)
48: call ISLoad(sizes,viewer,ierr);CHKERRA(ierr)
49: call ISGetIndicesF90(sizes,idx,ierr);CHKERRA(ierr)
50: call MatSetSizes(A,idx(1),idx(2),idx(3),idx(4),ierr);CHKERRA(ierr)
51: call ISRestoreIndicesF90(sizes,idx,ierr);CHKERRA(ierr)
52: call ISDestroy(sizes,ierr);CHKERRA(ierr)
53: call PetscViewerDestroy(viewer,ierr);CHKERRA(ierr)
54: call MatSetUp(A,ierr);CHKERRA(ierr)
55: write (name,'(a)')trim(dir)//'/A.dat'
56: call PetscViewerBinaryOpen(PETSC_COMM_WORLD,name,FILE_MODE_READ,viewer,ierr);CHKERRA(ierr)
57: call MatLoad(A,viewer,ierr);CHKERRA(ierr)
58: call PetscViewerDestroy(viewer,ierr);CHKERRA(ierr)
59: write (name,'(a)')trim(dir)//'/is_'//crank//'_'//csize//'.dat'
60: call PetscViewerBinaryOpen(PETSC_COMM_SELF,name,FILE_MODE_READ,viewer,ierr);CHKERRA(ierr)
61: call ISLoad(is,viewer,ierr);CHKERRA(ierr)
62: call PetscViewerDestroy(viewer,ierr);CHKERRA(ierr)
63: write (name,'(a)')trim(dir)//'/Neumann_'//crank//'_'//csize//'.dat'
64: call PetscViewerBinaryOpen(PETSC_COMM_SELF,name,FILE_MODE_READ,viewer,ierr);CHKERRA(ierr)
65: call MatSetBlockSizesFromMats(aux,A,A,ierr);CHKERRA(ierr)
66: call MatLoad(aux,viewer,ierr);CHKERRA(ierr)
67: call PetscViewerDestroy(viewer,ierr);CHKERRA(ierr)
68: call MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE,ierr);CHKERRA(ierr)
69: call MatSetOption(aux,MAT_SYMMETRIC,PETSC_TRUE,ierr);CHKERRA(ierr)
70: call KSPCreate(PETSC_COMM_WORLD,ksp,ierr);CHKERRA(ierr)
71: call KSPSetOperators(ksp,A,A,ierr);CHKERRA(ierr)
72: call KSPGetPC(ksp,pc,ierr);CHKERRA(ierr)
73: call PCSetType(pc,PCHPDDM,ierr);CHKERRA(ierr)
74: #if defined(PETSC_HAVE_HPDDM)
75: call PCHPDDMSetAuxiliaryMat(pc,is,aux,PETSC_NULL_FUNCTION,PETSC_NULL_INTEGER,ierr);CHKERRA(ierr)
76: call PCHPDDMHasNeumannMat(pc,PETSC_FALSE,ierr);CHKERRA(ierr)
77: #endif
78: call ISDestroy(is,ierr);CHKERRA(ierr)
79: call MatDestroy(aux,ierr);CHKERRA(ierr)
80: call KSPSetFromOptions(ksp,ierr);CHKERRA(ierr)
81: call MatCreateVecs(A,x,b,ierr);CHKERRA(ierr)
82: one = 1.0
83: call VecSet(b,one,ierr);CHKERRA(ierr);
84: call KSPSolve(ksp,b,x,ierr);CHKERRA(ierr)
85: call VecDestroy(x,ierr);CHKERRA(ierr)
86: call VecDestroy(b,ierr);CHKERRA(ierr)
87: call KSPDestroy(ksp,ierr);CHKERRA(ierr)
88: call MatDestroy(A,ierr);CHKERRA(ierr)
89: call PetscFinalize(ierr)
90: end
92: !/*TEST
93: !
94: ! test:
95: ! requires: hpddm slepc datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
96: ! output_file: output/ex76_1.out
97: ! nsize: 4
98: ! args: -ksp_rtol 1e-3 -ksp_converged_reason -pc_type {{bjacobi hpddm}shared output} -pc_hpddm_coarse_sub_pc_type lu -sub_pc_type lu -options_left no -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO
99: !
100: ! test:
101: ! requires: hpddm slepc datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
102: ! suffix: geneo
103: ! output_file: output/ex76_geneo_pc_hpddm_levels_1_eps_nev-5.out
104: ! nsize: 4
105: ! args: -ksp_converged_reason -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type cholesky -pc_hpddm_levels_1_eps_nev 5 -pc_hpddm_levels_1_st_pc_type cholesky -pc_hpddm_coarse_p {{1 2}shared output} -pc_hpddm_coarse_pc_type redundant -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO
106: !
107: ! test:
108: ! requires: hpddm slepc datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
109: ! suffix: fgmres_geneo_20_p_2
110: ! output_file: output/ex76_fgmres_geneo_20_p_2.out
111: ! nsize: 4
112: ! args: -ksp_converged_reason -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 20 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_type fgmres -pc_hpddm_coarse_mat_type {{baij sbaij}shared output} -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO
113: !
114: ! test:
115: ! requires: hpddm slepc datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
116: ! suffix: fgmres_geneo_20_p_2_geneo
117: ! output_file: output/ex76_fgmres_geneo_20_p_2.out
118: ! nsize: 4
119: ! args: -ksp_converged_reason -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type cholesky -pc_hpddm_levels_1_eps_nev 20 -pc_hpddm_levels_2_p 2 -pc_hpddm_levels_2_mat_type {{baij sbaij}shared output} -pc_hpddm_levels_2_eps_nev {{5 20}shared output} -pc_hpddm_levels_2_sub_pc_type cholesky -pc_hpddm_levels_2_ksp_type gmres -ksp_type fgmres -pc_hpddm_coarse_mat_type {{baij sbaij}shared output} -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO
120: !
121: !TEST*/