xref: /petsc/src/mat/tests/ex210.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Test MatCreateNest with block sizes.\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   Mat                    A, B, C, mats[2];
8*c4762a1bSJed Brown   ISLocalToGlobalMapping cmap, rmap;
9*c4762a1bSJed Brown   const PetscInt         indices[1] = {0};
10*c4762a1bSJed Brown   PetscMPIInt            size;
11*c4762a1bSJed Brown   PetscErrorCode         ierr;
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr;
14*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr);
15*c4762a1bSJed Brown   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Only coded for one process");
16*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD, &A);CHKERRQ(ierr);
17*c4762a1bSJed Brown   ierr = MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 1, 2);CHKERRQ(ierr);
18*c4762a1bSJed Brown   ierr = MatSetBlockSizes(A, 1, 2);CHKERRQ(ierr);
19*c4762a1bSJed Brown   ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr);
20*c4762a1bSJed Brown   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 2, 1, indices,PETSC_COPY_VALUES, &cmap);CHKERRQ(ierr);
21*c4762a1bSJed Brown   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 1, 1, indices,PETSC_COPY_VALUES, &rmap);CHKERRQ(ierr);
22*c4762a1bSJed Brown   ierr = MatSetLocalToGlobalMapping(A, rmap, cmap);CHKERRQ(ierr);
23*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD, &B);CHKERRQ(ierr);
24*c4762a1bSJed Brown   ierr = MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, 1, 1);CHKERRQ(ierr);
25*c4762a1bSJed Brown   ierr = MatSetBlockSizes(A, 1, 1);CHKERRQ(ierr);
26*c4762a1bSJed Brown   ierr = MatSetType(B,MATAIJ);CHKERRQ(ierr);
27*c4762a1bSJed Brown   ierr = MatSetLocalToGlobalMapping(B, rmap, rmap);CHKERRQ(ierr);
28*c4762a1bSJed Brown   ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
29*c4762a1bSJed Brown   ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
30*c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = MatSetUp(B);CHKERRQ(ierr);
32*c4762a1bSJed Brown   mats[0] = A;
33*c4762a1bSJed Brown   mats[1] = B;
34*c4762a1bSJed Brown   ierr = MatCreateNest(PETSC_COMM_WORLD, 1, NULL, 2, NULL, mats,&C);CHKERRQ(ierr);
35*c4762a1bSJed Brown   ierr = MatSetUp(C);CHKERRQ(ierr);
36*c4762a1bSJed Brown   ierr = MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
37*c4762a1bSJed Brown   ierr = MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = MatView(C, NULL);CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = PetscFinalize();
43*c4762a1bSJed Brown   return ierr;
44*c4762a1bSJed Brown }
45*c4762a1bSJed Brown 
46*c4762a1bSJed Brown /*TEST
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown     test:
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown TEST*/
51