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