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