Actual source code: ex195.c

petsc-3.13.1 2020-05-02
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";

 10:  #include <petscmat.h>


 13: int main(int argc,char **args)
 14: {
 15:   Mat                 A1,A2,A3,A4,A5,B,C,C1,nest;
 16:   Mat                 mata[4];
 17:   Mat                 aij;
 18:   MPI_Comm            comm;
 19:   PetscInt            m,M,n,istart,iend,ii,i,J,j,K=10;
 20:   PetscScalar         v;
 21:   PetscMPIInt         size;
 22:   PetscErrorCode      ierr;
 23:   PetscBool           equal;

 25:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 26:   comm = PETSC_COMM_WORLD;
 27:   MPI_Comm_size(comm,&size);

 29:   /*
 30:      Assemble the matrix for the five point stencil, YET AGAIN
 31:   */
 32:   MatCreate(comm,&A1);
 33:   m=2,n=2;
 34:   MatSetSizes(A1,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
 35:   MatSetFromOptions(A1);
 36:   MatSetUp(A1);
 37:   MatGetOwnershipRange(A1,&istart,&iend);
 38:   for (ii=istart; ii<iend; ii++) {
 39:     v = -1.0; i = ii/n; j = ii - i*n;
 40:     if (i>0)   {J = ii - n; MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);}
 41:     if (i<m-1) {J = ii + n; MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);}
 42:     if (j>0)   {J = ii - 1; MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);}
 43:     if (j<n-1) {J = ii + 1; MatSetValues(A1,1,&ii,1,&J,&v,INSERT_VALUES);}
 44:     v = 4.0; MatSetValues(A1,1,&ii,1,&ii,&v,INSERT_VALUES);
 45:   }
 46:   MatAssemblyBegin(A1,MAT_FINAL_ASSEMBLY);
 47:   MatAssemblyEnd(A1,MAT_FINAL_ASSEMBLY);
 48:   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);

 54:   /*create a nest matrix */
 55:   MatCreate(comm,&nest);
 56:   MatSetType(nest,MATNEST);
 57:   mata[0]=A1,mata[1]=A2,mata[2]=A3,mata[3]=A4;
 58:   MatNestSetSubMats(nest,2,NULL,2,NULL,mata);
 59:   MatSetUp(nest);
 60:   MatConvert(nest,MATAIJ,MAT_INITIAL_MATRIX,&aij);
 61:   MatView(aij,PETSC_VIEWER_STDOUT_WORLD);

 63:   /* create a dense matrix */
 64:   MatGetSize(nest,&M,NULL);
 65:   MatGetLocalSize(nest,&m,NULL);
 66:   MatCreateDense(comm,m,PETSC_DECIDE,M,K,NULL,&B);
 67:   MatSetRandom(B,PETSC_NULL);
 68:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 69:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 71:   /* C = nest*B_dense */
 72:   MatMatMult(nest,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C);
 73:   MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C);
 74:   MatMatMultEqual(nest,B,C,10,&equal);
 75:   if (!equal) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in C != nest*B_dense");
 76:   MatDestroy(&nest);

 78:   if (size > 1) { /* Do not know why this test fails for size = 1 */
 79:     MatCreateTranspose(A1,&A5); /* A1 is symmetric */
 80:     mata[0] = A5;
 81:     MatCreate(comm,&nest);
 82:     MatSetType(nest,MATNEST);
 83:     MatNestSetSubMats(nest,2,NULL,2,NULL,mata);
 84:     MatSetUp(nest);
 85:     MatMatMult(nest,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C1);
 86:     MatMatMult(nest,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C1);

 88:     MatEqual(C1,C,&equal);
 89:     if (!equal) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in C1 != C");
 90:     MatDestroy(&C1);
 91:     MatDestroy(&A5);
 92:   }

 94:   MatDestroy(&C);
 95:   MatDestroy(&B);
 96:   MatDestroy(&nest);
 97:   MatDestroy(&aij);
 98:   MatDestroy(&A1);
 99:   MatDestroy(&A2);
100:   MatDestroy(&A3);
101:   MatDestroy(&A4);
102:   PetscFinalize();
103:   return ierr;
104: }


107: /*TEST

109:    test:
110:       nsize: 2

112: TEST*/