1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatCreateComposite()\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown Include "petscmat.h" so that we can use matrices. 6c4762a1bSJed Brown automatically includes: 7c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 8c4762a1bSJed Brown petscmat.h - matrices 9c4762a1bSJed Brown petscis.h - index sets petscviewer.h - viewers 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown #include <petscmat.h> 12c4762a1bSJed Brown 13c4762a1bSJed Brown int main(int argc,char **args) 14c4762a1bSJed Brown { 15c4762a1bSJed Brown Mat *A,B; /* matrix */ 16c4762a1bSJed Brown Vec x,y,v,v2,z,z2; 17c4762a1bSJed Brown PetscReal rnorm; 18c4762a1bSJed Brown PetscInt n = 20; /* size of the matrix */ 19c4762a1bSJed Brown PetscInt nmat = 3; /* number of matrices */ 20c4762a1bSJed Brown PetscInt i; 21c4762a1bSJed Brown PetscRandom rctx; 22c4762a1bSJed Brown MatCompositeType type; 23c4762a1bSJed Brown PetscScalar scalings[5]={2,3,4,5,6}; 24c4762a1bSJed Brown 25*327415f7SBarry Smith PetscFunctionBeginUser; 269566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nmat",&nmat,NULL)); 29c4762a1bSJed Brown 30c4762a1bSJed Brown /* 31c4762a1bSJed Brown Create random matrices 32c4762a1bSJed Brown */ 339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nmat+3,&A)); 349566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rctx)); 359566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n/2,3,NULL,3,NULL,&A[0])); 36c4762a1bSJed Brown for (i = 1; i < nmat+1; i++) { 379566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,3,NULL,3,NULL,&A[i])); 38c4762a1bSJed Brown } 399566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n/2,n,3,NULL,3,NULL,&A[nmat+1])); 40c4762a1bSJed Brown for (i = 0; i < nmat+2; i++) { 419566063dSJacob Faibussowitsch PetscCall(MatSetRandom(A[i],rctx)); 42c4762a1bSJed Brown } 43c4762a1bSJed Brown 449566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A[1],&x,&y)); 459566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y,&z)); 469566063dSJacob Faibussowitsch PetscCall(VecDuplicate(z,&z2)); 479566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A[0],&v,NULL)); 489566063dSJacob Faibussowitsch PetscCall(VecDuplicate(v,&v2)); 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* Test MatMult of an ADDITIVE MatComposite B made up of A[1],A[2],A[3] with separate scalings */ 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Do MatMult with A[1],A[2],A[3] by hand and store the result in z */ 539566063dSJacob Faibussowitsch PetscCall(VecSet(x,1.0)); 549566063dSJacob Faibussowitsch PetscCall(MatMult(A[1],x,z)); 559566063dSJacob Faibussowitsch PetscCall(VecScale(z,scalings[1])); 56c4762a1bSJed Brown for (i = 2; i < nmat+1; i++) { 579566063dSJacob Faibussowitsch PetscCall(MatMult(A[i],x,z2)); 589566063dSJacob Faibussowitsch PetscCall(VecAXPY(z,scalings[i],z2)); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* Do MatMult using MatComposite and store the result in y */ 629566063dSJacob Faibussowitsch PetscCall(VecSet(y,0.0)); 639566063dSJacob Faibussowitsch PetscCall(MatCreateComposite(PETSC_COMM_WORLD,nmat,A+1,&B)); 649566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 659566063dSJacob Faibussowitsch PetscCall(MatCompositeSetScalings(B,&scalings[1])); 669566063dSJacob Faibussowitsch PetscCall(MatMultAdd(B,x,y,y)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* Diff y and z */ 699566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,-1.0,z)); 709566063dSJacob Faibussowitsch PetscCall(VecNorm(y,NORM_2,&rnorm)); 71c4762a1bSJed Brown if (rnorm > 10000.0*PETSC_MACHINE_EPSILON) { 729566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with composite add %g\n",(double)rnorm)); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* Test MatCompositeMerge on ADDITIVE MatComposite */ 769566063dSJacob Faibussowitsch PetscCall(MatCompositeSetMatStructure(B,DIFFERENT_NONZERO_PATTERN)); /* default */ 779566063dSJacob Faibussowitsch PetscCall(MatCompositeMerge(B)); 789566063dSJacob Faibussowitsch PetscCall(MatMult(B,x,y)); 799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 809566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,-1.0,z)); 819566063dSJacob Faibussowitsch PetscCall(VecNorm(y,NORM_2,&rnorm)); 82c4762a1bSJed Brown if (rnorm > 10000.0*PETSC_MACHINE_EPSILON) { 839566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with composite add after merge %g\n",(double)rnorm)); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* 87c4762a1bSJed Brown Test n x n/2 multiplicative composite B made up of A[0],A[1],A[2] with separate scalings 88c4762a1bSJed Brown */ 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* Do MatMult with A[0],A[1],A[2] by hand and store the result in z */ 919566063dSJacob Faibussowitsch PetscCall(VecSet(v,1.0)); 929566063dSJacob Faibussowitsch PetscCall(MatMult(A[0],v,z)); 939566063dSJacob Faibussowitsch PetscCall(VecScale(z,scalings[0])); 94c4762a1bSJed Brown for (i = 1; i < nmat; i++) { 959566063dSJacob Faibussowitsch PetscCall(MatMult(A[i],z,y)); 969566063dSJacob Faibussowitsch PetscCall(VecScale(y,scalings[i])); 979566063dSJacob Faibussowitsch PetscCall(VecCopy(y,z)); 98c4762a1bSJed Brown } 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* Do MatMult using MatComposite and store the result in y */ 1019566063dSJacob Faibussowitsch PetscCall(MatCreateComposite(PETSC_COMM_WORLD,nmat,A,&B)); 1029566063dSJacob Faibussowitsch PetscCall(MatCompositeSetType(B,MAT_COMPOSITE_MULTIPLICATIVE)); 1039566063dSJacob Faibussowitsch PetscCall(MatCompositeSetMergeType(B,MAT_COMPOSITE_MERGE_LEFT)); 1049566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 1059566063dSJacob Faibussowitsch PetscCall(MatCompositeSetScalings(B,&scalings[0])); 1069566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1079566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); /* do MatCompositeMerge() if -mat_composite_merge 1 */ 1089566063dSJacob Faibussowitsch PetscCall(MatMult(B,v,y)); 1099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* Diff y and z */ 1129566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,-1.0,z)); 1139566063dSJacob Faibussowitsch PetscCall(VecNorm(y,NORM_2,&rnorm)); 114c4762a1bSJed Brown if (rnorm > 10000.0*PETSC_MACHINE_EPSILON) { 1159566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with composite multiplicative %g\n",(double)rnorm)); 116c4762a1bSJed Brown } 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* 119c4762a1bSJed Brown Test n/2 x n multiplicative composite B made up of A[2], A[3], A[4] without separate scalings 120c4762a1bSJed Brown */ 1219566063dSJacob Faibussowitsch PetscCall(VecSet(x,1.0)); 1229566063dSJacob Faibussowitsch PetscCall(MatMult(A[2],x,z)); 123c4762a1bSJed Brown for (i = 3; i < nmat+1; i++) { 1249566063dSJacob Faibussowitsch PetscCall(MatMult(A[i],z,y)); 1259566063dSJacob Faibussowitsch PetscCall(VecCopy(y,z)); 126c4762a1bSJed Brown } 1279566063dSJacob Faibussowitsch PetscCall(MatMult(A[nmat+1],z,v)); 128c4762a1bSJed Brown 1299566063dSJacob Faibussowitsch PetscCall(MatCreateComposite(PETSC_COMM_WORLD,nmat,A+2,&B)); 1309566063dSJacob Faibussowitsch PetscCall(MatCompositeSetType(B,MAT_COMPOSITE_MULTIPLICATIVE)); 1319566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 1329566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1339566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); /* do MatCompositeMerge() if -mat_composite_merge 1 */ 1349566063dSJacob Faibussowitsch PetscCall(MatMult(B,x,v2)); 1359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 136c4762a1bSJed Brown 1379566063dSJacob Faibussowitsch PetscCall(VecAXPY(v2,-1.0,v)); 1389566063dSJacob Faibussowitsch PetscCall(VecNorm(v2,NORM_2,&rnorm)); 139c4762a1bSJed Brown if (rnorm > 10000.0*PETSC_MACHINE_EPSILON) { 1409566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with composite multiplicative %g\n",(double)rnorm)); 141c4762a1bSJed Brown } 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* 144c4762a1bSJed Brown Test get functions 145c4762a1bSJed Brown */ 1469566063dSJacob Faibussowitsch PetscCall(MatCreateComposite(PETSC_COMM_WORLD,nmat,A,&B)); 1479566063dSJacob Faibussowitsch PetscCall(MatCompositeGetNumberMat(B,&n)); 148c4762a1bSJed Brown if (nmat != n) { 1499566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with GetNumberMat %" PetscInt_FMT " != %" PetscInt_FMT "\n",nmat,n)); 150c4762a1bSJed Brown } 1519566063dSJacob Faibussowitsch PetscCall(MatCompositeGetMat(B,0,&A[nmat+2])); 152c4762a1bSJed Brown if (A[0] != A[nmat+2]) { 1539566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with GetMat\n")); 154c4762a1bSJed Brown } 1559566063dSJacob Faibussowitsch PetscCall(MatCompositeGetType(B,&type)); 156c4762a1bSJed Brown if (type != MAT_COMPOSITE_ADDITIVE) { 1579566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with GetType\n")); 158c4762a1bSJed Brown } 1599566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* 162c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 163c4762a1bSJed Brown are no longer needed. 164c4762a1bSJed Brown */ 1659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 1689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v2)); 1699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&z)); 1709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&z2)); 1719566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 172c4762a1bSJed Brown for (i = 0; i < nmat+2; i++) { 1739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A[i])); 174c4762a1bSJed Brown } 1759566063dSJacob Faibussowitsch PetscCall(PetscFree(A)); 176c4762a1bSJed Brown 1779566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 178b122ec5aSJacob Faibussowitsch return 0; 179c4762a1bSJed Brown } 180c4762a1bSJed Brown 181c4762a1bSJed Brown /*TEST 182c4762a1bSJed Brown 183c4762a1bSJed Brown test: 184c4762a1bSJed Brown nsize: 2 185c4762a1bSJed Brown requires: double 186c4762a1bSJed Brown args: -mat_composite_merge {{0 1}shared output} -mat_composite_merge_mvctx {{0 1}shared output} 187c4762a1bSJed Brown 188c4762a1bSJed Brown TEST*/ 189