1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatCopy() and MatStore/RetrieveValues().\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat C,A; 9c4762a1bSJed Brown PetscInt i, n = 10,midx[3],bs=1; 10c4762a1bSJed Brown PetscScalar v[3]; 11c4762a1bSJed Brown PetscBool flg,isAIJ; 12c4762a1bSJed Brown MatType type; 13c4762a1bSJed Brown PetscMPIInt size; 14c4762a1bSJed Brown 15*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 16*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 17*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 18*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL)); 19c4762a1bSJed Brown 20*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 21*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n)); 22*9566063dSJacob Faibussowitsch PetscCall(MatSetType(C,MATAIJ)); 23*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 24*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)C,"initial")); 25c4762a1bSJed Brown 26*9566063dSJacob Faibussowitsch PetscCall(MatGetType(C,&type)); 27c4762a1bSJed Brown if (size == 1) { 28*9566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)C,MATSEQAIJ,&isAIJ)); 29c4762a1bSJed Brown } else { 30*9566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)C,MATMPIAIJ,&isAIJ)); 31c4762a1bSJed Brown } 32*9566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(C,3,NULL)); 33*9566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(C,3,NULL,3,NULL)); 34*9566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(C,bs,3,NULL)); 35*9566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(C,bs,3,NULL,3,NULL)); 36*9566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(C,bs,3,NULL)); 37*9566063dSJacob Faibussowitsch PetscCall(MatMPISBAIJSetPreallocation(C,bs,3,NULL,3,NULL)); 38c4762a1bSJed Brown 39c4762a1bSJed Brown v[0] = -1.; v[1] = 2.; v[2] = -1.; 40c4762a1bSJed Brown for (i=1; i<n-1; i++) { 41c4762a1bSJed Brown midx[2] = i-1; midx[1] = i; midx[0] = i+1; 42*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&i,3,midx,v,INSERT_VALUES)); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown i = 0; midx[0] = 0; midx[1] = 1; 45c4762a1bSJed Brown v[0] = 2.0; v[1] = -1.; 46*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&i,2,midx,v,INSERT_VALUES)); 47c4762a1bSJed Brown i = n-1; midx[0] = n-2; midx[1] = n-1; 48c4762a1bSJed Brown v[0] = -1.0; v[1] = 2.; 49*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&i,2,midx,v,INSERT_VALUES)); 50c4762a1bSJed Brown 51*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 52*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 53*9566063dSJacob Faibussowitsch PetscCall(MatView(C,NULL)); 54*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(C,NULL,"-view")); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* test matduplicate */ 57*9566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C,MAT_COPY_VALUES,&A)); 58*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"duplicate_copy")); 59*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A,NULL,"-view")); 60*9566063dSJacob Faibussowitsch PetscCall(MatEqual(A,C,&flg)); 6128b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"MatDuplicate(C,MAT_COPY_VALUES,): Matrices are NOT equal"); 62*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* test matrices with different nonzero patterns - Note: A is created with different nonzero pattern of C! */ 65*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 66*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n)); 67*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 68*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 69*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 70*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 71c4762a1bSJed Brown 72*9566063dSJacob Faibussowitsch PetscCall(MatCopy(C,A,DIFFERENT_NONZERO_PATTERN)); 73*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"copy_diffnnz")); 74*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A,NULL,"-view")); 75*9566063dSJacob Faibussowitsch PetscCall(MatEqual(A,C,&flg)); 7628b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"MatCopy(C,A,DIFFERENT_NONZERO_PATTERN): Matrices are NOT equal"); 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* test matrices with same nonzero pattern */ 79*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 80*9566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&A)); 81*9566063dSJacob Faibussowitsch PetscCall(MatCopy(C,A,SAME_NONZERO_PATTERN)); 82*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"copy_samennz")); 83*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A,NULL,"-view")); 84*9566063dSJacob Faibussowitsch PetscCall(MatEqual(A,C,&flg)); 8528b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"MatCopy(C,A,SAME_NONZERO_PATTERN): Matrices are NOT equal"); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* test subset nonzero pattern */ 88*9566063dSJacob Faibussowitsch PetscCall(MatCopy(C,A,SUBSET_NONZERO_PATTERN)); 89*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"copy_subnnz")); 90*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A,NULL,"-view")); 91*9566063dSJacob Faibussowitsch PetscCall(MatEqual(A,C,&flg)); 9228b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"MatCopy(C,A,SUBSET_NONZERO_PATTERN): Matrices are NOT equal"); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* Test MatCopy on a matrix obtained after MatConvert from AIJ 95c4762a1bSJed Brown see https://lists.mcs.anl.gov/pipermail/petsc-dev/2019-April/024289.html */ 96*9566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(C,&flg)); 97c4762a1bSJed Brown if (flg) { 98c4762a1bSJed Brown Mat Cs,Cse; 99c4762a1bSJed Brown MatType Ctype,Cstype; 100c4762a1bSJed Brown 101*9566063dSJacob Faibussowitsch PetscCall(MatGetType(C,&Ctype)); 102*9566063dSJacob Faibussowitsch PetscCall(MatTranspose(C,MAT_INITIAL_MATRIX,&Cs)); 103*9566063dSJacob Faibussowitsch PetscCall(MatAXPY(Cs,1.0,C,DIFFERENT_NONZERO_PATTERN)); 104*9566063dSJacob Faibussowitsch PetscCall(MatConvert(Cs,MATAIJ,MAT_INPLACE_MATRIX,&Cs)); 105*9566063dSJacob Faibussowitsch PetscCall(MatSetOption(Cs,MAT_SYMMETRIC,PETSC_TRUE)); 106*9566063dSJacob Faibussowitsch PetscCall(MatGetType(Cs,&Cstype)); 107c4762a1bSJed Brown 108*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)Cs,"symm_initial")); 109*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Cs,NULL,"-view")); 110c4762a1bSJed Brown 111*9566063dSJacob Faibussowitsch PetscCall(MatConvert(Cs,Ctype,MAT_INITIAL_MATRIX,&Cse)); 112*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)Cse,"symm_conv_init")); 113*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Cse,NULL,"-view")); 114*9566063dSJacob Faibussowitsch PetscCall(MatMultEqual(Cs,Cse,5,&flg)); 11528b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatConvert MAT_INITIAL_MATRIX %s -> %s: Matrices are NOT multequal",Ctype,Cstype); 116c4762a1bSJed Brown 117*9566063dSJacob Faibussowitsch PetscCall(MatConvert(Cs,Ctype,MAT_REUSE_MATRIX,&Cse)); 118*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)Cse,"symm_conv_reuse")); 119*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Cse,NULL,"-view")); 120*9566063dSJacob Faibussowitsch PetscCall(MatMultEqual(Cs,Cse,5,&flg)); 12128b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatConvert MAT_REUSE_MATRIX %s -> %s: Matrices are NOT multequal",Ctype,Cstype); 122c4762a1bSJed Brown 123*9566063dSJacob Faibussowitsch PetscCall(MatCopy(Cs,Cse,SAME_NONZERO_PATTERN)); 124*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)Cse,"symm_conv_copy_samennz")); 125*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Cse,NULL,"-view")); 126*9566063dSJacob Faibussowitsch PetscCall(MatMultEqual(Cs,Cse,5,&flg)); 12728b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatCopy(...SAME_NONZERO_PATTERN) %s -> %s: Matrices are NOT multequal",Ctype,Cstype); 128c4762a1bSJed Brown 129*9566063dSJacob Faibussowitsch PetscCall(MatCopy(Cs,Cse,SUBSET_NONZERO_PATTERN)); 130*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)Cse,"symm_conv_copy_subnnz")); 131*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Cse,NULL,"-view")); 132*9566063dSJacob Faibussowitsch PetscCall(MatMultEqual(Cs,Cse,5,&flg)); 13328b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatCopy(...SUBSET_NONZERO_PATTERN) %s -> %s: Matrices are NOT multequal",Ctype,Cstype); 134c4762a1bSJed Brown 135*9566063dSJacob Faibussowitsch PetscCall(MatCopy(Cs,Cse,DIFFERENT_NONZERO_PATTERN)); 136*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)Cse,"symm_conv_copy_diffnnz")); 137*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(Cse,NULL,"-view")); 138*9566063dSJacob Faibussowitsch PetscCall(MatMultEqual(Cs,Cse,5,&flg)); 13928b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatCopy(...DIFFERENT_NONZERO_PATTERN) %s -> %s: Matrices are NOT multequal",Ctype,Cstype); 140c4762a1bSJed Brown 141*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cse)); 142*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cs)); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* test MatStore/RetrieveValues() */ 146c4762a1bSJed Brown if (isAIJ) { 147*9566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE)); 148*9566063dSJacob Faibussowitsch PetscCall(MatStoreValues(A)); 149*9566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(A)); 150*9566063dSJacob Faibussowitsch PetscCall(MatRetrieveValues(A)); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown 153*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 154*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 155*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 156b122ec5aSJacob Faibussowitsch return 0; 157c4762a1bSJed Brown } 158c4762a1bSJed Brown 159c4762a1bSJed Brown /*TEST 160c4762a1bSJed Brown 161c4762a1bSJed Brown testset: 162c4762a1bSJed Brown nsize: {{1 2}separate output} 163c4762a1bSJed Brown args: -view ::ascii_info -mat_type {{aij baij sbaij mpiaij mpibaij mpisbaij}separate output} -mat_block_size {{1 2}separate output} 164c4762a1bSJed Brown 165c4762a1bSJed Brown TEST*/ 166