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 DSGHEP.\n\n";
12 :
13 : #include <slepcds.h>
14 :
15 : /*
16 : Compute the norm of the j-th column of matrix mat in ds
17 : */
18 1 : PetscErrorCode ComputeNorm(DS ds,DSMatType mat,PetscInt j,PetscReal *onrm)
19 : {
20 1 : PetscScalar *X;
21 1 : PetscReal aux,nrm=0.0;
22 1 : PetscInt i,n,ld;
23 :
24 1 : PetscFunctionBeginUser;
25 1 : PetscCall(DSGetLeadingDimension(ds,&ld));
26 1 : PetscCall(DSGetDimensions(ds,&n,NULL,NULL,NULL));
27 1 : PetscCall(DSGetArray(ds,mat,&X));
28 11 : for (i=0;i<n;i++) {
29 10 : aux = PetscAbsScalar(X[i+j*ld]);
30 10 : nrm += aux*aux;
31 : }
32 1 : PetscCall(DSRestoreArray(ds,mat,&X));
33 1 : *onrm = PetscSqrtReal(nrm);
34 1 : PetscFunctionReturn(PETSC_SUCCESS);
35 : }
36 :
37 1 : int main(int argc,char **argv)
38 : {
39 1 : DS ds;
40 1 : SlepcSC sc;
41 1 : PetscReal re;
42 1 : PetscScalar *A,*B,*eig;
43 1 : PetscReal nrm;
44 1 : PetscInt i,j,n=10,ld;
45 1 : PetscViewer viewer;
46 1 : PetscBool verbose;
47 :
48 1 : PetscFunctionBeginUser;
49 1 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
50 1 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
51 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a System of type GHEP - dimension %" PetscInt_FMT ".\n",n));
52 1 : PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
53 :
54 : /* Create DS object */
55 1 : PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
56 1 : PetscCall(DSSetType(ds,DSGHEP));
57 1 : PetscCall(DSSetFromOptions(ds));
58 1 : ld = n+2; /* test leading dimension larger than n */
59 1 : PetscCall(DSAllocate(ds,ld));
60 1 : PetscCall(DSSetDimensions(ds,n,0,0));
61 :
62 : /* Set up viewer */
63 1 : PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
64 1 : PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
65 1 : PetscCall(DSView(ds,viewer));
66 1 : PetscCall(PetscViewerPopFormat(viewer));
67 1 : if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
68 :
69 : /* Fill with a symmetric Toeplitz matrix */
70 1 : PetscCall(DSGetArray(ds,DS_MAT_A,&A));
71 1 : PetscCall(DSGetArray(ds,DS_MAT_B,&B));
72 11 : for (i=0;i<n;i++) A[i+i*ld]=2.0;
73 3 : for (j=1;j<3;j++) {
74 19 : for (i=0;i<n-j;i++) { A[i+(i+j)*ld]=1.0; A[(i+j)+i*ld]=1.0; }
75 : }
76 3 : for (j=1;j<3;j++) { A[0+j*ld]=-1.0*(j+2); A[j+0*ld]=-1.0*(j+2); }
77 : /* Diagonal matrix */
78 11 : for (i=0;i<n;i++) B[i+i*ld]=0.1*(i+1);
79 1 : PetscCall(DSRestoreArray(ds,DS_MAT_A,&A));
80 1 : PetscCall(DSRestoreArray(ds,DS_MAT_B,&B));
81 1 : PetscCall(DSSetState(ds,DS_STATE_RAW));
82 1 : if (verbose) {
83 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
84 0 : PetscCall(DSView(ds,viewer));
85 : }
86 :
87 : /* Solve */
88 1 : PetscCall(PetscMalloc1(n,&eig));
89 1 : PetscCall(PetscNew(&sc));
90 1 : sc->comparison = SlepcCompareLargestMagnitude;
91 1 : sc->comparisonctx = NULL;
92 1 : sc->map = NULL;
93 1 : sc->mapobj = NULL;
94 1 : PetscCall(DSSetSlepcSC(ds,sc));
95 1 : PetscCall(DSSolve(ds,eig,NULL));
96 1 : PetscCall(DSSort(ds,eig,NULL,NULL,NULL,NULL));
97 1 : if (verbose) {
98 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
99 0 : PetscCall(DSView(ds,viewer));
100 : }
101 :
102 : /* Print eigenvalues */
103 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n"));
104 11 : for (i=0;i<n;i++) {
105 10 : re = PetscRealPart(eig[i]);
106 10 : PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re));
107 : }
108 :
109 : /* Eigenvectors */
110 1 : j = 0;
111 1 : PetscCall(DSVectors(ds,DS_MAT_X,&j,NULL)); /* all eigenvectors */
112 1 : PetscCall(ComputeNorm(ds,DS_MAT_X,0,&nrm));
113 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st vector = %.3f\n",(double)nrm));
114 1 : PetscCall(DSVectors(ds,DS_MAT_X,NULL,NULL)); /* all eigenvectors */
115 1 : if (verbose) {
116 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After vectors - - - - - - - - -\n"));
117 0 : PetscCall(DSView(ds,viewer));
118 : }
119 :
120 1 : PetscCall(PetscFree(eig));
121 1 : PetscCall(PetscFree(sc));
122 1 : PetscCall(DSDestroy(&ds));
123 1 : PetscCall(SlepcFinalize());
124 : return 0;
125 : }
126 :
127 : /*TEST
128 :
129 : test:
130 : suffix: 1
131 : requires: !single
132 :
133 : TEST*/
|