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 with complex eigenvalues.\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,*M,*wr,*wi;
22 1 : PetscReal re,im,nrm;
23 1 : PetscInt i,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 8 : for (i=0;i<n;i++) K[i+i*ld] = 2.0;
54 7 : for (i=1;i<n;i++) {
55 6 : K[i+(i-1)*ld] = -1.0;
56 6 : K[(i-1)+i*ld] = -1.0;
57 : }
58 1 : PetscCall(DSRestoreArray(ds,DS_MAT_E0,&K));
59 1 : PetscCall(DSGetArray(ds,DS_MAT_E2,&M));
60 8 : for (i=0;i<n;i++) M[i+i*ld] = 1.0;
61 1 : PetscCall(DSRestoreArray(ds,DS_MAT_E2,&M));
62 :
63 1 : if (verbose) {
64 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
65 0 : PetscCall(DSView(ds,viewer));
66 : }
67 :
68 : /* Solve */
69 1 : PetscCall(PetscMalloc2(d*n,&wr,d*n,&wi));
70 1 : PetscCall(DSGetSlepcSC(ds,&sc));
71 1 : sc->comparison = SlepcCompareLargestImaginary;
72 1 : sc->comparisonctx = NULL;
73 1 : sc->map = NULL;
74 1 : sc->mapobj = NULL;
75 1 : PetscCall(DSSolve(ds,wr,wi));
76 1 : PetscCall(DSSort(ds,wr,wi,NULL,NULL,NULL));
77 1 : if (verbose) {
78 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
79 0 : PetscCall(DSView(ds,viewer));
80 : }
81 :
82 : /* Print eigenvalues */
83 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n"));
84 15 : for (i=0;i<d*n;i++) {
85 : #if defined(PETSC_USE_COMPLEX)
86 : re = PetscRealPart(wr[i]);
87 : im = PetscImaginaryPart(wr[i]);
88 : #else
89 14 : re = wr[i];
90 14 : im = wi[i];
91 : #endif
92 14 : if (PetscAbs(im)<1e-10) PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re));
93 14 : else PetscCall(PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im));
94 : }
95 :
96 : /* Eigenvectors */
97 1 : PetscCall(DSVectors(ds,DS_MAT_X,NULL,NULL)); /* all eigenvectors */
98 1 : PetscCall(DSGetMat(ds,DS_MAT_X,&X));
99 1 : PetscCall(MatCreateVecs(X,NULL,&x0));
100 1 : PetscCall(MatGetColumnVector(X,x0,1));
101 1 : PetscCall(VecNorm(x0,NORM_2,&nrm));
102 1 : PetscCall(DSRestoreMat(ds,DS_MAT_X,&X));
103 1 : PetscCall(VecDestroy(&x0));
104 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of 2nd column of X = %.3f\n",(double)nrm));
105 1 : if (verbose) {
106 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After vectors - - - - - - - - -\n"));
107 0 : PetscCall(DSView(ds,viewer));
108 : }
109 :
110 1 : PetscCall(PetscFree2(wr,wi));
111 1 : PetscCall(DSDestroy(&ds));
112 1 : PetscCall(SlepcFinalize());
113 : return 0;
114 : }
115 :
116 : /*TEST
117 :
118 : test:
119 : suffix: 1
120 : args: -n 7
121 : requires: !complex
122 : filter: sed -e 's/-0.00000/0.00000/'
123 :
124 : test:
125 : suffix: 1_complex
126 : args: -n 7
127 : requires: complex
128 : filter: sed -e 's/-0.00000/0.00000/'
129 :
130 : TEST*/
|