xref: /petsc/src/mat/tests/ex160.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatMPIBAIJ format in sequential run \n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown int main(int argc,char **args)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   Mat            A,B;
8c4762a1bSJed Brown   PetscInt       i,rstart,rend;
9c4762a1bSJed Brown   PetscMPIInt    rank,size;
10c4762a1bSJed Brown   PetscErrorCode ierr;
11c4762a1bSJed Brown   PetscScalar    v;
12c4762a1bSJed Brown 
13c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
14*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
15*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   /* Create a MPIBAIJ matrix */
18*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,32,32));
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATMPIBAIJ));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqBAIJSetPreallocation(A,2,2,NULL));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIBAIJSetPreallocation(A,2,2,NULL,2,NULL));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   v    = 1.0;
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend));
26c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
27*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&i,1,&i,&v,INSERT_VALUES));
28c4762a1bSJed Brown   }
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   /* Convert A to AIJ format */
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B));
34c4762a1bSJed Brown 
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
37c4762a1bSJed Brown   ierr = PetscFinalize();
38c4762a1bSJed Brown   return ierr;
39c4762a1bSJed Brown }
40c4762a1bSJed Brown 
41c4762a1bSJed Brown /*TEST
42c4762a1bSJed Brown 
43c4762a1bSJed Brown    test:
44c4762a1bSJed Brown      output_file: output/ex160.out
45c4762a1bSJed Brown 
46c4762a1bSJed Brown TEST*/
47