1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Test MatCreate() with MAT_STRUCTURE_ONLY .\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **argv) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat mat; 9c4762a1bSJed Brown PetscInt m = 7,n,i,j,rstart,rend; 10c4762a1bSJed Brown PetscErrorCode ierr; 11c4762a1bSJed Brown PetscMPIInt size; 12c4762a1bSJed Brown PetscScalar v; 13c4762a1bSJed Brown PetscBool struct_only=PETSC_TRUE; 14c4762a1bSJed Brown 15c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 16*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 172c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 18c4762a1bSJed Brown 19*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON)); 20*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 21*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-struct_only",&struct_only,NULL)); 22c4762a1bSJed Brown n = m; 23c4762a1bSJed Brown 24c4762a1bSJed Brown /* ------- Assemble matrix, test MatValid() --------- */ 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&mat)); 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,m,n)); 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(mat)); 28c4762a1bSJed Brown if (struct_only) { 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(mat,MAT_STRUCTURE_ONLY,PETSC_TRUE)); 30c4762a1bSJed Brown } 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(mat)); 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(mat,&rstart,&rend)); 33c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 34c4762a1bSJed Brown for (j=0; j<n; j++) { 35c4762a1bSJed Brown v = 10.0*i+j; 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES)); 37c4762a1bSJed Brown } 38c4762a1bSJed Brown } 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY)); 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY)); 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(mat,PETSC_VIEWER_STDOUT_WORLD)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Free data structures */ 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mat)); 45c4762a1bSJed Brown ierr = PetscFinalize(); 46c4762a1bSJed Brown return ierr; 47c4762a1bSJed Brown } 48c4762a1bSJed Brown 49c4762a1bSJed Brown /*TEST 50c4762a1bSJed Brown 51c4762a1bSJed Brown test: 52c4762a1bSJed Brown output_file: output/ex107.out 53c4762a1bSJed Brown 54c4762a1bSJed Brown test: 55c4762a1bSJed Brown suffix: 2 56c4762a1bSJed Brown args: -mat_type baij -mat_block_size 2 -m 10 57c4762a1bSJed Brown 58c4762a1bSJed Brown TEST*/ 59