xref: /petsc/src/mat/tests/ex241.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
1c7a4214aSPierre Jolivet static char help[] = "Tests MATHTOOL\n\n";
2c7a4214aSPierre Jolivet 
3c7a4214aSPierre Jolivet #include <petscmat.h>
4c7a4214aSPierre Jolivet 
598e73e17SPierre Jolivet static PetscErrorCode GenEntries(PetscInt sdim,PetscInt M,PetscInt N,const PetscInt *J,const PetscInt *K,PetscScalar *ptr,void *ctx)
6c7a4214aSPierre Jolivet {
798e73e17SPierre Jolivet   PetscInt  d,j,k;
8c7a4214aSPierre Jolivet   PetscReal diff = 0.0,*coords = (PetscReal*)(ctx);
9c7a4214aSPierre Jolivet 
10c7a4214aSPierre Jolivet   PetscFunctionBeginUser;
1198e73e17SPierre Jolivet   for (j = 0; j < M; j++) {
1298e73e17SPierre Jolivet     for (k = 0; k < N; k++) {
1398e73e17SPierre Jolivet       diff = 0.0;
1498e73e17SPierre Jolivet       for (d = 0; d < sdim; d++) diff += (coords[J[j]*sdim+d] - coords[K[k]*sdim+d]) * (coords[J[j]*sdim+d] - coords[K[k]*sdim+d]);
1598e73e17SPierre Jolivet       ptr[j+M*k] = 1.0/(1.0e-2 + PetscSqrtReal(diff));
1698e73e17SPierre Jolivet     }
1798e73e17SPierre Jolivet   }
1898e73e17SPierre Jolivet   PetscFunctionReturn(0);
19c7a4214aSPierre Jolivet }
20c7a4214aSPierre Jolivet 
2198e73e17SPierre Jolivet static PetscErrorCode GenEntriesRectangular(PetscInt sdim,PetscInt M,PetscInt N,const PetscInt *J,const PetscInt *K,PetscScalar *ptr,void *ctx)
22c7a4214aSPierre Jolivet {
2398e73e17SPierre Jolivet   PetscInt  d,j,k;
24c7a4214aSPierre Jolivet   PetscReal diff = 0.0,**coords = (PetscReal**)(ctx);
25c7a4214aSPierre Jolivet 
26c7a4214aSPierre Jolivet   PetscFunctionBeginUser;
2798e73e17SPierre Jolivet   for (j = 0; j < M; j++) {
2898e73e17SPierre Jolivet     for (k = 0; k < N; k++) {
2998e73e17SPierre Jolivet       diff = 0.0;
3098e73e17SPierre Jolivet       for (d = 0; d < sdim; d++) diff += (coords[0][J[j]*sdim+d] - coords[1][K[k]*sdim+d]) * (coords[0][J[j]*sdim+d] - coords[1][K[k]*sdim+d]);
3198e73e17SPierre Jolivet       ptr[j+M*k] = 1.0/(1.0e-2 + PetscSqrtReal(diff));
3298e73e17SPierre Jolivet     }
3398e73e17SPierre Jolivet   }
3498e73e17SPierre Jolivet   PetscFunctionReturn(0);
35c7a4214aSPierre Jolivet }
36c7a4214aSPierre Jolivet 
37c7a4214aSPierre Jolivet int main(int argc,char **argv)
38c7a4214aSPierre Jolivet {
39c7a4214aSPierre Jolivet   Mat            A,AT,D,B,P,R,RT;
408f308287SPierre Jolivet   PetscInt       m = 100,dim = 3,M,K = 10,begin,n = 0,N,bs;
41c7a4214aSPierre Jolivet   PetscMPIInt    size;
42c7a4214aSPierre Jolivet   PetscScalar    *ptr;
43c7a4214aSPierre Jolivet   PetscReal      *coords,*gcoords,*scoords,*gscoords,*(ctx[2]),norm,epsilon;
4498e73e17SPierre Jolivet   MatHtoolKernel kernel = GenEntries;
45c7a4214aSPierre Jolivet   PetscBool      flg,sym = PETSC_FALSE;
46c7a4214aSPierre Jolivet   PetscRandom    rdm;
478f308287SPierre Jolivet   IS             iss,ist,is[2];
4898e73e17SPierre Jolivet   Vec            right,left,perm;
49c7a4214aSPierre Jolivet   PetscErrorCode ierr;
50c7a4214aSPierre Jolivet 
51c7a4214aSPierre Jolivet   ierr = PetscInitialize(&argc,&argv,(char*)NULL,help);if (ierr) return ierr;
52c7a4214aSPierre Jolivet   ierr = PetscOptionsGetInt(NULL,NULL,"-m_local",&m,NULL);CHKERRQ(ierr);
53c7a4214aSPierre Jolivet   ierr = PetscOptionsGetInt(NULL,NULL,"-n_local",&n,NULL);CHKERRQ(ierr);
54c7a4214aSPierre Jolivet   ierr = PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL);CHKERRQ(ierr);
55c7a4214aSPierre Jolivet   ierr = PetscOptionsGetInt(NULL,NULL,"-K",&K,NULL);CHKERRQ(ierr);
56c7a4214aSPierre Jolivet   ierr = PetscOptionsGetBool(NULL,NULL,"-symmetric",&sym,NULL);CHKERRQ(ierr);
57c7a4214aSPierre Jolivet   ierr = PetscOptionsGetReal(NULL,NULL,"-mat_htool_epsilon",&epsilon,NULL);CHKERRQ(ierr);
58c7a4214aSPierre Jolivet   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
59c7a4214aSPierre Jolivet   M = size*m;
60c7a4214aSPierre Jolivet   ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr);
61c7a4214aSPierre Jolivet   ierr = PetscMalloc1(m*dim,&coords);CHKERRQ(ierr);
62c7a4214aSPierre Jolivet   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rdm);CHKERRQ(ierr);
63c7a4214aSPierre Jolivet   ierr = PetscRandomGetValuesReal(rdm,m*dim,coords);CHKERRQ(ierr);
64c7a4214aSPierre Jolivet   ierr = PetscCalloc1(M*dim,&gcoords);CHKERRQ(ierr);
65c7a4214aSPierre Jolivet   ierr = MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,M,K,NULL,&B);CHKERRQ(ierr);
66c7a4214aSPierre Jolivet   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
67c7a4214aSPierre Jolivet   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68c7a4214aSPierre Jolivet   ierr = MatSetRandom(B,rdm);CHKERRQ(ierr);
69c7a4214aSPierre Jolivet   ierr = MatGetOwnershipRange(B,&begin,NULL);CHKERRQ(ierr);
70c7a4214aSPierre Jolivet   ierr = PetscArraycpy(gcoords+begin*dim,coords,m*dim);CHKERRQ(ierr);
71c7a4214aSPierre Jolivet   ierr = MPIU_Allreduce(MPI_IN_PLACE,gcoords,M*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr);
72c7a4214aSPierre Jolivet   ierr = MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&A);CHKERRQ(ierr);
73c7a4214aSPierre Jolivet   ierr = MatSetOption(A,MAT_SYMMETRIC,sym);CHKERRQ(ierr);
74c7a4214aSPierre Jolivet   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
75c7a4214aSPierre Jolivet   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
76c7a4214aSPierre Jolivet   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
77c7a4214aSPierre Jolivet   ierr = MatViewFromOptions(A,NULL,"-A_view");CHKERRQ(ierr);
788f308287SPierre Jolivet   ierr = MatGetOwnershipIS(A,is,NULL);CHKERRQ(ierr);
798f308287SPierre Jolivet   ierr = ISDuplicate(is[0],is+1);CHKERRQ(ierr);
808f308287SPierre Jolivet   ierr = MatIncreaseOverlap(A,1,is,2);CHKERRQ(ierr);
818f308287SPierre Jolivet   ierr = MatSetBlockSize(A,2);CHKERRQ(ierr);
828f308287SPierre Jolivet   ierr = MatIncreaseOverlap(A,1,is+1,1);CHKERRQ(ierr);
838f308287SPierre Jolivet   ierr = ISGetBlockSize(is[1],&bs);CHKERRQ(ierr);
84*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(bs != 2,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect block size %" PetscInt_FMT " != 2",bs);
858f308287SPierre Jolivet   ierr = MatSetBlockSize(A,1);CHKERRQ(ierr);
868f308287SPierre Jolivet   ierr = ISEqual(is[0],is[1],&flg);CHKERRQ(ierr);
87*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unequal index sets");
888f308287SPierre Jolivet   ierr = ISDestroy(is);CHKERRQ(ierr);
898f308287SPierre Jolivet   ierr = ISDestroy(is+1);CHKERRQ(ierr);
9098e73e17SPierre Jolivet   ierr = MatCreateVecs(A,&right,&left);CHKERRQ(ierr);
9198e73e17SPierre Jolivet   ierr = VecSetRandom(right,rdm);CHKERRQ(ierr);
9298e73e17SPierre Jolivet   ierr = MatMult(A,right,left);CHKERRQ(ierr);
9398e73e17SPierre Jolivet   ierr = MatHtoolGetPermutationSource(A,&iss);CHKERRQ(ierr);
9498e73e17SPierre Jolivet   ierr = MatHtoolGetPermutationTarget(A,&ist);CHKERRQ(ierr);
9598e73e17SPierre Jolivet   ierr = VecDuplicate(left,&perm);CHKERRQ(ierr);
9698e73e17SPierre Jolivet   ierr = VecCopy(left,perm);CHKERRQ(ierr);
9798e73e17SPierre Jolivet   ierr = VecPermute(perm,ist,PETSC_FALSE);CHKERRQ(ierr);
9898e73e17SPierre Jolivet   ierr = VecPermute(right,iss,PETSC_FALSE);CHKERRQ(ierr);
9998e73e17SPierre Jolivet   ierr = MatHtoolUsePermutation(A,PETSC_FALSE);CHKERRQ(ierr);
10098e73e17SPierre Jolivet   ierr = MatMult(A,right,left);CHKERRQ(ierr);
10198e73e17SPierre Jolivet   ierr = VecAXPY(left,-1.0,perm);CHKERRQ(ierr);
10298e73e17SPierre Jolivet   ierr = VecNorm(left,NORM_INFINITY,&norm);CHKERRQ(ierr);
103*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(PetscAbsReal(norm) > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||y(with permutation)-y(without permutation)|| = %g (> %g)",(double)PetscAbsReal(norm),(double)PETSC_SMALL);
10498e73e17SPierre Jolivet   ierr = MatHtoolUsePermutation(A,PETSC_TRUE);CHKERRQ(ierr);
10598e73e17SPierre Jolivet   ierr = VecDestroy(&perm);CHKERRQ(ierr);
10698e73e17SPierre Jolivet   ierr = VecDestroy(&left);CHKERRQ(ierr);
10798e73e17SPierre Jolivet   ierr = VecDestroy(&right);CHKERRQ(ierr);
10898e73e17SPierre Jolivet   ierr = ISDestroy(&ist);CHKERRQ(ierr);
10998e73e17SPierre Jolivet   ierr = ISDestroy(&iss);CHKERRQ(ierr);
110c7a4214aSPierre Jolivet   if (PetscAbsReal(epsilon) >= PETSC_SMALL) { /* when there is compression, it is more difficult to check against MATDENSE, so just compare symmetric and nonsymmetric assemblies */
111c7a4214aSPierre Jolivet     PetscReal relative;
112c7a4214aSPierre Jolivet     ierr = MatDestroy(&B);CHKERRQ(ierr);
113c7a4214aSPierre Jolivet     ierr = MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&B);CHKERRQ(ierr);
114c7a4214aSPierre Jolivet     ierr = MatSetOption(B,MAT_SYMMETRIC,(PetscBool)!sym);CHKERRQ(ierr);
115c7a4214aSPierre Jolivet     ierr = MatSetFromOptions(B);CHKERRQ(ierr);
116c7a4214aSPierre Jolivet     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
117c7a4214aSPierre Jolivet     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
118c7a4214aSPierre Jolivet     ierr = MatViewFromOptions(B,NULL,"-B_view");CHKERRQ(ierr);
119c7a4214aSPierre Jolivet     ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&P);CHKERRQ(ierr);
120c7a4214aSPierre Jolivet     ierr = MatNorm(P,NORM_FROBENIUS,&relative);CHKERRQ(ierr);
121c7a4214aSPierre Jolivet     ierr = MatConvert(B,MATDENSE,MAT_INITIAL_MATRIX,&R);CHKERRQ(ierr);
122c7a4214aSPierre Jolivet     ierr = MatAXPY(R,-1.0,P,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
123c7a4214aSPierre Jolivet     ierr = MatNorm(R,NORM_INFINITY,&norm);CHKERRQ(ierr);
124*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(PetscAbsReal(norm/relative) > epsilon,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A(!symmetric)-A(symmetric)|| = %g (> %g)",(double)PetscAbsReal(norm/relative),(double)epsilon);
125c7a4214aSPierre Jolivet     ierr = MatDestroy(&B);CHKERRQ(ierr);
126c7a4214aSPierre Jolivet     ierr = MatDestroy(&R);CHKERRQ(ierr);
127c7a4214aSPierre Jolivet     ierr = MatDestroy(&P);CHKERRQ(ierr);
128c7a4214aSPierre Jolivet   } else {
129c7a4214aSPierre Jolivet     ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&D);CHKERRQ(ierr);
130c7a4214aSPierre Jolivet     ierr = MatViewFromOptions(D,NULL,"-D_view");CHKERRQ(ierr);
131c7a4214aSPierre Jolivet     ierr = MatMultEqual(A,D,10,&flg);CHKERRQ(ierr);
132*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx");
1338b8ddd11SPierre Jolivet     ierr = MatMultTransposeEqual(A,D,10,&flg);CHKERRQ(ierr);
134*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx");
135c7a4214aSPierre Jolivet     ierr = MatMultAddEqual(A,D,10,&flg);CHKERRQ(ierr);
136*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"y+Ax != y+Dx");
137c7a4214aSPierre Jolivet     ierr = MatGetOwnershipRange(B,&begin,NULL);CHKERRQ(ierr);
138c7a4214aSPierre Jolivet     ierr = MatDenseGetArrayWrite(D,&ptr);CHKERRQ(ierr);
13998e73e17SPierre Jolivet     for (PetscInt i = begin; i < m+begin; ++i)
14098e73e17SPierre Jolivet       for (PetscInt j = 0; j < M; ++j) GenEntries(dim,1,1,&i,&j,ptr+i-begin+j*m,gcoords);
141c7a4214aSPierre Jolivet     ierr = MatDenseRestoreArrayWrite(D,&ptr);CHKERRQ(ierr);
142c7a4214aSPierre Jolivet     ierr = MatMultEqual(A,D,10,&flg);CHKERRQ(ierr);
143*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx");
144c7a4214aSPierre Jolivet     ierr = MatTranspose(D,MAT_INPLACE_MATRIX,&D);CHKERRQ(ierr);
145c7a4214aSPierre Jolivet     ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&AT);CHKERRQ(ierr);
146c7a4214aSPierre Jolivet     ierr = MatMultEqual(AT,D,10,&flg);CHKERRQ(ierr);
147*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx");
148c7a4214aSPierre Jolivet     ierr = MatTranspose(A,MAT_REUSE_MATRIX,&AT);CHKERRQ(ierr);
149c7a4214aSPierre Jolivet     ierr = MatMultEqual(AT,D,10,&flg);CHKERRQ(ierr);
150*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx");
151c7a4214aSPierre Jolivet     ierr = MatAXPY(D,-1.0,AT,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
152c7a4214aSPierre Jolivet     ierr = MatNorm(D,NORM_INFINITY,&norm);CHKERRQ(ierr);
153*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(PetscAbsReal(norm) > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A-D|| = %g (> %g)",(double)norm,(double)PETSC_SMALL);
154c7a4214aSPierre Jolivet     ierr = MatDestroy(&AT);CHKERRQ(ierr);
155c7a4214aSPierre Jolivet     ierr = MatDestroy(&D);CHKERRQ(ierr);
156c7a4214aSPierre Jolivet     ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P);CHKERRQ(ierr);
157c7a4214aSPierre Jolivet     ierr = MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
158c7a4214aSPierre Jolivet     ierr = MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
159c7a4214aSPierre Jolivet     ierr = MatMatMultEqual(A,B,P,10,&flg);CHKERRQ(ierr);
160*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ABx != Px");
1618b8ddd11SPierre Jolivet     ierr = MatTransposeMatMultEqual(A,B,P,10,&flg);CHKERRQ(ierr);
162*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^TBx != P^Tx");
163c7a4214aSPierre Jolivet     ierr = MatDestroy(&B);CHKERRQ(ierr);
164c7a4214aSPierre Jolivet     ierr = MatDestroy(&P);CHKERRQ(ierr);
165c7a4214aSPierre Jolivet     if (n) {
166c7a4214aSPierre Jolivet       ierr = PetscMalloc1(n*dim,&scoords);CHKERRQ(ierr);
167c7a4214aSPierre Jolivet       ierr = PetscRandomGetValuesReal(rdm,n*dim,scoords);CHKERRQ(ierr);
168c7a4214aSPierre Jolivet       N = n;
169c7a4214aSPierre Jolivet       ierr = MPIU_Allreduce(MPI_IN_PLACE,&N,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr);
170c7a4214aSPierre Jolivet       ierr = PetscCalloc1(N*dim,&gscoords);CHKERRQ(ierr);
171c7a4214aSPierre Jolivet       ierr = MPI_Exscan(&n,&begin,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr);
172c7a4214aSPierre Jolivet       ierr = PetscArraycpy(gscoords+begin*dim,scoords,n*dim);CHKERRQ(ierr);
173c7a4214aSPierre Jolivet       ierr = MPIU_Allreduce(MPI_IN_PLACE,gscoords,N*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr);
17498e73e17SPierre Jolivet       kernel = GenEntriesRectangular;
175c7a4214aSPierre Jolivet       ctx[0] = gcoords;
176c7a4214aSPierre Jolivet       ctx[1] = gscoords;
177c7a4214aSPierre Jolivet       ierr = MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,n,M,N,dim,coords,scoords,kernel,ctx,&R);CHKERRQ(ierr);
178c7a4214aSPierre Jolivet       ierr = MatSetFromOptions(R);CHKERRQ(ierr);
179c7a4214aSPierre Jolivet       ierr = MatAssemblyBegin(R,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
180c7a4214aSPierre Jolivet       ierr = MatAssemblyEnd(R,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
181c7a4214aSPierre Jolivet       ierr = MatViewFromOptions(R,NULL,"-R_view");CHKERRQ(ierr);
182c7a4214aSPierre Jolivet       ierr = MatConvert(R,MATDENSE,MAT_INITIAL_MATRIX,&D);CHKERRQ(ierr);
183c7a4214aSPierre Jolivet       ierr = MatViewFromOptions(D,NULL,"-D_view");CHKERRQ(ierr);
184c7a4214aSPierre Jolivet       ierr = MatMultEqual(R,D,10,&flg);CHKERRQ(ierr);
185*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Rx != Dx");
186c7a4214aSPierre Jolivet       ierr = MatTranspose(D,MAT_INPLACE_MATRIX,&D);CHKERRQ(ierr);
187c7a4214aSPierre Jolivet       ierr = MatTranspose(R,MAT_INITIAL_MATRIX,&RT);CHKERRQ(ierr);
188c7a4214aSPierre Jolivet       ierr = MatMultEqual(RT,D,10,&flg);CHKERRQ(ierr);
189*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx");
190c7a4214aSPierre Jolivet       ierr = MatTranspose(R,MAT_REUSE_MATRIX,&RT);CHKERRQ(ierr);
191c7a4214aSPierre Jolivet       ierr = MatMultEqual(RT,D,10,&flg);CHKERRQ(ierr);
192*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx");
193c7a4214aSPierre Jolivet       ierr = MatDestroy(&RT);CHKERRQ(ierr);
194c7a4214aSPierre Jolivet       ierr = MatDestroy(&D);CHKERRQ(ierr);
195c7a4214aSPierre Jolivet       ierr = MatCreateDense(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_DETERMINE,K,NULL,&B);CHKERRQ(ierr);
196c7a4214aSPierre Jolivet       ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
197c7a4214aSPierre Jolivet       ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
198c7a4214aSPierre Jolivet       ierr = MatSetRandom(B,rdm);CHKERRQ(ierr);
199c7a4214aSPierre Jolivet       ierr = MatMatMult(R,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P);CHKERRQ(ierr);
200c7a4214aSPierre Jolivet       ierr = MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
201c7a4214aSPierre Jolivet       ierr = MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
202c7a4214aSPierre Jolivet       ierr = MatMatMultEqual(R,B,P,10,&flg);CHKERRQ(ierr);
203*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"RBx != Px");
204c7a4214aSPierre Jolivet       ierr = MatDestroy(&B);CHKERRQ(ierr);
205c7a4214aSPierre Jolivet       ierr = MatDestroy(&P);CHKERRQ(ierr);
20698e73e17SPierre Jolivet       ierr = MatCreateVecs(R,&right,&left);CHKERRQ(ierr);
20798e73e17SPierre Jolivet       ierr = VecSetRandom(right,rdm);CHKERRQ(ierr);
20898e73e17SPierre Jolivet       ierr = MatMult(R,right,left);CHKERRQ(ierr);
20998e73e17SPierre Jolivet       ierr = MatHtoolGetPermutationSource(R,&iss);CHKERRQ(ierr);
21098e73e17SPierre Jolivet       ierr = MatHtoolGetPermutationTarget(R,&ist);CHKERRQ(ierr);
21198e73e17SPierre Jolivet       ierr = VecDuplicate(left,&perm);CHKERRQ(ierr);
21298e73e17SPierre Jolivet       ierr = VecCopy(left,perm);CHKERRQ(ierr);
21398e73e17SPierre Jolivet       ierr = VecPermute(perm,ist,PETSC_FALSE);CHKERRQ(ierr);
21498e73e17SPierre Jolivet       ierr = VecPermute(right,iss,PETSC_FALSE);CHKERRQ(ierr);
21598e73e17SPierre Jolivet       ierr = MatHtoolUsePermutation(R,PETSC_FALSE);CHKERRQ(ierr);
21698e73e17SPierre Jolivet       ierr = MatMult(R,right,left);CHKERRQ(ierr);
21798e73e17SPierre Jolivet       ierr = VecAXPY(left,-1.0,perm);CHKERRQ(ierr);
21898e73e17SPierre Jolivet       ierr = VecNorm(left,NORM_INFINITY,&norm);CHKERRQ(ierr);
219*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(PetscAbsReal(norm) > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||y(with permutation)-y(without permutation)|| = %g (> %g)",(double)PetscAbsReal(norm),(double)PETSC_SMALL);
22098e73e17SPierre Jolivet       ierr = MatHtoolUsePermutation(R,PETSC_TRUE);CHKERRQ(ierr);
22198e73e17SPierre Jolivet       ierr = VecDestroy(&perm);CHKERRQ(ierr);
22298e73e17SPierre Jolivet       ierr = VecDestroy(&left);CHKERRQ(ierr);
22398e73e17SPierre Jolivet       ierr = VecDestroy(&right);CHKERRQ(ierr);
22498e73e17SPierre Jolivet       ierr = ISDestroy(&ist);CHKERRQ(ierr);
22598e73e17SPierre Jolivet       ierr = ISDestroy(&iss);CHKERRQ(ierr);
226c7a4214aSPierre Jolivet       ierr = MatDestroy(&R);CHKERRQ(ierr);
227c7a4214aSPierre Jolivet       ierr = PetscFree(gscoords);CHKERRQ(ierr);
228c7a4214aSPierre Jolivet       ierr = PetscFree(scoords);CHKERRQ(ierr);
229c7a4214aSPierre Jolivet     }
230c7a4214aSPierre Jolivet   }
231c7a4214aSPierre Jolivet   ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr);
232c7a4214aSPierre Jolivet   ierr = MatDestroy(&A);CHKERRQ(ierr);
233c7a4214aSPierre Jolivet   ierr = PetscFree(gcoords);CHKERRQ(ierr);
234c7a4214aSPierre Jolivet   ierr = PetscFree(coords);CHKERRQ(ierr);
235c7a4214aSPierre Jolivet   ierr = PetscFinalize();
236c7a4214aSPierre Jolivet   return ierr;
237c7a4214aSPierre Jolivet }
238c7a4214aSPierre Jolivet 
239c7a4214aSPierre Jolivet /*TEST
240c7a4214aSPierre Jolivet 
241c7a4214aSPierre Jolivet    build:
242c7a4214aSPierre Jolivet       requires: htool
243c7a4214aSPierre Jolivet 
244c7a4214aSPierre Jolivet    test:
245c7a4214aSPierre Jolivet       requires: htool
246c7a4214aSPierre Jolivet       suffix: 1
247c7a4214aSPierre Jolivet       nsize: 4
248c7a4214aSPierre Jolivet       args: -m_local 80 -n_local 25 -mat_htool_epsilon 1.0e-11 -symmetric {{false true}shared output}
249c7a4214aSPierre Jolivet       output_file: output/ex101.out
250c7a4214aSPierre Jolivet 
251c7a4214aSPierre Jolivet    test:
252c7a4214aSPierre Jolivet       requires: htool
253c7a4214aSPierre Jolivet       suffix: 2
254c7a4214aSPierre Jolivet       nsize: 4
25598e73e17SPierre Jolivet       args: -m_local 120 -mat_htool_epsilon 1.0e-2 -mat_htool_compressor {{sympartialACA fullACA SVD}shared output} -mat_htool_clustering {{PCARegular PCAGeometric BoundingBox1Regular BoundingBox1Geometric}shared output}
256c7a4214aSPierre Jolivet       output_file: output/ex101.out
257c7a4214aSPierre Jolivet 
258c7a4214aSPierre Jolivet TEST*/
259