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 DSPEP.\n\n";
12 :
13 : #include <slepcds.h>
14 :
15 1 : int main(int argc,char **argv)
16 : {
17 1 : DS ds;
18 1 : SlepcSC sc;
19 1 : Mat X;
20 1 : Vec x0;
21 1 : PetscScalar *K,*C,*M,*wr,*wi,z=1.0;
22 1 : PetscReal re,im,nrm,*pbc;
23 1 : PetscInt i,j,n=10,d=2,ld;
24 1 : PetscViewer viewer;
25 1 : PetscBool verbose;
26 :
27 1 : PetscFunctionBeginUser;
28 1 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
29 1 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
30 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type PEP - n=%" PetscInt_FMT ".\n",n));
31 1 : PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
32 :
33 : /* Create DS object */
34 1 : PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
35 1 : PetscCall(DSSetType(ds,DSPEP));
36 1 : PetscCall(DSSetFromOptions(ds));
37 1 : PetscCall(DSPEPSetDegree(ds,d));
38 :
39 : /* Set dimensions */
40 1 : ld = n+2; /* test leading dimension larger than n */
41 1 : PetscCall(DSAllocate(ds,ld));
42 1 : PetscCall(DSSetDimensions(ds,n,0,0));
43 :
44 : /* Set up viewer */
45 1 : PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
46 1 : PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
47 1 : PetscCall(DSView(ds,viewer));
48 1 : PetscCall(PetscViewerPopFormat(viewer));
49 1 : if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
50 :
51 : /* Fill matrices */
52 1 : PetscCall(DSGetArray(ds,DS_MAT_E0,&K));
53 10 : for (i=0;i<n-1;i++) K[i+i*ld] = 2.0*n;
54 1 : K[n-1+(n-1)*ld] = 1.0*n;
55 10 : for (i=1;i<n;i++) {
56 9 : K[i+(i-1)*ld] = -1.0*n;
57 9 : K[(i-1)+i*ld] = -1.0*n;
58 : }
59 1 : PetscCall(DSRestoreArray(ds,DS_MAT_E0,&K));
60 1 : PetscCall(DSGetArray(ds,DS_MAT_E1,&C));
61 1 : C[n-1+(n-1)*ld] = 2.0*PETSC_PI/z;
62 1 : PetscCall(DSRestoreArray(ds,DS_MAT_E1,&C));
63 1 : PetscCall(DSGetArray(ds,DS_MAT_E2,&M));
64 10 : for (i=0;i<n-1;i++) M[i+i*ld] = -4.0*PETSC_PI*PETSC_PI/n;
65 1 : M[i+i*ld] = -2.0*PETSC_PI*PETSC_PI/n;
66 1 : PetscCall(DSRestoreArray(ds,DS_MAT_E2,&M));
67 :
68 1 : if (verbose) {
69 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
70 0 : PetscCall(DSView(ds,viewer));
71 : }
72 :
73 : /* Solve */
74 1 : PetscCall(PetscMalloc2(d*n,&wr,d*n,&wi));
75 1 : PetscCall(DSGetSlepcSC(ds,&sc));
76 1 : sc->comparison = SlepcCompareLargestReal;
77 1 : sc->comparisonctx = NULL;
78 1 : sc->map = NULL;
79 1 : sc->mapobj = NULL;
80 1 : PetscCall(DSSolve(ds,wr,wi));
81 1 : PetscCall(DSSort(ds,wr,wi,NULL,NULL,NULL));
82 1 : if (verbose) {
83 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
84 0 : PetscCall(DSView(ds,viewer));
85 : }
86 :
87 : /* Print polynomial coefficients */
88 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Polynomial coefficients (alpha,beta,gamma) =\n"));
89 1 : PetscCall(DSPEPGetCoefficients(ds,&pbc));
90 4 : for (j=0;j<3;j++) {
91 12 : for (i=0;i<d+1;i++) PetscCall(PetscViewerASCIIPrintf(viewer," %.5f",(double)pbc[j+3*i]));
92 3 : PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
93 : }
94 1 : PetscCall(PetscFree(pbc));
95 :
96 : /* Print eigenvalues */
97 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n"));
98 21 : for (i=0;i<d*n;i++) {
99 : #if defined(PETSC_USE_COMPLEX)
100 20 : re = PetscRealPart(wr[i]);
101 20 : im = PetscImaginaryPart(wr[i]);
102 : #else
103 : re = wr[i];
104 : im = wi[i];
105 : #endif
106 20 : if (PetscAbs(im)<1e-10) PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re));
107 20 : else PetscCall(PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im));
108 : }
109 :
110 : /* Eigenvectors */
111 1 : PetscCall(DSVectors(ds,DS_MAT_X,NULL,NULL)); /* all eigenvectors */
112 1 : PetscCall(DSGetMat(ds,DS_MAT_X,&X));
113 1 : PetscCall(MatCreateVecs(X,NULL,&x0));
114 1 : PetscCall(MatGetColumnVector(X,x0,0));
115 1 : PetscCall(VecNorm(x0,NORM_2,&nrm));
116 1 : PetscCall(DSRestoreMat(ds,DS_MAT_X,&X));
117 1 : PetscCall(VecDestroy(&x0));
118 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st vector = %.3f\n",(double)nrm));
119 1 : if (verbose) {
120 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After vectors - - - - - - - - -\n"));
121 0 : PetscCall(DSView(ds,viewer));
122 : }
123 :
124 1 : PetscCall(PetscFree2(wr,wi));
125 1 : PetscCall(DSDestroy(&ds));
126 1 : PetscCall(SlepcFinalize());
127 : return 0;
128 : }
129 :
130 : /*TEST
131 :
132 : test:
133 : suffix: 1
134 : requires: !single
135 :
136 : TEST*/
|