xref: /petsc/src/mat/tests/ex134.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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 bs,MatType mtype)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   const PetscInt    rc[]   = {0,1,2,3};
8c4762a1bSJed Brown   const PetscScalar vals[] = {100, 2, 3, 4, 5, 600, 7, 8,
9c4762a1bSJed Brown                               9,100,11,1200,13,14,15,1600,
10c4762a1bSJed Brown                               17,18,19,20,21,22,23,24,
11c4762a1bSJed Brown                               25,26,27,2800,29,30,31,32,
12c4762a1bSJed Brown                               33,34,35,36,37,38,39,40,
13c4762a1bSJed Brown                               41,42,43,44,45,46,47,48,
14c4762a1bSJed Brown                               49,50,51,52,53,54,55,56,
15c4762a1bSJed Brown                               57,58,49,60,61,62,63,64};
16c4762a1bSJed Brown   Mat               A;
17c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) || defined(PETSC_HAVE_MKL_CPARDISO)
18c4762a1bSJed Brown   Mat               F;
1983863137SStefano Zampini   MatSolverType     stype = MATSOLVERPETSC;
20c4762a1bSJed Brown   PetscRandom       rdm;
21c4762a1bSJed Brown   Vec               b,x,y;
22c4762a1bSJed Brown   PetscInt          i,j;
23560a203cSprj-   PetscReal         norm2,tol=100*PETSC_SMALL;
24c4762a1bSJed Brown   PetscBool         issbaij;
25c4762a1bSJed Brown #endif
26c4762a1bSJed Brown   PetscViewer       viewer;
27c4762a1bSJed Brown   PetscErrorCode    ierr;
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   PetscFunctionBegin;
30c4762a1bSJed Brown   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
31c4762a1bSJed Brown   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,4*bs,4*bs);CHKERRQ(ierr);
32c4762a1bSJed Brown   ierr = MatSetType(A,mtype);CHKERRQ(ierr);
33c4762a1bSJed Brown   ierr = MatMPIBAIJSetPreallocation(A,bs,2,NULL,2,NULL);CHKERRQ(ierr);
34c4762a1bSJed Brown   ierr = MatMPISBAIJSetPreallocation(A,bs,2,NULL,2,NULL);CHKERRQ(ierr);
35c4762a1bSJed Brown   ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
36c4762a1bSJed Brown   /* All processes contribute a global matrix */
37c4762a1bSJed Brown   ierr = MatSetValuesBlocked(A,4,rc,4,rc,vals,ADD_VALUES);CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
39c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
408fffc762SJacob Faibussowitsch   ierr = PetscPrintf(comm,"Matrix %s(%" PetscInt_FMT ")\n",mtype,bs);CHKERRQ(ierr);
41c4762a1bSJed Brown   ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr);
42c4762a1bSJed Brown   ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr);
43c4762a1bSJed Brown   ierr = MatView(A,viewer);CHKERRQ(ierr);
44c4762a1bSJed Brown   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
45c4762a1bSJed Brown   ierr = MatView(A,viewer);CHKERRQ(ierr);
46c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) || defined(PETSC_HAVE_MKL_CPARDISO)
47c4762a1bSJed Brown   ierr = PetscStrcmp(mtype,MATMPISBAIJ,&issbaij);CHKERRQ(ierr);
48c4762a1bSJed Brown   if (!issbaij) {
49c4762a1bSJed Brown     ierr = MatShift(A,10);CHKERRQ(ierr);
50c4762a1bSJed Brown   }
51c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rdm);CHKERRQ(ierr);
52c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr);
53c4762a1bSJed Brown   ierr = MatCreateVecs(A,&x,&y);CHKERRQ(ierr);
54c4762a1bSJed Brown   ierr = VecDuplicate(x,&b);CHKERRQ(ierr);
55c4762a1bSJed Brown   for (j=0; j<2; j++) {
56c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS)
57c4762a1bSJed Brown     if (j==0) stype = MATSOLVERMUMPS;
58c4762a1bSJed Brown #else
59c4762a1bSJed Brown     if (j==0) continue;
60c4762a1bSJed Brown #endif
61c4762a1bSJed Brown #if defined(PETSC_HAVE_MKL_CPARDISO)
62c4762a1bSJed Brown     if (j==1) stype = MATSOLVERMKL_CPARDISO;
63c4762a1bSJed Brown #else
64c4762a1bSJed Brown     if (j==1) continue;
65c4762a1bSJed Brown #endif
66c4762a1bSJed Brown     if (issbaij) {
67c4762a1bSJed Brown       ierr = MatGetFactor(A,stype,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
68c4762a1bSJed Brown       ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
69c4762a1bSJed Brown       ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
70c4762a1bSJed Brown     } else {
71c4762a1bSJed Brown       ierr = MatGetFactor(A,stype,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
72c4762a1bSJed Brown       ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
73c4762a1bSJed Brown       ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
74c4762a1bSJed Brown     }
75c4762a1bSJed Brown     for (i=0; i<10; i++) {
76c4762a1bSJed Brown       ierr = VecSetRandom(b,rdm);CHKERRQ(ierr);
77c4762a1bSJed Brown       ierr = MatSolve(F,b,y);CHKERRQ(ierr);
78c4762a1bSJed Brown       /* Check the error */
79c4762a1bSJed Brown       ierr = MatMult(A,y,x);CHKERRQ(ierr);
80c4762a1bSJed Brown       ierr = VecAXPY(x,-1.0,b);CHKERRQ(ierr);
81c4762a1bSJed Brown       ierr = VecNorm(x,NORM_2,&norm2);CHKERRQ(ierr);
82c4762a1bSJed Brown       if (norm2>tol) {
83c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_WORLD,"Error:MatSolve(), norm2: %g\n",(double)norm2);CHKERRQ(ierr);
84c4762a1bSJed Brown       }
85c4762a1bSJed Brown     }
86c4762a1bSJed Brown     ierr = MatDestroy(&F);CHKERRQ(ierr);
87c4762a1bSJed Brown   }
88c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
89c4762a1bSJed Brown   ierr = VecDestroy(&y);CHKERRQ(ierr);
90c4762a1bSJed Brown   ierr = VecDestroy(&b);CHKERRQ(ierr);
91c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr);
92c4762a1bSJed Brown #endif
93c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
94c4762a1bSJed Brown   PetscFunctionReturn(0);
95c4762a1bSJed Brown }
96c4762a1bSJed Brown 
97c4762a1bSJed Brown int main(int argc,char *argv[])
98c4762a1bSJed Brown {
99c4762a1bSJed Brown   PetscErrorCode ierr;
100c4762a1bSJed Brown   MPI_Comm       comm;
101c4762a1bSJed Brown   PetscMPIInt    size;
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
104c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
105ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
106*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 2,comm,PETSC_ERR_USER,"This example must be run with exactly two processes");
107c4762a1bSJed Brown   ierr = Assemble(comm,2,MATMPIBAIJ);CHKERRQ(ierr);
108c4762a1bSJed Brown   ierr = Assemble(comm,2,MATMPISBAIJ);CHKERRQ(ierr);
109c4762a1bSJed Brown   ierr = Assemble(comm,1,MATMPIBAIJ);CHKERRQ(ierr);
110c4762a1bSJed Brown   ierr = Assemble(comm,1,MATMPISBAIJ);CHKERRQ(ierr);
111c4762a1bSJed Brown   ierr = PetscFinalize();
112c4762a1bSJed Brown   return ierr;
113c4762a1bSJed Brown }
114c4762a1bSJed Brown 
115c4762a1bSJed Brown /*TEST
116c4762a1bSJed Brown 
117c4762a1bSJed Brown    test:
118c4762a1bSJed Brown       nsize: 2
11997929ea7SJunchao Zhang       args: -mat_ignore_lower_triangular
120c4762a1bSJed Brown       filter: sed -e "s~mem [0-9]*~mem~g"
121c4762a1bSJed Brown 
122c4762a1bSJed Brown TEST*/
123