xref: /petsc/src/mat/tests/ex135.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static const char help[] = "Test parallel assembly of SBAIJ matrices\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown PetscErrorCode Assemble(MPI_Comm comm,PetscInt n,MatType mtype)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   Mat            A;
8c4762a1bSJed Brown   PetscInt       first,last,i;
9c4762a1bSJed Brown   PetscMPIInt    rank,size;
10c4762a1bSJed Brown 
11c4762a1bSJed Brown   PetscFunctionBegin;
125f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
135f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A, PETSC_DECIDE,PETSC_DECIDE,n,n));
145f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATMPISBAIJ));
155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
165f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
175f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm,&rank));
18c4762a1bSJed Brown   if (rank < size-1) {
195f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMPISBAIJSetPreallocation(A,1,1,NULL,1,NULL));
20c4762a1bSJed Brown   } else {
215f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMPISBAIJSetPreallocation(A,1,2,NULL,0,NULL));
22c4762a1bSJed Brown   }
235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(A,&first,&last));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE));
25c4762a1bSJed Brown   last--;
26c4762a1bSJed Brown   for (i=first; i<=last; i++) {
275f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValue(A,i,i,2.,INSERT_VALUES));
285f80ce2aSJacob Faibussowitsch     if (i != n-1) CHKERRQ(MatSetValue(A,i,n-1,-1.,INSERT_VALUES));
29c4762a1bSJed Brown   }
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
33c4762a1bSJed Brown   PetscFunctionReturn(0);
34c4762a1bSJed Brown }
35c4762a1bSJed Brown 
36c4762a1bSJed Brown int main(int argc,char *argv[])
37c4762a1bSJed Brown {
38c4762a1bSJed Brown   MPI_Comm       comm;
39c4762a1bSJed Brown   PetscInt       n = 6;
40c4762a1bSJed Brown 
41*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
42c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(Assemble(comm,n,MATMPISBAIJ));
45*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
46*b122ec5aSJacob Faibussowitsch   return 0;
47c4762a1bSJed Brown }
48c4762a1bSJed Brown 
49c4762a1bSJed Brown /*TEST
50c4762a1bSJed Brown 
51c4762a1bSJed Brown    test:
52c4762a1bSJed Brown       nsize: 4
5397929ea7SJunchao Zhang       args: -n 1000 -mat_view ascii::ascii_info_detail
54dfd57a17SPierre Jolivet       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
55c4762a1bSJed Brown 
56c4762a1bSJed Brown TEST*/
57