Actual source code: test11.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: static char help[] = "Test BV block orthogonalization.\n\n";

 24: #include <slepcbv.h>

 28: int main(int argc,char **argv)
 29: {
 30:   PetscErrorCode    ierr;
 31:   BV                X,Y,Z,cached;
 32:   Mat               B,M;
 33:   Vec               v,t;
 34:   PetscInt          i,j,n=20,l=2,k=8,Istart,Iend;
 35:   PetscViewer       view;
 36:   PetscBool         verbose;
 37:   PetscReal         norm;
 38:   PetscScalar       alpha;
 39:   BVOrthogBlockType btype;

 41:   SlepcInitialize(&argc,&argv,(char*)0,help);
 42:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 43:   PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
 44:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 45:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 46:   PetscPrintf(PETSC_COMM_WORLD,"Test BV block orthogonalization (length %D, l=%D, k=%D).\n",n,l,k);

 48:   /* Create template vector */
 49:   VecCreate(PETSC_COMM_WORLD,&t);
 50:   VecSetSizes(t,PETSC_DECIDE,n);
 51:   VecSetFromOptions(t);

 53:   /* Create BV object X */
 54:   BVCreate(PETSC_COMM_WORLD,&X);
 55:   PetscObjectSetName((PetscObject)X,"X");
 56:   BVSetSizesFromVec(X,t,k);
 57:   BVSetFromOptions(X);
 58:   BVGetOrthogonalization(X,NULL,NULL,NULL,&btype);

 60:   /* Set up viewer */
 61:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 62:   if (verbose) {
 63:     PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
 64:   }

 66:   /* Fill X entries */
 67:   for (j=0;j<k;j++) {
 68:     BVGetColumn(X,j,&v);
 69:     VecSet(v,0.0);
 70:     for (i=0;i<=n/2;i++) {
 71:       if (i+j<n) {
 72:         alpha = (3.0*i+j-2)/(2*(i+j+1));
 73:         VecSetValue(v,i+j,alpha,INSERT_VALUES);
 74:       }
 75:     }
 76:     VecAssemblyBegin(v);
 77:     VecAssemblyEnd(v);
 78:     BVRestoreColumn(X,j,&v);
 79:   }
 80:   if (btype==BV_ORTHOG_BLOCK_GS) {  /* GS requires the leading columns to be orthogonal */
 81:     for (j=0;j<l;j++) {
 82:       BVOrthogonalizeColumn(X,j,NULL,&norm,NULL);
 83:       alpha = 1.0/norm;
 84:       BVScaleColumn(X,j,alpha);
 85:     }
 86:   }
 87:   if (verbose) {
 88:     BVView(X,view);
 89:   }

 91:   /* Create copy on Y */
 92:   BVDuplicate(X,&Y);
 93:   PetscObjectSetName((PetscObject)Y,"Y");
 94:   BVCopy(X,Y);
 95:   BVSetActiveColumns(Y,l,k);
 96:   BVSetActiveColumns(X,l,k);

 98:   /* Test BVOrthogonalize */
 99:   BVOrthogonalize(Y,NULL);
100:   if (verbose) {
101:     BVView(Y,view);
102:   }

104:   /* Check orthogonality */
105:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
106:   MatShift(M,1.0);   /* set leading part to identity */
107:   BVDot(Y,Y,M);
108:   MatShift(M,-1.0);
109:   MatNorm(M,NORM_1,&norm);
110:   if (norm<100*PETSC_MACHINE_EPSILON) {
111:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
112:   } else {
113:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm);
114:   }

116:   /* Create inner product matrix */
117:   MatCreate(PETSC_COMM_WORLD,&B);
118:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
119:   MatSetFromOptions(B);
120:   MatSetUp(B);
121:   PetscObjectSetName((PetscObject)B,"B");

123:   MatGetOwnershipRange(B,&Istart,&Iend);
124:   for (i=Istart;i<Iend;i++) {
125:     if (i>0) { MatSetValue(B,i,i-1,-1.0,INSERT_VALUES); }
126:     if (i<n-1) { MatSetValue(B,i,i+1,-1.0,INSERT_VALUES); }
127:     MatSetValue(B,i,i,2.0,INSERT_VALUES);
128:   }
129:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
130:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

132:   /* Prepare to repeat test, now with a non-standard inner product */
133:   BVCopy(X,Y);
134:   BVDuplicate(X,&Z);
135:   PetscObjectSetName((PetscObject)Z,"Z");
136:   BVSetActiveColumns(Z,l,k);
137:   BVSetMatrix(X,B,PETSC_FALSE);
138:   BVSetMatrix(Y,B,PETSC_FALSE);
139:   if (btype==BV_ORTHOG_BLOCK_GS) {  /* GS requires the leading columns to be orthogonal */
140:     for (j=0;j<l;j++) {
141:       BVOrthogonalizeColumn(Y,j,NULL,&norm,NULL);
142:       alpha = 1.0/norm;
143:       BVScaleColumn(Y,j,alpha);
144:     }
145:   }
146:   if (verbose) {
147:     BVView(X,view);
148:   }

150:   /* Test BVOrthogonalize */
151:   BVOrthogonalize(Y,NULL);
152:   if (verbose) {
153:     BVView(Y,view);
154:   }

156:   /* Extract cached BV and check it is equal to B*X */
157:   BVGetCachedBV(Y,&cached);
158:   BVMatMult(X,B,Z);
159:   BVMult(Z,-1.0,1.0,cached,NULL);
160:   BVNorm(Z,NORM_FROBENIUS,&norm);
161:   if (norm<100*PETSC_MACHINE_EPSILON) {
162:     PetscPrintf(PETSC_COMM_WORLD,"Residual ||cached-BX|| < 100*eps\n");
163:   } else {
164:     PetscPrintf(PETSC_COMM_WORLD,"Residual ||cached-BX||: %g\n",(double)norm);
165:   }

167:   /* Check orthogonality */
168:   MatZeroEntries(M);
169:   MatShift(M,1.0);   /* set leading part to identity */
170:   BVDot(Y,Y,M);
171:   MatShift(M,-1.0);
172:   MatNorm(M,NORM_1,&norm);
173:   if (norm<100*PETSC_MACHINE_EPSILON) {
174:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
175:   } else {
176:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm);
177:   }

179:   MatDestroy(&M);
180:   MatDestroy(&B);
181:   BVDestroy(&X);
182:   BVDestroy(&Y);
183:   BVDestroy(&Z);
184:   VecDestroy(&t);
185:   SlepcFinalize();
186:   return ierr;
187: }