Actual source code: test13.c

slepc-3.22.1 2024-10-28
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:   PetscFunctionBeginUser;
 25:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
 26:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 27:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
 28:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL));
 29:   PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
 30:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test BV with %" PetscInt_FMT " columns of dimension %" PetscInt_FMT ".\n",k,n));

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

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

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

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

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

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

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

 83: /*TEST

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

 99: TEST*/