1c4762a1bSJed Brown static char help[] = "Test MatTransposeColoring for SeqAIJ matrices. Used for '-matmattransmult_color' on MatMatTransposeMult \n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown #include <petsc/private/matimpl.h> /* Need struct _p_MatTransposeColoring for this test. */ 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **argv) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat A,R,C,C_dense,C_sparse,Rt_dense,P,PtAP; 9c4762a1bSJed Brown PetscInt row,col,m,n; 10c4762a1bSJed Brown PetscErrorCode ierr; 11c4762a1bSJed Brown MatScalar one =1.0,val; 12c4762a1bSJed Brown MatColoring mc; 13c4762a1bSJed Brown MatTransposeColoring matcoloring = 0; 14c4762a1bSJed Brown ISColoring iscoloring; 15c4762a1bSJed Brown PetscBool equal; 16c4762a1bSJed Brown PetscMPIInt size; 17c4762a1bSJed Brown 18c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 19*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 202c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* Create seqaij A */ 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&A)); 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,4,4,4,4)); 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATSEQAIJ)); 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 28*5f80ce2aSJacob Faibussowitsch row = 0; col=0; val=1.0; CHKERRQ(MatSetValues(A,1,&row,1,&col,&val,ADD_VALUES)); 29*5f80ce2aSJacob Faibussowitsch row = 1; col=3; val=2.0; CHKERRQ(MatSetValues(A,1,&row,1,&col,&val,ADD_VALUES)); 30*5f80ce2aSJacob Faibussowitsch row = 2; col=2; val=3.0; CHKERRQ(MatSetValues(A,1,&row,1,&col,&val,ADD_VALUES)); 31*5f80ce2aSJacob Faibussowitsch row = 3; col=0; val=4.0; CHKERRQ(MatSetValues(A,1,&row,1,&col,&val,ADD_VALUES)); 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(A,"A_")); 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n")); 37c4762a1bSJed Brown 38c4762a1bSJed Brown /* Create seqaij R */ 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&R)); 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(R,2,4,2,4)); 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(R,MATSEQAIJ)); 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(R)); 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(R)); 44*5f80ce2aSJacob Faibussowitsch row = 0; col=0; CHKERRQ(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES)); 45*5f80ce2aSJacob Faibussowitsch row = 0; col=1; CHKERRQ(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES)); 46c4762a1bSJed Brown 47*5f80ce2aSJacob Faibussowitsch row = 1; col=1; CHKERRQ(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES)); 48*5f80ce2aSJacob Faibussowitsch row = 1; col=2; CHKERRQ(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES)); 49*5f80ce2aSJacob Faibussowitsch row = 1; col=3; CHKERRQ(MatSetValues(R,1,&row,1,&col,&one,ADD_VALUES)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(R,MAT_FINAL_ASSEMBLY)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(R,MAT_FINAL_ASSEMBLY)); 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(R,"R_")); 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(R,PETSC_VIEWER_STDOUT_WORLD)); 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n")); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* C = A*R^T */ 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMult(A,R,MAT_INITIAL_MATRIX,2.0,&C)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(C,"ARt_")); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n")); 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* Create MatTransposeColoring from symbolic C=A*R^T */ 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatColoringCreate(C,&mc)); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatColoringSetDistance(mc,2)); 65*5f80ce2aSJacob Faibussowitsch /* CHKERRQ(MatColoringSetType(mc,MATCOLORINGSL)); */ 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatColoringSetFromOptions(mc)); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatColoringApply(mc,&iscoloring)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatColoringDestroy(&mc)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeColoringCreate(C,iscoloring,&matcoloring)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISColoringDestroy(&iscoloring)); 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* Create Rt_dense */ 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&Rt_dense)); 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(Rt_dense,4,matcoloring->ncolors,PETSC_DECIDE,PETSC_DECIDE)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(Rt_dense,MATDENSE)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqDenseSetPreallocation(Rt_dense,NULL)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(Rt_dense,MAT_FINAL_ASSEMBLY)); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(Rt_dense,MAT_FINAL_ASSEMBLY)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(Rt_dense,&m,&n)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Rt_dense: %" PetscInt_FMT ",%" PetscInt_FMT "\n",m,n)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* Get Rt_dense by Apply MatTransposeColoring to R */ 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransColoringApplySpToDen(matcoloring,R,Rt_dense)); 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* C_dense = A*Rt_dense */ 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A,Rt_dense,MAT_INITIAL_MATRIX,2.0,&C_dense)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(C_dense,"ARt_dense_")); 88*5f80ce2aSJacob Faibussowitsch /*CHKERRQ(MatView(C_dense,PETSC_VIEWER_STDOUT_WORLD)); */ 89*5f80ce2aSJacob Faibussowitsch /*CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n")); */ 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* Recover C from C_dense */ 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&C_sparse)); 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransColoringApplyDenToSp(matcoloring,C_dense,C_sparse)); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(C_sparse,"ARt_color_")); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C_sparse,PETSC_VIEWER_STDOUT_WORLD)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n")); 97c4762a1bSJed Brown 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C_dense)); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C_sparse)); 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Rt_dense)); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeColoringDestroy(&matcoloring)); 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* Test PtAP = P^T*A*P, P = R^T */ 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(R,MAT_INITIAL_MATRIX,&P)); 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A,P,MAT_INITIAL_MATRIX,2.0,&PtAP)); 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(PtAP,"PtAP_")); 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(PtAP,PETSC_VIEWER_STDOUT_WORLD)); 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&P)); 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* Test C = RARt */ 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatRARt(A,R,MAT_INITIAL_MATRIX,2.0,&C)); 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatRARt(A,R,MAT_REUSE_MATRIX,2.0,&C)); 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(C,PtAP,&equal)); 1152c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: PtAP != RARt"); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* Free spaces */ 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&R)); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&PtAP)); 122c4762a1bSJed Brown ierr = PetscFinalize(); 123c4762a1bSJed Brown return ierr; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126c4762a1bSJed Brown /*TEST 127c4762a1bSJed Brown 128c4762a1bSJed Brown test: 129c4762a1bSJed Brown output_file: output/ex161.out 130c4762a1bSJed Brown test: 131c4762a1bSJed Brown suffix: 2 132c20d7725SJed Brown args: -matmattransmult_via color 133c4762a1bSJed Brown output_file: output/ex161.out 134c4762a1bSJed Brown 135c4762a1bSJed Brown test: 136c4762a1bSJed Brown suffix: 3 137c20d7725SJed Brown args: -matmattransmult_via color -matden2sp_brows 3 138c4762a1bSJed Brown output_file: output/ex161.out 139c4762a1bSJed Brown 140c4762a1bSJed Brown test: 141c4762a1bSJed Brown suffix: 4 142c20d7725SJed Brown args: -matmattransmult_via color -matrart_via r*art 143c4762a1bSJed Brown output_file: output/ex161.out 144c4762a1bSJed Brown 145c4762a1bSJed Brown test: 146c4762a1bSJed Brown suffix: 5 147c20d7725SJed Brown args: -matmattransmult_via color -matrart_via coloring_rart 148c4762a1bSJed Brown output_file: output/ex161.out 149c4762a1bSJed Brown 150c4762a1bSJed Brown TEST*/ 151