Actual source code: test9.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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 BV matrix projection.\n\n";

 13: #include <slepcbv.h>

 15: int main(int argc,char **argv)
 16: {
 18:   Vec            t,v;
 19:   Mat            B,G,H0,H1;
 20:   BV             X,Y,Z;
 21:   PetscInt       i,j,n=20,kx=6,lx=3,ky=5,ly=2,Istart,Iend,col[5];
 22:   PetscScalar    alpha,value[] = { -1, 1, 1, 1, 1 };
 23:   PetscViewer    view;
 24:   PetscReal      norm;
 25:   PetscBool      verbose;

 27:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 28:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 29:   PetscOptionsGetInt(NULL,NULL,"-kx",&kx,NULL);
 30:   PetscOptionsGetInt(NULL,NULL,"-lx",&lx,NULL);
 31:   PetscOptionsGetInt(NULL,NULL,"-ky",&ky,NULL);
 32:   PetscOptionsGetInt(NULL,NULL,"-ly",&ly,NULL);
 33:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 34:   PetscPrintf(PETSC_COMM_WORLD,"Test BV projection (n=%D).\n",n);
 35:   PetscPrintf(PETSC_COMM_WORLD,"X has %D active columns (%D leading columns).\n",kx,lx);
 36:   PetscPrintf(PETSC_COMM_WORLD,"Y has %D active columns (%D leading columns).\n",ky,ly);

 38:   /* Set up viewer */
 39:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 40:   if (verbose) {
 41:     PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
 42:   }

 44:   /* Create non-symmetric matrix G (Toeplitz) */
 45:   MatCreate(PETSC_COMM_WORLD,&G);
 46:   MatSetSizes(G,PETSC_DECIDE,PETSC_DECIDE,n,n);
 47:   MatSetFromOptions(G);
 48:   MatSetUp(G);
 49:   PetscObjectSetName((PetscObject)G,"G");

 51:   MatGetOwnershipRange(G,&Istart,&Iend);
 52:   for (i=Istart;i<Iend;i++) {
 53:     col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
 54:     if (i==0) {
 55:       MatSetValues(G,1,&i,PetscMin(4,n-i),col+1,value+1,INSERT_VALUES);
 56:     } else {
 57:       MatSetValues(G,1,&i,PetscMin(5,n-i+1),col,value,INSERT_VALUES);
 58:     }
 59:   }
 60:   MatAssemblyBegin(G,MAT_FINAL_ASSEMBLY);
 61:   MatAssemblyEnd(G,MAT_FINAL_ASSEMBLY);
 62:   if (verbose) {
 63:     MatView(G,view);
 64:   }

 66:   /* Create symmetric matrix B (1-D Laplacian) */
 67:   MatCreate(PETSC_COMM_WORLD,&B);
 68:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
 69:   MatSetFromOptions(B);
 70:   MatSetUp(B);
 71:   PetscObjectSetName((PetscObject)B,"B");

 73:   MatGetOwnershipRange(B,&Istart,&Iend);
 74:   for (i=Istart;i<Iend;i++) {
 75:     if (i>0) { MatSetValue(B,i,i-1,-1.0,INSERT_VALUES); }
 76:     if (i<n-1) { MatSetValue(B,i,i+1,-1.0,INSERT_VALUES); }
 77:     MatSetValue(B,i,i,2.0,INSERT_VALUES);
 78:   }
 79:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 80:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 81:   MatCreateVecs(B,&t,NULL);
 82:   if (verbose) {
 83:     MatView(B,view);
 84:   }

 86:   /* Create BV object X */
 87:   BVCreate(PETSC_COMM_WORLD,&X);
 88:   PetscObjectSetName((PetscObject)X,"X");
 89:   BVSetSizesFromVec(X,t,kx+2);  /* two extra columns to test active columns */
 90:   BVSetFromOptions(X);

 92:   /* Fill X entries */
 93:   for (j=0;j<kx+2;j++) {
 94:     BVGetColumn(X,j,&v);
 95:     VecSet(v,0.0);
 96:     for (i=0;i<4;i++) {
 97:       if (i+j<n) {
 98: #if defined(PETSC_USE_COMPLEX)
 99:         alpha = (PetscReal)(3*i+j-2)+(PetscReal)(2*i)*PETSC_i;
100: #else
101:         alpha = (PetscReal)(3*i+j-2);
102: #endif
103:         VecSetValue(v,i+j,alpha,INSERT_VALUES);
104:       }
105:     }
106:     VecAssemblyBegin(v);
107:     VecAssemblyEnd(v);
108:     BVRestoreColumn(X,j,&v);
109:   }
110:   if (verbose) {
111:     BVView(X,view);
112:   }

114:   /* Duplicate BV object and store Z=G*X */
115:   BVDuplicate(X,&Z);
116:   PetscObjectSetName((PetscObject)Z,"Z");
117:   BVSetActiveColumns(X,0,kx);
118:   BVSetActiveColumns(Z,0,kx);
119:   BVMatMult(X,G,Z);
120:   BVSetActiveColumns(X,lx,kx);
121:   BVSetActiveColumns(Z,lx,kx);

123:   /* Create BV object Y */
124:   BVCreate(PETSC_COMM_WORLD,&Y);
125:   PetscObjectSetName((PetscObject)Y,"Y");
126:   BVSetSizesFromVec(Y,t,ky+1);
127:   BVSetFromOptions(Y);
128:   BVSetActiveColumns(Y,ly,ky);

130:   /* Fill Y entries */
131:   for (j=0;j<ky+1;j++) {
132:     BVGetColumn(Y,j,&v);
133: #if defined(PETSC_USE_COMPLEX)
134:     alpha = (PetscReal)(j+1)/4.0-(PetscReal)j*PETSC_i;
135: #else
136:     alpha = (PetscReal)(j+1)/4.0;
137: #endif
138:     VecSet(v,(PetscScalar)(j+1)/4.0);
139:     BVRestoreColumn(Y,j,&v);
140:   }
141:   if (verbose) {
142:     BVView(Y,view);
143:   }

145:   /* Test BVMatProject for non-symmetric matrix G */
146:   MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&H0);
147:   PetscObjectSetName((PetscObject)H0,"H0");
148:   BVMatProject(X,G,Y,H0);
149:   if (verbose) {
150:     MatView(H0,view);
151:   }

153:   /* Test BVMatProject with previously stored G*X */
154:   MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&H1);
155:   PetscObjectSetName((PetscObject)H1,"H1");
156:   BVMatProject(Z,NULL,Y,H1);
157:   if (verbose) {
158:     MatView(H1,view);
159:   }

161:   /* Check that H0 and H1 are equal */
162:   MatAXPY(H0,-1.0,H1,SAME_NONZERO_PATTERN);
163:   MatNorm(H0,NORM_1,&norm);
164:   if (norm<10*PETSC_MACHINE_EPSILON) {
165:     PetscPrintf(PETSC_COMM_WORLD,"||H0-H1|| < 10*eps\n");
166:   } else {
167:     PetscPrintf(PETSC_COMM_WORLD,"||H0-H1||=%g\n",(double)norm);
168:   }
169:   MatDestroy(&H0);
170:   MatDestroy(&H1);

172:   /* Test BVMatProject for symmetric matrix B with orthogonal projection */
173:   MatCreateSeqDense(PETSC_COMM_SELF,kx,kx,NULL,&H0);
174:   PetscObjectSetName((PetscObject)H0,"H0");
175:   BVMatProject(X,B,X,H0);
176:   if (verbose) {
177:     MatView(H0,view);
178:   }

180:   /* Repeat previous test with symmetry flag set */
181:   MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);
182:   MatCreateSeqDense(PETSC_COMM_SELF,kx,kx,NULL,&H1);
183:   PetscObjectSetName((PetscObject)H1,"H1");
184:   BVMatProject(X,B,X,H1);
185:   if (verbose) {
186:     MatView(H1,view);
187:   }

189:   /* Check that H0 and H1 are equal */
190:   MatAXPY(H0,-1.0,H1,SAME_NONZERO_PATTERN);
191:   MatNorm(H0,NORM_1,&norm);
192:   if (norm<10*PETSC_MACHINE_EPSILON) {
193:     PetscPrintf(PETSC_COMM_WORLD,"||H0-H1|| < 10*eps\n");
194:   } else {
195:     PetscPrintf(PETSC_COMM_WORLD,"||H0-H1||=%g\n",(double)norm);
196:   }
197:   MatDestroy(&H0);
198:   MatDestroy(&H1);

200:   BVDestroy(&X);
201:   BVDestroy(&Y);
202:   BVDestroy(&Z);
203:   MatDestroy(&B);
204:   MatDestroy(&G);
205:   VecDestroy(&t);
206:   SlepcFinalize();
207:   return ierr;
208: }

210: /*TEST

212:    test:
213:       suffix: 1
214:       nsize: 1
215:       args: -bv_type {{vecs contiguous svec mat}shared output}
216:       output_file: output/test9_1.out

218:    test:
219:       suffix: 1_svec_vecs
220:       nsize: 1
221:       args: -bv_type svec -bv_matmult vecs
222:       output_file: output/test9_1.out

224:    test:
225:       suffix: 1_cuda
226:       nsize: 1
227:       args: -bv_type svec -mat_type aijcusparse
228:       requires: cuda
229:       output_file: output/test9_1.out

231:    test:
232:       suffix: 2
233:       nsize: 2
234:       args: -bv_type {{vecs contiguous svec mat}shared output}
235:       output_file: output/test9_1.out

237:    test:
238:       suffix: 2_svec_vecs
239:       nsize: 2
240:       args: -bv_type svec -bv_matmult vecs
241:       output_file: output/test9_1.out

243: TEST*/