Actual source code: test6.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 orthogonalization functions with constraints.\n\n";

 13: #include <slepcbv.h>

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

 27:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 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 %D columns + %D constraints, of length %D.\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++) {
 51:         VecSetValue(C[j],i,1.0,INSERT_VALUES);
 52:       }
 53:       VecAssemblyBegin(C[j]);
 54:       VecAssemblyEnd(C[j]);
 55:     }
 56:     BVInsertConstraints(X,&nc,C);
 57:     VecDestroyVecs(nc,&C);
 58:   }

 60:   /* Set up viewer */
 61:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 62:   if (verbose) {
 63:     PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
 64:   }

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

 84:   /* Test BVOrthogonalizeColumn */
 85:   for (j=0;j<k;j++) {
 86:     BVOrthogonalizeColumn(X,j,NULL,&norm,NULL);
 87:     alpha = 1.0/norm;
 88:     BVScaleColumn(X,j,alpha);
 89:   }
 90:   if (verbose) {
 91:     BVView(X,view);
 92:   }

 94:   /* Check orthogonality */
 95:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
 96:   BVDot(X,X,M);
 97:   MatShift(M,-1.0);
 98:   MatNorm(M,NORM_1,&norm);
 99:   if (norm<100*PETSC_MACHINE_EPSILON) {
100:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
101:   } else {
102:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm);
103:   }

105:   MatDestroy(&M);
106:   BVDestroy(&X);
107:   VecDestroy(&t);
108:   SlepcFinalize();
109:   return ierr;
110: }

112: /*TEST

114:    test:
115:       suffix: 1
116:       nsize: 1
117:       args: -bv_type {{vecs contiguous svec mat}shared output}
118:       output_file: output/test6_1.out

120:    test:
121:       suffix: 1_cuda
122:       nsize: 1
123:       args: -bv_type svec -vec_type cuda
124:       requires: cuda
125:       output_file: output/test6_1.out

127:    test:
128:       suffix: 2
129:       nsize: 2
130:       args: -bv_type {{vecs contiguous svec mat}shared output}
131:       output_file: output/test6_1.out

133:    test:
134:       suffix: 3
135:       nsize: 1
136:       args: -bv_type {{vecs contiguous svec mat}shared output} -bv_orthog_type mgs
137:       output_file: output/test6_1.out

139:    test:
140:       suffix: 3_cuda
141:       nsize: 1
142:       args: -bv_type svec -vec_type cuda -bv_orthog_type mgs
143:       requires: cuda
144:       output_file: output/test6_1.out

146: TEST*/