Actual source code: ex26.c
slepc-3.8.0 2017-10-20
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2017, 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[] = "Computes the action of the square root of the 2-D Laplacian.\n\n"
12: "The command line options are:\n"
13: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
14: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";
16: #include <slepcmfn.h>
18: int main(int argc,char **argv)
19: {
20: Mat A; /* problem matrix */
21: MFN mfn;
22: FN f;
23: PetscReal norm,tol;
24: Vec v,y,z;
25: PetscInt N,n=10,m,Istart,Iend,i,j,II;
27: PetscBool flag,draw_sol;
29: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
31: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
32: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
33: if (!flag) m=n;
34: N = n*m;
35: PetscPrintf(PETSC_COMM_WORLD,"\nSquare root of Laplacian y=sqrt(A)*e_1, N=%D (%Dx%D grid)\n\n",N,n,m);
37: PetscOptionsHasName(NULL,NULL,"-draw_sol",&draw_sol);
39: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
40: Compute the discrete 2-D Laplacian, A
41: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
43: MatCreate(PETSC_COMM_WORLD,&A);
44: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
45: MatSetFromOptions(A);
46: MatSetUp(A);
48: MatGetOwnershipRange(A,&Istart,&Iend);
49: for (II=Istart;II<Iend;II++) {
50: i = II/n; j = II-i*n;
51: if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
52: if (i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
53: if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
54: if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
55: MatSetValue(A,II,II,4.0,INSERT_VALUES);
56: }
58: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
59: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
61: /* set symmetry flag so that solver can exploit it */
62: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
64: /* set v = e_1 */
65: MatCreateVecs(A,NULL,&v);
66: VecSetValue(v,0,1.0,INSERT_VALUES);
67: VecAssemblyBegin(v);
68: VecAssemblyEnd(v);
69: VecDuplicate(v,&y);
70: VecDuplicate(v,&z);
72: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73: Create the solver, set the matrix and the function
74: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75: MFNCreate(PETSC_COMM_WORLD,&mfn);
76: MFNSetOperator(mfn,A);
77: MFNGetFN(mfn,&f);
78: FNSetType(f,FNSQRT);
79: MFNSetErrorIfNotConverged(mfn,PETSC_TRUE);
80: MFNSetFromOptions(mfn);
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: First solve: y=sqrt(A)*v
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: MFNSolve(mfn,v,y);
87: VecNorm(y,NORM_2,&norm);
88: PetscPrintf(PETSC_COMM_WORLD," Intermediate vector has norm %g\n",(double)norm);
89: if (draw_sol) {
90: PetscViewerDrawSetPause(PETSC_VIEWER_DRAW_WORLD,-1);
91: VecView(y,PETSC_VIEWER_DRAW_WORLD);
92: }
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Second solve: z=sqrt(A)*y and compare against A*v
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: MFNSolve(mfn,y,z);
99: MFNGetTolerances(mfn,&tol,NULL);
101: MatMult(A,v,y); /* overwrite y */
102: VecAXPY(y,-1.0,z);
103: VecNorm(y,NORM_2,&norm);
105: if (norm<tol) {
106: PetscPrintf(PETSC_COMM_WORLD," Error norm is less than the requested tolerance\n\n");
107: } else {
108: PetscPrintf(PETSC_COMM_WORLD," Error norm larger than tolerance: %3.1e\n\n",(double)norm);
109: }
110: if (draw_sol) {
111: PetscViewerDrawSetPause(PETSC_VIEWER_DRAW_WORLD,-1);
112: VecView(z,PETSC_VIEWER_DRAW_WORLD);
113: }
115: /*
116: Free work space
117: */
118: MFNDestroy(&mfn);
119: MatDestroy(&A);
120: VecDestroy(&v);
121: VecDestroy(&y);
122: VecDestroy(&z);
123: SlepcFinalize();
124: return ierr;
125: }