xref: /petsc/src/mat/impls/shell/shellcnv.c (revision a5eb49655b3fdf389f9d65fc4214d6e1c240a941)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: shellcnv.c,v 1.3 1997/06/05 12:54:00 bsmith Exp balay $";
3 #endif
4 
5 
6 #include "petsc.h"
7 #include "src/mat/matimpl.h"        /*I "mat.h" I*/
8 #include "src/vec/vecimpl.h"
9 
10 #undef __FUNC__
11 #define __FUNC__ "MatConvert_Shell"
12 int MatConvert_Shell(Mat oldmat,MatType newtype, Mat *mat)
13 {
14   Vec      in,out;
15   int      ierr,i,M,m,size,*rows,start,end;
16   MPI_Comm comm;
17   Scalar   *array,zero = 0.0,one = 1.0;
18 
19   PetscValidHeaderSpecific(oldmat,MAT_COOKIE);
20   PetscValidPointer(mat);
21 
22   if (newtype != MATSEQDENSE || newtype != MATMPIDENSE) {
23     SETERRQ(PETSC_ERR_SUP,1,"Can only convert shell matrices to dense currently");
24   }
25   comm = oldmat->comm;
26 
27   MPI_Comm_size(comm,&size);
28 
29   ierr = MatGetOwnershipRange(oldmat,&start,&end); CHKERRQ(ierr);
30   ierr = VecCreateMPI(comm,end-start,PETSC_DECIDE,&in); CHKERRQ(ierr);
31   ierr = VecDuplicate(in,&out); CHKERRQ(ierr);
32   ierr = VecGetSize(in,&M); CHKERRQ(ierr);
33   ierr = VecGetLocalSize(in,&m); CHKERRQ(ierr);
34   rows = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(rows);
35   for ( i=0; i<m; i++ ) {rows[i] = start + i;}
36 
37   if (size == 1) {
38     ierr = MatCreateSeqDense(comm,M,M,PETSC_NULL,mat); CHKERRQ(ierr);
39   } else {
40     ierr = MatCreateMPIDense(comm,m,M,M,M,PETSC_NULL,mat); CHKERRQ(ierr);
41     /* ierr = MatCreateMPIAIJ(comm,m,m,M,M,0,0,0,0,mat); CHKERRQ(ierr); */
42   }
43 
44   for ( i=0; i<M; i++ ) {
45     ierr = VecSet(&zero,in); CHKERRQ(ierr);
46     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES); CHKERRQ(ierr);
47     ierr = VecAssemblyBegin(in); CHKERRQ(ierr);
48     ierr = VecAssemblyEnd(in); CHKERRQ(ierr);
49 
50     ierr = MatMult(oldmat,in,out); CHKERRQ(ierr);
51 
52     ierr = VecGetArray(out,&array); CHKERRQ(ierr);
53     ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
54     ierr = VecRestoreArray(out,&array); CHKERRQ(ierr);
55 
56   }
57   PetscFree(rows);
58   ierr = VecDestroy(in); CHKERRQ(ierr);
59   ierr = VecDestroy(out); CHKERRQ(ierr);
60   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
61   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
62   return 0;
63 }
64 
65 
66 
67