Actual source code: ex195.c
petsc-3.8.4 2018-03-24
1: /*
2: * ex195.c
3: *
4: * Created on: Aug 24, 2015
5: * Author: Fande Kong <fdkong.jd@gmail.com>
6: */
8: static char help[] = " Demonstrate the use of MatConvert_Nest_AIJ\n";
11: #include <petscmat.h>
14: int main(int argc,char **args)
15: {
16: Mat A1,A2,A3,A4,nest;
17: Mat mata[4];
18: Mat aij;
19: MPI_Comm comm;
20: PetscInt m,n,istart,iend,ii,i,J,j;
21: PetscScalar v;
22: PetscMPIInt size;
23: PetscErrorCode ierr;
25: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
26: comm = PETSC_COMM_WORLD;
27: MPI_Comm_size(comm,&size);
28: /*
29: Assemble the matrix for the five point stencil, YET AGAIN
30: */
31: MatCreate(comm,&A1);
32: m=2,n=2;
33: MatSetSizes(A1,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
34: MatSetFromOptions(A1);
35: MatSetUp(A1);
36: MatGetOwnershipRange(A1,&istart,&iend);
37: for (ii=istart; ii<iend; ii++) {
38: v = -1.0; i = ii/n; j = ii - i*n;
39: if (i>0) {J = ii - n; MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);}
40: if (i<m-1) {J = ii + n; MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);}
41: if (j>0) {J = ii - 1; MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);}
42: if (j<n-1) {J = ii + 1; MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);}
43: v = 4.0; MatSetValues(A1,1,&ii,1,&ii,&v,INSERT_VALUES);
44: }
45: MatAssemblyBegin(A1,MAT_FINAL_ASSEMBLY);
46: MatAssemblyEnd(A1,MAT_FINAL_ASSEMBLY);
47: MatView(A1,PETSC_VIEWER_STDOUT_WORLD);
48: MatDuplicate(A1,MAT_COPY_VALUES,&A2);
49: MatDuplicate(A1,MAT_COPY_VALUES,&A3);
50: MatDuplicate(A1,MAT_COPY_VALUES,&A4);
51: /*create a nest matrix */
52: MatCreate(comm,&nest);
53: MatSetType(nest,MATNEST);
54: mata[0]=A1,mata[1]=A2,mata[2]=A3,mata[3]=A4;
55: MatNestSetSubMats(nest,2,NULL,2,NULL,mata);
56: MatSetUp(nest);
57: MatConvert(nest,MATAIJ,MAT_INITIAL_MATRIX,&aij);
58: MatView(aij,PETSC_VIEWER_STDOUT_WORLD);
59: MatDestroy(&A1);
60: MatDestroy(&A2);
61: MatDestroy(&A3);
62: MatDestroy(&A4);
63: PetscFinalize();
64: return ierr;
65: }