Flow123d  release_2.2.0-914-gf1a3a4f
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  xA = NULL;
73  IAB = NULL;
74  IsB = NULL;
75  RHS1 = NULL;
76  RHS2 = NULL;
77  Sol1 = NULL;
78  Sol2 = NULL;
79 
80  // create A block index set
81  ISGetLocalSize(IsA, &loc_size_A);
82 
83  // create B block index set
85  ISCreateStride(PETSC_COMM_WORLD,loc_size_B,rows_ds_->begin()+loc_size_A,1,&IsB);
86 }
87 
88 
90 : LinSys_PETSC(other),
91  loc_size_A(other.loc_size_A), loc_size_B(other.loc_size_B), state(other.state),
92  Compl(other.Compl), ds_(other.ds_)
93 {
94  MatCopy(other.IA, IA, DIFFERENT_NONZERO_PATTERN);
95  MatCopy(other.IAB, IAB, DIFFERENT_NONZERO_PATTERN);
96  ISCopy(other.IsA, IsA);
97  ISCopy(other.IsB, IsB);
98  VecCopy(other.RHS1, RHS1);
99  VecCopy(other.RHS2, RHS2);
100  VecCopy(other.Sol1, Sol1);
101  VecCopy(other.Sol2, Sol2);
102 
103  B = NULL;
104  Bt = NULL;
105  xA = NULL;
106 }
107 
108 
109 
111 {
113 
114  ASSERT_PTR(Compl).error();
115  Compl->set_from_input( in_rec );
116 }
117 
118 /**
119  * COMPUTE A SCHUR COMPLEMENT OF A PETSC MATRIX
120  *
121  * given symmetric original matrix Orig has form
122  * A B x_1 RHS_1
123  * B' C * x_2 = RHS_2
124  * where the first block is given by index set IsA, and the second block by IsB
125  * user has to provide inverse IA of the A-block
126  * we suppose that original matrix have non-zero pattern for the schur complement
127  *
128  * we return: Shur - schur complement, ShurRHS - RHS of the complemented system:
129  * (B' * IA * B - C) * x_2 = (B' * IA * RHS_1 - RHS_2)
130  * IAB - a matrix to compute eliminated part of the solution:
131  * x_1 = IA * RHS_1 - IAB * x_2
132  *
133  * Actually as B' is stored separetly, the routine can be used also for nonsymetric original
134  * system
135  *
136  */
137 
139 {
140  START_TIMER("form schur complement");
141 
142 
143  PetscErrorCode ierr = 0;
144  MatReuse mat_reuse; // reuse structures after first computation of schur
145  MatStructure mat_subset_pattern;
146  PetscScalar *rhs_array, *sol_array;
147 
148  mat_reuse=MAT_REUSE_MATRIX;
149  mat_subset_pattern=SUBSET_NONZERO_PATTERN;
150  if (state==created) {
151  mat_reuse=MAT_INITIAL_MATRIX; // indicate first construction
152  mat_subset_pattern=DIFFERENT_NONZERO_PATTERN;
153 
154  // create complement system
155  // TODO: introduce LS as true object, clarify its internal states
156  // create RHS sub vecs RHS1, RHS2
157  VecGetArray(rhs_, &rhs_array);
158  VecCreateMPIWithArray(PETSC_COMM_WORLD,1,loc_size_A,PETSC_DETERMINE,rhs_array,&(RHS1));
159 
160  // create Solution sub vecs Sol1, Compl->solution
161  VecGetArray(solution_, &sol_array);
162  VecCreateMPIWithArray(PETSC_COMM_WORLD,1,loc_size_A,PETSC_DETERMINE,sol_array,&(Sol1));
163 
164  VecCreateMPIWithArray(PETSC_COMM_WORLD,1,loc_size_B,PETSC_DETERMINE,rhs_array+loc_size_A,&(RHS2));
165  VecCreateMPIWithArray(PETSC_COMM_WORLD,1,loc_size_B,PETSC_DETERMINE,sol_array+loc_size_A,&(Sol2));
166 
167  VecRestoreArray(rhs_, &rhs_array);
168  VecRestoreArray(solution_, &sol_array);
169 
170  VecGetArray( Sol2, &sol_array );
171  Compl->set_solution( sol_array );
172 
173  VecRestoreArray( Sol2, &sol_array );
174 
175 
176  }
177  DebugOut() << print_var(mat_reuse) << print_var(matrix_changed_) << print_var(state);
178 
179  // compose Schur complement
180  // Petsc need some fill estimate for results of multiplication in form nnz(A*B)/(nnz(A)+nnz(B))
181  // for the first Schur compl: IA*B is bounded by ( d*(d+1) )/( d*d+2*d ) <= 5/6 for d<=4
182  // B'*IA*B bounded by ( (d+1)*(d+1) )/ ( d*(d+1) + d ) ~ 1
183  // for the second Schur : IA*B have fill ratio ~ 1.
184  // B'*IA*B ... ( N/2 *(2*N-1) )/( 2 + 2*N ) <= 1.4
185  // nevertheless Petsc does not allows fill ratio below 1. so we use 1.1 for the first
186  // and 1.5 for the second multiplication
187 
188  // TODO:
189  // In order to let PETSC allocate structure of the complement we can not perform MatGetSubMatrix
190  // and MatAXPY on the same complement matrix, since one operation change the matrix structure for the other.
191  //
192  // Probably no way to make this optimal using high level methods. We should have our own
193  // format for schur complement matrix, store local systems and perform elimination localy.
194  // Or even better assembly the complement directly. (not compatible with raw P0 method)
195 
196  if (matrix_changed_) {
198 
199  // compute IAB=IA*B, loc_size_B removed
200  ierr+=MatGetSubMatrix(matrix_, IsA, IsB, mat_reuse, &B);
201  ierr+=MatMatMult(IA, B, mat_reuse, 1.0 ,&(IAB)); // 6/7 - fill estimate
202  // compute xA=Bt* IAB = Bt * IA * B, locSizeA removed
203  ierr+=MatGetSubMatrix(matrix_, IsB, IsA, mat_reuse, &(Bt));
204  ierr+=MatMatMult(Bt, IAB, mat_reuse, 1.9 ,&(xA)); // 1.1 - fill estimate (PETSC report values over 1.8)
205 
206  // get C block, loc_size_B removed
207  ierr+=MatGetSubMatrix( matrix_, IsB, IsB, mat_reuse, &C);
208 
209  if (state==created) MatDuplicate(C, MAT_DO_NOT_COPY_VALUES, const_cast<Mat *>( Compl->get_matrix() ) );
210  MatZeroEntries( *( Compl->get_matrix()) );
211 
212  // compute complement = (-1)cA+xA = Bt*IA*B - C
213  if ( is_negative_definite() ) {
214  ierr+=MatAXPY(*( Compl->get_matrix() ), 1, C, SUBSET_NONZERO_PATTERN);
215  ierr+=MatAXPY(*( Compl->get_matrix() ), -1, xA, mat_subset_pattern);
216  } else {
217  ierr+=MatAXPY(*( Compl->get_matrix() ), -1, C, SUBSET_NONZERO_PATTERN);
218  ierr+=MatAXPY(*( Compl->get_matrix() ), 1, xA, mat_subset_pattern);
219  }
221 
222  OLD_ASSERT( ierr == 0, "PETSC Error during calculation of Schur complement.\n");
223 
224  }
225 
226  form_rhs();
227 
228  matrix_changed_ = false;
229 
230  state=formed;
231 }
232 
234 {
235  START_TIMER("form rhs");
236  if (rhs_changed_ || matrix_changed_) {
237  MatMultTranspose(IAB, RHS1, *( Compl->get_rhs() ));
238  VecAXPY(*( Compl->get_rhs() ), -1, RHS2);
239  if ( is_negative_definite() ) {
240  VecScale(*( Compl->get_rhs() ), -1.0);
241  }
243  rhs_changed_ = false;
244  }
245 
246  state=formed;
247 }
248 
249 
250 
251 void SchurComplement::set_tolerances(double r_tol, double a_tol, unsigned int max_it)
252 {
253  LinSys_PETSC::set_tolerances(r_tol, a_tol, max_it);
254  if (Compl !=nullptr) Compl->set_tolerances(r_tol, a_tol, max_it);
255 }
256 
258 {
259  ASSERT_PTR(ls).error();
260  Compl = ls;
261 }
262 
263 
265 {
266  ds_ = new Distribution(loc_size_B, PETSC_COMM_WORLD);
267  return ds_;
268 }
269 
271 {
272  START_TIMER("create inversion matrix");
273  PetscInt ncols, pos_start, pos_start_IA;
274 
275  MatReuse mat_reuse=MAT_REUSE_MATRIX;
276  if (state==created) mat_reuse=MAT_INITIAL_MATRIX; // indicate first construction
277 
278  MatGetSubMatrix(matrix_, IsA, IsA, mat_reuse, &IA);
279  MatGetOwnershipRange(matrix_,&pos_start,PETSC_NULL);
280  MatGetOwnershipRange(IA,&pos_start_IA,PETSC_NULL);
281 
282  std::vector<PetscInt> submat_rows;
283  const PetscInt *cols;
284  const PetscScalar *vals;
285 
286  std::vector<unsigned int> processed_rows(loc_size_A,0);
287 
288  unsigned int mat_block=1; //actual processed block of matrix
289  for(unsigned int loc_row=0; loc_row < processed_rows.size(); loc_row++) {
290  if (processed_rows[loc_row] != 0) continue;
291 
292  PetscInt min=std::numeric_limits<int>::max(), max=-1, size_submat;
293  PetscInt b_vals = 0; // count of values stored in B-block of Orig system
294  submat_rows.clear();
295  MatGetRow(matrix_, loc_row + pos_start, &ncols, &cols, PETSC_NULL);
296  for (PetscInt i=0; i<ncols; i++) {
297  if (cols[i] < pos_start || cols[i] >= pos_start+loc_size_A) {
298  b_vals++;
299  } else {
300  if (cols[i] < min) {
301  min=cols[i];
302  }
303  if (cols[i] > max) {
304  max=cols[i];
305  }
306  }
307  }
308  size_submat = max - min + 1;
309  OLD_ASSERT(ncols-b_vals == size_submat, "Submatrix cannot contains empty values.\n");
310 
311  MatRestoreRow(matrix_, loc_row + pos_start, &ncols, &cols, PETSC_NULL);
312  arma::mat submat2(size_submat, size_submat);
313  submat2.zeros();
314  for (PetscInt i=0; i<size_submat; i++) {
315  processed_rows[ loc_row + i ] = mat_block;
316  submat_rows.push_back( i + loc_row + pos_start_IA );
317  MatGetRow(matrix_, i + loc_row + pos_start, &ncols, &cols, &vals);
318  for (PetscInt j=0; j<ncols; j++) {
319  if (cols[j] >= pos_start && cols[j] < pos_start+loc_size_A) {
320  submat2( i, cols[j] - loc_row - pos_start ) = vals[j];
321  }
322  }
323  MatRestoreRow(matrix_, i + loc_row + pos_start, &ncols, &cols, &vals);
324  }
325  // get inversion matrix
326  arma::mat invmat = submat2.i();
327  // stored to inversion IA matrix
328  const PetscInt* rows = &submat_rows[0];
329  MatSetValues(IA, submat_rows.size(), rows, submat_rows.size(), rows, invmat.memptr(), INSERT_VALUES);
330 
331  mat_block++;
332  }
333 
334  MatAssemblyBegin(IA, MAT_FINAL_ASSEMBLY);
335  MatAssemblyEnd(IA, MAT_FINAL_ASSEMBLY);
336 }
337 
338 
340 {
341  if (Compl != NULL) {
342  return Compl->get_solution_precision();
343  }
344  return std::numeric_limits<double>::infinity();
345 }
346 
347 
349  START_TIMER("SchurComplement::solve");
350  this->form_schur();
351 
352  //output schur complement in matlab file
353 // string output_file = FilePath("schur.m", FilePath::output_file);
354 // PetscViewer viewer;
355 // PetscViewerASCIIOpen(PETSC_COMM_WORLD, output_file.c_str(), &viewer);
356 // PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
357 // MatView( *const_cast<Mat*>(Compl->get_matrix()), viewer);
358 // VecView( *const_cast<Vec*>(Compl->get_rhs()), viewer);
359 
360  LinSys::SolveInfo si = Compl->solve();
361 // VecView(Compl->get_solution(), viewer);
362 
363  // TODO: Resolve step is not necessary inside of nonlinear solver. Can optimize.
364  this->resolve();
365  return si;
366 }
367 
368 
369 /**
370  * COMPUTE ELIMINATED PART OF THE ORIG. SYS. & RESTORE RHS and SOLUTION VECTORS
371  * x_1 = IA * RHS_1 - IAB * x_2
372  */
373 
375 {
376  this->form_schur();
377 
378  START_TIMER("SchurComplemet::resolve without form schur");
379 
380  chkerr(MatMult(IAB,Compl->get_solution(),Sol1));
381  chkerr(VecScale(Sol1,-1));
382  chkerr(MatMultAdd(IA,RHS1,Sol1,Sol1));
383 }
384 
385 
387 {
388  //DebugOut() << print_var(LinSys_PETSC::compute_residual());
389  resolve();
391 }
392 
393 
394 
395 /**
396  * SCHUR COMPLEMENT destructor
397  */
399 
400  if ( B != NULL ) MatDestroy(&B);
401  if ( Bt != NULL ) MatDestroy(&Bt);
402  if ( C != NULL ) MatDestroy(&C);
403  if ( xA != NULL ) MatDestroy(&xA);
404  if ( IA != NULL ) MatDestroy(&IA);
405  if ( IAB != NULL ) MatDestroy(&IAB);
406  if ( IsA != NULL ) ISDestroy(&IsA);
407  if ( IsB != NULL ) ISDestroy(&IsB);
408  if ( RHS1 != NULL ) VecDestroy(&RHS1);
409  if ( RHS2 != NULL ) VecDestroy(&RHS2);
410  if ( Sol1 != NULL ) VecDestroy(&Sol1);
411  if ( Sol2 != NULL ) VecDestroy(&Sol2);
412  if ( IA != NULL ) MatDestroy(&IA);
413 
414  if (Compl != NULL) delete Compl;
415 
416 }
bool is_negative_definite()
Definition: linsys.hh:560
void set_complement(LinSys_PETSC *ls)
Set complement LinSys object.
Definition: schur.cc:257
void resolve()
Definition: schur.cc:374
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:251
bool matrix_changed_
true if the matrix was changed since the last solve
Definition: linsys.hh:642
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:217
Vec rhs_
PETSc vector constructed with vx array.
LinSys_PETSC * Compl
Definition: schur.hh:159
LinSys::SolveInfo solve() override
Definition: schur.cc:348
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:643
Distribution * ds_
Definition: schur.hh:161
double get_solution_precision() override
~SchurComplement()
Definition: schur.cc:398
void form_rhs()
Definition: schur.cc:233
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:281
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:635
void set_from_input(const Input::Record in_rec) override
Definition: schur.cc:110
double compute_residual() override
Definition: schur.cc:386
void form_schur()
Definition: schur.cc:138
int loc_size_A
Definition: schur.hh:151
void set_solution(double *sol_array)
Definition: linsys.hh:289
double compute_residual() override
SchurComplement(IS ia, Distribution *ds)
Definition: schur.cc:61
void create_inversion_matrix()
create IA matrix
Definition: schur.cc:270
Support classes for parallel programing.
double get_solution_precision() override
get precision of solving
Definition: schur.cc:339
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:264
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:211
void set_tolerances(double r_tol, double a_tol, unsigned int max_it) override
Definition: linsys_PETSC.cc:86
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:645
Solver based on Multilevel BDDC - using corresponding class of OpenFTL package.
unsigned int lsize(int proc) const
get local size