Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
11 : static char help[] = "Test matrix logarithm.\n\n";
12 :
13 : #include <slepcfn.h>
14 :
15 : /*
16 : Compute matrix logarithm B = logm(A)
17 : */
18 2 : PetscErrorCode TestMatLog(FN fn,Mat A,PetscViewer viewer,PetscBool verbose,PetscBool inplace)
19 : {
20 2 : PetscBool set,flg;
21 2 : PetscScalar tau,eta;
22 2 : PetscInt n;
23 2 : Mat F,R;
24 2 : Vec v,f0;
25 2 : FN fnexp;
26 2 : PetscReal nrm;
27 :
28 2 : PetscFunctionBeginUser;
29 2 : PetscCall(MatGetSize(A,&n,NULL));
30 2 : PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&F));
31 2 : PetscCall(PetscObjectSetName((PetscObject)F,"F"));
32 2 : PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&R));
33 2 : PetscCall(PetscObjectSetName((PetscObject)R,"R"));
34 2 : PetscCall(FNGetScale(fn,&tau,&eta));
35 : /* compute matrix logarithm */
36 2 : if (inplace) {
37 0 : PetscCall(MatCopy(A,F,SAME_NONZERO_PATTERN));
38 0 : PetscCall(MatIsHermitianKnown(A,&set,&flg));
39 0 : if (set && flg) PetscCall(MatSetOption(F,MAT_HERMITIAN,PETSC_TRUE));
40 0 : PetscCall(FNEvaluateFunctionMat(fn,F,NULL));
41 2 : } else PetscCall(FNEvaluateFunctionMat(fn,A,F));
42 2 : if (verbose) {
43 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n"));
44 0 : PetscCall(MatView(A,viewer));
45 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed logm(A) - - - - - - -\n"));
46 0 : PetscCall(MatView(F,viewer));
47 : }
48 : /* check error ||expm(F)-A||_F */
49 2 : PetscCall(FNCreate(PETSC_COMM_WORLD,&fnexp));
50 2 : PetscCall(FNSetType(fnexp,FNEXP));
51 2 : PetscCall(MatCopy(F,R,SAME_NONZERO_PATTERN));
52 2 : if (eta!=1.0) PetscCall(MatScale(R,1.0/eta));
53 2 : PetscCall(FNEvaluateFunctionMat(fnexp,R,NULL));
54 2 : PetscCall(FNDestroy(&fnexp));
55 2 : PetscCall(MatAXPY(R,-tau,A,SAME_NONZERO_PATTERN));
56 2 : PetscCall(MatNorm(R,NORM_FROBENIUS,&nrm));
57 2 : if (nrm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"||expm(F)-A||_F < 100*eps\n"));
58 0 : else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"||expm(F)-A||_F = %g\n",(double)nrm));
59 : /* check FNEvaluateFunctionMatVec() */
60 2 : PetscCall(MatCreateVecs(A,&v,&f0));
61 2 : PetscCall(MatGetColumnVector(F,f0,0));
62 2 : PetscCall(FNEvaluateFunctionMatVec(fn,A,v));
63 2 : PetscCall(VecAXPY(v,-1.0,f0));
64 2 : PetscCall(VecNorm(v,NORM_2,&nrm));
65 2 : if (nrm>100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm));
66 2 : PetscCall(MatDestroy(&F));
67 2 : PetscCall(MatDestroy(&R));
68 2 : PetscCall(VecDestroy(&v));
69 2 : PetscCall(VecDestroy(&f0));
70 2 : PetscFunctionReturn(PETSC_SUCCESS);
71 : }
72 :
73 2 : int main(int argc,char **argv)
74 : {
75 2 : FN fn;
76 2 : Mat A;
77 2 : PetscInt i,j,n=10;
78 2 : PetscScalar *As;
79 2 : PetscViewer viewer;
80 2 : PetscBool verbose,inplace,random,triang;
81 :
82 2 : PetscFunctionBeginUser;
83 2 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
84 2 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
85 2 : PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
86 2 : PetscCall(PetscOptionsHasName(NULL,NULL,"-inplace",&inplace));
87 2 : PetscCall(PetscOptionsHasName(NULL,NULL,"-random",&random));
88 2 : PetscCall(PetscOptionsHasName(NULL,NULL,"-triang",&triang));
89 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Matrix logarithm, n=%" PetscInt_FMT ".\n",n));
90 :
91 : /* Create logarithm function object */
92 2 : PetscCall(FNCreate(PETSC_COMM_WORLD,&fn));
93 2 : PetscCall(FNSetType(fn,FNLOG));
94 2 : PetscCall(FNSetFromOptions(fn));
95 :
96 : /* Set up viewer */
97 2 : PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
98 2 : PetscCall(FNView(fn,viewer));
99 2 : if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
100 :
101 : /* Create matrices */
102 2 : PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A));
103 2 : PetscCall(PetscObjectSetName((PetscObject)A,"A"));
104 :
105 2 : if (random) PetscCall(MatSetRandom(A,NULL));
106 : else {
107 : /* Fill A with a non-symmetric Toeplitz matrix */
108 2 : PetscCall(MatDenseGetArray(A,&As));
109 152 : for (i=0;i<n;i++) As[i+i*n]=2.0;
110 6 : for (j=1;j<3;j++) {
111 298 : for (i=0;i<n-j;i++) {
112 294 : As[i+(i+j)*n]=1.0;
113 294 : if (!triang) As[(i+j)+i*n]=-1.0;
114 : }
115 : }
116 2 : As[(n-1)*n] = -5.0;
117 2 : As[0] = 2.01;
118 2 : PetscCall(MatDenseRestoreArray(A,&As));
119 : }
120 2 : PetscCall(TestMatLog(fn,A,viewer,verbose,inplace));
121 :
122 2 : PetscCall(MatDestroy(&A));
123 2 : PetscCall(FNDestroy(&fn));
124 2 : PetscCall(SlepcFinalize());
125 : return 0;
126 : }
127 :
128 : /*TEST
129 :
130 : testset:
131 : filter: grep -v "computing matrix functions"
132 : output_file: output/test13_1.out
133 : test:
134 : suffix: 1
135 : args: -fn_scale .04,2 -n 75
136 : requires: c99_complex !__float128
137 : test:
138 : suffix: 1_triang
139 : args: -fn_scale .04,2 -n 75 -triang
140 : requires: c99_complex !__float128
141 : test:
142 : suffix: 1_random
143 : args: -fn_scale .02,2 -n 75 -random
144 : requires: complex !__float128
145 : filter_output: sed -e 's/04/02/'
146 :
147 : TEST*/
|