Flow123d  jenkins-Flow123d-windows-release-multijob-285
schur.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @ingroup la
27  * @brief Assembly explicit Schur complement for the given linear system.
28  * Provides method for resolution of the full original vector of unknowns.
29  *
30  * Aim: Explicit schur should be faster then implicit, i.e.
31  *
32  * @todo
33  * - vyresit navaznost na lin sys - solve a export seq vektoru, redukce ... ?
34  * - inv_a - predava se pri konstrukci, ale neumoznuje jeji reuse - aktualizaci assemblace
35  * resp. nutno si na ni drzet ukazatel venku ... reseni ?
36  * - ? remove old_4_new - just for LSView
37  * - automatic preallocation
38  * - eliminated block given by IS
39  * - in place Schur
40  * - ? nemodifikovat puvodni system, leda skrze jeho metody
41  */
42 
43 #include <petscvec.h>
44 #include <algorithm>
45 #include <limits>
46 #include <petscmat.h>
47 #include <armadillo>
48 #include <petscis.h>
49 
50 #include "la/distribution.hh"
52 #include "system/system.hh"
53 #include "la/linsys.hh"
54 #include "la/linsys_BDDC.hh"
55 #include "la/schur.hh"
56 
57 /**
58  * Create Schur complement system.
59  * @param[in] orig : original system
60  * @param[in] inv_a : inversion of the A block
61  * @param[in] ia : index set of the A block,
62  * default continuous given by inv_a:
63  * proc 1 2 3
64  *
65  * Orig: ****** ****** ****
66  * IA : *** ** ***
67  *
68  *
69  */
70 
72 : LinSys_PETSC(ds), IsA(ia), state(created)
73 {
74  // check index set
75  ASSERT(IsA != NULL, "Index set IsA is not defined.\n" );
76 
77  // initialize variables
78  Compl = NULL;
79  IA = NULL;
80  B = NULL;
81  Bt = NULL;
82  xA = NULL;
83  IAB = NULL;
84  IsB = NULL;
85  RHS1 = NULL;
86  RHS2 = NULL;
87  Sol1 = NULL;
88  Sol2 = NULL;
89 
90  // create A block index set
91  ISGetLocalSize(IsA, &loc_size_A);
92 
93  // create B block index set
95  ISCreateStride(PETSC_COMM_WORLD,loc_size_B,rows_ds_->begin()+loc_size_A,1,&IsB);
96 }
97 
98 
100 : LinSys_PETSC(other),
101  loc_size_A(other.loc_size_A), loc_size_B(other.loc_size_B), state(other.state),
102  Compl(other.Compl), ds_(other.ds_)
103 {
104  MatCopy(other.IA, IA, DIFFERENT_NONZERO_PATTERN);
105  MatCopy(other.IAB, IAB, DIFFERENT_NONZERO_PATTERN);
106  ISCopy(other.IsA, IsA);
107  ISCopy(other.IsB, IsB);
108  VecCopy(other.RHS1, RHS1);
109  VecCopy(other.RHS2, RHS2);
110  VecCopy(other.Sol1, Sol1);
111  VecCopy(other.Sol2, Sol2);
112 
113  B = NULL;
114  Bt = NULL;
115  xA = NULL;
116 }
117 
118 
119 /**
120  * COMPUTE A SCHUR COMPLEMENT OF A PETSC MATRIX
121  *
122  * given symmetric original matrix Orig has form
123  * A B x_1 RHS_1
124  * B' C * x_2 = RHS_2
125  * where the first block is given by index set IsA, and the second block by IsB
126  * user has to provide inverse IA of the A-block
127  * we suppose that original matrix have non-zero pattern for the schur complement
128  *
129  * we return: Shur - schur complement, ShurRHS - RHS of the complemented system:
130  * (B' * IA * B - C) * x_2 = (B' * IA * RHS_1 - RHS_2)
131  * IAB - a matrix to compute eliminated part of the solution:
132  * x_1 = IA * RHS_1 - IAB * x_2
133  *
134  * Actually as B' is stored separetly, the routine can be used also for nonsymetric original
135  * system
136  *
137  */
138 
140 {
141  PetscErrorCode ierr = 0;
142  MatReuse mat_reuse; // reuse structures after first computation of schur
143  PetscScalar *rhs_array, *sol_array;
144 
145  mat_reuse=MAT_REUSE_MATRIX;
146  if (state==created) {
147  mat_reuse=MAT_INITIAL_MATRIX; // indicate first construction
148 
149  // create complement system
150  // TODO: introduce LS as true object, clarify its internal states
151  // create RHS sub vecs RHS1, RHS2
152  VecGetArray(rhs_, &rhs_array);
153  VecCreateMPIWithArray(PETSC_COMM_WORLD,1,loc_size_A,PETSC_DETERMINE,rhs_array,&(RHS1));
154 
155  // create Solution sub vecs Sol1, Compl->solution
156  VecGetArray(solution_, &sol_array);
157  VecCreateMPIWithArray(PETSC_COMM_WORLD,1,loc_size_A,PETSC_DETERMINE,sol_array,&(Sol1));
158 
159  VecCreateMPIWithArray(PETSC_COMM_WORLD,1,loc_size_B,PETSC_DETERMINE,rhs_array+loc_size_A,&(RHS2));
160  VecCreateMPIWithArray(PETSC_COMM_WORLD,1,loc_size_B,PETSC_DETERMINE,sol_array+loc_size_A,&(Sol2));
161 
162  VecRestoreArray(rhs_, &rhs_array);
163  VecRestoreArray(solution_, &sol_array);
164 
165  VecGetArray( Sol2, &sol_array );
166  Compl->set_solution( sol_array );
168  VecRestoreArray( Sol2, &sol_array );
169 
170  }
171 
172  // compose Schur complement
173  // Petsc need some fill estimate for results of multiplication in form nnz(A*B)/(nnz(A)+nnz(B))
174  // for the first Schur compl: IA*B is bounded by ( d*(d+1) )/( d*d+2*d ) <= 5/6 for d<=4
175  // B'*IA*B bounded by ( (d+1)*(d+1) )/ ( d*(d+1) + d ) ~ 1
176  // for the second Schur : IA*B have fill ratio ~ 1.
177  // B'*IA*B ... ( N/2 *(2*N-1) )/( 2 + 2*N ) <= 1.4
178  // nevertheless Petsc does not allows fill ratio below 1. so we use 1.1 for the first
179  // and 1.5 for the second multiplication
180 
181  if (matrix_changed_) {
183 
184  // compute IAB=IA*B, loc_size_B removed
185  ierr+=MatGetSubMatrix(matrix_, IsA, IsB, mat_reuse, &B);
186  ierr+=MatMatMult(IA, B, mat_reuse, 1.0 ,&(IAB)); // 6/7 - fill estimate
187  // compute xA=Bt* IAB = Bt * IA * B, locSizeA removed
188  ierr+=MatGetSubMatrix(matrix_, IsB, IsA, mat_reuse, &(Bt));
189  ierr+=MatMatMult(Bt, IAB, mat_reuse, 1.9 ,&(xA)); // 1.1 - fill estimate (PETSC report values over 1.8)
190 
191  // get C block, loc_size_B removed
192  ierr+=MatGetSubMatrix( matrix_, IsB, IsB, mat_reuse, const_cast<Mat *>( Compl->get_matrix() ) );
193  // compute complement = (-1)cA+xA = Bt*IA*B - C
194  if ( is_negative_definite() ) {
195  ierr+=MatAXPY(*( Compl->get_matrix() ), -1, xA, SUBSET_NONZERO_PATTERN);
196  } else {
197  ierr+=MatScale(*( Compl->get_matrix() ),-1.0);
198  ierr+=MatAXPY(*( Compl->get_matrix() ), 1, xA, SUBSET_NONZERO_PATTERN);
199  }
201 
202  ASSERT( ierr == 0, "PETSC Error during calculation of Schur complement.\n");
203 
204  }
205 
206  form_rhs();
207 
208  matrix_changed_ = false;
209 
210  state=formed;
211 }
212 
214 {
215  if (rhs_changed_ || matrix_changed_) {
216  MatMultTranspose(IAB, RHS1, *( Compl->get_rhs() ));
217  VecAXPY(*( Compl->get_rhs() ), -1, RHS2);
218  if ( is_negative_definite() ) {
219  VecScale(*( Compl->get_rhs() ), -1.0);
220  }
222  rhs_changed_ = false;
223  }
224 
225  state=formed;
226 }
227 
228 
229 /**
230  * COMPUTE ELIMINATED PART OF THE ORIG. SYS. & RESTORE RHS and SOLUTION VECTORS
231  * x_1 = IA * RHS_1 - IAB * x_2
232  */
233 
235 {
236  MatMult(IAB,Compl->get_solution(),Sol1);
237 
238  VecScale(Sol1,-1);
239 
240  MatMultAdd(IA,RHS1,Sol1,Sol1);
241 
242 }
243 
245 {
246  ASSERT(ls != nullptr, "NULL complement ls.\n");
247  Compl = ls;
248 }
249 
250 
252 {
253  ds_ = new Distribution(loc_size_B, PETSC_COMM_WORLD);
254  return ds_;
255 }
256 
258 {
259  PetscInt ncols, pos_start, pos_start_IA;
260 
261  MatReuse mat_reuse=MAT_REUSE_MATRIX;
262  if (state==created) mat_reuse=MAT_INITIAL_MATRIX; // indicate first construction
263 
264  MatGetSubMatrix(matrix_, IsA, IsA, mat_reuse, &IA);
265  MatGetOwnershipRange(matrix_,&pos_start,PETSC_NULL);
266  MatGetOwnershipRange(IA,&pos_start_IA,PETSC_NULL);
267 
268  std::vector<PetscInt> submat_rows;
269  const PetscInt *cols;
270  const PetscScalar *vals;
271 
272  std::vector<unsigned int> processed_rows(loc_size_A,0);
273 
274  unsigned int mat_block=1; //actual processed block of matrix
275  for(unsigned int loc_row=0; loc_row < processed_rows.size(); loc_row++) {
276  if (processed_rows[loc_row] != 0) continue;
277 
278  PetscInt min=std::numeric_limits<int>::max(), max=-1, size_submat;
279  PetscInt b_vals = 0; // count of values stored in B-block of Orig system
280  submat_rows.clear();
281  MatGetRow(matrix_, loc_row + pos_start, &ncols, &cols, PETSC_NULL);
282  for (PetscInt i=0; i<ncols; i++) {
283  if (cols[i] < pos_start || cols[i] >= pos_start+loc_size_A) {
284  b_vals++;
285  } else {
286  if (cols[i] < min) {
287  min=cols[i];
288  }
289  if (cols[i] > max) {
290  max=cols[i];
291  }
292  }
293  }
294  size_submat = max - min + 1;
295  ASSERT(ncols-b_vals == size_submat, "Submatrix cannot contains empty values.\n");
296 
297  MatRestoreRow(matrix_, loc_row + pos_start, &ncols, &cols, PETSC_NULL);
298  arma::mat submat2(size_submat, size_submat);
299  submat2.zeros();
300  for (PetscInt i=0; i<size_submat; i++) {
301  processed_rows[ loc_row + i ] = mat_block;
302  submat_rows.push_back( i + loc_row + pos_start_IA );
303  MatGetRow(matrix_, i + loc_row + pos_start, &ncols, &cols, &vals);
304  for (PetscInt j=0; j<ncols; j++) {
305  if (cols[j] >= pos_start && cols[j] < pos_start+loc_size_A) {
306  submat2( i, cols[j] - loc_row - pos_start ) = vals[j];
307  }
308  }
309  MatRestoreRow(matrix_, i + loc_row + pos_start, &ncols, &cols, &vals);
310  }
311  // get inversion matrix
312  arma::mat invmat = submat2.i();
313  // stored to inversion IA matrix
314  const PetscInt* rows = &submat_rows[0];
315  MatSetValues(IA, submat_rows.size(), rows, submat_rows.size(), rows, invmat.memptr(), INSERT_VALUES);
316 
317  mat_block++;
318  }
319 
320  MatAssemblyBegin(IA, MAT_FINAL_ASSEMBLY);
321  MatAssemblyEnd(IA, MAT_FINAL_ASSEMBLY);
322 }
323 
324 
326 {
327  if (Compl != NULL) {
328  return Compl->get_solution_precision();
329  }
330  return std::numeric_limits<double>::infinity();
331 }
332 
333 
335  this->form_schur();
336 
337  int converged_reason = Compl->solve();
338  this->resolve();
339 
340  return converged_reason;
341 }
342 
343 
344 /**
345  * SCHUR COMPLEMENT destructor
346  */
348 
349  if ( B != NULL ) MatDestroy(&B);
350  if ( Bt != NULL ) MatDestroy(&Bt);
351  if ( xA != NULL ) MatDestroy(&xA);
352  if ( IA != NULL ) MatDestroy(&IA);
353  if ( IAB != NULL ) MatDestroy(&IAB);
354  if ( IsA != NULL ) ISDestroy(&IsA);
355  if ( IsB != NULL ) ISDestroy(&IsB);
356  if ( RHS1 != NULL ) VecDestroy(&RHS1);
357  if ( RHS2 != NULL ) VecDestroy(&RHS2);
358  if ( Sol1 != NULL ) VecDestroy(&Sol1);
359  if ( Sol2 != NULL ) VecDestroy(&Sol2);
360  if ( IA != NULL ) MatDestroy(&IA);
361 
362  if (Compl != NULL) delete Compl;
363 
364 }
const Vec * get_rhs()
Definition: linsys_PETSC.hh:68
bool is_negative_definite()
Definition: linsys.hh:507
void set_complement(LinSys_PETSC *ls)
Set complement LinSys object.
Definition: schur.cc:244
void resolve()
Definition: schur.cc:234
bool matrix_changed_
true if the matrix was changed since the last solve
Definition: linsys.hh:588
Wrappers for linear systems based on MPIAIJ and MATIS format.
const Mat * get_matrix()
Definition: linsys_PETSC.hh:63
Definition: schur.hh:68
void set_rhs_changed()
Definition: linsys.hh:201
Vec rhs_
PETSc vector constructed with vx array.
LinSys_PETSC * Compl
Definition: schur.hh:144
void set_from_input(const Input::Record in_rec)
int loc_size_B
Definition: schur.hh:136
Input::Record in_rec_
Definition: linsys.hh:601
bool rhs_changed_
true if the right hand side was changed since the last solve
Definition: linsys.hh:589
Distribution * ds_
Definition: schur.hh:146
~SchurComplement()
Definition: schur.cc:347
void form_rhs()
Definition: schur.cc:213
Assembly explicit Schur complement for the given linear system. Provides method for resolution of the...
SchurState state
Definition: schur.hh:141
#define ASSERT(...)
Definition: global_defs.h:121
const Vec & get_solution()
Definition: linsys.hh:256
unsigned int begin(int proc) const
get starting local index
const Distribution * rows_ds_
final distribution of rows of MH matrix
Definition: linsys.hh:581
void form_schur()
Definition: schur.cc:139
int loc_size_A
Definition: schur.hh:136
int solve() override
Definition: schur.cc:334
void set_solution(double *sol_array)
Definition: linsys.hh:264
SchurComplement(IS ia, Distribution *ds)
Definition: schur.cc:71
void create_inversion_matrix()
create IA matrix
Definition: schur.cc:257
Support classes for parallel programing.
double get_solution_precision() override
get precision of solving
Definition: schur.cc:325
Distribution * make_complement_distribution()
get distribution of complement object if complement is defined
Definition: schur.cc:251
Definition: schur.hh:69
double get_solution_precision()
void set_matrix_changed()
Definition: linsys.hh:195
Mat matrix_
Petsc matrix of the problem.
Vec solution_
PETSc vector constructed with vb array.
Definition: linsys.hh:591
Solver based on Multilevel BDDC - using corresponding class of OpenFTL package.
unsigned int lsize(int proc) const
get local size