Actual source code: test17.c
slepc-3.22.1 2024-10-28
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 with complex eigenvalues.\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,*M,*wr,*wi;
22: PetscReal re,im,nrm;
23: PetscInt i,n=10,d=2,ld;
24: PetscViewer viewer;
25: PetscBool verbose;
27: PetscFunctionBeginUser;
28: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
29: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
30: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type PEP - n=%" PetscInt_FMT ".\n",n));
31: PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
33: /* Create DS object */
34: PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
35: PetscCall(DSSetType(ds,DSPEP));
36: PetscCall(DSSetFromOptions(ds));
37: PetscCall(DSPEPSetDegree(ds,d));
39: /* Set dimensions */
40: ld = n+2; /* test leading dimension larger than n */
41: PetscCall(DSAllocate(ds,ld));
42: PetscCall(DSSetDimensions(ds,n,0,0));
44: /* Set up viewer */
45: PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
46: PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
47: PetscCall(DSView(ds,viewer));
48: PetscCall(PetscViewerPopFormat(viewer));
49: if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
51: /* Fill matrices */
52: PetscCall(DSGetArray(ds,DS_MAT_E0,&K));
53: for (i=0;i<n;i++) K[i+i*ld] = 2.0;
54: for (i=1;i<n;i++) {
55: K[i+(i-1)*ld] = -1.0;
56: K[(i-1)+i*ld] = -1.0;
57: }
58: PetscCall(DSRestoreArray(ds,DS_MAT_E0,&K));
59: PetscCall(DSGetArray(ds,DS_MAT_E2,&M));
60: for (i=0;i<n;i++) M[i+i*ld] = 1.0;
61: PetscCall(DSRestoreArray(ds,DS_MAT_E2,&M));
63: if (verbose) {
64: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
65: PetscCall(DSView(ds,viewer));
66: }
68: /* Solve */
69: PetscCall(PetscMalloc2(d*n,&wr,d*n,&wi));
70: PetscCall(DSGetSlepcSC(ds,&sc));
71: sc->comparison = SlepcCompareLargestImaginary;
72: sc->comparisonctx = NULL;
73: sc->map = NULL;
74: sc->mapobj = NULL;
75: PetscCall(DSSolve(ds,wr,wi));
76: PetscCall(DSSort(ds,wr,wi,NULL,NULL,NULL));
77: if (verbose) {
78: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
79: PetscCall(DSView(ds,viewer));
80: }
82: /* Print eigenvalues */
83: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n"));
84: 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: re = wr[i];
90: im = wi[i];
91: #endif
92: if (PetscAbs(im)<1e-10) PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re));
93: else PetscCall(PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im));
94: }
96: /* Eigenvectors */
97: PetscCall(DSVectors(ds,DS_MAT_X,NULL,NULL)); /* all eigenvectors */
98: PetscCall(DSGetMat(ds,DS_MAT_X,&X));
99: PetscCall(MatCreateVecs(X,NULL,&x0));
100: PetscCall(MatGetColumnVector(X,x0,1));
101: PetscCall(VecNorm(x0,NORM_2,&nrm));
102: PetscCall(DSRestoreMat(ds,DS_MAT_X,&X));
103: PetscCall(VecDestroy(&x0));
104: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of 2nd column of X = %.3f\n",(double)nrm));
105: if (verbose) {
106: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After vectors - - - - - - - - -\n"));
107: PetscCall(DSView(ds,viewer));
108: }
110: PetscCall(PetscFree2(wr,wi));
111: PetscCall(DSDestroy(&ds));
112: PetscCall(SlepcFinalize());
113: return 0;
114: }
116: /*TEST
118: test:
119: suffix: 1
120: args: -n 7
121: requires: !complex
122: filter: sed -e 's/-0.00000/0.00000/'
124: test:
125: suffix: 1_complex
126: args: -n 7
127: requires: complex
128: filter: sed -e 's/-0.00000/0.00000/'
130: TEST*/