Actual source code: test1.c

slepc-3.18.3 2023-03-24
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.\n\n";

 13: #include <slepcbv.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Vec               t,v;
 18:   Mat               Q=NULL,M=NULL;
 19:   BV                X,Y;
 20:   PetscInt          i,j,n=10,k=5,l=3,nloc,lda;
 21:   PetscMPIInt       rank;
 22:   PetscScalar       *q,*z;
 23:   const PetscScalar *pX;
 24:   PetscReal         nrm;
 25:   PetscViewer       view;
 26:   PetscBool         verbose,matcuda,testlda=PETSC_FALSE;

 29:   SlepcInitialize(&argc,&argv,(char*)0,help);
 30:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 31:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 32:   PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
 33:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 34:   PetscOptionsHasName(NULL,NULL,"-matcuda",&matcuda);
 35:   PetscOptionsHasName(NULL,NULL,"-testlda",&testlda);
 36:   PetscPrintf(PETSC_COMM_WORLD,"Test BV with %" PetscInt_FMT " columns of dimension %" PetscInt_FMT ".\n",k,n);

 38:   /* Create template vector */
 39:   VecCreate(PETSC_COMM_WORLD,&t);
 40:   VecSetSizes(t,PETSC_DECIDE,n);
 41:   VecSetFromOptions(t);
 42:   VecGetLocalSize(t,&nloc);

 44:   /* Create BV object X */
 45:   BVCreate(PETSC_COMM_WORLD,&X);
 46:   PetscObjectSetName((PetscObject)X,"X");
 47:   BVSetSizesFromVec(X,t,k);
 48:   BVSetFromOptions(X);

 50:   /* Set up viewer */
 51:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 52:   PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_INFO_DETAIL);
 53:   BVView(X,view);
 54:   PetscViewerPopFormat(view);

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

 69:   /* Create BV object Y */
 70:   BVCreate(PETSC_COMM_WORLD,&Y);
 71:   PetscObjectSetName((PetscObject)Y,"Y");
 72:   BVSetSizesFromVec(Y,t,l);
 73:   BVSetFromOptions(Y);

 75:   /* Fill Y entries */
 76:   for (j=0;j<l;j++) {
 77:     BVGetColumn(Y,j,&v);
 78:     VecSet(v,(PetscScalar)(j+1)/4.0);
 79:     BVRestoreColumn(Y,j,&v);
 80:   }
 81:   if (verbose) BVView(Y,view);

 83:   /* Create Mat */
 84:   MatCreate(PETSC_COMM_SELF,&Q);
 85:   if (matcuda && PetscDefined(HAVE_CUDA)) MatSetType(Q,MATSEQDENSECUDA);
 86:   else MatSetType(Q,MATSEQDENSE);
 87:   MatSetSizes(Q,k,l,k,l);
 88:   if (testlda) MatDenseSetLDA(Q,k+2);
 89:   MatSeqDenseSetPreallocation(Q,NULL);
 90:   PetscObjectSetName((PetscObject)Q,"Q");
 91:   MatDenseGetArrayWrite(Q,&q);
 92:   MatDenseGetLDA(Q,&lda);
 93:   for (i=0;i<k;i++)
 94:     for (j=0;j<l;j++)
 95:       q[i+j*lda] = (i<j)? 2.0: -0.5;
 96:   MatDenseRestoreArrayWrite(Q,&q);
 97:   if (verbose) MatView(Q,NULL);

 99:   /* Test BVMult */
100:   BVMult(Y,2.0,1.0,X,Q);
101:   if (verbose) {
102:     PetscPrintf(PETSC_COMM_WORLD,"After BVMult - - - - - - - - -\n");
103:     BVView(Y,view);
104:   }

106:   /* Test BVMultVec */
107:   BVGetColumn(Y,0,&v);
108:   PetscMalloc1(k,&z);
109:   z[0] = 2.0;
110:   for (i=1;i<k;i++) z[i] = -0.5*z[i-1];
111:   BVMultVec(X,-1.0,1.0,v,z);
112:   PetscFree(z);
113:   BVRestoreColumn(Y,0,&v);
114:   if (verbose) {
115:     PetscPrintf(PETSC_COMM_WORLD,"After BVMultVec - - - - - - -\n");
116:     BVView(Y,view);
117:   }

119:   /* Test BVDot */
120:   MatCreate(PETSC_COMM_SELF,&M);
121:   if (matcuda && PetscDefined(HAVE_CUDA)) MatSetType(M,MATSEQDENSECUDA);
122:   else MatSetType(M,MATSEQDENSE);
123:   MatSetSizes(M,l,k,l,k);
124:   if (testlda) MatDenseSetLDA(M,l+2);
125:   MatSeqDenseSetPreallocation(M,NULL);
126:   PetscObjectSetName((PetscObject)M,"M");
127:   BVDot(X,Y,M);
128:   if (verbose) {
129:     PetscPrintf(PETSC_COMM_WORLD,"After BVDot - - - - - - - - -\n");
130:     MatView(M,NULL);
131:   }

133:   /* Test BVDotVec */
134:   BVGetColumn(Y,0,&v);
135:   PetscMalloc1(k,&z);
136:   BVDotVec(X,v,z);
137:   BVRestoreColumn(Y,0,&v);
138:   if (verbose) {
139:     PetscPrintf(PETSC_COMM_WORLD,"After BVDotVec - - - - - - -\n");
140:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,k,z,&v);
141:     PetscObjectSetName((PetscObject)v,"z");
142:     VecView(v,view);
143:     VecDestroy(&v);
144:   }
145:   PetscFree(z);

147:   /* Test BVMultInPlace and BVScale */
148:   BVMultInPlace(X,Q,1,l);
149:   BVScale(X,2.0);
150:   if (verbose) {
151:     PetscPrintf(PETSC_COMM_WORLD,"After BVMultInPlace - - - - -\n");
152:     BVView(X,view);
153:   }

155:   /* Test BVNorm */
156:   BVNormColumn(X,0,NORM_2,&nrm);
157:   PetscPrintf(PETSC_COMM_WORLD,"2-Norm of X[0] = %g\n",(double)nrm);
158:   BVNorm(X,NORM_FROBENIUS,&nrm);
159:   PetscPrintf(PETSC_COMM_WORLD,"Frobenius Norm of X = %g\n",(double)nrm);

161:   /* Test BVGetArrayRead */
162:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
163:   if (!rank) {
164:     PetscPrintf(PETSC_COMM_WORLD,"First row of X =\n");
165:     BVGetArrayRead(X,&pX);
166:     for (i=0;i<k;i++) PetscPrintf(PETSC_COMM_WORLD,"%g ",(double)PetscRealPart(pX[i*nloc]));
167:     PetscPrintf(PETSC_COMM_WORLD,"\n");
168:     BVRestoreArrayRead(X,&pX);
169:   }

171:   BVDestroy(&X);
172:   BVDestroy(&Y);
173:   MatDestroy(&Q);
174:   MatDestroy(&M);
175:   VecDestroy(&t);
176:   SlepcFinalize();
177:   return 0;
178: }

180: /*TEST

182:    test:
183:       args: -bv_type {{vecs contiguous svec mat}separate output} -verbose
184:       suffix: 1
185:       filter: sed -e 's/-0[.]/0./g'

187:    testset:
188:       args: -bv_type svec -vec_type cuda -verbose
189:       requires: cuda
190:       output_file: output/test1_1_cuda.out
191:       test:
192:          suffix: 1_cuda
193:       test:
194:          suffix: 1_cuda_mat
195:          args: -matcuda
196:          filter: sed -e "s/seqdensecuda/seqdense/"

198:    test:
199:       args: -bv_type {{vecs contiguous svec mat}separate output} -verbose -testlda
200:       suffix: 2
201:       filter: sed -e 's/-0[.]/0./g'

203:    testset:
204:       args: -bv_type svec -vec_type cuda -verbose -testlda
205:       requires: cuda
206:       output_file: output/test1_1_cuda.out
207:       test:
208:          suffix: 2_cuda
209:       test:
210:          suffix: 2_cuda_mat
211:          args: -matcuda
212:          filter: sed -e "s/seqdensecuda/seqdense/"

214: TEST*/