xref: /petsc/src/mat/tests/ex95.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
12c4762a1bSJed Brown 
13c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
14*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
15*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
16c4762a1bSJed Brown   prid = size;
17*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-prid",&prid,NULL));
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   m    = n = 10*size;
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&A));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DETERMINE,PETSC_DETERMINE,m,n));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATSEQAIJ));
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   a[0] = rank+1;
26c4762a1bSJed Brown   for (i=0; i<m-rank; i++) {
27c4762a1bSJed Brown     col  = i+rank;
28*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&i,1,&col,a,INSERT_VALUES));
29c4762a1bSJed Brown   }
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   if (rank == prid) {
34*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d] A: \n",rank));
35*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_SELF));
36c4762a1bSJed Brown   }
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   /* Test MatCreateMPIAIJSumSeqAIJ */
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateMPIAIJSumSeqAIJ(PETSC_COMM_WORLD,A,PETSC_DECIDE,PETSC_DECIDE,MAT_INITIAL_MATRIX,&B));
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   /* Test MAT_REUSE_MATRIX */
42c4762a1bSJed Brown   alpha = 0.1;
43c4762a1bSJed Brown   for (i=0; i<3; i++) {
44*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatScale(A,alpha));
45*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateMPIAIJSumSeqAIJ(PETSC_COMM_WORLD,A,PETSC_DECIDE,PETSC_DECIDE,MAT_REUSE_MATRIX,&B));
46c4762a1bSJed Brown   }
47*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B, PETSC_VIEWER_STDOUT_WORLD));
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
50c4762a1bSJed Brown   ierr = PetscFinalize();
51c4762a1bSJed Brown   return ierr;
52c4762a1bSJed Brown }
53c4762a1bSJed Brown 
54c4762a1bSJed Brown /*TEST
55c4762a1bSJed Brown 
56c4762a1bSJed Brown    test:
57c4762a1bSJed Brown       nsize: 3
58c4762a1bSJed Brown       filter: grep -v "MPI processes"
59c4762a1bSJed Brown 
60c4762a1bSJed Brown    test:
61c4762a1bSJed Brown       suffix: 2
62c4762a1bSJed Brown       filter: grep -v "MPI processes"
63c4762a1bSJed Brown 
64c4762a1bSJed Brown TEST*/
65