xref: /petsc/src/mat/impls/shell/shellcnv.c (revision 6098dc10f64a9429b4a6e4cb977ddd6c38617718)
1*6098dc10SBarry Smith #ifndef lint
2*6098dc10SBarry Smith static char vcid[] = "$Id: shell.c,v 1.43 1997/01/12 04:33:57 bsmith Exp bsmith $";
3*6098dc10SBarry Smith #endif
4*6098dc10SBarry Smith 
5*6098dc10SBarry Smith /*
6*6098dc10SBarry Smith    This provides a simple shell for Fortran (and C programmers) to
7*6098dc10SBarry Smith   create a very simple matrix class for use with KSP without coding
8*6098dc10SBarry Smith   much of anything.
9*6098dc10SBarry Smith */
10*6098dc10SBarry Smith 
11*6098dc10SBarry Smith #include "petsc.h"
12*6098dc10SBarry Smith #include "src/mat/matimpl.h"        /*I "mat.h" I*/
13*6098dc10SBarry Smith #include "src/vec/vecimpl.h"
14*6098dc10SBarry Smith 
15*6098dc10SBarry Smith typedef struct {
16*6098dc10SBarry Smith   int  M, N;                  /* number of global rows, columns */
17*6098dc10SBarry Smith   int  m, n;                  /* number of local rows, columns */
18*6098dc10SBarry Smith   int  (*destroy)(Mat);
19*6098dc10SBarry Smith   void *ctx;
20*6098dc10SBarry Smith } Mat_Shell;
21*6098dc10SBarry Smith 
22*6098dc10SBarry Smith #undef __FUNC__
23*6098dc10SBarry Smith #define __FUNC__ "MatShellGetContext"
24*6098dc10SBarry Smith /*@
25*6098dc10SBarry Smith     MatShellGetContext - Returns the user-provided context associated with a shell matrix.
26*6098dc10SBarry Smith 
27*6098dc10SBarry Smith     Input Parameter:
28*6098dc10SBarry Smith .   mat - the matrix, should have been created with MatCreateShell()
29*6098dc10SBarry Smith 
30*6098dc10SBarry Smith     Output Parameter:
31*6098dc10SBarry Smith .   ctx - the user provided context
32*6098dc10SBarry Smith 
33*6098dc10SBarry Smith     Notes:
34*6098dc10SBarry Smith     This routine is intended for use within various shell matrix routines,
35*6098dc10SBarry Smith     as set with MatShellSetOperation().
36*6098dc10SBarry Smith 
37*6098dc10SBarry Smith .keywords: matrix, shell, get, context
38*6098dc10SBarry Smith 
39*6098dc10SBarry Smith .seealso: MatCreateShell(), MatShellSetOperation()
40*6098dc10SBarry Smith @*/
41*6098dc10SBarry Smith int MatShellGetContext(Mat mat,void **ctx)
42*6098dc10SBarry Smith {
43*6098dc10SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
44*6098dc10SBarry Smith   if (mat->type != MATSHELL) *ctx = 0;
45*6098dc10SBarry Smith   else                       *ctx = ((Mat_Shell *) (mat->data))->ctx;
46*6098dc10SBarry Smith   return 0;
47*6098dc10SBarry Smith }
48*6098dc10SBarry Smith 
49*6098dc10SBarry Smith #undef __FUNC__
50*6098dc10SBarry Smith #define __FUNC__ "MatGetSize_Shell"
51*6098dc10SBarry Smith static int MatGetSize_Shell(Mat mat,int *M,int *N)
52*6098dc10SBarry Smith {
53*6098dc10SBarry Smith   Mat_Shell *shell = (Mat_Shell *) mat->data;
54*6098dc10SBarry Smith   *M = shell->M; *N = shell->N;
55*6098dc10SBarry Smith   return 0;
56*6098dc10SBarry Smith }
57*6098dc10SBarry Smith 
58*6098dc10SBarry Smith #undef __FUNC__
59*6098dc10SBarry Smith #define __FUNC__ "MatGetLocalSize_Shell"
60*6098dc10SBarry Smith static int MatGetLocalSize_Shell(Mat mat,int *m,int *n)
61*6098dc10SBarry Smith {
62*6098dc10SBarry Smith   Mat_Shell *shell = (Mat_Shell *) mat->data;
63*6098dc10SBarry Smith   *m = shell->n; *n = shell->n;
64*6098dc10SBarry Smith   return 0;
65*6098dc10SBarry Smith }
66*6098dc10SBarry Smith 
67*6098dc10SBarry Smith #undef __FUNC__
68*6098dc10SBarry Smith #define __FUNC__ "MatDestroy_Shell"
69*6098dc10SBarry Smith static int MatDestroy_Shell(PetscObject obj)
70*6098dc10SBarry Smith {
71*6098dc10SBarry Smith   int       ierr;
72*6098dc10SBarry Smith   Mat       mat = (Mat) obj;
73*6098dc10SBarry Smith   Mat_Shell *shell;
74*6098dc10SBarry Smith 
75*6098dc10SBarry Smith   shell = (Mat_Shell *) mat->data;
76*6098dc10SBarry Smith   if (shell->destroy) {ierr = (*shell->destroy)(mat);CHKERRQ(ierr);}
77*6098dc10SBarry Smith   PetscFree(shell);
78*6098dc10SBarry Smith   PLogObjectDestroy(mat);
79*6098dc10SBarry Smith   PetscHeaderDestroy(mat);
80*6098dc10SBarry Smith   return 0;
81*6098dc10SBarry Smith }
82*6098dc10SBarry Smith 
83*6098dc10SBarry Smith #undef __FUNC__
84*6098dc10SBarry Smith #define __FUNC__ "MatConvert_Shell"
85*6098dc10SBarry Smith int MatConvert_Shell(Mat oldmat,MatType newtype, Mat *mat)
86*6098dc10SBarry Smith {
87*6098dc10SBarry Smith   Vec      in,out;
88*6098dc10SBarry Smith   int      ierr,i,M,m,size,*rows,start,end;
89*6098dc10SBarry Smith   MPI_Comm comm;
90*6098dc10SBarry Smith   Scalar   *array,zero = 0.0,one = 1.0;
91*6098dc10SBarry Smith 
92*6098dc10SBarry Smith   PetscValidHeaderSpecific(oldmat,MAT_COOKIE);
93*6098dc10SBarry Smith   PetscValidPointer(mat);
94*6098dc10SBarry Smith 
95*6098dc10SBarry Smith   if (newtype != MATSEQDENSE || newtype != MATMPIDENSE) {
96*6098dc10SBarry Smith     SETERRQ(PETSC_ERR_SUP,1,"Can only convert shell matrices to dense currently");
97*6098dc10SBarry Smith   }
98*6098dc10SBarry Smith   comm = oldmat->comm;
99*6098dc10SBarry Smith 
100*6098dc10SBarry Smith   MPI_Comm_size(comm,&size);
101*6098dc10SBarry Smith 
102*6098dc10SBarry Smith   ierr = MatGetOwnershipRange(oldmat,&start,&end); CHKERRQ(ierr);
103*6098dc10SBarry Smith   ierr = VecCreateMPI(comm,end-start,PETSC_DECIDE,&in); CHKERRQ(ierr);
104*6098dc10SBarry Smith   ierr = VecDuplicate(in,&out); CHKERRQ(ierr);
105*6098dc10SBarry Smith   ierr = VecGetSize(in,&M); CHKERRQ(ierr);
106*6098dc10SBarry Smith   ierr = VecGetLocalSize(in,&m); CHKERRQ(ierr);
107*6098dc10SBarry Smith   rows = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(rows);
108*6098dc10SBarry Smith   for ( i=0; i<m; i++ ) {rows[i] = start + i;}
109*6098dc10SBarry Smith 
110*6098dc10SBarry Smith   if (size == 1) {
111*6098dc10SBarry Smith     ierr = MatCreateSeqDense(comm,M,M,PETSC_NULL,mat); CHKERRQ(ierr);
112*6098dc10SBarry Smith   } else {
113*6098dc10SBarry Smith     ierr = MatCreateMPIDense(comm,m,M,M,M,PETSC_NULL,mat); CHKERRQ(ierr);
114*6098dc10SBarry Smith     /* ierr = MatCreateMPIAIJ(comm,m,m,M,M,0,0,0,0,mat); CHKERRQ(ierr); */
115*6098dc10SBarry Smith   }
116*6098dc10SBarry Smith 
117*6098dc10SBarry Smith   for ( i=0; i<M; i++ ) {
118*6098dc10SBarry Smith 
119*6098dc10SBarry Smith     ierr = VecSet(&zero,in); CHKERRQ(ierr);
120*6098dc10SBarry Smith     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES); CHKERRQ(ierr);
121*6098dc10SBarry Smith     ierr = VecAssemblyBegin(in); CHKERRQ(ierr);
122*6098dc10SBarry Smith     ierr = VecAssemblyEnd(in); CHKERRQ(ierr);
123*6098dc10SBarry Smith 
124*6098dc10SBarry Smith     ierr = MatMult(oldmat,in,out); CHKERRQ(ierr);
125*6098dc10SBarry Smith 
126*6098dc10SBarry Smith     ierr = VecGetArray(out,&array); CHKERRQ(ierr);
127*6098dc10SBarry Smith     ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
128*6098dc10SBarry Smith     ierr = VecRestoreArray(out,&array); CHKERRQ(ierr);
129*6098dc10SBarry Smith 
130*6098dc10SBarry Smith   }
131*6098dc10SBarry Smith   PetscFree(rows);
132*6098dc10SBarry Smith   ierr = VecDestroy(in); CHKERRQ(ierr);
133*6098dc10SBarry Smith   ierr = VecDestroy(out); CHKERRQ(ierr);
134*6098dc10SBarry Smith   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
135*6098dc10SBarry Smith   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
136*6098dc10SBarry Smith   return 0;
137*6098dc10SBarry Smith }
138*6098dc10SBarry Smith 
139*6098dc10SBarry Smith static int MatGetOwnershipRange_Shell(Mat mat, int *rstart,int *rend)
140*6098dc10SBarry Smith {
141*6098dc10SBarry Smith   MPI_Scan(&mat->m,rend,1,MPI_INT,MPI_SUM,mat->comm);
142*6098dc10SBarry Smith   *rstart = *rend - mat->m;
143*6098dc10SBarry Smith   return 0;
144*6098dc10SBarry Smith }
145*6098dc10SBarry Smith 
146*6098dc10SBarry Smith 
147*6098dc10SBarry Smith 
148*6098dc10SBarry Smith 
149*6098dc10SBarry Smith static struct _MatOps MatOps = {0,
150*6098dc10SBarry Smith        0,
151*6098dc10SBarry Smith        0,
152*6098dc10SBarry Smith        0,
153*6098dc10SBarry Smith        0,
154*6098dc10SBarry Smith        0,
155*6098dc10SBarry Smith        0,
156*6098dc10SBarry Smith        0,
157*6098dc10SBarry Smith        0,
158*6098dc10SBarry Smith        0,
159*6098dc10SBarry Smith        0,
160*6098dc10SBarry Smith        0,
161*6098dc10SBarry Smith        0,
162*6098dc10SBarry Smith        0,
163*6098dc10SBarry Smith        0,
164*6098dc10SBarry Smith        0,
165*6098dc10SBarry Smith        0,
166*6098dc10SBarry Smith        0,
167*6098dc10SBarry Smith        0,
168*6098dc10SBarry Smith        0,
169*6098dc10SBarry Smith        0,
170*6098dc10SBarry Smith        0,
171*6098dc10SBarry Smith        0,
172*6098dc10SBarry Smith        0,
173*6098dc10SBarry Smith        0,
174*6098dc10SBarry Smith        0,
175*6098dc10SBarry Smith        0,
176*6098dc10SBarry Smith        0,
177*6098dc10SBarry Smith        0,
178*6098dc10SBarry Smith        0,
179*6098dc10SBarry Smith        MatGetSize_Shell,
180*6098dc10SBarry Smith        MatGetLocalSize_Shell,
181*6098dc10SBarry Smith        MatGetOwnershipRange_Shell,
182*6098dc10SBarry Smith        0,
183*6098dc10SBarry Smith        0,
184*6098dc10SBarry Smith        0,
185*6098dc10SBarry Smith        0,
186*6098dc10SBarry Smith        0 };
187*6098dc10SBarry Smith 
188*6098dc10SBarry Smith #undef __FUNC__
189*6098dc10SBarry Smith #define __FUNC__ "MatCreateShell"
190*6098dc10SBarry Smith /*@C
191*6098dc10SBarry Smith    MatCreateShell - Creates a new matrix class for use with a user-defined
192*6098dc10SBarry Smith    private data storage format.
193*6098dc10SBarry Smith 
194*6098dc10SBarry Smith    Input Parameters:
195*6098dc10SBarry Smith .  comm - MPI communicator
196*6098dc10SBarry Smith .  m - number of local rows
197*6098dc10SBarry Smith .  n - number of local columns
198*6098dc10SBarry Smith .  M - number of global rows
199*6098dc10SBarry Smith .  N - number of global columns
200*6098dc10SBarry Smith .  ctx - pointer to data needed by the shell matrix routines
201*6098dc10SBarry Smith 
202*6098dc10SBarry Smith    Output Parameter:
203*6098dc10SBarry Smith .  A - the matrix
204*6098dc10SBarry Smith 
205*6098dc10SBarry Smith    Usage:
206*6098dc10SBarry Smith $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
207*6098dc10SBarry Smith $    MatShellSetOperation(mat,MATOP_MULT,mult);
208*6098dc10SBarry Smith $    [ Use matrix for operations that have been set ]
209*6098dc10SBarry Smith $    MatDestroy(mat);
210*6098dc10SBarry Smith 
211*6098dc10SBarry Smith    Notes:
212*6098dc10SBarry Smith    The shell matrix type is intended to provide a simple class to use
213*6098dc10SBarry Smith    with KSP (such as, for use with matrix-free methods). You should not
214*6098dc10SBarry Smith    use the shell type if you plan to define a complete matrix class.
215*6098dc10SBarry Smith 
216*6098dc10SBarry Smith    PETSc requires that matrices and vectors being used for certain
217*6098dc10SBarry Smith    operations are partitioned accordingly.  For example, when
218*6098dc10SBarry Smith    creating a shell matrix, A, that supports parallel matrix-vector
219*6098dc10SBarry Smith    products using MatMult(A,x,y) the user should set the number
220*6098dc10SBarry Smith    of local matrix rows to be the number of local elements of the
221*6098dc10SBarry Smith    corresponding result vector, y. Note that this is information is
222*6098dc10SBarry Smith    required for use of the matrix interface routines, even though
223*6098dc10SBarry Smith    the shell matrix may not actually be physically partitioned.
224*6098dc10SBarry Smith    For example,
225*6098dc10SBarry Smith 
226*6098dc10SBarry Smith $
227*6098dc10SBarry Smith $     Vec x, y
228*6098dc10SBarry Smith $     Mat A
229*6098dc10SBarry Smith $
230*6098dc10SBarry Smith $     VecCreate(comm,M,&y);
231*6098dc10SBarry Smith $     VecCreate(comm,N,&x);
232*6098dc10SBarry Smith $     VecGetLocalSize(y,&m);
233*6098dc10SBarry Smith $     MatCreateShell(comm,m,N,M,N,ctx,&A);
234*6098dc10SBarry Smith $     MatShellSetOperation(mat,MATOP_MULT,mult);
235*6098dc10SBarry Smith $     MatMult(A,x,y);
236*6098dc10SBarry Smith $     MatDestroy(A);
237*6098dc10SBarry Smith $     VecDestroy(y); VecDestroy(x);
238*6098dc10SBarry Smith $
239*6098dc10SBarry Smith 
240*6098dc10SBarry Smith .keywords: matrix, shell, create
241*6098dc10SBarry Smith 
242*6098dc10SBarry Smith .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext()
243*6098dc10SBarry Smith @*/
244*6098dc10SBarry Smith int MatCreateShell(MPI_Comm comm,int m,int n,int M,int N,void *ctx,Mat *A)
245*6098dc10SBarry Smith {
246*6098dc10SBarry Smith   Mat       B;
247*6098dc10SBarry Smith   Mat_Shell *b;
248*6098dc10SBarry Smith 
249*6098dc10SBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSHELL,comm);
250*6098dc10SBarry Smith   PLogObjectCreate(B);
251*6098dc10SBarry Smith   B->factor    = 0;
252*6098dc10SBarry Smith   B->destroy   = MatDestroy_Shell;
253*6098dc10SBarry Smith   B->assembled = PETSC_TRUE;
254*6098dc10SBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
255*6098dc10SBarry Smith 
256*6098dc10SBarry Smith   b          = PetscNew(Mat_Shell); CHKPTRQ(b);
257*6098dc10SBarry Smith   PetscMemzero(b,sizeof(Mat_Shell));
258*6098dc10SBarry Smith   B->data   = (void *) b;
259*6098dc10SBarry Smith   b->M = M; B->M = M;
260*6098dc10SBarry Smith   b->N = N; B->N = N;
261*6098dc10SBarry Smith   b->m = m; B->m = m;
262*6098dc10SBarry Smith   b->n = n; B->n = n;
263*6098dc10SBarry Smith   b->ctx     = ctx;
264*6098dc10SBarry Smith   *A = B;
265*6098dc10SBarry Smith   return 0;
266*6098dc10SBarry Smith }
267*6098dc10SBarry Smith 
268*6098dc10SBarry Smith #undef __FUNC__
269*6098dc10SBarry Smith #define __FUNC__ "MatShellSetOperation"
270*6098dc10SBarry Smith /*@C
271*6098dc10SBarry Smith     MatShellSetOperation - Allows user to set a matrix operation for
272*6098dc10SBarry Smith                            a shell matrix.
273*6098dc10SBarry Smith 
274*6098dc10SBarry Smith     Input Parameters:
275*6098dc10SBarry Smith .   mat - the shell matrix
276*6098dc10SBarry Smith .   op - the name of the operation
277*6098dc10SBarry Smith .   f - the function that provides the operation.
278*6098dc10SBarry Smith 
279*6098dc10SBarry Smith     Usage:
280*6098dc10SBarry Smith $      extern int usermult(Mat,Vec,Vec);
281*6098dc10SBarry Smith $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
282*6098dc10SBarry Smith $      ierr = MatShellSetOperation(A,MATOP_MULT,usermult);
283*6098dc10SBarry Smith 
284*6098dc10SBarry Smith     Notes:
285*6098dc10SBarry Smith     See the file petsc/include/mat.h for a complete list of matrix
286*6098dc10SBarry Smith     operations, which all have the form MATOP_<OPERATION>, where
287*6098dc10SBarry Smith     <OPERATION> is the name (in all capital letters) of the
288*6098dc10SBarry Smith     user interface routine (e.g., MatMult() -> MATOP_MULT).
289*6098dc10SBarry Smith 
290*6098dc10SBarry Smith     All user-provided functions should have the same calling
291*6098dc10SBarry Smith     sequence as the usual matrix interface routines, since they
292*6098dc10SBarry Smith     are intended to be accessed via the usual matrix interface
293*6098dc10SBarry Smith     routines, e.g.,
294*6098dc10SBarry Smith $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
295*6098dc10SBarry Smith 
296*6098dc10SBarry Smith     Within each user-defined routine, the user should call
297*6098dc10SBarry Smith     MatShellGetContext() to obtain the user-defined context that was
298*6098dc10SBarry Smith     set by MatCreateShell().
299*6098dc10SBarry Smith 
300*6098dc10SBarry Smith .keywords: matrix, shell, set, operation
301*6098dc10SBarry Smith 
302*6098dc10SBarry Smith .seealso: MatCreateShell(), MatShellGetContext()
303*6098dc10SBarry Smith @*/
304*6098dc10SBarry Smith int MatShellSetOperation(Mat mat,MatOperation op, void *f)
305*6098dc10SBarry Smith {
306*6098dc10SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
307*6098dc10SBarry Smith 
308*6098dc10SBarry Smith   if (op == MATOP_DESTROY) {
309*6098dc10SBarry Smith     if (mat->type == MATSHELL) {
310*6098dc10SBarry Smith        Mat_Shell *shell = (Mat_Shell *) mat->data;
311*6098dc10SBarry Smith        shell->destroy                 = (int (*)(Mat)) f;
312*6098dc10SBarry Smith     }
313*6098dc10SBarry Smith     else mat->destroy                 = (int (*)(PetscObject)) f;
314*6098dc10SBarry Smith   }
315*6098dc10SBarry Smith   else if (op == MATOP_VIEW) mat->view  = (int (*)(PetscObject,Viewer)) f;
316*6098dc10SBarry Smith   else      (((void**)&mat->ops)[op]) = f;
317*6098dc10SBarry Smith 
318*6098dc10SBarry Smith   return 0;
319*6098dc10SBarry Smith }
320*6098dc10SBarry Smith 
321*6098dc10SBarry Smith 
322*6098dc10SBarry Smith 
323