Actual source code: test16.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 pseudo-orthogonalization.\n\n";
13: #include <slepcds.h>
15: int main(int argc,char **argv)
16: {
18: DS ds;
19: PetscReal *s,*ns;
20: PetscScalar *A;
21: PetscInt i,j,n=10;
22: PetscViewer viewer;
23: PetscBool verbose;
25: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
26: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
27: PetscPrintf(PETSC_COMM_WORLD,"Test pseudo-orthogonalization for GHIEP - dimension %D.\n",n);
28: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
30: /* Create DS object */
31: DSCreate(PETSC_COMM_WORLD,&ds);
32: DSSetType(ds,DSGHIEP);
33: DSSetFromOptions(ds);
34: DSAllocate(ds,n);
35: DSSetDimensions(ds,n,0,0,0);
37: /* Set up viewer */
38: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
39: if (verbose) {
40: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
41: }
43: /* Fill with a symmetric Toeplitz matrix */
44: DSGetArray(ds,DS_MAT_A,&A);
45: for (i=0;i<n;i++) A[i+i*n]=2.0;
46: for (j=1;j<3;j++) {
47: for (i=0;i<n-j;i++) { A[i+(i+j)*n]=1.0; A[(i+j)+i*n]=1.0; }
48: }
49: for (j=1;j<3;j++) { A[0+j*n]=-1.0*(j+2); A[j+0*n]=-1.0*(j+2); }
50: DSRestoreArray(ds,DS_MAT_A,&A);
51: DSSetState(ds,DS_STATE_RAW);
52: if (verbose) {
53: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
54: DSView(ds,viewer);
55: }
57: /* Signature matrix */
58: PetscMalloc2(n,&s,n,&ns);
59: s[0] = -1.0;
60: for (i=1;i<n-1;i++) s[i]=1.0;
61: s[n-1] = -1.0;
63: /* Orthogonalize and show signature */
64: DSPseudoOrthogonalize(ds,DS_MAT_A,n,s,NULL,ns);
65: if (verbose) {
66: PetscPrintf(PETSC_COMM_WORLD,"After pseudo-orthogonalize - - - - - - - - -\n");
67: DSView(ds,viewer);
68: }
69: PetscPrintf(PETSC_COMM_WORLD,"Resulting signature:\n");
70: for (i=0;i<n;i++) {
71: PetscPrintf(PETSC_COMM_WORLD,"%g\n",(double)ns[i]);
72: }
73: PetscFree2(s,ns);
75: DSDestroy(&ds);
76: SlepcFinalize();
77: return ierr;
78: }
80: /*TEST
82: test:
83: suffix: 1
85: TEST*/