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 PetscErrorCode ierr; 11c4762a1bSJed Brown PetscMPIInt rank,size; 12c4762a1bSJed Brown PetscScalar v; 13c4762a1bSJed Brown char mtype[256]; 14c4762a1bSJed Brown 15c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 16ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 17c4762a1bSJed Brown 18c4762a1bSJed Brown /* This example does not work correctly for np > 2 */ 19ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 20*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(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 */ 23c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); 24c4762a1bSJed Brown ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr); 25c4762a1bSJed Brown ierr = MatSetFromOptions(C);CHKERRQ(ierr); 26c4762a1bSJed Brown ierr = MatSetUp(C);CHKERRQ(ierr); 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; 30c4762a1bSJed Brown if (i>0) {J = Ii - n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 31c4762a1bSJed Brown if (i<m-1) {J = Ii + n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 32c4762a1bSJed Brown if (j>0) {J = Ii - 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 33c4762a1bSJed Brown if (j<n-1) {J = Ii + 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 34c4762a1bSJed Brown v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr); 35c4762a1bSJed Brown } 36c4762a1bSJed Brown } 37c4762a1bSJed Brown 38c4762a1bSJed Brown ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 39c4762a1bSJed Brown ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 40c4762a1bSJed Brown ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 42c4762a1bSJed Brown ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 43c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 44c4762a1bSJed Brown 45c4762a1bSJed Brown ierr = PetscStrcpy(mtype,MATSAME);CHKERRQ(ierr); 46589a23caSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-conv_mat_type",mtype,sizeof(mtype),NULL);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = MatConvert(C,mtype,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_IMPL);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Free data structures */ 56c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 57c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 58c4762a1bSJed Brown 59c4762a1bSJed Brown ierr = PetscFinalize(); 60c4762a1bSJed Brown return ierr; 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