xref: /petsc/src/mat/tests/ex233.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Tests MatMPI{AIJ,BAIJ,SBAIJ}SetPreallocationCSR\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petscmat.h>
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown int main(int argc,char **argv)
6*c4762a1bSJed Brown {
7*c4762a1bSJed Brown   PetscErrorCode ierr;
8*c4762a1bSJed Brown   PetscInt       ia[3]={0,2,4};
9*c4762a1bSJed Brown   PetscInt       ja[4]={0,1,0,1};
10*c4762a1bSJed Brown   PetscScalar    c[16]={0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15};
11*c4762a1bSJed Brown   PetscInt       ia2[5]={0,4,8,12,16};
12*c4762a1bSJed Brown   PetscInt       ja2[16]={0,1,2,3,0,1,2,3,0,1,2,3,0,1,2,3};
13*c4762a1bSJed Brown   PetscScalar    c2[16]={0,1,4,5,2,3,6,7,8,9,12,13,10,11,14,15};
14*c4762a1bSJed Brown   PetscMPIInt    size,rank;
15*c4762a1bSJed Brown   Mat            ssbaij;
16*c4762a1bSJed Brown   PetscBool      rect = PETSC_FALSE;
17*c4762a1bSJed Brown 
18*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
19*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
20*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
21*c4762a1bSJed Brown   if (size < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is an example with more then one processors");
22*c4762a1bSJed Brown   if (rank) {
23*c4762a1bSJed Brown     PetscInt i;
24*c4762a1bSJed Brown     for (i = 0; i < 3; i++) ia[i] = 0;
25*c4762a1bSJed Brown     for (i = 0; i < 5; i++) ia2[i] = 0;
26*c4762a1bSJed Brown   }
27*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-rect",&rect,NULL);CHKERRQ(ierr);
28*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&ssbaij);CHKERRQ(ierr);
29*c4762a1bSJed Brown   ierr = MatSetBlockSize(ssbaij,2);CHKERRQ(ierr);
30*c4762a1bSJed Brown   if (rect) {
31*c4762a1bSJed Brown     ierr = MatSetType(ssbaij,MATMPIBAIJ);CHKERRQ(ierr);
32*c4762a1bSJed Brown     ierr = MatSetSizes(ssbaij,4,6,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
33*c4762a1bSJed Brown   } else {
34*c4762a1bSJed Brown     ierr = MatSetType(ssbaij,MATMPISBAIJ);CHKERRQ(ierr);
35*c4762a1bSJed Brown     ierr = MatSetSizes(ssbaij,4,4,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
36*c4762a1bSJed Brown   }
37*c4762a1bSJed Brown   ierr = MatSetFromOptions(ssbaij);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = MatMPIAIJSetPreallocationCSR(ssbaij,ia2,ja2,c2);CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = MatMPIBAIJSetPreallocationCSR(ssbaij,2,ia,ja,c);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = MatMPISBAIJSetPreallocationCSR(ssbaij,2,ia,ja,c);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = MatViewFromOptions(ssbaij,NULL,"-view");CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = MatDestroy(&ssbaij);CHKERRQ(ierr);
43*c4762a1bSJed Brown   ierr = PetscFinalize();
44*c4762a1bSJed Brown   return ierr;
45*c4762a1bSJed Brown }
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown /*TEST
48*c4762a1bSJed Brown 
49*c4762a1bSJed Brown   test:
50*c4762a1bSJed Brown     filter: grep -v type | sed -e "s/\.//g"
51*c4762a1bSJed Brown     suffix: aijbaij_csr
52*c4762a1bSJed Brown     nsize: 2
53*c4762a1bSJed Brown     args: -mat_type {{aij baij}} -view -rect {{0 1}}
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown   test:
56*c4762a1bSJed Brown     filter: sed -e "s/\.//g"
57*c4762a1bSJed Brown     suffix: sbaij_csr
58*c4762a1bSJed Brown     nsize: 2
59*c4762a1bSJed Brown     args: -mat_type sbaij -view -rect {{0 1}}
60*c4762a1bSJed Brown 
61*c4762a1bSJed Brown TEST*/
62