Actual source code: mpisbaijspooles.c

  1: /*$Id: mpisbaijspooles.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/
  2: /* 
  3:    Provides an interface to the Spooles parallel sparse solver (MPI SPOOLES)
  4: */

 6:  #include src/mat/impls/aij/seq/spooles/spooles.h
 7:  #include src/mat/impls/sbaij/mpi/mpisbaij.h

 11: int MatDestroy_MPISBAIJSpooles(Mat A) {
 13: 
 15:   /* MPISBAIJ_Spooles isn't really the matrix that USES spooles, */
 16:   /* rather it is a factory class for creating a symmetric matrix that can */
 17:   /* invoke Spooles' parallel cholesky solver. */
 18:   /* As a result, we don't have to clean up the stuff set for use in spooles */
 19:   /* as in MatDestroy_MPIAIJ_Spooles. */
 20:   MatConvert_Spooles_Base(A,MATMPISBAIJ,&A);
 21:   (*A->ops->destroy)(A);
 22:   return(0);
 23: }

 27: int MatAssemblyEnd_MPISBAIJSpooles(Mat A,MatAssemblyType mode) {
 28:   int         ierr,bs;
 29:   Mat_Spooles *lu=(Mat_Spooles *)(A->spptr);

 32:   (*lu->MatAssemblyEnd)(A,mode);
 33:   MatGetBlockSize(A,&bs);
 34:   if (bs > 1) SETERRQ1(1,"Block size %d not supported by Spooles",bs);
 35:   lu->MatCholeskyFactorSymbolic  = A->ops->choleskyfactorsymbolic;
 36:   A->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MPISBAIJSpooles;
 37:   return(0);
 38: }

 40: /* 
 41:   input:
 42:    F:                 numeric factor
 43:   output:
 44:    nneg, nzero, npos: global matrix inertia in all processors
 45: */

 49: int MatGetInertia_MPISBAIJSpooles(Mat F,int *nneg,int *nzero,int *npos)
 50: {
 51:   Mat_Spooles *lu = (Mat_Spooles*)F->spptr;
 52:   int         ierr,neg,zero,pos,sbuf[3],rbuf[3];

 55:   FrontMtx_inertia(lu->frontmtx, &neg, &zero, &pos) ;
 56:   sbuf[0] = neg; sbuf[1] = zero; sbuf[2] = pos;
 57:   MPI_Allreduce(sbuf,rbuf,3,MPI_INT,MPI_SUM,F->comm);
 58:   *nneg  = rbuf[0]; *nzero = rbuf[1]; *npos  = rbuf[2];
 59:   return(0);
 60: }

 62: /* Note the Petsc r permutation is ignored */
 65: int MatCholeskyFactorSymbolic_MPISBAIJSpooles(Mat A,IS r,MatFactorInfo *info,Mat *F)
 66: {
 67:   Mat           B;
 68:   Mat_Spooles   *lu;
 69:   int           ierr;
 70: 

 73:   /* Create the factorization matrix */
 74:   MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);
 75:   MatSetType(B,A->type_name);
 76:   MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);
 77: 
 78:   B->ops->choleskyfactornumeric = MatFactorNumeric_MPIAIJSpooles;
 79:   B->ops->getinertia            = MatGetInertia_MPISBAIJSpooles;
 80:   B->factor                     = FACTOR_CHOLESKY;

 82:   lu                       = (Mat_Spooles*)(B->spptr);
 83:   lu->options.pivotingflag = SPOOLES_NO_PIVOTING;
 84:   lu->flg                  = DIFFERENT_NONZERO_PATTERN;
 85:   lu->options.useQR        = PETSC_FALSE;
 86:   lu->options.symflag      = SPOOLES_SYMMETRIC;  /* default */

 88:   MPI_Comm_dup(A->comm,&(lu->comm_spooles));
 89:   *F = B;
 90:   return(0);
 91: }

 93: EXTERN_C_BEGIN
 96: int MatMPISBAIJSetPreallocation_MPISBAIJSpooles(Mat  B,int bs,int d_nz,int *d_nnz,int o_nz,int *o_nnz)
 97: {
 98:   Mat         A;
 99:   Mat_Spooles *lu = (Mat_Spooles*)B->spptr;
100:   int         ierr;

103:   /*
104:     After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix
105:     into Spooles type so that the block jacobi preconditioner (for example) can use Spooles.  I would
106:     like this to be done in the MatCreate routine, but the creation of this inner matrix requires
107:     block size info so that PETSc can determine the local size properly.  The block size info is set
108:     in the preallocation routine.
109:   */
110:   (*lu->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz);
111:   A    = ((Mat_MPISBAIJ *)B->data)->A;
112:   MatConvert_SeqSBAIJ_SeqSBAIJSpooles(A,MATSEQSBAIJSPOOLES,&A);
113:   return(0);
114: }
115: EXTERN_C_END

117: EXTERN_C_BEGIN
120: int MatConvert_MPISBAIJ_MPISBAIJSpooles(Mat A,const MatType type,Mat *newmat) {
121:   /* This routine is only called to convert a MATMPISBAIJ matrix */
122:   /* to a MATMPISBAIJSPOOLES matrix, so we will ignore 'MatType type'. */
123:   int         ierr;
124:   Mat         B=*newmat;
125:   Mat_Spooles *lu;
126:   void        (*f)(void);

129:   if (B != A) {
130:     /* This routine is inherited, so we know the type is correct. */
131:     MatDuplicate(A,MAT_COPY_VALUES,&B);
132:   }

134:   PetscNew(Mat_Spooles,&lu);
135:   B->spptr                       = (void*)lu;

137:   lu->basetype                   = MATMPISBAIJ;
138:   lu->MatDuplicate               = A->ops->duplicate;
139:   lu->MatCholeskyFactorSymbolic  = A->ops->choleskyfactorsymbolic;
140:   lu->MatLUFactorSymbolic        = A->ops->lufactorsymbolic;
141:   lu->MatView                    = A->ops->view;
142:   lu->MatAssemblyEnd             = A->ops->assemblyend;
143:   lu->MatDestroy                 = A->ops->destroy;
144:   B->ops->duplicate              = MatDuplicate_MPISBAIJSpooles;
145:   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MPISBAIJSpooles;
146:   B->ops->assemblyend            = MatAssemblyEnd_MPISBAIJSpooles;
147:   B->ops->destroy                = MatDestroy_MPISBAIJSpooles;

149:   /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */
150:   PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",&f);
151:   if (f) {
152:     lu->MatPreallocate = (int (*)(Mat,int,int,int*,int,int*))f;
153:     PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
154:                                              "MatMPISBAIJSetPreallocation_MPISBAIJSpooles",
155:                                              MatMPISBAIJSetPreallocation_MPISBAIJSpooles);
156:   }

158:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaijspooles_mpisbaij_C",
159:                                            "MatConvert_Spooles_Base",MatConvert_Spooles_Base);
160:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_mpisbaijspooles_C",
161:                                            "MatConvert_MPISBAIJ_MPISBAIJSpooles",
162:                                            MatConvert_MPISBAIJ_MPISBAIJSpooles);

164:   PetscObjectChangeTypeName((PetscObject)B,MATMPISBAIJSPOOLES);
165:   *newmat = B;
166:   return(0);
167: }
168: EXTERN_C_END

172: int MatDuplicate_MPISBAIJSpooles(Mat A, MatDuplicateOption op, Mat *M) {
173:   int         ierr;
174:   Mat_Spooles *lu=(Mat_Spooles *)A->spptr;

177:   (*lu->MatDuplicate)(A,op,M);
178:   PetscMemcpy((*M)->spptr,lu,sizeof(Mat_Spooles));
179:   return(0);
180: }

182: /*MC
183:   MATMPISBAIJSPOOLES - MATMPISBAIJSPOOLES = "mpisbaijspooles" - a matrix type providing direct solvers (Cholesky) for distributed symmetric
184:   matrices via the external package Spooles.

186:   If Spooles is installed (see the manual for
187:   instructions on how to declare the existence of external packages),
188:   a matrix type can be constructed which invokes Spooles solvers.
189:   After calling MatCreate(...,A), simply call MatSetType(A,MATMPISBAIJSPOOLES).
190:   This matrix type is only supported for double precision real.

192:   This matrix inherits from MATMPISBAIJ.  As a result, MatMPISBAIJSetPreallocation is 
193:   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from 
194:   the MATMPISBAIJ type without data copy.

196:   Options Database Keys:
197: + -mat_type mpisbaijspooles - sets the matrix type to mpisbaijspooles during a call to MatSetFromOptions()
198: . -mat_spooles_tau <tau> - upper bound on the magnitude of the largest element in L or U
199: . -mat_spooles_seed <seed> - random number seed used for ordering
200: . -mat_spooles_msglvl <msglvl> - message output level
201: . -mat_spooles_ordering <BestOfNDandMS,MMD,MS,ND> - ordering used
202: . -mat_spooles_maxdomainsize <n> - maximum subgraph size used by Spooles orderings
203: . -mat_spooles_maxzeros <n> - maximum number of zeros inside a supernode
204: . -mat_spooles_maxsize <n> - maximum size of a supernode
205: . -mat_spooles_FrontMtxInfo <true,fase> - print Spooles information about the computed factorization
206: . -mat_spooles_symmetryflag <0,1,2> - 0: SPOOLES_SYMMETRIC, 1: SPOOLES_HERMITIAN, 2: SPOOLES_NONSYMMETRIC
207: . -mat_spooles_patchAndGoFlag <0,1,2> - 0: no patch, 1: use PatchAndGo strategy 1, 2: use PatchAndGo strategy 2
208: . -mat_spooles_toosmall <dt> - drop tolerance for PatchAndGo strategy 1
209: . -mat_spooles_storeids <bool integer> - if nonzero, stores row and col numbers where patches were applied in an IV object
210: . -mat_spooles_fudge <delta> - fudge factor for rescaling diagonals with PatchAndGo strategy 2
211: - -mat_spooles_storevalues <bool integer> - if nonzero and PatchAndGo strategy 2 is used, store change in diagonal value in a DV object

213:    Level: beginner

215: .seealso: MATSEQSBAIJSPOOLES, MATSEQAIJSPOOLES, MATMPIAIJSPOOLES, PCCHOLESKY
216: M*/

218: EXTERN_C_BEGIN
221: int MatCreate_MPISBAIJSpooles(Mat A) {

225:   /* Change type name before calling MatSetType to force proper construction of MPISBAIJ */
226:   /*   and MPISBAIJSpooles types */
227:   PetscObjectChangeTypeName((PetscObject)A,MATMPISBAIJSPOOLES);
228:   MatSetType(A,MATMPISBAIJ);
229:   MatConvert_MPISBAIJ_MPISBAIJSpooles(A,MATMPISBAIJSPOOLES,&A);
230:   return(0);
231: }
232: EXTERN_C_END