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