Actual source code: ex195.c
petsc-3.10.5 2019-03-28
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";
10: #include <petscmat.h>
13: int main(int argc,char **args)
14: {
15: Mat A1,A2,A3,A4,nest;
16: Mat mata[4];
17: Mat aij;
18: MPI_Comm comm;
19: PetscInt m,n,istart,iend,ii,i,J,j;
20: PetscScalar v;
21: PetscMPIInt size;
22: PetscErrorCode ierr;
24: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
25: comm = PETSC_COMM_WORLD;
26: MPI_Comm_size(comm,&size);
27: /*
28: Assemble the matrix for the five point stencil, YET AGAIN
29: */
30: MatCreate(comm,&A1);
31: m=2,n=2;
32: MatSetSizes(A1,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
33: MatSetFromOptions(A1);
34: MatSetUp(A1);
35: MatGetOwnershipRange(A1,&istart,&iend);
36: for (ii=istart; ii<iend; ii++) {
37: v = -1.0; i = ii/n; j = ii - i*n;
38: if (i>0) {J = ii - n; MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);}
39: if (i<m-1) {J = ii + n; MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);}
40: if (j>0) {J = ii - 1; MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);}
41: if (j<n-1) {J = ii + 1; MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);}
42: v = 4.0; MatSetValues(A1,1,&ii,1,&ii,&v,INSERT_VALUES);
43: }
44: MatAssemblyBegin(A1,MAT_FINAL_ASSEMBLY);
45: MatAssemblyEnd(A1,MAT_FINAL_ASSEMBLY);
46: MatView(A1,PETSC_VIEWER_STDOUT_WORLD);
47: MatDuplicate(A1,MAT_COPY_VALUES,&A2);
48: MatDuplicate(A1,MAT_COPY_VALUES,&A3);
49: MatDuplicate(A1,MAT_COPY_VALUES,&A4);
50: /*create a nest matrix */
51: MatCreate(comm,&nest);
52: MatSetType(nest,MATNEST);
53: mata[0]=A1,mata[1]=A2,mata[2]=A3,mata[3]=A4;
54: MatNestSetSubMats(nest,2,NULL,2,NULL,mata);
55: MatSetUp(nest);
56: MatConvert(nest,MATAIJ,MAT_INITIAL_MATRIX,&aij);
57: MatView(aij,PETSC_VIEWER_STDOUT_WORLD);
58: MatDestroy(&nest);
59: MatDestroy(&aij);
60: MatDestroy(&A1);
61: MatDestroy(&A2);
62: MatDestroy(&A3);
63: MatDestroy(&A4);
64: PetscFinalize();
65: return ierr;
66: }
69: /*TEST
71: test:
72: nsize: 2
74: TEST*/