Actual source code: ex195.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  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>


 16: int main(int argc,char **args)
 17: {
 18:   Mat                 A1,A2,A3,A4,nest;
 19:   Mat                 mata[4];
 20:   Mat                 aij;
 21:   MPI_Comm            comm;
 22:   PetscInt            m,n,istart,iend,ii,i,J,j;
 23:   PetscScalar         v;
 24:   PetscMPIInt         size;
 25:   PetscErrorCode      ierr;

 27:   PetscInitialize(&argc,&args,(char*)0,help);
 28:   comm = PETSC_COMM_WORLD;
 29:   MPI_Comm_size(comm,&size);
 30:   /*
 31:      Assemble the matrix for the five point stencil, YET AGAIN
 32:   */
 33:   MatCreate(comm,&A1);
 34:   m=2,n=2;
 35:   MatSetSizes(A1,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
 36:   MatSetFromOptions(A1);
 37:   MatSetUp(A1);
 38:   MatGetOwnershipRange(A1,&istart,&iend);
 39:   for (ii=istart; ii<iend; ii++) {
 40:     v = -1.0; i = ii/n; j = ii - i*n;
 41:     if (i>0)   {J = ii - n; MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);}
 42:     if (i<m-1) {J = ii + n; MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);}
 43:     if (j>0)   {J = ii - 1; MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);}
 44:     if (j<n-1) {J = ii + 1; MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);}
 45:     v = 4.0; MatSetValues(A1,1,&ii,1,&ii,&v,INSERT_VALUES);
 46:   }
 47:   MatAssemblyBegin(A1,MAT_FINAL_ASSEMBLY);
 48:   MatAssemblyEnd(A1,MAT_FINAL_ASSEMBLY);
 49:   MatView(A1,PETSC_VIEWER_STDOUT_WORLD);
 50:   MatDuplicate(A1,MAT_COPY_VALUES,&A2);
 51:   MatDuplicate(A1,MAT_COPY_VALUES,&A3);
 52:   MatDuplicate(A1,MAT_COPY_VALUES,&A4);
 53:   /*create a nest matrix */
 54:   MatCreate(comm,&nest);
 55:   MatSetType(nest,MATNEST);
 56:   mata[0]=A1,mata[1]=A2,mata[2]=A3,mata[3]=A4;
 57:   MatNestSetSubMats(nest,2,NULL,2,NULL,mata);
 58:   MatSetUp(nest);
 59:   MatConvert(nest,MATAIJ,MAT_INITIAL_MATRIX,&aij);
 60:   MatView(aij,PETSC_VIEWER_STDOUT_WORLD);
 61:   MatDestroy(&A1);
 62:   MatDestroy(&A2);
 63:   MatDestroy(&A3);
 64:   MatDestroy(&A4);
 65:   PetscFinalize();
 66:   return 0;
 67: }