Actual source code: test7.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2014, 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,col[3];
 35:   PetscScalar    value[3],*pZ;
 36:   PetscReal      norm;
 37:   PetscViewer    view;
 38:   PetscBool      verbose,FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;

 40:   SlepcInitialize(&argc,&argv,(char*)0,help);
 41:   PetscOptionsGetInt(NULL,"-n",&n,NULL);
 42:   PetscOptionsGetInt(NULL,"-k",&k,NULL);
 43:   PetscOptionsHasName(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:   if (Istart==0) FirstBlock=PETSC_TRUE;
 55:   if (Iend==n) LastBlock=PETSC_TRUE;
 56:   value[0]=-1.0; value[1]=2.0; value[2]=-1.0;
 57:   for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
 58:     col[0]=i-1; col[1]=i; col[2]=i+1;
 59:     MatSetValues(B,1,&i,3,col,value,INSERT_VALUES);
 60:   }
 61:   if (LastBlock) {
 62:     i=n-1; col[0]=n-2; col[1]=n-1;
 63:     MatSetValues(B,1,&i,2,col,value,INSERT_VALUES);
 64:   }
 65:   if (FirstBlock) {
 66:     i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
 67:     MatSetValues(B,1,&i,2,col,value,INSERT_VALUES);
 68:   }
 69:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 70:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 71:   MatGetVecs(B,&t,NULL);

 73:   /* Create BV object X */
 74:   BVCreate(PETSC_COMM_WORLD,&X);
 75:   PetscObjectSetName((PetscObject)X,"X");
 76:   BVSetSizesFromVec(X,t,k);
 77:   BVSetFromOptions(X);

 79:   /* Set up viewer */
 80:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 81:   if (verbose) {
 82:     PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
 83:   }

 85:   /* Fill X entries */
 86:   for (j=0;j<k;j++) {
 87:     BVGetColumn(X,j,&v);
 88:     VecZeroEntries(v);
 89:     for (i=Istart;i<PetscMin(j+1,Iend);i++) {
 90:       VecSetValue(v,i,1.0,INSERT_VALUES);
 91:     }
 92:     VecAssemblyBegin(v);
 93:     VecAssemblyEnd(v);
 94:     BVRestoreColumn(X,j,&v);
 95:   }
 96:   if (verbose) {
 97:     MatView(B,view);
 98:     BVView(X,view);
 99:   }

101:   /* Create BV object Y */
102:   BVDuplicateResize(X,k+4,&Y);
103:   PetscObjectSetName((PetscObject)Y,"Y");
104:   BVSetActiveColumns(Y,2,k+2);

106:   /* Test BVMatMult */
107:   BVMatMult(X,B,Y);
108:   if (verbose) {
109:     BVView(Y,view);
110:   }

112:   /* Create BV object Z */
113:   BVDuplicate(X,&Z);
114:   PetscObjectSetName((PetscObject)Z,"Z");

116:   /* Fill Z entries */
117:   for (j=0;j<k;j++) {
118:     BVGetColumn(Z,j,&v);
119:     VecZeroEntries(v);
120:     if (!Istart) { VecSetValue(v,0,1.0,ADD_VALUES); }
121:     if (j<n && j>=Istart && j<Iend) { VecSetValue(v,j,1.0,ADD_VALUES); }
122:     if (j+1<n && j>=Istart && j<Iend) { VecSetValue(v,j+1,-1.0,ADD_VALUES); }
123:     VecAssemblyBegin(v);
124:     VecAssemblyEnd(v);
125:     BVRestoreColumn(Z,j,&v);
126:   }
127:   if (verbose) {
128:     BVView(Z,view);
129:   }

131:   /* Save a copy of Z */
132:   BVDuplicate(Z,&Zcopy);
133:   BVCopy(Z,Zcopy);

135:   /* Test BVAXPY, check result of previous operations */
136:   BVAXPY(Z,-1.0,Y);
137:   BVNorm(Z,NORM_FROBENIUS,&norm);
138:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error: %g\n",(double)norm);

140:   /* Test BVMatMultColumn, multiply Y(:,2), result in Y(:,3) */
141:   BVMatMultColumn(Y,B,2);
142:   if (verbose) {
143:     BVView(Y,view);
144:   }

146:   /* Test BVGetArray, modify Z to match Y */
147:   BVCopy(Zcopy,Z);
148:   BVGetArray(Z,&pZ);
149:   if (Istart==0) {
150:     if (Iend<3) SETERRQ(PETSC_COMM_SELF,1,"First process must have at least 3 rows");
151:     pZ[Iend]   = 5.0;   /* modify 3 first entries of second column */
152:     pZ[Iend+1] = -4.0;
153:     pZ[Iend+2] = 1.0;
154:   }
155:   BVRestoreArray(Z,&pZ);
156:   if (verbose) {
157:     BVView(Z,view);
158:   }

160:   /* Check result again with BVAXPY */
161:   BVAXPY(Z,-1.0,Y);
162:   BVNorm(Z,NORM_FROBENIUS,&norm);
163:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error: %g\n",(double)norm);

165:   BVDestroy(&X);
166:   BVDestroy(&Y);
167:   BVDestroy(&Z);
168:   BVDestroy(&Zcopy);
169:   MatDestroy(&B);
170:   VecDestroy(&t);
171:   SlepcFinalize();
172:   return 0;
173: }