Actual source code: test12.c
slepc-3.21.1 2024-04-26
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: PetscFunctionBeginUser;
26: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
27: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
28: PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
29: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Matrix function of symmetric/Hermitian matrix, n=%" PetscInt_FMT ".\n",n));
31: /* Create function object */
32: PetscCall(FNCreate(PETSC_COMM_WORLD,&fn));
33: PetscCall(FNSetType(fn,FNEXP)); /* default to exponential */
34: #if defined(PETSC_USE_COMPLEX)
35: alpha = PetscCMPLX(0.3,0.8);
36: beta = PetscCMPLX(1.1,-0.1);
37: #else
38: alpha = 0.3;
39: beta = 1.1;
40: #endif
41: PetscCall(FNSetScale(fn,alpha,beta));
42: PetscCall(FNSetFromOptions(fn));
44: /* Set up viewer */
45: PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
46: if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
48: /* Create a symmetric/Hermitian Toeplitz matrix */
49: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A));
50: PetscCall(PetscObjectSetName((PetscObject)A,"A"));
51: PetscCall(MatDenseGetArray(A,&As));
52: for (i=0;i<n;i++) As[i+i*n]=2.0;
53: for (j=1;j<3;j++) {
54: for (i=0;i<n-j;i++) {
55: #if defined(PETSC_USE_COMPLEX)
56: As[i+(i+j)*n]=PetscCMPLX(1.0,0.1); As[(i+j)+i*n]=PetscCMPLX(1.0,-0.1);
57: #else
58: As[i+(i+j)*n]=0.5; As[(i+j)+i*n]=0.5;
59: #endif
60: }
61: }
62: PetscCall(MatDenseRestoreArray(A,&As));
63: if (verbose) {
64: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n"));
65: PetscCall(MatView(A,viewer));
66: }
68: /* compute matrix function */
69: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&F));
70: PetscCall(PetscObjectSetName((PetscObject)F,"F"));
71: PetscCall(FNEvaluateFunctionMat(fn,A,F));
72: if (verbose) {
73: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed f(A) - - - - - - -\n"));
74: PetscCall(MatView(F,viewer));
75: }
77: /* Repeat with MAT_HERMITIAN flag set */
78: PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE));
79: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&G));
80: PetscCall(PetscObjectSetName((PetscObject)G,"G"));
81: PetscCall(FNEvaluateFunctionMat(fn,A,G));
82: if (verbose) {
83: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed f(A) symm - - - - - - -\n"));
84: PetscCall(MatView(G,viewer));
85: }
87: /* compare the two results */
88: PetscCall(MatAXPY(F,-1.0,G,SAME_NONZERO_PATTERN));
89: PetscCall(MatNorm(F,NORM_FROBENIUS,&nrm));
90: if (nrm>100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of F-G is %g\n",(double)nrm));
91: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed results match.\n"));
93: PetscCall(MatDestroy(&A));
94: PetscCall(MatDestroy(&F));
95: PetscCall(MatDestroy(&G));
96: PetscCall(FNDestroy(&fn));
97: PetscCall(SlepcFinalize());
98: return 0;
99: }
101: /*TEST
103: test:
104: suffix: 1
105: nsize: 1
106: args: -fn_type {{exp sqrt}shared output}
107: output_file: output/test12_1.out
109: test:
110: suffix: 1_rational
111: nsize: 1
112: args: -fn_type rational -fn_rational_numerator 2,-1.5 -fn_rational_denominator 1,0.8
113: output_file: output/test12_1.out
115: TEST*/