Actual source code: ex26.c
slepc-3.7.4 2017-05-17
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: static char help[] = "Computes the action of the square root of the 2-D Laplacian.\n\n"
23: "The command line options are:\n"
24: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
25: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";
27: #include <slepcmfn.h>
31: int main(int argc,char **argv)
32: {
33: Mat A; /* problem matrix */
34: MFN mfn;
35: FN f;
36: PetscReal norm;
37: Vec v,y,z;
38: PetscInt N,n=10,m,Istart,Iend,i,j,II;
39: PetscErrorCode ierr;
40: PetscBool flag,draw_sol;
41: MFNConvergedReason reason;
43: SlepcInitialize(&argc,&argv,(char*)0,help);
45: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
46: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
47: if (!flag) m=n;
48: N = n*m;
49: PetscPrintf(PETSC_COMM_WORLD,"\nSquare root of Laplacian y=sqrt(A)*e_1, N=%D (%Dx%D grid)\n\n",N,n,m);
51: PetscOptionsHasName(NULL,NULL,"-draw_sol",&draw_sol);
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Compute the discrete 2-D Laplacian, A
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57: MatCreate(PETSC_COMM_WORLD,&A);
58: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
59: MatSetFromOptions(A);
60: MatSetUp(A);
62: MatGetOwnershipRange(A,&Istart,&Iend);
63: for (II=Istart;II<Iend;II++) {
64: i = II/n; j = II-i*n;
65: if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
66: if (i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
67: if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
68: if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
69: MatSetValue(A,II,II,4.0,INSERT_VALUES);
70: }
72: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
73: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
75: /* set symmetry flag so that solver can exploit it */
76: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
78: /* set v = e_1 */
79: MatCreateVecs(A,NULL,&v);
80: VecSetValue(v,0,1.0,INSERT_VALUES);
81: VecAssemblyBegin(v);
82: VecAssemblyEnd(v);
83: VecDuplicate(v,&y);
84: VecDuplicate(v,&z);
86: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: Create the solver, set the matrix and the function
88: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
89: MFNCreate(PETSC_COMM_WORLD,&mfn);
90: MFNSetOperator(mfn,A);
91: MFNGetFN(mfn,&f);
92: FNSetType(f,FNSQRT);
93: MFNSetFromOptions(mfn);
95: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96: First solve: y=sqrt(A)*v
97: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99: MFNSolve(mfn,v,y);
100: MFNGetConvergedReason(mfn,&reason);
101: if (reason<0) SETERRQ(PETSC_COMM_WORLD,1,"Solver did not converge");
102: VecNorm(y,NORM_2,&norm);
103:
104: PetscPrintf(PETSC_COMM_WORLD," Intermediate vector has norm %g\n",(double)norm);
105: if (draw_sol) {
106: PetscViewerDrawSetPause(PETSC_VIEWER_DRAW_WORLD,-1);
107: VecView(y,PETSC_VIEWER_DRAW_WORLD);
108: }
110: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111: Second solve: z=sqrt(A)*y and compare against A*v
112: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114: MFNSolve(mfn,y,z);
115: MFNGetConvergedReason(mfn,&reason);
116: if (reason<0) SETERRQ(PETSC_COMM_WORLD,1,"Solver did not converge");
118: MatMult(A,v,y); /* overwrite y */
119: VecAXPY(y,-1.0,z);
120: VecNorm(y,NORM_2,&norm);
121:
122: if (norm<100*PETSC_MACHINE_EPSILON) {
123: PetscPrintf(PETSC_COMM_WORLD," Error norm is less than 100*epsilon\n\n");
124: } else {
125: PetscPrintf(PETSC_COMM_WORLD," Error norm %3.1e\n\n",(double)norm);
126: }
127: if (draw_sol) {
128: PetscViewerDrawSetPause(PETSC_VIEWER_DRAW_WORLD,-1);
129: VecView(z,PETSC_VIEWER_DRAW_WORLD);
130: }
132: /*
133: Free work space
134: */
135: MFNDestroy(&mfn);
136: MatDestroy(&A);
137: VecDestroy(&v);
138: VecDestroy(&y);
139: VecDestroy(&z);
140: SlepcFinalize();
141: return ierr;
142: }