Actual source code: test3.c
slepc-3.6.1 2015-09-03
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: static char help[] = "Test matrix exponential.\n\n";
24: #include <slepcfn.h>
28: int main(int argc,char **argv)
29: {
31: FN fn;
32: Mat A,B;
33: PetscInt i,j,n=10;
34: PetscReal nrm;
35: PetscScalar *As,tau=1.0,eta=1.0;
36: PetscViewer viewer;
37: PetscBool verbose;
39: SlepcInitialize(&argc,&argv,(char*)0,help);
40: PetscOptionsGetInt(NULL,"-n",&n,NULL);
41: PetscOptionsGetScalar(NULL,"-tau",&tau,NULL);
42: PetscOptionsGetScalar(NULL,"-eta",&eta,NULL);
43: PetscOptionsHasName(NULL,"-verbose",&verbose);
44: PetscPrintf(PETSC_COMM_WORLD,"Matrix exponential, n=%D.\n",n);
46: /* Create exponential function eta*exp(tau*x) */
47: FNCreate(PETSC_COMM_WORLD,&fn);
48: FNSetType(fn,FNEXP);
49: FNSetScale(fn,tau,eta);
51: /* Set up viewer */
52: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
53: FNView(fn,viewer);
54: if (verbose) {
55: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
56: }
58: /* Create matrices */
59: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
60: PetscObjectSetName((PetscObject)A,"A");
61: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&B);
62: PetscObjectSetName((PetscObject)B,"B");
64: /* Fill A with a symmetric Toeplitz matrix */
65: MatDenseGetArray(A,&As);
66: for (i=0;i<n;i++) As[i+i*n]=2.0;
67: for (j=1;j<3;j++) {
68: for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
69: }
70: MatDenseRestoreArray(A,&As);
71: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
72: if (verbose) {
73: PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
74: MatView(A,viewer);
75: }
77: /* Compute matrix exponential */
78: FNEvaluateFunctionMat(fn,A,B);
79: if (verbose) {
80: PetscPrintf(PETSC_COMM_WORLD,"Computed f(A) - - - - - - -\n");
81: MatView(B,viewer);
82: }
83: MatNorm(B,NORM_1,&nrm);
84: PetscPrintf(PETSC_COMM_WORLD,"The 1-norm of f(A) is %g\n",(double)nrm);
86: /* Repeat with non-symmetric A */
87: MatDenseGetArray(A,&As);
88: for (j=1;j<3;j++) {
89: for (i=0;i<n-j;i++) { As[(i+j)+i*n]=-1.0; }
90: }
91: MatDenseRestoreArray(A,&As);
92: MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
93: if (verbose) {
94: PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
95: MatView(A,viewer);
96: }
98: /* Compute matrix exponential */
99: FNEvaluateFunctionMat(fn,A,B);
100: if (verbose) {
101: PetscPrintf(PETSC_COMM_WORLD,"Computed f(A) - - - - - - -\n");
102: MatView(B,viewer);
103: }
104: MatNorm(B,NORM_1,&nrm);
105: PetscPrintf(PETSC_COMM_WORLD,"The 1-norm of f(A) is %g\n",(double)nrm);
107: MatDestroy(&A);
108: MatDestroy(&B);
109: FNDestroy(&fn);
110: SlepcFinalize();
111: return 0;
112: }