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