Actual source code: test13.c

slepc-3.17.2 2022-08-09
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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: {
 17:   Vec            t,v,z;
 18:   BV             X;
 19:   PetscInt       i,j,n=10,k=5,l=3;
 20:   PetscReal      nrm;
 21:   PetscViewer    view;
 22:   PetscBool      verbose;

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

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

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

 42:   /* Set up viewer */
 43:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 44:   if (verbose) PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);

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

 59:   /* Test BVDotColumn */
 60:   BVDotColumn(X,2,NULL);
 61:   if (verbose) {
 62:     PetscPrintf(PETSC_COMM_WORLD,"After BVDotColumn - - - - - - -\n");
 63:     BVGetBufferVec(X,&z);
 64:     VecView(z,view);
 65:   }
 66:   /* Test BVMultColumn */
 67:   BVMultColumn(X,-1.0,1.0,2,NULL);
 68:   if (verbose) {
 69:     PetscPrintf(PETSC_COMM_WORLD,"After BVMultColumn - - - - - - - - -\n");
 70:     BVView(X,view);
 71:   }

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

 76:   BVDestroy(&X);
 77:   VecDestroy(&t);
 78:   SlepcFinalize();
 79:   return 0;
 80: }

 82: /*TEST

 84:    testset:
 85:       output_file: output/test13_1.out
 86:       test:
 87:          suffix: 1
 88:          args: -bv_type {{vecs contiguous svec mat}shared output}
 89:       test:
 90:          suffix: 1_cuda
 91:          args: -bv_type svec -vec_type cuda
 92:          requires: cuda

 94: TEST*/