Actual source code: ex26.c
slepc-3.22.1 2024-10-28
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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"
15: "To draw the solution run with -mfn_view_solution draw -draw_pause -1\n\n";
17: #include <slepcmfn.h>
19: int main(int argc,char **argv)
20: {
21: Mat A; /* problem matrix */
22: MFN mfn;
23: FN f;
24: PetscReal norm,tol;
25: Vec v,y,z;
26: PetscInt N,n=10,m,Istart,Iend,i,j,II;
27: PetscBool flag;
29: PetscFunctionBeginUser;
30: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
32: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
33: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
34: if (!flag) m=n;
35: N = n*m;
36: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSquare root of Laplacian y=sqrt(A)*e_1, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
38: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
39: Compute the discrete 2-D Laplacian, A
40: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
42: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
43: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
44: PetscCall(MatSetFromOptions(A));
46: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
47: for (II=Istart;II<Iend;II++) {
48: i = II/n; j = II-i*n;
49: if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
50: if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
51: if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
52: if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
53: PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
54: }
56: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
57: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
59: /* set symmetry flag so that solver can exploit it */
60: PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE));
62: /* set v = e_1 */
63: PetscCall(MatCreateVecs(A,NULL,&v));
64: PetscCall(VecSetValue(v,0,1.0,INSERT_VALUES));
65: PetscCall(VecAssemblyBegin(v));
66: PetscCall(VecAssemblyEnd(v));
67: PetscCall(VecDuplicate(v,&y));
68: PetscCall(VecDuplicate(v,&z));
70: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71: Create the solver, set the matrix and the function
72: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73: PetscCall(MFNCreate(PETSC_COMM_WORLD,&mfn));
74: PetscCall(MFNSetOperator(mfn,A));
75: PetscCall(MFNGetFN(mfn,&f));
76: PetscCall(FNSetType(f,FNSQRT));
77: PetscCall(MFNSetErrorIfNotConverged(mfn,PETSC_TRUE));
78: PetscCall(MFNSetFromOptions(mfn));
80: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81: First solve: y=sqrt(A)*v
82: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
84: PetscCall(MFNSolve(mfn,v,y));
85: PetscCall(VecNorm(y,NORM_2,&norm));
86: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Intermediate vector has norm %g\n",(double)norm));
88: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89: Second solve: z=sqrt(A)*y and compare against A*v
90: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
92: PetscCall(MFNSolve(mfn,y,z));
93: PetscCall(MFNGetTolerances(mfn,&tol,NULL));
95: PetscCall(MatMult(A,v,y)); /* overwrite y */
96: PetscCall(VecAXPY(y,-1.0,z));
97: PetscCall(VecNorm(y,NORM_2,&norm));
99: if (norm<tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Error norm is less than the requested tolerance\n\n"));
100: else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Error norm larger than tolerance: %3.1e\n\n",(double)norm));
102: /*
103: Free work space
104: */
105: PetscCall(MFNDestroy(&mfn));
106: PetscCall(MatDestroy(&A));
107: PetscCall(VecDestroy(&v));
108: PetscCall(VecDestroy(&y));
109: PetscCall(VecDestroy(&z));
110: PetscCall(SlepcFinalize());
111: return 0;
112: }
114: /*TEST
116: test:
117: suffix: 1
118: args: -mfn_tol 1e-4
120: TEST*/