xref: /petsc/src/mat/impls/is/matis.c (revision b7ce53b6f5b49c6ff22bb6f08f6eff23909c7dbf)
1be1d678aSKris Buschelman 
2b4319ba4SBarry Smith /*
3b4319ba4SBarry Smith     Creates a matrix class for using the Neumann-Neumann type preconditioners.
4b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
5b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
6b4319ba4SBarry Smith    product is handled "implicitly".
7b4319ba4SBarry Smith 
8b4319ba4SBarry Smith      We provide:
9b4319ba4SBarry Smith          MatMult()
10b4319ba4SBarry Smith 
11b4319ba4SBarry Smith     Currently this allows for only one subdomain per processor.
12b4319ba4SBarry Smith 
13b4319ba4SBarry Smith */
14b4319ba4SBarry Smith 
15c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h>      /*I "petscmat.h" I*/
16b4319ba4SBarry Smith 
17b4319ba4SBarry Smith #undef __FUNCT__
18*b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS"
19*b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M)
20*b7ce53b6SStefano Zampini {
21*b7ce53b6SStefano Zampini   Mat_IS                 *matis = (Mat_IS*)(mat->data);
22*b7ce53b6SStefano Zampini   /* info on mat */
23*b7ce53b6SStefano Zampini   /* ISLocalToGlobalMapping rmapping,cmapping; */
24*b7ce53b6SStefano Zampini   PetscInt               bs,rows,cols;
25*b7ce53b6SStefano Zampini   PetscInt               lrows,lcols;
26*b7ce53b6SStefano Zampini   PetscInt               local_rows,local_cols;
27*b7ce53b6SStefano Zampini   PetscBool              isdense,issbaij,issbaij_red;
28*b7ce53b6SStefano Zampini   /* values insertion */
29*b7ce53b6SStefano Zampini   PetscScalar            *array;
30*b7ce53b6SStefano Zampini   PetscInt               *local_indices,*global_indices;
31*b7ce53b6SStefano Zampini   /* work */
32*b7ce53b6SStefano Zampini   PetscInt               i,j,index_row;
33*b7ce53b6SStefano Zampini   PetscErrorCode         ierr;
34*b7ce53b6SStefano Zampini 
35*b7ce53b6SStefano Zampini   PetscFunctionBegin;
36*b7ce53b6SStefano Zampini   /* MISSING CHECKS
37*b7ce53b6SStefano Zampini     - rectangular case not covered (it is not allowed by MATIS)
38*b7ce53b6SStefano Zampini   */
39*b7ce53b6SStefano Zampini   /* get info from mat */
40*b7ce53b6SStefano Zampini   /* ierr = MatGetLocalToGlobalMapping(mat,&rmapping,&cmapping);CHKERRQ(ierr); */
41*b7ce53b6SStefano Zampini   ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
42*b7ce53b6SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
43*b7ce53b6SStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
44*b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
45*b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
46*b7ce53b6SStefano Zampini 
47*b7ce53b6SStefano Zampini   /* work */
48*b7ce53b6SStefano Zampini   ierr = PetscMalloc1(local_rows,&local_indices);CHKERRQ(ierr);
49*b7ce53b6SStefano Zampini   for (i=0;i<local_rows;i++) local_indices[i]=i;
50*b7ce53b6SStefano Zampini   /* map indices of local mat to global values */
51*b7ce53b6SStefano Zampini   ierr = PetscMalloc(PetscMax(local_cols,local_rows)*sizeof(*global_indices),&global_indices);CHKERRQ(ierr);
52*b7ce53b6SStefano Zampini   /* ierr = ISLocalToGlobalMappingApply(rmapping,local_rows,local_indices,global_indices);CHKERRQ(ierr); */
53*b7ce53b6SStefano Zampini   ierr = ISLocalToGlobalMappingApply(matis->mapping,local_rows,local_indices,global_indices);CHKERRQ(ierr);
54*b7ce53b6SStefano Zampini 
55*b7ce53b6SStefano Zampini   if (issbaij) {
56*b7ce53b6SStefano Zampini     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
57*b7ce53b6SStefano Zampini   }
58*b7ce53b6SStefano Zampini 
59*b7ce53b6SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
60*b7ce53b6SStefano Zampini     Mat         new_mat;
61*b7ce53b6SStefano Zampini     MatType     new_mat_type;
62*b7ce53b6SStefano Zampini     Vec         vec_dnz,vec_onz;
63*b7ce53b6SStefano Zampini     PetscScalar *my_dnz,*my_onz;
64*b7ce53b6SStefano Zampini     PetscInt    *dnz,*onz,*mat_ranges,*row_ownership;
65*b7ce53b6SStefano Zampini     PetscInt    index_col,owner;
66*b7ce53b6SStefano Zampini     PetscMPIInt nsubdomains;
67*b7ce53b6SStefano Zampini 
68*b7ce53b6SStefano Zampini     /* determining new matrix type */
69*b7ce53b6SStefano Zampini     ierr = MPI_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
70*b7ce53b6SStefano Zampini     if (issbaij_red) {
71*b7ce53b6SStefano Zampini       new_mat_type = MATSBAIJ;
72*b7ce53b6SStefano Zampini     } else {
73*b7ce53b6SStefano Zampini       if (bs>1) {
74*b7ce53b6SStefano Zampini         new_mat_type = MATBAIJ;
75*b7ce53b6SStefano Zampini       } else {
76*b7ce53b6SStefano Zampini         new_mat_type = MATAIJ;
77*b7ce53b6SStefano Zampini       }
78*b7ce53b6SStefano Zampini     }
79*b7ce53b6SStefano Zampini 
80*b7ce53b6SStefano Zampini     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr);
81*b7ce53b6SStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)mat),&new_mat);CHKERRQ(ierr);
82*b7ce53b6SStefano Zampini     ierr = MatSetSizes(new_mat,PETSC_DECIDE,PETSC_DECIDE,rows,cols);CHKERRQ(ierr);
83*b7ce53b6SStefano Zampini     ierr = MatSetBlockSize(new_mat,bs);CHKERRQ(ierr);
84*b7ce53b6SStefano Zampini     ierr = MatSetType(new_mat,new_mat_type);CHKERRQ(ierr);
85*b7ce53b6SStefano Zampini     ierr = MatSetUp(new_mat);CHKERRQ(ierr);
86*b7ce53b6SStefano Zampini     ierr = MatGetLocalSize(new_mat,&lrows,&lcols);CHKERRQ(ierr);
87*b7ce53b6SStefano Zampini 
88*b7ce53b6SStefano Zampini     /*
89*b7ce53b6SStefano Zampini       preallocation
90*b7ce53b6SStefano Zampini     */
91*b7ce53b6SStefano Zampini 
92*b7ce53b6SStefano Zampini     ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)new_mat),lrows,lcols,dnz,onz);CHKERRQ(ierr);
93*b7ce53b6SStefano Zampini     /*
94*b7ce53b6SStefano Zampini        Some vectors are needed to sum up properly on shared interface dofs.
95*b7ce53b6SStefano Zampini        Preallocation macros cannot do the job.
96*b7ce53b6SStefano Zampini        Note that preallocation is not exact, since it overestimates nonzeros
97*b7ce53b6SStefano Zampini     */
98*b7ce53b6SStefano Zampini     ierr = MatGetVecs(new_mat,NULL,&vec_dnz);CHKERRQ(ierr);
99*b7ce53b6SStefano Zampini     /* ierr = VecSetLocalToGlobalMapping(vec_dnz,rmapping);CHKERRQ(ierr); */
100*b7ce53b6SStefano Zampini     ierr = VecSetLocalToGlobalMapping(vec_dnz,matis->mapping);CHKERRQ(ierr);
101*b7ce53b6SStefano Zampini     ierr = VecDuplicate(vec_dnz,&vec_onz);CHKERRQ(ierr);
102*b7ce53b6SStefano Zampini     /* All processes need to compute entire row ownership */
103*b7ce53b6SStefano Zampini     ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
104*b7ce53b6SStefano Zampini     ierr = MatGetOwnershipRanges(new_mat,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
105*b7ce53b6SStefano Zampini     for (i=0;i<nsubdomains;i++) {
106*b7ce53b6SStefano Zampini       for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
107*b7ce53b6SStefano Zampini         row_ownership[j]=i;
108*b7ce53b6SStefano Zampini       }
109*b7ce53b6SStefano Zampini     }
110*b7ce53b6SStefano Zampini 
111*b7ce53b6SStefano Zampini     /*
112*b7ce53b6SStefano Zampini        my_dnz and my_onz contains exact contribution to preallocation from each local mat
113*b7ce53b6SStefano Zampini        then, they will be summed up properly. This way, preallocation is always sufficient
114*b7ce53b6SStefano Zampini     */
115*b7ce53b6SStefano Zampini     ierr = PetscMalloc1(local_rows,&my_dnz);CHKERRQ(ierr);
116*b7ce53b6SStefano Zampini     ierr = PetscMalloc1(local_rows,&my_onz);CHKERRQ(ierr);
117*b7ce53b6SStefano Zampini     ierr = PetscMemzero(my_dnz,local_rows*sizeof(*my_dnz));CHKERRQ(ierr);
118*b7ce53b6SStefano Zampini     ierr = PetscMemzero(my_onz,local_rows*sizeof(*my_onz));CHKERRQ(ierr);
119*b7ce53b6SStefano Zampini     /* preallocation as a MATAIJ */
120*b7ce53b6SStefano Zampini     if (isdense) { /* special case for dense local matrices */
121*b7ce53b6SStefano Zampini       for (i=0;i<local_rows;i++) {
122*b7ce53b6SStefano Zampini         index_row = global_indices[i];
123*b7ce53b6SStefano Zampini         for (j=i;j<local_rows;j++) {
124*b7ce53b6SStefano Zampini           owner = row_ownership[index_row];
125*b7ce53b6SStefano Zampini           index_col = global_indices[j];
126*b7ce53b6SStefano Zampini           if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
127*b7ce53b6SStefano Zampini             my_dnz[i] += 1.0;
128*b7ce53b6SStefano Zampini           } else { /* offdiag block */
129*b7ce53b6SStefano Zampini             my_onz[i] += 1.0;
130*b7ce53b6SStefano Zampini           }
131*b7ce53b6SStefano Zampini           /* same as before, interchanging rows and cols */
132*b7ce53b6SStefano Zampini           if (i != j) {
133*b7ce53b6SStefano Zampini             owner = row_ownership[index_col];
134*b7ce53b6SStefano Zampini             if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
135*b7ce53b6SStefano Zampini               my_dnz[j] += 1.0;
136*b7ce53b6SStefano Zampini             } else {
137*b7ce53b6SStefano Zampini               my_onz[j] += 1.0;
138*b7ce53b6SStefano Zampini             }
139*b7ce53b6SStefano Zampini           }
140*b7ce53b6SStefano Zampini         }
141*b7ce53b6SStefano Zampini       }
142*b7ce53b6SStefano Zampini     } else {
143*b7ce53b6SStefano Zampini       for (i=0;i<local_rows;i++) {
144*b7ce53b6SStefano Zampini         PetscInt ncols;
145*b7ce53b6SStefano Zampini         const PetscInt *cols;
146*b7ce53b6SStefano Zampini         index_row = global_indices[i];
147*b7ce53b6SStefano Zampini         ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
148*b7ce53b6SStefano Zampini         for (j=0;j<ncols;j++) {
149*b7ce53b6SStefano Zampini           owner = row_ownership[index_row];
150*b7ce53b6SStefano Zampini           index_col = global_indices[cols[j]];
151*b7ce53b6SStefano Zampini           if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
152*b7ce53b6SStefano Zampini             my_dnz[i] += 1.0;
153*b7ce53b6SStefano Zampini           } else { /* offdiag block */
154*b7ce53b6SStefano Zampini             my_onz[i] += 1.0;
155*b7ce53b6SStefano Zampini           }
156*b7ce53b6SStefano Zampini           /* same as before, interchanging rows and cols */
157*b7ce53b6SStefano Zampini           if (issbaij) {
158*b7ce53b6SStefano Zampini             owner = row_ownership[index_col];
159*b7ce53b6SStefano Zampini             if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
160*b7ce53b6SStefano Zampini               my_dnz[j] += 1.0;
161*b7ce53b6SStefano Zampini             } else {
162*b7ce53b6SStefano Zampini               my_onz[j] += 1.0;
163*b7ce53b6SStefano Zampini             }
164*b7ce53b6SStefano Zampini           }
165*b7ce53b6SStefano Zampini         }
166*b7ce53b6SStefano Zampini         ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
167*b7ce53b6SStefano Zampini       }
168*b7ce53b6SStefano Zampini     }
169*b7ce53b6SStefano Zampini     ierr = VecSet(vec_dnz,0.0);CHKERRQ(ierr);
170*b7ce53b6SStefano Zampini     ierr = VecSet(vec_onz,0.0);CHKERRQ(ierr);
171*b7ce53b6SStefano Zampini     if (local_rows) { /* multilevel guard */
172*b7ce53b6SStefano Zampini       ierr = VecSetValuesLocal(vec_dnz,local_rows,local_indices,my_dnz,ADD_VALUES);CHKERRQ(ierr);
173*b7ce53b6SStefano Zampini       ierr = VecSetValuesLocal(vec_onz,local_rows,local_indices,my_onz,ADD_VALUES);CHKERRQ(ierr);
174*b7ce53b6SStefano Zampini     }
175*b7ce53b6SStefano Zampini     ierr = VecAssemblyBegin(vec_dnz);CHKERRQ(ierr);
176*b7ce53b6SStefano Zampini     ierr = VecAssemblyBegin(vec_onz);CHKERRQ(ierr);
177*b7ce53b6SStefano Zampini     ierr = VecAssemblyEnd(vec_dnz);CHKERRQ(ierr);
178*b7ce53b6SStefano Zampini     ierr = VecAssemblyEnd(vec_onz);CHKERRQ(ierr);
179*b7ce53b6SStefano Zampini     ierr = PetscFree(my_dnz);CHKERRQ(ierr);
180*b7ce53b6SStefano Zampini     ierr = PetscFree(my_onz);CHKERRQ(ierr);
181*b7ce53b6SStefano Zampini     ierr = PetscFree(row_ownership);CHKERRQ(ierr);
182*b7ce53b6SStefano Zampini 
183*b7ce53b6SStefano Zampini     /* set computed preallocation in dnz and onz */
184*b7ce53b6SStefano Zampini     ierr = VecGetArray(vec_dnz,&array);CHKERRQ(ierr);
185*b7ce53b6SStefano Zampini     for (i=0; i<lrows; i++) dnz[i] = (PetscInt)PetscRealPart(array[i]);
186*b7ce53b6SStefano Zampini     ierr = VecRestoreArray(vec_dnz,&array);CHKERRQ(ierr);
187*b7ce53b6SStefano Zampini     ierr = VecGetArray(vec_onz,&array);CHKERRQ(ierr);
188*b7ce53b6SStefano Zampini     for (i=0;i<lrows;i++) onz[i] = (PetscInt)PetscRealPart(array[i]);
189*b7ce53b6SStefano Zampini     ierr = VecRestoreArray(vec_onz,&array);CHKERRQ(ierr);
190*b7ce53b6SStefano Zampini     ierr = VecDestroy(&vec_dnz);CHKERRQ(ierr);
191*b7ce53b6SStefano Zampini     ierr = VecDestroy(&vec_onz);CHKERRQ(ierr);
192*b7ce53b6SStefano Zampini 
193*b7ce53b6SStefano Zampini     /* Resize preallocation if overestimated */
194*b7ce53b6SStefano Zampini     for (i=0;i<lrows;i++) {
195*b7ce53b6SStefano Zampini       dnz[i] = PetscMin(dnz[i],lcols);
196*b7ce53b6SStefano Zampini       onz[i] = PetscMin(onz[i],cols-lcols);
197*b7ce53b6SStefano Zampini     }
198*b7ce53b6SStefano Zampini     /* set preallocation */
199*b7ce53b6SStefano Zampini     ierr = MatMPIAIJSetPreallocation(new_mat,0,dnz,0,onz);CHKERRQ(ierr);
200*b7ce53b6SStefano Zampini     for (i=0;i<lrows/bs;i++) {
201*b7ce53b6SStefano Zampini       dnz[i] = dnz[i*bs]/bs;
202*b7ce53b6SStefano Zampini       onz[i] = onz[i*bs]/bs;
203*b7ce53b6SStefano Zampini     }
204*b7ce53b6SStefano Zampini     ierr = MatMPIBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr);
205*b7ce53b6SStefano Zampini     for (i=0;i<lrows/bs;i++) {
206*b7ce53b6SStefano Zampini       dnz[i] = dnz[i]-i;
207*b7ce53b6SStefano Zampini     }
208*b7ce53b6SStefano Zampini     ierr = MatMPISBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr);
209*b7ce53b6SStefano Zampini     ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
210*b7ce53b6SStefano Zampini     *M = new_mat;
211*b7ce53b6SStefano Zampini   } else {
212*b7ce53b6SStefano Zampini     PetscInt mbs,mrows,mcols;
213*b7ce53b6SStefano Zampini     /* some checks */
214*b7ce53b6SStefano Zampini     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
215*b7ce53b6SStefano Zampini     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
216*b7ce53b6SStefano Zampini     if (mrows != rows) {
217*b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
218*b7ce53b6SStefano Zampini     }
219*b7ce53b6SStefano Zampini     if (mrows != rows) {
220*b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
221*b7ce53b6SStefano Zampini     }
222*b7ce53b6SStefano Zampini     if (mbs != bs) {
223*b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
224*b7ce53b6SStefano Zampini     }
225*b7ce53b6SStefano Zampini     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
226*b7ce53b6SStefano Zampini   }
227*b7ce53b6SStefano Zampini   /* set local to global mappings */
228*b7ce53b6SStefano Zampini   /* ierr = MatSetLocalToGlobalMapping(*M,rmapping,cmapping);CHKERRQ(ierr); */
229*b7ce53b6SStefano Zampini   /* Set values */
230*b7ce53b6SStefano Zampini   if (isdense) { /* special case for dense local matrices */
231*b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
232*b7ce53b6SStefano Zampini     ierr = MatDenseGetArray(matis->A,&array);CHKERRQ(ierr);
233*b7ce53b6SStefano Zampini     ierr = MatSetValues(*M,local_rows,global_indices,local_cols,global_indices,array,ADD_VALUES);CHKERRQ(ierr);
234*b7ce53b6SStefano Zampini     ierr = MatDenseRestoreArray(matis->A,&array);CHKERRQ(ierr);
235*b7ce53b6SStefano Zampini     ierr = PetscFree(local_indices);CHKERRQ(ierr);
236*b7ce53b6SStefano Zampini     ierr = PetscFree(global_indices);CHKERRQ(ierr);
237*b7ce53b6SStefano Zampini   } else { /* very basic values insertion for all other matrix types */
238*b7ce53b6SStefano Zampini     ierr = PetscFree(local_indices);CHKERRQ(ierr);
239*b7ce53b6SStefano Zampini     for (i=0;i<local_rows;i++) {
240*b7ce53b6SStefano Zampini       ierr = MatGetRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr);
241*b7ce53b6SStefano Zampini       /* ierr = MatSetValuesLocal(*M,1,&i,j,local_indices,array,ADD_VALUES);CHKERRQ(ierr); */
242*b7ce53b6SStefano Zampini       ierr = ISLocalToGlobalMappingApply(matis->mapping,j,local_indices,global_indices);CHKERRQ(ierr);
243*b7ce53b6SStefano Zampini       ierr = ISLocalToGlobalMappingApply(matis->mapping,1,&i,&index_row);CHKERRQ(ierr);
244*b7ce53b6SStefano Zampini       ierr = MatSetValues(*M,1,&index_row,j,global_indices,array,ADD_VALUES);CHKERRQ(ierr);
245*b7ce53b6SStefano Zampini       ierr = MatRestoreRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr);
246*b7ce53b6SStefano Zampini     }
247*b7ce53b6SStefano Zampini     ierr = PetscFree(global_indices);CHKERRQ(ierr);
248*b7ce53b6SStefano Zampini   }
249*b7ce53b6SStefano Zampini   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
250*b7ce53b6SStefano Zampini   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
251*b7ce53b6SStefano Zampini   if (isdense) {
252*b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
253*b7ce53b6SStefano Zampini   }
254*b7ce53b6SStefano Zampini   if (issbaij) {
255*b7ce53b6SStefano Zampini     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
256*b7ce53b6SStefano Zampini   }
257*b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
258*b7ce53b6SStefano Zampini }
259*b7ce53b6SStefano Zampini 
260*b7ce53b6SStefano Zampini #undef __FUNCT__
261*b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ"
262*b7ce53b6SStefano Zampini /*@
263*b7ce53b6SStefano Zampini     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
264*b7ce53b6SStefano Zampini 
265*b7ce53b6SStefano Zampini   Input Parameter:
266*b7ce53b6SStefano Zampini .  mat - the matrix (should be of type MATIS)
267*b7ce53b6SStefano Zampini .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
268*b7ce53b6SStefano Zampini 
269*b7ce53b6SStefano Zampini   Output Parameter:
270*b7ce53b6SStefano Zampini .  newmat - the matrix in AIJ format
271*b7ce53b6SStefano Zampini 
272*b7ce53b6SStefano Zampini   Level: developer
273*b7ce53b6SStefano Zampini 
274*b7ce53b6SStefano Zampini   Notes:
275*b7ce53b6SStefano Zampini 
276*b7ce53b6SStefano Zampini .seealso: MATIS
277*b7ce53b6SStefano Zampini @*/
278*b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
279*b7ce53b6SStefano Zampini {
280*b7ce53b6SStefano Zampini   PetscErrorCode ierr;
281*b7ce53b6SStefano Zampini 
282*b7ce53b6SStefano Zampini   PetscFunctionBegin;
283*b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
284*b7ce53b6SStefano Zampini   PetscValidLogicalCollectiveEnum(mat,reuse,2);
285*b7ce53b6SStefano Zampini   PetscValidPointer(newmat,3);
286*b7ce53b6SStefano Zampini   if (reuse != MAT_INITIAL_MATRIX) {
287*b7ce53b6SStefano Zampini     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
288*b7ce53b6SStefano Zampini     PetscCheckSameComm(mat,1,*newmat,3);
289*b7ce53b6SStefano Zampini   }
290*b7ce53b6SStefano Zampini   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
291*b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
292*b7ce53b6SStefano Zampini }
293*b7ce53b6SStefano Zampini 
294*b7ce53b6SStefano Zampini #undef __FUNCT__
295ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS"
296ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
297ad6194a2SStefano Zampini {
298ad6194a2SStefano Zampini   PetscErrorCode ierr;
299ad6194a2SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
300ad6194a2SStefano Zampini   PetscInt       bs,m,n,M,N;
301ad6194a2SStefano Zampini   Mat            B,localmat;
302ad6194a2SStefano Zampini 
303ad6194a2SStefano Zampini   PetscFunctionBegin;
304ad6194a2SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
305ad6194a2SStefano Zampini   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
306ad6194a2SStefano Zampini   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
307ad6194a2SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,matis->mapping,&B);CHKERRQ(ierr);
308ad6194a2SStefano Zampini   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
309ad6194a2SStefano Zampini   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
310ad6194a2SStefano Zampini   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
311ad6194a2SStefano Zampini   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
312ad6194a2SStefano Zampini   *newmat = B;
313ad6194a2SStefano Zampini   PetscFunctionReturn(0);
314ad6194a2SStefano Zampini }
315ad6194a2SStefano Zampini 
316ad6194a2SStefano Zampini #undef __FUNCT__
31769796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS"
31869796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
31969796d55SStefano Zampini {
32069796d55SStefano Zampini   PetscErrorCode ierr;
32169796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
32269796d55SStefano Zampini   PetscBool      local_sym;
32369796d55SStefano Zampini 
32469796d55SStefano Zampini   PetscFunctionBegin;
32569796d55SStefano Zampini   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
32669796d55SStefano Zampini   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
32769796d55SStefano Zampini   PetscFunctionReturn(0);
32869796d55SStefano Zampini }
32969796d55SStefano Zampini 
33069796d55SStefano Zampini #undef __FUNCT__
33169796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS"
33269796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
33369796d55SStefano Zampini {
33469796d55SStefano Zampini   PetscErrorCode ierr;
33569796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
33669796d55SStefano Zampini   PetscBool      local_sym;
33769796d55SStefano Zampini 
33869796d55SStefano Zampini   PetscFunctionBegin;
33969796d55SStefano Zampini   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
34069796d55SStefano Zampini   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
34169796d55SStefano Zampini   PetscFunctionReturn(0);
34269796d55SStefano Zampini }
34369796d55SStefano Zampini 
34469796d55SStefano Zampini #undef __FUNCT__
345b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
346dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A)
347b4319ba4SBarry Smith {
348dfbe8321SBarry Smith   PetscErrorCode ierr;
349b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
350b4319ba4SBarry Smith 
351b4319ba4SBarry Smith   PetscFunctionBegin;
3526bf464f9SBarry Smith   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
3536bf464f9SBarry Smith   ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr);
3546bf464f9SBarry Smith   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
3556bf464f9SBarry Smith   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
3566bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr);
357bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
358dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
359bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
360*b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
361*b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
362b4319ba4SBarry Smith   PetscFunctionReturn(0);
363b4319ba4SBarry Smith }
364b4319ba4SBarry Smith 
365b4319ba4SBarry Smith #undef __FUNCT__
366b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
367dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
368b4319ba4SBarry Smith {
369dfbe8321SBarry Smith   PetscErrorCode ierr;
370b4319ba4SBarry Smith   Mat_IS         *is  = (Mat_IS*)A->data;
371b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
372b4319ba4SBarry Smith 
373b4319ba4SBarry Smith   PetscFunctionBegin;
374b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
375ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
376ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
377b4319ba4SBarry Smith 
378b4319ba4SBarry Smith   /* multiply the local matrix */
379b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
380b4319ba4SBarry Smith 
381b4319ba4SBarry Smith   /* scatter product back into global memory */
3822dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
383ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
384ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
385b4319ba4SBarry Smith   PetscFunctionReturn(0);
386b4319ba4SBarry Smith }
387b4319ba4SBarry Smith 
388b4319ba4SBarry Smith #undef __FUNCT__
3892e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
390*b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
3912e74eeadSLisandro Dalcin {
392650997f4SStefano Zampini   Vec            temp_vec;
3932e74eeadSLisandro Dalcin   PetscErrorCode ierr;
3942e74eeadSLisandro Dalcin 
3952e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
396650997f4SStefano Zampini   if (v3 != v2) {
397650997f4SStefano Zampini     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
398650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
399650997f4SStefano Zampini   } else {
400650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
401650997f4SStefano Zampini     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
402650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
403650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
404650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
405650997f4SStefano Zampini   }
4062e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
4072e74eeadSLisandro Dalcin }
4082e74eeadSLisandro Dalcin 
4092e74eeadSLisandro Dalcin #undef __FUNCT__
4102e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
4112e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y)
4122e74eeadSLisandro Dalcin {
4132e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
4142e74eeadSLisandro Dalcin   PetscErrorCode ierr;
4152e74eeadSLisandro Dalcin 
4162e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  y = A' * x */
4172e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
418ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
419ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4202e74eeadSLisandro Dalcin 
4212e74eeadSLisandro Dalcin   /* multiply the local matrix */
4222e74eeadSLisandro Dalcin   ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr);
4232e74eeadSLisandro Dalcin 
4242e74eeadSLisandro Dalcin   /* scatter product back into global vector */
4252e74eeadSLisandro Dalcin   ierr = VecSet(y,0);CHKERRQ(ierr);
426ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
427ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4282e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
4292e74eeadSLisandro Dalcin }
4302e74eeadSLisandro Dalcin 
4312e74eeadSLisandro Dalcin #undef __FUNCT__
4322e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
4332e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
4342e74eeadSLisandro Dalcin {
435650997f4SStefano Zampini   Vec            temp_vec;
4362e74eeadSLisandro Dalcin   PetscErrorCode ierr;
4372e74eeadSLisandro Dalcin 
4382e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
439650997f4SStefano Zampini   if (v3 != v2) {
440650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
441650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
442650997f4SStefano Zampini   } else {
443650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
444650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
445650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
446650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
447650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
448650997f4SStefano Zampini   }
4492e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
4502e74eeadSLisandro Dalcin }
4512e74eeadSLisandro Dalcin 
4522e74eeadSLisandro Dalcin #undef __FUNCT__
453b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
454dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
455b4319ba4SBarry Smith {
456b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
457dfbe8321SBarry Smith   PetscErrorCode ierr;
458b4319ba4SBarry Smith   PetscViewer    sviewer;
459b4319ba4SBarry Smith 
460b4319ba4SBarry Smith   PetscFunctionBegin;
461dd2fa690SBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
462b4319ba4SBarry Smith   ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
463b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
464b4319ba4SBarry Smith   ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
465b4319ba4SBarry Smith   PetscFunctionReturn(0);
466b4319ba4SBarry Smith }
467b4319ba4SBarry Smith 
468b4319ba4SBarry Smith #undef __FUNCT__
469b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
470784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
471b4319ba4SBarry Smith {
472dfbe8321SBarry Smith   PetscErrorCode ierr;
4734e4c7dbeSStefano Zampini   PetscInt       n,bs;
474b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
475b4319ba4SBarry Smith   IS             from,to;
476b4319ba4SBarry Smith   Vec            global;
477b4319ba4SBarry Smith 
478b4319ba4SBarry Smith   PetscFunctionBegin;
479784ac674SJed Brown   PetscCheckSameComm(A,1,rmapping,2);
480ce94432eSBarry Smith   if (rmapping != cmapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical");
48170cf5478SStefano Zampini   if (is->mapping) { /* Currenly destroys the objects that will be created by this routine. Is there anything else that should be checked? */
48270cf5478SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
48370cf5478SStefano Zampini     ierr = VecDestroy(&is->x);CHKERRQ(ierr);
48470cf5478SStefano Zampini     ierr = VecDestroy(&is->y);CHKERRQ(ierr);
48570cf5478SStefano Zampini     ierr = VecScatterDestroy(&is->ctx);CHKERRQ(ierr);
4861c47cb0fSStefano Zampini     ierr = MatDestroy(&is->A);CHKERRQ(ierr);
48770cf5478SStefano Zampini   }
488784ac674SJed Brown   ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr);
4896bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
490784ac674SJed Brown   is->mapping = rmapping;
491fa7f1dd8SStefano Zampini /*
492fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
493fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
494fa7f1dd8SStefano Zampini */
495b4319ba4SBarry Smith 
496b4319ba4SBarry Smith   /* Create the local matrix A */
497784ac674SJed Brown   ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr);
4984e4c7dbeSStefano Zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
499f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
500f69a0ea3SMatthew Knepley   ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr);
5014e4c7dbeSStefano Zampini   ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr);
502ff130e51SJed Brown   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
503ff130e51SJed Brown   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
504b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
505b4319ba4SBarry Smith 
506b4319ba4SBarry Smith   /* Create the local work vectors */
5074e4c7dbeSStefano Zampini   ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr);
5084e4c7dbeSStefano Zampini   ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr);
5094e4c7dbeSStefano Zampini   ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr);
510ff130e51SJed Brown   ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr);
511ff130e51SJed Brown   ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr);
5124e4c7dbeSStefano Zampini   ierr = VecSetFromOptions(is->x);CHKERRQ(ierr);
513b4319ba4SBarry Smith   ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr);
514b4319ba4SBarry Smith 
515b4319ba4SBarry Smith   /* setup the global to local scatter */
516b4319ba4SBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr);
517784ac674SJed Brown   ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
5180298fd71SBarry Smith   ierr = MatGetVecs(A,&global,NULL);CHKERRQ(ierr);
519b4319ba4SBarry Smith   ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr);
5206bf464f9SBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
5216bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
5226bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
523b4319ba4SBarry Smith   PetscFunctionReturn(0);
524b4319ba4SBarry Smith }
525b4319ba4SBarry Smith 
5262e74eeadSLisandro Dalcin #undef __FUNCT__
5272e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
5282e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
5292e74eeadSLisandro Dalcin {
5302e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
5312e74eeadSLisandro Dalcin   PetscInt       rows_l[2048],cols_l[2048];
5322e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5332e74eeadSLisandro Dalcin 
5342e74eeadSLisandro Dalcin   PetscFunctionBegin;
5352e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
536f23aa3ddSBarry Smith   if (m > 2048 || n > 2048) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n);
5372e74eeadSLisandro Dalcin #endif
5382e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr);
5392e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr);
5402e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
5412e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5422e74eeadSLisandro Dalcin }
5432e74eeadSLisandro Dalcin 
5442e74eeadSLisandro Dalcin #undef ISG2LMapSetUp
5452e74eeadSLisandro Dalcin #undef ISG2LMapApply
546b4319ba4SBarry Smith 
547b4319ba4SBarry Smith #undef __FUNCT__
548b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
54913f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
550b4319ba4SBarry Smith {
551dfbe8321SBarry Smith   PetscErrorCode ierr;
552b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
553b4319ba4SBarry Smith 
554b4319ba4SBarry Smith   PetscFunctionBegin;
555b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
556b4319ba4SBarry Smith   PetscFunctionReturn(0);
557b4319ba4SBarry Smith }
558b4319ba4SBarry Smith 
559b4319ba4SBarry Smith #undef __FUNCT__
560f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
561f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
562f0006bf2SLisandro Dalcin {
563f0006bf2SLisandro Dalcin   PetscErrorCode ierr;
564f0006bf2SLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
565f0006bf2SLisandro Dalcin 
566f0006bf2SLisandro Dalcin   PetscFunctionBegin;
567f0006bf2SLisandro Dalcin   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
568f0006bf2SLisandro Dalcin   PetscFunctionReturn(0);
569f0006bf2SLisandro Dalcin }
570f0006bf2SLisandro Dalcin 
571f0006bf2SLisandro Dalcin #undef __FUNCT__
5722e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS"
5732b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
5742e74eeadSLisandro Dalcin {
5752e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
5760298fd71SBarry Smith   PetscInt       n_l = 0, *rows_l = NULL;
5772e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5782e74eeadSLisandro Dalcin 
5792e74eeadSLisandro Dalcin   PetscFunctionBegin;
580ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
5812e74eeadSLisandro Dalcin   if (n) {
582785e854fSJed Brown     ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr);
5832e74eeadSLisandro Dalcin     ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr);
5842e74eeadSLisandro Dalcin   }
5852b40b63fSBarry Smith   ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr);
586c31cb41cSBarry Smith   ierr = PetscFree(rows_l);CHKERRQ(ierr);
5872e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5882e74eeadSLisandro Dalcin }
5892e74eeadSLisandro Dalcin 
5902e74eeadSLisandro Dalcin #undef __FUNCT__
591b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS"
5922b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
593b4319ba4SBarry Smith {
594b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
595dfbe8321SBarry Smith   PetscErrorCode ierr;
596f4df32b1SMatthew Knepley   PetscInt       i;
597b4319ba4SBarry Smith   PetscScalar    *array;
598b4319ba4SBarry Smith 
599b4319ba4SBarry Smith   PetscFunctionBegin;
600ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
601b4319ba4SBarry Smith   {
602b4319ba4SBarry Smith     /*
603b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
604b4319ba4SBarry Smith        work properly in the interface nodes.
605b4319ba4SBarry Smith     */
606b4319ba4SBarry Smith     Vec         counter;
607b4319ba4SBarry Smith     PetscScalar one=1.0, zero=0.0;
6080298fd71SBarry Smith     ierr = MatGetVecs(A,&counter,NULL);CHKERRQ(ierr);
6092dcb1b2aSMatthew Knepley     ierr = VecSet(counter,zero);CHKERRQ(ierr);
6102dcb1b2aSMatthew Knepley     ierr = VecSet(is->x,one);CHKERRQ(ierr);
611ca9f406cSSatish Balay     ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
612ca9f406cSSatish Balay     ierr = VecScatterEnd  (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
613ca9f406cSSatish Balay     ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
614ca9f406cSSatish Balay     ierr = VecScatterEnd  (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6156bf464f9SBarry Smith     ierr = VecDestroy(&counter);CHKERRQ(ierr);
616b4319ba4SBarry Smith   }
617958c9bccSBarry Smith   if (!n) {
618b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
619b4319ba4SBarry Smith   } else {
620b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
6212205254eSKarl Rupp 
622b4319ba4SBarry Smith     ierr = VecGetArray(is->x,&array);CHKERRQ(ierr);
6232b40b63fSBarry Smith     ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
624b4319ba4SBarry Smith     for (i=0; i<n; i++) {
625f4df32b1SMatthew Knepley       ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
626b4319ba4SBarry Smith     }
627b4319ba4SBarry Smith     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
628b4319ba4SBarry Smith     ierr = MatAssemblyEnd  (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
629b4319ba4SBarry Smith     ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr);
630b4319ba4SBarry Smith   }
631b4319ba4SBarry Smith   PetscFunctionReturn(0);
632b4319ba4SBarry Smith }
633b4319ba4SBarry Smith 
634b4319ba4SBarry Smith #undef __FUNCT__
635b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
636dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
637b4319ba4SBarry Smith {
638b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
639dfbe8321SBarry Smith   PetscErrorCode ierr;
640b4319ba4SBarry Smith 
641b4319ba4SBarry Smith   PetscFunctionBegin;
642b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
643b4319ba4SBarry Smith   PetscFunctionReturn(0);
644b4319ba4SBarry Smith }
645b4319ba4SBarry Smith 
646b4319ba4SBarry Smith #undef __FUNCT__
647b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
648dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
649b4319ba4SBarry Smith {
650b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
651dfbe8321SBarry Smith   PetscErrorCode ierr;
652b4319ba4SBarry Smith 
653b4319ba4SBarry Smith   PetscFunctionBegin;
654b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
655b4319ba4SBarry Smith   PetscFunctionReturn(0);
656b4319ba4SBarry Smith }
657b4319ba4SBarry Smith 
658b4319ba4SBarry Smith #undef __FUNCT__
659b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
6607087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
661b4319ba4SBarry Smith {
662b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)mat->data;
663b4319ba4SBarry Smith 
664b4319ba4SBarry Smith   PetscFunctionBegin;
665b4319ba4SBarry Smith   *local = is->A;
666b4319ba4SBarry Smith   PetscFunctionReturn(0);
667b4319ba4SBarry Smith }
668b4319ba4SBarry Smith 
669b4319ba4SBarry Smith #undef __FUNCT__
670b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
671b4319ba4SBarry Smith /*@
672b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
673b4319ba4SBarry Smith 
674b4319ba4SBarry Smith   Input Parameter:
675b4319ba4SBarry Smith .  mat - the matrix
676b4319ba4SBarry Smith 
677b4319ba4SBarry Smith   Output Parameter:
678b4319ba4SBarry Smith .  local - the local matrix usually MATSEQAIJ
679b4319ba4SBarry Smith 
680b4319ba4SBarry Smith   Level: advanced
681b4319ba4SBarry Smith 
682b4319ba4SBarry Smith   Notes:
683b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
684b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
685b4319ba4SBarry Smith   of the MatSetValues() operation.
686b4319ba4SBarry Smith 
687b4319ba4SBarry Smith .seealso: MATIS
688b4319ba4SBarry Smith @*/
6897087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
690b4319ba4SBarry Smith {
6914ac538c5SBarry Smith   PetscErrorCode ierr;
692b4319ba4SBarry Smith 
693b4319ba4SBarry Smith   PetscFunctionBegin;
6940700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
695b4319ba4SBarry Smith   PetscValidPointer(local,2);
6964ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
697b4319ba4SBarry Smith   PetscFunctionReturn(0);
698b4319ba4SBarry Smith }
699b4319ba4SBarry Smith 
7003b03a366Sstefano_zampini #undef __FUNCT__
7013b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
7023b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
7033b03a366Sstefano_zampini {
7043b03a366Sstefano_zampini   Mat_IS         *is = (Mat_IS*)mat->data;
7053b03a366Sstefano_zampini   PetscInt       nrows,ncols,orows,ocols;
7063b03a366Sstefano_zampini   PetscErrorCode ierr;
7073b03a366Sstefano_zampini 
7083b03a366Sstefano_zampini   PetscFunctionBegin;
7094e4c7dbeSStefano Zampini   if (is->A) {
7103b03a366Sstefano_zampini     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
7113b03a366Sstefano_zampini     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
712f23aa3ddSBarry Smith     if (orows != nrows || ocols != ncols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local MATIS matrix should be of size %dx%d (you passed a %dx%d matrix)\n",orows,ocols,nrows,ncols);
7134e4c7dbeSStefano Zampini   }
7143b03a366Sstefano_zampini   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
7153b03a366Sstefano_zampini   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
7163b03a366Sstefano_zampini   is->A = local;
7173b03a366Sstefano_zampini   PetscFunctionReturn(0);
7183b03a366Sstefano_zampini }
7193b03a366Sstefano_zampini 
7203b03a366Sstefano_zampini #undef __FUNCT__
7213b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
7223b03a366Sstefano_zampini /*@
7233b03a366Sstefano_zampini     MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix.
7243b03a366Sstefano_zampini 
7253b03a366Sstefano_zampini   Input Parameter:
7263b03a366Sstefano_zampini .  mat - the matrix
7273b03a366Sstefano_zampini .  local - the local matrix usually MATSEQAIJ
7283b03a366Sstefano_zampini 
7293b03a366Sstefano_zampini   Output Parameter:
7303b03a366Sstefano_zampini 
7313b03a366Sstefano_zampini   Level: advanced
7323b03a366Sstefano_zampini 
7333b03a366Sstefano_zampini   Notes:
7343b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
7353b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
7363b03a366Sstefano_zampini 
7373b03a366Sstefano_zampini .seealso: MATIS
7383b03a366Sstefano_zampini @*/
7393b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
7403b03a366Sstefano_zampini {
7413b03a366Sstefano_zampini   PetscErrorCode ierr;
7423b03a366Sstefano_zampini 
7433b03a366Sstefano_zampini   PetscFunctionBegin;
7443b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
745*b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
7463b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
7473b03a366Sstefano_zampini   PetscFunctionReturn(0);
7483b03a366Sstefano_zampini }
7493b03a366Sstefano_zampini 
7506726f965SBarry Smith #undef __FUNCT__
7516726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
7526726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A)
7536726f965SBarry Smith {
7546726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
7556726f965SBarry Smith   PetscErrorCode ierr;
7566726f965SBarry Smith 
7576726f965SBarry Smith   PetscFunctionBegin;
7586726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
7596726f965SBarry Smith   PetscFunctionReturn(0);
7606726f965SBarry Smith }
7616726f965SBarry Smith 
7626726f965SBarry Smith #undef __FUNCT__
7632e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
7642e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
7652e74eeadSLisandro Dalcin {
7662e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
7672e74eeadSLisandro Dalcin   PetscErrorCode ierr;
7682e74eeadSLisandro Dalcin 
7692e74eeadSLisandro Dalcin   PetscFunctionBegin;
7702e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
7712e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
7722e74eeadSLisandro Dalcin }
7732e74eeadSLisandro Dalcin 
7742e74eeadSLisandro Dalcin #undef __FUNCT__
7752e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
7762e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
7772e74eeadSLisandro Dalcin {
7782e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
7792e74eeadSLisandro Dalcin   PetscErrorCode ierr;
7802e74eeadSLisandro Dalcin 
7812e74eeadSLisandro Dalcin   PetscFunctionBegin;
7822e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
7832e74eeadSLisandro Dalcin   ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr);
7842e74eeadSLisandro Dalcin 
7852e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
7862e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
787ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
788ca9f406cSSatish Balay   ierr = VecScatterEnd  (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7892e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
7902e74eeadSLisandro Dalcin }
7912e74eeadSLisandro Dalcin 
7922e74eeadSLisandro Dalcin #undef __FUNCT__
7936726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
794ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
7956726f965SBarry Smith {
7966726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
7976726f965SBarry Smith   PetscErrorCode ierr;
7986726f965SBarry Smith 
7996726f965SBarry Smith   PetscFunctionBegin;
8004e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
8016726f965SBarry Smith   PetscFunctionReturn(0);
8026726f965SBarry Smith }
8036726f965SBarry Smith 
804284134d9SBarry Smith #undef __FUNCT__
805284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
806284134d9SBarry Smith /*@
807284134d9SBarry Smith     MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
808284134d9SBarry Smith        process but not across processes.
809284134d9SBarry Smith 
810284134d9SBarry Smith    Input Parameters:
811284134d9SBarry Smith +     comm - MPI communicator that will share the matrix
8124e4c7dbeSStefano Zampini .     bs - local and global block size of the matrix
813284134d9SBarry Smith .     m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products
814284134d9SBarry Smith -     map - mapping that defines the global number for each local number
815284134d9SBarry Smith 
816284134d9SBarry Smith    Output Parameter:
817284134d9SBarry Smith .    A - the resulting matrix
818284134d9SBarry Smith 
8198e6c10adSSatish Balay    Level: advanced
8208e6c10adSSatish Balay 
821284134d9SBarry Smith    Notes: See MATIS for more details
8228cda6cd7SBarry Smith           m and n are NOT related to the size of the map, they are the size of the part of the vector owned
8238cda6cd7SBarry Smith           by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points
8248cda6cd7SBarry Smith           plus the ghost points to global indices.
825284134d9SBarry Smith 
826284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
827284134d9SBarry Smith @*/
8284e4c7dbeSStefano Zampini PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A)
829284134d9SBarry Smith {
830284134d9SBarry Smith   PetscErrorCode ierr;
831284134d9SBarry Smith 
832284134d9SBarry Smith   PetscFunctionBegin;
833284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
834d3ee2243SStefano Zampini   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
835284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
836284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
8373b03a366Sstefano_zampini   ierr = MatSetUp(*A);CHKERRQ(ierr);
838784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr);
839284134d9SBarry Smith   PetscFunctionReturn(0);
840284134d9SBarry Smith }
841284134d9SBarry Smith 
842b4319ba4SBarry Smith /*MC
8436a818285SBarry Smith    MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners, see PCBDDC.
844b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
845b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
846b4319ba4SBarry Smith    product is handled "implicitly".
847b4319ba4SBarry Smith 
848b4319ba4SBarry Smith    Operations Provided:
8496726f965SBarry Smith +  MatMult()
8502e74eeadSLisandro Dalcin .  MatMultAdd()
8512e74eeadSLisandro Dalcin .  MatMultTranspose()
8522e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
8536726f965SBarry Smith .  MatZeroEntries()
8546726f965SBarry Smith .  MatSetOption()
8552e74eeadSLisandro Dalcin .  MatZeroRows()
8566726f965SBarry Smith .  MatZeroRowsLocal()
8572e74eeadSLisandro Dalcin .  MatSetValues()
8586726f965SBarry Smith .  MatSetValuesLocal()
8592e74eeadSLisandro Dalcin .  MatScale()
8602e74eeadSLisandro Dalcin .  MatGetDiagonal()
8616726f965SBarry Smith -  MatSetLocalToGlobalMapping()
862b4319ba4SBarry Smith 
863b4319ba4SBarry Smith    Options Database Keys:
864b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
865b4319ba4SBarry Smith 
866b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
867b4319ba4SBarry Smith 
868b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
869b4319ba4SBarry Smith 
870b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
871b4319ba4SBarry Smith           MatISGetLocalMat()
872b4319ba4SBarry Smith 
873b4319ba4SBarry Smith   Level: advanced
874b4319ba4SBarry Smith 
8756a818285SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), PCBDDC
876b4319ba4SBarry Smith 
877b4319ba4SBarry Smith M*/
878b4319ba4SBarry Smith 
879b4319ba4SBarry Smith #undef __FUNCT__
880b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
8818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
882b4319ba4SBarry Smith {
883dfbe8321SBarry Smith   PetscErrorCode ierr;
884b4319ba4SBarry Smith   Mat_IS         *b;
885b4319ba4SBarry Smith 
886b4319ba4SBarry Smith   PetscFunctionBegin;
887b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
888b4319ba4SBarry Smith   A->data = (void*)b;
889b4319ba4SBarry Smith   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
890b4319ba4SBarry Smith 
891b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
8922e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
8932e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
8942e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
895b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
896b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
8972e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
898b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
899f0006bf2SLisandro Dalcin   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
9002e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
901b4319ba4SBarry Smith   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
902b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
903b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
904b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
9056726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
9062e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
9072e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
9086726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
90969796d55SStefano Zampini   A->ops->ishermitian             = MatIsHermitian_IS;
91069796d55SStefano Zampini   A->ops->issymmetric             = MatIsSymmetric_IS;
911ad6194a2SStefano Zampini   A->ops->duplicate               = MatDuplicate_IS;
912b4319ba4SBarry Smith 
91326283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
91426283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
915b4319ba4SBarry Smith 
916*b7ce53b6SStefano Zampini   /* zeros pointer */
917b4319ba4SBarry Smith   b->A       = 0;
918b4319ba4SBarry Smith   b->ctx     = 0;
919b4319ba4SBarry Smith   b->x       = 0;
920b4319ba4SBarry Smith   b->y       = 0;
921b09366fdSStefano Zampini   b->mapping = 0;
922*b7ce53b6SStefano Zampini 
923*b7ce53b6SStefano Zampini   /* special MATIS functions */
924bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
925bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
926*b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
92717667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
928b4319ba4SBarry Smith   PetscFunctionReturn(0);
929b4319ba4SBarry Smith }
930