Flow123d  release_3.0.0-958-ga07bd24
schur.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file schur.cc
15  * @ingroup la
16  * @brief Assembly explicit Schur complement for the given linear system.
17  * Provides method for resolution of the full original vector of unknowns.
18  *
19  * Aim: Explicit schur should be faster then implicit, i.e.
20  *
21  * @todo
22  * - vyresit navaznost na lin sys - solve a export seq vektoru, redukce ... ?
23  * - inv_a - predava se pri konstrukci, ale neumoznuje jeji reuse - aktualizaci assemblace
24  * resp. nutno si na ni drzet ukazatel venku ... reseni ?
25  * - ? remove old_4_new - just for LSView
26  * - automatic preallocation
27  * - eliminated block given by IS
28  * - in place Schur
29  * - ? nemodifikovat puvodni system, leda skrze jeho metody
30  */
31 
32 #include <petscvec.h>
33 #include <algorithm>
34 #include <limits>
35 #include <petscmat.h>
36 #include <armadillo>
37 #include <petscis.h>
38 
39 #include "system/sys_profiler.hh"
40 #include "la/distribution.hh"
42 #include "system/system.hh"
43 #include "la/linsys.hh"
44 #include "la/linsys_BDDC.hh"
45 #include "la/schur.hh"
46 
47 /**
48  * Create Schur complement system.
49  * @param[in] orig : original system
50  * @param[in] inv_a : inversion of the A block
51  * @param[in] ia : index set of the A block,
52  * default continuous given by inv_a:
53  * proc 1 2 3
54  *
55  * Orig: ****** ****** ****
56  * IA : *** ** ***
57  *
58  *
59  */
60 
62 : LinSys_PETSC(ds), IsA(ia), state(created)
63 {
64  // check index set
65  OLD_ASSERT(IsA != NULL, "Index set IsA is not defined.\n" );
66 
67  // initialize variables
68  Compl = NULL;
69  IA = NULL;
70  B = NULL;
71  Bt = NULL;
72  C = NULL;
73  xA = NULL;
74  IAB = NULL;
75  IsB = NULL;
76  RHS1 = NULL;
77  RHS2 = NULL;
78  Sol1 = NULL;
79  Sol2 = NULL;
80 
81  // create A block index set
82  ISGetLocalSize(IsA, &loc_size_A);
83 
84  // create B block index set
86  ISCreateStride(PETSC_COMM_WORLD,loc_size_B,rows_ds_->begin()+loc_size_A,1,&IsB);
87 }
88 
89 
91 : LinSys_PETSC(other),
92  loc_size_A(other.loc_size_A), loc_size_B(other.loc_size_B), state(other.state),
93  Compl(other.Compl), ds_(other.ds_)
94 {
95  MatCopy(other.IA, IA, DIFFERENT_NONZERO_PATTERN);
96  MatCopy(other.IAB, IAB, DIFFERENT_NONZERO_PATTERN);
97  ISCopy(other.IsA, IsA);
98  ISCopy(other.IsB, IsB);
99  VecCopy(other.RHS1, RHS1);
100  VecCopy(other.RHS2, RHS2);
101  VecCopy(other.Sol1, Sol1);
102  VecCopy(other.Sol2, Sol2);
103 
104  B = NULL;
105  Bt = NULL;
106  xA = NULL;
107  C = NULL;
108 }
109 
110 
111 
113 {
115 
116  ASSERT_PTR(Compl).error();
117  Compl->set_from_input( in_rec );
118 }
119 
120 /**
121  * COMPUTE A SCHUR COMPLEMENT OF A PETSC MATRIX
122  *
123  * given symmetric original matrix Orig has form
124  * A B x_1 RHS_1
125  * B' C * x_2 = RHS_2
126  * where the first block is given by index set IsA, and the second block by IsB
127  * user has to provide inverse IA of the A-block
128  * we suppose that original matrix have non-zero pattern for the schur complement
129  *
130  * we return: Shur - schur complement, ShurRHS - RHS of the complemented system:
131  * (B' * IA * B - C) * x_2 = (B' * IA * RHS_1 - RHS_2)
132  * IAB - a matrix to compute eliminated part of the solution:
133  * x_1 = IA * RHS_1 - IAB * x_2
134  *
135  * Actually as B' is stored separetly, the routine can be used also for nonsymetric original
136  * system
137  *
138  */
139 
141 {
142  START_TIMER("form schur complement");
143 
144 
145  PetscErrorCode ierr = 0;
146  MatReuse mat_reuse; // reuse structures after first computation of schur
147  MatStructure mat_subset_pattern;
148  PetscScalar *rhs_array, *sol_array;
149 
150  mat_reuse=MAT_REUSE_MATRIX;
151  mat_subset_pattern=SUBSET_NONZERO_PATTERN;
152  if (state==created) {
153  mat_reuse=MAT_INITIAL_MATRIX; // indicate first construction
154  mat_subset_pattern=DIFFERENT_NONZERO_PATTERN;
155 
156  // create complement system
157  // TODO: introduce LS as true object, clarify its internal states
158  // create RHS sub vecs RHS1, RHS2
159  VecGetArray(rhs_, &rhs_array);
160  VecCreateMPIWithArray(PETSC_COMM_WORLD,1,loc_size_A,PETSC_DETERMINE,rhs_array,&(RHS1));
161 
162  // create Solution sub vecs Sol1, Compl->solution
163  VecGetArray(solution_, &sol_array);
164  VecCreateMPIWithArray(PETSC_COMM_WORLD,1,loc_size_A,PETSC_DETERMINE,sol_array,&(Sol1));
165 
166  VecCreateMPIWithArray(PETSC_COMM_WORLD,1,loc_size_B,PETSC_DETERMINE,rhs_array+loc_size_A,&(RHS2));
167  VecCreateMPIWithArray(PETSC_COMM_WORLD,1,loc_size_B,PETSC_DETERMINE,sol_array+loc_size_A,&(Sol2));
168 
169  VecRestoreArray(rhs_, &rhs_array);
170  VecRestoreArray(solution_, &sol_array);
171 
172  VecGetArray( Sol2, &sol_array );
173  Compl->set_solution( sol_array );
174 
175  VecRestoreArray( Sol2, &sol_array );
176 
177 
178  }
179  DebugOut() << print_var(mat_reuse) << print_var(matrix_changed_) << print_var(state);
180 
181  // compose Schur complement
182  // Petsc need some fill estimate for results of multiplication in form nnz(A*B)/(nnz(A)+nnz(B))
183  // for the first Schur compl: IA*B is bounded by ( d*(d+1) )/( d*d+2*d ) <= 5/6 for d<=4
184  // B'*IA*B bounded by ( (d+1)*(d+1) )/ ( d*(d+1) + d ) ~ 1
185  // for the second Schur : IA*B have fill ratio ~ 1.
186  // B'*IA*B ... ( N/2 *(2*N-1) )/( 2 + 2*N ) <= 1.4
187  // nevertheless Petsc does not allows fill ratio below 1. so we use 1.1 for the first
188  // and 1.5 for the second multiplication
189 
190  // TODO:
191  // In order to let PETSC allocate structure of the complement we can not perform MatGetSubMatrix
192  // and MatAXPY on the same complement matrix, since one operation change the matrix structure for the other.
193  //
194  // Probably no way to make this optimal using high level methods. We should have our own
195  // format for schur complement matrix, store local systems and perform elimination localy.
196  // Or even better assembly the complement directly. (not compatible with raw P0 method)
197 
198  if (matrix_changed_) {
200 
201  // compute IAB=IA*B, loc_size_B removed
202  ierr+=MatGetSubMatrix(matrix_, IsA, IsB, mat_reuse, &B);
203  ierr+=MatMatMult(IA, B, mat_reuse, 1.0 ,&(IAB)); // 6/7 - fill estimate
204  // compute xA=Bt* IAB = Bt * IA * B, locSizeA removed
205  ierr+=MatGetSubMatrix(matrix_, IsB, IsA, mat_reuse, &(Bt));
206  ierr+=MatMatMult(Bt, IAB, mat_reuse, 1.9 ,&(xA)); // 1.1 - fill estimate (PETSC report values over 1.8)
207 
208  // get C block, loc_size_B removed
209  ierr+=MatGetSubMatrix( matrix_, IsB, IsB, mat_reuse, &C);
210 
211  if (state==created) MatDuplicate(C, MAT_DO_NOT_COPY_VALUES, const_cast<Mat *>( Compl->get_matrix() ) );
212  MatZeroEntries( *( Compl->get_matrix()) );
213 
214  // compute complement = (-1)cA+xA = Bt*IA*B - C
215  if ( is_negative_definite() ) {
216  ierr+=MatAXPY(*( Compl->get_matrix() ), 1, C, SUBSET_NONZERO_PATTERN);
217  ierr+=MatAXPY(*( Compl->get_matrix() ), -1, xA, mat_subset_pattern);
218  } else {
219  ierr+=MatAXPY(*( Compl->get_matrix() ), -1, C, SUBSET_NONZERO_PATTERN);
220  ierr+=MatAXPY(*( Compl->get_matrix() ), 1, xA, mat_subset_pattern);
221  }
223 
224  OLD_ASSERT( ierr == 0, "PETSC Error during calculation of Schur complement.\n");
225 
226  }
227 
228  form_rhs();
229 
230  matrix_changed_ = false;
231 
232  state=formed;
233 }
234 
236 {
237  START_TIMER("form rhs");
238  if (rhs_changed_ || matrix_changed_) {
239  MatMultTranspose(IAB, RHS1, *( Compl->get_rhs() ));
240  VecAXPY(*( Compl->get_rhs() ), -1, RHS2);
241  if ( is_negative_definite() ) {
242  VecScale(*( Compl->get_rhs() ), -1.0);
243  }
245  rhs_changed_ = false;
246  }
247 
248  state=formed;
249 }
250 
251 
252 
253 void SchurComplement::set_tolerances(double r_tol, double a_tol, unsigned int max_it)
254 {
255  LinSys_PETSC::set_tolerances(r_tol, a_tol, max_it);
256  if (Compl !=nullptr) Compl->set_tolerances(r_tol, a_tol, max_it);
257 }
258 
260 {
261  ASSERT_PTR(ls).error();
262  Compl = ls;
263 }
264 
265 
267 {
268  ds_ = new Distribution(loc_size_B, PETSC_COMM_WORLD);
269  return ds_;
270 }
271 
273 {
274  START_TIMER("create inversion matrix");
275  PetscInt ncols, pos_start, pos_start_IA;
276 
277  MatReuse mat_reuse=MAT_REUSE_MATRIX;
278  if (state==created) mat_reuse=MAT_INITIAL_MATRIX; // indicate first construction
279 
280  MatGetSubMatrix(matrix_, IsA, IsA, mat_reuse, &IA);
281  MatGetOwnershipRange(matrix_,&pos_start,PETSC_NULL);
282  MatGetOwnershipRange(IA,&pos_start_IA,PETSC_NULL);
283 
284  std::vector<PetscInt> submat_rows;
285  const PetscInt *cols;
286  const PetscScalar *vals;
287 
288  std::vector<unsigned int> processed_rows(loc_size_A,0);
289 
290  unsigned int mat_block=1; //actual processed block of matrix
291  for(unsigned int loc_row=0; loc_row < processed_rows.size(); loc_row++) {
292  if (processed_rows[loc_row] != 0) continue;
293 
294  PetscInt min=std::numeric_limits<int>::max(), max=-1, size_submat;
295  PetscInt b_vals = 0; // count of values stored in B-block of Orig system
296  submat_rows.clear();
297  MatGetRow(matrix_, loc_row + pos_start, &ncols, &cols, PETSC_NULL);
298  for (PetscInt i=0; i<ncols; i++) {
299  if (cols[i] < pos_start || cols[i] >= pos_start+loc_size_A) {
300  b_vals++;
301  } else {
302  if (cols[i] < min) {
303  min=cols[i];
304  }
305  if (cols[i] > max) {
306  max=cols[i];
307  }
308  }
309  }
310  size_submat = max - min + 1;
311  OLD_ASSERT(ncols-b_vals == size_submat, "Submatrix cannot contains empty values.\n");
312 
313  MatRestoreRow(matrix_, loc_row + pos_start, &ncols, &cols, PETSC_NULL);
314  arma::mat submat2(size_submat, size_submat);
315  submat2.zeros();
316  for (PetscInt i=0; i<size_submat; i++) {
317  processed_rows[ loc_row + i ] = mat_block;
318  submat_rows.push_back( i + loc_row + pos_start_IA );
319  MatGetRow(matrix_, i + loc_row + pos_start, &ncols, &cols, &vals);
320  for (PetscInt j=0; j<ncols; j++) {
321  if (cols[j] >= pos_start && cols[j] < pos_start+loc_size_A) {
322  submat2( i, cols[j] - loc_row - pos_start ) = vals[j];
323  }
324  }
325  MatRestoreRow(matrix_, i + loc_row + pos_start, &ncols, &cols, &vals);
326  }
327  // get inversion matrix
328  arma::mat invmat = submat2.i();
329  // stored to inversion IA matrix
330  const PetscInt* rows = &submat_rows[0];
331  MatSetValues(IA, submat_rows.size(), rows, submat_rows.size(), rows, invmat.memptr(), INSERT_VALUES);
332 
333  mat_block++;
334  }
335 
336  MatAssemblyBegin(IA, MAT_FINAL_ASSEMBLY);
337  MatAssemblyEnd(IA, MAT_FINAL_ASSEMBLY);
338 }
339 
340 
342 {
343  if (Compl != NULL) {
344  return Compl->get_solution_precision();
345  }
346  return std::numeric_limits<double>::infinity();
347 }
348 
349 
351  START_TIMER("SchurComplement::solve");
352  this->form_schur();
353 
354  //output schur complement in matlab file
355 // string output_file = FilePath("schur.m", FilePath::output_file);
356 // PetscViewer viewer;
357 // PetscViewerASCIIOpen(PETSC_COMM_WORLD, output_file.c_str(), &viewer);
358 // PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
359 // MatView( *const_cast<Mat*>(Compl->get_matrix()), viewer);
360 // VecView( *const_cast<Vec*>(Compl->get_rhs()), viewer);
361 
362  LinSys::SolveInfo si = Compl->solve();
363 // VecView(Compl->get_solution(), viewer);
364 
365  // TODO: Resolve step is not necessary inside of nonlinear solver. Can optimize.
366  this->resolve();
367  return si;
368 }
369 
370 
371 /**
372  * COMPUTE ELIMINATED PART OF THE ORIG. SYS. & RESTORE RHS and SOLUTION VECTORS
373  * x_1 = IA * RHS_1 - IAB * x_2
374  */
375 
377 {
378  this->form_schur();
379 
380  START_TIMER("SchurComplemet::resolve without form schur");
381 
382  chkerr(MatMult(IAB,Compl->get_solution(),Sol1));
383  chkerr(VecScale(Sol1,-1));
384  chkerr(MatMultAdd(IA,RHS1,Sol1,Sol1));
385 }
386 
387 
389 {
390  //DebugOut() << print_var(LinSys_PETSC::compute_residual());
391  resolve();
393 }
394 
395 
396 
397 /**
398  * SCHUR COMPLEMENT destructor
399  */
401 
402  if ( B != NULL ) chkerr(MatDestroy(&B));
403  if ( Bt != NULL ) chkerr(MatDestroy(&Bt));
404  if ( C != NULL ) chkerr(MatDestroy(&C));
405  if ( xA != NULL ) chkerr(MatDestroy(&xA));
406  if ( IA != NULL ) chkerr(MatDestroy(&IA));
407  if ( IAB != NULL ) chkerr(MatDestroy(&IAB));
408  if ( IsA != NULL ) chkerr(ISDestroy(&IsA));
409  if ( IsB != NULL ) chkerr(ISDestroy(&IsB));
410  if ( RHS1 != NULL ) chkerr(VecDestroy(&RHS1));
411  if ( RHS2 != NULL ) chkerr(VecDestroy(&RHS2));
412  if ( Sol1 != NULL ) chkerr(VecDestroy(&Sol1));
413  if ( Sol2 != NULL ) chkerr(VecDestroy(&Sol2));
414  if ( IA != NULL ) chkerr(MatDestroy(&IA));
415 
416  if (Compl != NULL) delete Compl;
417 
418 }
bool is_negative_definite()
Definition: linsys.hh:580
void set_complement(LinSys_PETSC *ls)
Set complement LinSys object.
Definition: schur.cc:259
void resolve()
Definition: schur.cc:376
const Mat * get_matrix() override
Definition: linsys_PETSC.hh:69
void set_tolerances(double r_tol, double a_tol, unsigned int max_it) override
Definition: schur.cc:253
bool matrix_changed_
true if the matrix was changed since the last solve
Definition: linsys.hh:661
void set_from_input(const Input::Record in_rec) override
LinSys::SolveInfo solve() override
Wrappers for linear systems based on MPIAIJ and MATIS format.
Definition: schur.hh:59
void set_rhs_changed()
Definition: linsys.hh:218
Vec rhs_
PETSc vector constructed with vx array.
LinSys_PETSC * Compl
Definition: schur.hh:159
LinSys::SolveInfo solve() override
Definition: schur.cc:350
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
Definition: system.hh:147
int loc_size_B
Definition: schur.hh:151
bool rhs_changed_
true if the right hand side was changed since the last solve
Definition: linsys.hh:662
Distribution * ds_
Definition: schur.hh:161
double get_solution_precision() override
~SchurComplement()
Definition: schur.cc:400
void form_rhs()
Definition: schur.cc:235
Assembly explicit Schur complement for the given linear system. Provides method for resolution of the...
SchurState state
Definition: schur.hh:156
#define OLD_ASSERT(...)
Definition: global_defs.h:131
const Vec & get_solution()
Definition: linsys.hh:282
unsigned int begin(int proc) const
get starting local index
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
#define START_TIMER(tag)
Starts a timer with specified tag.
const Distribution * rows_ds_
final distribution of rows of MH matrix
Definition: linsys.hh:654
void set_from_input(const Input::Record in_rec) override
Definition: schur.cc:112
double compute_residual() override
Definition: schur.cc:388
void form_schur()
Definition: schur.cc:140
int loc_size_A
Definition: schur.hh:151
double compute_residual() override
SchurComplement(IS ia, Distribution *ds)
Definition: schur.cc:61
void set_solution(Vec sol_vec)
Definition: linsys.hh:290
void create_inversion_matrix()
create IA matrix
Definition: schur.cc:272
Support classes for parallel programing.
double get_solution_precision() override
get precision of solving
Definition: schur.cc:341
const Vec * get_rhs() override
Definition: linsys_PETSC.hh:74
Distribution * make_complement_distribution()
get distribution of complement object if complement is defined
Definition: schur.cc:266
Definition: schur.hh:60
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
Definition: asserts.hh:335
void set_matrix_changed()
Definition: linsys.hh:212
void set_tolerances(double r_tol, double a_tol, unsigned int max_it) override
Definition: linsys_PETSC.cc:90
Mat matrix_
Petsc matrix of the problem.
#define print_var(var)
Definition: logger.hh:260
#define DebugOut()
Macro defining &#39;debug&#39; record of log.
Definition: logger.hh:252
Vec solution_
PETSc vector constructed with vb array.
Definition: linsys.hh:664
Solver based on Multilevel BDDC - using corresponding class of OpenFTL package.
unsigned int lsize(int proc) const
get local size