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