xref: /petsc/src/mat/tests/ex95.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "Testing MatCreateMPIAIJSumSeqAIJ().\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc,char **argv)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   Mat            A,B;
8c4762a1bSJed Brown   MatScalar      a[1],alpha;
9c4762a1bSJed Brown   PetscMPIInt    size,rank;
10c4762a1bSJed Brown   PetscInt       m,n,i,col, prid;
11c4762a1bSJed Brown 
12*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
135f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
145f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
15c4762a1bSJed Brown   prid = size;
165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-prid",&prid,NULL));
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   m    = n = 10*size;
195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&A));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DETERMINE,PETSC_DETERMINE,m,n));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATSEQAIJ));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   a[0] = rank+1;
25c4762a1bSJed Brown   for (i=0; i<m-rank; i++) {
26c4762a1bSJed Brown     col  = i+rank;
275f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&i,1,&col,a,INSERT_VALUES));
28c4762a1bSJed Brown   }
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   if (rank == prid) {
335f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d] A: \n",rank));
345f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_SELF));
35c4762a1bSJed Brown   }
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /* Test MatCreateMPIAIJSumSeqAIJ */
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateMPIAIJSumSeqAIJ(PETSC_COMM_WORLD,A,PETSC_DECIDE,PETSC_DECIDE,MAT_INITIAL_MATRIX,&B));
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   /* Test MAT_REUSE_MATRIX */
41c4762a1bSJed Brown   alpha = 0.1;
42c4762a1bSJed Brown   for (i=0; i<3; i++) {
435f80ce2aSJacob Faibussowitsch     CHKERRQ(MatScale(A,alpha));
445f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateMPIAIJSumSeqAIJ(PETSC_COMM_WORLD,A,PETSC_DECIDE,PETSC_DECIDE,MAT_REUSE_MATRIX,&B));
45c4762a1bSJed Brown   }
465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B, PETSC_VIEWER_STDOUT_WORLD));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
49*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
50*b122ec5aSJacob Faibussowitsch   return 0;
51c4762a1bSJed Brown }
52c4762a1bSJed Brown 
53c4762a1bSJed Brown /*TEST
54c4762a1bSJed Brown 
55c4762a1bSJed Brown    test:
56c4762a1bSJed Brown       nsize: 3
57c4762a1bSJed Brown       filter: grep -v "MPI processes"
58c4762a1bSJed Brown 
59c4762a1bSJed Brown    test:
60c4762a1bSJed Brown       suffix: 2
61c4762a1bSJed Brown       filter: grep -v "MPI processes"
62c4762a1bSJed Brown 
63c4762a1bSJed Brown TEST*/
64