Actual source code: test2.c
slepc-3.17.2 2022-08-09
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[] = "Tests the case when both arguments of MFNSolve() are the same Vec.\n\n"
12: "The command line options are:\n"
13: " -t <sval>, where <sval> = scalar value that multiplies the argument.\n"
14: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
15: " -m <m>, where <m> = number of grid subdivisions in y dimension.\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;
25: PetscScalar t=0.3;
26: PetscInt N,n=25,m,Istart,Iend,II,i,j;
27: PetscBool flag;
28: Vec v,y;
30: SlepcInitialize(&argc,&argv,(char*)0,help);
32: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
33: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
34: if (!flag) m=n;
35: N = n*m;
36: PetscOptionsGetScalar(NULL,NULL,"-t",&t,NULL);
37: PetscPrintf(PETSC_COMM_WORLD,"\nMatrix exponential y=exp(t*A)*e, of the 2-D Laplacian, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m);
39: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
40: Build the 2-D Laplacian
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 v = ones(n,1) */
62: MatCreateVecs(A,NULL,&y);
63: MatCreateVecs(A,NULL,&v);
64: VecSet(v,1.0);
66: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67: Create the solver and set various options
68: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
70: FNCreate(PETSC_COMM_WORLD,&f);
71: FNSetType(f,FNEXP);
73: MFNCreate(PETSC_COMM_WORLD,&mfn);
74: MFNSetOperator(mfn,A);
75: MFNSetFN(mfn,f);
76: MFNSetErrorIfNotConverged(mfn,PETSC_TRUE);
77: MFNSetFromOptions(mfn);
79: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80: Solve the problem, y=exp(t*A)*v
81: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83: FNSetScale(f,t,1.0);
84: MFNSolve(mfn,v,y);
85: VecNorm(y,NORM_2,&norm);
86: PetscPrintf(PETSC_COMM_WORLD," Computed vector at time t=%.4g has norm %g\n\n",(double)PetscRealPart(t),(double)norm);
88: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89: Repeat the computation in two steps, overwriting v:
90: v=exp(0.5*t*A)*v, v=exp(0.5*t*A)*v
91: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93: FNSetScale(f,0.5*t,1.0);
94: MFNSolve(mfn,v,v);
95: MFNSolve(mfn,v,v);
96: /* compute norm of difference */
97: VecAXPY(y,-1.0,v);
98: VecNorm(y,NORM_2,&norm);
99: if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD," The norm of the difference is <100*eps\n\n");
100: else PetscPrintf(PETSC_COMM_WORLD," The norm of the difference is %g\n\n",(double)norm);
102: /*
103: Free work space
104: */
105: MFNDestroy(&mfn);
106: FNDestroy(&f);
107: MatDestroy(&A);
108: VecDestroy(&v);
109: VecDestroy(&y);
110: SlepcFinalize();
111: return 0;
112: }
114: /*TEST
116: testset:
117: args: -mfn_type {{krylov expokit}}
118: output_file: output/test2_1.out
119: test:
120: suffix: 1
121: test:
122: suffix: 1_cuda
123: args: -mat_type aijcusparse
124: requires: cuda
126: test:
127: suffix: 3
128: args: -mfn_type expokit -t 0.6 -mfn_ncv 35
129: requires: !__float128
131: TEST*/