Actual source code: test3.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: static char help[] = "Test MFN interface functions.\n\n"
12: "The command line options are:\n"
13: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n\n";
15: #include <slepcmfn.h>
17: int main(int argc,char **argv)
18: {
19: Mat A,B;
20: MFN mfn;
21: FN f;
22: MFNConvergedReason reason;
23: PetscReal norm,tol;
24: Vec v,y;
25: PetscInt N,n=4,Istart,Iend,i,j,II,ncv,its,maxit;
26: PetscBool flg;
27: PetscErrorCode ierr;
28: PetscViewerAndFormat *vf;
30: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
31: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
32: N = n*n;
33: PetscPrintf(PETSC_COMM_WORLD,"\nSquare root of Laplacian y=sqrt(A)*e_1, N=%D (%Dx%D grid)\n\n",N,n,n);
35: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36: Compute the discrete 2-D Laplacian, A
37: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
39: MatCreate(PETSC_COMM_WORLD,&A);
40: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
41: MatSetFromOptions(A);
42: MatSetUp(A);
44: MatGetOwnershipRange(A,&Istart,&Iend);
45: for (II=Istart;II<Iend;II++) {
46: i = II/n; j = II-i*n;
47: if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
48: if (i<n-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
49: if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
50: if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
51: MatSetValue(A,II,II,4.0,INSERT_VALUES);
52: }
54: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
55: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
56: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
58: MatCreateVecs(A,NULL,&v);
59: VecSetValue(v,0,1.0,INSERT_VALUES);
60: VecAssemblyBegin(v);
61: VecAssemblyEnd(v);
62: VecDuplicate(v,&y);
64: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
65: Create the solver, set the matrix and the function
66: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
67: MFNCreate(PETSC_COMM_WORLD,&mfn);
68: MFNSetOperator(mfn,A);
69: MFNGetFN(mfn,&f);
70: FNSetType(f,FNSQRT);
72: /* test some interface functions */
73: MFNGetOperator(mfn,&B);
74: MatView(B,PETSC_VIEWER_STDOUT_WORLD);
75: MFNSetOptionsPrefix(mfn,"myprefix_");
76: MFNSetTolerances(mfn,1e-4,500);
77: MFNSetDimensions(mfn,6);
78: MFNSetErrorIfNotConverged(mfn,PETSC_TRUE);
79: /* test monitors */
80: PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf);
81: MFNMonitorSet(mfn,(PetscErrorCode (*)(MFN,PetscInt,PetscReal,void*))MFNMonitorDefault,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
82: /* SVDMonitorCancel(svd); */
83: MFNSetFromOptions(mfn);
85: /* query properties and print them */
86: MFNGetTolerances(mfn,&tol,&maxit);
87: PetscPrintf(PETSC_COMM_WORLD," Tolerance: %g, max iterations: %D\n",(double)tol,maxit);
88: MFNGetDimensions(mfn,&ncv);
89: PetscPrintf(PETSC_COMM_WORLD," Subspace dimension: %D\n",ncv);
90: MFNGetErrorIfNotConverged(mfn,&flg);
91: if (flg) { PetscPrintf(PETSC_COMM_WORLD," Erroring out if convergence fails\n"); }
93: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94: Solve y=sqrt(A)*v
95: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97: MFNSolve(mfn,v,y);
98: MFNGetConvergedReason(mfn,&reason);
99: MFNGetIterationNumber(mfn,&its);
100: PetscPrintf(PETSC_COMM_WORLD," Finished - converged reason = %d\n",(int)reason);
101: /* PetscPrintf(PETSC_COMM_WORLD," its = %D\n",its); */
102: VecNorm(y,NORM_2,&norm);
103: PetscPrintf(PETSC_COMM_WORLD," sqrt(A)*v has norm %g\n",(double)norm);
105: /*
106: Free work space
107: */
108: MFNDestroy(&mfn);
109: MatDestroy(&A);
110: VecDestroy(&v);
111: VecDestroy(&y);
112: SlepcFinalize();
113: return ierr;
114: }
116: /*TEST
118: test:
119: suffix: 1
120: args: -myprefix_mfn_monitor_cancel -myprefix_mfn_converged_reason -myprefix_mfn_view
122: TEST*/