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