Actual source code: ex64.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

  2: static char help[] = "Saves 4by4 block matrix.\n\n";

  4:  #include <petscmat.h>

  6: int main(int argc,char **args)
  7: {
  8:   Mat            A;
  9:   PetscInt       i,j;
 11:   PetscMPIInt    size;
 12:   PetscViewer    fd;
 13:   PetscScalar    values[16],one = 1.0;
 14:   Vec            x;

 16:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 17:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 18:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,1,"Can only run on one processor");

 20:   /*
 21:      Open binary file.  Note that we use FILE_MODE_WRITE to indicate
 22:      writing to this file.
 23:   */
 24:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"4by4",FILE_MODE_WRITE,&fd);

 26:   MatCreateSeqBAIJ(PETSC_COMM_WORLD,4,12,12,0,0,&A);

 28:   for (i=0; i<16; i++) values[i] = i;
 29:   for (i=0; i<4; i++) values[4*i+i] += 5;
 30:   i    = 0; j = 0;
 31:   MatSetValuesBlocked(A,1,&i,1,&j,values,INSERT_VALUES);

 33:   for (i=0; i<16; i++) values[i] = i;
 34:   i    = 0; j = 2;
 35:   MatSetValuesBlocked(A,1,&i,1,&j,values,INSERT_VALUES);

 37:   for (i=0; i<16; i++) values[i] = i;
 38:   i    = 1; j = 0;
 39:   MatSetValuesBlocked(A,1,&i,1,&j,values,INSERT_VALUES);

 41:   for (i=0; i<16; i++) values[i] = i;for (i=0; i<4; i++) values[4*i+i] += 6;
 42:   i    = 1; j = 1;
 43:   MatSetValuesBlocked(A,1,&i,1,&j,values,INSERT_VALUES);

 45:   for (i=0; i<16; i++) values[i] = i;
 46:   i    = 2; j = 0;
 47:   MatSetValuesBlocked(A,1,&i,1,&j,values,INSERT_VALUES);

 49:   for (i=0; i<16; i++) values[i] = i;for (i=0; i<4; i++) values[4*i+i] += 7;
 50:   i    = 2; j = 2;
 51:   MatSetValuesBlocked(A,1,&i,1,&j,values,INSERT_VALUES);

 53:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 54:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 55:   MatView(A,fd);
 56:   MatDestroy(&A);

 58:   VecCreateSeq(PETSC_COMM_WORLD,12,&x);
 59:   VecSet(x,one);
 60:   VecView(x,fd);
 61:   VecDestroy(&x);

 63:   PetscViewerDestroy(&fd);
 64:   PetscFinalize();
 65:   return ierr;
 66: }