Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
11 : static char help[] = "Test pseudo-orthogonalization.\n\n";
12 :
13 : #include <slepcds.h>
14 :
15 1 : int main(int argc,char **argv)
16 : {
17 1 : DS ds;
18 1 : PetscReal *s,*ns;
19 1 : PetscScalar *A;
20 1 : PetscInt i,j,n=10;
21 1 : PetscViewer viewer;
22 1 : PetscBool verbose;
23 :
24 1 : PetscFunctionBeginUser;
25 1 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
26 1 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
27 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test pseudo-orthogonalization for GHIEP - dimension %" PetscInt_FMT ".\n",n));
28 1 : PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
29 :
30 : /* Create DS object */
31 1 : PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
32 1 : PetscCall(DSSetType(ds,DSGHIEP));
33 1 : PetscCall(DSSetFromOptions(ds));
34 1 : PetscCall(DSAllocate(ds,n));
35 1 : PetscCall(DSSetDimensions(ds,n,0,0));
36 :
37 : /* Set up viewer */
38 1 : PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
39 1 : if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
40 :
41 : /* Fill with a symmetric Toeplitz matrix */
42 1 : PetscCall(DSGetArray(ds,DS_MAT_A,&A));
43 11 : for (i=0;i<n;i++) A[i+i*n]=2.0;
44 3 : for (j=1;j<3;j++) {
45 19 : for (i=0;i<n-j;i++) { A[i+(i+j)*n]=1.0; A[(i+j)+i*n]=1.0; }
46 : }
47 3 : for (j=1;j<3;j++) { A[0+j*n]=-1.0*(j+2); A[j+0*n]=-1.0*(j+2); }
48 1 : PetscCall(DSRestoreArray(ds,DS_MAT_A,&A));
49 1 : PetscCall(DSSetState(ds,DS_STATE_RAW));
50 1 : if (verbose) {
51 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
52 0 : PetscCall(DSView(ds,viewer));
53 : }
54 :
55 : /* Signature matrix */
56 1 : PetscCall(PetscMalloc2(n,&s,n,&ns));
57 1 : s[0] = -1.0;
58 9 : for (i=1;i<n-1;i++) s[i]=1.0;
59 1 : s[n-1] = -1.0;
60 :
61 : /* Orthogonalize and show signature */
62 1 : PetscCall(DSPseudoOrthogonalize(ds,DS_MAT_A,n,s,NULL,ns));
63 1 : if (verbose) {
64 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After pseudo-orthogonalize - - - - - - - - -\n"));
65 0 : PetscCall(DSView(ds,viewer));
66 : }
67 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Resulting signature:\n"));
68 11 : for (i=0;i<n;i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%g\n",(double)ns[i]));
69 1 : PetscCall(PetscFree2(s,ns));
70 :
71 1 : PetscCall(DSDestroy(&ds));
72 1 : PetscCall(SlepcFinalize());
73 : return 0;
74 : }
75 :
76 : /*TEST
77 :
78 : test:
79 : suffix: 1
80 :
81 : TEST*/
|