Actual source code: ex32.c

  2: static char help[] = "Reads in a matrix and vector in ASCII slap format. Writes\n\
  3: them using the PETSc sparse format. Input parameters are:\n\
  4:   -fin <filename> : input file\n\
  5:   -fout <filename> : output file\n\n";

 7:  #include petscmat.h

 11: int main(int argc,char **args)
 12: {
 13:   Mat            A;
 14:   Vec            b;
 15:   char           filein[PETSC_MAX_PATH_LEN],fileout[PETSC_MAX_PATH_LEN];
 16:   PetscInt       i,j,m,n,nnz,start,end,*col,*row,*brow,length;
 18:   PetscMPIInt    size,rank;
 19:   PetscScalar    *val,*bval;
 20:   FILE*          file;
 21:   PetscViewer    view;
 22:   PetscTruth     opt;

 24:   PetscInitialize(&argc,&args,(char *)0,help);

 26:   /* Read in matrix and RHS */
 27:   PetscOptionsGetString(PETSC_NULL,"-fin",filein,PETSC_MAX_PATH_LEN-1,&opt);
 28:   if (!opt) {
 29:     SETERRQ(PETSC_ERR_ARG_WRONG, "No filename was specified for this test");
 30:   }
 31:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 32:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 34:   PetscFOpen(PETSC_COMM_SELF,filein,"r",&file);

 36:   fscanf(file,"  NUNKNS =%d  NCOEFF =%d\n",&n,&nnz);
 37:   fscanf(file,"  JA POINTER IN SLAPSV\n");

 39:   MatCreateSeqAIJ(PETSC_COMM_WORLD,n,n,20,0,&A);
 40:   VecCreate(PETSC_COMM_WORLD,&b);
 41:   VecSetSizes(b,PETSC_DECIDE,n);
 42:   VecSetFromOptions(b);

 44:   PetscMalloc((n+1)*sizeof(PetscInt),&col);
 45:   for (i=0; i<n+1; i++)
 46:     fscanf(file,"     I=%d%d\n",&j,&col[i]);
 47:   fscanf(file,"  EOD JA\n");

 49:   PetscMalloc(nnz*sizeof(PetscScalar),&val);
 50:   PetscMalloc(nnz*sizeof(PetscInt),&row);
 51:   fscanf(file,"  COEFFICIENT MATRIX IN SLAPSV: I, IA, A\n");
 52:   for (i=0; i<nnz; i++) {
 53:     fscanf(file,"    %d%d%le\n",&j,&row[i],(double*)&val[i]);
 54:     row[i]--;
 55:   }
 56:   fscanf(file,"  EOD IA\n");

 58:   PetscMalloc(n*sizeof(PetscScalar),&bval);
 59:   PetscMalloc(n*sizeof(PetscInt),&brow);
 60:   fscanf(file,"  RESIDUAL IN SLAPSV ;IRHS=%d\n",&j);
 61:   for (i=0; i<n; i++) {
 62:     fscanf(file,"      %d%le%d\n",&j,(double*)(bval+i),&j);
 63:     brow[i] = i;
 64:   }
 65:   fscanf(file,"  EOD RESIDUAL");
 66:   fclose(file);

 68:   m = n/size+1;
 69:   start = rank*m;
 70:   end = (rank+1)*m; if (end > n) end = n;
 71:   for (j=start; j<end; j++) {
 72:     length = col[j+1]-col[j];
 73:     MatSetValues(A,length,&row[col[j]-1],1,&j,&val[col[j]-1],INSERT_VALUES);
 74:   }
 75:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 76:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 78:   VecGetOwnershipRange(b,&start,&end);
 79:   VecSetValues(b,end-start,brow+start,bval+start,INSERT_VALUES);
 80:   VecAssemblyBegin(b);
 81:   VecAssemblyEnd(b);

 83:   PetscFree(col);
 84:   PetscFree(val);
 85:   PetscFree(row);
 86:   PetscFree(bval);
 87:   PetscFree(brow);

 89:   PetscPrintf(PETSC_COMM_SELF,"Reading matrix completes.\n");
 90:   PetscOptionsGetString(PETSC_NULL,"-fout",fileout,PETSC_MAX_PATH_LEN-1,PETSC_NULL);
 91:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,fileout,FILE_MODE_WRITE,&view);
 92:   MatView(A,view);
 93:   VecView(b,view);
 94:   PetscViewerDestroy(view);

 96:   VecDestroy(b);
 97:   MatDestroy(A);

 99:   PetscFinalize();
100:   return 0;
101: }