Actual source code: test1.c

slepc-3.21.2 2024-09-25
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,ldx,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:   PetscFunctionBeginUser;
 29:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 30:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 31:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
 32:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL));
 33:   PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
 34:   PetscCall(PetscOptionsHasName(NULL,NULL,"-matcuda",&matcuda));
 35:   PetscCall(PetscOptionsHasName(NULL,NULL,"-testlda",&testlda));
 36:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test BV with %" PetscInt_FMT " columns of dimension %" PetscInt_FMT ".\n",k,n));

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

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

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

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

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

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

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

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

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

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

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

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

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

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

171:   PetscCall(BVDestroy(&X));
172:   PetscCall(BVDestroy(&Y));
173:   PetscCall(MatDestroy(&Q));
174:   PetscCall(MatDestroy(&M));
175:   PetscCall(VecDestroy(&t));
176:   PetscCall(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_svec_cuda.out
191:       test:
192:          suffix: 1_svec_cuda
193:       test:
194:          suffix: 1_svec_cuda_mat
195:          args: -matcuda
196:          filter: sed -e "s/seqdensecuda/seqdense/"

198:    testset:
199:       args: -bv_type mat -vec_type cuda -verbose
200:       requires: cuda
201:       output_file: output/test1_1_mat_cuda.out
202:       filter: sed -e "s/seqdensecuda/seqdense/"
203:       test:
204:          suffix: 1_mat_cuda
205:       test:
206:          suffix: 1_mat_cuda_mat
207:          args: -matcuda

209:    test:
210:       args: -bv_type {{vecs contiguous svec mat}separate output} -verbose -testlda
211:       suffix: 2
212:       filter: sed -e 's/-0[.]/0./g'

214:    testset:
215:       args: -bv_type svec -vec_type cuda -verbose -testlda
216:       requires: cuda
217:       output_file: output/test1_1_svec_cuda.out
218:       test:
219:          suffix: 2_svec_cuda
220:       test:
221:          suffix: 2_svec_cuda_mat
222:          args: -matcuda
223:          filter: sed -e "s/seqdensecuda/seqdense/"

225:    testset:
226:       args: -bv_type mat -vec_type cuda -verbose -testlda
227:       requires: cuda
228:       output_file: output/test1_1_mat_cuda.out
229:       filter: sed -e "s/seqdensecuda/seqdense/"
230:       test:
231:          suffix: 2_mat_cuda
232:       test:
233:          suffix: 2_mat_cuda_mat
234:          args: -matcuda

236: TEST*/