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; ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL);CHKERRQ(ierr); ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr); ierr = MatSetType(C,MATAIJ);CHKERRQ(ierr); ierr = MatSetFromOptions(C);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject)C,"initial");CHKERRQ(ierr); ierr = MatGetType(C,&type);CHKERRQ(ierr); if (size == 1) { ierr = PetscObjectTypeCompare((PetscObject)C,MATSEQAIJ,&isAIJ);CHKERRQ(ierr); } else { ierr = PetscObjectTypeCompare((PetscObject)C,MATMPIAIJ,&isAIJ);CHKERRQ(ierr); } ierr = MatSeqAIJSetPreallocation(C,3,NULL);CHKERRQ(ierr); ierr = MatMPIAIJSetPreallocation(C,3,NULL,3,NULL);CHKERRQ(ierr); ierr = MatSeqBAIJSetPreallocation(C,bs,3,NULL);CHKERRQ(ierr); ierr = MatMPIBAIJSetPreallocation(C,bs,3,NULL,3,NULL);CHKERRQ(ierr); ierr = MatSeqSBAIJSetPreallocation(C,bs,3,NULL);CHKERRQ(ierr); ierr = MatMPISBAIJSetPreallocation(C,bs,3,NULL,3,NULL);CHKERRQ(ierr); v[0] = -1.; v[1] = 2.; v[2] = -1.; for (i=1; i %s: Matrices are NOT multequal",Ctype,Cstype); ierr = MatConvert(Cs,Ctype,MAT_REUSE_MATRIX,&Cse);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject)Cse,"symm_conv_reuse");CHKERRQ(ierr); ierr = MatViewFromOptions(Cse,NULL,"-view");CHKERRQ(ierr); ierr = MatMultEqual(Cs,Cse,5,&flg);CHKERRQ(ierr); if (!flg) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatConvert MAT_REUSE_MATRIX %s -> %s: Matrices are NOT multequal",Ctype,Cstype); ierr = MatCopy(Cs,Cse,SAME_NONZERO_PATTERN);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject)Cse,"symm_conv_copy_samennz");CHKERRQ(ierr); ierr = MatViewFromOptions(Cse,NULL,"-view");CHKERRQ(ierr); ierr = MatMultEqual(Cs,Cse,5,&flg);CHKERRQ(ierr); if (!flg) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatCopy(...SAME_NONZERO_PATTERN) %s -> %s: Matrices are NOT multequal",Ctype,Cstype); ierr = MatCopy(Cs,Cse,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject)Cse,"symm_conv_copy_subnnz");CHKERRQ(ierr); ierr = MatViewFromOptions(Cse,NULL,"-view");CHKERRQ(ierr); ierr = MatMultEqual(Cs,Cse,5,&flg);CHKERRQ(ierr); if (!flg) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatCopy(...SUBSET_NONZERO_PATTERN) %s -> %s: Matrices are NOT multequal",Ctype,Cstype); ierr = MatCopy(Cs,Cse,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject)Cse,"symm_conv_copy_diffnnz");CHKERRQ(ierr); ierr = MatViewFromOptions(Cse,NULL,"-view");CHKERRQ(ierr); ierr = MatMultEqual(Cs,Cse,5,&flg);CHKERRQ(ierr); if (!flg) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatCopy(...DIFFERENT_NONZERO_PATTERN) %s -> %s: Matrices are NOT multequal",Ctype,Cstype); ierr = MatDestroy(&Cse);CHKERRQ(ierr); ierr = MatDestroy(&Cs);CHKERRQ(ierr); } /* test MatStore/RetrieveValues() */ if (isAIJ) { ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);CHKERRQ(ierr); ierr = MatStoreValues(A);CHKERRQ(ierr); ierr = MatZeroEntries(A);CHKERRQ(ierr); ierr = MatRetrieveValues(A);CHKERRQ(ierr); } ierr = MatDestroy(&C);CHKERRQ(ierr); ierr = MatDestroy(&A);CHKERRQ(ierr); 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*/