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