Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
11 : static char help[] = "Test BV operations.\n\n";
12 :
13 : #include <slepcbv.h>
14 :
15 8 : int main(int argc,char **argv)
16 : {
17 8 : Vec t,v;
18 8 : Mat Q=NULL,M=NULL;
19 8 : BV X,Y;
20 8 : PetscInt i,j,n=10,k=5,l=3,ldx,lda;
21 8 : PetscMPIInt rank;
22 8 : PetscScalar *q,*z;
23 8 : const PetscScalar *pX;
24 8 : PetscReal nrm;
25 8 : PetscViewer view;
26 8 : PetscBool verbose,matcuda,mathip,testlda=PETSC_FALSE;
27 :
28 8 : PetscFunctionBeginUser;
29 8 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
30 8 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
31 8 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
32 8 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL));
33 8 : PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
34 8 : PetscCall(PetscOptionsHasName(NULL,NULL,"-matcuda",&matcuda));
35 8 : PetscCall(PetscOptionsHasName(NULL,NULL,"-mathip",&mathip));
36 8 : PetscCall(PetscOptionsHasName(NULL,NULL,"-testlda",&testlda));
37 8 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test BV with %" PetscInt_FMT " columns of dimension %" PetscInt_FMT ".\n",k,n));
38 :
39 : /* Create template vector */
40 8 : PetscCall(VecCreate(PETSC_COMM_WORLD,&t));
41 8 : PetscCall(VecSetSizes(t,PETSC_DECIDE,n));
42 8 : PetscCall(VecSetFromOptions(t));
43 :
44 : /* Create BV object X */
45 8 : PetscCall(BVCreate(PETSC_COMM_WORLD,&X));
46 8 : PetscCall(PetscObjectSetName((PetscObject)X,"X"));
47 8 : PetscCall(BVSetSizesFromVec(X,t,k));
48 8 : PetscCall(BVSetFromOptions(X));
49 :
50 : /* Set up viewer */
51 8 : PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view));
52 8 : PetscCall(PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_INFO_DETAIL));
53 8 : PetscCall(BVView(X,view));
54 8 : PetscCall(PetscViewerPopFormat(view));
55 :
56 : /* Fill X entries */
57 48 : for (j=0;j<k;j++) {
58 40 : PetscCall(BVGetColumn(X,j,&v));
59 40 : PetscCall(VecSet(v,0.0));
60 200 : for (i=0;i<4;i++) {
61 160 : if (i+j<n) PetscCall(VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES));
62 : }
63 40 : PetscCall(VecAssemblyBegin(v));
64 40 : PetscCall(VecAssemblyEnd(v));
65 40 : PetscCall(BVRestoreColumn(X,j,&v));
66 : }
67 8 : if (verbose) PetscCall(BVView(X,view));
68 :
69 : /* Create BV object Y */
70 8 : PetscCall(BVCreate(PETSC_COMM_WORLD,&Y));
71 8 : PetscCall(PetscObjectSetName((PetscObject)Y,"Y"));
72 8 : PetscCall(BVSetSizesFromVec(Y,t,l));
73 8 : PetscCall(BVSetFromOptions(Y));
74 :
75 : /* Fill Y entries */
76 32 : for (j=0;j<l;j++) {
77 24 : PetscCall(BVGetColumn(Y,j,&v));
78 24 : PetscCall(VecSet(v,(PetscScalar)(j+1)/4.0));
79 24 : PetscCall(BVRestoreColumn(Y,j,&v));
80 : }
81 8 : if (verbose) PetscCall(BVView(Y,view));
82 :
83 : /* Create Mat */
84 8 : PetscCall(MatCreate(PETSC_COMM_SELF,&Q));
85 8 : if (matcuda && PetscDefined(HAVE_CUDA)) PetscCall(MatSetType(Q,MATSEQDENSECUDA));
86 8 : else if (mathip && PetscDefined(HAVE_HIP)) PetscCall(MatSetType(Q,MATSEQDENSEHIP));
87 8 : else PetscCall(MatSetType(Q,MATSEQDENSE));
88 8 : PetscCall(MatSetSizes(Q,k,l,k,l));
89 8 : if (testlda) PetscCall(MatDenseSetLDA(Q,k+2));
90 8 : PetscCall(MatSeqDenseSetPreallocation(Q,NULL));
91 8 : PetscCall(PetscObjectSetName((PetscObject)Q,"Q"));
92 8 : PetscCall(MatDenseGetArrayWrite(Q,&q));
93 8 : PetscCall(MatDenseGetLDA(Q,&lda));
94 48 : for (i=0;i<k;i++)
95 160 : for (j=0;j<l;j++)
96 216 : q[i+j*lda] = (i<j)? 2.0: -0.5;
97 8 : PetscCall(MatDenseRestoreArrayWrite(Q,&q));
98 8 : if (verbose) PetscCall(MatView(Q,NULL));
99 :
100 : /* Test BVMult */
101 8 : PetscCall(BVMult(Y,2.0,1.0,X,Q));
102 8 : if (verbose) {
103 8 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After BVMult - - - - - - - - -\n"));
104 8 : PetscCall(BVView(Y,view));
105 : }
106 :
107 : /* Test BVMultVec */
108 8 : PetscCall(BVGetColumn(Y,0,&v));
109 8 : PetscCall(PetscMalloc1(k,&z));
110 8 : z[0] = 2.0;
111 40 : for (i=1;i<k;i++) z[i] = -0.5*z[i-1];
112 8 : PetscCall(BVMultVec(X,-1.0,1.0,v,z));
113 8 : PetscCall(PetscFree(z));
114 8 : PetscCall(BVRestoreColumn(Y,0,&v));
115 8 : if (verbose) {
116 8 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After BVMultVec - - - - - - -\n"));
117 8 : PetscCall(BVView(Y,view));
118 : }
119 :
120 : /* Test BVDot */
121 8 : PetscCall(MatCreate(PETSC_COMM_SELF,&M));
122 8 : if (matcuda && PetscDefined(HAVE_CUDA)) PetscCall(MatSetType(M,MATSEQDENSECUDA));
123 8 : else if (mathip && PetscDefined(HAVE_HIP)) PetscCall(MatSetType(M,MATSEQDENSEHIP));
124 8 : else PetscCall(MatSetType(M,MATSEQDENSE));
125 8 : PetscCall(MatSetSizes(M,l,k,l,k));
126 8 : if (testlda) PetscCall(MatDenseSetLDA(M,l+2));
127 8 : PetscCall(MatSeqDenseSetPreallocation(M,NULL));
128 8 : PetscCall(PetscObjectSetName((PetscObject)M,"M"));
129 8 : PetscCall(BVDot(X,Y,M));
130 8 : if (verbose) {
131 8 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After BVDot - - - - - - - - -\n"));
132 8 : PetscCall(MatView(M,NULL));
133 : }
134 :
135 : /* Test BVDotVec */
136 8 : PetscCall(BVGetColumn(Y,0,&v));
137 8 : PetscCall(PetscMalloc1(k,&z));
138 8 : PetscCall(BVDotVec(X,v,z));
139 8 : PetscCall(BVRestoreColumn(Y,0,&v));
140 8 : if (verbose) {
141 8 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After BVDotVec - - - - - - -\n"));
142 8 : PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,k,z,&v));
143 8 : PetscCall(PetscObjectSetName((PetscObject)v,"z"));
144 8 : PetscCall(VecView(v,view));
145 8 : PetscCall(VecDestroy(&v));
146 : }
147 8 : PetscCall(PetscFree(z));
148 :
149 : /* Test BVMultInPlace and BVScale */
150 8 : PetscCall(BVMultInPlace(X,Q,1,l));
151 8 : PetscCall(BVScale(X,2.0));
152 8 : if (verbose) {
153 8 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After BVMultInPlace - - - - -\n"));
154 8 : PetscCall(BVView(X,view));
155 : }
156 :
157 : /* Test BVNorm */
158 8 : PetscCall(BVNormColumn(X,0,NORM_2,&nrm));
159 8 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"2-Norm of X[0] = %g\n",(double)nrm));
160 8 : PetscCall(BVNorm(X,NORM_FROBENIUS,&nrm));
161 8 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Frobenius Norm of X = %g\n",(double)nrm));
162 :
163 : /* Test BVGetArrayRead */
164 8 : PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
165 8 : if (!rank) {
166 8 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"First row of X =\n"));
167 8 : PetscCall(BVGetLeadingDimension(X,&ldx));
168 8 : PetscCall(BVGetArrayRead(X,&pX));
169 48 : for (i=0;i<k;i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%g ",(double)PetscRealPart(pX[i*ldx])));
170 8 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
171 8 : PetscCall(BVRestoreArrayRead(X,&pX));
172 : }
173 :
174 8 : PetscCall(BVDestroy(&X));
175 8 : PetscCall(BVDestroy(&Y));
176 8 : PetscCall(MatDestroy(&Q));
177 8 : PetscCall(MatDestroy(&M));
178 8 : PetscCall(VecDestroy(&t));
179 8 : PetscCall(SlepcFinalize());
180 : return 0;
181 : }
182 :
183 : /*TEST
184 :
185 : test:
186 : args: -bv_type {{vecs contiguous svec mat}separate output} -verbose
187 : suffix: 1
188 : filter: sed -e 's/-0[.]/0./g'
189 :
190 : testset:
191 : args: -bv_type svec -vec_type cuda -verbose
192 : requires: cuda
193 : output_file: output/test1_1_svec_gpu.out
194 : test:
195 : suffix: 1_svec_cuda
196 : test:
197 : suffix: 1_svec_cuda_mat
198 : args: -matcuda
199 : filter: sed -e "s/seqdensecuda/seqdense/"
200 :
201 : testset:
202 : args: -bv_type mat -vec_type cuda -verbose
203 : requires: cuda
204 : output_file: output/test1_1_mat_gpu.out
205 : filter: sed -e "s/seqdensecuda/seqdense/"
206 : test:
207 : suffix: 1_mat_cuda
208 : test:
209 : suffix: 1_mat_cuda_mat
210 : args: -matcuda
211 :
212 : testset:
213 : args: -bv_type svec -vec_type hip -verbose
214 : requires: hip
215 : output_file: output/test1_1_svec_gpu.out
216 : filter: sed -e "s/seqhip/seqcuda/" | sed -e "s/HIP/CUDA/"
217 : test:
218 : suffix: 1_svec_hip
219 : test:
220 : suffix: 1_svec_hip_mat
221 : args: -mathip
222 : filter: sed -e "s/seqhip/seqcuda/" | sed -e "s/seqdensehip/seqdense/" | sed -e "s/HIP/CUDA/"
223 :
224 : testset:
225 : args: -bv_type mat -vec_type hip -verbose
226 : requires: hip
227 : output_file: output/test1_1_mat_gpu.out
228 : filter: sed -e "s/seqdensehip/seqdense/" | sed -e "s/HIP/CUDA/"
229 : test:
230 : suffix: 1_mat_hip
231 : test:
232 : suffix: 1_mat_hip_mat
233 : args: -mathip
234 :
235 : test:
236 : args: -bv_type {{vecs contiguous svec mat}separate output} -verbose -testlda
237 : suffix: 2
238 : filter: sed -e 's/-0[.]/0./g'
239 :
240 : testset:
241 : args: -bv_type svec -vec_type cuda -verbose -testlda
242 : requires: cuda
243 : output_file: output/test1_1_svec_gpu.out
244 : test:
245 : suffix: 2_svec_cuda
246 : test:
247 : suffix: 2_svec_cuda_mat
248 : args: -matcuda
249 : filter: sed -e "s/seqdensecuda/seqdense/"
250 :
251 : testset:
252 : args: -bv_type mat -vec_type cuda -verbose -testlda
253 : requires: cuda
254 : output_file: output/test1_1_mat_gpu.out
255 : filter: sed -e "s/seqdensecuda/seqdense/"
256 : test:
257 : suffix: 2_mat_cuda
258 : test:
259 : suffix: 2_mat_cuda_mat
260 : args: -matcuda
261 :
262 : testset:
263 : args: -bv_type svec -vec_type hip -verbose -testlda
264 : requires: hip
265 : output_file: output/test1_1_svec_gpu.out
266 : test:
267 : suffix: 2_svec_hip
268 : filter: sed -e "s/seqhip/seqcuda/" | sed -e "s/HIP/CUDA/"
269 : test:
270 : suffix: 2_svec_hip_mat
271 : args: -mathip
272 : filter: sed -e "s/seqhip/seqcuda/" | sed -e "s/seqdensehip/seqdense/" | sed -e "s/HIP/CUDA/"
273 :
274 : testset:
275 : args: -bv_type mat -vec_type hip -verbose -testlda
276 : requires: hip
277 : output_file: output/test1_1_mat_gpu.out
278 : filter: sed -e "s/seqdensehip/seqdense/" | sed -e "s/HIP/CUDA/"
279 : test:
280 : suffix: 2_mat_hip
281 : test:
282 : suffix: 2_mat_hip_mat
283 : args: -mathip
284 :
285 : TEST*/
|