Actual source code: clique.cxx

petsc-dev 2014-02-02
Report Typos and Errors
  1: #include <../src/mat/impls/aij/mpi/clique/matcliqueimpl.h> /*I "petscmat.h" I*/
  2: /*
  3:  Provides an interface to the Clique sparse solver (http://poulson.github.com/Clique/)
  4: */

  8: PetscErrorCode PetscCliqueFinalizePackage(void)
  9: {
 11:   cliq::Finalize();
 12:   return(0);
 13: }

 17: PetscErrorCode PetscCliqueInitializePackage(void)
 18: {

 22:   if (cliq::Initialized()) return(0);
 23:   { /* We have already initialized MPI, so this song and dance is just to pass these variables (which won't be used by Clique) through the interface that needs references */
 24:     int zero = 0;
 25:     char **nothing = 0;
 26:     cliq::Initialize(zero,nothing);
 27:   }
 28:   PetscRegisterFinalize(PetscCliqueFinalizePackage);
 29:   return(0);
 30: }

 32: /*
 33:   MatConvertToClique: Convert Petsc aij matrix to Clique matrix

 35:   input:
 36: +   A     - matrix in seqaij or mpiaij format
 37: -   reuse - denotes if the destination matrix is to be created or reused. Currently
 38:             MAT_REUSE_MATRIX is only supported for inplace conversion, otherwise use MAT_INITIAL_MATRIX.

 40:   output:
 41: .   cliq - Clique context
 42: */
 45: PetscErrorCode MatConvertToClique(Mat A,MatReuse reuse,Mat_Clique *cliq)
 46: {
 47:   PetscErrorCode        ierr;
 48:   PetscInt              i,j,rstart,rend,ncols;
 49:   const PetscInt        *cols;
 50:   const PetscCliqScalar *vals;
 51:   cliq::DistSparseMatrix<PetscCliqScalar> *cmat;

 54:   if (reuse == MAT_INITIAL_MATRIX){
 55:     /* create Clique matrix */
 56:     cmat = new cliq::DistSparseMatrix<PetscCliqScalar>(A->rmap->N,cliq->cliq_comm);
 57:     cliq->cmat = cmat;
 58:   } else {
 59:     cmat = cliq->cmat;
 60:   }
 61:   /* fill matrix values */
 62:   MatGetOwnershipRange(A,&rstart,&rend);
 63:   const int firstLocalRow = cmat->FirstLocalRow();
 64:   const int localHeight = cmat->LocalHeight();
 65:   if (rstart != firstLocalRow || rend-rstart != localHeight) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"matrix rowblock distribution does not match");

 67:   cmat->StartAssembly();
 68:   //cmat->Reserve( 7*localHeight ); ???
 69:   for (i=rstart; i<rend; i++){
 70:     MatGetRow(A,i,&ncols,&cols,&vals);
 71:     for (j=0; j<ncols; j++){
 72:       cmat->Update(i,cols[j],vals[j]);
 73:     }
 74:     MatRestoreRow(A,i,&ncols,&cols,&vals);
 75:   }
 76:   cmat->StopAssembly();
 77:   return(0);
 78: }

 82: static PetscErrorCode MatMult_Clique(Mat A,Vec X,Vec Y)
 83: {
 84:   PetscErrorCode        ierr;
 85:   PetscInt              i;
 86:   const PetscCliqScalar *x;
 87:   Mat_Clique            *cliq=(Mat_Clique*)A->spptr;
 88:   cliq::DistSparseMatrix<PetscCliqScalar> *cmat=cliq->cmat;
 89:   cliq::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A));

 92:   if (!cmat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Clique matrix cmat is not created yet");
 93:   VecGetArrayRead(X,(const PetscScalar **)&x);

 95:   cliq::DistMultiVec<PetscCliqScalar> xc(A->cmap->N,1,cxxcomm);
 96:   cliq::DistMultiVec<PetscCliqScalar> yc(A->rmap->N,1,cxxcomm);
 97:   for (i=0; i<A->cmap->n; i++) {
 98:     xc.SetLocal(i,0,x[i]);
 99:   }
100:   cliq::Multiply(1.0,*cmat,xc,0.0,yc);
101:   VecRestoreArrayRead(X,(const PetscScalar **)&x);

103:   for (i=0; i<A->cmap->n; i++) {
104:     VecSetValueLocal(Y,i,yc.GetLocal(i,0),INSERT_VALUES);
105:   }
106:   VecAssemblyBegin(Y);
107:   VecAssemblyEnd(Y);
108:   return(0);
109: }

113: PetscErrorCode MatView_Clique(Mat A,PetscViewer viewer)
114: {
116:   PetscBool      iascii;

119:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
120:   if (iascii) {
121:     PetscViewerFormat format;
122:     PetscViewerGetFormat(viewer,&format);
123:     if (format == PETSC_VIEWER_ASCII_INFO) {
124:       PetscViewerASCIIPrintf(viewer,"Clique run parameters:\n");
125:     } else if (format == PETSC_VIEWER_DEFAULT) { /* matrix A is factored matrix, remove this block */
126:       Mat Aaij;
127:       PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
128:       PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer);
129:       PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
130:       PetscPrintf(PetscObjectComm((PetscObject)viewer),"Clique matrix\n");
131:       MatComputeExplicitOperator(A,&Aaij);
132:       MatView(Aaij,PETSC_VIEWER_STDOUT_WORLD);
133:       MatDestroy(&Aaij);
134:     } else SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Format");
135:   }
136:   return(0);
137: }

141: PetscErrorCode MatDestroy_Clique(Mat A)
142: {
144:   Mat_Clique     *cliq=(Mat_Clique*)A->spptr;

147:   printf("MatDestroy_Clique ...\n");
148:   if (cliq && cliq->CleanUpClique) {
149:     /* Terminate instance, deallocate memories */
150:     printf("MatDestroy_Clique ... destroy clique struct \n");
151:     PetscCommDestroy(&(cliq->cliq_comm));
152:     // free cmat here
153:     delete cliq->cmat;
154:     delete cliq->frontTree;
155:     delete cliq->rhs;
156:     delete cliq->xNodal;
157:     delete cliq->info;
158:     delete cliq->inverseMap;
159:   }
160:   if (cliq && cliq->Destroy) {
161:     cliq->Destroy(A);
162:   }
163:   PetscFree(A->spptr);

165:   /* clear composed functions */
166:   PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);

168:   return(0);
169: }

173: PetscErrorCode MatSolve_Clique(Mat A,Vec B,Vec X)
174: {
175:   PetscErrorCode        ierr;
176:   PetscInt              i,rank;
177:   const PetscCliqScalar *b;
178:   Mat_Clique            *cliq=(Mat_Clique*)A->spptr;
179:   cliq::DistMultiVec<PetscCliqScalar> *bc=cliq->rhs;
180:   cliq::DistNodalMultiVec<PetscCliqScalar> *xNodal=cliq->xNodal;

183:   VecGetArrayRead(B,(const PetscScalar **)&b);
184:   for (i=0; i<A->rmap->n; i++) {
185:     bc->SetLocal(i,0,b[i]);
186:   }
187:   VecRestoreArrayRead(B,(const PetscScalar **)&b);

189:   xNodal->Pull( *cliq->inverseMap, *cliq->info, *bc );
190:   cliq::Solve( *cliq->info, *cliq->frontTree, *xNodal);
191:   xNodal->Push( *cliq->inverseMap, *cliq->info, *bc );

193:   MPI_Comm_rank(cliq->cliq_comm,&rank);
194:   for (i=0; i<bc->LocalHeight(); i++) {
195:     VecSetValue(X,rank*bc->Blocksize()+i,bc->GetLocal(i,0),INSERT_VALUES);
196:   }
197:   VecAssemblyBegin(X);
198:   VecAssemblyEnd(X);
199:   return(0);
200: }

204: PetscErrorCode MatCholeskyFactorNumeric_Clique(Mat F,Mat A,const MatFactorInfo *info)
205: {
206:   PetscErrorCode    ierr;
207:   Mat_Clique        *cliq=(Mat_Clique*)F->spptr;
208:   PETSC_UNUSED
209:   cliq::DistSparseMatrix<PetscCliqScalar> *cmat;

212:   cmat = cliq->cmat;
213:   if (cliq->matstruc == SAME_NONZERO_PATTERN){ /* successing numerical factorization */
214:     /* Update cmat */
215:     MatConvertToClique(A,MAT_REUSE_MATRIX,cliq);
216:   }

218:   /* Numeric factorization */
219:   cliq::LDL( *cliq->info, *cliq->frontTree, cliq::LDL_1D);
220:   //L.frontType = cliq::SYMM_2D;

222:   // refactor
223:   //cliq::ChangeFrontType( *cliq->frontTree, cliq::LDL_2D );
224:   //*(cliq->frontTree.frontType) = cliq::LDL_2D;
225:   //cliq::LDL( *cliq->info, *cliq->frontTree, cliq::LDL_2D );

227:   cliq->matstruc = SAME_NONZERO_PATTERN;
228:   F->assembled   = PETSC_TRUE;
229:   return(0);
230: }

234: PetscErrorCode MatCholeskyFactorSymbolic_Clique(Mat F,Mat A,IS r,const MatFactorInfo *info)
235: {
236:   PetscErrorCode    ierr;
237:   Mat_Clique        *Acliq=(Mat_Clique*)F->spptr;
238:   cliq::DistSparseMatrix<PetscCliqScalar> *cmat;
239:   cliq::DistSeparatorTree                 sepTree;
240:   cliq::DistMap                           map;

243:   /* Convert A to Aclique */
244:   MatConvertToClique(A,MAT_INITIAL_MATRIX,Acliq);
245:   cmat = Acliq->cmat;

247:   cliq::NestedDissection( cmat->DistGraph(), map, sepTree, *Acliq->info, PETSC_TRUE, Acliq->numDistSeps, Acliq->numSeqSeps, Acliq->cutoff);
248:   map.FormInverse( *Acliq->inverseMap );
249:   Acliq->frontTree = new cliq::DistSymmFrontTree<PetscCliqScalar>( cliq::TRANSPOSE, *cmat, map, sepTree, *Acliq->info );

251:   Acliq->matstruc      = DIFFERENT_NONZERO_PATTERN;
252:   Acliq->CleanUpClique = PETSC_TRUE;
253:   return(0);
254: }

256: /*MC
257:      MATSOLVERCLIQUE  - A solver package providing direct solvers for distributed
258:   and sequential matrices via the external package Clique.

260:   Use ./configure --download-clique to have PETSc installed with Clique

262:   Options Database Keys:
263: + -mat_clique_    -
264: - -mat_clique_ <integer> -

266:   Level: beginner

268: .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage

270: M*/

274: PetscErrorCode MatFactorGetSolverPackage_Clique(Mat A,const MatSolverPackage *type)
275: {
277:   *type = MATSOLVERCLIQUE;
278:   return(0);
279: }

283: PETSC_EXTERN PetscErrorCode MatGetFactor_aij_clique(Mat A,MatFactorType ftype,Mat *F)
284: {
285:   Mat            B;
286:   Mat_Clique     *cliq;

290:   PetscCliqueInitializePackage();
291:   MatCreate(PetscObjectComm((PetscObject)A),&B);
292:   MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);
293:   MatSetType(B,((PetscObject)A)->type_name);
294:   MatSetUp(B);

296:   if (ftype == MAT_FACTOR_CHOLESKY){
297:     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_Clique;
298:     B->ops->choleskyfactornumeric  = MatCholeskyFactorNumeric_Clique;
299:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");

301:   PetscNewLog(B,&cliq);
302:   B->spptr            = (void*)cliq;
303:   cliq::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A));
304:   PetscCommDuplicate(cxxcomm,&(cliq->cliq_comm),NULL);
305:   cliq->rhs           = new cliq::DistMultiVec<PetscCliqScalar>(A->rmap->N,1,cliq->cliq_comm);
306:   cliq->xNodal        = new cliq::DistNodalMultiVec<PetscCliqScalar>();
307:   cliq->info          = new cliq::DistSymmInfo;
308:   cliq->inverseMap    = new cliq::DistMap;
309:   cliq->CleanUpClique = PETSC_FALSE;
310:   cliq->Destroy       = B->ops->destroy;

312:   B->ops->view    = MatView_Clique;
313:   B->ops->mult    = MatMult_Clique; /* for cliq->cmat */
314:   B->ops->solve   = MatSolve_Clique;

316:   B->ops->destroy = MatDestroy_Clique;
317:   B->factortype   = ftype;
318:   B->assembled    = PETSC_FALSE;

320:   /* Set Clique options */
321:   PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"Clique Options","Mat");
322:   cliq->cutoff      = 128;  /* maximum size of leaf node */
323:   cliq->numDistSeps = 1;    /* number of distributed separators to try */
324:   cliq->numSeqSeps  = 1;    /* number of sequential separators to try */
325:   PetscOptionsEnd();

327:   *F = B;
328:   return(0);
329: }