Actual source code: test1.c

slepc-3.16.1 2021-11-17
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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=NULL,M=NULL;
 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,matcuda;

 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:   PetscOptionsHasName(NULL,NULL,"-matcuda",&matcuda);
 35:   PetscPrintf(PETSC_COMM_WORLD,"Test BV with %D columns of dimension %D.\n",k,n);

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

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

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

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

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

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

 88:   /* Create Mat */
 89:   if (matcuda) {
 90: #if defined(PETSC_HAVE_CUDA)
 91:     MatCreateSeqDenseCUDA(PETSC_COMM_SELF,k,l,NULL,&Q);
 92: #endif
 93:   } else {
 94:     MatCreateSeqDense(PETSC_COMM_SELF,k,l,NULL,&Q);
 95:   }
 96:   PetscObjectSetName((PetscObject)Q,"Q");
 97:   MatDenseGetArray(Q,&q);
 98:   for (i=0;i<k;i++)
 99:     for (j=0;j<l;j++)
100:       q[i+j*k] = (i<j)? 2.0: -0.5;
101:   MatDenseRestoreArray(Q,&q);
102:   if (verbose) {
103:     MatView(Q,NULL);
104:   }

106:   /* Test BVMult */
107:   BVMult(Y,2.0,1.0,X,Q);
108:   if (verbose) {
109:     PetscPrintf(PETSC_COMM_WORLD,"After BVMult - - - - - - - - -\n");
110:     BVView(Y,view);
111:   }

113:   /* Test BVMultVec */
114:   BVGetColumn(Y,0,&v);
115:   PetscMalloc1(k,&z);
116:   z[0] = 2.0;
117:   for (i=1;i<k;i++) z[i] = -0.5*z[i-1];
118:   BVMultVec(X,-1.0,1.0,v,z);
119:   PetscFree(z);
120:   BVRestoreColumn(Y,0,&v);
121:   if (verbose) {
122:     PetscPrintf(PETSC_COMM_WORLD,"After BVMultVec - - - - - - -\n");
123:     BVView(Y,view);
124:   }

126:   /* Test BVDot */
127:   if (matcuda) {
128: #if defined(PETSC_HAVE_CUDA)
129:     MatCreateSeqDenseCUDA(PETSC_COMM_SELF,l,k,NULL,&M);
130: #endif
131:   } else {
132:     MatCreateSeqDense(PETSC_COMM_SELF,l,k,NULL,&M);
133:   }
134:   PetscObjectSetName((PetscObject)M,"M");
135:   BVDot(X,Y,M);
136:   if (verbose) {
137:     PetscPrintf(PETSC_COMM_WORLD,"After BVDot - - - - - - - - -\n");
138:     MatView(M,NULL);
139:   }

141:   /* Test BVDotVec */
142:   BVGetColumn(Y,0,&v);
143:   PetscMalloc1(k,&z);
144:   BVDotVec(X,v,z);
145:   BVRestoreColumn(Y,0,&v);
146:   if (verbose) {
147:     PetscPrintf(PETSC_COMM_WORLD,"After BVDotVec - - - - - - -\n");
148:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,k,z,&v);
149:     PetscObjectSetName((PetscObject)v,"z");
150:     VecView(v,view);
151:     VecDestroy(&v);
152:   }
153:   PetscFree(z);

155:   /* Test BVMultInPlace and BVScale */
156:   BVMultInPlace(X,Q,1,l);
157:   BVScale(X,2.0);
158:   if (verbose) {
159:     PetscPrintf(PETSC_COMM_WORLD,"After BVMultInPlace - - - - -\n");
160:     BVView(X,view);
161:   }

163:   /* Test BVNorm */
164:   BVNormColumn(X,0,NORM_2,&nrm);
165:   PetscPrintf(PETSC_COMM_WORLD,"2-Norm of X[0] = %g\n",(double)nrm);
166:   BVNorm(X,NORM_FROBENIUS,&nrm);
167:   PetscPrintf(PETSC_COMM_WORLD,"Frobenius Norm of X = %g\n",(double)nrm);

169:   /* Test BVGetArrayRead */
170:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
171:   if (!rank) {
172:     PetscPrintf(PETSC_COMM_WORLD,"First row of X =\n");
173:     BVGetArrayRead(X,&pX);
174:     for (i=0;i<k;i++) {
175:       PetscPrintf(PETSC_COMM_WORLD,"%g ",(double)PetscRealPart(pX[i*nloc]));
176:     }
177:     PetscPrintf(PETSC_COMM_WORLD,"\n");
178:     BVRestoreArrayRead(X,&pX);
179:   }

181:   BVDestroy(&X);
182:   BVDestroy(&Y);
183:   MatDestroy(&Q);
184:   MatDestroy(&M);
185:   VecDestroy(&t);
186:   SlepcFinalize();
187:   return ierr;
188: }

190: /*TEST

192:    test:
193:       suffix: 1
194:       nsize: 1
195:       args: -bv_type {{vecs contiguous svec mat}separate output} -verbose

197:    testset:
198:       args: -bv_type svec -vec_type cuda -verbose
199:       requires: cuda
200:       output_file: output/test1_1_cuda.out
201:       test:
202:          suffix: 1_cuda
203:       test:
204:          suffix: 1_cuda_mat
205:          args: -matcuda
206:          filter: sed -e "s/seqdensecuda/seqdense/"

208: TEST*/