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