Actual source code: test2.c
slepc-3.18.2 2023-01-26
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 dense LME functions.\n\n";
13: #include <slepclme.h>
15: int main(int argc,char **argv)
16: {
17: LME lme;
18: Mat A,B,C,X;
19: PetscInt i,j,n=10,k=2;
20: PetscScalar *As,*Bs,*Cs,*Xs;
21: PetscViewer viewer;
22: PetscBool verbose;
25: SlepcInitialize(&argc,&argv,(char*)0,help);
26: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
27: PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
28: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
29: PetscPrintf(PETSC_COMM_WORLD,"Dense matrix equations, n=%" PetscInt_FMT ".\n",n);
31: /* Create LME object */
32: LMECreate(PETSC_COMM_WORLD,&lme);
34: /* Set up viewer */
35: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
36: if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
38: /* Create matrices */
39: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
40: PetscObjectSetName((PetscObject)A,"A");
41: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&B);
42: PetscObjectSetName((PetscObject)B,"B");
43: MatCreateSeqDense(PETSC_COMM_SELF,n,k,NULL,&C);
44: PetscObjectSetName((PetscObject)C,"C");
45: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&X);
46: PetscObjectSetName((PetscObject)X,"X");
48: /* Fill A with an upper Hessenberg Toeplitz matrix */
49: MatDenseGetArray(A,&As);
50: for (i=0;i<n;i++) As[i+i*n]=3.0-(PetscReal)n/2;
51: for (i=0;i<n-1;i++) As[i+1+i*n]=0.5;
52: for (j=1;j<3;j++) {
53: for (i=0;i<n-j;i++) As[i+(i+j)*n]=1.0;
54: }
55: MatDenseRestoreArray(A,&As);
57: if (verbose) {
58: PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
59: MatView(A,viewer);
60: }
62: /* Fill B with the 1-D Laplacian matrix */
63: MatDenseGetArray(B,&Bs);
64: for (i=0;i<n;i++) Bs[i+i*n]=2.0;
65: for (i=0;i<n-1;i++) { Bs[i+1+i*n]=-1; Bs[i+(i+1)*n]=-1; }
66: MatDenseRestoreArray(B,&Bs);
68: if (verbose) {
69: PetscPrintf(PETSC_COMM_WORLD,"Matrix B - - - - - - - -\n");
70: MatView(B,viewer);
71: }
73: /* Solve Lyapunov equation A*X+X*A'= -B */
74: PetscPrintf(PETSC_COMM_WORLD,"Solving Lyapunov equation for B\n");
75: MatDenseGetArray(A,&As);
76: MatDenseGetArray(B,&Bs);
77: MatDenseGetArray(X,&Xs);
78: LMEDenseLyapunov(lme,n,As,n,Bs,n,Xs,n);
79: MatDenseRestoreArray(A,&As);
80: MatDenseRestoreArray(B,&Bs);
81: MatDenseRestoreArray(X,&Xs);
82: if (verbose) {
83: PetscPrintf(PETSC_COMM_WORLD,"Solution X - - - - - - - -\n");
84: MatView(X,viewer);
85: }
87: /* Fill C with a full-rank nx2 matrix */
88: MatDenseGetArray(C,&Cs);
89: for (i=0;i<k;i++) Cs[i+i*n] = (i%2)? -1.0: 1.0;
90: MatDenseRestoreArray(C,&Cs);
92: if (verbose) {
93: PetscPrintf(PETSC_COMM_WORLD,"Matrix C - - - - - - - -\n");
94: MatView(C,viewer);
95: }
97: /* Solve Lyapunov equation A*X+X*A'= -C*C' */
98: PetscPrintf(PETSC_COMM_WORLD,"Solving Lyapunov equation for C (Cholesky)\n");
99: MatDenseGetArray(A,&As);
100: MatDenseGetArray(C,&Cs);
101: MatDenseGetArray(X,&Xs);
102: LMEDenseHessLyapunovChol(lme,n,As,n,2,Cs,n,Xs,n,NULL);
103: MatDenseRestoreArray(A,&As);
104: MatDenseRestoreArray(C,&Cs);
105: MatDenseRestoreArray(X,&Xs);
106: if (verbose) {
107: PetscPrintf(PETSC_COMM_WORLD,"Solution X - - - - - - - -\n");
108: MatView(X,viewer);
109: }
111: MatDestroy(&A);
112: MatDestroy(&B);
113: MatDestroy(&C);
114: MatDestroy(&X);
115: LMEDestroy(&lme);
116: SlepcFinalize();
117: return 0;
118: }
120: /*TEST
122: test:
123: args: -info :lme
124: requires: double
125: filter: sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/1e-\\1/g"
127: TEST*/