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