xref: /petsc/src/mat/tests/ex182.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests using MatShift() to create a constant diagonal matrix\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown int main(int argc,char **argv)
7*c4762a1bSJed Brown {
8*c4762a1bSJed Brown   Mat            A,F;
9*c4762a1bSJed Brown   MatFactorInfo  info;
10*c4762a1bSJed Brown   PetscErrorCode ierr;
11*c4762a1bSJed Brown   PetscInt       m = 10;
12*c4762a1bSJed Brown   IS             perm;
13*c4762a1bSJed Brown   PetscMPIInt    size;
14*c4762a1bSJed Brown   PetscBool      issbaij;
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*) 0,help);if (ierr) return ierr;
17*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
20*c4762a1bSJed Brown   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);CHKERRQ(ierr);
21*c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
22*c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
23*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
24*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
25*c4762a1bSJed Brown 
26*c4762a1bSJed Brown   ierr = MatShift(A,1.0);CHKERRQ(ierr);
27*c4762a1bSJed Brown 
28*c4762a1bSJed Brown   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
29*c4762a1bSJed Brown   if (size == 1 && !issbaij) {
30*c4762a1bSJed Brown     ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
31*c4762a1bSJed Brown     ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
32*c4762a1bSJed Brown     ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&perm);CHKERRQ(ierr);
33*c4762a1bSJed Brown     ierr = MatLUFactorSymbolic(F,A,perm,perm,&info);CHKERRQ(ierr);
34*c4762a1bSJed Brown     ierr = MatLUFactorNumeric(F,A,&info);CHKERRQ(ierr);
35*c4762a1bSJed Brown     ierr = MatDestroy(&F);CHKERRQ(ierr);
36*c4762a1bSJed Brown     ierr = ISDestroy(&perm);CHKERRQ(ierr);
37*c4762a1bSJed Brown   }
38*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = PetscFinalize();
40*c4762a1bSJed Brown   return ierr;
41*c4762a1bSJed Brown }
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown /*TEST
44*c4762a1bSJed Brown 
45*c4762a1bSJed Brown    test:
46*c4762a1bSJed Brown       requires: define(PETSC_USE_INFO)
47*c4762a1bSJed Brown       args: -info
48*c4762a1bSJed Brown       filter: grep malloc | sort -b
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown    test:
51*c4762a1bSJed Brown       suffix: 2
52*c4762a1bSJed Brown       nsize: 2
53*c4762a1bSJed Brown       requires: define(PETSC_USE_INFO)
54*c4762a1bSJed Brown       args: -info ex182info
55*c4762a1bSJed Brown       filter: grep -h malloc "ex182info.0" | sort -b
56*c4762a1bSJed Brown 
57*c4762a1bSJed Brown    test:
58*c4762a1bSJed Brown       suffix: 3
59*c4762a1bSJed Brown       requires: define(PETSC_USE_INFO)
60*c4762a1bSJed Brown       args: -info -mat_type baij
61*c4762a1bSJed Brown       filter: grep malloc | sort -b
62*c4762a1bSJed Brown 
63*c4762a1bSJed Brown    test:
64*c4762a1bSJed Brown       suffix: 4
65*c4762a1bSJed Brown       nsize: 2
66*c4762a1bSJed Brown       requires: define(PETSC_USE_INFO)
67*c4762a1bSJed Brown       args: -info ex182info -mat_type baij
68*c4762a1bSJed Brown       filter: grep -h malloc "ex182info.1" | sort -b
69*c4762a1bSJed Brown 
70*c4762a1bSJed Brown    test:
71*c4762a1bSJed Brown       suffix: 5
72*c4762a1bSJed Brown       requires: define(PETSC_USE_INFO)
73*c4762a1bSJed Brown       args: -info -mat_type sbaij
74*c4762a1bSJed Brown       filter: grep malloc | sort  -b
75*c4762a1bSJed Brown 
76*c4762a1bSJed Brown    test:
77*c4762a1bSJed Brown       suffix: 6
78*c4762a1bSJed Brown       nsize: 2
79*c4762a1bSJed Brown       requires: define(PETSC_USE_INFO)
80*c4762a1bSJed Brown       args: -info ex182info -mat_type sbaij
81*c4762a1bSJed Brown       filter: grep -h malloc "ex182info.0" | sort -b
82*c4762a1bSJed Brown 
83*c4762a1bSJed Brown TEST*/
84