xref: /petsc/src/mat/impls/shell/shellcnv.c (revision e090d5668ba2b2ea997ebb925e3a05be0dc5d9ab)
1*e090d566SSatish Balay /*$Id: shellcnv.c,v 1.13 2000/04/12 04:23:17 bsmith Exp balay $*/
26098dc10SBarry Smith 
3*e090d566SSatish Balay #include "src/mat/matimpl.h"        /*I "petscmat.h" I*/
46098dc10SBarry Smith #include "src/vec/vecimpl.h"
56098dc10SBarry Smith 
66098dc10SBarry Smith #undef __FUNC__
7b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatConvert_Shell"
86098dc10SBarry Smith int MatConvert_Shell(Mat oldmat,MatType newtype,Mat *mat)
96098dc10SBarry Smith {
106098dc10SBarry Smith   Vec      in,out;
116098dc10SBarry Smith   int      ierr,i,M,m,size,*rows,start,end;
126098dc10SBarry Smith   MPI_Comm comm;
136098dc10SBarry Smith   Scalar   *array,zero = 0.0,one = 1.0;
146098dc10SBarry Smith 
153a40ed3dSBarry Smith   PetscFunctionBegin;
166098dc10SBarry Smith   PetscValidHeaderSpecific(oldmat,MAT_COOKIE);
176098dc10SBarry Smith   PetscValidPointer(mat);
186098dc10SBarry Smith 
19f1af5d2fSBarry Smith   if (newtype != MATSEQDENSE && newtype != MATMPIDENSE) {
206098dc10SBarry Smith     SETERRQ(PETSC_ERR_SUP,1,"Can only convert shell matrices to dense currently");
216098dc10SBarry Smith   }
226098dc10SBarry Smith   comm = oldmat->comm;
236098dc10SBarry Smith 
24d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
256098dc10SBarry Smith 
266098dc10SBarry Smith   ierr = MatGetOwnershipRange(oldmat,&start,&end);CHKERRQ(ierr);
276098dc10SBarry Smith   ierr = VecCreateMPI(comm,end-start,PETSC_DECIDE,&in);CHKERRQ(ierr);
286098dc10SBarry Smith   ierr = VecDuplicate(in,&out);CHKERRQ(ierr);
296098dc10SBarry Smith   ierr = VecGetSize(in,&M);CHKERRQ(ierr);
306098dc10SBarry Smith   ierr = VecGetLocalSize(in,&m);CHKERRQ(ierr);
316098dc10SBarry Smith   rows = (int*)PetscMalloc((m+1)*sizeof(int));CHKPTRQ(rows);
326098dc10SBarry Smith   for (i=0; i<m; i++) {rows[i] = start + i;}
336098dc10SBarry Smith 
346098dc10SBarry Smith   if (size == 1) {
356098dc10SBarry Smith     ierr = MatCreateSeqDense(comm,M,M,PETSC_NULL,mat);CHKERRQ(ierr);
366098dc10SBarry Smith   } else {
376098dc10SBarry Smith     ierr = MatCreateMPIDense(comm,m,M,M,M,PETSC_NULL,mat);CHKERRQ(ierr);
386098dc10SBarry Smith     /* ierr = MatCreateMPIAIJ(comm,m,m,M,M,0,0,0,0,mat);CHKERRQ(ierr); */
396098dc10SBarry Smith   }
406098dc10SBarry Smith 
416098dc10SBarry Smith   for (i=0; i<M; i++) {
426098dc10SBarry Smith     ierr = VecSet(&zero,in);CHKERRQ(ierr);
436098dc10SBarry Smith     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr);
446098dc10SBarry Smith     ierr = VecAssemblyBegin(in);CHKERRQ(ierr);
456098dc10SBarry Smith     ierr = VecAssemblyEnd(in);CHKERRQ(ierr);
466098dc10SBarry Smith 
476098dc10SBarry Smith     ierr = MatMult(oldmat,in,out);CHKERRQ(ierr);
486098dc10SBarry Smith 
496098dc10SBarry Smith     ierr = VecGetArray(out,&array);CHKERRQ(ierr);
506098dc10SBarry Smith     ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr);
516098dc10SBarry Smith     ierr = VecRestoreArray(out,&array);CHKERRQ(ierr);
526098dc10SBarry Smith 
536098dc10SBarry Smith   }
54606d414cSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
556098dc10SBarry Smith   ierr = VecDestroy(in);CHKERRQ(ierr);
566098dc10SBarry Smith   ierr = VecDestroy(out);CHKERRQ(ierr);
576098dc10SBarry Smith   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
586098dc10SBarry Smith   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
593a40ed3dSBarry Smith   PetscFunctionReturn(0);
606098dc10SBarry Smith }
616098dc10SBarry Smith 
626098dc10SBarry Smith 
636098dc10SBarry Smith 
64