Actual source code: test7.c

slepc-3.8.2 2017-12-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, 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 multiplication of a Mat times a BV.\n\n";

 13: #include <slepcbv.h>

 15: int main(int argc,char **argv)
 16: {
 18:   Vec            t,v;
 19:   Mat            B;
 20:   BV             X,Y,Z,Zcopy;
 21:   PetscInt       i,j,n=10,k=5,Istart,Iend;
 22:   PetscScalar    *pZ;
 23:   PetscReal      norm;
 24:   PetscViewer    view;
 25:   PetscBool      verbose;
 26:   BVMatMultType  vmm;

 28:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 29:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 30:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 31:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 32:   PetscPrintf(PETSC_COMM_WORLD,"Test BVMatMult (n=%D, k=%D).\n",n,k);

 34:   /* Create Laplacian matrix */
 35:   MatCreate(PETSC_COMM_WORLD,&B);
 36:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
 37:   MatSetFromOptions(B);
 38:   MatSetUp(B);
 39:   PetscObjectSetName((PetscObject)B,"B");

 41:   MatGetOwnershipRange(B,&Istart,&Iend);
 42:   for (i=Istart;i<Iend;i++) {
 43:     if (i>0) { MatSetValue(B,i,i-1,-1.0,INSERT_VALUES); }
 44:     if (i<n-1) { MatSetValue(B,i,i+1,-1.0,INSERT_VALUES); }
 45:     MatSetValue(B,i,i,2.0,INSERT_VALUES);
 46:   }
 47:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 48:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 49:   MatCreateVecs(B,&t,NULL);

 51:   /* Create BV object X */
 52:   BVCreate(PETSC_COMM_WORLD,&X);
 53:   PetscObjectSetName((PetscObject)X,"X");
 54:   BVSetSizesFromVec(X,t,k);
 55:   BVSetMatMultMethod(X,BV_MATMULT_VECS);
 56:   BVSetFromOptions(X);
 57:   BVGetMatMultMethod(X,&vmm);
 58:   PetscPrintf(PETSC_COMM_WORLD,"Using method: %s\n",BVMatMultTypes[vmm]);

 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=Istart;i<PetscMin(j+1,Iend);i++) {
 71:       VecSetValue(v,i,1.0,INSERT_VALUES);
 72:     }
 73:     VecAssemblyBegin(v);
 74:     VecAssemblyEnd(v);
 75:     BVRestoreColumn(X,j,&v);
 76:   }
 77:   if (verbose) {
 78:     MatView(B,view);
 79:     BVView(X,view);
 80:   }

 82:   /* Create BV object Y */
 83:   BVDuplicateResize(X,k+4,&Y);
 84:   PetscObjectSetName((PetscObject)Y,"Y");
 85:   BVSetActiveColumns(Y,2,k+2);

 87:   /* Test BVMatMult */
 88:   BVMatMult(X,B,Y);
 89:   if (verbose) {
 90:     BVView(Y,view);
 91:   }

 93:   /* Create BV object Z */
 94:   BVDuplicate(X,&Z);
 95:   PetscObjectSetName((PetscObject)Z,"Z");

 97:   /* Fill Z entries */
 98:   for (j=0;j<k;j++) {
 99:     BVGetColumn(Z,j,&v);
100:     VecSet(v,0.0);
101:     if (!Istart) { VecSetValue(v,0,1.0,ADD_VALUES); }
102:     if (j<n && j>=Istart && j<Iend) { VecSetValue(v,j,1.0,ADD_VALUES); }
103:     if (j+1<n && j>=Istart && j<Iend) { VecSetValue(v,j+1,-1.0,ADD_VALUES); }
104:     VecAssemblyBegin(v);
105:     VecAssemblyEnd(v);
106:     BVRestoreColumn(Z,j,&v);
107:   }
108:   if (verbose) {
109:     BVView(Z,view);
110:   }

112:   /* Save a copy of Z */
113:   BVDuplicate(Z,&Zcopy);
114:   BVCopy(Z,Zcopy);

116:   /* Test BVMult, check result of previous operations */
117:   BVMult(Z,-1.0,1.0,Y,NULL);
118:   BVNorm(Z,NORM_FROBENIUS,&norm);
119:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error: %g\n",(double)norm);

121:   /* Test BVMatMultColumn, multiply Y(:,2), result in Y(:,3) */
122:   BVMatMultColumn(Y,B,2);
123:   if (verbose) {
124:     BVView(Y,view);
125:   }

127:   /* Test BVGetArray, modify Z to match Y */
128:   BVCopy(Zcopy,Z);
129:   BVGetArray(Z,&pZ);
130:   if (Istart==0) {
131:     if (Iend<3) SETERRQ(PETSC_COMM_WORLD,1,"First process must have at least 3 rows");
132:     pZ[Iend]   = 5.0;   /* modify 3 first entries of second column */
133:     pZ[Iend+1] = -4.0;
134:     pZ[Iend+2] = 1.0;
135:   }
136:   BVRestoreArray(Z,&pZ);
137:   if (verbose) {
138:     BVView(Z,view);
139:   }

141:   /* Check result again with BVMult */
142:   BVMult(Z,-1.0,1.0,Y,NULL);
143:   BVNorm(Z,NORM_FROBENIUS,&norm);
144:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error: %g\n",(double)norm);

146:   BVDestroy(&X);
147:   BVDestroy(&Y);
148:   BVDestroy(&Z);
149:   BVDestroy(&Zcopy);
150:   MatDestroy(&B);
151:   VecDestroy(&t);
152:   SlepcFinalize();
153:   return ierr;
154: }