Actual source code: ex45.c

slepc-3.15.2 2021-09-20
Report Typos and Errors
  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[] = "Computes a partial generalized singular value decomposition (GSVD).\n"
 12:   "The command line options are:\n"
 13:   "  -m <m>, where <m> = number of rows of A.\n"
 14:   "  -n <n>, where <n> = number of columns of A.\n"
 15:   "  -p <p>, where <p> = number of rows of B.\n\n";

 17: #include <slepcsvd.h>

 19: int main(int argc,char **argv)
 20: {
 21:   Mat            A,B;             /* operator matrices */
 22:   Vec            u,v,x;           /* singular vectors */
 23:   SVD            svd;             /* singular value problem solver context */
 24:   SVDType        type;
 25:   Vec            uv,aux[2];
 26:   PetscReal      error,tol,sigma;
 27:   PetscInt       m=100,n,p=14,i,j,Istart,Iend,nsv,maxit,its,nconv;
 28:   PetscBool      flg;

 31:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 33:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 34:   PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg);
 35:   if (!flg) n = m;
 36:   PetscOptionsGetInt(NULL,NULL,"-p",&p,&flg);
 37:   PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized singular value decomposition, (%D+%D)x%D\n\n",m,p,n);

 39:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 40:                           Build the matrices
 41:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 43:   MatCreate(PETSC_COMM_WORLD,&A);
 44:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n);
 45:   MatSetFromOptions(A);
 46:   MatSetUp(A);

 48:   MatGetOwnershipRange(A,&Istart,&Iend);
 49:   for (i=Istart;i<Iend;i++) {
 50:     if (i>0 && i-1<n) { MatSetValue(A,i,i-1,-1.0,INSERT_VALUES); }
 51:     if (i+1<n) { MatSetValue(A,i,i+1,-1.0,INSERT_VALUES); }
 52:     if (i<n) { MatSetValue(A,i,i,2.0,INSERT_VALUES); }
 53:     if (i>n) { MatSetValue(A,i,n-1,1.0,INSERT_VALUES); }
 54:   }
 55:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 56:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 58:   MatCreate(PETSC_COMM_WORLD,&B);
 59:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,p,n);
 60:   MatSetFromOptions(B);
 61:   MatSetUp(B);

 63:   MatGetOwnershipRange(B,&Istart,&Iend);
 64:   for (i=Istart;i<Iend;i++) {
 65:     for (j=0;j<=PetscMin(i,n-1);j++) {
 66:       MatSetValue(B,i,j,1.0,INSERT_VALUES);
 67:     }
 68:   }
 69:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 70:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 72:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 73:           Create the singular value solver and set various options
 74:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 76:   /*
 77:      Create singular value solver context
 78:   */
 79:   SVDCreate(PETSC_COMM_WORLD,&svd);

 81:   /*
 82:      Set operators and problem type
 83:   */
 84:   SVDSetOperators(svd,A,B);
 85:   SVDSetProblemType(svd,SVD_GENERALIZED);

 87:   /*
 88:      Set solver parameters at runtime
 89:   */
 90:   SVDSetFromOptions(svd);

 92:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 93:                       Solve the singular value system
 94:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 96:   SVDSolve(svd);
 97:   SVDGetIterationNumber(svd,&its);
 98:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);

100:   /*
101:      Optional: Get some information from the solver and display it
102:   */
103:   SVDGetType(svd,&type);
104:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
105:   SVDGetDimensions(svd,&nsv,NULL,NULL);
106:   PetscPrintf(PETSC_COMM_WORLD," Number of requested singular values: %D\n",nsv);
107:   SVDGetTolerances(svd,&tol,&maxit);
108:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%D\n",(double)tol,maxit);

110:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111:                     Display solution and clean up
112:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

114:   /*
115:      Get number of converged singular triplets
116:   */
117:   SVDGetConverged(svd,&nconv);
118:   PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate singular triplets: %D\n\n",nconv);

120:   if (nconv>0) {
121:     /*
122:        Create vectors. The interface returns u and v as stacked on top of each other
123:        [u;v] so need to create a special vector (VecNest) to extract them
124:     */
125:     MatCreateVecs(A,&x,&u);
126:     MatCreateVecs(B,NULL,&v);
127:     aux[0] = u;
128:     aux[1] = v;
129:     VecCreateNest(PETSC_COMM_WORLD,2,NULL,aux,&uv);

131:     /*
132:        Display singular values and relative errors
133:     */
134:     PetscPrintf(PETSC_COMM_WORLD,
135:          "          sigma           relative error\n"
136:          "  --------------------- ------------------\n");
137:     for (i=0;i<nconv;i++) {
138:       /*
139:          Get converged singular triplets: i-th singular value is stored in sigma
140:       */
141:       SVDGetSingularTriplet(svd,i,&sigma,uv,x);
142:       /* at this point, u and v can be used normally as individual vectors */

144:       /*
145:          Compute the error associated to each singular triplet
146:       */
147:       SVDComputeError(svd,i,SVD_ERROR_RELATIVE,&error);

149:       PetscPrintf(PETSC_COMM_WORLD,"       % 6f      ",(double)sigma);
150:       PetscPrintf(PETSC_COMM_WORLD,"   % 12g\n",(double)error);
151:     }
152:     PetscPrintf(PETSC_COMM_WORLD,"\n");
153:     VecDestroy(&x);
154:     VecDestroy(&u);
155:     VecDestroy(&v);
156:     VecDestroy(&uv);
157:   }

159:   /*
160:      Free work space
161:   */
162:   SVDDestroy(&svd);
163:   MatDestroy(&A);
164:   MatDestroy(&B);
165:   SlepcFinalize();
166:   return ierr;
167: }

169: /*TEST

171:    test:
172:       args: -m 20 -n 10 -p 6 -svd_type lapack
173:       suffix: 1
174:       filter: sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
175:       requires: double

177:    test:
178:       args: -m 15 -n 20 -p 20 -svd_type lapack
179:       suffix: 2
180:       filter: sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
181:       requires: double

183: TEST*/