Actual source code: test11.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2013, 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 in DSNHEP.\n\n";
24: #include <slepcds.h>
25: #include <slepc-private/dsimpl.h> /* for DSViewMat_Private() */
29: int main(int argc,char **argv)
30: {
32: DS ds;
33: PetscScalar *A;
34: PetscInt i,j,n=10,ld;
35: PetscViewer viewer;
36: PetscBool verbose;
38: SlepcInitialize(&argc,&argv,(char*)0,help);
39: PetscOptionsGetInt(NULL,"-n",&n,NULL);
40: PetscPrintf(PETSC_COMM_WORLD,"Compute non-symmetric matrix exponential - dimension %D.\n",n);
41: PetscOptionsHasName(NULL,"-verbose",&verbose);
43: /* Create DS object */
44: DSCreate(PETSC_COMM_WORLD,&ds);
45: DSSetType(ds,DSNHEP);
46: DSSetFromOptions(ds);
47: ld = n+2; /* test leading dimension larger than n */
48: DSAllocate(ds,ld);
49: DSSetDimensions(ds,n,0,0,0);
51: /* Set up viewer */
52: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
53: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
54: DSView(ds,viewer);
55: PetscViewerPopFormat(viewer);
56: if (verbose) {
57: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
58: }
60: /* Fill with a symmetric Toeplitz matrix */
61: DSGetArray(ds,DS_MAT_A,&A);
62: for (i=1;i<n;i++) A[i+(i-1)*ld]=-1.0;
63: for (j=0;j<4;j++) {
64: for (i=0;i<n-j;i++) A[i+(i+j)*ld]=1.0;
65: }
66: DSRestoreArray(ds,DS_MAT_A,&A);
67: DSSetState(ds,DS_STATE_INTERMEDIATE);
68: if (verbose) {
69: PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
70: DSView(ds,viewer);
71: }
73: /* Compute matrix exponential */
74: DSComputeFunction(ds,SLEPC_FUNCTION_EXP);
75: if (verbose) {
76: PetscPrintf(PETSC_COMM_WORLD,"Computed f(A) - - - - - - -\n");
77: DSViewMat_Private(ds,viewer,DS_MAT_F);
78: }
80: DSDestroy(&ds);
81: SlepcFinalize();
82: return 0;
83: }