Actual source code: ex49.c

  1: /*$Id: ex49.c,v 1.24 2001/08/07 21:30:08 bsmith Exp $*/

  3: static char help[] = "Tests MatTranspose(), MatNorm(), MatValid(), and MatAXPY().nn";

 5:  #include petscmat.h

  7: int main(int argc,char **argv)
  8: {
  9:   Mat          mat,tmat = 0;
 10:   int          m = 4,n,i,j,ierr,size,rank;
 11:   int          rstart,rend,rect = 0,nd,bs,*diag,*bdlen;
 12:   PetscTruth   flg,isbdiag;
 13:   PetscScalar  v,**diagv;
 14:   PetscReal    normf,normi,norm1;
 15:   MatInfo      info;
 16: 
 17:   PetscInitialize(&argc,&argv,(char*)0,help);
 18:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 19:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 20:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 21:   n = m;
 22:   PetscOptionsHasName(PETSC_NULL,"-rect1",&flg);
 23:   if (flg) {n += 2; rect = 1;}
 24:   PetscOptionsHasName(PETSC_NULL,"-rect2",&flg);
 25:   if (flg) {n -= 2; rect = 1;}

 27:   /* Create and assemble matrix */
 28:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,n,&mat);
 29:   MatSetFromOptions(mat);
 30:   MatGetOwnershipRange(mat,&rstart,&rend);
 31:   for (i=rstart; i<rend; i++) {
 32:     for (j=0; j<n; j++) {
 33:       v=10*i+j;
 34:       MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES);
 35:     }
 36:   }
 37:   MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
 38:   MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);

 40:   /* Test whether matrix has been corrupted (just to demonstrate this
 41:      routine) not needed in most application codes. */
 42:   MatValid(mat,(PetscTruth*)&flg);
 43:   if (!flg) SETERRQ(1,"Corrupted matrix.");

 45:   /* Print info about original matrix */
 46:   MatGetInfo(mat,MAT_GLOBAL_SUM,&info);
 47:   PetscPrintf(PETSC_COMM_WORLD,"original matrix nonzeros = %d, allocated nonzeros = %dn",
 48:                     (int)info.nz_used,(int)info.nz_allocated);
 49:   MatNorm(mat,NORM_FROBENIUS,&normf);
 50:   MatNorm(mat,NORM_1,&norm1);
 51:   MatNorm(mat,NORM_INFINITY,&normi);
 52:   PetscPrintf(PETSC_COMM_WORLD,"original: Frobenious norm = %g, one norm = %g, infinity norm = %gn",
 53:                      normf,norm1,normi);
 54:   MatView(mat,PETSC_VIEWER_STDOUT_WORLD);

 56:   PetscTypeCompare((PetscObject)mat,MATSEQBDIAG,&isbdiag);
 57:   if (!isbdiag) {
 58:     PetscTypeCompare((PetscObject)mat,MATMPIBDIAG,&isbdiag);
 59:   }
 60:   if (isbdiag) {
 61:     MatBDiagGetData(mat,&nd,&bs,&diag,&bdlen,&diagv);
 62:     PetscPrintf(PETSC_COMM_WORLD,"original matrix: block diag format: %d diagonals, block size = %dn",nd,bs);
 63:     for (i=0; i<nd; i++) {
 64:       PetscPrintf(PETSC_COMM_WORLD," diag=%d, bdlength=%dn",diag[i],bdlen[i]);
 65:     }
 66:   }

 68:   /* Form matrix transpose */
 69:   PetscOptionsHasName(PETSC_NULL,"-in_place",&flg);
 70:   if (!rect && flg) {
 71:     MatTranspose(mat,0);   /* in-place transpose */
 72:     tmat = mat; mat = 0;
 73:   } else {      /* out-of-place transpose */
 74:     MatTranspose(mat,&tmat);
 75:   }

 77:   /* Print info about transpose matrix */
 78:   MatGetInfo(tmat,MAT_GLOBAL_SUM,&info);
 79:   PetscPrintf(PETSC_COMM_WORLD,"transpose matrix nonzeros = %d, allocated nonzeros = %dn",
 80:                      (int)info.nz_used,(int)info.nz_allocated);
 81:   MatNorm(tmat,NORM_FROBENIUS,&normf);
 82:   MatNorm(tmat,NORM_1,&norm1);
 83:   MatNorm(tmat,NORM_INFINITY,&normi);
 84:   PetscPrintf(PETSC_COMM_WORLD,"transpose: Frobenious norm = %g, one norm = %g, infinity norm = %gn",
 85:                      normf,norm1,normi);
 86:   MatView(tmat,PETSC_VIEWER_STDOUT_WORLD);

 88:   if (isbdiag) {
 89:     MatBDiagGetData(tmat,&nd,&bs,&diag,&bdlen,&diagv);
 90:     PetscPrintf(PETSC_COMM_WORLD,"transpose matrix: block diag format: %d diagonals, block size = %dn",nd,bs);
 91:     for (i=0; i<nd; i++) {
 92:       PetscPrintf(PETSC_COMM_WORLD," diag=%d, bdlength=%dn",diag[i],bdlen[i]);
 93:     }
 94:   }

 96:   /* Test MatAXPY */
 97:   if (mat && !rect) {
 98:     PetscScalar alpha = 1.0;
 99:     PetscOptionsGetScalar(PETSC_NULL,"-alpha",&alpha,PETSC_NULL);
100:     PetscPrintf(PETSC_COMM_WORLD,"matrix addition:  B = B + alpha * An");
101:     MatAXPY(&alpha,mat,tmat,DIFFERENT_NONZERO_PATTERN);
102:     MatView(tmat,PETSC_VIEWER_STDOUT_WORLD);
103:   }

105:   /* Free data structures */
106:   MatDestroy(tmat);
107:   if (mat) {MatDestroy(mat);}

109:   PetscFinalize();
110:   return 0;
111: }
112: