1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests converting a matrix to another format with MatConvert().\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,j,m = 5,n = 4,Ii,J; 10c4762a1bSJed Brown PetscMPIInt rank,size; 11c4762a1bSJed Brown PetscScalar v; 12c4762a1bSJed Brown char mtype[256]; 13c4762a1bSJed Brown 14*327415f7SBarry Smith PetscFunctionBeginUser; 159566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 17c4762a1bSJed Brown 18c4762a1bSJed Brown /* This example does not work correctly for np > 2 */ 199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 20be096a46SBarry Smith PetscCheck(size <= 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Use np <= 2"); 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* Create the matrix for the five point stencil, YET AGAIN */ 239566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 249566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 259566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 269566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 27c4762a1bSJed Brown for (i=0; i<m; i++) { 28c4762a1bSJed Brown for (j=2*rank; j<2*rank+2; j++) { 29c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 309566063dSJacob Faibussowitsch if (i>0) {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 319566063dSJacob Faibussowitsch if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 329566063dSJacob Faibussowitsch if (j>0) {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 339566063dSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 349566063dSJacob Faibussowitsch v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 35c4762a1bSJed Brown } 36c4762a1bSJed Brown } 37c4762a1bSJed Brown 389566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 399566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 409566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO)); 419566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 429566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 439566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 44c4762a1bSJed Brown 459566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(mtype,MATSAME)); 469566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-conv_mat_type",mtype,sizeof(mtype),NULL)); 479566063dSJacob Faibussowitsch PetscCall(MatConvert(C,mtype,MAT_INITIAL_MATRIX,&A)); 489566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO)); 499566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 509566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 519566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_IMPL)); 529566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 539566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Free data structures */ 569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 58c4762a1bSJed Brown 599566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 60b122ec5aSJacob Faibussowitsch return 0; 61c4762a1bSJed Brown } 62c4762a1bSJed Brown 63c4762a1bSJed Brown /*TEST 64c4762a1bSJed Brown 65c4762a1bSJed Brown test: 66c4762a1bSJed Brown args: -conv_mat_type seqaij 67c4762a1bSJed Brown 68c4762a1bSJed Brown TEST*/ 69