Actual source code: ex5.c

  1: /*$Id: ex5.c,v 1.25 2001/04/10 19:35:44 bsmith Exp $*/
  2: 
  3: static char help[] = "Tests MatMult(), MatMultAdd(), MatMultTranspose().n
  4: Also MatMultTransposeAdd(), MatScale(), MatGetDiagonal(), and MatDiagonalScale().nn";

 6:  #include petscmat.h

  8: int main(int argc,char **args)
  9: {
 10:   Mat        C;
 11:   Vec        s,u,w,x,y,z;
 12:   int        ierr,i,j,m = 8,n,rstart,rend,vstart,vend;
 13:   Scalar     one = 1.0,negone = -1.0,v,alpha=0.1;
 14:   double     norm;
 15:   PetscTruth flg;

 17:   PetscInitialize(&argc,&args,(char *)0,help);
 18:   PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON);
 19:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 20:   n = m;
 21:   PetscOptionsHasName(PETSC_NULL,"-rectA",&flg);
 22:   if (flg) n += 2;
 23:   PetscOptionsHasName(PETSC_NULL,"-rectB",&flg);
 24:   if (flg) n -= 2;

 26:   /* ---------- Assemble matrix and vectors ----------- */

 28:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,n,&C);
 29:   MatSetFromOptions(C);
 30:   MatGetOwnershipRange(C,&rstart,&rend);
 31:   VecCreate(PETSC_COMM_WORLD,PETSC_DECIDE,m,&x);
 32:   VecSetFromOptions(x);
 33:   VecDuplicate(x,&z);
 34:   VecDuplicate(x,&w);
 35:   VecCreate(PETSC_COMM_WORLD,PETSC_DECIDE,n,&y);
 36:   VecSetFromOptions(y);
 37:   VecDuplicate(y,&u);
 38:   VecDuplicate(y,&s);
 39:   VecGetOwnershipRange(y,&vstart,&vend);

 41:   /* Assembly */
 42:   for (i=rstart; i<rend; i++) {
 43:     v = 100*(i+1);
 44:     VecSetValues(z,1,&i,&v,INSERT_VALUES);
 45:     for (j=0; j<n; j++) {
 46:       v=10*(i+1)+j+1;
 47:       MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES);
 48:     }
 49:   }

 51:   /* Flush off proc Vec values and do more assembly */
 52:   VecAssemblyBegin(z);
 53:   for (i=vstart; i<vend; i++) {
 54:     v = one*((double)i);
 55:     VecSetValues(y,1,&i,&v,INSERT_VALUES);
 56:     v = 100.0*i;
 57:     VecSetValues(u,1,&i,&v,INSERT_VALUES);
 58:   }

 60:   /* Flush off proc Mat values and do more assembly */
 61:   MatAssemblyBegin(C,MAT_FLUSH_ASSEMBLY);
 62:   for (i=rstart; i<rend; i++) {
 63:     for (j=0; j<n; j++) {
 64:       v=10*(i+1)+j+1;
 65:       MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES);
 66:     }
 67:   }
 68:   /* Try overlap Coomunication with the next stage XXXSetValues */
 69:   VecAssemblyEnd(z);
 70:   MatAssemblyEnd(C,MAT_FLUSH_ASSEMBLY);

 72:   /* The Assembly for the second Stage */
 73:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 74:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
 75:   VecAssemblyBegin(y);
 76:   VecAssemblyEnd(y);
 77:   MatScale(&alpha,C);
 78:   VecAssemblyBegin(u);
 79:   VecAssemblyEnd(u);

 81:   /* ------------ Test MatMult(), MatMultAdd()  ---------- */

 83:   PetscPrintf(PETSC_COMM_WORLD,"testing MatMult()n");
 84:   MatMult(C,y,x);
 85:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);
 86:   PetscPrintf(PETSC_COMM_WORLD,"testing MatMultAdd()n");
 87:   MatMultAdd(C,y,z,w);
 88:   VecAXPY(&one,z,x);
 89:   VecAXPY(&negone,w,x);
 90:   VecNorm(x,NORM_2,&norm);
 91:   if (norm > 1.e-8){
 92:     PetscPrintf(PETSC_COMM_WORLD,"Norm of error difference = %gn",norm);
 93:   }

 95:   /* ------- Test MatMultTranspose(), MatMultTransposeAdd() ------- */

 97:   for (i=rstart; i<rend; i++) {
 98:     v = one*((double)i);
 99:     VecSetValues(x,1,&i,&v,INSERT_VALUES);
100:   }
101:   VecAssemblyBegin(x);
102:   VecAssemblyEnd(x);
103:   PetscPrintf(PETSC_COMM_WORLD,"testing MatMultTranspose()n");
104:   MatMultTranspose(C,x,y);
105:   VecView(y,PETSC_VIEWER_STDOUT_WORLD);

107:   PetscPrintf(PETSC_COMM_WORLD,"testing MatMultTransposeAdd()n");
108:   MatMultTransposeAdd(C,x,u,s);
109:   VecAXPY(&one,u,y);
110:   VecAXPY(&negone,s,y);
111:   VecNorm(y,NORM_2,&norm);
112:   if (norm > 1.e-8){
113:     PetscPrintf(PETSC_COMM_WORLD,"Norm of error difference = %gn",norm);
114:   }

116:   /* -------------------- Test MatGetDiagonal() ------------------ */

118:   PetscPrintf(PETSC_COMM_WORLD,"testing MatGetDiagonal(), MatDiagonalScale()n");
119:   MatView(C,PETSC_VIEWER_STDOUT_WORLD);
120:   VecSet(&one,x);
121:   MatGetDiagonal(C,x);
122:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);
123:   for (i=vstart; i<vend; i++) {
124:     v = one*((double)(i+1));
125:     VecSetValues(y,1,&i,&v,INSERT_VALUES);
126:   }

128:   /* -------------------- Test () MatDiagonalScale ------------------ */
129:   PetscOptionsHasName(PETSC_NULL,"-test_diagonalscale",&flg);
130:   if (flg) {
131:     MatDiagonalScale(C,x,y);
132:     MatView(C,PETSC_VIEWER_STDOUT_WORLD);
133:   }
134:   /* Free data structures */
135:   VecDestroy(u); VecDestroy(s);
136:   VecDestroy(w); VecDestroy(x);
137:   VecDestroy(y); VecDestroy(z);
138:   MatDestroy(C);

140:   PetscFinalize();
141:   return 0;
142: }