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