Actual source code: mumps.c
1: #define PETSCMAT_DLL
3: /*
4: Provides an interface to the MUMPS sparse solver
5: */
6: #include src/mat/impls/aij/seq/aij.h
7: #include src/mat/impls/aij/mpi/mpiaij.h
8: #include src/mat/impls/sbaij/seq/sbaij.h
9: #include src/mat/impls/sbaij/mpi/mpisbaij.h
12: #if defined(PETSC_USE_COMPLEX)
13: #include "zmumps_c.h"
14: #else
15: #include "dmumps_c.h"
16: #endif
18: #define JOB_INIT -1
19: #define JOB_END -2
20: /* macros s.t. indices match MUMPS documentation */
21: #define ICNTL(I) icntl[(I)-1]
22: #define CNTL(I) cntl[(I)-1]
23: #define INFOG(I) infog[(I)-1]
24: #define INFO(I) info[(I)-1]
25: #define RINFOG(I) rinfog[(I)-1]
26: #define RINFO(I) rinfo[(I)-1]
28: typedef struct {
29: #if defined(PETSC_USE_COMPLEX)
30: ZMUMPS_STRUC_C id;
31: #else
32: DMUMPS_STRUC_C id;
33: #endif
34: MatStructure matstruc;
35: PetscMPIInt myid,size;
36: PetscInt *irn,*jcn,sym,nSolve;
37: PetscScalar *val;
38: MPI_Comm comm_mumps;
39: VecScatter scat_rhs, scat_sol;
40: PetscTruth isAIJ,CleanUpMUMPS;
41: Vec b_seq,x_seq;
42: PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
43: PetscErrorCode (*MatView)(Mat,PetscViewer);
44: PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
45: PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
46: PetscErrorCode (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*);
47: PetscErrorCode (*MatDestroy)(Mat);
48: PetscErrorCode (*specialdestroy)(Mat);
49: PetscErrorCode (*MatPreallocate)(Mat,int,int,int*,int,int*);
50: } Mat_MUMPS;
52: EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
54: PetscErrorCode MatConvert_SBAIJ_SBAIJMUMPS(Mat,MatType,MatReuse,Mat*);
56: /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */
57: /*
58: input:
59: A - matrix in mpiaij or mpisbaij (bs=1) format
60: shift - 0: C style output triple; 1: Fortran style output triple.
61: valOnly - FALSE: spaces are allocated and values are set for the triple
62: TRUE: only the values in v array are updated
63: output:
64: nnz - dim of r, c, and v (number of local nonzero entries of A)
65: r, c, v - row and col index, matrix values (matrix triples)
66: */
67: PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) {
68: PetscInt *ai, *aj, *bi, *bj, rstart,nz, *garray;
70: PetscInt i,j,jj,jB,irow,m=A->rmap.n,*ajj,*bjj,countA,countB,colA_start,jcol;
71: PetscInt *row,*col;
72: PetscScalar *av, *bv,*val;
73: Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
76: if (mumps->isAIJ){
77: Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data;
78: Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(mat->A)->data;
79: Mat_SeqAIJ *bb=(Mat_SeqAIJ*)(mat->B)->data;
80: nz = aa->nz + bb->nz;
81: ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart;
82: garray = mat->garray;
83: av=aa->a; bv=bb->a;
84:
85: } else {
86: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data;
87: Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
88: Mat_SeqBAIJ *bb=(Mat_SeqBAIJ*)(mat->B)->data;
89: if (A->rmap.bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", A->rmap.bs);
90: nz = aa->nz + bb->nz;
91: ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart;
92: garray = mat->garray;
93: av=aa->a; bv=bb->a;
94: }
96: if (!valOnly){
97: PetscMalloc(nz*sizeof(PetscInt) ,&row);
98: PetscMalloc(nz*sizeof(PetscInt),&col);
99: PetscMalloc(nz*sizeof(PetscScalar),&val);
100: *r = row; *c = col; *v = val;
101: } else {
102: row = *r; col = *c; val = *v;
103: }
104: *nnz = nz;
106: jj = 0; irow = rstart;
107: for ( i=0; i<m; i++ ) {
108: ajj = aj + ai[i]; /* ptr to the beginning of this row */
109: countA = ai[i+1] - ai[i];
110: countB = bi[i+1] - bi[i];
111: bjj = bj + bi[i];
113: /* get jB, the starting local col index for the 2nd B-part */
114: colA_start = rstart + ajj[0]; /* the smallest col index for A */
115: j=-1;
116: do {
117: j++;
118: if (j == countB) break;
119: jcol = garray[bjj[j]];
120: } while (jcol < colA_start);
121: jB = j;
122:
123: /* B-part, smaller col index */
124: colA_start = rstart + ajj[0]; /* the smallest col index for A */
125: for (j=0; j<jB; j++){
126: jcol = garray[bjj[j]];
127: if (!valOnly){
128: row[jj] = irow + shift; col[jj] = jcol + shift;
130: }
131: val[jj++] = *bv++;
132: }
133: /* A-part */
134: for (j=0; j<countA; j++){
135: if (!valOnly){
136: row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
137: }
138: val[jj++] = *av++;
139: }
140: /* B-part, larger col index */
141: for (j=jB; j<countB; j++){
142: if (!valOnly){
143: row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
144: }
145: val[jj++] = *bv++;
146: }
147: irow++;
148: }
149:
150: return(0);
151: }
156: PetscErrorCode MatConvert_MUMPS_Base(Mat A,MatType type,MatReuse reuse,Mat *newmat)
157: {
159: Mat B=*newmat;
160: Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
161: void (*f)(void);
164: if (reuse == MAT_INITIAL_MATRIX) {
165: MatDuplicate(A,MAT_COPY_VALUES,&B);
166: }
167: B->ops->duplicate = mumps->MatDuplicate;
168: B->ops->view = mumps->MatView;
169: B->ops->assemblyend = mumps->MatAssemblyEnd;
170: B->ops->lufactorsymbolic = mumps->MatLUFactorSymbolic;
171: B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic;
172: B->ops->destroy = mumps->MatDestroy;
174: /* put back original composed preallocation function */
175: PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);
176: if (f) {
177: PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C","",(PetscVoidFunction)mumps->MatPreallocate);
178: }
180: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_aijmumps_C","",PETSC_NULL);
181: PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_seqaij_C","",PETSC_NULL);
182: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpiaij_aijmumps_C","",PETSC_NULL);
183: PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_mpiaij_C","",PETSC_NULL);
184: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C","",PETSC_NULL);
185: PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C","",PETSC_NULL);
186: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C","",PETSC_NULL);
187: PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C","",PETSC_NULL);
189: PetscObjectChangeTypeName((PetscObject)B,type);
190: *newmat = B;
191: return(0);
192: }
197: PetscErrorCode MatDestroy_MUMPS(Mat A)
198: {
199: Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
201: PetscMPIInt size=lu->size;
202: PetscErrorCode (*specialdestroy)(Mat);
204: if (lu->CleanUpMUMPS) {
205: /* Terminate instance, deallocate memories */
206: if (size > 1){
207: PetscFree(lu->id.sol_loc);
208: VecScatterDestroy(lu->scat_rhs);
209: VecDestroy(lu->b_seq);
210: VecScatterDestroy(lu->scat_sol);
211: VecDestroy(lu->x_seq);
212: PetscFree(lu->val);
213: }
214: lu->id.job=JOB_END;
215: #if defined(PETSC_USE_COMPLEX)
216: zmumps_c(&lu->id);
217: #else
218: dmumps_c(&lu->id);
219: #endif
220: PetscFree(lu->irn);
221: PetscFree(lu->jcn);
222: MPI_Comm_free(&(lu->comm_mumps));
223: }
224: specialdestroy = lu->specialdestroy;
225: (*specialdestroy)(A);
226: (*A->ops->destroy)(A);
227: return(0);
228: }
232: PetscErrorCode MatDestroy_AIJMUMPS(Mat A)
233: {
235: PetscMPIInt size;
238: MPI_Comm_size(A->comm,&size);
239: if (size==1) {
240: MatConvert_MUMPS_Base(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);
241: } else {
242: MatConvert_MUMPS_Base(A,MATMPIAIJ,MAT_REUSE_MATRIX,&A);
243: }
244: return(0);
245: }
249: PetscErrorCode MatDestroy_SBAIJMUMPS(Mat A)
250: {
252: PetscMPIInt size;
255: MPI_Comm_size(A->comm,&size);
256: if (size==1) {
257: MatConvert_MUMPS_Base(A,MATSEQSBAIJ,MAT_REUSE_MATRIX,&A);
258: } else {
259: MatConvert_MUMPS_Base(A,MATMPISBAIJ,MAT_REUSE_MATRIX,&A);
260: }
261: return(0);
262: }
266: PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
267: Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
271: PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");
272: PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",lu->id.sym);
273: PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",lu->id.par);
274: PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",lu->id.ICNTL(1));
275: PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));
276: PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",lu->id.ICNTL(3));
277: PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",lu->id.ICNTL(4));
278: PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",lu->id.ICNTL(5));
279: PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",lu->id.ICNTL(6));
280: PetscViewerASCIIPrintf(viewer," ICNTL(7) (matrix ordering): %d \n",lu->id.ICNTL(7));
281: PetscViewerASCIIPrintf(viewer," ICNTL(8) (scalling strategy): %d \n",lu->id.ICNTL(8));
282: PetscViewerASCIIPrintf(viewer," ICNTL(9) (A/A^T x=b is solved): %d \n",lu->id.ICNTL(9));
283: PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));
284: PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",lu->id.ICNTL(11));
285: if (!lu->myid && lu->id.ICNTL(11)>0) {
286: PetscPrintf(PETSC_COMM_SELF," RINFOG(4) (inf norm of input mat): %g\n",lu->id.RINFOG(4));
287: PetscPrintf(PETSC_COMM_SELF," RINFOG(5) (inf norm of solution): %g\n",lu->id.RINFOG(5));
288: PetscPrintf(PETSC_COMM_SELF," RINFOG(6) (inf norm of residual): %g\n",lu->id.RINFOG(6));
289: PetscPrintf(PETSC_COMM_SELF," RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",lu->id.RINFOG(7),lu->id.RINFOG(8));
290: PetscPrintf(PETSC_COMM_SELF," RINFOG(9) (error estimate): %g \n",lu->id.RINFOG(9));
291: PetscPrintf(PETSC_COMM_SELF," RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));
292:
293: }
294: PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",lu->id.ICNTL(12));
295: PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",lu->id.ICNTL(13));
296: PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));
297: /* ICNTL(15-17) not used */
298: PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",lu->id.ICNTL(18));
299: PetscViewerASCIIPrintf(viewer," ICNTL(19) (Shur complement info): %d \n",lu->id.ICNTL(19));
300: PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",lu->id.ICNTL(20));
301: PetscViewerASCIIPrintf(viewer," ICNTL(21) (solution struct): %d \n",lu->id.ICNTL(21));
303: PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",lu->id.CNTL(1));
304: PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));
305: PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",lu->id.CNTL(3));
306: PetscViewerASCIIPrintf(viewer," CNTL(4) (value of static pivoting): %g \n",lu->id.CNTL(4));
308: /* infomation local to each processor */
309: if (!lu->myid) {PetscPrintf(PETSC_COMM_SELF, " RINFO(1) (local estimated flops for the elimination after analysis): \n");}
310: PetscSynchronizedPrintf(A->comm," [%d] %g \n",lu->myid,lu->id.RINFO(1));
311: PetscSynchronizedFlush(A->comm);
312: if (!lu->myid) {PetscPrintf(PETSC_COMM_SELF, " RINFO(2) (local estimated flops for the assembly after factorization): \n");}
313: PetscSynchronizedPrintf(A->comm," [%d] %g \n",lu->myid,lu->id.RINFO(2));
314: PetscSynchronizedFlush(A->comm);
315: if (!lu->myid) {PetscPrintf(PETSC_COMM_SELF, " RINFO(3) (local estimated flops for the elimination after factorization): \n");}
316: PetscSynchronizedPrintf(A->comm," [%d] %g \n",lu->myid,lu->id.RINFO(3));
317: PetscSynchronizedFlush(A->comm);
318: /*
319: if (!lu->myid) {PetscPrintf(PETSC_COMM_SELF, " INFO(2) (info about error or warning ): \n");}
320: PetscSynchronizedPrintf(A->comm," [%d] %d \n",lu->myid,lu->id.INFO(2));
321: PetscSynchronizedFlush(A->comm);
322: */
324: if (!lu->myid) {PetscPrintf(PETSC_COMM_SELF, " INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");}
325: PetscSynchronizedPrintf(A->comm," [%d] %d \n",lu->myid,lu->id.INFO(15));
326: PetscSynchronizedFlush(A->comm);
328: if (!lu->myid) {PetscPrintf(PETSC_COMM_SELF, " INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");}
329: PetscSynchronizedPrintf(A->comm," [%d] %d \n",lu->myid,lu->id.INFO(16));
330: PetscSynchronizedFlush(A->comm);
332: if (!lu->myid) {PetscPrintf(PETSC_COMM_SELF, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");}
333: PetscSynchronizedPrintf(A->comm," [%d] %d \n",lu->myid,lu->id.INFO(23));
334: PetscSynchronizedFlush(A->comm);
336: if (!lu->myid){ /* information from the host */
337: PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));
338: PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));
339: PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));
341: PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));
342: PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));
343: PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));
344: PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));
345: PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));
346: PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));
347: PetscViewerASCIIPrintf(viewer," INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",lu->id.INFOG(9));
348: PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));
349: PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));
350: PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));
351: PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));
352: PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));
353: PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));
354: PetscViewerASCIIPrintf(viewer," INFOG(16) (estimated size (in MB) of all MUMPS internal data for factorization after analysis: value on the most memory consuming processor): %d \n",lu->id.INFOG(16));
355: PetscViewerASCIIPrintf(viewer," INFOG(17) (estimated size of all MUMPS internal data for factorization after analysis: sum over all processors): %d \n",lu->id.INFOG(17));
356: PetscViewerASCIIPrintf(viewer," INFOG(18) (size of all MUMPS internal data allocated during factorization: value on the most memory consuming processor): %d \n",lu->id.INFOG(18));
357: PetscViewerASCIIPrintf(viewer," INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",lu->id.INFOG(19));
358: PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));
359: PetscViewerASCIIPrintf(viewer," INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",lu->id.INFOG(21));
360: PetscViewerASCIIPrintf(viewer," INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",lu->id.INFOG(22));
361: PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));
362: PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));
363: PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));
364: }
366: return(0);
367: }
371: PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) {
372: PetscErrorCode ierr;
373: PetscTruth iascii;
374: PetscViewerFormat format;
375: Mat_MUMPS *mumps=(Mat_MUMPS*)(A->spptr);
378: (*mumps->MatView)(A,viewer);
380: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
381: if (iascii) {
382: PetscViewerGetFormat(viewer,&format);
383: if (format == PETSC_VIEWER_ASCII_INFO){
384: MatFactorInfo_MUMPS(A,viewer);
385: }
386: }
387: return(0);
388: }
392: PetscErrorCode MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) {
393: Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
394: PetscScalar *array;
395: Vec x_seq;
396: IS is_iden,is_petsc;
397: VecScatter scat_rhs=lu->scat_rhs,scat_sol=lu->scat_sol;
399: PetscInt i;
402: lu->id.nrhs = 1;
403: x_seq = lu->b_seq;
404: if (lu->size > 1){
405: /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
406: VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat_rhs);
407: VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat_rhs);
408: if (!lu->myid) {VecGetArray(x_seq,&array);}
409: } else { /* size == 1 */
410: VecCopy(b,x);
411: VecGetArray(x,&array);
412: }
413: if (!lu->myid) { /* define rhs on the host */
414: #if defined(PETSC_USE_COMPLEX)
415: lu->id.rhs = (mumps_double_complex*)array;
416: #else
417: lu->id.rhs = array;
418: #endif
419: }
420: if (lu->size == 1){
421: VecRestoreArray(x,&array);
422: } else if (!lu->myid){
423: VecRestoreArray(x_seq,&array);
424: }
426: if (lu->size > 1){
427: /* distributed solution */
428: lu->id.ICNTL(21) = 1;
429: if (!lu->nSolve){
430: /* Create x_seq=sol_loc for repeated use */
431: PetscInt lsol_loc;
432: PetscScalar *sol_loc;
433: lsol_loc = lu->id.INFO(23); /* length of sol_loc */
434: PetscMalloc((1+lsol_loc)*(sizeof(PetscScalar)+sizeof(PetscInt)),&sol_loc);
435: lu->id.isol_loc = (PetscInt *)(sol_loc + lsol_loc);
436: lu->id.lsol_loc = lsol_loc;
437: lu->id.sol_loc = (F_DOUBLE *)sol_loc;
438: VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);
439: }
440: }
442: /* solve phase */
443: /*-------------*/
444: lu->id.job = 3;
445: #if defined(PETSC_USE_COMPLEX)
446: zmumps_c(&lu->id);
447: #else
448: dmumps_c(&lu->id);
449: #endif
450: if (lu->id.INFOG(1) < 0) {
451: SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
452: }
454: if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */
455: if (!lu->nSolve){ /* create scatter scat_sol */
456: ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden); /* from */
457: for (i=0; i<lu->id.lsol_loc; i++){
458: lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */
459: }
460: ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc); /* to */
461: VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);
462: ISDestroy(is_iden);
463: ISDestroy(is_petsc);
464: }
465: VecScatterBegin(lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD,lu->scat_sol);
466: VecScatterEnd(lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD,lu->scat_sol);
467: }
468: lu->nSolve++;
469: return(0);
470: }
472: /*
473: input:
474: F: numeric factor
475: output:
476: nneg: total number of negative pivots
477: nzero: 0
478: npos: (global dimension of F) - nneg
479: */
483: PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
484: {
485: Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr;
487: PetscMPIInt size;
490: MPI_Comm_size(F->comm,&size);
491: /* MUMPS 4.3.1 calls ScaLAPACK when ICNTL(13)=0 (default), which does not offer the possibility to compute the inertia of a dense matrix. Set ICNTL(13)=1 to skip ScaLAPACK */
492: if (size > 1 && lu->id.ICNTL(13) != 1){
493: SETERRQ1(PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",lu->id.INFOG(13));
494: }
495: if (nneg){
496: if (!lu->myid){
497: *nneg = lu->id.INFOG(12);
498: }
499: MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);
500: }
501: if (nzero) *nzero = 0;
502: if (npos) *npos = F->rmap.N - (*nneg);
503: return(0);
504: }
508: PetscErrorCode MatFactorNumeric_AIJMUMPS(Mat A,MatFactorInfo *info,Mat *F)
509: {
510: Mat_MUMPS *lu =(Mat_MUMPS*)(*F)->spptr;
511: Mat_MUMPS *lua=(Mat_MUMPS*)(A)->spptr;
513: PetscInt rnz,nnz,nz=0,i,M=A->rmap.N,*ai,*aj,icntl;
514: PetscTruth valOnly,flg;
515: Mat F_diag;
518: if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
519: (*F)->ops->solve = MatSolve_AIJMUMPS;
521: /* Initialize a MUMPS instance */
522: MPI_Comm_rank(A->comm, &lu->myid);
523: MPI_Comm_size(A->comm,&lu->size);
524: lua->myid = lu->myid; lua->size = lu->size;
525: lu->id.job = JOB_INIT;
526: MPI_Comm_dup(A->comm,&(lu->comm_mumps));
527: MPICCommToFortranComm(lu->comm_mumps,&(lu->id.comm_fortran));
529: /* Set mumps options */
530: PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");
531: lu->id.par=1; /* host participates factorizaton and solve */
532: lu->id.sym=lu->sym;
533: if (lu->sym == 2){
534: PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);
535: if (flg && icntl == 1) lu->id.sym=icntl; /* matrix is spd */
536: }
537: #if defined(PETSC_USE_COMPLEX)
538: zmumps_c(&lu->id);
539: #else
540: dmumps_c(&lu->id);
541: #endif
542:
543: if (lu->size == 1){
544: lu->id.ICNTL(18) = 0; /* centralized assembled matrix input */
545: } else {
546: lu->id.ICNTL(18) = 3; /* distributed assembled matrix input */
547: }
549: icntl=-1;
550: PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);
551: if ((flg && icntl > 0) || PetscLogPrintInfo) {
552: lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
553: } else { /* no output */
554: lu->id.ICNTL(1) = 0; /* error message, default= 6 */
555: lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */
556: lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */
557: lu->id.ICNTL(4) = 0; /* level of printing, 0,1,2,3,4, default=2 */
558: }
559: PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): matrix prescaling (0 to 7)","None",lu->id.ICNTL(6),&lu->id.ICNTL(6),PETSC_NULL);
560: icntl=-1;
561: PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);
562: if (flg) {
563: if (icntl== 1){
564: SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
565: } else {
566: lu->id.ICNTL(7) = icntl;
567: }
568: }
569: PetscOptionsInt("-mat_mumps_icntl_9","ICNTL(9): A or A^T x=b to be solved. 1: A; otherwise: A^T","None",lu->id.ICNTL(9),&lu->id.ICNTL(9),PETSC_NULL);
570: PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",lu->id.ICNTL(10),&lu->id.ICNTL(10),PETSC_NULL);
571: PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): error analysis, a positive value returns statistics (by -ksp_view)","None",lu->id.ICNTL(11),&lu->id.ICNTL(11),PETSC_NULL);
572: PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);
573: PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);
574: PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage of estimated workspace increase","None",lu->id.ICNTL(14),&lu->id.ICNTL(14),PETSC_NULL);
575: PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);
577: PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",lu->id.CNTL(1),&lu->id.CNTL(1),PETSC_NULL);
578: PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",lu->id.CNTL(2),&lu->id.CNTL(2),PETSC_NULL);
579: PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",lu->id.CNTL(3),&lu->id.CNTL(3),PETSC_NULL);
580: PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",lu->id.CNTL(4),&lu->id.CNTL(4),PETSC_NULL);
581: PetscOptionsEnd();
582: }
584: /* define matrix A */
585: switch (lu->id.ICNTL(18)){
586: case 0: /* centralized assembled matrix input (size=1) */
587: if (!lu->myid) {
588: if (lua->isAIJ){
589: Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data;
590: nz = aa->nz;
591: ai = aa->i; aj = aa->j; lu->val = aa->a;
592: } else {
593: Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
594: nz = aa->nz;
595: ai = aa->i; aj = aa->j; lu->val = aa->a;
596: }
597: if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
598: PetscMalloc(nz*sizeof(PetscInt),&lu->irn);
599: PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);
600: nz = 0;
601: for (i=0; i<M; i++){
602: rnz = ai[i+1] - ai[i];
603: while (rnz--) { /* Fortran row/col index! */
604: lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
605: }
606: }
607: }
608: }
609: break;
610: case 3: /* distributed assembled matrix input (size>1) */
611: if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
612: valOnly = PETSC_FALSE;
613: } else {
614: valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
615: }
616: MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);
617: break;
618: default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
619: }
621: /* analysis phase */
622: /*----------------*/
623: if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
624: lu->id.job = 1;
626: lu->id.n = M;
627: switch (lu->id.ICNTL(18)){
628: case 0: /* centralized assembled matrix input */
629: if (!lu->myid) {
630: lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
631: if (lu->id.ICNTL(6)>1){
632: #if defined(PETSC_USE_COMPLEX)
633: lu->id.a = (mumps_double_complex*)lu->val;
634: #else
635: lu->id.a = lu->val;
636: #endif
637: }
638: }
639: break;
640: case 3: /* distributed assembled matrix input (size>1) */
641: lu->id.nz_loc = nnz;
642: lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
643: if (lu->id.ICNTL(6)>1) {
644: #if defined(PETSC_USE_COMPLEX)
645: lu->id.a_loc = (mumps_double_complex*)lu->val;
646: #else
647: lu->id.a_loc = lu->val;
648: #endif
649: }
650: /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
651: IS is_iden;
652: Vec b;
653: if (!lu->myid){
654: VecCreateSeq(PETSC_COMM_SELF,A->cmap.N,&lu->b_seq);
655: ISCreateStride(PETSC_COMM_SELF,A->cmap.N,0,1,&is_iden);
656: } else {
657: VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);
658: ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
659: }
660: VecCreate(A->comm,&b);
661: VecSetSizes(b,A->rmap.n,PETSC_DECIDE);
662: VecSetFromOptions(b);
664: VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);
665: ISDestroy(is_iden);
666: VecDestroy(b);
667: break;
668: }
669: #if defined(PETSC_USE_COMPLEX)
670: zmumps_c(&lu->id);
671: #else
672: dmumps_c(&lu->id);
673: #endif
674: if (lu->id.INFOG(1) < 0) {
675: SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
676: }
677: }
679: /* numerical factorization phase */
680: /*-------------------------------*/
681: lu->id.job = 2;
682: if(!lu->id.ICNTL(18)) {
683: if (!lu->myid) {
684: #if defined(PETSC_USE_COMPLEX)
685: lu->id.a = (mumps_double_complex*)lu->val;
686: #else
687: lu->id.a = lu->val;
688: #endif
689: }
690: } else {
691: #if defined(PETSC_USE_COMPLEX)
692: lu->id.a_loc = (mumps_double_complex*)lu->val;
693: #else
694: lu->id.a_loc = lu->val;
695: #endif
696: }
697: #if defined(PETSC_USE_COMPLEX)
698: zmumps_c(&lu->id);
699: #else
700: dmumps_c(&lu->id);
701: #endif
702: if (lu->id.INFOG(1) < 0) {
703: if (lu->id.INFO(1) == -13) {
704: SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2));
705: } else {
706: SETERRQ2(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",lu->id.INFO(1),lu->id.INFO(2));
707: }
708: }
710: if (!lu->myid && lu->id.ICNTL(16) > 0){
711: SETERRQ1(PETSC_ERR_LIB," lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
712: }
714: if (lu->size > 1){
715: if ((*F)->factor == FACTOR_LU){
716: F_diag = ((Mat_MPIAIJ *)(*F)->data)->A;
717: } else {
718: F_diag = ((Mat_MPISBAIJ *)(*F)->data)->A;
719: }
720: F_diag->assembled = PETSC_TRUE;
721: if (lu->nSolve){
722: VecScatterDestroy(lu->scat_sol);
723: PetscFree(lu->id.sol_loc);
724: VecDestroy(lu->x_seq);
725: }
726: }
727: (*F)->assembled = PETSC_TRUE;
728: lu->matstruc = SAME_NONZERO_PATTERN;
729: lu->CleanUpMUMPS = PETSC_TRUE;
730: lu->nSolve = 0;
731: return(0);
732: }
734: /* Note the Petsc r and c permutations are ignored */
737: PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
738: Mat B;
739: Mat_MUMPS *lu;
743: /* Create the factorization matrix */
744: MatCreate(A->comm,&B);
745: MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);
746: MatSetType(B,A->type_name);
747: MatSeqAIJSetPreallocation(B,0,PETSC_NULL);
748: MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);
750: B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS;
751: B->factor = FACTOR_LU;
752: lu = (Mat_MUMPS*)B->spptr;
753: lu->sym = 0;
754: lu->matstruc = DIFFERENT_NONZERO_PATTERN;
756: *F = B;
757: return(0);
758: }
760: /* Note the Petsc r permutation is ignored */
763: PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) {
764: Mat B;
765: Mat_MUMPS *lu;
769: /* Create the factorization matrix */
770: MatCreate(A->comm,&B);
771: MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);
772: MatSetType(B,A->type_name);
773: MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);
774: MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);
776: B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS;
777: B->ops->getinertia = MatGetInertia_SBAIJMUMPS;
778: B->factor = FACTOR_CHOLESKY;
779: lu = (Mat_MUMPS*)B->spptr;
780: lu->sym = 2;
781: lu->matstruc = DIFFERENT_NONZERO_PATTERN;
783: *F = B;
784: return(0);
785: }
789: PetscErrorCode MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) {
791: Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
794: (*mumps->MatAssemblyEnd)(A,mode);
796: mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
797: mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
798: A->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
799: return(0);
800: }
805: PetscErrorCode MatConvert_AIJ_AIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
806: {
808: PetscMPIInt size;
809: MPI_Comm comm;
810: Mat B=*newmat;
811: Mat_MUMPS *mumps;
814: PetscObjectGetComm((PetscObject)A,&comm);
815: PetscNew(Mat_MUMPS,&mumps);
817: if (reuse == MAT_INITIAL_MATRIX) {
818: MatDuplicate(A,MAT_COPY_VALUES,&B);
819: /* A may have special container that is not duplicated,
820: e.g., A is obtainted from MatMatMult(,&A). Save B->ops instead */
821: mumps->MatDuplicate = B->ops->duplicate;
822: mumps->MatView = B->ops->view;
823: mumps->MatAssemblyEnd = B->ops->assemblyend;
824: mumps->MatLUFactorSymbolic = B->ops->lufactorsymbolic;
825: mumps->MatCholeskyFactorSymbolic = B->ops->choleskyfactorsymbolic;
826: mumps->MatDestroy = B->ops->destroy;
827: } else {
828: mumps->MatDuplicate = A->ops->duplicate;
829: mumps->MatView = A->ops->view;
830: mumps->MatAssemblyEnd = A->ops->assemblyend;
831: mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
832: mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
833: mumps->MatDestroy = A->ops->destroy;
834: }
835: mumps->specialdestroy = MatDestroy_AIJMUMPS;
836: mumps->CleanUpMUMPS = PETSC_FALSE;
837: mumps->isAIJ = PETSC_TRUE;
839: B->spptr = (void*)mumps;
840: B->ops->duplicate = MatDuplicate_MUMPS;
841: B->ops->view = MatView_MUMPS;
842: B->ops->assemblyend = MatAssemblyEnd_AIJMUMPS;
843: B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
844: B->ops->destroy = MatDestroy_MUMPS;
846: MPI_Comm_size(comm,&size);
847: if (size == 1) {
848: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C",
849: "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);
850: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C",
851: "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);
852: } else {
853: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C",
854: "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);
855: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C",
856: "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);
857: }
859: PetscInfo(0,"Using MUMPS for LU factorization and solves.\n");
860: PetscObjectChangeTypeName((PetscObject)B,newtype);
861: *newmat = B;
862: return(0);
863: }
866: /*MC
867: MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed
868: and sequential matrices via the external package MUMPS.
870: If MUMPS is installed (see the manual for instructions
871: on how to declare the existence of external packages),
872: a matrix type can be constructed which invokes MUMPS solvers.
873: After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS).
875: If created with a single process communicator, this matrix type inherits from MATSEQAIJ.
876: Otherwise, this matrix type inherits from MATMPIAIJ. Hence for single process communicators,
877: MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
878: for communicators controlling multiple processes. It is recommended that you call both of
879: the above preallocation routines for simplicity. One can also call MatConvert for an inplace
880: conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
881: without data copy.
883: Options Database Keys:
884: + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions()
885: . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
886: . -mat_mumps_icntl_4 <0,1,2,3,4> - print level
887: . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
888: . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
889: . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
890: . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
891: . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
892: . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
893: . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
894: . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
895: . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
896: . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
897: . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
898: - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
900: Level: beginner
902: .seealso: MATSBAIJMUMPS
903: M*/
908: PetscErrorCode MatCreate_AIJMUMPS(Mat A)
909: {
911: PetscMPIInt size;
912:
914: MPI_Comm_size(A->comm,&size);
915: if (size == 1) {
916: MatSetType(A,MATSEQAIJ);
917: } else {
918: MatSetType(A,MATMPIAIJ);
919: /*
920: Mat A_diag = ((Mat_MPIAIJ *)A->data)->A;
921: MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,MAT_REUSE_MATRIX,&A_diag);
922: */
923: }
924: MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,MAT_REUSE_MATRIX,&A);
925: return(0);
926: }
931: PetscErrorCode MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode)
932: {
934: Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
937: (*mumps->MatAssemblyEnd)(A,mode);
938: mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
939: A->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS;
940: return(0);
941: }
946: PetscErrorCode MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat B,PetscInt bs,PetscInt d_nz,PetscInt *d_nnz,PetscInt o_nz,PetscInt *o_nnz)
947: {
948: Mat A;
949: Mat_MUMPS *mumps=(Mat_MUMPS*)B->spptr;
953: /*
954: After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix
955: into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS. I would
956: like this to be done in the MatCreate routine, but the creation of this inner matrix requires
957: block size info so that PETSc can determine the local size properly. The block size info is set
958: in the preallocation routine.
959: */
960: (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz);
961: A = ((Mat_MPISBAIJ *)B->data)->A;
962: MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);
963: return(0);
964: }
970: PetscErrorCode MatConvert_SBAIJ_SBAIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
971: {
973: PetscMPIInt size;
974: MPI_Comm comm;
975: Mat B=*newmat;
976: Mat_MUMPS *mumps;
977: void (*f)(void);
980: if (reuse == MAT_INITIAL_MATRIX) {
981: MatDuplicate(A,MAT_COPY_VALUES,&B);
982: }
984: PetscObjectGetComm((PetscObject)A,&comm);
985: PetscNew(Mat_MUMPS,&mumps);
987: mumps->MatDuplicate = A->ops->duplicate;
988: mumps->MatView = A->ops->view;
989: mumps->MatAssemblyEnd = A->ops->assemblyend;
990: mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
991: mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
992: mumps->MatDestroy = A->ops->destroy;
993: mumps->specialdestroy = MatDestroy_SBAIJMUMPS;
994: mumps->CleanUpMUMPS = PETSC_FALSE;
995: mumps->isAIJ = PETSC_FALSE;
996:
997: B->spptr = (void*)mumps;
998: B->ops->duplicate = MatDuplicate_MUMPS;
999: B->ops->view = MatView_MUMPS;
1000: B->ops->assemblyend = MatAssemblyEnd_SBAIJMUMPS;
1001: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS;
1002: B->ops->destroy = MatDestroy_MUMPS;
1004: MPI_Comm_size(comm,&size);
1005: if (size == 1) {
1006: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C",
1007: "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);
1008: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C",
1009: "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);
1010: } else {
1011: /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */
1012: PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);
1013: if (f) { /* This case should always be true when this routine is called */
1014: mumps->MatPreallocate = (PetscErrorCode (*)(Mat,int,int,int*,int,int*))f;
1015: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
1016: "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS",
1017: MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);
1018: }
1019: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C",
1020: "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);
1021: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C",
1022: "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);
1023: }
1025: PetscInfo(0,"Using MUMPS for Cholesky factorization and solves.\n");
1026: PetscObjectChangeTypeName((PetscObject)B,newtype);
1027: *newmat = B;
1028: return(0);
1029: }
1034: PetscErrorCode MatDuplicate_MUMPS(Mat A, MatDuplicateOption op, Mat *M) {
1036: Mat_MUMPS *lu=(Mat_MUMPS *)A->spptr;
1039: (*lu->MatDuplicate)(A,op,M);
1040: PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));
1041: return(0);
1042: }
1044: /*MC
1045: MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for
1046: distributed and sequential matrices via the external package MUMPS.
1048: If MUMPS is installed (see the manual for instructions
1049: on how to declare the existence of external packages),
1050: a matrix type can be constructed which invokes MUMPS solvers.
1051: After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS).
1053: If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ.
1054: Otherwise, this matrix type inherits from MATMPISBAIJ. Hence for single process communicators,
1055: MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported
1056: for communicators controlling multiple processes. It is recommended that you call both of
1057: the above preallocation routines for simplicity. One can also call MatConvert for an inplace
1058: conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size)
1059: without data copy.
1061: Options Database Keys:
1062: + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions()
1063: . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
1064: . -mat_mumps_icntl_4 <0,...,4> - print level
1065: . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
1066: . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
1067: . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
1068: . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
1069: . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
1070: . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
1071: . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
1072: . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
1073: . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
1074: . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
1075: . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
1076: - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
1078: Level: beginner
1080: .seealso: MATAIJMUMPS
1081: M*/
1086: PetscErrorCode MatCreate_SBAIJMUMPS(Mat A)
1087: {
1089: PetscMPIInt size;
1092: MPI_Comm_size(A->comm,&size);
1093: if (size == 1) {
1094: MatSetType(A,MATSEQSBAIJ);
1095: } else {
1096: MatSetType(A,MATMPISBAIJ);
1097: }
1098: MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);
1099: return(0);
1100: }