QR_MUMPS
 All Classes Files Functions Variables Enumerations Enumerator Pages
qrm_analysis_mod.F90
Go to the documentation of this file.
1 !! ##############################################################################################
2 !!
3 !! Copyright 2012 CNRS, INPT
4 !!
5 !! This file is part of qr_mumps.
6 !!
7 !! qr_mumps is free software: you can redistribute it and/or modify
8 !! it under the terms of the GNU Lesser General Public License as
9 !! published by the Free Software Foundation, either version 3 of
10 !! the License, or (at your option) any later version.
11 !!
12 !! qr_mumps is distributed in the hope that it will be useful,
13 !! but WITHOUT ANY WARRANTY; without even the implied warranty of
14 !! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 !! GNU Lesser General Public License for more details.
16 !!
17 !! You can find a copy of the GNU Lesser General Public License
18 !! in the qr_mumps/doc directory.
19 !!
20 !! ##############################################################################################
21 
22 
23 !! ##############################################################################################
33 
34 
37 
39  interface qrm_analyse
40  subroutine _qrm_analyse(qrm_mat, transp)
41  use _qrm_spmat_mod
42  type(_qrm_spmat_type) :: qrm_mat
43  character, optional, intent(in) :: transp
44  end subroutine _qrm_analyse
45  end interface
46 
49  subroutine _qrm_compute_graph(qrm_mat, graph, work)
50  use _qrm_spmat_mod
51  type(_qrm_spmat_type) :: qrm_mat
52  type(_qrm_spmat_type), intent(out) :: graph
53  integer, target, optional :: work(:)
54  end subroutine _qrm_compute_graph
55  end interface
56 
59  subroutine _qrm_detect_singletons(graph, scol, srow, mrperm, mcperm, nrsing, ncsing)
60  use _qrm_spmat_mod
61  type(_qrm_spmat_type), intent(in) :: graph
62  integer, allocatable :: scol(:), srow(:)
63  integer, allocatable :: mcperm(:), mrperm(:)
64  integer :: nrsing, ncsing
65  end subroutine _qrm_detect_singletons
66  end interface
67 
69  interface qrm_do_colamd
70  subroutine _qrm_do_colamd(graph, cperm)
71  use _qrm_spmat_mod
72  type(_qrm_spmat_type) :: graph
73  integer :: cperm(:)
74  end subroutine _qrm_do_colamd
75  end interface
76 
78  interface qrm_do_metis
79  subroutine _qrm_do_metis(graph, cperm)
80  use _qrm_spmat_mod
81  type(_qrm_spmat_type) :: graph
82  integer :: cperm(:)
83  end subroutine _qrm_do_metis
84  end interface
85 
87  interface qrm_do_scotch
88  subroutine _qrm_do_scotch(graph, cperm)
89  use _qrm_spmat_mod
90  type(_qrm_spmat_type) :: graph
91  integer :: cperm(:)
92  end subroutine _qrm_do_scotch
93  end interface
94 
96  interface qrm_ata_graph
97  subroutine _qrm_ata_graph(g_csc, ata_graph)
98  use _qrm_spmat_mod
99  type(_qrm_spmat_type), intent(in) :: g_csc
100  type(_qrm_spmat_type), intent(out) :: ata_graph
101  end subroutine _qrm_ata_graph
102  end interface
103 
105  interface qrm_do_ordering
106  subroutine _qrm_do_ordering(graph, cperm, cperm_in)
107  use _qrm_spmat_mod
108  type(_qrm_spmat_type) :: graph
109  integer :: cperm(:)
110  integer, pointer :: cperm_in(:)
111  end subroutine _qrm_do_ordering
112  end interface
113 
115  interface qrm_elim_tree
116  subroutine _qrm_elim_tree(graph, cperm, parent)
117  use _qrm_spmat_mod
118  type(_qrm_spmat_type), intent(in) :: graph
119  integer, intent(in) :: cperm(:)
120  integer, allocatable :: parent(:)
121  end subroutine _qrm_elim_tree
122  end interface
123 
125  interface qrm_rowcount
126  subroutine _qrm_rowcount(graph, parent, porder, rc)
127  use _qrm_spmat_mod
128  type(_qrm_spmat_type) :: graph
129  integer :: parent(:), porder(:), rc(:)
130  end subroutine _qrm_rowcount
131  end interface
132 
134  interface qrm_symbolic
135  subroutine _qrm_symbolic(graph)
136  use _qrm_spmat_mod
137  type(_qrm_spmat_type), intent(inout) :: graph
138  end subroutine _qrm_symbolic
139  end interface
140 
142  interface qrm_rowperm
143  subroutine _qrm_rowperm(qrm_mat, cperm, rperm, nvar, stair)
144  use _qrm_spmat_mod
145  type(_qrm_spmat_type) :: qrm_mat
146  integer :: cperm(:), rperm(:), nvar(:), stair(:)
147  end subroutine _qrm_rowperm
148  end interface
149 
152  subroutine _qrm_attach_singletons(qrm_mat, scol, srow, mrperm, mcperm, nrsing, ncsing)
153  use _qrm_spmat_mod
154  type(_qrm_spmat_type):: qrm_mat
155  integer :: scol(:), srow(:), mcperm(:), mrperm(:)
156  integer :: nrsing, ncsing
157  end subroutine _qrm_attach_singletons
158  end interface
159 
160 end module _qrm_analysis_mod
161 
Generic interface for the ::_qrm_ata_graph routine.
Generic interface for the ::_qrm_do_scotch routine.
subroutine _qrm_do_colamd(graph, cperm)
This subroutine computes the fill reducing ordering using COLAMD.
Generic interface for the ::_qrm_do_metis routine.
Generic interface for the ::_qrm_symbolic routine.
Generic interface for the ::_qrm_do_ordering routine.
This module contains the generic interfaces for all the analysis routines.
subroutine _qrm_do_ordering(graph, cperm, cperm_in)
This routine computes (through different methods) a column permutation of the input matrix in order t...
Generic interface for the ::_qrm_do_colamd routine.
subroutine _qrm_attach_singletons(qrm_mat, scol, srow, mrperm, mcperm, nrsing, ncsing)
This subroutine merges the results of the singletons detection into the results of the analysis phase...
subroutine _qrm_rowperm(graph, cperm, rperm, nvar, stair)
This routine computes a row permutation that puts the input matrix into a staircase format...
Definition: qrm_rowperm.F90:68
subroutine _qrm_detect_singletons(graph, scol, srow, mrperm, mcperm, nrsing, ncsing)
This subroutine detects singletons in a matrix.
Generic interface for the ::_qrm_compute_graph routine.
Generic interface for the ::_qrm_detect_singletons routine.
Generic interface for the ::_qrm_rowcount routine.
subroutine _qrm_analyse(qrm_mat, transp)
This is the driver routine for the analysis phase.
Definition: qrm_analyse.F90:64
Generic interface for the ::_qrm_rowperm routine.
subroutine _qrm_compute_graph(qrm_mat, graph)
Computes the adjacency graph of a matrix.
Generic interface for the ::_qrm_attach_singletons routine.
subroutine _qrm_do_metis(graph, cperm)
Please refer to:
Generic interface for the ::_qrm_elim_tree routine.
subroutine _qrm_do_scotch(graph, cperm)
Please refer to:
subroutine _qrm_rowcount(graph, parent, porder, rc)
This subroutine computes the rowcount of the R factor.
This type defines the data structure used to store a matrix.
This module contains the definition of the basic sparse matrix type and of the associated methods...
Generic interface for the ::_qrm_analyse routine.
subroutine _qrm_elim_tree(graph, cperm, parent)
This subroutine builds the elimination tree for A'A.
subroutine _qrm_symbolic(graph)
This subroutine computes the symbolic QR factorization of a matrix.
subroutine _qrm_ata_graph(g_csc, ata_graph)
This subroutine computes the fill reducing ordering using METIS.