Actual source code: ex76f.F90
1: !
2: ! Description: Solves a linear systems using PCHPDDM.
3: !
5: program main
6: #include <petsc/finclude/petscksp.h>
7: use petscis
8: use petscvec
9: use petscmat
10: use petscksp
11: implicit none
12: Vec x,b
13: Mat A,aux,Y,C
14: KSP ksp
15: PC pc
16: IS is,sizes
17: PetscScalar one
18: PetscInt, pointer :: idx(:)
19: PetscMPIInt rank,size
20: PetscInt m,N
21: PetscViewer viewer
22: character*(PETSC_MAX_PATH_LEN) dir,name
23: character*(8) fmt
24: character(1) crank,csize
25: PetscBool flg
26: PetscErrorCode ierr
28: PetscCallA(PetscInitialize(ierr))
30: PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD,size,ierr))
31: N = 1
32: PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-rhs',N,flg,ierr))
33: PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
34: PetscCallA(MatCreate(PETSC_COMM_WORLD,A,ierr))
35: PetscCallA(MatCreate(PETSC_COMM_SELF,aux,ierr))
36: PetscCallA(ISCreate(PETSC_COMM_SELF,is,ierr))
37: dir = '.'
38: PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-load_dir',dir,flg,ierr))
39: fmt = '(I1)'
40: write (crank,fmt) rank
41: write (csize,fmt) size
42: write (name,'(a)')trim(dir)//'/sizes_'//crank//'_'//csize//'.dat'
43: PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_SELF,name,FILE_MODE_READ, viewer,ierr))
44: PetscCallA(ISCreate(PETSC_COMM_SELF,sizes,ierr))
45: PetscCallA(ISLoad(sizes,viewer,ierr))
46: PetscCallA(ISGetIndices(sizes,idx,ierr))
47: PetscCallA(MatSetSizes(A,idx(1),idx(2),idx(3),idx(4),ierr))
48: PetscCallA(ISRestoreIndices(sizes,idx,ierr))
49: PetscCallA(ISDestroy(sizes,ierr))
50: PetscCallA(PetscViewerDestroy(viewer,ierr))
51: write (name,'(a)')trim(dir)//'/A.dat'
52: PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_WORLD,name,FILE_MODE_READ,viewer,ierr))
53: PetscCallA(MatLoad(A,viewer,ierr))
54: PetscCallA(PetscViewerDestroy(viewer,ierr))
55: write (name,'(a)')trim(dir)//'/is_'//crank//'_'//csize//'.dat'
56: PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_SELF,name,FILE_MODE_READ,viewer,ierr))
57: PetscCallA(ISLoad(is,viewer,ierr))
58: PetscCallA(PetscViewerDestroy(viewer,ierr))
59: write (name,'(a)')trim(dir)//'/Neumann_'//crank//'_'//csize//'.dat'
60: PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_SELF,name,FILE_MODE_READ,viewer,ierr))
61: PetscCallA(MatSetBlockSizesFromMats(aux,A,A,ierr))
62: PetscCallA(MatLoad(aux,viewer,ierr))
63: PetscCallA(PetscViewerDestroy(viewer,ierr))
64: PetscCallA(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE,ierr))
65: PetscCallA(MatSetOption(aux,MAT_SYMMETRIC,PETSC_TRUE,ierr))
66: PetscCallA(KSPCreate(PETSC_COMM_WORLD,ksp,ierr))
67: PetscCallA(KSPSetOperators(ksp,A,A,ierr))
68: PetscCallA(KSPGetPC(ksp,pc,ierr))
69: PetscCallA(PCSetType(pc,PCHPDDM,ierr))
70: #if defined(PETSC_HAVE_HPDDM) && defined(PETSC_HAVE_DYNAMIC_LIBRARIES) && defined(PETSC_USE_SHARED_LIBRARIES)
71: PetscCallA(PCHPDDMSetAuxiliaryMat(pc,is,aux,PETSC_NULL_FUNCTION,PETSC_NULL_INTEGER,ierr))
72: PetscCallA(PCHPDDMHasNeumannMat(pc,PETSC_FALSE,ierr))
73: #endif
74: PetscCallA(ISDestroy(is,ierr))
75: PetscCallA(MatDestroy(aux,ierr))
76: PetscCallA(KSPSetFromOptions(ksp,ierr))
77: PetscCallA(MatCreateVecs(A,x,b,ierr))
78: one = 1.0
79: PetscCallA(VecSet(b,one,ierr))
80: PetscCallA(KSPSolve(ksp,b,x,ierr))
81: PetscCallA(VecGetLocalSize(x,m,ierr))
82: PetscCallA(VecDestroy(x,ierr))
83: PetscCallA(VecDestroy(b,ierr))
84: if (N .gt. 1) then
85: PetscCallA(PetscOptionsClearValue(PETSC_NULL_OPTIONS,'-ksp_converged_reason',ierr))
86: PetscCallA(KSPSetFromOptions(ksp,ierr))
87: PetscCallA(MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,PETSC_DECIDE,N,PETSC_NULL_SCALAR_ARRAY,C,ierr))
88: PetscCallA(MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,PETSC_DECIDE,N,PETSC_NULL_SCALAR_ARRAY,Y,ierr))
89: PetscCallA(MatSetRandom(C,PETSC_NULL_RANDOM,ierr))
90: PetscCallA(KSPMatSolve(ksp,C,Y,ierr))
91: PetscCallA(MatDestroy(Y,ierr))
92: PetscCallA(MatDestroy(C,ierr))
93: endif
94: PetscCallA(KSPDestroy(ksp,ierr))
95: PetscCallA(MatDestroy(A,ierr))
96: PetscCallA(PetscFinalize(ierr))
97: end
99: !/*TEST
100: !
101: ! test:
102: ! requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
103: ! output_file: output/ex76_1.out
104: ! nsize: 4
105: ! 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
106: !
107: ! test:
108: ! requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
109: ! suffix: geneo
110: ! output_file: output/ex76_geneo_pc_hpddm_levels_1_eps_nev-5.out
111: ! nsize: 4
112: ! 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
113: !
114: ! test:
115: ! requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
116: ! suffix: fgmres_geneo_20_p_2
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 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
120: !
121: ! test:
122: ! requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
123: ! suffix: fgmres_geneo_20_p_2_geneo
124: ! output_file: output/ex76_fgmres_geneo_20_p_2.out
125: ! nsize: 4
126: ! 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
127: ! # PCHPDDM + KSPHPDDM test to exercise multilevel + multiple RHS in one go
128: ! test:
129: ! requires: hpddm slepc datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
130: ! suffix: fgmres_geneo_20_p_2_geneo_rhs
131: ! output_file: output/ex76_fgmres_geneo_20_p_2.out
132: ! nsize: 4
133: ! 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 -pc_hpddm_levels_2_eps_nev 5 -pc_hpddm_levels_2_sub_pc_type cholesky -pc_hpddm_levels_2_ksp_max_it 10 -pc_hpddm_levels_2_ksp_type hpddm -pc_hpddm_levels_2_ksp_hpddm_type gmres -ksp_type hpddm -ksp_hpddm_variant flexible -pc_hpddm_coarse_mat_type baij -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO -rhs 4
134: !
135: !TEST*/