Flow123d  jenkins-Flow123d-windows32-release-multijob-28
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  //ISView(IsA, PETSC_VIEWER_STDOUT_WORLD);
93 
94  // create B block index set
96  ISCreateStride(PETSC_COMM_WORLD,loc_size_B,rows_ds_->begin()+loc_size_A,1,&IsB);
97  //ISView(IsB, PETSC_VIEWER_STDOUT_WORLD);
98 }
99 
100 
102 : LinSys_PETSC(other),
103  loc_size_A(other.loc_size_A), loc_size_B(other.loc_size_B), state(other.state),
104  Compl(other.Compl), ds_(other.ds_)
105 {
106  MatCopy(other.IA, IA, DIFFERENT_NONZERO_PATTERN);
107  MatCopy(other.IAB, IAB, DIFFERENT_NONZERO_PATTERN);
108  ISCopy(other.IsA, IsA);
109  ISCopy(other.IsB, IsB);
110  VecCopy(other.RHS1, RHS1);
111  VecCopy(other.RHS2, RHS2);
112  VecCopy(other.Sol1, Sol1);
113  VecCopy(other.Sol2, Sol2);
114 
115  B = NULL;
116  Bt = NULL;
117  xA = NULL;
118 }
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  PetscErrorCode ierr = 0;
144  MatReuse mat_reuse; // reuse structures after first computation of schur
145  PetscScalar *rhs_array, *sol_array;
146 
147  mat_reuse=MAT_REUSE_MATRIX;
148  if (state==created) {
149  mat_reuse=MAT_INITIAL_MATRIX; // indicate first construction
150 
151  // create complement system
152  // TODO: introduce LS as true object, clarify its internal states
153  // create RHS sub vecs RHS1, RHS2
154  VecGetArray(rhs_, &rhs_array);
155  VecCreateMPIWithArray(PETSC_COMM_WORLD,1,loc_size_A,PETSC_DETERMINE,rhs_array,&(RHS1));
156 
157  // create Solution sub vecs Sol1, Compl->solution
158  VecGetArray(solution_, &sol_array);
159  VecCreateMPIWithArray(PETSC_COMM_WORLD,1,loc_size_A,PETSC_DETERMINE,sol_array,&(Sol1));
160 
161  VecCreateMPIWithArray(PETSC_COMM_WORLD,1,loc_size_B,PETSC_DETERMINE,rhs_array+loc_size_A,&(RHS2));
162  VecCreateMPIWithArray(PETSC_COMM_WORLD,1,loc_size_B,PETSC_DETERMINE,sol_array+loc_size_A,&(Sol2));
163 
164  VecRestoreArray(rhs_, &rhs_array);
165  VecRestoreArray(solution_, &sol_array);
166 
167  VecGetArray( Sol2, &sol_array );
168  Compl->set_solution( sol_array );
170  VecRestoreArray( Sol2, &sol_array );
171 
172  }
173 
174  //DBGMSG("Compute Schur complement of\n");
175  //MatView(matrix_,PETSC_VIEWER_STDOUT_WORLD);
176  //DBGMSG("inverse IA:\n");
177  //MatView(Schur->IA,PETSC_VIEWER_STDOUT_WORLD);
178  // compose Schur complement
179  // Petsc need some fill estimate for results of multiplication in form nnz(A*B)/(nnz(A)+nnz(B))
180  // for the first Schur compl: IA*B is bounded by ( d*(d+1) )/( d*d+2*d ) <= 5/6 for d<=4
181  // B'*IA*B bounded by ( (d+1)*(d+1) )/ ( d*(d+1) + d ) ~ 1
182  // for the second Schur : IA*B have fill ratio ~ 1.
183  // B'*IA*B ... ( N/2 *(2*N-1) )/( 2 + 2*N ) <= 1.4
184  // nevertheless Petsc does not allows fill ratio below 1. so we use 1.1 for the first
185  // and 1.5 for the second multiplication
186 
187  if (matrix_changed_) {
189 
190  // compute IAB=IA*B, loc_size_B removed
191  ierr+=MatGetSubMatrix(matrix_, IsA, IsB, mat_reuse, &B);
192  //DBGMSG(" B:\n");
193  //MatView(B,PETSC_VIEWER_STDOUT_WORLD);
194  ierr+=MatMatMult(IA, B, mat_reuse, 1.0 ,&(IAB)); // 6/7 - fill estimate
195  //DBGMSG(" IAB:\n");
196  //MatView(IAB,PETSC_VIEWER_STDOUT_WORLD);
197  // compute xA=Bt* IAB = Bt * IA * B, locSizeA removed
198  ierr+=MatGetSubMatrix(matrix_, IsB, IsA, mat_reuse, &(Bt));
199  ierr+=MatMatMult(Bt, IAB, mat_reuse, 1.9 ,&(xA)); // 1.1 - fill estimate (PETSC report values over 1.8)
200  //DBGMSG("xA:\n");
201  //MatView(xA,PETSC_VIEWER_STDOUT_WORLD);
202 
203  // get C block, loc_size_B removed
204  ierr+=MatGetSubMatrix( matrix_, IsB, IsB, mat_reuse, const_cast<Mat *>( Compl->get_matrix() ) );
205  // compute complement = (-1)cA+xA = Bt*IA*B - C
206  if ( is_negative_definite() ) {
207  ierr+=MatAXPY(*( Compl->get_matrix() ), -1, xA, SUBSET_NONZERO_PATTERN);
208  } else {
209  ierr+=MatScale(*( Compl->get_matrix() ),-1.0);
210  ierr+=MatAXPY(*( Compl->get_matrix() ), 1, xA, SUBSET_NONZERO_PATTERN);
211  }
213  //DBGMSG("C block:\n");
214 
215  //MatView(Compl->get_matrix(),PETSC_VIEWER_STDOUT_WORLD);
216  //DBGMSG("C block:\n");
217  //MatView(Schur->Compl->A,PETSC_VIEWER_STDOUT_WORLD);
218  //
219  //TODO: MatAXPY - umoznuje nasobit -1, t.j. bylo by lepe vytvorit konvencni Schuruv doplnek zde,
220  // a ve funkci schur1 (a ne v schur2) uzit metodu "scale" z tohoto objektu - kvuli negativni definitnosti
221  // usetri se tim jeden MatScale
222 
223  ASSERT( ierr == 0, "PETSC Error during calculation of Schur complement.\n");
224 
225  }
226 
227  /*
228  PetscViewerASCIIOpen(PETSC_COMM_WORLD,"matAinv.output",&myViewer);
229  PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
230  MatView( IA, myViewer );
231  PetscViewerDestroy(myViewer);
232 
233  PetscViewerASCIIOpen(PETSC_COMM_WORLD,"matSchur.output",&myViewer);
234  PetscViewerSetFormat(myViewer,PETSC_VIEWER_ASCII_MATLAB);
235  MatView( Compl->get_matrix( ), myViewer );
236  PetscViewerDestroy(myViewer);
237 */
238 
239  form_rhs();
240 
241  matrix_changed_ = false;
242 
243  state=formed;
244 }
245 
247 {
248  if (rhs_changed_ || matrix_changed_) {
249  MatMultTranspose(IAB, RHS1, *( Compl->get_rhs() ));
250  VecAXPY(*( Compl->get_rhs() ), -1, RHS2);
251  if ( is_negative_definite() ) {
252  VecScale(*( Compl->get_rhs() ), -1.0);
253  }
255  rhs_changed_ = false;
256  }
257 
258  state=formed;
259 }
260 
261 
262 /**
263  * COMPUTE ELIMINATED PART OF THE ORIG. SYS. & RESTORE RHS and SOLUTION VECTORS
264  * x_1 = IA * RHS_1 - IAB * x_2
265  */
266 
268 {
269  MatMult(IAB,Compl->get_solution(),Sol1);
270 
271  VecScale(Sol1,-1);
272 
273  MatMultAdd(IA,RHS1,Sol1,Sol1);
274 
275 }
276 
278 {
279  ASSERT(ls != nullptr, "NULL complement ls.\n");
280  Compl = ls;
281 }
282 
283 
285 {
286  ds_ = new Distribution(loc_size_B, PETSC_COMM_WORLD);
287  return ds_;
288 }
289 
291 {
292  PetscInt ncols, pos_start, pos_start_IA;
293 
294  MatReuse mat_reuse=MAT_REUSE_MATRIX;
295  if (state==created) mat_reuse=MAT_INITIAL_MATRIX; // indicate first construction
296 
297  MatGetSubMatrix(matrix_, IsA, IsA, mat_reuse, &IA);
298  //MatView(IA,PETSC_VIEWER_STDOUT_WORLD);
299  MatGetOwnershipRange(matrix_,&pos_start,PETSC_NULL);
300  MatGetOwnershipRange(IA,&pos_start_IA,PETSC_NULL);
301 
302  std::vector<PetscInt> submat_rows;
303  const PetscInt *cols;
304  const PetscScalar *vals;
305 
306  std::vector<unsigned int> processed_rows(loc_size_A,0);
307 
308  unsigned int mat_block=1; //actual processed block of matrix
309  for(unsigned int loc_row=0; loc_row < processed_rows.size(); loc_row++) {
310  if (processed_rows[loc_row] != 0) continue;
311 
312  PetscInt min=std::numeric_limits<int>::max(), max=-1, size_submat;
313  PetscInt b_vals = 0; // count of values stored in B-block of Orig system
314  submat_rows.clear();
315  MatGetRow(matrix_, loc_row + pos_start, &ncols, &cols, PETSC_NULL);
316  for (PetscInt i=0; i<ncols; i++) {
317  if (cols[i] < pos_start || cols[i] >= pos_start+loc_size_A) {
318  b_vals++;
319  } else {
320  if (cols[i] < min) {
321  min=cols[i];
322  }
323  if (cols[i] > max) {
324  max=cols[i];
325  }
326  }
327  }
328  size_submat = max - min + 1;
329  ASSERT(ncols-b_vals == size_submat, "Submatrix cannot contains empty values.\n");
330 
331  MatRestoreRow(matrix_, loc_row + pos_start, &ncols, &cols, PETSC_NULL);
332  arma::mat submat2(size_submat, size_submat);
333  submat2.zeros();
334  for (PetscInt i=0; i<size_submat; i++) {
335  processed_rows[ loc_row + i ] = mat_block;
336  submat_rows.push_back( i + loc_row + pos_start_IA );
337  MatGetRow(matrix_, i + loc_row + pos_start, &ncols, &cols, &vals);
338  for (PetscInt j=0; j<ncols; j++) {
339  if (cols[j] >= pos_start && cols[j] < pos_start+loc_size_A) {
340  submat2( i, cols[j] - loc_row - pos_start ) = vals[j];
341  }
342  }
343  MatRestoreRow(matrix_, i + loc_row + pos_start, &ncols, &cols, &vals);
344  }
345  // test output
346 // xprintf(Msg, "__ Get submat: rank %d, MIN-MAX %d %d, size %d\n", rank, min, max, size_submat);
347 // for (int i=0; i<size_submat; i++) {
348 // for (int j=0; j<size_submat; j++) {
349 // xprintf(Msg, "%2.0f ", submat2(i,j));
350 // }
351 // xprintf(Msg, "\n");
352 // }
353 // xprintf(Msg, "\n");
354  // get inversion matrix
355  arma::mat invmat = submat2.i();
356  // stored to inversion IA matrix
357  const PetscInt* rows = &submat_rows[0];
358  MatSetValues(IA, submat_rows.size(), rows, submat_rows.size(), rows, invmat.memptr(), INSERT_VALUES);
359 
360  mat_block++;
361  }
362 
363  MatAssemblyBegin(IA, MAT_FINAL_ASSEMBLY);
364  MatAssemblyEnd(IA, MAT_FINAL_ASSEMBLY);
365 }
366 
367 
369 {
370  if (Compl != NULL) {
371  return Compl->get_solution_precision();
372  }
373  return std::numeric_limits<double>::infinity();
374 }
375 
376 
378  this->form_schur();
379 
380  int converged_reason = Compl->solve();
381  this->resolve();
382 
383  return converged_reason;
384 }
385 
386 
387 /**
388  * SCHUR COMPLEMENT destructor
389  */
391 
392  if ( B != NULL ) MatDestroy(&B);
393  if ( Bt != NULL ) MatDestroy(&Bt);
394  if ( xA != NULL ) MatDestroy(&xA);
395  if ( IA != NULL ) MatDestroy(&IA);
396  if ( IAB != NULL ) MatDestroy(&IAB);
397  if ( IsA != NULL ) ISDestroy(&IsA);
398  if ( IsB != NULL ) ISDestroy(&IsB);
399  if ( RHS1 != NULL ) VecDestroy(&RHS1);
400  if ( RHS2 != NULL ) VecDestroy(&RHS2);
401  if ( Sol1 != NULL ) VecDestroy(&Sol1);
402  if ( Sol2 != NULL ) VecDestroy(&Sol2);
403  if ( IA != NULL ) MatDestroy(&IA);
404 
405  if (Compl != NULL) delete Compl;
406 
407 }
const Vec * get_rhs()
Definition: linsys_PETSC.hh:68
bool is_negative_definite()
Definition: linsys.hh:528
void set_complement(LinSys_PETSC *ls)
Set complement LinSys object.
Definition: schur.cc:277
void resolve()
Definition: schur.cc:267
bool matrix_changed_
true if the matrix was changed since the last solve
Definition: linsys.hh:610
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:211
Vec rhs_
PETSc vector constructed with vx array.
LinSys_PETSC * Compl
Definition: schur.hh:146
void set_from_input(const Input::Record in_rec)
int loc_size_B
Definition: schur.hh:136
Input::Record in_rec_
Definition: linsys.hh:623
bool rhs_changed_
true if the right hand side was changed since the last solve
Definition: linsys.hh:611
Distribution * ds_
Definition: schur.hh:148
~SchurComplement()
Definition: schur.cc:390
void form_rhs()
Definition: schur.cc:246
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:120
const Vec & get_solution()
Definition: linsys.hh:266
unsigned int begin(int proc) const
get starting local index
const Distribution * rows_ds_
final distribution of rows of MH matrix
Definition: linsys.hh:603
void form_schur()
Definition: schur.cc:141
int loc_size_A
Definition: schur.hh:136
int solve() override
Definition: schur.cc:377
void set_solution(double *sol_array)
Definition: linsys.hh:274
SchurComplement(IS ia, Distribution *ds)
Definition: schur.cc:71
void create_inversion_matrix()
create IA matrix
Definition: schur.cc:290
Support classes for parallel programing.
double get_solution_precision() override
get precision of solving
Definition: schur.cc:368
Distribution * make_complement_distribution()
get distribution of complement object if complement is defined
Definition: schur.cc:284
Definition: schur.hh:69
double get_solution_precision()
void set_matrix_changed()
Definition: linsys.hh:205
Mat matrix_
Petsc matrix of the problem.
Vec solution_
PETSc vector constructed with vb array.
Definition: linsys.hh:613
Solver based on Multilevel BDDC - using corresponding class of OpenFTL package.
unsigned int lsize(int proc) const
get local size