static char help[] = "Tests MatCopy() and MatStore/RetrieveValues().\n\n"; #include int main(int argc,char **args) { Mat C,A; PetscInt i, n = 10,midx[3],bs=1; PetscErrorCode ierr; PetscScalar v[3]; PetscBool flg,isAIJ; MatType type; PetscMPIInt size; ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL)); CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n)); CHKERRQ(MatSetType(C,MATAIJ)); CHKERRQ(MatSetFromOptions(C)); CHKERRQ(PetscObjectSetName((PetscObject)C,"initial")); CHKERRQ(MatGetType(C,&type)); if (size == 1) { CHKERRQ(PetscObjectTypeCompare((PetscObject)C,MATSEQAIJ,&isAIJ)); } else { CHKERRQ(PetscObjectTypeCompare((PetscObject)C,MATMPIAIJ,&isAIJ)); } CHKERRQ(MatSeqAIJSetPreallocation(C,3,NULL)); CHKERRQ(MatMPIAIJSetPreallocation(C,3,NULL,3,NULL)); CHKERRQ(MatSeqBAIJSetPreallocation(C,bs,3,NULL)); CHKERRQ(MatMPIBAIJSetPreallocation(C,bs,3,NULL,3,NULL)); CHKERRQ(MatSeqSBAIJSetPreallocation(C,bs,3,NULL)); CHKERRQ(MatMPISBAIJSetPreallocation(C,bs,3,NULL,3,NULL)); v[0] = -1.; v[1] = 2.; v[2] = -1.; for (i=1; i %s: Matrices are NOT multequal",Ctype,Cstype); CHKERRQ(MatConvert(Cs,Ctype,MAT_REUSE_MATRIX,&Cse)); CHKERRQ(PetscObjectSetName((PetscObject)Cse,"symm_conv_reuse")); CHKERRQ(MatViewFromOptions(Cse,NULL,"-view")); CHKERRQ(MatMultEqual(Cs,Cse,5,&flg)); PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatConvert MAT_REUSE_MATRIX %s -> %s: Matrices are NOT multequal",Ctype,Cstype); CHKERRQ(MatCopy(Cs,Cse,SAME_NONZERO_PATTERN)); CHKERRQ(PetscObjectSetName((PetscObject)Cse,"symm_conv_copy_samennz")); CHKERRQ(MatViewFromOptions(Cse,NULL,"-view")); CHKERRQ(MatMultEqual(Cs,Cse,5,&flg)); PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatCopy(...SAME_NONZERO_PATTERN) %s -> %s: Matrices are NOT multequal",Ctype,Cstype); CHKERRQ(MatCopy(Cs,Cse,SUBSET_NONZERO_PATTERN)); CHKERRQ(PetscObjectSetName((PetscObject)Cse,"symm_conv_copy_subnnz")); CHKERRQ(MatViewFromOptions(Cse,NULL,"-view")); CHKERRQ(MatMultEqual(Cs,Cse,5,&flg)); PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatCopy(...SUBSET_NONZERO_PATTERN) %s -> %s: Matrices are NOT multequal",Ctype,Cstype); CHKERRQ(MatCopy(Cs,Cse,DIFFERENT_NONZERO_PATTERN)); CHKERRQ(PetscObjectSetName((PetscObject)Cse,"symm_conv_copy_diffnnz")); CHKERRQ(MatViewFromOptions(Cse,NULL,"-view")); CHKERRQ(MatMultEqual(Cs,Cse,5,&flg)); PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatCopy(...DIFFERENT_NONZERO_PATTERN) %s -> %s: Matrices are NOT multequal",Ctype,Cstype); CHKERRQ(MatDestroy(&Cse)); CHKERRQ(MatDestroy(&Cs)); } /* test MatStore/RetrieveValues() */ if (isAIJ) { CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE)); CHKERRQ(MatStoreValues(A)); CHKERRQ(MatZeroEntries(A)); CHKERRQ(MatRetrieveValues(A)); } CHKERRQ(MatDestroy(&C)); CHKERRQ(MatDestroy(&A)); ierr = PetscFinalize(); return ierr; } /*TEST testset: nsize: {{1 2}separate output} args: -view ::ascii_info -mat_type {{aij baij sbaij mpiaij mpibaij mpisbaij}separate output} -mat_block_size {{1 2}separate output} TEST*/