Actual source code: test12.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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 function evaluation via diagonalization.\n\n";
13: #include <slepcfn.h>
15: int main(int argc,char **argv)
16: {
18: FN fn;
19: Mat A,F,G;
20: PetscInt i,j,n=10;
21: PetscReal nrm;
22: PetscScalar *As,alpha,beta;
23: PetscViewer viewer;
24: PetscBool verbose;
26: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
27: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
28: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
29: PetscPrintf(PETSC_COMM_WORLD,"Matrix function of symmetric/Hermitian matrix, n=%D.\n",n);
31: /* Create function object */
32: FNCreate(PETSC_COMM_WORLD,&fn);
33: FNSetType(fn,FNEXP); /* default to exponential */
34: #if defined(PETSC_USE_COMPLEX)
35: alpha = 0.3+0.8*PETSC_i;
36: beta = 1.1-0.1*PETSC_i;
37: #else
38: alpha = 0.3;
39: beta = 1.1;
40: #endif
41: FNSetScale(fn,alpha,beta);
42: FNSetFromOptions(fn);
44: /* Set up viewer */
45: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
46: if (verbose) {
47: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
48: }
50: /* Create a symmetric/Hermitian Toeplitz matrix */
51: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
52: PetscObjectSetName((PetscObject)A,"A");
53: MatDenseGetArray(A,&As);
54: for (i=0;i<n;i++) As[i+i*n]=2.0;
55: for (j=1;j<3;j++) {
56: for (i=0;i<n-j;i++) {
57: #if defined(PETSC_USE_COMPLEX)
58: As[i+(i+j)*n]=1.0+0.1*PETSC_i; As[(i+j)+i*n]=1.0-0.1*PETSC_i;
59: #else
60: As[i+(i+j)*n]=0.5; As[(i+j)+i*n]=0.5;
61: #endif
62: }
63: }
64: MatDenseRestoreArray(A,&As);
65: if (verbose) {
66: PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
67: MatView(A,viewer);
68: }
70: /* compute matrix function */
71: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&F);
72: PetscObjectSetName((PetscObject)F,"F");
73: FNEvaluateFunctionMat(fn,A,F);
74: if (verbose) {
75: PetscPrintf(PETSC_COMM_WORLD,"Computed f(A) - - - - - - -\n");
76: MatView(F,viewer);
77: }
79: /* Repeat with MAT_HERMITIAN flag set */
80: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
81: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&G);
82: PetscObjectSetName((PetscObject)G,"G");
83: FNEvaluateFunctionMat(fn,A,G);
84: if (verbose) {
85: PetscPrintf(PETSC_COMM_WORLD,"Computed f(A) symm - - - - - - -\n");
86: MatView(G,viewer);
87: }
89: /* compare the two results */
90: MatAXPY(F,-1.0,G,SAME_NONZERO_PATTERN);
91: MatNorm(F,NORM_FROBENIUS,&nrm);
92: if (nrm>100*PETSC_MACHINE_EPSILON) {
93: PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of F-G is %g\n",(double)nrm);
94: } else {
95: PetscPrintf(PETSC_COMM_WORLD,"Computed results match.\n");
96: }
98: MatDestroy(&A);
99: MatDestroy(&F);
100: MatDestroy(&G);
101: FNDestroy(&fn);
102: SlepcFinalize();
103: return ierr;
104: }
106: /*TEST
108: test:
109: suffix: 1
110: nsize: 1
111: args: -fn_type {{exp sqrt}shared output}
112: output_file: output/test12_1.out
114: test:
115: suffix: 1_rational
116: nsize: 1
117: args: -fn_type rational -fn_rational_numerator 2,-1.5 -fn_rational_denominator 1,0.8
118: output_file: output/test12_1.out
120: TEST*/