Actual source code: test1.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.\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;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

174:   BVDestroy(&X);
175:   BVDestroy(&Y);
176:   MatDestroy(&Q);
177:   MatDestroy(&M);
178:   VecDestroy(&t);
179:   SlepcFinalize();
180:   return 0;
181: }

183: /*TEST

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

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

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

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

217: TEST*/