Actual source code: mumps.c
petsc-3.7.7 2017-09-25
2: /*
3: Provides an interface to the MUMPS sparse solver
4: */
6: #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I "petscmat.h" I*/
7: #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
8: #include <petscblaslapack.h>
10: EXTERN_C_BEGIN
11: #if defined(PETSC_USE_COMPLEX)
12: #if defined(PETSC_USE_REAL_SINGLE)
13: #include <cmumps_c.h>
14: #else
15: #include <zmumps_c.h>
16: #endif
17: #else
18: #if defined(PETSC_USE_REAL_SINGLE)
19: #include <smumps_c.h>
20: #else
21: #include <dmumps_c.h>
22: #endif
23: #endif
24: EXTERN_C_END
25: #define JOB_INIT -1
26: #define JOB_FACTSYMBOLIC 1
27: #define JOB_FACTNUMERIC 2
28: #define JOB_SOLVE 3
29: #define JOB_END -2
31: /* calls to MUMPS */
32: #if defined(PETSC_USE_COMPLEX)
33: #if defined(PETSC_USE_REAL_SINGLE)
34: #define PetscMUMPS_c cmumps_c
35: #else
36: #define PetscMUMPS_c zmumps_c
37: #endif
38: #else
39: #if defined(PETSC_USE_REAL_SINGLE)
40: #define PetscMUMPS_c smumps_c
41: #else
42: #define PetscMUMPS_c dmumps_c
43: #endif
44: #endif
46: /* declare MumpsScalar */
47: #if defined(PETSC_USE_COMPLEX)
48: #if defined(PETSC_USE_REAL_SINGLE)
49: #define MumpsScalar mumps_complex
50: #else
51: #define MumpsScalar mumps_double_complex
52: #endif
53: #else
54: #define MumpsScalar PetscScalar
55: #endif
57: /* macros s.t. indices match MUMPS documentation */
58: #define ICNTL(I) icntl[(I)-1]
59: #define CNTL(I) cntl[(I)-1]
60: #define INFOG(I) infog[(I)-1]
61: #define INFO(I) info[(I)-1]
62: #define RINFOG(I) rinfog[(I)-1]
63: #define RINFO(I) rinfo[(I)-1]
65: typedef struct {
66: #if defined(PETSC_USE_COMPLEX)
67: #if defined(PETSC_USE_REAL_SINGLE)
68: CMUMPS_STRUC_C id;
69: #else
70: ZMUMPS_STRUC_C id;
71: #endif
72: #else
73: #if defined(PETSC_USE_REAL_SINGLE)
74: SMUMPS_STRUC_C id;
75: #else
76: DMUMPS_STRUC_C id;
77: #endif
78: #endif
80: MatStructure matstruc;
81: PetscMPIInt myid,size;
82: PetscInt *irn,*jcn,nz,sym;
83: PetscScalar *val;
84: MPI_Comm comm_mumps;
85: PetscBool isAIJ;
86: PetscInt ICNTL9_pre; /* check if ICNTL(9) is changed from previous MatSolve */
87: VecScatter scat_rhs, scat_sol; /* used by MatSolve() */
88: Vec b_seq,x_seq;
89: PetscInt ninfo,*info; /* display INFO */
90: PetscInt sizeredrhs;
91: PetscInt *schur_pivots;
92: PetscInt schur_B_lwork;
93: PetscScalar *schur_work;
94: PetscScalar *schur_sol;
95: PetscInt schur_sizesol;
96: PetscBool schur_factored;
97: PetscBool schur_inverted;
99: PetscErrorCode (*Destroy)(Mat);
100: PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
101: } Mat_MUMPS;
103: extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
107: static PetscErrorCode MatMumpsResetSchur_Private(Mat_MUMPS* mumps)
108: {
112: PetscFree2(mumps->id.listvar_schur,mumps->id.schur);
113: PetscFree(mumps->id.redrhs);
114: PetscFree(mumps->schur_sol);
115: PetscFree(mumps->schur_pivots);
116: PetscFree(mumps->schur_work);
117: mumps->id.size_schur = 0;
118: mumps->id.ICNTL(19) = 0;
119: return(0);
120: }
124: static PetscErrorCode MatMumpsFactorSchur_Private(Mat_MUMPS* mumps)
125: {
126: PetscBLASInt B_N,B_ierr,B_slda;
130: if (mumps->schur_factored) {
131: return(0);
132: }
133: PetscBLASIntCast(mumps->id.size_schur,&B_N);
134: PetscBLASIntCast(mumps->id.schur_lld,&B_slda);
135: if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
136: if (!mumps->schur_pivots) {
137: PetscMalloc1(B_N,&mumps->schur_pivots);
138: }
139: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
140: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,&B_ierr));
141: PetscFPTrapPop();
142: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
143: } else { /* either full or lower-triangular (not packed) */
144: char ord[2];
145: if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
146: sprintf(ord,"L");
147: } else { /* ICNTL(19) == 1 lower triangular stored by rows */
148: sprintf(ord,"U");
149: }
150: if (mumps->id.sym == 2) {
151: if (!mumps->schur_pivots) {
152: PetscScalar lwork;
154: PetscMalloc1(B_N,&mumps->schur_pivots);
155: mumps->schur_B_lwork=-1;
156: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
157: PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,&lwork,&mumps->schur_B_lwork,&B_ierr));
158: PetscFPTrapPop();
159: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYTRF Lapack routine %d",(int)B_ierr);
160: PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&mumps->schur_B_lwork);
161: PetscMalloc1(mumps->schur_B_lwork,&mumps->schur_work);
162: }
163: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
164: PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,mumps->schur_work,&mumps->schur_B_lwork,&B_ierr));
165: PetscFPTrapPop();
166: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
167: } else {
168: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
169: PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_slda,&B_ierr));
170: PetscFPTrapPop();
171: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
172: }
173: }
174: mumps->schur_factored = PETSC_TRUE;
175: return(0);
176: }
180: static PetscErrorCode MatMumpsInvertSchur_Private(Mat_MUMPS* mumps)
181: {
182: PetscBLASInt B_N,B_ierr,B_slda;
186: MatMumpsFactorSchur_Private(mumps);
187: PetscBLASIntCast(mumps->id.size_schur,&B_N);
188: PetscBLASIntCast(mumps->id.schur_lld,&B_slda);
189: if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
190: if (!mumps->schur_work) {
191: PetscScalar lwork;
193: mumps->schur_B_lwork = -1;
194: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
195: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,(PetscScalar*)mumps->id.schur,&B_N,mumps->schur_pivots,&lwork,&mumps->schur_B_lwork,&B_ierr));
196: PetscFPTrapPop();
197: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
198: PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&mumps->schur_B_lwork);
199: PetscMalloc1(mumps->schur_B_lwork,&mumps->schur_work);
200: }
201: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
202: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,(PetscScalar*)mumps->id.schur,&B_N,mumps->schur_pivots,mumps->schur_work,&mumps->schur_B_lwork,&B_ierr));
203: PetscFPTrapPop();
204: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
205: } else { /* either full or lower-triangular (not packed) */
206: char ord[2];
207: if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
208: sprintf(ord,"L");
209: } else { /* ICNTL(19) == 1 lower triangular stored by rows */
210: sprintf(ord,"U");
211: }
212: if (mumps->id.sym == 2) {
213: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
214: PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_N,mumps->schur_pivots,mumps->schur_work,&B_ierr));
215: PetscFPTrapPop();
216: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
217: } else {
218: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
219: PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_(ord,&B_N,(PetscScalar*)mumps->id.schur,&B_N,&B_ierr));
220: PetscFPTrapPop();
221: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
222: }
223: }
224: mumps->schur_inverted = PETSC_TRUE;
225: return(0);
226: }
230: static PetscErrorCode MatMumpsSolveSchur_Private(Mat_MUMPS* mumps, PetscBool sol_in_redrhs)
231: {
232: PetscBLASInt B_N,B_Nrhs,B_ierr,B_slda,B_rlda;
233: PetscScalar one=1.,zero=0.;
237: MatMumpsFactorSchur_Private(mumps);
238: PetscBLASIntCast(mumps->id.size_schur,&B_N);
239: PetscBLASIntCast(mumps->id.schur_lld,&B_slda);
240: PetscBLASIntCast(mumps->id.nrhs,&B_Nrhs);
241: PetscBLASIntCast(mumps->id.lredrhs,&B_rlda);
242: if (mumps->schur_inverted) {
243: PetscInt sizesol = B_Nrhs*B_N;
244: if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
245: PetscFree(mumps->schur_sol);
246: PetscMalloc1(sizesol,&mumps->schur_sol);
247: mumps->schur_sizesol = sizesol;
248: }
249: if (!mumps->sym) {
250: char type[2];
251: if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
252: if (!mumps->id.ICNTL(9)) { /* transpose solve */
253: sprintf(type,"N");
254: } else {
255: sprintf(type,"T");
256: }
257: } else { /* stored by columns */
258: if (!mumps->id.ICNTL(9)) { /* transpose solve */
259: sprintf(type,"T");
260: } else {
261: sprintf(type,"N");
262: }
263: }
264: PetscStackCallBLAS("BLASgemm",BLASgemm_(type,"N",&B_N,&B_Nrhs,&B_N,&one,(PetscScalar*)mumps->id.schur,&B_slda,(PetscScalar*)mumps->id.redrhs,&B_rlda,&zero,mumps->schur_sol,&B_rlda));
265: } else {
266: char ord[2];
267: if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
268: sprintf(ord,"L");
269: } else { /* ICNTL(19) == 1 lower triangular stored by rows */
270: sprintf(ord,"U");
271: }
272: PetscStackCallBLAS("BLASsymm",BLASsymm_("L",ord,&B_N,&B_Nrhs,&one,(PetscScalar*)mumps->id.schur,&B_slda,(PetscScalar*)mumps->id.redrhs,&B_rlda,&zero,mumps->schur_sol,&B_rlda));
273: }
274: if (sol_in_redrhs) {
275: PetscMemcpy(mumps->id.redrhs,mumps->schur_sol,sizesol*sizeof(PetscScalar));
276: }
277: } else { /* Schur complement has not been inverted */
278: MumpsScalar *orhs=NULL;
280: if (!sol_in_redrhs) {
281: PetscInt sizesol = B_Nrhs*B_N;
282: if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
283: PetscFree(mumps->schur_sol);
284: PetscMalloc1(sizesol,&mumps->schur_sol);
285: mumps->schur_sizesol = sizesol;
286: }
287: orhs = mumps->id.redrhs;
288: PetscMemcpy(mumps->schur_sol,mumps->id.redrhs,sizesol*sizeof(PetscScalar));
289: mumps->id.redrhs = (MumpsScalar*)mumps->schur_sol;
290: }
291: if (!mumps->sym) { /* MUMPS always return a full Schur matrix */
292: char type[2];
293: if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
294: if (!mumps->id.ICNTL(9)) { /* transpose solve */
295: sprintf(type,"N");
296: } else {
297: sprintf(type,"T");
298: }
299: } else { /* stored by columns */
300: if (!mumps->id.ICNTL(9)) { /* transpose solve */
301: sprintf(type,"T");
302: } else {
303: sprintf(type,"N");
304: }
305: }
306: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
307: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_(type,&B_N,&B_Nrhs,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,(PetscScalar*)mumps->id.redrhs,&B_rlda,&B_ierr));
308: PetscFPTrapPop();
309: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRS Lapack routine %d",(int)B_ierr);
310: } else { /* either full or lower-triangular (not packed) */
311: char ord[2];
312: if (mumps->id.ICNTL(19) == 2 || mumps->id.ICNTL(19) == 3) { /* lower triangular stored by columns or full matrix */
313: sprintf(ord,"L");
314: } else { /* ICNTL(19) == 1 lower triangular stored by rows */
315: sprintf(ord,"U");
316: }
317: if (mumps->id.sym == 2) {
318: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
319: PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_(ord,&B_N,&B_Nrhs,(PetscScalar*)mumps->id.schur,&B_slda,mumps->schur_pivots,(PetscScalar*)mumps->id.redrhs,&B_rlda,&B_ierr));
320: PetscFPTrapPop();
321: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRS Lapack routine %d",(int)B_ierr);
322: } else {
323: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
324: PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_(ord,&B_N,&B_Nrhs,(PetscScalar*)mumps->id.schur,&B_slda,(PetscScalar*)mumps->id.redrhs,&B_rlda,&B_ierr));
325: PetscFPTrapPop();
326: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRS Lapack routine %d",(int)B_ierr);
327: }
328: }
329: if (!sol_in_redrhs) {
330: mumps->id.redrhs = orhs;
331: }
332: }
333: return(0);
334: }
338: static PetscErrorCode MatMumpsHandleSchur_Private(Mat_MUMPS* mumps, PetscBool expansion)
339: {
343: if (!mumps->id.ICNTL(19)) { /* do nothing when Schur complement has not been computed */
344: return(0);
345: }
346: if (!expansion) { /* prepare for the condensation step */
347: PetscInt sizeredrhs = mumps->id.nrhs*mumps->id.size_schur;
348: /* allocate MUMPS internal array to store reduced right-hand sides */
349: if (!mumps->id.redrhs || sizeredrhs > mumps->sizeredrhs) {
350: PetscFree(mumps->id.redrhs);
351: mumps->id.lredrhs = mumps->id.size_schur;
352: PetscMalloc1(mumps->id.nrhs*mumps->id.lredrhs,&mumps->id.redrhs);
353: mumps->sizeredrhs = mumps->id.nrhs*mumps->id.lredrhs;
354: }
355: mumps->id.ICNTL(26) = 1; /* condensation phase */
356: } else { /* prepare for the expansion step */
357: /* solve Schur complement (this has to be done by the MUMPS user, so basically us) */
358: MatMumpsSolveSchur_Private(mumps,PETSC_TRUE);
359: mumps->id.ICNTL(26) = 2; /* expansion phase */
360: PetscMUMPS_c(&mumps->id);
361: if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
362: /* restore defaults */
363: mumps->id.ICNTL(26) = -1;
364: }
365: return(0);
366: }
368: /*
369: MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz]
371: input:
372: A - matrix in aij,baij or sbaij (bs=1) format
373: shift - 0: C style output triple; 1: Fortran style output triple.
374: reuse - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
375: MAT_REUSE_MATRIX: only the values in v array are updated
376: output:
377: nnz - dim of r, c, and v (number of local nonzero entries of A)
378: r, c, v - row and col index, matrix values (matrix triples)
380: The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is
381: freed with PetscFree((mumps->irn); This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means
382: that the PetscMalloc() cannot easily be replaced with a PetscMalloc3().
384: */
388: PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
389: {
390: const PetscInt *ai,*aj,*ajj,M=A->rmap->n;
391: PetscInt nz,rnz,i,j;
393: PetscInt *row,*col;
394: Mat_SeqAIJ *aa=(Mat_SeqAIJ*)A->data;
397: *v=aa->a;
398: if (reuse == MAT_INITIAL_MATRIX) {
399: nz = aa->nz;
400: ai = aa->i;
401: aj = aa->j;
402: *nnz = nz;
403: PetscMalloc1(2*nz, &row);
404: col = row + nz;
406: nz = 0;
407: for (i=0; i<M; i++) {
408: rnz = ai[i+1] - ai[i];
409: ajj = aj + ai[i];
410: for (j=0; j<rnz; j++) {
411: row[nz] = i+shift; col[nz++] = ajj[j] + shift;
412: }
413: }
414: *r = row; *c = col;
415: }
416: return(0);
417: }
421: PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
422: {
423: Mat_SeqBAIJ *aa=(Mat_SeqBAIJ*)A->data;
424: const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2;
425: PetscInt bs,M,nz,idx=0,rnz,i,j,k,m;
427: PetscInt *row,*col;
430: MatGetBlockSize(A,&bs);
431: M = A->rmap->N/bs;
432: *v = aa->a;
433: if (reuse == MAT_INITIAL_MATRIX) {
434: ai = aa->i; aj = aa->j;
435: nz = bs2*aa->nz;
436: *nnz = nz;
437: PetscMalloc1(2*nz, &row);
438: col = row + nz;
440: for (i=0; i<M; i++) {
441: ajj = aj + ai[i];
442: rnz = ai[i+1] - ai[i];
443: for (k=0; k<rnz; k++) {
444: for (j=0; j<bs; j++) {
445: for (m=0; m<bs; m++) {
446: row[idx] = i*bs + m + shift;
447: col[idx++] = bs*(ajj[k]) + j + shift;
448: }
449: }
450: }
451: }
452: *r = row; *c = col;
453: }
454: return(0);
455: }
459: PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
460: {
461: const PetscInt *ai, *aj,*ajj,M=A->rmap->n;
462: PetscInt nz,rnz,i,j;
464: PetscInt *row,*col;
465: Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)A->data;
468: *v = aa->a;
469: if (reuse == MAT_INITIAL_MATRIX) {
470: nz = aa->nz;
471: ai = aa->i;
472: aj = aa->j;
473: *v = aa->a;
474: *nnz = nz;
475: PetscMalloc1(2*nz, &row);
476: col = row + nz;
478: nz = 0;
479: for (i=0; i<M; i++) {
480: rnz = ai[i+1] - ai[i];
481: ajj = aj + ai[i];
482: for (j=0; j<rnz; j++) {
483: row[nz] = i+shift; col[nz++] = ajj[j] + shift;
484: }
485: }
486: *r = row; *c = col;
487: }
488: return(0);
489: }
493: PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
494: {
495: const PetscInt *ai,*aj,*ajj,*adiag,M=A->rmap->n;
496: PetscInt nz,rnz,i,j;
497: const PetscScalar *av,*v1;
498: PetscScalar *val;
499: PetscErrorCode ierr;
500: PetscInt *row,*col;
501: Mat_SeqAIJ *aa=(Mat_SeqAIJ*)A->data;
502: PetscBool missing;
505: ai =aa->i; aj=aa->j;av=aa->a;
506: adiag=aa->diag;
507: MatMissingDiagonal_SeqAIJ(A,&missing,&i);
508: if (reuse == MAT_INITIAL_MATRIX) {
509: /* count nz in the uppper triangular part of A */
510: nz = 0;
511: if (missing) {
512: for (i=0; i<M; i++) {
513: if (PetscUnlikely(adiag[i] >= ai[i+1])) {
514: for (j=ai[i];j<ai[i+1];j++) {
515: if (aj[j] < i) continue;
516: nz++;
517: }
518: } else {
519: nz += ai[i+1] - adiag[i];
520: }
521: }
522: } else {
523: for (i=0; i<M; i++) nz += ai[i+1] - adiag[i];
524: }
525: *nnz = nz;
527: PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
528: col = row + nz;
529: val = (PetscScalar*)(col + nz);
531: nz = 0;
532: if (missing) {
533: for (i=0; i<M; i++) {
534: if (PetscUnlikely(adiag[i] >= ai[i+1])) {
535: for (j=ai[i];j<ai[i+1];j++) {
536: if (aj[j] < i) continue;
537: row[nz] = i+shift;
538: col[nz] = aj[j]+shift;
539: val[nz] = av[j];
540: nz++;
541: }
542: } else {
543: rnz = ai[i+1] - adiag[i];
544: ajj = aj + adiag[i];
545: v1 = av + adiag[i];
546: for (j=0; j<rnz; j++) {
547: row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
548: }
549: }
550: }
551: } else {
552: for (i=0; i<M; i++) {
553: rnz = ai[i+1] - adiag[i];
554: ajj = aj + adiag[i];
555: v1 = av + adiag[i];
556: for (j=0; j<rnz; j++) {
557: row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
558: }
559: }
560: }
561: *r = row; *c = col; *v = val;
562: } else {
563: nz = 0; val = *v;
564: if (missing) {
565: for (i=0; i <M; i++) {
566: if (PetscUnlikely(adiag[i] >= ai[i+1])) {
567: for (j=ai[i];j<ai[i+1];j++) {
568: if (aj[j] < i) continue;
569: val[nz++] = av[j];
570: }
571: } else {
572: rnz = ai[i+1] - adiag[i];
573: v1 = av + adiag[i];
574: for (j=0; j<rnz; j++) {
575: val[nz++] = v1[j];
576: }
577: }
578: }
579: } else {
580: for (i=0; i <M; i++) {
581: rnz = ai[i+1] - adiag[i];
582: v1 = av + adiag[i];
583: for (j=0; j<rnz; j++) {
584: val[nz++] = v1[j];
585: }
586: }
587: }
588: }
589: return(0);
590: }
594: PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
595: {
596: const PetscInt *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
597: PetscErrorCode ierr;
598: PetscInt rstart,nz,i,j,jj,irow,countA,countB;
599: PetscInt *row,*col;
600: const PetscScalar *av, *bv,*v1,*v2;
601: PetscScalar *val;
602: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data;
603: Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)(mat->A)->data;
604: Mat_SeqBAIJ *bb = (Mat_SeqBAIJ*)(mat->B)->data;
607: ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
608: av=aa->a; bv=bb->a;
610: garray = mat->garray;
612: if (reuse == MAT_INITIAL_MATRIX) {
613: nz = aa->nz + bb->nz;
614: *nnz = nz;
615: PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
616: col = row + nz;
617: val = (PetscScalar*)(col + nz);
619: *r = row; *c = col; *v = val;
620: } else {
621: row = *r; col = *c; val = *v;
622: }
624: jj = 0; irow = rstart;
625: for (i=0; i<m; i++) {
626: ajj = aj + ai[i]; /* ptr to the beginning of this row */
627: countA = ai[i+1] - ai[i];
628: countB = bi[i+1] - bi[i];
629: bjj = bj + bi[i];
630: v1 = av + ai[i];
631: v2 = bv + bi[i];
633: /* A-part */
634: for (j=0; j<countA; j++) {
635: if (reuse == MAT_INITIAL_MATRIX) {
636: row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
637: }
638: val[jj++] = v1[j];
639: }
641: /* B-part */
642: for (j=0; j < countB; j++) {
643: if (reuse == MAT_INITIAL_MATRIX) {
644: row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
645: }
646: val[jj++] = v2[j];
647: }
648: irow++;
649: }
650: return(0);
651: }
655: PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
656: {
657: const PetscInt *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
658: PetscErrorCode ierr;
659: PetscInt rstart,nz,i,j,jj,irow,countA,countB;
660: PetscInt *row,*col;
661: const PetscScalar *av, *bv,*v1,*v2;
662: PetscScalar *val;
663: Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data;
664: Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(mat->A)->data;
665: Mat_SeqAIJ *bb = (Mat_SeqAIJ*)(mat->B)->data;
668: ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
669: av=aa->a; bv=bb->a;
671: garray = mat->garray;
673: if (reuse == MAT_INITIAL_MATRIX) {
674: nz = aa->nz + bb->nz;
675: *nnz = nz;
676: PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
677: col = row + nz;
678: val = (PetscScalar*)(col + nz);
680: *r = row; *c = col; *v = val;
681: } else {
682: row = *r; col = *c; val = *v;
683: }
685: jj = 0; irow = rstart;
686: for (i=0; i<m; i++) {
687: ajj = aj + ai[i]; /* ptr to the beginning of this row */
688: countA = ai[i+1] - ai[i];
689: countB = bi[i+1] - bi[i];
690: bjj = bj + bi[i];
691: v1 = av + ai[i];
692: v2 = bv + bi[i];
694: /* A-part */
695: for (j=0; j<countA; j++) {
696: if (reuse == MAT_INITIAL_MATRIX) {
697: row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
698: }
699: val[jj++] = v1[j];
700: }
702: /* B-part */
703: for (j=0; j < countB; j++) {
704: if (reuse == MAT_INITIAL_MATRIX) {
705: row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
706: }
707: val[jj++] = v2[j];
708: }
709: irow++;
710: }
711: return(0);
712: }
716: PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
717: {
718: Mat_MPIBAIJ *mat = (Mat_MPIBAIJ*)A->data;
719: Mat_SeqBAIJ *aa = (Mat_SeqBAIJ*)(mat->A)->data;
720: Mat_SeqBAIJ *bb = (Mat_SeqBAIJ*)(mat->B)->data;
721: const PetscInt *ai = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
722: const PetscInt *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
723: const PetscInt bs2=mat->bs2;
724: PetscErrorCode ierr;
725: PetscInt bs,nz,i,j,k,n,jj,irow,countA,countB,idx;
726: PetscInt *row,*col;
727: const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2;
728: PetscScalar *val;
731: MatGetBlockSize(A,&bs);
732: if (reuse == MAT_INITIAL_MATRIX) {
733: nz = bs2*(aa->nz + bb->nz);
734: *nnz = nz;
735: PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
736: col = row + nz;
737: val = (PetscScalar*)(col + nz);
739: *r = row; *c = col; *v = val;
740: } else {
741: row = *r; col = *c; val = *v;
742: }
744: jj = 0; irow = rstart;
745: for (i=0; i<mbs; i++) {
746: countA = ai[i+1] - ai[i];
747: countB = bi[i+1] - bi[i];
748: ajj = aj + ai[i];
749: bjj = bj + bi[i];
750: v1 = av + bs2*ai[i];
751: v2 = bv + bs2*bi[i];
753: idx = 0;
754: /* A-part */
755: for (k=0; k<countA; k++) {
756: for (j=0; j<bs; j++) {
757: for (n=0; n<bs; n++) {
758: if (reuse == MAT_INITIAL_MATRIX) {
759: row[jj] = irow + n + shift;
760: col[jj] = rstart + bs*ajj[k] + j + shift;
761: }
762: val[jj++] = v1[idx++];
763: }
764: }
765: }
767: idx = 0;
768: /* B-part */
769: for (k=0; k<countB; k++) {
770: for (j=0; j<bs; j++) {
771: for (n=0; n<bs; n++) {
772: if (reuse == MAT_INITIAL_MATRIX) {
773: row[jj] = irow + n + shift;
774: col[jj] = bs*garray[bjj[k]] + j + shift;
775: }
776: val[jj++] = v2[idx++];
777: }
778: }
779: }
780: irow += bs;
781: }
782: return(0);
783: }
787: PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
788: {
789: const PetscInt *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
790: PetscErrorCode ierr;
791: PetscInt rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
792: PetscInt *row,*col;
793: const PetscScalar *av, *bv,*v1,*v2;
794: PetscScalar *val;
795: Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data;
796: Mat_SeqAIJ *aa =(Mat_SeqAIJ*)(mat->A)->data;
797: Mat_SeqAIJ *bb =(Mat_SeqAIJ*)(mat->B)->data;
800: ai=aa->i; aj=aa->j; adiag=aa->diag;
801: bi=bb->i; bj=bb->j; garray = mat->garray;
802: av=aa->a; bv=bb->a;
804: rstart = A->rmap->rstart;
806: if (reuse == MAT_INITIAL_MATRIX) {
807: nza = 0; /* num of upper triangular entries in mat->A, including diagonals */
808: nzb = 0; /* num of upper triangular entries in mat->B */
809: for (i=0; i<m; i++) {
810: nza += (ai[i+1] - adiag[i]);
811: countB = bi[i+1] - bi[i];
812: bjj = bj + bi[i];
813: for (j=0; j<countB; j++) {
814: if (garray[bjj[j]] > rstart) nzb++;
815: }
816: }
818: nz = nza + nzb; /* total nz of upper triangular part of mat */
819: *nnz = nz;
820: PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);
821: col = row + nz;
822: val = (PetscScalar*)(col + nz);
824: *r = row; *c = col; *v = val;
825: } else {
826: row = *r; col = *c; val = *v;
827: }
829: jj = 0; irow = rstart;
830: for (i=0; i<m; i++) {
831: ajj = aj + adiag[i]; /* ptr to the beginning of the diagonal of this row */
832: v1 = av + adiag[i];
833: countA = ai[i+1] - adiag[i];
834: countB = bi[i+1] - bi[i];
835: bjj = bj + bi[i];
836: v2 = bv + bi[i];
838: /* A-part */
839: for (j=0; j<countA; j++) {
840: if (reuse == MAT_INITIAL_MATRIX) {
841: row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
842: }
843: val[jj++] = v1[j];
844: }
846: /* B-part */
847: for (j=0; j < countB; j++) {
848: if (garray[bjj[j]] > rstart) {
849: if (reuse == MAT_INITIAL_MATRIX) {
850: row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
851: }
852: val[jj++] = v2[j];
853: }
854: }
855: irow++;
856: }
857: return(0);
858: }
862: PetscErrorCode MatGetDiagonal_MUMPS(Mat A,Vec v)
863: {
865: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: MUMPS factor");
866: return(0);
867: }
871: PetscErrorCode MatDestroy_MUMPS(Mat A)
872: {
873: Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
877: PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);
878: VecScatterDestroy(&mumps->scat_rhs);
879: VecScatterDestroy(&mumps->scat_sol);
880: VecDestroy(&mumps->b_seq);
881: VecDestroy(&mumps->x_seq);
882: PetscFree(mumps->id.perm_in);
883: PetscFree(mumps->irn);
884: PetscFree(mumps->info);
885: MatMumpsResetSchur_Private(mumps);
886: mumps->id.job = JOB_END;
887: PetscMUMPS_c(&mumps->id);
888: MPI_Comm_free(&mumps->comm_mumps);
889: if (mumps->Destroy) {
890: (mumps->Destroy)(A);
891: }
892: PetscFree(A->spptr);
894: /* clear composed functions */
895: PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);
896: PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);
897: PetscObjectComposeFunction((PetscObject)A,"MatFactorInvertSchurComplement_C",NULL);
898: PetscObjectComposeFunction((PetscObject)A,"MatFactorCreateSchurComplement_C",NULL);
899: PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSchurComplement_C",NULL);
900: PetscObjectComposeFunction((PetscObject)A,"MatFactorSolveSchurComplement_C",NULL);
901: PetscObjectComposeFunction((PetscObject)A,"MatFactorSolveSchurComplementTranspose_C",NULL);
902: PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);
903: PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);
904: PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);
905: PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);
906: PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);
907: PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);
908: PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);
909: PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);
910: return(0);
911: }
915: PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
916: {
917: Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
918: PetscScalar *array;
919: Vec b_seq;
920: IS is_iden,is_petsc;
921: PetscErrorCode ierr;
922: PetscInt i;
923: PetscBool second_solve = PETSC_FALSE;
924: static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE;
927: PetscCitationsRegister("@article{MUMPS01,\n author = {P.~R. Amestoy and I.~S. Duff and J.-Y. L'Excellent and J. Koster},\n title = {A fully asynchronous multifrontal solver using distributed dynamic scheduling},\n journal = {SIAM Journal on Matrix Analysis and Applications},\n volume = {23},\n number = {1},\n pages = {15--41},\n year = {2001}\n}\n",&cite1);
928: PetscCitationsRegister("@article{MUMPS02,\n author = {P.~R. Amestoy and A. Guermouche and J.-Y. L'Excellent and S. Pralet},\n title = {Hybrid scheduling for the parallel solution of linear systems},\n journal = {Parallel Computing},\n volume = {32},\n number = {2},\n pages = {136--156},\n year = {2006}\n}\n",&cite2);
930: if (A->errortype) {
931: PetscInfo2(A,"MatSolve is called with singular matrix factor, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
932: VecSetInf(x);
933: return(0);
934: }
936: mumps->id.nrhs = 1;
937: b_seq = mumps->b_seq;
938: if (mumps->size > 1) {
939: /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
940: VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);
941: VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);
942: if (!mumps->myid) {VecGetArray(b_seq,&array);}
943: } else { /* size == 1 */
944: VecCopy(b,x);
945: VecGetArray(x,&array);
946: }
947: if (!mumps->myid) { /* define rhs on the host */
948: mumps->id.nrhs = 1;
949: mumps->id.rhs = (MumpsScalar*)array;
950: }
952: /*
953: handle condensation step of Schur complement (if any)
954: We set by default ICNTL(26) == -1 when Schur indices have been provided by the user.
955: According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase
956: Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system.
957: This requires an extra call to PetscMUMPS_c and the computation of the factors for S
958: */
959: if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) {
960: second_solve = PETSC_TRUE;
961: MatMumpsHandleSchur_Private(mumps,PETSC_FALSE);
962: }
963: /* solve phase */
964: /*-------------*/
965: mumps->id.job = JOB_SOLVE;
966: PetscMUMPS_c(&mumps->id);
967: if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
969: /* handle expansion step of Schur complement (if any) */
970: if (second_solve) {
971: MatMumpsHandleSchur_Private(mumps,PETSC_TRUE);
972: }
974: if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */
975: if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
976: /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
977: VecScatterDestroy(&mumps->scat_sol);
978: }
979: if (!mumps->scat_sol) { /* create scatter scat_sol */
980: ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden); /* from */
981: for (i=0; i<mumps->id.lsol_loc; i++) {
982: mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */
983: }
984: ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc); /* to */
985: VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);
986: ISDestroy(&is_iden);
987: ISDestroy(&is_petsc);
989: mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */
990: }
992: VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
993: VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);
994: }
995: return(0);
996: }
1000: PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
1001: {
1002: Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
1006: mumps->id.ICNTL(9) = 0;
1007: MatSolve_MUMPS(A,b,x);
1008: mumps->id.ICNTL(9) = 1;
1009: return(0);
1010: }
1014: PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
1015: {
1017: PetscBool flg;
1018: Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
1019: PetscInt i,nrhs,M;
1020: PetscScalar *array,*bray;
1023: PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
1024: if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
1025: PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
1026: if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
1027: if (B->rmap->n != X->rmap->n) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution");
1029: MatGetSize(B,&M,&nrhs);
1030: mumps->id.nrhs = nrhs;
1031: mumps->id.lrhs = M;
1033: if (mumps->size == 1) {
1034: /* copy B to X */
1035: MatDenseGetArray(B,&bray);
1036: MatDenseGetArray(X,&array);
1037: PetscMemcpy(array,bray,M*nrhs*sizeof(PetscScalar));
1038: MatDenseRestoreArray(B,&bray);
1039: mumps->id.rhs = (MumpsScalar*)array;
1040: /* handle condensation step of Schur complement (if any) */
1041: MatMumpsHandleSchur_Private(mumps,PETSC_FALSE);
1043: /* solve phase */
1044: /*-------------*/
1045: mumps->id.job = JOB_SOLVE;
1046: PetscMUMPS_c(&mumps->id);
1047: if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
1049: /* handle expansion step of Schur complement (if any) */
1050: MatMumpsHandleSchur_Private(mumps,PETSC_TRUE);
1051: MatDenseRestoreArray(X,&array);
1052: } else { /*--------- parallel case --------*/
1053: PetscInt lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save;
1054: MumpsScalar *sol_loc,*sol_loc_save;
1055: IS is_to,is_from;
1056: PetscInt k,proc,j,m;
1057: const PetscInt *rstart;
1058: Vec v_mpi,b_seq,x_seq;
1059: VecScatter scat_rhs,scat_sol;
1061: /* create x_seq to hold local solution */
1062: isol_loc_save = mumps->id.isol_loc; /* save it for MatSovle() */
1063: sol_loc_save = mumps->id.sol_loc;
1065: lsol_loc = mumps->id.INFO(23);
1066: nlsol_loc = nrhs*lsol_loc; /* length of sol_loc */
1067: PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);
1068: mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1069: mumps->id.isol_loc = isol_loc;
1071: VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,(PetscScalar*)sol_loc,&x_seq);
1073: /* copy rhs matrix B into vector v_mpi */
1074: MatGetLocalSize(B,&m,NULL);
1075: MatDenseGetArray(B,&bray);
1076: VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);
1077: MatDenseRestoreArray(B,&bray);
1079: /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */
1080: /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B;
1081: iidx: inverse of idx, will be used by scattering xx_seq -> X */
1082: PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);
1083: MatGetOwnershipRanges(B,&rstart);
1084: k = 0;
1085: for (proc=0; proc<mumps->size; proc++){
1086: for (j=0; j<nrhs; j++){
1087: for (i=rstart[proc]; i<rstart[proc+1]; i++){
1088: iidx[j*M + i] = k;
1089: idx[k++] = j*M + i;
1090: }
1091: }
1092: }
1094: if (!mumps->myid) {
1095: VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);
1096: ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);
1097: ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);
1098: } else {
1099: VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);
1100: ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);
1101: ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);
1102: }
1103: VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);
1104: VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);
1105: ISDestroy(&is_to);
1106: ISDestroy(&is_from);
1107: VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);
1109: if (!mumps->myid) { /* define rhs on the host */
1110: VecGetArray(b_seq,&bray);
1111: mumps->id.rhs = (MumpsScalar*)bray;
1112: VecRestoreArray(b_seq,&bray);
1113: }
1115: /* solve phase */
1116: /*-------------*/
1117: mumps->id.job = JOB_SOLVE;
1118: PetscMUMPS_c(&mumps->id);
1119: if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
1121: /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */
1122: MatDenseGetArray(X,&array);
1123: VecPlaceArray(v_mpi,array);
1124:
1125: /* create scatter scat_sol */
1126: PetscMalloc1(nlsol_loc,&idxx);
1127: ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);
1128: for (i=0; i<lsol_loc; i++) {
1129: isol_loc[i] -= 1; /* change Fortran style to C style */
1130: idxx[i] = iidx[isol_loc[i]];
1131: for (j=1; j<nrhs; j++){
1132: idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M];
1133: }
1134: }
1135: ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);
1136: VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);
1137: VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);
1138: ISDestroy(&is_from);
1139: ISDestroy(&is_to);
1140: VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);
1141: MatDenseRestoreArray(X,&array);
1143: /* free spaces */
1144: mumps->id.sol_loc = sol_loc_save;
1145: mumps->id.isol_loc = isol_loc_save;
1147: PetscFree2(sol_loc,isol_loc);
1148: PetscFree2(idx,iidx);
1149: PetscFree(idxx);
1150: VecDestroy(&x_seq);
1151: VecDestroy(&v_mpi);
1152: VecDestroy(&b_seq);
1153: VecScatterDestroy(&scat_rhs);
1154: VecScatterDestroy(&scat_sol);
1155: }
1156: return(0);
1157: }
1159: #if !defined(PETSC_USE_COMPLEX)
1160: /*
1161: input:
1162: F: numeric factor
1163: output:
1164: nneg: total number of negative pivots
1165: nzero: 0
1166: npos: (global dimension of F) - nneg
1167: */
1171: PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
1172: {
1173: Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1175: PetscMPIInt size;
1178: MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);
1179: /* 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 */
1180: if (size > 1 && mumps->id.ICNTL(13) != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",mumps->id.INFOG(13));
1182: if (nneg) *nneg = mumps->id.INFOG(12);
1183: if (nzero || npos) {
1184: if (mumps->id.ICNTL(24) != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"-mat_mumps_icntl_24 must be set as 1 for null pivot row detection");
1185: if (nzero) *nzero = mumps->id.INFOG(28);
1186: if (npos) *npos = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
1187: }
1188: return(0);
1189: }
1190: #endif /* !defined(PETSC_USE_COMPLEX) */
1194: PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
1195: {
1196: Mat_MUMPS *mumps =(Mat_MUMPS*)(F)->spptr;
1198: Mat F_diag;
1199: PetscBool isMPIAIJ;
1202: if (mumps->id.INFOG(1) < 0) {
1203: if (mumps->id.INFOG(1) == -6) {
1204: PetscInfo2(A,"MatFactorNumeric is called with singular matrix structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1205: }
1206: PetscInfo2(A,"MatFactorNumeric is called after analysis phase fails, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1207: return(0);
1208: }
1210: (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);
1212: /* numerical factorization phase */
1213: /*-------------------------------*/
1214: mumps->id.job = JOB_FACTNUMERIC;
1215: if (!mumps->id.ICNTL(18)) { /* A is centralized */
1216: if (!mumps->myid) {
1217: mumps->id.a = (MumpsScalar*)mumps->val;
1218: }
1219: } else {
1220: mumps->id.a_loc = (MumpsScalar*)mumps->val;
1221: }
1222: PetscMUMPS_c(&mumps->id);
1223: if (mumps->id.INFOG(1) < 0) {
1224: if (A->erroriffailure) {
1225: SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1226: } else {
1227: if (mumps->id.INFOG(1) == -10) { /* numerically singular matrix */
1228: PetscInfo2(F,"matrix is numerically singular, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1229: F->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1230: } else if (mumps->id.INFOG(1) == -13) {
1231: PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, cannot allocate required memory %d megabytes\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1232: F->errortype = MAT_FACTOR_OUTMEMORY;
1233: } else if (mumps->id.INFOG(1) == -8 || mumps->id.INFOG(1) == -9 || (-16 < mumps->id.INFOG(1) && mumps->id.INFOG(1) < -10) ) {
1234: PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d, problem with workarray \n",mumps->id.INFOG(1),mumps->id.INFO(2));
1235: F->errortype = MAT_FACTOR_OUTMEMORY;
1236: } else {
1237: PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1238: F->errortype = MAT_FACTOR_OTHER;
1239: }
1240: }
1241: }
1242: if (!mumps->myid && mumps->id.ICNTL(16) > 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB," mumps->id.ICNTL(16):=%d\n",mumps->id.INFOG(16));
1244: (F)->assembled = PETSC_TRUE;
1245: mumps->matstruc = SAME_NONZERO_PATTERN;
1246: mumps->schur_factored = PETSC_FALSE;
1247: mumps->schur_inverted = PETSC_FALSE;
1249: /* just to be sure that ICNTL(19) value returned by a call from MatMumpsGetIcntl is always consistent */
1250: if (!mumps->sym && mumps->id.ICNTL(19) && mumps->id.ICNTL(19) != 1) mumps->id.ICNTL(19) = 3;
1252: if (mumps->size > 1) {
1253: PetscInt lsol_loc;
1254: PetscScalar *sol_loc;
1256: PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);
1257: if (isMPIAIJ) F_diag = ((Mat_MPIAIJ*)(F)->data)->A;
1258: else F_diag = ((Mat_MPISBAIJ*)(F)->data)->A;
1259: F_diag->assembled = PETSC_TRUE;
1261: /* distributed solution; Create x_seq=sol_loc for repeated use */
1262: if (mumps->x_seq) {
1263: VecScatterDestroy(&mumps->scat_sol);
1264: PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);
1265: VecDestroy(&mumps->x_seq);
1266: }
1267: lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
1268: PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);
1269: mumps->id.lsol_loc = lsol_loc;
1270: mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1271: VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);
1272: }
1273: return(0);
1274: }
1276: /* Sets MUMPS options from the options database */
1279: PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
1280: {
1281: Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr;
1283: PetscInt icntl,info[40],i,ninfo=40;
1284: PetscBool flg;
1287: PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");
1288: PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);
1289: if (flg) mumps->id.ICNTL(1) = icntl;
1290: PetscOptionsInt("-mat_mumps_icntl_2","ICNTL(2): output stream for diagnostic printing, statistics, and warning","None",mumps->id.ICNTL(2),&icntl,&flg);
1291: if (flg) mumps->id.ICNTL(2) = icntl;
1292: PetscOptionsInt("-mat_mumps_icntl_3","ICNTL(3): output stream for global information, collected on the host","None",mumps->id.ICNTL(3),&icntl,&flg);
1293: if (flg) mumps->id.ICNTL(3) = icntl;
1295: PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);
1296: if (flg) mumps->id.ICNTL(4) = icntl;
1297: if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
1299: PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)","None",mumps->id.ICNTL(6),&icntl,&flg);
1300: if (flg) mumps->id.ICNTL(6) = icntl;
1302: PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis","None",mumps->id.ICNTL(7),&icntl,&flg);
1303: if (flg) {
1304: if (icntl== 1 && mumps->size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
1305: else mumps->id.ICNTL(7) = icntl;
1306: }
1308: PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),NULL);
1309: /* PetscOptionsInt("-mat_mumps_icntl_9","ICNTL(9): computes the solution using A or A^T","None",mumps->id.ICNTL(9),&mumps->id.ICNTL(9),NULL); handled by MatSolveTranspose_MUMPS() */
1310: PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);
1311: PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to an error analysis (via -ksp_view)","None",mumps->id.ICNTL(11),&mumps->id.ICNTL(11),NULL);
1312: PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)","None",mumps->id.ICNTL(12),&mumps->id.ICNTL(12),NULL);
1313: PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting","None",mumps->id.ICNTL(13),&mumps->id.ICNTL(13),NULL);
1314: PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage increase in the estimated working space","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),NULL);
1315: PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);
1316: if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */
1317: MatMumpsResetSchur_Private(mumps);
1318: }
1319: /* PetscOptionsInt("-mat_mumps_icntl_20","ICNTL(20): the format (dense or sparse) of the right-hand sides","None",mumps->id.ICNTL(20),&mumps->id.ICNTL(20),NULL); -- sparse rhs is not supported in PETSc API */
1320: /* PetscOptionsInt("-mat_mumps_icntl_21","ICNTL(21): the distribution (centralized or distributed) of the solution vectors","None",mumps->id.ICNTL(21),&mumps->id.ICNTL(21),NULL); we only use distributed solution vector */
1322: PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)","None",mumps->id.ICNTL(22),&mumps->id.ICNTL(22),NULL);
1323: PetscOptionsInt("-mat_mumps_icntl_23","ICNTL(23): max size of the working memory (MB) that can allocate per processor","None",mumps->id.ICNTL(23),&mumps->id.ICNTL(23),NULL);
1324: PetscOptionsInt("-mat_mumps_icntl_24","ICNTL(24): detection of null pivot rows (0 or 1)","None",mumps->id.ICNTL(24),&mumps->id.ICNTL(24),NULL);
1325: if (mumps->id.ICNTL(24)) {
1326: mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
1327: }
1329: PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): compute a solution of a deficient matrix and a null space basis","None",mumps->id.ICNTL(25),&mumps->id.ICNTL(25),NULL);
1330: PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): drives the solution phase if a Schur complement matrix","None",mumps->id.ICNTL(26),&mumps->id.ICNTL(26),NULL);
1331: PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): the blocking size for multiple right-hand sides","None",mumps->id.ICNTL(27),&mumps->id.ICNTL(27),NULL);
1332: PetscOptionsInt("-mat_mumps_icntl_28","ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering","None",mumps->id.ICNTL(28),&mumps->id.ICNTL(28),NULL);
1333: PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),NULL);
1334: PetscOptionsInt("-mat_mumps_icntl_30","ICNTL(30): compute user-specified set of entries in inv(A)","None",mumps->id.ICNTL(30),&mumps->id.ICNTL(30),NULL);
1335: PetscOptionsInt("-mat_mumps_icntl_31","ICNTL(31): indicates which factors may be discarded during factorization","None",mumps->id.ICNTL(31),&mumps->id.ICNTL(31),NULL);
1336: /* PetscOptionsInt("-mat_mumps_icntl_32","ICNTL(32): performs the forward elemination of the right-hand sides during factorization","None",mumps->id.ICNTL(32),&mumps->id.ICNTL(32),NULL); -- not supported by PETSc API */
1337: PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);
1339: PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);
1340: PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL);
1341: PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);
1342: PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL);
1343: PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL);
1345: PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL);
1347: PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);
1348: if (ninfo) {
1349: if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo);
1350: PetscMalloc1(ninfo,&mumps->info);
1351: mumps->ninfo = ninfo;
1352: for (i=0; i<ninfo; i++) {
1353: if (info[i] < 0 || info[i]>40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index of INFO %d must between 1 and 40\n",ninfo);
1354: else {
1355: mumps->info[i] = info[i];
1356: }
1357: }
1358: }
1360: PetscOptionsEnd();
1361: return(0);
1362: }
1366: PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
1367: {
1371: MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);
1372: MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);
1373: MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));
1375: mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);
1377: mumps->id.job = JOB_INIT;
1378: mumps->id.par = 1; /* host participates factorizaton and solve */
1379: mumps->id.sym = mumps->sym;
1380: PetscMUMPS_c(&mumps->id);
1382: mumps->scat_rhs = NULL;
1383: mumps->scat_sol = NULL;
1385: /* set PETSc-MUMPS default options - override MUMPS default */
1386: mumps->id.ICNTL(3) = 0;
1387: mumps->id.ICNTL(4) = 0;
1388: if (mumps->size == 1) {
1389: mumps->id.ICNTL(18) = 0; /* centralized assembled matrix input */
1390: } else {
1391: mumps->id.ICNTL(18) = 3; /* distributed assembled matrix input */
1392: mumps->id.ICNTL(20) = 0; /* rhs is in dense format */
1393: mumps->id.ICNTL(21) = 1; /* distributed solution */
1394: }
1396: /* schur */
1397: mumps->id.size_schur = 0;
1398: mumps->id.listvar_schur = NULL;
1399: mumps->id.schur = NULL;
1400: mumps->sizeredrhs = 0;
1401: mumps->schur_pivots = NULL;
1402: mumps->schur_work = NULL;
1403: mumps->schur_sol = NULL;
1404: mumps->schur_sizesol = 0;
1405: mumps->schur_factored = PETSC_FALSE;
1406: mumps->schur_inverted = PETSC_FALSE;
1407: return(0);
1408: }
1412: PetscErrorCode MatFactorSymbolic_MUMPS_ReportIfError(Mat F,Mat A,const MatFactorInfo *info,Mat_MUMPS *mumps)
1413: {
1417: if (mumps->id.INFOG(1) < 0) {
1418: if (A->erroriffailure) {
1419: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
1420: } else {
1421: if (mumps->id.INFOG(1) == -6) {
1422: PetscInfo2(F,"matrix is singular in structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1423: F->errortype = MAT_FACTOR_STRUCT_ZEROPIVOT;
1424: } else if (mumps->id.INFOG(1) == -5 || mumps->id.INFOG(1) == -7) {
1425: PetscInfo2(F,"problem of workspace, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1426: F->errortype = MAT_FACTOR_OUTMEMORY;
1427: } else {
1428: PetscInfo2(F,"Error reported by MUMPS in analysis phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1429: F->errortype = MAT_FACTOR_OTHER;
1430: }
1431: }
1432: }
1433: return(0);
1434: }
1436: /* Note Petsc r(=c) permutation is used when mumps->id.ICNTL(7)==1 with centralized assembled matrix input; otherwise r and c are ignored */
1439: PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1440: {
1441: Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr;
1443: Vec b;
1444: IS is_iden;
1445: const PetscInt M = A->rmap->N;
1448: mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1450: /* Set MUMPS options from the options database */
1451: PetscSetMUMPSFromOptions(F,A);
1453: (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);
1455: /* analysis phase */
1456: /*----------------*/
1457: mumps->id.job = JOB_FACTSYMBOLIC;
1458: mumps->id.n = M;
1459: switch (mumps->id.ICNTL(18)) {
1460: case 0: /* centralized assembled matrix input */
1461: if (!mumps->myid) {
1462: mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1463: if (mumps->id.ICNTL(6)>1) {
1464: mumps->id.a = (MumpsScalar*)mumps->val;
1465: }
1466: if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
1467: /*
1468: PetscBool flag;
1469: ISEqual(r,c,&flag);
1470: if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
1471: ISView(r,PETSC_VIEWER_STDOUT_SELF);
1472: */
1473: if (!mumps->myid) {
1474: const PetscInt *idx;
1475: PetscInt i,*perm_in;
1477: PetscMalloc1(M,&perm_in);
1478: ISGetIndices(r,&idx);
1480: mumps->id.perm_in = perm_in;
1481: for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
1482: ISRestoreIndices(r,&idx);
1483: }
1484: }
1485: }
1486: break;
1487: case 3: /* distributed assembled matrix input (size>1) */
1488: mumps->id.nz_loc = mumps->nz;
1489: mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1490: if (mumps->id.ICNTL(6)>1) {
1491: mumps->id.a_loc = (MumpsScalar*)mumps->val;
1492: }
1493: /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1494: if (!mumps->myid) {
1495: VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);
1496: ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);
1497: } else {
1498: VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);
1499: ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
1500: }
1501: MatCreateVecs(A,NULL,&b);
1502: VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);
1503: ISDestroy(&is_iden);
1504: VecDestroy(&b);
1505: break;
1506: }
1507: PetscMUMPS_c(&mumps->id);
1508: MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);
1510: F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1511: F->ops->solve = MatSolve_MUMPS;
1512: F->ops->solvetranspose = MatSolveTranspose_MUMPS;
1513: F->ops->matsolve = MatMatSolve_MUMPS;
1514: return(0);
1515: }
1517: /* Note the Petsc r and c permutations are ignored */
1520: PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1521: {
1522: Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr;
1524: Vec b;
1525: IS is_iden;
1526: const PetscInt M = A->rmap->N;
1529: mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1531: /* Set MUMPS options from the options database */
1532: PetscSetMUMPSFromOptions(F,A);
1534: (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);
1536: /* analysis phase */
1537: /*----------------*/
1538: mumps->id.job = JOB_FACTSYMBOLIC;
1539: mumps->id.n = M;
1540: switch (mumps->id.ICNTL(18)) {
1541: case 0: /* centralized assembled matrix input */
1542: if (!mumps->myid) {
1543: mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1544: if (mumps->id.ICNTL(6)>1) {
1545: mumps->id.a = (MumpsScalar*)mumps->val;
1546: }
1547: }
1548: break;
1549: case 3: /* distributed assembled matrix input (size>1) */
1550: mumps->id.nz_loc = mumps->nz;
1551: mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1552: if (mumps->id.ICNTL(6)>1) {
1553: mumps->id.a_loc = (MumpsScalar*)mumps->val;
1554: }
1555: /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1556: if (!mumps->myid) {
1557: VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);
1558: ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);
1559: } else {
1560: VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);
1561: ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
1562: }
1563: MatCreateVecs(A,NULL,&b);
1564: VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);
1565: ISDestroy(&is_iden);
1566: VecDestroy(&b);
1567: break;
1568: }
1569: PetscMUMPS_c(&mumps->id);
1570: MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);
1572: F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1573: F->ops->solve = MatSolve_MUMPS;
1574: F->ops->solvetranspose = MatSolveTranspose_MUMPS;
1575: return(0);
1576: }
1578: /* Note the Petsc r permutation and factor info are ignored */
1581: PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1582: {
1583: Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr;
1585: Vec b;
1586: IS is_iden;
1587: const PetscInt M = A->rmap->N;
1590: mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1592: /* Set MUMPS options from the options database */
1593: PetscSetMUMPSFromOptions(F,A);
1595: (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);
1597: /* analysis phase */
1598: /*----------------*/
1599: mumps->id.job = JOB_FACTSYMBOLIC;
1600: mumps->id.n = M;
1601: switch (mumps->id.ICNTL(18)) {
1602: case 0: /* centralized assembled matrix input */
1603: if (!mumps->myid) {
1604: mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1605: if (mumps->id.ICNTL(6)>1) {
1606: mumps->id.a = (MumpsScalar*)mumps->val;
1607: }
1608: }
1609: break;
1610: case 3: /* distributed assembled matrix input (size>1) */
1611: mumps->id.nz_loc = mumps->nz;
1612: mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1613: if (mumps->id.ICNTL(6)>1) {
1614: mumps->id.a_loc = (MumpsScalar*)mumps->val;
1615: }
1616: /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1617: if (!mumps->myid) {
1618: VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);
1619: ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);
1620: } else {
1621: VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);
1622: ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);
1623: }
1624: MatCreateVecs(A,NULL,&b);
1625: VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);
1626: ISDestroy(&is_iden);
1627: VecDestroy(&b);
1628: break;
1629: }
1630: PetscMUMPS_c(&mumps->id);
1631: MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);
1634: F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1635: F->ops->solve = MatSolve_MUMPS;
1636: F->ops->solvetranspose = MatSolve_MUMPS;
1637: F->ops->matsolve = MatMatSolve_MUMPS;
1638: #if defined(PETSC_USE_COMPLEX)
1639: F->ops->getinertia = NULL;
1640: #else
1641: F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
1642: #endif
1643: return(0);
1644: }
1648: PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
1649: {
1650: PetscErrorCode ierr;
1651: PetscBool iascii;
1652: PetscViewerFormat format;
1653: Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
1656: /* check if matrix is mumps type */
1657: if (A->ops->solve != MatSolve_MUMPS) return(0);
1659: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1660: if (iascii) {
1661: PetscViewerGetFormat(viewer,&format);
1662: if (format == PETSC_VIEWER_ASCII_INFO) {
1663: PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");
1664: PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",mumps->id.sym);
1665: PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",mumps->id.par);
1666: PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",mumps->id.ICNTL(1));
1667: PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));
1668: PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",mumps->id.ICNTL(3));
1669: PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",mumps->id.ICNTL(4));
1670: PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",mumps->id.ICNTL(5));
1671: PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",mumps->id.ICNTL(6));
1672: PetscViewerASCIIPrintf(viewer," ICNTL(7) (sequential matrix ordering):%d \n",mumps->id.ICNTL(7));
1673: PetscViewerASCIIPrintf(viewer," ICNTL(8) (scaling strategy): %d \n",mumps->id.ICNTL(8));
1674: PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",mumps->id.ICNTL(10));
1675: PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",mumps->id.ICNTL(11));
1676: if (mumps->id.ICNTL(11)>0) {
1677: PetscViewerASCIIPrintf(viewer," RINFOG(4) (inf norm of input mat): %g\n",mumps->id.RINFOG(4));
1678: PetscViewerASCIIPrintf(viewer," RINFOG(5) (inf norm of solution): %g\n",mumps->id.RINFOG(5));
1679: PetscViewerASCIIPrintf(viewer," RINFOG(6) (inf norm of residual): %g\n",mumps->id.RINFOG(6));
1680: PetscViewerASCIIPrintf(viewer," RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));
1681: PetscViewerASCIIPrintf(viewer," RINFOG(9) (error estimate): %g \n",mumps->id.RINFOG(9));
1682: PetscViewerASCIIPrintf(viewer," RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));
1683: }
1684: PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",mumps->id.ICNTL(12));
1685: PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",mumps->id.ICNTL(13));
1686: PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));
1687: /* ICNTL(15-17) not used */
1688: PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",mumps->id.ICNTL(18));
1689: PetscViewerASCIIPrintf(viewer," ICNTL(19) (Schur complement info): %d \n",mumps->id.ICNTL(19));
1690: PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",mumps->id.ICNTL(20));
1691: PetscViewerASCIIPrintf(viewer," ICNTL(21) (solution struct): %d \n",mumps->id.ICNTL(21));
1692: PetscViewerASCIIPrintf(viewer," ICNTL(22) (in-core/out-of-core facility): %d \n",mumps->id.ICNTL(22));
1693: PetscViewerASCIIPrintf(viewer," ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));
1695: PetscViewerASCIIPrintf(viewer," ICNTL(24) (detection of null pivot rows): %d \n",mumps->id.ICNTL(24));
1696: PetscViewerASCIIPrintf(viewer," ICNTL(25) (computation of a null space basis): %d \n",mumps->id.ICNTL(25));
1697: PetscViewerASCIIPrintf(viewer," ICNTL(26) (Schur options for rhs or solution): %d \n",mumps->id.ICNTL(26));
1698: PetscViewerASCIIPrintf(viewer," ICNTL(27) (experimental parameter): %d \n",mumps->id.ICNTL(27));
1699: PetscViewerASCIIPrintf(viewer," ICNTL(28) (use parallel or sequential ordering): %d \n",mumps->id.ICNTL(28));
1700: PetscViewerASCIIPrintf(viewer," ICNTL(29) (parallel ordering): %d \n",mumps->id.ICNTL(29));
1702: PetscViewerASCIIPrintf(viewer," ICNTL(30) (user-specified set of entries in inv(A)): %d \n",mumps->id.ICNTL(30));
1703: PetscViewerASCIIPrintf(viewer," ICNTL(31) (factors is discarded in the solve phase): %d \n",mumps->id.ICNTL(31));
1704: PetscViewerASCIIPrintf(viewer," ICNTL(33) (compute determinant): %d \n",mumps->id.ICNTL(33));
1706: PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",mumps->id.CNTL(1));
1707: PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));
1708: PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",mumps->id.CNTL(3));
1709: PetscViewerASCIIPrintf(viewer," CNTL(4) (value of static pivoting): %g \n",mumps->id.CNTL(4));
1710: PetscViewerASCIIPrintf(viewer," CNTL(5) (fixation for null pivots): %g \n",mumps->id.CNTL(5));
1712: /* infomation local to each processor */
1713: PetscViewerASCIIPrintf(viewer, " RINFO(1) (local estimated flops for the elimination after analysis): \n");
1714: PetscViewerASCIIPushSynchronized(viewer);
1715: PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(1));
1716: PetscViewerFlush(viewer);
1717: PetscViewerASCIIPrintf(viewer, " RINFO(2) (local estimated flops for the assembly after factorization): \n");
1718: PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(2));
1719: PetscViewerFlush(viewer);
1720: PetscViewerASCIIPrintf(viewer, " RINFO(3) (local estimated flops for the elimination after factorization): \n");
1721: PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(3));
1722: PetscViewerFlush(viewer);
1724: PetscViewerASCIIPrintf(viewer, " INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");
1725: PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(15));
1726: PetscViewerFlush(viewer);
1728: PetscViewerASCIIPrintf(viewer, " INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");
1729: PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(16));
1730: PetscViewerFlush(viewer);
1732: PetscViewerASCIIPrintf(viewer, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");
1733: PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(23));
1734: PetscViewerFlush(viewer);
1736: if (mumps->ninfo && mumps->ninfo <= 40){
1737: PetscInt i;
1738: for (i=0; i<mumps->ninfo; i++){
1739: PetscViewerASCIIPrintf(viewer, " INFO(%d): \n",mumps->info[i]);
1740: PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));
1741: PetscViewerFlush(viewer);
1742: }
1743: }
1746: PetscViewerASCIIPopSynchronized(viewer);
1748: if (!mumps->myid) { /* information from the host */
1749: PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));
1750: PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));
1751: PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));
1752: PetscViewerASCIIPrintf(viewer," (RINFOG(12) RINFOG(13))*2^INFOG(34) (determinant): (%g,%g)*(2^%d)\n",mumps->id.RINFOG(12),mumps->id.RINFOG(13),mumps->id.INFOG(34));
1754: PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));
1755: PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));
1756: PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));
1757: PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));
1758: PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));
1759: PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));
1760: PetscViewerASCIIPrintf(viewer," INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",mumps->id.INFOG(9));
1761: PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));
1762: PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));
1763: PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));
1764: PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));
1765: PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));
1766: PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));
1767: 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",mumps->id.INFOG(16));
1768: PetscViewerASCIIPrintf(viewer," INFOG(17) (estimated size of all MUMPS internal data for factorization after analysis: sum over all processors): %d \n",mumps->id.INFOG(17));
1769: PetscViewerASCIIPrintf(viewer," INFOG(18) (size of all MUMPS internal data allocated during factorization: value on the most memory consuming processor): %d \n",mumps->id.INFOG(18));
1770: PetscViewerASCIIPrintf(viewer," INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",mumps->id.INFOG(19));
1771: PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));
1772: PetscViewerASCIIPrintf(viewer," INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",mumps->id.INFOG(21));
1773: PetscViewerASCIIPrintf(viewer," INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",mumps->id.INFOG(22));
1774: PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));
1775: PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));
1776: PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));
1777: PetscViewerASCIIPrintf(viewer," INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));
1778: PetscViewerASCIIPrintf(viewer," INFOG(29) (after factorization: effective number of entries in the factors (sum over all processors)): %d\n",mumps->id.INFOG(29));
1779: PetscViewerASCIIPrintf(viewer," INFOG(30, 31) (after solution: size in Mbytes of memory used during solution phase): %d, %d\n",mumps->id.INFOG(30),mumps->id.INFOG(31));
1780: PetscViewerASCIIPrintf(viewer," INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));
1781: PetscViewerASCIIPrintf(viewer," INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));
1782: PetscViewerASCIIPrintf(viewer," INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));
1783: }
1784: }
1785: }
1786: return(0);
1787: }
1791: PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
1792: {
1793: Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr;
1796: info->block_size = 1.0;
1797: info->nz_allocated = mumps->id.INFOG(20);
1798: info->nz_used = mumps->id.INFOG(20);
1799: info->nz_unneeded = 0.0;
1800: info->assemblies = 0.0;
1801: info->mallocs = 0.0;
1802: info->memory = 0.0;
1803: info->fill_ratio_given = 0;
1804: info->fill_ratio_needed = 0;
1805: info->factor_mallocs = 0;
1806: return(0);
1807: }
1809: /* -------------------------------------------------------------------------------------------*/
1812: PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is)
1813: {
1814: Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1815: const PetscInt *idxs;
1816: PetscInt size,i;
1820: if (mumps->size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS parallel Schur complements not yet supported from PETSc\n");
1821: ISGetLocalSize(is,&size);
1822: if (mumps->id.size_schur != size) {
1823: PetscFree2(mumps->id.listvar_schur,mumps->id.schur);
1824: mumps->id.size_schur = size;
1825: mumps->id.schur_lld = size;
1826: PetscMalloc2(size,&mumps->id.listvar_schur,size*size,&mumps->id.schur);
1827: }
1828: ISGetIndices(is,&idxs);
1829: PetscMemcpy(mumps->id.listvar_schur,idxs,size*sizeof(PetscInt));
1830: /* MUMPS expects Fortran style indices */
1831: for (i=0;i<size;i++) mumps->id.listvar_schur[i]++;
1832: ISRestoreIndices(is,&idxs);
1833: if (size) { /* turn on Schur switch if we the set of indices is not empty */
1834: if (F->factortype == MAT_FACTOR_LU) {
1835: mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */
1836: } else {
1837: mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */
1838: }
1839: /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */
1840: mumps->id.ICNTL(26) = -1;
1841: }
1842: return(0);
1843: }
1845: /* -------------------------------------------------------------------------------------------*/
1848: PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S)
1849: {
1850: Mat St;
1851: Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1852: PetscScalar *array;
1853: #if defined(PETSC_USE_COMPLEX)
1854: PetscScalar im = PetscSqrtScalar((PetscScalar)-1.0);
1855: #endif
1859: if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
1860: else if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
1862: MatCreate(PetscObjectComm((PetscObject)F),&St);
1863: MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);
1864: MatSetType(St,MATDENSE);
1865: MatSetUp(St);
1866: MatDenseGetArray(St,&array);
1867: if (!mumps->sym) { /* MUMPS always return a full matrix */
1868: if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1869: PetscInt i,j,N=mumps->id.size_schur;
1870: for (i=0;i<N;i++) {
1871: for (j=0;j<N;j++) {
1872: #if !defined(PETSC_USE_COMPLEX)
1873: PetscScalar val = mumps->id.schur[i*N+j];
1874: #else
1875: PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1876: #endif
1877: array[j*N+i] = val;
1878: }
1879: }
1880: } else { /* stored by columns */
1881: PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));
1882: }
1883: } else { /* either full or lower-triangular (not packed) */
1884: if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */
1885: PetscInt i,j,N=mumps->id.size_schur;
1886: for (i=0;i<N;i++) {
1887: for (j=i;j<N;j++) {
1888: #if !defined(PETSC_USE_COMPLEX)
1889: PetscScalar val = mumps->id.schur[i*N+j];
1890: #else
1891: PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1892: #endif
1893: array[i*N+j] = val;
1894: array[j*N+i] = val;
1895: }
1896: }
1897: } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */
1898: PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));
1899: } else { /* ICNTL(19) == 1 lower triangular stored by rows */
1900: PetscInt i,j,N=mumps->id.size_schur;
1901: for (i=0;i<N;i++) {
1902: for (j=0;j<i+1;j++) {
1903: #if !defined(PETSC_USE_COMPLEX)
1904: PetscScalar val = mumps->id.schur[i*N+j];
1905: #else
1906: PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
1907: #endif
1908: array[i*N+j] = val;
1909: array[j*N+i] = val;
1910: }
1911: }
1912: }
1913: }
1914: MatDenseRestoreArray(St,&array);
1915: *S = St;
1916: return(0);
1917: }
1919: /* -------------------------------------------------------------------------------------------*/
1922: PetscErrorCode MatFactorGetSchurComplement_MUMPS(Mat F,Mat* S)
1923: {
1924: Mat St;
1925: Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1929: if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
1930: else if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
1932: /* It should be the responsibility of the user to handle different ICNTL(19) cases and factorization stages if they want to work with the raw data */
1933: MatCreateSeqDense(PetscObjectComm((PetscObject)F),mumps->id.size_schur,mumps->id.size_schur,(PetscScalar*)mumps->id.schur,&St);
1934: *S = St;
1935: return(0);
1936: }
1938: /* -------------------------------------------------------------------------------------------*/
1941: PetscErrorCode MatFactorInvertSchurComplement_MUMPS(Mat F)
1942: {
1943: Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1947: if (!mumps->id.ICNTL(19)) { /* do nothing */
1948: return(0);
1949: }
1950: if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
1951: MatMumpsInvertSchur_Private(mumps);
1952: return(0);
1953: }
1955: /* -------------------------------------------------------------------------------------------*/
1958: PetscErrorCode MatFactorSolveSchurComplement_MUMPS(Mat F, Vec rhs, Vec sol)
1959: {
1960: Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
1961: MumpsScalar *orhs;
1962: PetscScalar *osol,*nrhs,*nsol;
1963: PetscInt orhs_size,osol_size,olrhs_size;
1967: if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
1968: if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
1970: /* swap pointers */
1971: orhs = mumps->id.redrhs;
1972: olrhs_size = mumps->id.lredrhs;
1973: orhs_size = mumps->sizeredrhs;
1974: osol = mumps->schur_sol;
1975: osol_size = mumps->schur_sizesol;
1976: VecGetArray(rhs,&nrhs);
1977: VecGetArray(sol,&nsol);
1978: mumps->id.redrhs = (MumpsScalar*)nrhs;
1979: VecGetLocalSize(rhs,&mumps->sizeredrhs);
1980: mumps->id.lredrhs = mumps->sizeredrhs;
1981: mumps->schur_sol = nsol;
1982: VecGetLocalSize(sol,&mumps->schur_sizesol);
1984: /* solve Schur complement */
1985: mumps->id.nrhs = 1;
1986: MatMumpsSolveSchur_Private(mumps,PETSC_FALSE);
1987: /* restore pointers */
1988: VecRestoreArray(rhs,&nrhs);
1989: VecRestoreArray(sol,&nsol);
1990: mumps->id.redrhs = orhs;
1991: mumps->id.lredrhs = olrhs_size;
1992: mumps->sizeredrhs = orhs_size;
1993: mumps->schur_sol = osol;
1994: mumps->schur_sizesol = osol_size;
1995: return(0);
1996: }
1998: /* -------------------------------------------------------------------------------------------*/
2001: PetscErrorCode MatFactorSolveSchurComplementTranspose_MUMPS(Mat F, Vec rhs, Vec sol)
2002: {
2003: Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2004: MumpsScalar *orhs;
2005: PetscScalar *osol,*nrhs,*nsol;
2006: PetscInt orhs_size,osol_size;
2010: if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
2011: else if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
2013: /* swap pointers */
2014: orhs = mumps->id.redrhs;
2015: orhs_size = mumps->sizeredrhs;
2016: osol = mumps->schur_sol;
2017: osol_size = mumps->schur_sizesol;
2018: VecGetArray(rhs,&nrhs);
2019: VecGetArray(sol,&nsol);
2020: mumps->id.redrhs = (MumpsScalar*)nrhs;
2021: VecGetLocalSize(rhs,&mumps->sizeredrhs);
2022: mumps->schur_sol = nsol;
2023: VecGetLocalSize(sol,&mumps->schur_sizesol);
2025: /* solve Schur complement */
2026: mumps->id.nrhs = 1;
2027: mumps->id.ICNTL(9) = 0;
2028: MatMumpsSolveSchur_Private(mumps,PETSC_FALSE);
2029: mumps->id.ICNTL(9) = 1;
2030: /* restore pointers */
2031: VecRestoreArray(rhs,&nrhs);
2032: VecRestoreArray(sol,&nsol);
2033: mumps->id.redrhs = orhs;
2034: mumps->sizeredrhs = orhs_size;
2035: mumps->schur_sol = osol;
2036: mumps->schur_sizesol = osol_size;
2037: return(0);
2038: }
2040: /* -------------------------------------------------------------------------------------------*/
2043: PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
2044: {
2045: Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2048: mumps->id.ICNTL(icntl) = ival;
2049: return(0);
2050: }
2054: PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
2055: {
2056: Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2059: *ival = mumps->id.ICNTL(icntl);
2060: return(0);
2061: }
2065: /*@
2066: MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
2068: Logically Collective on Mat
2070: Input Parameters:
2071: + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2072: . icntl - index of MUMPS parameter array ICNTL()
2073: - ival - value of MUMPS ICNTL(icntl)
2075: Options Database:
2076: . -mat_mumps_icntl_<icntl> <ival>
2078: Level: beginner
2080: References:
2081: . MUMPS Users' Guide
2083: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2084: @*/
2085: PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
2086: {
2088:
2091: if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2094: PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));
2095: return(0);
2096: }
2100: /*@
2101: MatMumpsGetIcntl - Get MUMPS parameter ICNTL()
2103: Logically Collective on Mat
2105: Input Parameters:
2106: + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2107: - icntl - index of MUMPS parameter array ICNTL()
2109: Output Parameter:
2110: . ival - value of MUMPS ICNTL(icntl)
2112: Level: beginner
2114: References:
2115: . MUMPS Users' Guide
2117: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2118: @*/
2119: PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
2120: {
2125: if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2128: PetscUseMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2129: return(0);
2130: }
2132: /* -------------------------------------------------------------------------------------------*/
2135: PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
2136: {
2137: Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2140: mumps->id.CNTL(icntl) = val;
2141: return(0);
2142: }
2146: PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
2147: {
2148: Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2151: *val = mumps->id.CNTL(icntl);
2152: return(0);
2153: }
2157: /*@
2158: MatMumpsSetCntl - Set MUMPS parameter CNTL()
2160: Logically Collective on Mat
2162: Input Parameters:
2163: + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2164: . icntl - index of MUMPS parameter array CNTL()
2165: - val - value of MUMPS CNTL(icntl)
2167: Options Database:
2168: . -mat_mumps_cntl_<icntl> <val>
2170: Level: beginner
2172: References:
2173: . MUMPS Users' Guide
2175: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2176: @*/
2177: PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
2178: {
2183: if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2186: PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));
2187: return(0);
2188: }
2192: /*@
2193: MatMumpsGetCntl - Get MUMPS parameter CNTL()
2195: Logically Collective on Mat
2197: Input Parameters:
2198: + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2199: - icntl - index of MUMPS parameter array CNTL()
2201: Output Parameter:
2202: . val - value of MUMPS CNTL(icntl)
2204: Level: beginner
2206: References:
2207: . MUMPS Users' Guide
2209: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2210: @*/
2211: PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
2212: {
2217: if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2220: PetscUseMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2221: return(0);
2222: }
2226: PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
2227: {
2228: Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2231: *info = mumps->id.INFO(icntl);
2232: return(0);
2233: }
2237: PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
2238: {
2239: Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2242: *infog = mumps->id.INFOG(icntl);
2243: return(0);
2244: }
2248: PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
2249: {
2250: Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2253: *rinfo = mumps->id.RINFO(icntl);
2254: return(0);
2255: }
2259: PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
2260: {
2261: Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
2264: *rinfog = mumps->id.RINFOG(icntl);
2265: return(0);
2266: }
2270: /*@
2271: MatMumpsGetInfo - Get MUMPS parameter INFO()
2273: Logically Collective on Mat
2275: Input Parameters:
2276: + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2277: - icntl - index of MUMPS parameter array INFO()
2279: Output Parameter:
2280: . ival - value of MUMPS INFO(icntl)
2282: Level: beginner
2284: References:
2285: . MUMPS Users' Guide
2287: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2288: @*/
2289: PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
2290: {
2295: if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2297: PetscUseMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2298: return(0);
2299: }
2303: /*@
2304: MatMumpsGetInfog - Get MUMPS parameter INFOG()
2306: Logically Collective on Mat
2308: Input Parameters:
2309: + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2310: - icntl - index of MUMPS parameter array INFOG()
2312: Output Parameter:
2313: . ival - value of MUMPS INFOG(icntl)
2315: Level: beginner
2317: References:
2318: . MUMPS Users' Guide
2320: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2321: @*/
2322: PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
2323: {
2328: if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2330: PetscUseMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));
2331: return(0);
2332: }
2336: /*@
2337: MatMumpsGetRinfo - Get MUMPS parameter RINFO()
2339: Logically Collective on Mat
2341: Input Parameters:
2342: + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2343: - icntl - index of MUMPS parameter array RINFO()
2345: Output Parameter:
2346: . val - value of MUMPS RINFO(icntl)
2348: Level: beginner
2350: References:
2351: . MUMPS Users' Guide
2353: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2354: @*/
2355: PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
2356: {
2361: if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2363: PetscUseMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2364: return(0);
2365: }
2369: /*@
2370: MatMumpsGetRinfog - Get MUMPS parameter RINFOG()
2372: Logically Collective on Mat
2374: Input Parameters:
2375: + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2376: - icntl - index of MUMPS parameter array RINFOG()
2378: Output Parameter:
2379: . val - value of MUMPS RINFOG(icntl)
2381: Level: beginner
2383: References:
2384: . MUMPS Users' Guide
2386: .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2387: @*/
2388: PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
2389: {
2394: if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2396: PetscUseMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));
2397: return(0);
2398: }
2400: /*MC
2401: MATSOLVERMUMPS - A matrix type providing direct solvers (LU and Cholesky) for
2402: distributed and sequential matrices via the external package MUMPS.
2404: Works with MATAIJ and MATSBAIJ matrices
2406: Use ./configure --download-mumps --download-scalapack --download-parmetis --download-metis --download-ptscotch to have PETSc installed with MUMPS
2408: Use -pc_type cholesky or lu -pc_factor_mat_solver_package mumps to us this direct solver
2410: Options Database Keys:
2411: + -mat_mumps_icntl_1 - ICNTL(1): output stream for error messages
2412: . -mat_mumps_icntl_2 - ICNTL(2): output stream for diagnostic printing, statistics, and warning
2413: . -mat_mumps_icntl_3 - ICNTL(3): output stream for global information, collected on the host
2414: . -mat_mumps_icntl_4 - ICNTL(4): level of printing (0 to 4)
2415: . -mat_mumps_icntl_6 - ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)
2416: . -mat_mumps_icntl_7 - ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis
2417: . -mat_mumps_icntl_8 - ICNTL(8): scaling strategy (-2 to 8 or 77)
2418: . -mat_mumps_icntl_10 - ICNTL(10): max num of refinements
2419: . -mat_mumps_icntl_11 - ICNTL(11): statistics related to an error analysis (via -ksp_view)
2420: . -mat_mumps_icntl_12 - ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)
2421: . -mat_mumps_icntl_13 - ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting
2422: . -mat_mumps_icntl_14 - ICNTL(14): percentage increase in the estimated working space
2423: . -mat_mumps_icntl_19 - ICNTL(19): computes the Schur complement
2424: . -mat_mumps_icntl_22 - ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)
2425: . -mat_mumps_icntl_23 - ICNTL(23): max size of the working memory (MB) that can allocate per processor
2426: . -mat_mumps_icntl_24 - ICNTL(24): detection of null pivot rows (0 or 1)
2427: . -mat_mumps_icntl_25 - ICNTL(25): compute a solution of a deficient matrix and a null space basis
2428: . -mat_mumps_icntl_26 - ICNTL(26): drives the solution phase if a Schur complement matrix
2429: . -mat_mumps_icntl_28 - ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering
2430: . -mat_mumps_icntl_29 - ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis
2431: . -mat_mumps_icntl_30 - ICNTL(30): compute user-specified set of entries in inv(A)
2432: . -mat_mumps_icntl_31 - ICNTL(31): indicates which factors may be discarded during factorization
2433: . -mat_mumps_icntl_33 - ICNTL(33): compute determinant
2434: . -mat_mumps_cntl_1 - CNTL(1): relative pivoting threshold
2435: . -mat_mumps_cntl_2 - CNTL(2): stopping criterion of refinement
2436: . -mat_mumps_cntl_3 - CNTL(3): absolute pivoting threshold
2437: . -mat_mumps_cntl_4 - CNTL(4): value for static pivoting
2438: - -mat_mumps_cntl_5 - CNTL(5): fixation for null pivots
2440: Level: beginner
2442: Notes: When a MUMPS factorization fails inside a KSP solve, for example with a KSP_DIVERGED_PCSETUP_FAILED, one can find the MUMPS information about the failure by calling
2443: $ KSPGetPC(ksp,&pc);
2444: $ PCFactorGetMatrix(pc,&mat);
2445: $ MatMumpsGetInfo(mat,....);
2446: $ MatMumpsGetInfog(mat,....); etc.
2447: Or you can run with -ksp_error_if_not_converged and the program will be stopped and the information printed in the error message.
2449: .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage, MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog(), KSPGetPC(), PCGetFactor(), PCFactorGetMatrix()
2451: M*/
2455: static PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
2456: {
2458: *type = MATSOLVERMUMPS;
2459: return(0);
2460: }
2462: /* MatGetFactor for Seq and MPI AIJ matrices */
2465: static PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
2466: {
2467: Mat B;
2469: Mat_MUMPS *mumps;
2470: PetscBool isSeqAIJ;
2473: /* Create the factorization matrix */
2474: PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);
2475: MatCreate(PetscObjectComm((PetscObject)A),&B);
2476: MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2477: MatSetType(B,((PetscObject)A)->type_name);
2478: if (isSeqAIJ) {
2479: MatSeqAIJSetPreallocation(B,0,NULL);
2480: } else {
2481: MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);
2482: }
2484: PetscNewLog(B,&mumps);
2486: B->ops->view = MatView_MUMPS;
2487: B->ops->getinfo = MatGetInfo_MUMPS;
2488: B->ops->getdiagonal = MatGetDiagonal_MUMPS;
2490: PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);
2491: PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2492: PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);
2493: PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2494: PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);
2495: PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);
2496: PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);
2497: PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2498: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2499: PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2500: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2501: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2502: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2503: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2504: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);
2506: if (ftype == MAT_FACTOR_LU) {
2507: B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
2508: B->factortype = MAT_FACTOR_LU;
2509: if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
2510: else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
2511: mumps->sym = 0;
2512: } else {
2513: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2514: B->factortype = MAT_FACTOR_CHOLESKY;
2515: if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
2516: else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
2517: #if defined(PETSC_USE_COMPLEX)
2518: mumps->sym = 2;
2519: #else
2520: if (A->spd_set && A->spd) mumps->sym = 1;
2521: else mumps->sym = 2;
2522: #endif
2523: }
2525: /* set solvertype */
2526: PetscFree(B->solvertype);
2527: PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);
2529: mumps->isAIJ = PETSC_TRUE;
2530: mumps->Destroy = B->ops->destroy;
2531: B->ops->destroy = MatDestroy_MUMPS;
2532: B->spptr = (void*)mumps;
2534: PetscInitializeMUMPS(A,mumps);
2536: *F = B;
2537: return(0);
2538: }
2540: /* MatGetFactor for Seq and MPI SBAIJ matrices */
2543: static PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
2544: {
2545: Mat B;
2547: Mat_MUMPS *mumps;
2548: PetscBool isSeqSBAIJ;
2551: if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
2552: if (A->rmap->bs > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with block size > 1 with MUMPS Cholesky, use AIJ matrix instead");
2553: PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);
2554: /* Create the factorization matrix */
2555: MatCreate(PetscObjectComm((PetscObject)A),&B);
2556: MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2557: MatSetType(B,((PetscObject)A)->type_name);
2558: PetscNewLog(B,&mumps);
2559: if (isSeqSBAIJ) {
2560: MatSeqSBAIJSetPreallocation(B,1,0,NULL);
2562: mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
2563: } else {
2564: MatMPISBAIJSetPreallocation(B,1,0,NULL,0,NULL);
2566: mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
2567: }
2569: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2570: B->ops->view = MatView_MUMPS;
2571: B->ops->getdiagonal = MatGetDiagonal_MUMPS;
2573: PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);
2574: PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2575: PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);
2576: PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2577: PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);
2578: PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);
2579: PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);
2580: PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2581: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2582: PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2583: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2584: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2585: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2586: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2587: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);
2589: B->factortype = MAT_FACTOR_CHOLESKY;
2590: #if defined(PETSC_USE_COMPLEX)
2591: mumps->sym = 2;
2592: #else
2593: if (A->spd_set && A->spd) mumps->sym = 1;
2594: else mumps->sym = 2;
2595: #endif
2597: /* set solvertype */
2598: PetscFree(B->solvertype);
2599: PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);
2601: mumps->isAIJ = PETSC_FALSE;
2602: mumps->Destroy = B->ops->destroy;
2603: B->ops->destroy = MatDestroy_MUMPS;
2604: B->spptr = (void*)mumps;
2606: PetscInitializeMUMPS(A,mumps);
2608: *F = B;
2609: return(0);
2610: }
2614: static PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
2615: {
2616: Mat B;
2618: Mat_MUMPS *mumps;
2619: PetscBool isSeqBAIJ;
2622: /* Create the factorization matrix */
2623: PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);
2624: MatCreate(PetscObjectComm((PetscObject)A),&B);
2625: MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2626: MatSetType(B,((PetscObject)A)->type_name);
2627: if (isSeqBAIJ) {
2628: MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,NULL);
2629: } else {
2630: MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,NULL,0,NULL);
2631: }
2633: PetscNewLog(B,&mumps);
2634: if (ftype == MAT_FACTOR_LU) {
2635: B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
2636: B->factortype = MAT_FACTOR_LU;
2637: if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
2638: else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
2639: mumps->sym = 0;
2640: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
2642: B->ops->view = MatView_MUMPS;
2643: B->ops->getdiagonal = MatGetDiagonal_MUMPS;
2645: PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);
2646: PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);
2647: PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);
2648: PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);
2649: PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);
2650: PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);
2651: PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);
2652: PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);
2653: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);
2654: PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);
2655: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);
2656: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);
2657: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);
2658: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);
2659: PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);
2661: /* set solvertype */
2662: PetscFree(B->solvertype);
2663: PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);
2665: mumps->isAIJ = PETSC_TRUE;
2666: mumps->Destroy = B->ops->destroy;
2667: B->ops->destroy = MatDestroy_MUMPS;
2668: B->spptr = (void*)mumps;
2670: PetscInitializeMUMPS(A,mumps);
2672: *F = B;
2673: return(0);
2674: }
2678: PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MUMPS(void)
2679: {
2683: MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);
2684: MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);
2685: MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);
2686: MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);
2687: MatSolverPackageRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);
2688: MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);
2689: MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);
2690: MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);
2691: MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);
2692: MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);
2693: return(0);
2694: }