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