PageRenderTime 56ms CodeModel.GetById 28ms RepoModel.GetById 1ms app.codeStats 0ms

/src/mat/color/impls/minpack/color.c

https://bitbucket.org/petsc/petsc
C | 375 lines | 242 code | 55 blank | 78 comment | 31 complexity | f7035dd58092759356bfaeaa7e8e70a5 MD5 | raw file
Possible License(s): BSD-3-Clause, LGPL-2.1
  1. /*
  2. Routines that call the kernel minpack coloring subroutines
  3. */
  4. #include <petsc/private/matimpl.h>
  5. #include <petsc/private/isimpl.h>
  6. #include <../src/mat/color/impls/minpack/color.h>
  7. /*
  8. MatFDColoringDegreeSequence_Minpack - Calls the MINPACK routine seqr() that
  9. computes the degree sequence required by MINPACK coloring routines.
  10. */
  11. PETSC_INTERN PetscErrorCode MatFDColoringDegreeSequence_Minpack(PetscInt m,const PetscInt *cja,const PetscInt *cia,const PetscInt *rja,const PetscInt *ria,PetscInt **seq)
  12. {
  13. PetscInt *work;
  14. PetscErrorCode ierr;
  15. PetscFunctionBegin;
  16. ierr = PetscMalloc1(m,&work);CHKERRQ(ierr);
  17. ierr = PetscMalloc1(m,seq);CHKERRQ(ierr);
  18. MINPACKdegr(&m,cja,cia,rja,ria,*seq,work);
  19. ierr = PetscFree(work);CHKERRQ(ierr);
  20. PetscFunctionReturn(0);
  21. }
  22. /*
  23. MatFDColoringMinimumNumberofColors_Private - For a given sparse
  24. matrix computes the minimum number of colors needed.
  25. */
  26. PetscErrorCode MatFDColoringMinimumNumberofColors_Private(PetscInt m,PetscInt *ia,PetscInt *minc)
  27. {
  28. PetscInt i,c = 0;
  29. PetscFunctionBegin;
  30. for (i=0; i<m; i++) c = PetscMax(c,ia[i+1]-ia[i]);
  31. *minc = c;
  32. PetscFunctionReturn(0);
  33. }
  34. static PetscErrorCode MatColoringApply_SL(MatColoring mc,ISColoring *iscoloring)
  35. {
  36. PetscErrorCode ierr;
  37. PetscInt *list,*work,clique,*seq,*coloring,n;
  38. const PetscInt *ria,*rja,*cia,*cja;
  39. PetscInt ncolors,i;
  40. PetscBool done;
  41. Mat mat = mc->mat;
  42. Mat mat_seq = mat;
  43. PetscMPIInt size;
  44. MPI_Comm comm;
  45. ISColoring iscoloring_seq;
  46. PetscInt bs = 1,rstart,rend,N_loc,nc;
  47. ISColoringValue *colors_loc;
  48. PetscBool flg1,flg2;
  49. PetscFunctionBegin;
  50. if (mc->dist != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SL may only do distance 2 coloring");
  51. /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
  52. ierr = PetscObjectBaseTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr);
  53. ierr = PetscObjectBaseTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr);
  54. if (flg1 || flg2) {
  55. ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
  56. }
  57. ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
  58. ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
  59. if (size > 1) {
  60. /* create a sequential iscoloring on all processors */
  61. ierr = MatGetSeqNonzeroStructure(mat,&mat_seq);CHKERRQ(ierr);
  62. }
  63. ierr = MatGetRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);CHKERRQ(ierr);
  64. ierr = MatGetColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);CHKERRQ(ierr);
  65. if (!done) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Ordering requires IJ");
  66. ierr = MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);CHKERRQ(ierr);
  67. ierr = PetscMalloc2(n,&list,4*n,&work);CHKERRQ(ierr);
  68. MINPACKslo(&n,cja,cia,rja,ria,seq,list,&clique,work,work+n,work+2*n,work+3*n);
  69. ierr = PetscMalloc1(n,&coloring);CHKERRQ(ierr);
  70. MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);
  71. ierr = PetscFree2(list,work);CHKERRQ(ierr);
  72. ierr = PetscFree(seq);CHKERRQ(ierr);
  73. ierr = MatRestoreRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&ria,&rja,&done);CHKERRQ(ierr);
  74. ierr = MatRestoreColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&cia,&cja,&done);CHKERRQ(ierr);
  75. /* shift coloring numbers to start at zero and shorten */
  76. if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Maximum color size exceeded");
  77. {
  78. ISColoringValue *s = (ISColoringValue*) coloring;
  79. for (i=0; i<n; i++) {
  80. s[i] = (ISColoringValue) (coloring[i]-1);
  81. }
  82. ierr = MatColoringPatch(mat_seq,ncolors,n,s,iscoloring);CHKERRQ(ierr);
  83. }
  84. if (size > 1) {
  85. ierr = MatDestroySeqNonzeroStructure(&mat_seq);CHKERRQ(ierr);
  86. /* convert iscoloring_seq to a parallel iscoloring */
  87. iscoloring_seq = *iscoloring;
  88. rstart = mat->rmap->rstart/bs;
  89. rend = mat->rmap->rend/bs;
  90. N_loc = rend - rstart; /* number of local nodes */
  91. /* get local colors for each local node */
  92. ierr = PetscMalloc1(N_loc+1,&colors_loc);CHKERRQ(ierr);
  93. for (i=rstart; i<rend; i++) {
  94. colors_loc[i-rstart] = iscoloring_seq->colors[i];
  95. }
  96. /* create a parallel iscoloring */
  97. nc = iscoloring_seq->n;
  98. ierr = ISColoringCreate(comm,nc,N_loc,colors_loc,PETSC_OWN_POINTER,iscoloring);CHKERRQ(ierr);
  99. ierr = ISColoringDestroy(&iscoloring_seq);CHKERRQ(ierr);
  100. }
  101. PetscFunctionReturn(0);
  102. }
  103. /*MC
  104. MATCOLORINGSL - implements the SL (smallest last) coloring routine
  105. Level: beginner
  106. Notes:
  107. Supports only distance two colorings (for computation of Jacobians)
  108. This is a sequential algorithm
  109. References:
  110. . 1. - TF Coleman and J More, "Estimation of sparse Jacobian matrices and graph coloring," SIAM Journal on Numerical Analysis, vol. 20, no. 1,
  111. pp. 187-209, 1983.
  112. .seealso: MatColoringCreate(), MatColoring, MatColoringSetType(), MATCOLORINGGREEDY, MatColoringType
  113. M*/
  114. PETSC_EXTERN PetscErrorCode MatColoringCreate_SL(MatColoring mc)
  115. {
  116. PetscFunctionBegin;
  117. mc->dist = 2;
  118. mc->data = NULL;
  119. mc->ops->apply = MatColoringApply_SL;
  120. mc->ops->view = NULL;
  121. mc->ops->destroy = NULL;
  122. mc->ops->setfromoptions = NULL;
  123. PetscFunctionReturn(0);
  124. }
  125. static PetscErrorCode MatColoringApply_LF(MatColoring mc,ISColoring *iscoloring)
  126. {
  127. PetscErrorCode ierr;
  128. PetscInt *list,*work,*seq,*coloring,n;
  129. const PetscInt *ria,*rja,*cia,*cja;
  130. PetscInt n1, none,ncolors,i;
  131. PetscBool done;
  132. Mat mat = mc->mat;
  133. Mat mat_seq = mat;
  134. PetscMPIInt size;
  135. MPI_Comm comm;
  136. ISColoring iscoloring_seq;
  137. PetscInt bs = 1,rstart,rend,N_loc,nc;
  138. ISColoringValue *colors_loc;
  139. PetscBool flg1,flg2;
  140. PetscFunctionBegin;
  141. if (mc->dist != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LF may only do distance 2 coloring");
  142. /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
  143. ierr = PetscObjectBaseTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr);
  144. ierr = PetscObjectBaseTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr);
  145. if (flg1 || flg2) {
  146. ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
  147. }
  148. ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
  149. ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
  150. if (size > 1) {
  151. /* create a sequential iscoloring on all processors */
  152. ierr = MatGetSeqNonzeroStructure(mat,&mat_seq);CHKERRQ(ierr);
  153. }
  154. ierr = MatGetRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);CHKERRQ(ierr);
  155. ierr = MatGetColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);CHKERRQ(ierr);
  156. if (!done) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Ordering requires IJ");
  157. ierr = MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);CHKERRQ(ierr);
  158. ierr = PetscMalloc2(n,&list,4*n,&work);CHKERRQ(ierr);
  159. n1 = n - 1;
  160. none = -1;
  161. MINPACKnumsrt(&n,&n1,seq,&none,list,work+2*n,work+n);
  162. ierr = PetscMalloc1(n,&coloring);CHKERRQ(ierr);
  163. MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);
  164. ierr = PetscFree2(list,work);CHKERRQ(ierr);
  165. ierr = PetscFree(seq);CHKERRQ(ierr);
  166. ierr = MatRestoreRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&ria,&rja,&done);CHKERRQ(ierr);
  167. ierr = MatRestoreColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&cia,&cja,&done);CHKERRQ(ierr);
  168. /* shift coloring numbers to start at zero and shorten */
  169. if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Maximum color size exceeded");
  170. {
  171. ISColoringValue *s = (ISColoringValue*) coloring;
  172. for (i=0; i<n; i++) s[i] = (ISColoringValue) (coloring[i]-1);
  173. ierr = MatColoringPatch(mat_seq,ncolors,n,s,iscoloring);CHKERRQ(ierr);
  174. }
  175. if (size > 1) {
  176. ierr = MatDestroySeqNonzeroStructure(&mat_seq);CHKERRQ(ierr);
  177. /* convert iscoloring_seq to a parallel iscoloring */
  178. iscoloring_seq = *iscoloring;
  179. rstart = mat->rmap->rstart/bs;
  180. rend = mat->rmap->rend/bs;
  181. N_loc = rend - rstart; /* number of local nodes */
  182. /* get local colors for each local node */
  183. ierr = PetscMalloc1(N_loc+1,&colors_loc);CHKERRQ(ierr);
  184. for (i=rstart; i<rend; i++) colors_loc[i-rstart] = iscoloring_seq->colors[i];
  185. /* create a parallel iscoloring */
  186. nc = iscoloring_seq->n;
  187. ierr = ISColoringCreate(comm,nc,N_loc,colors_loc,PETSC_OWN_POINTER,iscoloring);CHKERRQ(ierr);
  188. ierr = ISColoringDestroy(&iscoloring_seq);CHKERRQ(ierr);
  189. }
  190. PetscFunctionReturn(0);
  191. }
  192. /*MC
  193. MATCOLORINGLF - implements the LF (largest first) coloring routine
  194. Level: beginner
  195. Notes:
  196. Supports only distance two colorings (for computation of Jacobians)
  197. This is a sequential algorithm
  198. References:
  199. . 1. - TF Coleman and J More, "Estimation of sparse Jacobian matrices and graph coloring," SIAM Journal on Numerical Analysis, vol. 20, no. 1,
  200. pp. 187-209, 1983.
  201. .seealso: MatColoringCreate(), MatColoring, MatColoringSetType(), MATCOLORINGGREEDY, MatColoringType
  202. M*/
  203. PETSC_EXTERN PetscErrorCode MatColoringCreate_LF(MatColoring mc)
  204. {
  205. PetscFunctionBegin;
  206. mc->dist = 2;
  207. mc->data = NULL;
  208. mc->ops->apply = MatColoringApply_LF;
  209. mc->ops->view = NULL;
  210. mc->ops->destroy = NULL;
  211. mc->ops->setfromoptions = NULL;
  212. PetscFunctionReturn(0);
  213. }
  214. static PetscErrorCode MatColoringApply_ID(MatColoring mc,ISColoring *iscoloring)
  215. {
  216. PetscErrorCode ierr;
  217. PetscInt *list,*work,clique,*seq,*coloring,n;
  218. const PetscInt *ria,*rja,*cia,*cja;
  219. PetscInt ncolors,i;
  220. PetscBool done;
  221. Mat mat = mc->mat;
  222. Mat mat_seq = mat;
  223. PetscMPIInt size;
  224. MPI_Comm comm;
  225. ISColoring iscoloring_seq;
  226. PetscInt bs = 1,rstart,rend,N_loc,nc;
  227. ISColoringValue *colors_loc;
  228. PetscBool flg1,flg2;
  229. PetscFunctionBegin;
  230. if (mc->dist != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"IDO may only do distance 2 coloring");
  231. /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
  232. ierr = PetscObjectBaseTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr);
  233. ierr = PetscObjectBaseTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr);
  234. if (flg1 || flg2) {
  235. ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
  236. }
  237. ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
  238. ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
  239. if (size > 1) {
  240. /* create a sequential iscoloring on all processors */
  241. ierr = MatGetSeqNonzeroStructure(mat,&mat_seq);CHKERRQ(ierr);
  242. }
  243. ierr = MatGetRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&ria,&rja,&done);CHKERRQ(ierr);
  244. ierr = MatGetColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,&n,&cia,&cja,&done);CHKERRQ(ierr);
  245. if (!done) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Ordering requires IJ");
  246. ierr = MatFDColoringDegreeSequence_Minpack(n,cja,cia,rja,ria,&seq);CHKERRQ(ierr);
  247. ierr = PetscMalloc2(n,&list,4*n,&work);CHKERRQ(ierr);
  248. MINPACKido(&n,&n,cja,cia,rja,ria,seq,list,&clique,work,work+n,work+2*n,work+3*n);
  249. ierr = PetscMalloc1(n,&coloring);CHKERRQ(ierr);
  250. MINPACKseq(&n,cja,cia,rja,ria,list,coloring,&ncolors,work);
  251. ierr = PetscFree2(list,work);CHKERRQ(ierr);
  252. ierr = PetscFree(seq);CHKERRQ(ierr);
  253. ierr = MatRestoreRowIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&ria,&rja,&done);CHKERRQ(ierr);
  254. ierr = MatRestoreColumnIJ(mat_seq,1,PETSC_FALSE,PETSC_TRUE,NULL,&cia,&cja,&done);CHKERRQ(ierr);
  255. /* shift coloring numbers to start at zero and shorten */
  256. if (ncolors > IS_COLORING_MAX-1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Maximum color size exceeded");
  257. {
  258. ISColoringValue *s = (ISColoringValue*) coloring;
  259. for (i=0; i<n; i++) {
  260. s[i] = (ISColoringValue) (coloring[i]-1);
  261. }
  262. ierr = MatColoringPatch(mat_seq,ncolors,n,s,iscoloring);CHKERRQ(ierr);
  263. }
  264. if (size > 1) {
  265. ierr = MatDestroySeqNonzeroStructure(&mat_seq);CHKERRQ(ierr);
  266. /* convert iscoloring_seq to a parallel iscoloring */
  267. iscoloring_seq = *iscoloring;
  268. rstart = mat->rmap->rstart/bs;
  269. rend = mat->rmap->rend/bs;
  270. N_loc = rend - rstart; /* number of local nodes */
  271. /* get local colors for each local node */
  272. ierr = PetscMalloc1(N_loc+1,&colors_loc);CHKERRQ(ierr);
  273. for (i=rstart; i<rend; i++) {
  274. colors_loc[i-rstart] = iscoloring_seq->colors[i];
  275. }
  276. /* create a parallel iscoloring */
  277. nc = iscoloring_seq->n;
  278. ierr = ISColoringCreate(comm,nc,N_loc,colors_loc,PETSC_OWN_POINTER,iscoloring);CHKERRQ(ierr);
  279. ierr = ISColoringDestroy(&iscoloring_seq);CHKERRQ(ierr);
  280. }
  281. PetscFunctionReturn(0);
  282. }
  283. /*MC
  284. MATCOLORINGID - implements the ID (incidence degree) coloring routine
  285. Level: beginner
  286. Notes:
  287. Supports only distance two colorings (for computation of Jacobians)
  288. This is a sequential algorithm
  289. References:
  290. . 1. - TF Coleman and J More, "Estimation of sparse Jacobian matrices and graph coloring," SIAM Journal on Numerical Analysis, vol. 20, no. 1,
  291. pp. 187-209, 1983.
  292. .seealso: MatColoringCreate(), MatColoring, MatColoringSetType(), MATCOLORINGGREEDY, MatColoringType
  293. M*/
  294. PETSC_EXTERN PetscErrorCode MatColoringCreate_ID(MatColoring mc)
  295. {
  296. PetscFunctionBegin;
  297. mc->dist = 2;
  298. mc->data = NULL;
  299. mc->ops->apply = MatColoringApply_ID;
  300. mc->ops->view = NULL;
  301. mc->ops->destroy = NULL;
  302. mc->ops->setfromoptions = NULL;
  303. PetscFunctionReturn(0);
  304. }