Actual source code: test2.c
slepc-3.18.2 2023-01-26
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;
31: SlepcInitialize(&argc,&argv,(char*)0,help);
33: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
34: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
35: if (!flag) m=n;
36: N = n*m;
37: PetscOptionsGetScalar(NULL,NULL,"-t",&t,NULL);
38: 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);
40: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
41: Build the 2-D Laplacian
42: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
44: MatCreate(PETSC_COMM_WORLD,&A);
45: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
46: MatSetFromOptions(A);
47: MatSetUp(A);
49: MatGetOwnershipRange(A,&Istart,&Iend);
50: for (II=Istart;II<Iend;II++) {
51: i = II/n; j = II-i*n;
52: if (i>0) MatSetValue(A,II,II-n,-1.0,INSERT_VALUES);
53: if (i<m-1) MatSetValue(A,II,II+n,-1.0,INSERT_VALUES);
54: if (j>0) MatSetValue(A,II,II-1,-1.0,INSERT_VALUES);
55: if (j<n-1) MatSetValue(A,II,II+1,-1.0,INSERT_VALUES);
56: MatSetValue(A,II,II,4.0,INSERT_VALUES);
57: }
59: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
60: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
62: /* set v = ones(n,1) */
63: MatCreateVecs(A,NULL,&y);
64: MatCreateVecs(A,NULL,&v);
65: VecSet(v,1.0);
67: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68: Create the solver and set various options
69: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71: FNCreate(PETSC_COMM_WORLD,&f);
72: FNSetType(f,FNEXP);
74: MFNCreate(PETSC_COMM_WORLD,&mfn);
75: MFNSetOperator(mfn,A);
76: MFNSetFN(mfn,f);
77: MFNSetErrorIfNotConverged(mfn,PETSC_TRUE);
78: MFNSetFromOptions(mfn);
80: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81: Solve the problem, y=exp(t*A)*v
82: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
84: FNSetScale(f,t,1.0);
85: MFNSolve(mfn,v,y);
86: VecNorm(y,NORM_2,&norm);
87: PetscPrintf(PETSC_COMM_WORLD," Computed vector at time t=%.4g has norm %g\n\n",(double)PetscRealPart(t),(double)norm);
89: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
90: Repeat the computation in two steps, overwriting v:
91: v=exp(0.5*t*A)*v, v=exp(0.5*t*A)*v
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94: FNSetScale(f,0.5*t,1.0);
95: MFNSolve(mfn,v,v);
96: MFNSolve(mfn,v,v);
97: /* compute norm of difference */
98: VecAXPY(y,-1.0,v);
99: VecNorm(y,NORM_2,&norm);
100: if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD," The norm of the difference is <100*eps\n\n");
101: else PetscPrintf(PETSC_COMM_WORLD," The norm of the difference is %g\n\n",(double)norm);
103: /*
104: Free work space
105: */
106: MFNDestroy(&mfn);
107: FNDestroy(&f);
108: MatDestroy(&A);
109: VecDestroy(&v);
110: VecDestroy(&y);
111: SlepcFinalize();
112: return 0;
113: }
115: /*TEST
117: testset:
118: args: -mfn_type {{krylov expokit}}
119: output_file: output/test2_1.out
120: test:
121: suffix: 1
122: test:
123: suffix: 1_cuda
124: args: -mat_type aijcusparse
125: requires: cuda
127: test:
128: suffix: 3
129: args: -mfn_type expokit -t 0.6 -mfn_ncv 35
130: requires: !__float128
132: TEST*/