1c4762a1bSJed Brown static char help[] = "Test MatMatMult(), MatTranspose(), MatTransposeMatMult() for Dense and Elemental matrices.\n\n"; 2c4762a1bSJed Brown /* 3c4762a1bSJed Brown Example: 4c4762a1bSJed Brown mpiexec -n <np> ./ex104 -mat_type elemental 5c4762a1bSJed Brown */ 6c4762a1bSJed Brown 7c4762a1bSJed Brown #include <petscmat.h> 8c4762a1bSJed Brown 9c4762a1bSJed Brown int main(int argc,char **argv) 10c4762a1bSJed Brown { 11c4762a1bSJed Brown Mat A,B,C,D; 12c4762a1bSJed Brown PetscInt i,M=10,N=5,j,nrows,ncols,am,an,rstart,rend; 13c4762a1bSJed Brown PetscErrorCode ierr; 14c4762a1bSJed Brown PetscRandom r; 15c20d7725SJed Brown PetscBool equal,Aiselemental; 16c4762a1bSJed Brown PetscReal fill = 1.0; 17c4762a1bSJed Brown IS isrows,iscols; 18c4762a1bSJed Brown const PetscInt *rows,*cols; 19c4762a1bSJed Brown PetscScalar *v,rval; 20c4762a1bSJed Brown #if defined(PETSC_HAVE_ELEMENTAL) 21c4762a1bSJed Brown PetscBool Test_MatMatMult=PETSC_TRUE; 22c4762a1bSJed Brown #else 23c4762a1bSJed Brown PetscBool Test_MatMatMult=PETSC_FALSE; 24c4762a1bSJed Brown #endif 25c4762a1bSJed Brown PetscMPIInt size; 26c4762a1bSJed Brown 27c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 28*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 29c4762a1bSJed Brown 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N)); 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATDENSE)); 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&r)); 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(r)); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Set local matrix entries */ 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipIS(A,&isrows,&iscols)); 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(isrows,&nrows)); 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isrows,&rows)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(iscols,&ncols)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(iscols,&cols)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nrows*ncols,&v)); 47c4762a1bSJed Brown for (i=0; i<nrows; i++) { 48c4762a1bSJed Brown for (j=0; j<ncols; j++) { 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(r,&rval)); 50c4762a1bSJed Brown v[i*ncols+j] = rval; 51c4762a1bSJed Brown } 52c4762a1bSJed Brown } 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES)); 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isrows,&rows)); 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(iscols,&cols)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrows)); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscols)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&r)); 61c4762a1bSJed Brown 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATELEMENTAL,&Aiselemental)); 63c20d7725SJed Brown 64c20d7725SJed Brown /* Test MatCreateTranspose() and MatTranspose() */ 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateTranspose(A,&C)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX,&B)); /* B = A^T */ 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(C,B,10,&equal)); 682c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"A^T*x != (x^T*A)^T"); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 70c20d7725SJed Brown 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&B)); 72c20d7725SJed Brown if (!Aiselemental) { 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(B,MAT_INPLACE_MATRIX,&B)); 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(C,B,10,&equal)); 752c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"C*x != B*x"); 76c20d7725SJed Brown } 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 78c20d7725SJed Brown 79c20d7725SJed Brown /* Test B = C*A for matrix type transpose and seqdense */ 80c20d7725SJed Brown if (size == 1 && !Aiselemental) { 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&B)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(C,A,B,10,&equal)); 832c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B != C*A for matrix type transpose and seqdense"); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 85c20d7725SJed Brown } 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* Test MatMatMult() */ 89c4762a1bSJed Brown if (Test_MatMatMult) { 90c4762a1bSJed Brown #if !defined(PETSC_HAVE_ELEMENTAL) 912c71b3e2SJacob Faibussowitsch PetscCheckFalse(size > 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"This test requires ELEMENTAL"); 92c4762a1bSJed Brown #endif 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX,&B)); /* B = A^T */ 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C)); /* C = B*A = A^T*A */ 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(B,A,C,10,&equal)); 972c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B*A*x != C*x"); 98c4762a1bSJed Brown 99c4762a1bSJed Brown /* Test MatDuplicate for matrix product */ 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&D)); 101c20d7725SJed Brown 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&D)); 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 105c4762a1bSJed Brown } 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* Test MatTransposeMatMult() */ 108c20d7725SJed Brown if (!Aiselemental) { 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&D)); /* D = A^T*A */ 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(A,A,MAT_REUSE_MATRIX,fill,&D)); 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMultEqual(A,A,D,10,&equal)); 1122c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x"); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* Test MatDuplicate for matrix product */ 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(D,MAT_COPY_VALUES,&C)); 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&D)); 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* Test D*x = A^T*C*A*x, where C is in AIJ format */ 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&am,&an)); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 122c4762a1bSJed Brown if (size == 1) { 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,am,am)); 124c4762a1bSJed Brown } else { 125*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,am,am,PETSC_DECIDE,PETSC_DECIDE)); 126c4762a1bSJed Brown } 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(C,&rstart,&rend)); 130c4762a1bSJed Brown v[0] = 1.0; 131c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 132*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&i,1,&i,v,INSERT_VALUES)); 133c4762a1bSJed Brown } 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* B = C*A, D = A^T*B */ 138*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(C,A,MAT_INITIAL_MATRIX,1.0,&B)); 139*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(A,B,MAT_INITIAL_MATRIX,fill,&D)); 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMultEqual(A,B,D,10,&equal)); 1412c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*B*x"); 142c4762a1bSJed Brown 143*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&D)); 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 146c4762a1bSJed Brown } 147c4762a1bSJed Brown 148c4762a1bSJed Brown /* Test MatMatTransposeMult() */ 149c20d7725SJed Brown if (!Aiselemental) { 150c4762a1bSJed Brown PetscReal diff, scale; 151c4762a1bSJed Brown PetscInt am, an, aM, aN; 152c4762a1bSJed Brown 153*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A, &am, &an)); 154*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A, &aM, &aN)); 155*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PetscObjectComm((PetscObject)A),PETSC_DECIDE, an, aM + 10, aN, NULL, &B)); 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(B, NULL)); 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 158*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMult(A,B,MAT_INITIAL_MATRIX,fill,&D)); /* D = A*A^T */ 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* Test MatDuplicate for matrix product */ 162*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(D,MAT_COPY_VALUES,&C)); 163c4762a1bSJed Brown 164*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMult(A,B,MAT_REUSE_MATRIX,fill,&D)); 165*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(C, -1., D, SAME_NONZERO_PATTERN)); 166c4762a1bSJed Brown 167*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C, NORM_FROBENIUS, &diff)); 168*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(D, NORM_FROBENIUS, &scale)); 1692c71b3e2SJacob Faibussowitsch PetscCheckFalse(diff > PETSC_SMALL * scale,PetscObjectComm((PetscObject)D), PETSC_ERR_PLIB, "MatMatTransposeMult() differs between MAT_INITIAL_MATRIX and MAT_REUSE_MATRIX"); 170*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 171c4762a1bSJed Brown 172*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMultEqual(A,B,D,10,&equal)); 1732c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x"); 174*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&D)); 175*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 176c4762a1bSJed Brown 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 180*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 181c4762a1bSJed Brown ierr = PetscFinalize(); 182c4762a1bSJed Brown return ierr; 183c4762a1bSJed Brown } 184c4762a1bSJed Brown 185c4762a1bSJed Brown /*TEST 186c4762a1bSJed Brown 187c4762a1bSJed Brown test: 188c4762a1bSJed Brown output_file: output/ex104.out 189c4762a1bSJed Brown 190c4762a1bSJed Brown test: 191c4762a1bSJed Brown suffix: 2 192c4762a1bSJed Brown nsize: 2 193c4762a1bSJed Brown output_file: output/ex104.out 194c4762a1bSJed Brown 195c4762a1bSJed Brown test: 196c4762a1bSJed Brown suffix: 3 197c4762a1bSJed Brown nsize: 4 198c4762a1bSJed Brown output_file: output/ex104.out 199c4762a1bSJed Brown args: -M 23 -N 31 200c4762a1bSJed Brown 201c4762a1bSJed Brown test: 202c4762a1bSJed Brown suffix: 4 203c4762a1bSJed Brown nsize: 4 204c4762a1bSJed Brown output_file: output/ex104.out 205c4762a1bSJed Brown args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via cyclic 206c4762a1bSJed Brown 207c4762a1bSJed Brown test: 208c4762a1bSJed Brown suffix: 5 209c4762a1bSJed Brown nsize: 4 210c4762a1bSJed Brown output_file: output/ex104.out 211c4762a1bSJed Brown args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via allgatherv 212c4762a1bSJed Brown 213c20d7725SJed Brown test: 214c20d7725SJed Brown suffix: 6 215c20d7725SJed Brown args: -mat_type elemental 216c20d7725SJed Brown requires: elemental 217c20d7725SJed Brown output_file: output/ex104.out 218c20d7725SJed Brown 219c20d7725SJed Brown test: 220c20d7725SJed Brown suffix: 7 221c20d7725SJed Brown nsize: 2 222c20d7725SJed Brown args: -mat_type elemental 223c20d7725SJed Brown requires: elemental 224c20d7725SJed Brown output_file: output/ex104.out 225c20d7725SJed Brown 226c4762a1bSJed Brown TEST*/ 227