xref: /petsc/src/mat/tests/ex20.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests converting a matrix to another format with MatConvert().\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown int main(int argc,char **args)
7*c4762a1bSJed Brown {
8*c4762a1bSJed Brown   Mat            C,A;
9*c4762a1bSJed Brown   PetscInt       i,j,m = 5,n = 4,Ii,J;
10*c4762a1bSJed Brown   PetscErrorCode ierr;
11*c4762a1bSJed Brown   PetscMPIInt    rank,size;
12*c4762a1bSJed Brown   PetscScalar    v;
13*c4762a1bSJed Brown   char           mtype[256];
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
16*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
17*c4762a1bSJed Brown 
18*c4762a1bSJed Brown   /* This example does not work correctly for np > 2 */
19*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
20*c4762a1bSJed Brown   if (size > 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Use np <= 2");
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown   /* Create the matrix for the five point stencil, YET AGAIN */
23*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
24*c4762a1bSJed Brown   ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr);
25*c4762a1bSJed Brown   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
26*c4762a1bSJed Brown   ierr = MatSetUp(C);CHKERRQ(ierr);
27*c4762a1bSJed Brown   for (i=0; i<m; i++) {
28*c4762a1bSJed Brown     for (j=2*rank; j<2*rank+2; j++) {
29*c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
30*c4762a1bSJed Brown       if (i>0)   {J = Ii - n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
31*c4762a1bSJed Brown       if (i<m-1) {J = Ii + n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
32*c4762a1bSJed Brown       if (j>0)   {J = Ii - 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
33*c4762a1bSJed Brown       if (j<n-1) {J = Ii + 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
34*c4762a1bSJed Brown       v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr);
35*c4762a1bSJed Brown     }
36*c4762a1bSJed Brown   }
37*c4762a1bSJed Brown 
38*c4762a1bSJed Brown   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
43*c4762a1bSJed Brown   ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
44*c4762a1bSJed Brown 
45*c4762a1bSJed Brown   ierr = PetscStrcpy(mtype,MATSAME);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = PetscOptionsGetString(NULL,NULL,"-conv_mat_type",mtype,256,NULL);CHKERRQ(ierr);
47*c4762a1bSJed Brown   ierr = MatConvert(C,mtype,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
48*c4762a1bSJed Brown   ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
49*c4762a1bSJed Brown   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
50*c4762a1bSJed Brown   ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
51*c4762a1bSJed Brown   ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_IMPL);CHKERRQ(ierr);
52*c4762a1bSJed Brown   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
53*c4762a1bSJed Brown   ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown   /* Free data structures */
56*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
58*c4762a1bSJed Brown 
59*c4762a1bSJed Brown   ierr = PetscFinalize();
60*c4762a1bSJed Brown   return ierr;
61*c4762a1bSJed Brown }
62*c4762a1bSJed Brown 
63*c4762a1bSJed Brown 
64*c4762a1bSJed Brown /*TEST
65*c4762a1bSJed Brown 
66*c4762a1bSJed Brown    test:
67*c4762a1bSJed Brown       args: -conv_mat_type seqaij
68*c4762a1bSJed Brown 
69*c4762a1bSJed Brown TEST*/
70