Actual source code: bddc.c
petsc-3.6.1 2015-08-06
1: /* TODOLIST
3: Solvers
4: - Add support for cholesky for coarse solver (similar to local solvers)
5: - Propagate ksp prefixes for solvers to mat objects?
7: User interface
8: - ** DM attached to pc?
10: Debugging output
11: - * Better management of verbosity levels of debugging output
13: Extra
14: - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
15: - BDDC with MG framework?
17: FETIDP
18: - Move FETIDP code to its own classes
20: MATIS related operations contained in BDDC code
21: - Provide general case for subassembling
23: */
25: #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/ /* includes for fortran wrappers */
26: #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
27: #include <petscblaslapack.h>
29: /* -------------------------------------------------------------------------- */
32: PetscErrorCode PCSetFromOptions_BDDC(PetscOptions *PetscOptionsObject,PC pc)
33: {
34: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
38: PetscOptionsHead(PetscOptionsObject,"BDDC options");
39: /* Verbose debugging */
40: PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);
41: /* Primal space cumstomization */
42: PetscOptionsBool("-pc_bddc_use_local_mat_graph","Use or not adjacency graph of local mat for interface analysis","none",pcbddc->use_local_adj,&pcbddc->use_local_adj,NULL);
43: PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);
44: PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);
45: PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);
46: PetscOptionsBool("-pc_bddc_use_true_nnsp","Use near null space attached to the matrix without modifications","none",pcbddc->use_nnsp_true,&pcbddc->use_nnsp_true,NULL);
47: PetscOptionsBool("-pc_bddc_use_qr_single","Use QR factorization for single constraints on cc (QR is used when multiple constraints are present)","none",pcbddc->use_qr_single,&pcbddc->use_qr_single,NULL);
48: /* Change of basis */
49: PetscOptionsBool("-pc_bddc_use_change_of_basis","Use or not internal change of basis on local edge nodes","none",pcbddc->use_change_of_basis,&pcbddc->use_change_of_basis,NULL);
50: PetscOptionsBool("-pc_bddc_use_change_on_faces","Use or not internal change of basis on local face nodes","none",pcbddc->use_change_on_faces,&pcbddc->use_change_on_faces,NULL);
51: if (!pcbddc->use_change_of_basis) {
52: pcbddc->use_change_on_faces = PETSC_FALSE;
53: }
54: /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
55: PetscOptionsBool("-pc_bddc_switch_static","Switch on static condensation ops around the interface preconditioner","none",pcbddc->switch_static,&pcbddc->switch_static,NULL);
56: PetscOptionsInt("-pc_bddc_coarse_redistribute","Number of procs where to redistribute coarse problem","none",pcbddc->redistribute_coarse,&pcbddc->redistribute_coarse,NULL);
57: PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);
58: PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);
59: PetscOptionsBool("-pc_bddc_use_coarse_estimates","Use estimated eigenvalues for coarse problem","none",pcbddc->use_coarse_estimates,&pcbddc->use_coarse_estimates,NULL);
60: PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);
61: PetscOptionsBool("-pc_bddc_schur_rebuild","Whether or not the interface graph for Schur principal minors has to be rebuilt (i.e. define the interface without any adjacency)","none",pcbddc->sub_schurs_rebuild,&pcbddc->sub_schurs_rebuild,NULL);
62: PetscOptionsInt("-pc_bddc_schur_layers","Number of dofs' layers for the computation of principal minors (i.e. -1 uses all dofs)","none",pcbddc->sub_schurs_layers,&pcbddc->sub_schurs_layers,NULL);
63: PetscOptionsBool("-pc_bddc_schur_use_useradj","Whether or not the CSR graph specified by the user should be used for computing successive layers (default is to use adj of local mat)","none",pcbddc->sub_schurs_use_useradj,&pcbddc->sub_schurs_use_useradj,NULL);
64: PetscOptionsBool("-pc_bddc_deluxe_faster","Faster application of deluxe scaling (requires extra work during setup)","none",pcbddc->faster_deluxe,&pcbddc->faster_deluxe,NULL);
65: PetscOptionsReal("-pc_bddc_adaptive_threshold","Threshold to be used for adaptive selection of constraints","none",pcbddc->adaptive_threshold,&pcbddc->adaptive_threshold,NULL);
66: PetscOptionsInt("-pc_bddc_adaptive_nmin","Minimum number of constraints per connected components","none",pcbddc->adaptive_nmin,&pcbddc->adaptive_nmin,NULL);
67: PetscOptionsInt("-pc_bddc_adaptive_nmax","Maximum number of constraints per connected components","none",pcbddc->adaptive_nmax,&pcbddc->adaptive_nmax,NULL);
68: PetscOptionsBool("-pc_bddc_symmetric","Symmetric computation of primal basis functions","none",pcbddc->symmetric_primal,&pcbddc->symmetric_primal,NULL);
69: PetscOptionsInt("-pc_bddc_coarse_adj","Number of processors where to map the coarse adjacency list","none",pcbddc->coarse_adj_red,&pcbddc->coarse_adj_red,NULL);
70: PetscOptionsTail();
71: return(0);
72: }
73: /* -------------------------------------------------------------------------- */
76: static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change)
77: {
78: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
82: MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);
83: PetscObjectReference((PetscObject)change);
84: pcbddc->user_ChangeOfBasisMatrix = change;
85: return(0);
86: }
89: /*@
90: PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs
92: Collective on PC
94: Input Parameters:
95: + pc - the preconditioning context
96: - change - the change of basis matrix
98: Level: intermediate
100: Notes:
102: .seealso: PCBDDC
103: @*/
104: PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change)
105: {
112: if (pc->mat) {
113: PetscInt rows_c,cols_c,rows,cols;
114: MatGetSize(pc->mat,&rows,&cols);
115: MatGetSize(change,&rows_c,&cols_c);
116: if (rows_c != rows) {
117: SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of rows for change of basis matrix! %d != %d",rows_c,rows);
118: }
119: if (cols_c != cols) {
120: SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of columns for change of basis matrix! %d != %d",cols_c,cols);
121: }
122: MatGetLocalSize(pc->mat,&rows,&cols);
123: MatGetLocalSize(change,&rows_c,&cols_c);
124: if (rows_c != rows) {
125: SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local rows for change of basis matrix! %d != %d",rows_c,rows);
126: }
127: if (cols_c != cols) {
128: SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local columns for change of basis matrix! %d != %d",cols_c,cols);
129: }
130: }
131: PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat),(pc,change));
132: return(0);
133: }
134: /* -------------------------------------------------------------------------- */
137: static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
138: {
139: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
143: ISDestroy(&pcbddc->user_primal_vertices);
144: PetscObjectReference((PetscObject)PrimalVertices);
145: pcbddc->user_primal_vertices = PrimalVertices;
146: return(0);
147: }
150: /*@
151: PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
153: Collective
155: Input Parameters:
156: + pc - the preconditioning context
157: - PrimalVertices - index set of primal vertices in local numbering (can be empty)
159: Level: intermediate
161: Notes:
163: .seealso: PCBDDC
164: @*/
165: PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
166: {
173: PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));
174: return(0);
175: }
176: /* -------------------------------------------------------------------------- */
179: static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
180: {
181: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
184: pcbddc->coarsening_ratio = k;
185: return(0);
186: }
190: /*@
191: PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
193: Logically collective on PC
195: Input Parameters:
196: + pc - the preconditioning context
197: - k - coarsening ratio (H/h at the coarser level)
199: Options Database Keys:
200: . -pc_bddc_coarsening_ratio
202: Level: intermediate
204: Notes:
205: Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
207: .seealso: PCBDDC, PCBDDCSetLevels()
208: @*/
209: PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
210: {
216: PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));
217: return(0);
218: }
220: /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
223: static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
224: {
225: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
228: pcbddc->use_exact_dirichlet_trick = flg;
229: return(0);
230: }
234: PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
235: {
241: PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));
242: return(0);
243: }
247: static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
248: {
249: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
252: pcbddc->current_level = level;
253: return(0);
254: }
258: PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
259: {
265: PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));
266: return(0);
267: }
271: static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
272: {
273: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
276: pcbddc->max_levels = levels;
277: return(0);
278: }
282: /*@
283: PCBDDCSetLevels - Sets the maximum number of levels for multilevel
285: Logically collective on PC
287: Input Parameters:
288: + pc - the preconditioning context
289: - levels - the maximum number of levels (max 9)
291: Options Database Keys:
292: . -pc_bddc_levels
294: Level: intermediate
296: Notes:
297: Default value is 0, i.e. traditional one-level BDDC
299: .seealso: PCBDDC, PCBDDCSetCoarseningRatio()
300: @*/
301: PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
302: {
308: if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
309: PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));
310: return(0);
311: }
312: /* -------------------------------------------------------------------------- */
316: static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
317: {
318: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
322: PetscObjectReference((PetscObject)NullSpace);
323: MatNullSpaceDestroy(&pcbddc->NullSpace);
324: pcbddc->NullSpace = NullSpace;
325: return(0);
326: }
330: /*@
331: PCBDDCSetNullSpace - Set nullspace for BDDC operator
333: Logically collective on PC and MatNullSpace
335: Input Parameters:
336: + pc - the preconditioning context
337: - NullSpace - Null space of the linear operator to be preconditioned (Pmat)
339: Level: intermediate
341: Notes:
343: .seealso: PCBDDC
344: @*/
345: PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
346: {
353: PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));
354: return(0);
355: }
356: /* -------------------------------------------------------------------------- */
360: static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
361: {
362: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
366: /* last user setting takes precendence -> destroy any other customization */
367: ISDestroy(&pcbddc->DirichletBoundariesLocal);
368: ISDestroy(&pcbddc->DirichletBoundaries);
369: PetscObjectReference((PetscObject)DirichletBoundaries);
370: pcbddc->DirichletBoundaries = DirichletBoundaries;
371: pcbddc->recompute_topography = PETSC_TRUE;
372: return(0);
373: }
377: /*@
378: PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
380: Collective
382: Input Parameters:
383: + pc - the preconditioning context
384: - DirichletBoundaries - parallel IS defining the Dirichlet boundaries
386: Level: intermediate
388: Notes:
389: Provide the information if you used MatZeroRows/Columns routines. Any process can list any global node
391: .seealso: PCBDDC, PCBDDCSetDirichletBoundariesLocal(), MatZeroRows(), MatZeroRowsColumns()
392: @*/
393: PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
394: {
401: PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));
402: return(0);
403: }
404: /* -------------------------------------------------------------------------- */
408: static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
409: {
410: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
414: /* last user setting takes precendence -> destroy any other customization */
415: ISDestroy(&pcbddc->DirichletBoundariesLocal);
416: ISDestroy(&pcbddc->DirichletBoundaries);
417: PetscObjectReference((PetscObject)DirichletBoundaries);
418: pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
419: pcbddc->recompute_topography = PETSC_TRUE;
420: return(0);
421: }
425: /*@
426: PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
428: Collective
430: Input Parameters:
431: + pc - the preconditioning context
432: - DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
434: Level: intermediate
436: Notes:
438: .seealso: PCBDDC, PCBDDCSetDirichletBoundaries(), MatZeroRows(), MatZeroRowsColumns()
439: @*/
440: PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
441: {
448: PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));
449: return(0);
450: }
451: /* -------------------------------------------------------------------------- */
455: static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
456: {
457: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
461: /* last user setting takes precendence -> destroy any other customization */
462: ISDestroy(&pcbddc->NeumannBoundariesLocal);
463: ISDestroy(&pcbddc->NeumannBoundaries);
464: PetscObjectReference((PetscObject)NeumannBoundaries);
465: pcbddc->NeumannBoundaries = NeumannBoundaries;
466: pcbddc->recompute_topography = PETSC_TRUE;
467: return(0);
468: }
472: /*@
473: PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
475: Collective
477: Input Parameters:
478: + pc - the preconditioning context
479: - NeumannBoundaries - parallel IS defining the Neumann boundaries
481: Level: intermediate
483: Notes:
484: Any process can list any global node
486: .seealso: PCBDDC, PCBDDCSetNeumannBoundariesLocal()
487: @*/
488: PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
489: {
496: PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));
497: return(0);
498: }
499: /* -------------------------------------------------------------------------- */
503: static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
504: {
505: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
509: /* last user setting takes precendence -> destroy any other customization */
510: ISDestroy(&pcbddc->NeumannBoundariesLocal);
511: ISDestroy(&pcbddc->NeumannBoundaries);
512: PetscObjectReference((PetscObject)NeumannBoundaries);
513: pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
514: pcbddc->recompute_topography = PETSC_TRUE;
515: return(0);
516: }
520: /*@
521: PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
523: Collective
525: Input Parameters:
526: + pc - the preconditioning context
527: - NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
529: Level: intermediate
531: Notes:
533: .seealso: PCBDDC, PCBDDCSetNeumannBoundaries()
534: @*/
535: PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
536: {
543: PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));
544: return(0);
545: }
546: /* -------------------------------------------------------------------------- */
550: static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
551: {
552: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
555: *DirichletBoundaries = pcbddc->DirichletBoundaries;
556: return(0);
557: }
561: /*@
562: PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
564: Collective
566: Input Parameters:
567: . pc - the preconditioning context
569: Output Parameters:
570: . DirichletBoundaries - index set defining the Dirichlet boundaries
572: Level: intermediate
574: Notes:
575: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
577: .seealso: PCBDDC
578: @*/
579: PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
580: {
585: PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));
586: return(0);
587: }
588: /* -------------------------------------------------------------------------- */
592: static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
593: {
594: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
597: *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
598: return(0);
599: }
603: /*@
604: PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
606: Collective
608: Input Parameters:
609: . pc - the preconditioning context
611: Output Parameters:
612: . DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
614: Level: intermediate
616: Notes:
617: The IS returned could be the same passed in earlier by the user (if provided with PCBDDCSetDirichletBoundariesLocal) or a global-to-local map of the global IS (if provided with PCBDDCSetDirichletBoundaries).
618: In the latter case, the IS will be available after PCSetUp.
620: .seealso: PCBDDC
621: @*/
622: PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
623: {
628: PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));
629: return(0);
630: }
631: /* -------------------------------------------------------------------------- */
635: static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
636: {
637: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
640: *NeumannBoundaries = pcbddc->NeumannBoundaries;
641: return(0);
642: }
646: /*@
647: PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
649: Collective
651: Input Parameters:
652: . pc - the preconditioning context
654: Output Parameters:
655: . NeumannBoundaries - index set defining the Neumann boundaries
657: Level: intermediate
659: Notes:
660: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
662: .seealso: PCBDDC
663: @*/
664: PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
665: {
670: PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));
671: return(0);
672: }
673: /* -------------------------------------------------------------------------- */
677: static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
678: {
679: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
682: *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
683: return(0);
684: }
688: /*@
689: PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
691: Collective
693: Input Parameters:
694: . pc - the preconditioning context
696: Output Parameters:
697: . NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
699: Level: intermediate
701: Notes:
702: The IS returned could be the same passed in earlier by the user (if provided with PCBDDCSetNeumannBoundariesLocal) or a global-to-local map of the global IS (if provided with PCBDDCSetNeumannBoundaries).
703: In the latter case, the IS will be available after PCSetUp.
705: .seealso: PCBDDC
706: @*/
707: PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
708: {
713: PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));
714: return(0);
715: }
716: /* -------------------------------------------------------------------------- */
720: static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
721: {
722: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
723: PCBDDCGraph mat_graph = pcbddc->mat_graph;
727: /* free old CSR */
728: PCBDDCGraphResetCSR(mat_graph);
729: /* TODO: PCBDDCGraphSetAdjacency */
730: /* get CSR into graph structure */
731: if (copymode == PETSC_COPY_VALUES) {
732: PetscMalloc1(nvtxs+1,&mat_graph->xadj);
733: PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);
734: PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));
735: PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));
736: } else if (copymode == PETSC_OWN_POINTER) {
737: mat_graph->xadj = (PetscInt*)xadj;
738: mat_graph->adjncy = (PetscInt*)adjncy;
739: }
740: mat_graph->nvtxs_csr = nvtxs;
741: return(0);
742: }
746: /*@
747: PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local matrix
749: Not collective
751: Input Parameters:
752: + pc - the preconditioning context
753: . nvtxs - number of local vertices of the graph (i.e., the size of the local problem)
754: . xadj, adjncy - the CSR graph
755: - copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
757: Level: intermediate
759: Notes:
761: .seealso: PCBDDC,PetscCopyMode
762: @*/
763: PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
764: {
765: void (*f)(void) = 0;
772: if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
773: SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
774: }
775: PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));
776: /* free arrays if PCBDDC is not the PC type */
777: PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);
778: if (!f && copymode == PETSC_OWN_POINTER) {
779: PetscFree(xadj);
780: PetscFree(adjncy);
781: }
782: return(0);
783: }
784: /* -------------------------------------------------------------------------- */
788: static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
789: {
790: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
791: PetscInt i;
795: /* Destroy ISes if they were already set */
796: for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
797: ISDestroy(&pcbddc->ISForDofsLocal[i]);
798: }
799: PetscFree(pcbddc->ISForDofsLocal);
800: /* last user setting takes precendence -> destroy any other customization */
801: for (i=0;i<pcbddc->n_ISForDofs;i++) {
802: ISDestroy(&pcbddc->ISForDofs[i]);
803: }
804: PetscFree(pcbddc->ISForDofs);
805: pcbddc->n_ISForDofs = 0;
806: /* allocate space then set */
807: if (n_is) {
808: PetscMalloc1(n_is,&pcbddc->ISForDofsLocal);
809: }
810: for (i=0;i<n_is;i++) {
811: PetscObjectReference((PetscObject)ISForDofs[i]);
812: pcbddc->ISForDofsLocal[i]=ISForDofs[i];
813: }
814: pcbddc->n_ISForDofsLocal=n_is;
815: if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
816: pcbddc->recompute_topography = PETSC_TRUE;
817: return(0);
818: }
822: /*@
823: PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
825: Collective
827: Input Parameters:
828: + pc - the preconditioning context
829: . n_is - number of index sets defining the fields
830: - ISForDofs - array of IS describing the fields in local ordering
832: Level: intermediate
834: Notes:
835: n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
837: .seealso: PCBDDC
838: @*/
839: PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
840: {
841: PetscInt i;
847: for (i=0;i<n_is;i++) {
850: }
851: PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));
852: return(0);
853: }
854: /* -------------------------------------------------------------------------- */
858: static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
859: {
860: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
861: PetscInt i;
865: /* Destroy ISes if they were already set */
866: for (i=0;i<pcbddc->n_ISForDofs;i++) {
867: ISDestroy(&pcbddc->ISForDofs[i]);
868: }
869: PetscFree(pcbddc->ISForDofs);
870: /* last user setting takes precendence -> destroy any other customization */
871: for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
872: ISDestroy(&pcbddc->ISForDofsLocal[i]);
873: }
874: PetscFree(pcbddc->ISForDofsLocal);
875: pcbddc->n_ISForDofsLocal = 0;
876: /* allocate space then set */
877: if (n_is) {
878: PetscMalloc1(n_is,&pcbddc->ISForDofs);
879: }
880: for (i=0;i<n_is;i++) {
881: PetscObjectReference((PetscObject)ISForDofs[i]);
882: pcbddc->ISForDofs[i]=ISForDofs[i];
883: }
884: pcbddc->n_ISForDofs=n_is;
885: if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
886: pcbddc->recompute_topography = PETSC_TRUE;
887: return(0);
888: }
892: /*@
893: PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
895: Collective
897: Input Parameters:
898: + pc - the preconditioning context
899: . n_is - number of index sets defining the fields
900: - ISForDofs - array of IS describing the fields in global ordering
902: Level: intermediate
904: Notes:
905: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
907: .seealso: PCBDDC
908: @*/
909: PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
910: {
911: PetscInt i;
917: for (i=0;i<n_is;i++) {
920: }
921: PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));
922: return(0);
923: }
925: /* -------------------------------------------------------------------------- */
928: /* -------------------------------------------------------------------------- */
929: /*
930: PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
931: guess if a transformation of basis approach has been selected.
933: Input Parameter:
934: + pc - the preconditioner contex
936: Application Interface Routine: PCPreSolve()
938: Notes:
939: The interface routine PCPreSolve() is not usually called directly by
940: the user, but instead is called by KSPSolve().
941: */
942: static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
943: {
945: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
946: PC_IS *pcis = (PC_IS*)(pc->data);
947: Vec used_vec;
948: PetscBool copy_rhs = PETSC_TRUE;
951: /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
952: if (ksp) {
953: PetscBool iscg;
954: PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);
955: if (!iscg) {
956: PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);
957: }
958: }
959: /* Creates parallel work vectors used in presolve */
960: if (!pcbddc->original_rhs) {
961: VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);
962: }
963: if (!pcbddc->temp_solution) {
964: VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);
965: }
967: if (x) {
968: PetscObjectReference((PetscObject)x);
969: used_vec = x;
970: } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */
971: PetscObjectReference((PetscObject)pcbddc->temp_solution);
972: used_vec = pcbddc->temp_solution;
973: VecSet(used_vec,0.0);
974: }
976: /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */
977: if (ksp) {
978: /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */
979: KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);
980: if (!pcbddc->ksp_guess_nonzero) {
981: VecSet(used_vec,0.0);
982: }
983: }
985: pcbddc->rhs_change = PETSC_FALSE;
987: /* Take into account zeroed rows -> change rhs and store solution removed */
988: if (rhs) {
989: IS dirIS = NULL;
991: /* DirichletBoundariesLocal may not be consistent among neighbours; gets a dirichlet dofs IS from graph (may be cached) */
992: PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);
993: if (dirIS) {
994: Mat_IS *matis = (Mat_IS*)pc->pmat->data;
995: PetscInt dirsize,i,*is_indices;
996: PetscScalar *array_x;
997: const PetscScalar *array_diagonal;
999: MatGetDiagonal(pc->pmat,pcis->vec1_global);
1000: VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);
1001: VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);
1002: VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);
1003: VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
1004: VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
1005: ISGetLocalSize(dirIS,&dirsize);
1006: VecGetArray(pcis->vec1_N,&array_x);
1007: VecGetArrayRead(pcis->vec2_N,&array_diagonal);
1008: ISGetIndices(dirIS,(const PetscInt**)&is_indices);
1009: for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
1010: ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);
1011: VecRestoreArrayRead(pcis->vec2_N,&array_diagonal);
1012: VecRestoreArray(pcis->vec1_N,&array_x);
1013: VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);
1014: VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);
1015: pcbddc->rhs_change = PETSC_TRUE;
1016: ISDestroy(&dirIS);
1017: }
1018: }
1020: /* remove the computed solution or the initial guess from the rhs */
1021: if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) {
1022: /* store the original rhs */
1023: if (copy_rhs) {
1024: VecCopy(rhs,pcbddc->original_rhs);
1025: copy_rhs = PETSC_FALSE;
1026: }
1027: pcbddc->rhs_change = PETSC_TRUE;
1028: VecScale(used_vec,-1.0);
1029: MatMultAdd(pc->pmat,used_vec,rhs,rhs);
1030: VecScale(used_vec,-1.0);
1031: VecCopy(used_vec,pcbddc->temp_solution);
1032: if (ksp) {
1033: KSPSetInitialGuessNonzero(ksp,PETSC_FALSE);
1034: }
1035: }
1036: VecDestroy(&used_vec);
1038: /* store partially computed solution and set initial guess */
1039: if (x && pcbddc->use_exact_dirichlet_trick) {
1040: VecSet(x,0.0);
1041: VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1042: VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1043: KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);
1044: VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);
1045: VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);
1046: if (ksp) {
1047: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
1048: }
1049: }
1051: if (pcbddc->ChangeOfBasisMatrix) {
1052: PCBDDCChange_ctx change_ctx;
1054: /* get change ctx */
1055: MatShellGetContext(pcbddc->new_global_mat,&change_ctx);
1057: /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */
1058: MatDestroy(&change_ctx->original_mat);
1059: PetscObjectReference((PetscObject)pc->mat);
1060: change_ctx->original_mat = pc->mat;
1062: /* change iteration matrix */
1063: MatDestroy(&pc->mat);
1064: PetscObjectReference((PetscObject)pcbddc->new_global_mat);
1065: pc->mat = pcbddc->new_global_mat;
1067: /* store the original rhs */
1068: if (copy_rhs) {
1069: VecCopy(rhs,pcbddc->original_rhs);
1070: copy_rhs = PETSC_FALSE;
1071: }
1073: /* change rhs */
1074: MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);
1075: VecCopy(pcis->vec1_global,rhs);
1076: pcbddc->rhs_change = PETSC_TRUE;
1077: }
1079: /* remove nullspace if present */
1080: if (ksp && x && pcbddc->NullSpace) {
1081: MatNullSpaceRemove(pcbddc->NullSpace,x);
1082: /* store the original rhs */
1083: if (copy_rhs) {
1084: VecCopy(rhs,pcbddc->original_rhs);
1085: copy_rhs = PETSC_FALSE;
1086: }
1087: pcbddc->rhs_change = PETSC_TRUE;
1088: MatNullSpaceRemove(pcbddc->NullSpace,rhs);
1089: }
1090: return(0);
1091: }
1093: /* -------------------------------------------------------------------------- */
1096: /* -------------------------------------------------------------------------- */
1097: /*
1098: PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1099: approach has been selected. Also, restores rhs to its original state.
1101: Input Parameter:
1102: + pc - the preconditioner contex
1104: Application Interface Routine: PCPostSolve()
1106: Notes:
1107: The interface routine PCPostSolve() is not usually called directly by
1108: the user, but instead is called by KSPSolve().
1109: */
1110: static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1111: {
1113: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
1116: if (pcbddc->ChangeOfBasisMatrix) {
1117: PCBDDCChange_ctx change_ctx;
1119: /* get change ctx */
1120: MatShellGetContext(pcbddc->new_global_mat,&change_ctx);
1122: /* restore iteration matrix */
1123: MatDestroy(&pc->mat);
1124: PetscObjectReference((PetscObject)change_ctx->original_mat);
1125: pc->mat = change_ctx->original_mat;
1127: /* get solution in original basis */
1128: if (x) {
1129: PC_IS *pcis = (PC_IS*)(pc->data);
1130: MatMult(change_ctx->global_change,x,pcis->vec1_global);
1131: VecCopy(pcis->vec1_global,x);
1132: }
1133: }
1135: /* add solution removed in presolve */
1136: if (x && pcbddc->rhs_change) {
1137: VecAXPY(x,1.0,pcbddc->temp_solution);
1138: }
1140: /* restore rhs to its original state */
1141: if (rhs && pcbddc->rhs_change) {
1142: VecCopy(pcbddc->original_rhs,rhs);
1143: }
1144: pcbddc->rhs_change = PETSC_FALSE;
1146: /* restore ksp guess state */
1147: if (ksp) {
1148: KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);
1149: }
1150: return(0);
1151: }
1152: /* -------------------------------------------------------------------------- */
1155: /* -------------------------------------------------------------------------- */
1156: /*
1157: PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
1158: by setting data structures and options.
1160: Input Parameter:
1161: + pc - the preconditioner context
1163: Application Interface Routine: PCSetUp()
1165: Notes:
1166: The interface routine PCSetUp() is not usually called directly by
1167: the user, but instead is called by PCApply() if necessary.
1168: */
1169: PetscErrorCode PCSetUp_BDDC(PC pc)
1170: {
1172: PC_BDDC* pcbddc = (PC_BDDC*)pc->data;
1173: Mat_IS* matis;
1174: MatNullSpace nearnullspace;
1175: PetscBool computetopography,computesolvers,computesubschurs;
1176: PetscBool computeconstraintsmatrix;
1177: PetscBool new_nearnullspace_provided,ismatis;
1180: PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);
1181: if (!ismatis) {
1182: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS");
1183: }
1184: matis = (Mat_IS*)pc->pmat->data;
1185: /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
1186: /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1187: Also, BDDC directly build the Dirichlet problem */
1188: /* split work */
1189: if (pc->setupcalled) {
1190: if (pc->flag == SAME_NONZERO_PATTERN) {
1191: computetopography = PETSC_FALSE;
1192: computesolvers = PETSC_TRUE;
1193: } else { /* DIFFERENT_NONZERO_PATTERN */
1194: computetopography = PETSC_TRUE;
1195: computesolvers = PETSC_TRUE;
1196: }
1197: } else {
1198: computetopography = PETSC_TRUE;
1199: computesolvers = PETSC_TRUE;
1200: }
1201: if (pcbddc->recompute_topography) {
1202: computetopography = PETSC_TRUE;
1203: }
1204: computeconstraintsmatrix = PETSC_FALSE;
1205: if (pcbddc->adaptive_threshold > 0.0 && !pcbddc->use_deluxe_scaling) {
1206: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute adaptive constraints without deluxe scaling. Rerun with -pc_bddc_use_deluxe_scaling");
1207: }
1208: pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold > 0.0 && pcbddc->use_deluxe_scaling);
1209: if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE;
1211: computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling);
1212: if (pcbddc->faster_deluxe && pcbddc->adaptive_selection && pcbddc->use_change_of_basis) {
1213: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute faster deluxe if adaptivity and change of basis are both requested. Rerun with -pc_bddc_deluxe_faster false");
1214: }
1215: /* Get stdout for dbg */
1216: if (pcbddc->dbg_flag) {
1217: if (!pcbddc->dbg_viewer) {
1218: pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1219: PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);
1220: }
1221: PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);
1222: }
1224: if (pcbddc->user_ChangeOfBasisMatrix) {
1225: /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
1226: pcbddc->use_change_of_basis = PETSC_FALSE;
1227: PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);
1228: } else {
1229: MatDestroy(&pcbddc->local_mat);
1230: PetscObjectReference((PetscObject)matis->A);
1231: pcbddc->local_mat = matis->A;
1232: }
1234: /* workaround for reals */
1235: #if !defined(PETSC_USE_COMPLEX)
1236: if (matis->A->symmetric_set) {
1237: MatSetOption(pcbddc->local_mat,MAT_HERMITIAN,matis->A->symmetric);
1238: }
1239: #endif
1241: /* Set up all the "iterative substructuring" common block without computing solvers */
1242: {
1243: Mat temp_mat;
1245: temp_mat = matis->A;
1246: matis->A = pcbddc->local_mat;
1247: PCISSetUp(pc,PETSC_FALSE);
1248: pcbddc->local_mat = matis->A;
1249: matis->A = temp_mat;
1250: }
1252: /* Analyze interface and setup sub_schurs data */
1253: if (computetopography) {
1254: PCBDDCAnalyzeInterface(pc);
1255: computeconstraintsmatrix = PETSC_TRUE;
1256: }
1258: /* Setup local dirichlet solver ksp_D and sub_schurs solvers */
1259: if (computesolvers) {
1260: PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs;
1262: if (computesubschurs && computetopography) {
1263: PCBDDCInitSubSchurs(pc);
1264: }
1265: if (sub_schurs->use_mumps) {
1266: if (computesubschurs) {
1267: PCBDDCSetUpSubSchurs(pc);
1268: }
1269: PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);
1270: } else {
1271: PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);
1272: if (computesubschurs) {
1273: PCBDDCSetUpSubSchurs(pc);
1274: }
1275: }
1276: if (pcbddc->adaptive_selection) {
1277: PCBDDCAdaptiveSelection(pc);
1278: computeconstraintsmatrix = PETSC_TRUE;
1279: }
1280: }
1282: /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1283: new_nearnullspace_provided = PETSC_FALSE;
1284: MatGetNearNullSpace(pc->pmat,&nearnullspace);
1285: if (pcbddc->onearnullspace) { /* already used nearnullspace */
1286: if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1287: new_nearnullspace_provided = PETSC_TRUE;
1288: } else {
1289: /* determine if the two nullspaces are different (should be lightweight) */
1290: if (nearnullspace != pcbddc->onearnullspace) {
1291: new_nearnullspace_provided = PETSC_TRUE;
1292: } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1293: PetscInt i;
1294: const Vec *nearnullvecs;
1295: PetscObjectState state;
1296: PetscInt nnsp_size;
1297: MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);
1298: for (i=0;i<nnsp_size;i++) {
1299: PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);
1300: if (pcbddc->onearnullvecs_state[i] != state) {
1301: new_nearnullspace_provided = PETSC_TRUE;
1302: break;
1303: }
1304: }
1305: }
1306: }
1307: } else {
1308: if (!nearnullspace) { /* both nearnullspaces are null */
1309: new_nearnullspace_provided = PETSC_FALSE;
1310: } else { /* nearnullspace attached later */
1311: new_nearnullspace_provided = PETSC_TRUE;
1312: }
1313: }
1315: /* Setup constraints and related work vectors */
1316: /* reset primal space flags */
1317: pcbddc->new_primal_space = PETSC_FALSE;
1318: pcbddc->new_primal_space_local = PETSC_FALSE;
1319: if (computeconstraintsmatrix || new_nearnullspace_provided) {
1320: /* It also sets the primal space flags */
1321: PCBDDCConstraintsSetUp(pc);
1322: /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1323: PCBDDCSetUpLocalWorkVectors(pc);
1324: }
1326: if (computesolvers || pcbddc->new_primal_space) {
1327: if (pcbddc->use_change_of_basis) {
1328: PC_IS *pcis = (PC_IS*)(pc->data);
1330: PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);
1331: /* get submatrices */
1332: MatDestroy(&pcis->A_IB);
1333: MatDestroy(&pcis->A_BI);
1334: MatDestroy(&pcis->A_BB);
1335: MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);
1336: MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);
1337: MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);
1338: /* set flag in pcis to not reuse submatrices during PCISCreate */
1339: pcis->reusesubmatrices = PETSC_FALSE;
1340: } else if (!pcbddc->user_ChangeOfBasisMatrix) {
1341: MatDestroy(&pcbddc->local_mat);
1342: PetscObjectReference((PetscObject)matis->A);
1343: pcbddc->local_mat = matis->A;
1344: }
1345: /* SetUp coarse and local Neumann solvers */
1346: PCBDDCSetUpSolvers(pc);
1347: /* SetUp Scaling operator */
1348: PCBDDCScalingSetUp(pc);
1349: }
1351: if (pcbddc->dbg_flag) {
1352: PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);
1353: }
1354: return(0);
1355: }
1357: /* -------------------------------------------------------------------------- */
1358: /*
1359: PCApply_BDDC - Applies the BDDC operator to a vector.
1361: Input Parameters:
1362: + pc - the preconditioner context
1363: - r - input vector (global)
1365: Output Parameter:
1366: . z - output vector (global)
1368: Application Interface Routine: PCApply()
1369: */
1372: PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
1373: {
1374: PC_IS *pcis = (PC_IS*)(pc->data);
1375: PC_BDDC *pcbddc = (PC_BDDC*)(pc->data);
1376: PetscInt n_B = pcis->n_B, n_D = pcis->n - n_B;
1377: PetscErrorCode ierr;
1378: const PetscScalar one = 1.0;
1379: const PetscScalar m_one = -1.0;
1380: const PetscScalar zero = 0.0;
1382: /* This code is similar to that provided in nn.c for PCNN
1383: NN interface preconditioner changed to BDDC
1384: Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */
1387: if (!pcbddc->use_exact_dirichlet_trick) {
1388: VecCopy(r,z);
1389: /* First Dirichlet solve */
1390: VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1391: VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1392: /*
1393: Assembling right hand side for BDDC operator
1394: - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1395: - pcis->vec1_B the interface part of the global vector z
1396: */
1397: if (n_D) {
1398: KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);
1399: VecScale(pcis->vec2_D,m_one);
1400: if (pcbddc->switch_static) { MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D); }
1401: MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);
1402: } else {
1403: VecSet(pcis->vec1_B,zero);
1404: }
1405: VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1406: VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1407: PCBDDCScalingRestriction(pc,z,pcis->vec1_B);
1408: } else {
1409: if (pcbddc->switch_static) {
1410: VecSet(pcis->vec1_D,zero);
1411: }
1412: PCBDDCScalingRestriction(pc,r,pcis->vec1_B);
1413: }
1415: /* Apply interface preconditioner
1416: input/output vecs: pcis->vec1_B and pcis->vec1_D */
1417: PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);
1419: /* Apply transpose of partition of unity operator */
1420: PCBDDCScalingExtension(pc,pcis->vec1_B,z);
1422: /* Second Dirichlet solve and assembling of output */
1423: VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1424: VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1425: if (n_B) {
1426: MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);
1427: if (pcbddc->switch_static) { MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D); }
1428: } else if (pcbddc->switch_static) {
1429: MatMult(pcis->A_II,pcis->vec1_D,pcis->vec3_D);
1430: }
1431: KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);
1432: if (!pcbddc->use_exact_dirichlet_trick) {
1433: if (pcbddc->switch_static) {
1434: VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);
1435: } else {
1436: VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);
1437: }
1438: VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
1439: VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
1440: } else {
1441: if (pcbddc->switch_static) {
1442: VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);
1443: } else {
1444: VecScale(pcis->vec4_D,m_one);
1445: }
1446: VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);
1447: VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);
1448: }
1449: return(0);
1450: }
1452: /* -------------------------------------------------------------------------- */
1453: /*
1454: PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
1456: Input Parameters:
1457: + pc - the preconditioner context
1458: - r - input vector (global)
1460: Output Parameter:
1461: . z - output vector (global)
1463: Application Interface Routine: PCApplyTranspose()
1464: */
1467: PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
1468: {
1469: PC_IS *pcis = (PC_IS*)(pc->data);
1470: PC_BDDC *pcbddc = (PC_BDDC*)(pc->data);
1471: PetscInt n_B = pcis->n_B, n_D = pcis->n - n_B;
1472: PetscErrorCode ierr;
1473: const PetscScalar one = 1.0;
1474: const PetscScalar m_one = -1.0;
1475: const PetscScalar zero = 0.0;
1478: if (!pcbddc->use_exact_dirichlet_trick) {
1479: VecCopy(r,z);
1480: /* First Dirichlet solve */
1481: VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1482: VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1483: /*
1484: Assembling right hand side for BDDC operator
1485: - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1486: - pcis->vec1_B the interface part of the global vector z
1487: */
1488: if (n_D) {
1489: KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);
1490: VecScale(pcis->vec2_D,m_one);
1491: if (pcbddc->switch_static) { MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D); }
1492: MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);
1493: } else {
1494: VecSet(pcis->vec1_B,zero);
1495: }
1496: VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1497: VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1498: PCBDDCScalingRestriction(pc,z,pcis->vec1_B);
1499: } else {
1500: if (pcbddc->switch_static) {
1501: VecSet(pcis->vec1_D,zero);
1502: }
1503: PCBDDCScalingRestriction(pc,r,pcis->vec1_B);
1504: }
1506: /* Apply interface preconditioner
1507: input/output vecs: pcis->vec1_B and pcis->vec1_D */
1508: PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);
1510: /* Apply transpose of partition of unity operator */
1511: PCBDDCScalingExtension(pc,pcis->vec1_B,z);
1513: /* Second Dirichlet solve and assembling of output */
1514: VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1515: VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1516: if (n_B) {
1517: MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);
1518: if (pcbddc->switch_static) { MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D); }
1519: } else if (pcbddc->switch_static) {
1520: MatMultTranspose(pcis->A_II,pcis->vec1_D,pcis->vec3_D);
1521: }
1522: KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);
1523: if (!pcbddc->use_exact_dirichlet_trick) {
1524: if (pcbddc->switch_static) {
1525: VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);
1526: } else {
1527: VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);
1528: }
1529: VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
1530: VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
1531: } else {
1532: if (pcbddc->switch_static) {
1533: VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);
1534: } else {
1535: VecScale(pcis->vec4_D,m_one);
1536: }
1537: VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);
1538: VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);
1539: }
1540: return(0);
1541: }
1542: /* -------------------------------------------------------------------------- */
1546: PetscErrorCode PCDestroy_BDDC(PC pc)
1547: {
1548: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
1552: /* free data created by PCIS */
1553: PCISDestroy(pc);
1554: /* free BDDC custom data */
1555: PCBDDCResetCustomization(pc);
1556: /* destroy objects related to topography */
1557: PCBDDCResetTopography(pc);
1558: /* free allocated graph structure */
1559: PetscFree(pcbddc->mat_graph);
1560: /* free allocated sub schurs structure */
1561: PetscFree(pcbddc->sub_schurs);
1562: /* destroy objects for scaling operator */
1563: PCBDDCScalingDestroy(pc);
1564: PetscFree(pcbddc->deluxe_ctx);
1565: /* free solvers stuff */
1566: PCBDDCResetSolvers(pc);
1567: /* free global vectors needed in presolve */
1568: VecDestroy(&pcbddc->temp_solution);
1569: VecDestroy(&pcbddc->original_rhs);
1570: /* free stuff for change of basis hooks */
1571: if (pcbddc->new_global_mat) {
1572: PCBDDCChange_ctx change_ctx;
1573: MatShellGetContext(pcbddc->new_global_mat,&change_ctx);
1574: MatDestroy(&change_ctx->original_mat);
1575: MatDestroy(&change_ctx->global_change);
1576: VecDestroyVecs(2,&change_ctx->work);
1577: PetscFree(change_ctx);
1578: }
1579: MatDestroy(&pcbddc->new_global_mat);
1580: /* remove functions */
1581: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);
1582: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);
1583: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);
1584: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);
1585: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);
1586: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);
1587: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);
1588: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);
1589: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);
1590: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);
1591: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);
1592: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);
1593: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);
1594: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);
1595: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);
1596: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);
1597: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);
1598: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);
1599: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);
1600: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);
1601: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);
1602: /* Free the private data structure */
1603: PetscFree(pc->data);
1604: return(0);
1605: }
1606: /* -------------------------------------------------------------------------- */
1610: static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1611: {
1612: FETIDPMat_ctx mat_ctx;
1613: Vec copy_standard_rhs;
1614: PC_IS* pcis;
1615: PC_BDDC* pcbddc;
1619: MatShellGetContext(fetidp_mat,&mat_ctx);
1620: pcis = (PC_IS*)mat_ctx->pc->data;
1621: pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1623: /*
1624: change of basis for physical rhs if needed
1625: It also changes the rhs in case of dirichlet boundaries
1626: TODO: better management when FETIDP will have its own class
1627: */
1628: VecDuplicate(standard_rhs,©_standard_rhs);
1629: VecCopy(standard_rhs,copy_standard_rhs);
1630: PCPreSolve_BDDC(mat_ctx->pc,NULL,copy_standard_rhs,NULL);
1631: /* store vectors for computation of fetidp final solution */
1632: VecScatterBegin(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);
1633: VecScatterEnd(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);
1634: /* scale rhs since it should be unassembled */
1635: /* TODO use counter scaling? (also below) */
1636: VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
1637: VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
1638: /* Apply partition of unity */
1639: VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);
1640: /* PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B); */
1641: if (!pcbddc->switch_static) {
1642: /* compute partially subassembled Schur complement right-hand side */
1643: KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);
1644: MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);
1645: VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);
1646: VecSet(copy_standard_rhs,0.0);
1647: VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);
1648: VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);
1649: /* PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B); */
1650: VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
1651: VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
1652: VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);
1653: }
1654: VecDestroy(©_standard_rhs);
1655: /* BDDC rhs */
1656: VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);
1657: if (pcbddc->switch_static) {
1658: VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);
1659: }
1660: /* apply BDDC */
1661: PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);
1662: /* Application of B_delta and assembling of rhs for fetidp fluxes */
1663: VecSet(fetidp_flux_rhs,0.0);
1664: MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);
1665: VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);
1666: VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);
1667: return(0);
1668: }
1672: /*@
1673: PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETI-DP linear system using the physical right-hand side
1675: Collective
1677: Input Parameters:
1678: + fetidp_mat - the FETI-DP matrix object obtained by a call to PCBDDCCreateFETIDPOperators
1679: - standard_rhs - the right-hand side of the original linear system
1681: Output Parameters:
1682: . fetidp_flux_rhs - the right-hand side for the FETI-DP linear system
1684: Level: developer
1686: Notes:
1688: .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetSolution
1689: @*/
1690: PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
1691: {
1692: FETIDPMat_ctx mat_ctx;
1696: MatShellGetContext(fetidp_mat,&mat_ctx);
1697: PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));
1698: return(0);
1699: }
1700: /* -------------------------------------------------------------------------- */
1704: static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1705: {
1706: FETIDPMat_ctx mat_ctx;
1707: PC_IS* pcis;
1708: PC_BDDC* pcbddc;
1712: MatShellGetContext(fetidp_mat,&mat_ctx);
1713: pcis = (PC_IS*)mat_ctx->pc->data;
1714: pcbddc = (PC_BDDC*)mat_ctx->pc->data;
1716: /* apply B_delta^T */
1717: VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
1718: VecScatterEnd (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
1719: MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);
1720: /* compute rhs for BDDC application */
1721: VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);
1722: if (pcbddc->switch_static) {
1723: VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);
1724: }
1725: /* apply BDDC */
1726: PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);
1727: /* put values into standard global vector */
1728: VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);
1729: VecScatterEnd (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);
1730: if (!pcbddc->switch_static) {
1731: /* compute values into the interior if solved for the partially subassembled Schur complement */
1732: MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);
1733: VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);
1734: KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);
1735: }
1736: VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);
1737: VecScatterEnd (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);
1738: /* final change of basis if needed
1739: Is also sums the dirichlet part removed during RHS assembling */
1740: PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);
1741: return(0);
1742: }
1746: /*@
1747: PCBDDCMatFETIDPGetSolution - Compute the physical solution using the solution of the FETI-DP linear system
1749: Collective
1751: Input Parameters:
1752: + fetidp_mat - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators
1753: - fetidp_flux_sol - the solution of the FETI-DP linear system
1755: Output Parameters:
1756: . standard_sol - the solution defined on the physical domain
1758: Level: developer
1760: Notes:
1762: .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetRHS
1763: @*/
1764: PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
1765: {
1766: FETIDPMat_ctx mat_ctx;
1770: MatShellGetContext(fetidp_mat,&mat_ctx);
1771: PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));
1772: return(0);
1773: }
1774: /* -------------------------------------------------------------------------- */
1776: extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1777: extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
1778: extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1779: extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1780: extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
1781: extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1785: static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1786: {
1788: FETIDPMat_ctx fetidpmat_ctx;
1789: Mat newmat;
1790: FETIDPPC_ctx fetidppc_ctx;
1791: PC newpc;
1792: MPI_Comm comm;
1796: PetscObjectGetComm((PetscObject)pc,&comm);
1797: /* FETIDP linear matrix */
1798: PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);
1799: PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);
1800: MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);
1801: MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);
1802: MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);
1803: MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);
1804: MatSetUp(newmat);
1805: /* FETIDP preconditioner */
1806: PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);
1807: PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);
1808: PCCreate(comm,&newpc);
1809: PCSetType(newpc,PCSHELL);
1810: PCShellSetContext(newpc,fetidppc_ctx);
1811: PCShellSetApply(newpc,FETIDPPCApply);
1812: PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);
1813: PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);
1814: PCSetOperators(newpc,newmat,newmat);
1815: PCSetUp(newpc);
1816: /* return pointers for objects created */
1817: *fetidp_mat=newmat;
1818: *fetidp_pc=newpc;
1819: return(0);
1820: }
1824: /*@
1825: PCBDDCCreateFETIDPOperators - Create FETI-DP operators
1827: Collective
1829: Input Parameters:
1830: . pc - the BDDC preconditioning context (setup should have been called before)
1832: Output Parameters:
1833: + fetidp_mat - shell FETI-DP matrix object
1834: - fetidp_pc - shell Dirichlet preconditioner for FETI-DP matrix
1836: Options Database Keys:
1837: . -fetidp_fullyredundant <false> - use or not a fully redundant set of Lagrange multipliers
1839: Level: developer
1841: Notes:
1842: Currently the only operations provided for FETI-DP matrix are MatMult and MatMultTranspose
1844: .seealso: PCBDDC, PCBDDCMatFETIDPGetRHS, PCBDDCMatFETIDPGetSolution
1845: @*/
1846: PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
1847: {
1852: if (pc->setupcalled) {
1853: PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));
1854: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
1855: return(0);
1856: }
1857: /* -------------------------------------------------------------------------- */
1858: /*MC
1859: PCBDDC - Balancing Domain Decomposition by Constraints.
1861: An implementation of the BDDC preconditioner based on
1863: .vb
1864: [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
1865: [2] A. Klawonn and O. B. Widlund. "Dual-Primal FETI Methods for Linear Elasticity", http://cs.nyu.edu/csweb/Research/TechReports/TR2004-855/TR2004-855.pdf
1866: [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
1867: [4] C. Pechstein and C. R. Dohrmann. "Modern domain decomposition methods BDDC, deluxe scaling, and an algebraic approach", Seminar talk, Linz, December 2013, http://people.ricam.oeaw.ac.at/c.pechstein/pechstein-bddc2013.pdf
1868: .ve
1870: The matrix to be preconditioned (Pmat) must be of type MATIS.
1872: Currently works with MATIS matrices with local matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
1874: It also works with unsymmetric and indefinite problems.
1876: Unlike 'conventional' interface preconditioners, PCBDDC iterates over all degrees of freedom, not just those on the interface. This allows the use of approximate solvers on the subdomains.
1878: Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace()
1880: Boundary nodes are split in vertices, edges and faces classes using information from the local to global mapping of dofs and the local connectivity graph of nodes. The latter can be customized by using PCBDDCSetLocalAdjacencyGraph()
1881: Additional information on dofs can be provided by using PCBDDCSetDofsSplitting(), PCBDDCSetDirichletBoundaries(), PCBDDCSetNeumannBoundaries(), and PCBDDCSetPrimalVerticesLocalIS()
1883: Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). Non-singular modes are retained via SVD.
1885: Change of basis is performed similarly to [2] when requested. When more than one constraint is present on a single connected component (i.e. an edge or a face), a robust method based on local QR factorizations is used.
1886: User defined change of basis can be passed to PCBDDC by using PCBDDCSetChangeOfBasisMat()
1888: The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using a MatPartitioning object.
1890: Adaptive selection of primal constraints [4] is supported for SPD systems with high-contrast in the coefficients if MUMPS is present. Future versions of the code will also consider using MKL_PARDISO or PASTIX.
1892: An experimental interface to the FETI-DP method is available. FETI-DP operators could be created using PCBDDCCreateFETIDPOperators(). A stand-alone class for the FETI-DP method will be provided in the next releases.
1893: Deluxe scaling is not supported yet for FETI-DP.
1895: Options Database Keys (some of them, run with -h for a complete list):
1897: . -pc_bddc_use_vertices <true> - use or not vertices in primal space
1898: . -pc_bddc_use_edges <true> - use or not edges in primal space
1899: . -pc_bddc_use_faces <false> - use or not faces in primal space
1900: . -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems
1901: . -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only)
1902: . -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested
1903: . -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1])
1904: . -pc_bddc_levels <0> - maximum number of levels for multilevel
1905: . -pc_bddc_coarsening_ratio <8> - number of subdomains which will be aggregated together at the coarser level (e.g. H/h ratio at the coarser level, significative only in the multilevel case)
1906: . -pc_bddc_redistribute <0> - size of a subset of processors where the coarse problem will be remapped (the value is ignored if not at the coarsest level)
1907: . -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling
1908: . -pc_bddc_schur_layers <-1> - select the economic version of deluxe scaling by specifying the number of layers (-1 corresponds to the original deluxe scaling)
1909: . -pc_bddc_adaptive_threshold <0.0> - when a value greater than one is specified, adaptive selection of constraints is performed on edges and faces (requires deluxe scaling and MUMPS installed)
1910: - -pc_bddc_check_level <0> - set verbosity level of debugging output
1912: Options for Dirichlet, Neumann or coarse solver can be set with
1913: .vb
1914: -pc_bddc_dirichlet_
1915: -pc_bddc_neumann_
1916: -pc_bddc_coarse_
1917: .ve
1918: e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. PCBDDC uses by default KPSPREONLY and PCLU.
1920: When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified as
1921: .vb
1922: -pc_bddc_dirichlet_lN_
1923: -pc_bddc_neumann_lN_
1924: -pc_bddc_coarse_lN_
1925: .ve
1926: Note that level number ranges from the finest (0) to the coarsest (N).
1927: In order to specify options for the BDDC operators at the coarser levels (and not for the solvers), prepend -pc_bddc_coarse_ or -pc_bddc_coarse_l to the option, e.g.
1928: .vb
1929: -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3
1930: .ve
1931: will use a threshold of 5 for constraints' selection at the first coarse level and will redistribute the coarse problem of the first coarse level on 3 processors
1933: Level: intermediate
1935: Developer notes:
1937: Contributed by Stefano Zampini
1939: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS
1940: M*/
1944: PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1945: {
1946: PetscErrorCode ierr;
1947: PC_BDDC *pcbddc;
1950: /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1951: PetscNewLog(pc,&pcbddc);
1952: pc->data = (void*)pcbddc;
1954: /* create PCIS data structure */
1955: PCISCreate(pc);
1957: /* BDDC customization */
1958: pcbddc->use_local_adj = PETSC_TRUE;
1959: pcbddc->use_vertices = PETSC_TRUE;
1960: pcbddc->use_edges = PETSC_TRUE;
1961: pcbddc->use_faces = PETSC_FALSE;
1962: pcbddc->use_change_of_basis = PETSC_FALSE;
1963: pcbddc->use_change_on_faces = PETSC_FALSE;
1964: pcbddc->switch_static = PETSC_FALSE;
1965: pcbddc->use_nnsp_true = PETSC_FALSE;
1966: pcbddc->use_qr_single = PETSC_FALSE;
1967: pcbddc->symmetric_primal = PETSC_TRUE;
1968: pcbddc->dbg_flag = 0;
1969: /* private */
1970: pcbddc->local_primal_size = 0;
1971: pcbddc->local_primal_size_cc = 0;
1972: pcbddc->local_primal_ref_node = 0;
1973: pcbddc->local_primal_ref_mult = 0;
1974: pcbddc->n_vertices = 0;
1975: pcbddc->primal_indices_local_idxs = 0;
1976: pcbddc->recompute_topography = PETSC_FALSE;
1977: pcbddc->coarse_size = -1;
1978: pcbddc->new_primal_space = PETSC_FALSE;
1979: pcbddc->new_primal_space_local = PETSC_FALSE;
1980: pcbddc->global_primal_indices = 0;
1981: pcbddc->onearnullspace = 0;
1982: pcbddc->onearnullvecs_state = 0;
1983: pcbddc->user_primal_vertices = 0;
1984: pcbddc->NullSpace = 0;
1985: pcbddc->temp_solution = 0;
1986: pcbddc->original_rhs = 0;
1987: pcbddc->local_mat = 0;
1988: pcbddc->ChangeOfBasisMatrix = 0;
1989: pcbddc->user_ChangeOfBasisMatrix = 0;
1990: pcbddc->new_global_mat = 0;
1991: pcbddc->coarse_vec = 0;
1992: pcbddc->coarse_ksp = 0;
1993: pcbddc->coarse_phi_B = 0;
1994: pcbddc->coarse_phi_D = 0;
1995: pcbddc->coarse_psi_B = 0;
1996: pcbddc->coarse_psi_D = 0;
1997: pcbddc->vec1_P = 0;
1998: pcbddc->vec1_R = 0;
1999: pcbddc->vec2_R = 0;
2000: pcbddc->local_auxmat1 = 0;
2001: pcbddc->local_auxmat2 = 0;
2002: pcbddc->R_to_B = 0;
2003: pcbddc->R_to_D = 0;
2004: pcbddc->ksp_D = 0;
2005: pcbddc->ksp_R = 0;
2006: pcbddc->NeumannBoundaries = 0;
2007: pcbddc->NeumannBoundariesLocal = 0;
2008: pcbddc->DirichletBoundaries = 0;
2009: pcbddc->DirichletBoundariesLocal = 0;
2010: pcbddc->user_provided_isfordofs = PETSC_FALSE;
2011: pcbddc->n_ISForDofs = 0;
2012: pcbddc->n_ISForDofsLocal = 0;
2013: pcbddc->ISForDofs = 0;
2014: pcbddc->ISForDofsLocal = 0;
2015: pcbddc->ConstraintMatrix = 0;
2016: pcbddc->use_exact_dirichlet_trick = PETSC_TRUE;
2017: pcbddc->coarse_loc_to_glob = 0;
2018: pcbddc->coarsening_ratio = 8;
2019: pcbddc->coarse_adj_red = 0;
2020: pcbddc->current_level = 0;
2021: pcbddc->max_levels = 0;
2022: pcbddc->use_coarse_estimates = PETSC_FALSE;
2023: pcbddc->redistribute_coarse = 0;
2024: pcbddc->coarse_subassembling = 0;
2025: pcbddc->coarse_subassembling_init = 0;
2027: /* create local graph structure */
2028: PCBDDCGraphCreate(&pcbddc->mat_graph);
2030: /* scaling */
2031: pcbddc->work_scaling = 0;
2032: pcbddc->use_deluxe_scaling = PETSC_FALSE;
2033: pcbddc->faster_deluxe = PETSC_FALSE;
2035: /* create sub schurs structure */
2036: PCBDDCSubSchursCreate(&pcbddc->sub_schurs);
2037: pcbddc->sub_schurs_rebuild = PETSC_FALSE;
2038: pcbddc->sub_schurs_layers = -1;
2039: pcbddc->sub_schurs_use_useradj = PETSC_FALSE;
2041: pcbddc->computed_rowadj = PETSC_FALSE;
2043: /* adaptivity */
2044: pcbddc->adaptive_threshold = 0.0;
2045: pcbddc->adaptive_nmax = 0;
2046: pcbddc->adaptive_nmin = 0;
2048: /* function pointers */
2049: pc->ops->apply = PCApply_BDDC;
2050: pc->ops->applytranspose = PCApplyTranspose_BDDC;
2051: pc->ops->setup = PCSetUp_BDDC;
2052: pc->ops->destroy = PCDestroy_BDDC;
2053: pc->ops->setfromoptions = PCSetFromOptions_BDDC;
2054: pc->ops->view = 0;
2055: pc->ops->applyrichardson = 0;
2056: pc->ops->applysymmetricleft = 0;
2057: pc->ops->applysymmetricright = 0;
2058: pc->ops->presolve = PCPreSolve_BDDC;
2059: pc->ops->postsolve = PCPostSolve_BDDC;
2061: /* composing function */
2062: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);
2063: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);
2064: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);
2065: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);
2066: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);
2067: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);
2068: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);
2069: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);
2070: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);
2071: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);
2072: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);
2073: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);
2074: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);
2075: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);
2076: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);
2077: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);
2078: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);
2079: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);
2080: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);
2081: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);
2082: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);
2083: return(0);
2084: }