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; 28ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 29c4762a1bSJed Brown 30c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 33c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 34c4762a1bSJed Brown ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr); 35c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 36c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 37c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&r);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Set local matrix entries */ 41c4762a1bSJed Brown ierr = MatGetOwnershipIS(A,&isrows,&iscols);CHKERRQ(ierr); 42c4762a1bSJed Brown ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr); 43c4762a1bSJed Brown ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr); 44c4762a1bSJed Brown ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr); 45c4762a1bSJed Brown ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = PetscMalloc1(nrows*ncols,&v);CHKERRQ(ierr); 47c4762a1bSJed Brown for (i=0; i<nrows; i++) { 48c4762a1bSJed Brown for (j=0; j<ncols; j++) { 49c4762a1bSJed Brown ierr = PetscRandomGetValue(r,&rval);CHKERRQ(ierr); 50c4762a1bSJed Brown v[i*ncols+j] = rval; 51c4762a1bSJed Brown } 52c4762a1bSJed Brown } 53c4762a1bSJed Brown ierr = MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 56c4762a1bSJed Brown ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr); 57c4762a1bSJed Brown ierr = ISRestoreIndices(iscols,&cols);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = ISDestroy(&isrows);CHKERRQ(ierr); 59c4762a1bSJed Brown ierr = ISDestroy(&iscols);CHKERRQ(ierr); 60c4762a1bSJed Brown ierr = PetscRandomDestroy(&r);CHKERRQ(ierr); 61c4762a1bSJed Brown 62c20d7725SJed Brown ierr = PetscObjectTypeCompare((PetscObject)A,MATELEMENTAL,&Aiselemental);CHKERRQ(ierr); 63c20d7725SJed Brown 64c20d7725SJed Brown /* Test MatCreateTranspose() and MatTranspose() */ 65c4762a1bSJed Brown ierr = MatCreateTranspose(A,&C);CHKERRQ(ierr); 66c4762a1bSJed Brown ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); /* B = A^T */ 67c4762a1bSJed Brown ierr = MatMultEqual(C,B,10,&equal);CHKERRQ(ierr); 68*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"A^T*x != (x^T*A)^T"); 69c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 70c20d7725SJed Brown 71c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 72c20d7725SJed Brown if (!Aiselemental) { 73c4762a1bSJed Brown ierr = MatTranspose(B,MAT_INPLACE_MATRIX,&B);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = MatMultEqual(C,B,10,&equal);CHKERRQ(ierr); 75*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"C*x != B*x"); 76c20d7725SJed Brown } 77c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 78c20d7725SJed Brown 79c20d7725SJed Brown /* Test B = C*A for matrix type transpose and seqdense */ 80c20d7725SJed Brown if (size == 1 && !Aiselemental) { 81c20d7725SJed Brown ierr = MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&B);CHKERRQ(ierr); 82c20d7725SJed Brown ierr = MatMatMultEqual(C,A,B,10,&equal);CHKERRQ(ierr); 83*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B != C*A for matrix type transpose and seqdense"); 84c20d7725SJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 85c20d7725SJed Brown } 86c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* Test MatMatMult() */ 89c4762a1bSJed Brown if (Test_MatMatMult) { 90c4762a1bSJed Brown #if !defined(PETSC_HAVE_ELEMENTAL) 91*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(size > 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"This test requires ELEMENTAL"); 92c4762a1bSJed Brown #endif 93c4762a1bSJed Brown ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); /* B = A^T */ 94c4762a1bSJed Brown ierr = MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); /* C = B*A = A^T*A */ 95c4762a1bSJed Brown ierr = MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr); 96c20d7725SJed Brown ierr = MatMatMultEqual(B,A,C,10,&equal);CHKERRQ(ierr); 97*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B*A*x != C*x"); 98c4762a1bSJed Brown 99c4762a1bSJed Brown /* Test MatDuplicate for matrix product */ 100c4762a1bSJed Brown ierr = MatDuplicate(C,MAT_COPY_VALUES,&D);CHKERRQ(ierr); 101c20d7725SJed Brown 102c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 103c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 104c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 105c4762a1bSJed Brown } 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* Test MatTransposeMatMult() */ 108c20d7725SJed Brown if (!Aiselemental) { 109c4762a1bSJed Brown ierr = MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); /* D = A^T*A */ 110c4762a1bSJed Brown ierr = MatTransposeMatMult(A,A,MAT_REUSE_MATRIX,fill,&D);CHKERRQ(ierr); 111c20d7725SJed Brown ierr = MatTransposeMatMultEqual(A,A,D,10,&equal);CHKERRQ(ierr); 112*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x"); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* Test MatDuplicate for matrix product */ 115c4762a1bSJed Brown ierr = MatDuplicate(D,MAT_COPY_VALUES,&C);CHKERRQ(ierr); 116c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 117c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* Test D*x = A^T*C*A*x, where C is in AIJ format */ 120c4762a1bSJed Brown ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); 122c4762a1bSJed Brown if (size == 1) { 123c4762a1bSJed Brown ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,am,am);CHKERRQ(ierr); 124c4762a1bSJed Brown } else { 125c4762a1bSJed Brown ierr = MatSetSizes(C,am,am,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown ierr = MatSetFromOptions(C);CHKERRQ(ierr); 128c4762a1bSJed Brown ierr = MatSetUp(C);CHKERRQ(ierr); 129c4762a1bSJed Brown ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 130c4762a1bSJed Brown v[0] = 1.0; 131c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 132c4762a1bSJed Brown ierr = MatSetValues(C,1,&i,1,&i,v,INSERT_VALUES);CHKERRQ(ierr); 133c4762a1bSJed Brown } 134c4762a1bSJed Brown ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 135c4762a1bSJed Brown ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* B = C*A, D = A^T*B */ 138c4762a1bSJed Brown ierr = MatMatMult(C,A,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 139c4762a1bSJed Brown ierr = MatTransposeMatMult(A,B,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = MatTransposeMatMultEqual(A,B,D,10,&equal);CHKERRQ(ierr); 141*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*B*x"); 142c4762a1bSJed Brown 143c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 144c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 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 153c4762a1bSJed Brown ierr = MatGetLocalSize(A, &am, &an);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = MatGetSize(A, &aM, &aN);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = MatCreateDense(PetscObjectComm((PetscObject)A),PETSC_DECIDE, an, aM + 10, aN, NULL, &B);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = MatSetRandom(B, NULL);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 159c4762a1bSJed Brown ierr = MatMatTransposeMult(A,B,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); /* D = A*A^T */ 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* Test MatDuplicate for matrix product */ 162c4762a1bSJed Brown ierr = MatDuplicate(D,MAT_COPY_VALUES,&C);CHKERRQ(ierr); 163c4762a1bSJed Brown 164c4762a1bSJed Brown ierr = MatMatTransposeMult(A,B,MAT_REUSE_MATRIX,fill,&D);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = MatAXPY(C, -1., D, SAME_NONZERO_PATTERN);CHKERRQ(ierr); 166c4762a1bSJed Brown 167c4762a1bSJed Brown ierr = MatNorm(C, NORM_FROBENIUS, &diff);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = MatNorm(D, NORM_FROBENIUS, &scale);CHKERRQ(ierr); 169*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(diff > PETSC_SMALL * scale,PetscObjectComm((PetscObject)D), PETSC_ERR_PLIB, "MatMatTransposeMult() differs between MAT_INITIAL_MATRIX and MAT_REUSE_MATRIX"); 170c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 171c4762a1bSJed Brown 172c4762a1bSJed Brown ierr = MatMatTransposeMultEqual(A,B,D,10,&equal);CHKERRQ(ierr); 173*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x"); 174c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 175c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 176c4762a1bSJed Brown 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 180c4762a1bSJed Brown ierr = PetscFree(v);CHKERRQ(ierr); 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