Actual source code: test6.c

slepc-3.18.0 2022-10-01
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 orthogonalization functions with constraints.\n\n";

 13: #include <slepcbv.h>

 15: int main(int argc,char **argv)
 16: {
 17:   BV             X;
 18:   Mat            M;
 19:   Vec            v,t,*C;
 20:   PetscInt       i,j,n=20,k=8,nc=2,nloc;
 21:   PetscViewer    view;
 22:   PetscBool      verbose;
 23:   PetscReal      norm;
 24:   PetscScalar    alpha;

 27:   SlepcInitialize(&argc,&argv,(char*)0,help);
 28:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 29:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 30:   PetscOptionsGetInt(NULL,NULL,"-nc",&nc,NULL);
 31:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 32:   PetscPrintf(PETSC_COMM_WORLD,"Test BV orthogonalization with %" PetscInt_FMT " columns + %" PetscInt_FMT " constraints, of length %" PetscInt_FMT ".\n",k,nc,n);

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

 40:   /* Create BV object X */
 41:   BVCreate(PETSC_COMM_WORLD,&X);
 42:   PetscObjectSetName((PetscObject)X,"X");
 43:   BVSetSizes(X,nloc,n,k);
 44:   BVSetFromOptions(X);

 46:   /* Generate constraints and attach them to X */
 47:   if (nc>0) {
 48:     VecDuplicateVecs(t,nc,&C);
 49:     for (j=0;j<nc;j++) {
 50:       for (i=0;i<=j;i++) VecSetValue(C[j],i,1.0,INSERT_VALUES);
 51:       VecAssemblyBegin(C[j]);
 52:       VecAssemblyEnd(C[j]);
 53:     }
 54:     BVInsertConstraints(X,&nc,C);
 55:     VecDestroyVecs(nc,&C);
 56:   }

 58:   /* Set up viewer */
 59:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 60:   if (verbose) PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);

 62:   /* Fill X entries */
 63:   for (j=0;j<k;j++) {
 64:     BVGetColumn(X,j,&v);
 65:     VecSet(v,0.0);
 66:     for (i=0;i<=n/2;i++) {
 67:       if (i+j<n) {
 68:         alpha = (3.0*i+j-2)/(2*(i+j+1));
 69:         VecSetValue(v,i+j,alpha,INSERT_VALUES);
 70:       }
 71:     }
 72:     VecAssemblyBegin(v);
 73:     VecAssemblyEnd(v);
 74:     BVRestoreColumn(X,j,&v);
 75:   }
 76:   if (verbose) BVView(X,view);

 78:   /* Test BVOrthogonalizeColumn */
 79:   for (j=0;j<k;j++) {
 80:     BVOrthogonalizeColumn(X,j,NULL,&norm,NULL);
 81:     alpha = 1.0/norm;
 82:     BVScaleColumn(X,j,alpha);
 83:   }
 84:   if (verbose) BVView(X,view);

 86:   /* Check orthogonality */
 87:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
 88:   BVDot(X,X,M);
 89:   MatShift(M,-1.0);
 90:   MatNorm(M,NORM_1,&norm);
 91:   if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
 92:   else PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm);

 94:   MatDestroy(&M);
 95:   BVDestroy(&X);
 96:   VecDestroy(&t);
 97:   SlepcFinalize();
 98:   return 0;
 99: }

101: /*TEST

103:    testset:
104:       output_file: output/test6_1.out
105:       test:
106:          suffix: 1
107:          args: -bv_type {{vecs contiguous svec mat}shared output}
108:       test:
109:          suffix: 1_cuda
110:          args: -bv_type svec -vec_type cuda
111:          requires: cuda
112:       test:
113:          suffix: 2
114:          nsize: 2
115:          args: -bv_type {{vecs contiguous svec mat}shared output}
116:       test:
117:          suffix: 3
118:          args: -bv_type {{vecs contiguous svec mat}shared output} -bv_orthog_type mgs
119:       test:
120:          suffix: 3_cuda
121:          args: -bv_type svec -vec_type cuda -bv_orthog_type mgs
122:          requires: cuda

124: TEST*/