Actual source code: test9.c

slepc-3.17.0 2022-03-31
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 BV matrix projection.\n\n";

 13: #include <slepcbv.h>

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

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

 37:   /* Set up viewer */
 38:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 39:   if (verbose) PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);

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

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

 58:   /* Create symmetric matrix B (1-D Laplacian) */
 59:   MatCreate(PETSC_COMM_WORLD,&B);
 60:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
 61:   MatSetFromOptions(B);
 62:   MatSetUp(B);
 63:   PetscObjectSetName((PetscObject)B,"B");

 65:   MatGetOwnershipRange(B,&Istart,&Iend);
 66:   for (i=Istart;i<Iend;i++) {
 67:     if (i>0) MatSetValue(B,i,i-1,-1.0,INSERT_VALUES);
 68:     if (i<n-1) MatSetValue(B,i,i+1,-1.0,INSERT_VALUES);
 69:     MatSetValue(B,i,i,2.0,INSERT_VALUES);
 70:   }
 71:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 72:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 73:   MatCreateVecs(B,&t,NULL);
 74:   if (verbose) MatView(B,view);

 76:   /* Create BV object X */
 77:   BVCreate(PETSC_COMM_WORLD,&X);
 78:   PetscObjectSetName((PetscObject)X,"X");
 79:   BVSetSizesFromVec(X,t,kx+2);  /* two extra columns to test active columns */
 80:   BVSetFromOptions(X);

 82:   /* Fill X entries */
 83:   for (j=0;j<kx+2;j++) {
 84:     BVGetColumn(X,j,&v);
 85:     VecSet(v,0.0);
 86:     for (i=0;i<4;i++) {
 87:       if (i+j<n) {
 88: #if defined(PETSC_USE_COMPLEX)
 89:         alpha = PetscCMPLX((PetscReal)(3*i+j-2),(PetscReal)(2*i));
 90: #else
 91:         alpha = (PetscReal)(3*i+j-2);
 92: #endif
 93:         VecSetValue(v,i+j,alpha,INSERT_VALUES);
 94:       }
 95:     }
 96:     VecAssemblyBegin(v);
 97:     VecAssemblyEnd(v);
 98:     BVRestoreColumn(X,j,&v);
 99:   }
100:   if (verbose) BVView(X,view);

102:   /* Duplicate BV object and store Z=G*X */
103:   BVDuplicate(X,&Z);
104:   PetscObjectSetName((PetscObject)Z,"Z");
105:   BVSetActiveColumns(X,0,kx);
106:   BVSetActiveColumns(Z,0,kx);
107:   BVMatMult(X,G,Z);
108:   BVSetActiveColumns(X,lx,kx);
109:   BVSetActiveColumns(Z,lx,kx);

111:   /* Create BV object Y */
112:   BVCreate(PETSC_COMM_WORLD,&Y);
113:   PetscObjectSetName((PetscObject)Y,"Y");
114:   BVSetSizesFromVec(Y,t,ky+1);
115:   BVSetFromOptions(Y);
116:   BVSetActiveColumns(Y,ly,ky);

118:   /* Fill Y entries */
119:   for (j=0;j<ky+1;j++) {
120:     BVGetColumn(Y,j,&v);
121: #if defined(PETSC_USE_COMPLEX)
122:     alpha = PetscCMPLX((PetscReal)(j+1)/4.0,-(PetscReal)j);
123: #else
124:     alpha = (PetscReal)(j+1)/4.0;
125: #endif
126:     VecSet(v,(PetscScalar)(j+1)/4.0);
127:     BVRestoreColumn(Y,j,&v);
128:   }
129:   if (verbose) BVView(Y,view);

131:   /* Test BVMatProject for non-symmetric matrix G */
132:   MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&H0);
133:   PetscObjectSetName((PetscObject)H0,"H0");
134:   BVMatProject(X,G,Y,H0);
135:   if (verbose) MatView(H0,view);

137:   /* Test BVMatProject with previously stored G*X */
138:   MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&H1);
139:   PetscObjectSetName((PetscObject)H1,"H1");
140:   BVMatProject(Z,NULL,Y,H1);
141:   if (verbose) MatView(H1,view);

143:   /* Check that H0 and H1 are equal */
144:   MatAXPY(H0,-1.0,H1,SAME_NONZERO_PATTERN);
145:   MatNorm(H0,NORM_1,&norm);
146:   if (norm<10*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"||H0-H1|| < 10*eps\n");
147:   else PetscPrintf(PETSC_COMM_WORLD,"||H0-H1||=%g\n",(double)norm);
148:   MatDestroy(&H0);
149:   MatDestroy(&H1);

151:   /* Test BVMatProject for symmetric matrix B with orthogonal projection */
152:   MatCreateSeqDense(PETSC_COMM_SELF,kx,kx,NULL,&H0);
153:   PetscObjectSetName((PetscObject)H0,"H0");
154:   BVMatProject(X,B,X,H0);
155:   if (verbose) MatView(H0,view);

157:   /* Repeat previous test with symmetry flag set */
158:   MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);
159:   MatCreateSeqDense(PETSC_COMM_SELF,kx,kx,NULL,&H1);
160:   PetscObjectSetName((PetscObject)H1,"H1");
161:   BVMatProject(X,B,X,H1);
162:   if (verbose) MatView(H1,view);

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

172:   BVDestroy(&X);
173:   BVDestroy(&Y);
174:   BVDestroy(&Z);
175:   MatDestroy(&B);
176:   MatDestroy(&G);
177:   VecDestroy(&t);
178:   SlepcFinalize();
179:   return 0;
180: }

182: /*TEST

184:    testset:
185:       output_file: output/test9_1.out
186:       test:
187:          suffix: 1
188:          args: -bv_type {{vecs contiguous svec mat}shared output}
189:       test:
190:          suffix: 1_svec_vecs
191:          args: -bv_type svec -bv_matmult vecs
192:       test:
193:          suffix: 1_cuda
194:          args: -bv_type svec -mat_type aijcusparse
195:          requires: cuda
196:       test:
197:          suffix: 2
198:          nsize: 2
199:          args: -bv_type {{vecs contiguous svec mat}shared output}
200:       test:
201:          suffix: 2_svec_vecs
202:          nsize: 2
203:          args: -bv_type svec -bv_matmult vecs

205: TEST*/