Actual source code: test24.c

slepc-3.22.1 2024-10-28
Report Typos and Errors
  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 DSGSVD with compact storage and rectangular matrix A.\n\n";

 13: #include <slepcds.h>

 15: int main(int argc,char **argv)
 16: {
 17:   DS             ds;
 18:   Mat            X;
 19:   Vec            x0;
 20:   SlepcSC        sc;
 21:   PetscReal      *T,*D,sigma,rnorm,aux,cond;
 22:   PetscScalar    *U,*V,*w,d;
 23:   PetscInt       i,n=10,m,l=0,k=0,ld;
 24:   PetscViewer    viewer;
 25:   PetscBool      verbose,extrarow;

 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 GSVD with compact storage - dimension %" PetscInt_FMT "x%" PetscInt_FMT ".\n",n+1,n));
 31:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL));
 32:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
 33:   PetscCheck(l<=n && k<=n && l<=k,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Wrong value of dimensions");
 34:   PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
 35:   PetscCall(PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow));
 36:   m = n+1;

 38:   /* Create DS object */
 39:   PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
 40:   PetscCall(DSSetType(ds,DSGSVD));
 41:   PetscCall(DSSetFromOptions(ds));
 42:   ld = n+2;  /* test leading dimension larger than n */
 43:   PetscCall(DSAllocate(ds,ld));
 44:   PetscCall(DSSetDimensions(ds,m,l,k));
 45:   PetscCall(DSGSVDSetDimensions(ds,n,n));
 46:   PetscCall(DSSetCompact(ds,PETSC_TRUE));
 47:   PetscCall(DSSetExtraRow(ds,extrarow));

 49:   /* Set up viewer */
 50:   PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
 51:   PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
 52:   PetscCall(DSView(ds,viewer));
 53:   PetscCall(PetscViewerPopFormat(viewer));

 55:   /* Fill A and B with lower/upper arrow-bidiagonal matrices
 56:      verifying that [A;B] has orthonormal columns */
 57:   PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T));
 58:   PetscCall(DSGetArrayReal(ds,DS_MAT_D,&D));
 59:   for (i=0;i<n;i++) T[i] = (PetscReal)(i+1)/(n+1); /* diagonal of matrix A */
 60:   for (i=0;i<k;i++) D[i] = PetscSqrtReal(1.0-T[i]*T[i]);
 61:   for (i=l;i<k;i++) {
 62:     T[i+ld] = PetscSqrtReal((1.0-T[k]*T[k])/(1.0+T[i]*T[i]/(D[i]*D[i])))*0.5*(1.0/k); /* upper diagonal of matrix A */
 63:     T[i+2*ld] = -T[i+ld]*T[i]/D[i]; /* upper diagonal of matrix B */
 64:   }
 65:   aux = 1.0-T[k]*T[k];
 66:   for (i=l;i<k;i++) aux -= T[i+ld]*T[i+ld]+T[i+2*ld]*T[i+2*ld];
 67:   T[k+ld] = PetscSqrtReal((1.0-aux)*.1);
 68:   aux -= T[k+ld]*T[k+ld];
 69:   D[k] = PetscSqrtReal(aux);
 70:   for (i=k+1;i<n;i++) {
 71:     T[i-1+2*ld] = -T[i-1+ld]*T[i]/D[i-1]; /* upper diagonal of matrix B */
 72:     aux = 1.0-T[i]*T[i]-T[2*ld+i-1]*T[2*ld+i-1];
 73:     T[i+ld] = PetscSqrtReal(aux)*.1; /* upper diagonal of matrix A */
 74:     D[i] = PetscSqrtReal(aux-T[i+ld]*T[i+ld]);
 75:   }
 76:   if (extrarow) { T[n]=-1.0; T[n-1+2*ld]=1.0; }
 77:   /* Fill locked eigenvalues */
 78:   PetscCall(PetscMalloc1(n,&w));
 79:   for (i=0;i<l;i++) w[i] = T[i]/D[i];
 80:   PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T));
 81:   PetscCall(DSRestoreArrayReal(ds,DS_MAT_D,&D));
 82:   if (l==0 && k==0) PetscCall(DSSetState(ds,DS_STATE_INTERMEDIATE));
 83:   else PetscCall(DSSetState(ds,DS_STATE_RAW));
 84:   if (verbose) {
 85:     PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
 86:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
 87:     PetscCall(DSView(ds,viewer));
 88:   }

 90:   /* Condition number */
 91:   PetscCall(DSCond(ds,&cond));
 92:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Condition number = %.3f\n",(double)cond));

 94:   /* Solve */
 95:   PetscCall(DSGetSlepcSC(ds,&sc));
 96:   sc->comparison    = SlepcCompareLargestReal;
 97:   sc->comparisonctx = NULL;
 98:   sc->map           = NULL;
 99:   sc->mapobj        = NULL;
100:   PetscCall(DSSolve(ds,w,NULL));
101:   PetscCall(DSSort(ds,w,NULL,NULL,NULL,NULL));
102:   if (extrarow) PetscCall(DSUpdateExtraRow(ds));
103:   PetscCall(DSSynchronize(ds,w,NULL));
104:   if (verbose) {
105:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
106:     PetscCall(DSView(ds,viewer));
107:   }

109:   /* Print singular values */
110:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed singular values =\n"));
111:   for (i=0;i<n;i++) {
112:     sigma = PetscRealPart(w[i]);
113:     PetscCall(PetscViewerASCIIPrintf(viewer,"  %.5f\n",(double)sigma));
114:   }

116:   if (extrarow) {
117:     /* Check that extra row is correct */
118:     PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T));
119:     PetscCall(DSGetArray(ds,DS_MAT_U,&U));
120:     PetscCall(DSGetArray(ds,DS_MAT_V,&V));
121:     d = 0.0;
122:     for (i=0;i<n;i++) d += T[i+ld]+U[n+i*ld];
123:     if (PetscAbsScalar(d)>10*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: there is a mismatch in A's extra row of %g\n",(double)PetscAbsScalar(d)));
124:     d = 0.0;
125:     for (i=0;i<n;i++) d += T[i+2*ld]-V[n-1+i*ld];
126:     if (PetscAbsScalar(d)>10*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: there is a mismatch in B's extra row of %g\n",(double)PetscAbsScalar(d)));
127:     PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T));
128:     PetscCall(DSRestoreArray(ds,DS_MAT_U,&U));
129:     PetscCall(DSRestoreArray(ds,DS_MAT_V,&V));
130:   }

132:   /* Singular vectors */
133:   PetscCall(DSVectors(ds,DS_MAT_X,NULL,NULL));  /* all singular vectors */
134:   PetscCall(DSGetMat(ds,DS_MAT_X,&X));
135:   PetscCall(MatCreateVecs(X,NULL,&x0));
136:   PetscCall(MatGetColumnVector(X,x0,0));
137:   PetscCall(VecNorm(x0,NORM_2,&rnorm));
138:   PetscCall(DSRestoreMat(ds,DS_MAT_X,&X));
139:   PetscCall(VecDestroy(&x0));
140:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st X vector = %.3f\n",(double)rnorm));

142:   PetscCall(DSGetMat(ds,DS_MAT_U,&X));
143:   PetscCall(MatCreateVecs(X,NULL,&x0));
144:   PetscCall(MatGetColumnVector(X,x0,0));
145:   PetscCall(VecNorm(x0,NORM_2,&rnorm));
146:   PetscCall(DSRestoreMat(ds,DS_MAT_U,&X));
147:   PetscCall(VecDestroy(&x0));
148:   if (PetscAbs(rnorm-1.0)>10*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: the 1st U vector has norm %g\n",(double)rnorm));

150:   PetscCall(DSGetMat(ds,DS_MAT_V,&X));
151:   PetscCall(MatCreateVecs(X,NULL,&x0));
152:   PetscCall(MatGetColumnVector(X,x0,0));
153:   PetscCall(VecNorm(x0,NORM_2,&rnorm));
154:   PetscCall(DSRestoreMat(ds,DS_MAT_V,&X));
155:   PetscCall(VecDestroy(&x0));
156:   if (PetscAbs(rnorm-1.0)>10*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: the 1st V vector has norm %g\n",(double)rnorm));

158:   PetscCall(PetscFree(w));
159:   PetscCall(DSDestroy(&ds));
160:   PetscCall(SlepcFinalize());
161:   return 0;
162: }

164: /*TEST

166:    testset:
167:       requires: double
168:       output_file: output/test24_1.out
169:       test:
170:          suffix: 1
171:          args: -l 1 -k 4
172:       test:
173:          suffix: 1_extrarow
174:          filter: sed -e "s/extrarow//"
175:          args: -l 1 -k 4 -extrarow

177: TEST*/