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