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