Actual source code: test1.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 operations.\n\n";

 13: #include <slepcbv.h>

 15: int main(int argc,char **argv)
 16: {
 17:   PetscErrorCode    ierr;
 18:   Vec               t,v;
 19:   Mat               Q,M;
 20:   BV                X,Y;
 21:   PetscInt          i,j,n=10,k=5,l=3,nloc;
 22:   PetscMPIInt       rank;
 23:   PetscScalar       *q,*z;
 24:   const PetscScalar *pX;
 25:   PetscReal         nrm;
 26:   PetscViewer       view;
 27:   PetscBool         verbose;

 29:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 30:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 31:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 32:   PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
 33:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 34:   PetscPrintf(PETSC_COMM_WORLD,"Test BV with %D columns of dimension %D.\n",k,n);

 36:   /* Create template vector */
 37:   VecCreate(PETSC_COMM_WORLD,&t);
 38:   VecSetSizes(t,PETSC_DECIDE,n);
 39:   VecSetFromOptions(t);
 40:   VecGetLocalSize(t,&nloc);

 42:   /* Create BV object X */
 43:   BVCreate(PETSC_COMM_WORLD,&X);
 44:   PetscObjectSetName((PetscObject)X,"X");
 45:   BVSetSizesFromVec(X,t,k);
 46:   BVSetFromOptions(X);

 48:   /* Set up viewer */
 49:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 50:   PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_INFO_DETAIL);
 51:   BVView(X,view);
 52:   PetscViewerPopFormat(view);

 54:   /* Fill X entries */
 55:   for (j=0;j<k;j++) {
 56:     BVGetColumn(X,j,&v);
 57:     VecSet(v,0.0);
 58:     for (i=0;i<4;i++) {
 59:       if (i+j<n) {
 60:         VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);
 61:       }
 62:     }
 63:     VecAssemblyBegin(v);
 64:     VecAssemblyEnd(v);
 65:     BVRestoreColumn(X,j,&v);
 66:   }
 67:   if (verbose) {
 68:     BVView(X,view);
 69:   }

 71:   /* Create BV object Y */
 72:   BVCreate(PETSC_COMM_WORLD,&Y);
 73:   PetscObjectSetName((PetscObject)Y,"Y");
 74:   BVSetSizesFromVec(Y,t,l);
 75:   BVSetFromOptions(Y);

 77:   /* Fill Y entries */
 78:   for (j=0;j<l;j++) {
 79:     BVGetColumn(Y,j,&v);
 80:     VecSet(v,(PetscScalar)(j+1)/4.0);
 81:     BVRestoreColumn(Y,j,&v);
 82:   }
 83:   if (verbose) {
 84:     BVView(Y,view);
 85:   }

 87:   /* Create Mat */
 88:   MatCreateSeqDense(PETSC_COMM_SELF,k,l,NULL,&Q);
 89:   PetscObjectSetName((PetscObject)Q,"Q");
 90:   MatDenseGetArray(Q,&q);
 91:   for (i=0;i<k;i++)
 92:     for (j=0;j<l;j++)
 93:       q[i+j*k] = (i<j)? 2.0: -0.5;
 94:   MatDenseRestoreArray(Q,&q);
 95:   if (verbose) {
 96:     MatView(Q,NULL);
 97:   }

 99:   /* Test BVMult */
100:   BVMult(Y,2.0,1.0,X,Q);
101:   if (verbose) {
102:     PetscPrintf(PETSC_COMM_WORLD,"After BVMult - - - - - - - - -\n");
103:     BVView(Y,view);
104:   }

106:   /* Test BVMultVec */
107:   BVGetColumn(Y,0,&v);
108:   PetscMalloc1(k,&z);
109:   z[0] = 2.0;
110:   for (i=1;i<k;i++) z[i] = -0.5*z[i-1];
111:   BVMultVec(X,-1.0,1.0,v,z);
112:   PetscFree(z);
113:   BVRestoreColumn(Y,0,&v);
114:   if (verbose) {
115:     PetscPrintf(PETSC_COMM_WORLD,"After BVMultVec - - - - - - -\n");
116:     BVView(Y,view);
117:   }

119:   /* Test BVDot */
120:   MatCreateSeqDense(PETSC_COMM_SELF,l,k,NULL,&M);
121:   PetscObjectSetName((PetscObject)M,"M");
122:   BVDot(X,Y,M);
123:   if (verbose) {
124:     PetscPrintf(PETSC_COMM_WORLD,"After BVDot - - - - - - - - -\n");
125:     MatView(M,NULL);
126:   }

128:   /* Test BVDotVec */
129:   BVGetColumn(Y,0,&v);
130:   PetscMalloc1(k,&z);
131:   BVDotVec(X,v,z);
132:   BVRestoreColumn(Y,0,&v);
133:   if (verbose) {
134:     PetscPrintf(PETSC_COMM_WORLD,"After BVDotVec - - - - - - -\n");
135:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,k,z,&v);
136:     PetscObjectSetName((PetscObject)v,"z");
137:     VecView(v,view);
138:     VecDestroy(&v);
139:   }
140:   PetscFree(z);

142:   /* Test BVMultInPlace and BVScale */
143:   BVMultInPlace(X,Q,1,l);
144:   BVScale(X,2.0);
145:   if (verbose) {
146:     PetscPrintf(PETSC_COMM_WORLD,"After BVMultInPlace - - - - -\n");
147:     BVView(X,view);
148:   }

150:   /* Test BVNorm */
151:   BVNormColumn(X,0,NORM_2,&nrm);
152:   PetscPrintf(PETSC_COMM_WORLD,"2-Norm of X[0] = %g\n",(double)nrm);
153:   BVNorm(X,NORM_FROBENIUS,&nrm);
154:   PetscPrintf(PETSC_COMM_WORLD,"Frobenius Norm of X = %g\n",(double)nrm);

156:   /* Test BVGetArrayRead */
157:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
158:   if (!rank) {
159:     PetscPrintf(PETSC_COMM_WORLD,"First row of X =\n");
160:     BVGetArrayRead(X,&pX);
161:     for (i=0;i<k;i++) {
162:       PetscPrintf(PETSC_COMM_WORLD,"%g ",(double)PetscRealPart(pX[i*nloc]));
163:     }
164:     PetscPrintf(PETSC_COMM_WORLD,"\n");
165:     BVRestoreArrayRead(X,&pX);
166:   }

168:   BVDestroy(&X);
169:   BVDestroy(&Y);
170:   MatDestroy(&Q);
171:   MatDestroy(&M);
172:   VecDestroy(&t);
173:   SlepcFinalize();
174:   return ierr;
175: }

177: /*TEST

179:    test:
180:       suffix: 1
181:       nsize: 1
182:       args: -bv_type {{vecs contiguous svec mat}separate output} -verbose

184:    test:
185:       suffix: 1_cuda
186:       nsize: 1
187:       args: -bv_type svec -vec_type cuda -verbose
188:       requires: cuda
189:       filter: sed -e "s/type: seqcuda/type: seq/"

191: TEST*/