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 tensor BV.\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 S,M,Q;
19 8 : BV U,V,UU;
20 8 : PetscInt i,ii,j,jj,n=10,k=6,l=3,d=3,deg,id,lds;
21 8 : PetscScalar *pS,*q;
22 8 : PetscReal norm;
23 8 : PetscViewer view;
24 8 : PetscBool verbose;
25 :
26 8 : PetscFunctionBeginUser;
27 8 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
28 8 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
29 8 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
30 8 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL));
31 8 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-d",&d,NULL));
32 8 : PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
33 8 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test tensor BV of degree %" PetscInt_FMT " with %" PetscInt_FMT " columns of dimension %" PetscInt_FMT "*d.\n",d,k,n));
34 :
35 : /* Create template vector */
36 8 : PetscCall(VecCreate(PETSC_COMM_WORLD,&t));
37 8 : PetscCall(VecSetSizes(t,PETSC_DECIDE,n));
38 8 : PetscCall(VecSetFromOptions(t));
39 :
40 : /* Create BV object U */
41 8 : PetscCall(BVCreate(PETSC_COMM_WORLD,&U));
42 8 : PetscCall(PetscObjectSetName((PetscObject)U,"U"));
43 8 : PetscCall(BVSetSizesFromVec(U,t,k+d-1));
44 8 : PetscCall(BVSetFromOptions(U));
45 8 : PetscCall(PetscObjectSetName((PetscObject)U,"U"));
46 :
47 : /* Fill first d columns of U */
48 32 : for (j=0;j<d;j++) {
49 24 : PetscCall(BVGetColumn(U,j,&v));
50 24 : PetscCall(VecSet(v,0.0));
51 120 : for (i=0;i<4;i++) {
52 96 : if (i+j<n) PetscCall(VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES));
53 : }
54 24 : PetscCall(VecAssemblyBegin(v));
55 24 : PetscCall(VecAssemblyEnd(v));
56 24 : PetscCall(BVRestoreColumn(U,j,&v));
57 : }
58 :
59 : /* Create tensor BV */
60 8 : PetscCall(BVCreateTensor(U,d,&V));
61 8 : PetscCall(PetscObjectSetName((PetscObject)V,"V"));
62 8 : PetscCall(BVTensorGetDegree(V,°));
63 8 : PetscCheck(deg==d,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Wrong degree");
64 :
65 : /* Set up viewer */
66 8 : PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view));
67 8 : PetscCall(PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_INFO_DETAIL));
68 8 : PetscCall(BVView(V,view));
69 8 : PetscCall(PetscViewerPopFormat(view));
70 8 : if (verbose) {
71 0 : PetscCall(PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB));
72 0 : PetscCall(BVView(V,view));
73 : }
74 :
75 : /* Build first column from previously introduced coefficients */
76 8 : PetscCall(BVTensorBuildFirstColumn(V,d));
77 8 : if (verbose) {
78 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After building the first column - - - - -\n"));
79 0 : PetscCall(BVView(V,view));
80 : }
81 :
82 : /* Test orthogonalization */
83 8 : PetscCall(BVTensorGetFactors(V,&UU,&S));
84 8 : PetscCall(BVGetActiveColumns(UU,NULL,&j));
85 8 : PetscCall(BVGetSizes(UU,NULL,NULL,&id));
86 8 : PetscCheck(id==k+d-1,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Wrong dimensions");
87 8 : lds = id*d;
88 48 : for (jj=1;jj<k;jj++) {
89 : /* set new orthogonal column in U */
90 40 : PetscCall(BVGetColumn(UU,j,&v));
91 40 : PetscCall(VecSet(v,0.0));
92 200 : for (i=0;i<4;i++) {
93 160 : if (i+j<n) PetscCall(VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES));
94 : }
95 40 : PetscCall(VecAssemblyBegin(v));
96 40 : PetscCall(VecAssemblyEnd(v));
97 40 : PetscCall(BVRestoreColumn(UU,j,&v));
98 40 : PetscCall(BVOrthonormalizeColumn(UU,j,PETSC_TRUE,NULL,NULL));
99 40 : j++;
100 40 : PetscCall(BVSetActiveColumns(UU,0,j));
101 : /* set new column of S */
102 40 : PetscCall(MatDenseGetArray(S,&pS));
103 160 : for (ii=0;ii<d;ii++) {
104 720 : for (i=0;i<ii+jj+1;i++) {
105 600 : pS[i+ii*id+jj*lds] = (PetscScalar)(2*ii+i+0.5*jj);
106 : }
107 : }
108 40 : PetscCall(MatDenseRestoreArray(S,&pS));
109 40 : PetscCall(BVOrthonormalizeColumn(V,jj,PETSC_TRUE,NULL,NULL));
110 : }
111 8 : PetscCall(BVTensorRestoreFactors(V,&UU,&S));
112 8 : if (verbose) {
113 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After orthogonalization - - - - -\n"));
114 0 : PetscCall(BVView(V,view));
115 : }
116 :
117 : /* Check orthogonality */
118 8 : PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M));
119 8 : PetscCall(BVDot(V,V,M));
120 8 : PetscCall(MatShift(M,-1.0));
121 8 : PetscCall(MatNorm(M,NORM_1,&norm));
122 8 : if (norm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n"));
123 0 : else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm));
124 :
125 : /* Test BVTensorCompress */
126 8 : PetscCall(BVSetActiveColumns(V,0,l));
127 8 : PetscCall(BVTensorCompress(V,0));
128 8 : if (verbose) {
129 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After BVTensorCompress - - - - -\n"));
130 0 : PetscCall(BVView(V,view));
131 : }
132 :
133 : /* Create Mat */
134 8 : PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,k,l,NULL,&Q));
135 8 : PetscCall(PetscObjectSetName((PetscObject)Q,"Q"));
136 8 : PetscCall(MatDenseGetArray(Q,&q));
137 56 : for (i=0;i<k;i++)
138 192 : for (j=0;j<l;j++)
139 264 : q[i+j*k] = (i<j)? 2.0: -0.5;
140 8 : PetscCall(MatDenseRestoreArray(Q,&q));
141 8 : if (verbose) PetscCall(MatView(Q,NULL));
142 :
143 : /* Test BVMultInPlace */
144 8 : PetscCall(BVMultInPlace(V,Q,1,l));
145 8 : if (verbose) {
146 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After BVMultInPlace - - - - -\n"));
147 0 : PetscCall(BVView(V,view));
148 : }
149 :
150 : /* Test BVNorm */
151 8 : PetscCall(BVNorm(V,NORM_1,&norm));
152 8 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm: %g\n",(double)norm));
153 :
154 8 : PetscCall(BVDestroy(&U));
155 8 : PetscCall(BVDestroy(&V));
156 8 : PetscCall(MatDestroy(&Q));
157 8 : PetscCall(MatDestroy(&M));
158 8 : PetscCall(VecDestroy(&t));
159 8 : PetscCall(SlepcFinalize());
160 : return 0;
161 : }
162 :
163 : /*TEST
164 :
165 : testset:
166 : nsize: 2
167 : output_file: output/test16_1.out
168 : filter: grep -v "doing matmult"
169 : test:
170 : suffix: 1
171 : args: -bv_type {{vecs contiguous svec mat}}
172 : test:
173 : suffix: 1_cuda
174 : args: -bv_type {{vecs svec mat}} -vec_type cuda
175 : requires: cuda
176 : test:
177 : suffix: 1_hip
178 : args: -bv_type {{vecs svec mat}} -vec_type hip
179 : requires: hip
180 :
181 : TEST*/
|