Actual source code: test13.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 using internal buffer instead of array arguments.\n\n";

 13: #include <slepcbv.h>

 15: int main(int argc,char **argv)
 16: {
 18:   Vec            t,v,z;
 19:   BV             X;
 20:   PetscInt       i,j,n=10,k=5,l=3;
 21:   PetscReal      nrm;
 22:   PetscViewer    view;
 23:   PetscBool      verbose;

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

 32:   /* Create template vector */
 33:   VecCreate(PETSC_COMM_WORLD,&t);
 34:   VecSetSizes(t,PETSC_DECIDE,n);
 35:   VecSetFromOptions(t);

 37:   /* Create BV object X */
 38:   BVCreate(PETSC_COMM_WORLD,&X);
 39:   PetscObjectSetName((PetscObject)X,"X");
 40:   BVSetSizesFromVec(X,t,k);
 41:   BVSetFromOptions(X);

 43:   /* Set up viewer */
 44:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 45:   if (verbose) {
 46:     PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
 47:   }    

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

 66:   /* Test BVDotColumn */
 67:   BVDotColumn(X,2,NULL);
 68:   if (verbose) {
 69:     PetscPrintf(PETSC_COMM_WORLD,"After BVDotColumn - - - - - - -\n");
 70:     BVGetBufferVec(X,&z);
 71:     VecView(z,view);
 72:   }
 73:   /* Test BVMultColumn */
 74:   BVMultColumn(X,-1.0,1.0,2,NULL);
 75:   if (verbose) {
 76:     PetscPrintf(PETSC_COMM_WORLD,"After BVMultColumn - - - - - - - - -\n");
 77:     BVView(X,view);
 78:   }

 80:   BVNorm(X,NORM_FROBENIUS,&nrm);
 81:   PetscPrintf(PETSC_COMM_WORLD,"Frobenius Norm or X = %g\n",(double)nrm);

 83:   BVDestroy(&X);
 84:   VecDestroy(&t);
 85:   SlepcFinalize();
 86:   return ierr;
 87: }

 89: /*TEST

 91:    test:
 92:       suffix: 1
 93:       nsize: 1
 94:       args: -bv_type {{vecs contiguous svec mat}shared output}
 95:       output_file: output/test13_1.out

 97:    test:
 98:       suffix: 1_cuda
 99:       nsize: 1
100:       args: -bv_type svec -vec_type cuda
101:       requires: cuda
102:       output_file: output/test13_1.out

104: TEST*/