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