1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatCreateComposite()\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*T 5c4762a1bSJed Brown Concepts: Mat^composite matrices 6c4762a1bSJed Brown Processors: n 7c4762a1bSJed Brown T*/ 8c4762a1bSJed Brown 9c4762a1bSJed Brown /* 10c4762a1bSJed Brown Include "petscmat.h" so that we can use matrices. 11c4762a1bSJed Brown automatically includes: 12c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 13c4762a1bSJed Brown petscmat.h - matrices 14c4762a1bSJed Brown petscis.h - index sets petscviewer.h - viewers 15c4762a1bSJed Brown */ 16c4762a1bSJed Brown #include <petscmat.h> 17c4762a1bSJed Brown 18c4762a1bSJed Brown int main(int argc,char **args) 19c4762a1bSJed Brown { 20c4762a1bSJed Brown Mat *A,B; /* matrix */ 21c4762a1bSJed Brown Vec x,y,v,v2,z,z2; 22c4762a1bSJed Brown PetscReal rnorm; 23c4762a1bSJed Brown PetscInt n = 20; /* size of the matrix */ 24c4762a1bSJed Brown PetscInt nmat = 3; /* number of matrices */ 25c4762a1bSJed Brown PetscInt i; 26c4762a1bSJed Brown PetscRandom rctx; 27c4762a1bSJed Brown MatCompositeType type; 28c4762a1bSJed Brown PetscScalar scalings[5]={2,3,4,5,6}; 29c4762a1bSJed Brown 30*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nmat",&nmat,NULL)); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* 35c4762a1bSJed Brown Create random matrices 36c4762a1bSJed Brown */ 375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nmat+3,&A)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rctx)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n/2,3,NULL,3,NULL,&A[0])); 40c4762a1bSJed Brown for (i = 1; i < nmat+1; i++) { 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,3,NULL,3,NULL,&A[i])); 42c4762a1bSJed Brown } 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n/2,n,3,NULL,3,NULL,&A[nmat+1])); 44c4762a1bSJed Brown for (i = 0; i < nmat+2; i++) { 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(A[i],rctx)); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A[1],&x,&y)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(y,&z)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(z,&z2)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A[0],&v,NULL)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(v,&v2)); 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* Test MatMult of an ADDITIVE MatComposite B made up of A[1],A[2],A[3] with separate scalings */ 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Do MatMult with A[1],A[2],A[3] by hand and store the result in z */ 575f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(x,1.0)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A[1],x,z)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(z,scalings[1])); 60c4762a1bSJed Brown for (i = 2; i < nmat+1; i++) { 615f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A[i],x,z2)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(z,scalings[i],z2)); 63c4762a1bSJed Brown } 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* Do MatMult using MatComposite and store the result in y */ 665f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(y,0.0)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateComposite(PETSC_COMM_WORLD,nmat,A+1,&B)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeSetScalings(B,&scalings[1])); 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(B,x,y,y)); 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* Diff y and z */ 735f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,-1.0,z)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&rnorm)); 75c4762a1bSJed Brown if (rnorm > 10000.0*PETSC_MACHINE_EPSILON) { 765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error with composite add %g\n",(double)rnorm)); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Test MatCompositeMerge on ADDITIVE MatComposite */ 805f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeSetMatStructure(B,DIFFERENT_NONZERO_PATTERN)); /* default */ 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeMerge(B)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(B,x,y)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,-1.0,z)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&rnorm)); 86c4762a1bSJed Brown if (rnorm > 10000.0*PETSC_MACHINE_EPSILON) { 875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error with composite add after merge %g\n",(double)rnorm)); 88c4762a1bSJed Brown } 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* 91c4762a1bSJed Brown Test n x n/2 multiplicative composite B made up of A[0],A[1],A[2] with separate scalings 92c4762a1bSJed Brown */ 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* Do MatMult with A[0],A[1],A[2] by hand and store the result in z */ 955f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v,1.0)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A[0],v,z)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(z,scalings[0])); 98c4762a1bSJed Brown for (i = 1; i < nmat; i++) { 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A[i],z,y)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(y,scalings[i])); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(y,z)); 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* Do MatMult using MatComposite and store the result in y */ 1055f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateComposite(PETSC_COMM_WORLD,nmat,A,&B)); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeSetType(B,MAT_COMPOSITE_MULTIPLICATIVE)); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeSetMergeType(B,MAT_COMPOSITE_MERGE_LEFT)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeSetScalings(B,&scalings[0])); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); /* do MatCompositeMerge() if -mat_composite_merge 1 */ 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(B,v,y)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* Diff y and z */ 1165f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,-1.0,z)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&rnorm)); 118c4762a1bSJed Brown if (rnorm > 10000.0*PETSC_MACHINE_EPSILON) { 1195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error with composite multiplicative %g\n",(double)rnorm)); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* 123c4762a1bSJed Brown Test n/2 x n multiplicative composite B made up of A[2], A[3], A[4] without separate scalings 124c4762a1bSJed Brown */ 1255f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(x,1.0)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A[2],x,z)); 127c4762a1bSJed Brown for (i = 3; i < nmat+1; i++) { 1285f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A[i],z,y)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(y,z)); 130c4762a1bSJed Brown } 1315f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A[nmat+1],z,v)); 132c4762a1bSJed Brown 1335f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateComposite(PETSC_COMM_WORLD,nmat,A+2,&B)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeSetType(B,MAT_COMPOSITE_MULTIPLICATIVE)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); /* do MatCompositeMerge() if -mat_composite_merge 1 */ 1385f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(B,x,v2)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 140c4762a1bSJed Brown 1415f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(v2,-1.0,v)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(v2,NORM_2,&rnorm)); 143c4762a1bSJed Brown if (rnorm > 10000.0*PETSC_MACHINE_EPSILON) { 1445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error with composite multiplicative %g\n",(double)rnorm)); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* 148c4762a1bSJed Brown Test get functions 149c4762a1bSJed Brown */ 1505f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateComposite(PETSC_COMM_WORLD,nmat,A,&B)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeGetNumberMat(B,&n)); 152c4762a1bSJed Brown if (nmat != n) { 1535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error with GetNumberMat %" PetscInt_FMT " != %" PetscInt_FMT "\n",nmat,n)); 154c4762a1bSJed Brown } 1555f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeGetMat(B,0,&A[nmat+2])); 156c4762a1bSJed Brown if (A[0] != A[nmat+2]) { 1575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error with GetMat\n")); 158c4762a1bSJed Brown } 1595f80ce2aSJacob Faibussowitsch CHKERRQ(MatCompositeGetType(B,&type)); 160c4762a1bSJed Brown if (type != MAT_COMPOSITE_ADDITIVE) { 1615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error with GetType\n")); 162c4762a1bSJed Brown } 1635f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* 166c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 167c4762a1bSJed Brown are no longer needed. 168c4762a1bSJed Brown */ 1695f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 1715f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v)); 1725f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v2)); 1735f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&z)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&z2)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rctx)); 176c4762a1bSJed Brown for (i = 0; i < nmat+2; i++) { 1775f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A[i])); 178c4762a1bSJed Brown } 1795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(A)); 180c4762a1bSJed Brown 181*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 182*b122ec5aSJacob Faibussowitsch return 0; 183c4762a1bSJed Brown } 184c4762a1bSJed Brown 185c4762a1bSJed Brown /*TEST 186c4762a1bSJed Brown 187c4762a1bSJed Brown test: 188c4762a1bSJed Brown nsize: 2 189c4762a1bSJed Brown requires: double 190c4762a1bSJed Brown args: -mat_composite_merge {{0 1}shared output} -mat_composite_merge_mvctx {{0 1}shared output} 191c4762a1bSJed Brown 192c4762a1bSJed Brown TEST*/ 193