Actual source code: test7.c

slepc-3.7.4 2017-05-17
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 multiplication of a Mat times a BV.\n\n";

 24: #include <slepcbv.h>

 28: int main(int argc,char **argv)
 29: {
 31:   Vec            t,v;
 32:   Mat            B;
 33:   BV             X,Y,Z,Zcopy;
 34:   PetscInt       i,j,n=10,k=5,Istart,Iend;
 35:   PetscScalar    *pZ;
 36:   PetscReal      norm;
 37:   PetscViewer    view;
 38:   PetscBool      verbose;

 40:   SlepcInitialize(&argc,&argv,(char*)0,help);
 41:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 42:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 43:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 44:   PetscPrintf(PETSC_COMM_WORLD,"Test BVMatMult (n=%D, k=%D).\n",n,k);

 46:   /* Create Laplacian matrix */
 47:   MatCreate(PETSC_COMM_WORLD,&B);
 48:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
 49:   MatSetFromOptions(B);
 50:   MatSetUp(B);
 51:   PetscObjectSetName((PetscObject)B,"B");

 53:   MatGetOwnershipRange(B,&Istart,&Iend);
 54:   for (i=Istart;i<Iend;i++) {
 55:     if (i>0) { MatSetValue(B,i,i-1,-1.0,INSERT_VALUES); }
 56:     if (i<n-1) { MatSetValue(B,i,i+1,-1.0,INSERT_VALUES); }
 57:     MatSetValue(B,i,i,2.0,INSERT_VALUES);
 58:   }
 59:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 60:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 61:   MatCreateVecs(B,&t,NULL);

 63:   /* Create BV object X */
 64:   BVCreate(PETSC_COMM_WORLD,&X);
 65:   PetscObjectSetName((PetscObject)X,"X");
 66:   BVSetSizesFromVec(X,t,k);
 67:   BVSetFromOptions(X);

 69:   /* Set up viewer */
 70:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 71:   if (verbose) {
 72:     PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
 73:   }

 75:   /* Fill X entries */
 76:   for (j=0;j<k;j++) {
 77:     BVGetColumn(X,j,&v);
 78:     VecSet(v,0.0);
 79:     for (i=Istart;i<PetscMin(j+1,Iend);i++) {
 80:       VecSetValue(v,i,1.0,INSERT_VALUES);
 81:     }
 82:     VecAssemblyBegin(v);
 83:     VecAssemblyEnd(v);
 84:     BVRestoreColumn(X,j,&v);
 85:   }
 86:   if (verbose) {
 87:     MatView(B,view);
 88:     BVView(X,view);
 89:   }

 91:   /* Create BV object Y */
 92:   BVDuplicateResize(X,k+4,&Y);
 93:   PetscObjectSetName((PetscObject)Y,"Y");
 94:   BVSetActiveColumns(Y,2,k+2);

 96:   /* Test BVMatMult */
 97:   BVMatMult(X,B,Y);
 98:   if (verbose) {
 99:     BVView(Y,view);
100:   }

102:   /* Create BV object Z */
103:   BVDuplicate(X,&Z);
104:   PetscObjectSetName((PetscObject)Z,"Z");

106:   /* Fill Z entries */
107:   for (j=0;j<k;j++) {
108:     BVGetColumn(Z,j,&v);
109:     VecSet(v,0.0);
110:     if (!Istart) { VecSetValue(v,0,1.0,ADD_VALUES); }
111:     if (j<n && j>=Istart && j<Iend) { VecSetValue(v,j,1.0,ADD_VALUES); }
112:     if (j+1<n && j>=Istart && j<Iend) { VecSetValue(v,j+1,-1.0,ADD_VALUES); }
113:     VecAssemblyBegin(v);
114:     VecAssemblyEnd(v);
115:     BVRestoreColumn(Z,j,&v);
116:   }
117:   if (verbose) {
118:     BVView(Z,view);
119:   }

121:   /* Save a copy of Z */
122:   BVDuplicate(Z,&Zcopy);
123:   BVCopy(Z,Zcopy);

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

130:   /* Test BVMatMultColumn, multiply Y(:,2), result in Y(:,3) */
131:   BVMatMultColumn(Y,B,2);
132:   if (verbose) {
133:     BVView(Y,view);
134:   }

136:   /* Test BVGetArray, modify Z to match Y */
137:   BVCopy(Zcopy,Z);
138:   BVGetArray(Z,&pZ);
139:   if (Istart==0) {
140:     if (Iend<3) SETERRQ(PETSC_COMM_SELF,1,"First process must have at least 3 rows");
141:     pZ[Iend]   = 5.0;   /* modify 3 first entries of second column */
142:     pZ[Iend+1] = -4.0;
143:     pZ[Iend+2] = 1.0;
144:   }
145:   BVRestoreArray(Z,&pZ);
146:   if (verbose) {
147:     BVView(Z,view);
148:   }

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

155:   BVDestroy(&X);
156:   BVDestroy(&Y);
157:   BVDestroy(&Z);
158:   BVDestroy(&Zcopy);
159:   MatDestroy(&B);
160:   VecDestroy(&t);
161:   SlepcFinalize();
162:   return ierr;
163: }