Actual source code: test3.c
slepc-3.16.3 2022-04-11
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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[] = "Test matrix exponential.\n\n";
13: #include <slepcfn.h>
15: /*
16: Compute matrix exponential B = expm(A)
17: */
18: PetscErrorCode TestMatExp(FN fn,Mat A,PetscViewer viewer,PetscBool verbose,PetscBool inplace,PetscBool checkerror)
19: {
21: PetscScalar tau,eta;
22: PetscBool set,flg;
23: PetscInt n;
24: Mat F,R,Finv,Acopy;
25: Vec v,f0;
26: FN finv;
27: PetscReal nrm,nrmf;
30: MatGetSize(A,&n,NULL);
31: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&F);
32: PetscObjectSetName((PetscObject)F,"F");
33: /* compute matrix exponential */
34: if (inplace) {
35: MatCopy(A,F,SAME_NONZERO_PATTERN);
36: MatIsHermitianKnown(A,&set,&flg);
37: if (set && flg) { MatSetOption(F,MAT_HERMITIAN,PETSC_TRUE); }
38: FNEvaluateFunctionMat(fn,F,NULL);
39: } else {
40: MatDuplicate(A,MAT_COPY_VALUES,&Acopy);
41: FNEvaluateFunctionMat(fn,A,F);
42: /* check that A has not been modified */
43: MatAXPY(Acopy,-1.0,A,SAME_NONZERO_PATTERN);
44: MatNorm(Acopy,NORM_1,&nrm);
45: if (nrm>100*PETSC_MACHINE_EPSILON) {
46: PetscPrintf(PETSC_COMM_WORLD,"Warning: the input matrix has changed by %g\n",(double)nrm);
47: }
48: MatDestroy(&Acopy);
49: }
50: if (verbose) {
51: PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
52: MatView(A,viewer);
53: PetscPrintf(PETSC_COMM_WORLD,"Computed expm(A) - - - - - - -\n");
54: MatView(F,viewer);
55: }
56: /* print matrix norm for checking */
57: MatNorm(F,NORM_1,&nrmf);
58: PetscPrintf(PETSC_COMM_WORLD,"The 1-norm of f(A) is %g\n",(double)nrmf);
59: if (checkerror) {
60: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Finv);
61: PetscObjectSetName((PetscObject)Finv,"Finv");
62: FNGetScale(fn,&tau,&eta);
63: /* compute inverse exp(-tau*A)/eta */
64: FNCreate(PETSC_COMM_WORLD,&finv);
65: FNSetType(finv,FNEXP);
66: FNSetFromOptions(finv);
67: FNSetScale(finv,-tau,1.0/eta);
68: if (inplace) {
69: MatCopy(A,Finv,SAME_NONZERO_PATTERN);
70: MatIsHermitianKnown(A,&set,&flg);
71: if (set && flg) { MatSetOption(Finv,MAT_HERMITIAN,PETSC_TRUE); }
72: FNEvaluateFunctionMat(finv,Finv,NULL);
73: } else {
74: FNEvaluateFunctionMat(finv,A,Finv);
75: }
76: FNDestroy(&finv);
77: /* check error ||F*Finv-I||_F */
78: MatMatMult(F,Finv,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&R);
79: MatShift(R,-1.0);
80: MatNorm(R,NORM_FROBENIUS,&nrm);
81: if (nrm<100*PETSC_MACHINE_EPSILON) {
82: PetscPrintf(PETSC_COMM_WORLD,"||exp(A)*exp(-A)-I||_F < 100*eps\n");
83: } else {
84: PetscPrintf(PETSC_COMM_WORLD,"||exp(A)*exp(-A)-I||_F = %g\n",(double)nrm);
85: }
86: MatDestroy(&R);
87: MatDestroy(&Finv);
88: }
89: /* check FNEvaluateFunctionMatVec() */
90: MatCreateVecs(A,&v,&f0);
91: MatGetColumnVector(F,f0,0);
92: FNEvaluateFunctionMatVec(fn,A,v);
93: VecAXPY(v,-1.0,f0);
94: VecNorm(v,NORM_2,&nrm);
95: if (nrm/nrmf>100*PETSC_MACHINE_EPSILON) {
96: PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm);
97: }
98: MatDestroy(&F);
99: VecDestroy(&v);
100: VecDestroy(&f0);
101: return(0);
102: }
104: int main(int argc,char **argv)
105: {
107: FN fn;
108: Mat A;
109: PetscInt i,j,n=10;
110: PetscScalar *As;
111: PetscViewer viewer;
112: PetscBool verbose,inplace,checkerror;
114: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
115: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
116: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
117: PetscOptionsHasName(NULL,NULL,"-inplace",&inplace);
118: PetscOptionsHasName(NULL,NULL,"-checkerror",&checkerror);
119: PetscPrintf(PETSC_COMM_WORLD,"Matrix exponential, n=%D.\n",n);
121: /* Create exponential function object */
122: FNCreate(PETSC_COMM_WORLD,&fn);
123: FNSetType(fn,FNEXP);
124: FNSetFromOptions(fn);
126: /* Set up viewer */
127: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
128: FNView(fn,viewer);
129: if (verbose) {
130: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
131: }
133: /* Create matrices */
134: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
135: PetscObjectSetName((PetscObject)A,"A");
137: /* Fill A with a symmetric Toeplitz matrix */
138: MatDenseGetArray(A,&As);
139: for (i=0;i<n;i++) As[i+i*n]=2.0;
140: for (j=1;j<3;j++) {
141: for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
142: }
143: MatDenseRestoreArray(A,&As);
144: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
145: TestMatExp(fn,A,viewer,verbose,inplace,checkerror);
147: /* Repeat with non-symmetric A */
148: MatDenseGetArray(A,&As);
149: for (j=1;j<3;j++) {
150: for (i=0;i<n-j;i++) { As[(i+j)+i*n]=-1.0; }
151: }
152: MatDenseRestoreArray(A,&As);
153: MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
154: TestMatExp(fn,A,viewer,verbose,inplace,checkerror);
156: MatDestroy(&A);
157: FNDestroy(&fn);
158: SlepcFinalize();
159: return ierr;
160: }
162: /*TEST
164: testset:
165: filter: grep -v "computing matrix functions"
166: output_file: output/test3_1.out
167: test:
168: suffix: 1
169: args: -fn_method {{0 1}}
170: test:
171: suffix: 1_subdiagonalpade
172: args: -fn_method {{2 3}}
173: requires: c99_complex !single
174: test:
175: suffix: 1_cuda
176: args: -fn_method 4
177: requires: cuda
178: test:
179: suffix: 1_magma
180: args: -fn_method {{5 6 7 8}}
181: requires: cuda magma
182: test:
183: suffix: 2
184: args: -inplace -fn_method{{0 1}}
185: test:
186: suffix: 2_subdiagonalpade
187: args: -inplace -fn_method{{2 3}}
188: requires: c99_complex !single
189: test:
190: suffix: 2_cuda
191: args: -inplace -fn_method 4
192: requires: cuda
193: test:
194: suffix: 2_magma
195: args: -inplace -fn_method {{5 6 7 8}}
196: requires: cuda magma
198: testset:
199: args: -fn_scale 0.1
200: filter: grep -v "computing matrix functions"
201: output_file: output/test3_3.out
202: test:
203: suffix: 3
204: args: -fn_method {{0 1}}
205: test:
206: suffix: 3_subdiagonalpade
207: args: -fn_method {{2 3}}
208: requires: c99_complex !single
210: testset:
211: args: -n 120 -fn_scale 0.6,1.5
212: filter: grep -v "computing matrix functions"
213: output_file: output/test3_4.out
214: test:
215: suffix: 4
216: args: -fn_method {{0 1}}
217: requires: !single
218: test:
219: suffix: 4_subdiagonalpade
220: args: -fn_method {{2 3}}
221: requires: c99_complex !single
223: test:
224: suffix: 5
225: args: -fn_scale 30 -fn_method {{2 3}}
226: filter: grep -v "computing matrix functions"
227: requires: c99_complex !single
228: output_file: output/test3_5.out
230: test:
231: suffix: 6
232: args: -fn_scale 1e-9 -fn_method {{2 3}}
233: filter: grep -v "computing matrix functions"
234: requires: c99_complex !single
235: output_file: output/test3_6.out
237: TEST*/