MAST
complex_solver_base.cpp
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2019 Manav Bhatia
4  *
5  * This library is free software; you can redistribute it and/or
6  * modify it under the terms of the GNU Lesser General Public
7  * License as published by the Free Software Foundation; either
8  * version 2.1 of the License, or (at your option) any later version.
9  *
10  * This library is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13  * Lesser General Public License for more details.
14  *
15  * You should have received a copy of the GNU Lesser General Public
16  * License along with this library; if not, write to the Free Software
17  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
18  */
19 
20 
21 // MAST includes
24 #include "base/nonlinear_system.h"
25 
26 
27 // libMesh includes
28 #include "libmesh/numeric_vector.h"
29 #include "libmesh/system.h"
30 #include "libmesh/petsc_matrix.h"
31 #include "libmesh/petsc_vector.h"
32 #include "libmesh/libmesh_common.h"
33 #include "libmesh/dof_map.h"
34 
35 // PETSc includes
36 #include <petscmat.h>
37 
38 
40 tol (1.0e-3),
41 max_iters (20),
42 _assembly (nullptr) {
43 
44 }
45 
46 
47 
49 
50 }
51 
52 
53 
54 void
56 
57  libmesh_assert(!_assembly);
58 
59  _assembly = &assemble;
60 }
61 
62 
63 
64 void
66 
68 
69  // remove the real part of the vector
70  std::string nm;
71  nm = sys.name();
72  nm += "real_sol";
73 
74  if (sys.have_vector(nm))
75  sys.remove_vector(nm);
76 
77  // and its sensitivity
78  nm += "_sens";
79 
80  if (sys.have_vector(nm))
81  sys.remove_vector(nm);
82 
83  // remove the imaginary part of the vector
84  nm = sys.name();
85  nm += "imag_sol";
86 
87  if (sys.have_vector(nm))
88  sys.remove_vector(nm);
89 
90  // and its sensitivity
91  nm += "_sens";
92 
93  if (sys.have_vector(nm))
94  sys.remove_vector(nm);
95 
96 
97  _assembly = nullptr;
98 }
99 
100 
101 
102 libMesh::NumericVector<Real>&
104 
105  libmesh_assert(_assembly);
106 
108 
109  std::string nm;
110  nm = sys.name();
111  nm += "real_sol";
112  if (if_sens)
113  nm += "_sens";
114 
115  if (!sys.have_vector(nm))
116  sys.add_vector(nm);
117 
118  return sys.get_vector(nm);
119 }
120 
121 
122 
123 const libMesh::NumericVector<Real>&
125 
126  libmesh_assert(_assembly);
127 
129 
130  std::string nm;
131  nm = sys.name();
132  nm += "real_sol";
133  if (if_sens)
134  nm += "_sens";
135 
136  if (!sys.have_vector(nm))
137  sys.add_vector(nm);
138 
139  return sys.get_vector(nm);
140 }
141 
142 
143 libMesh::NumericVector<Real>&
145 
146  libmesh_assert(_assembly);
147 
149 
150  std::string nm;
151  nm = sys.name();
152  nm += "imag_sol";
153  if (if_sens)
154  nm += "_sens";
155 
156  if (!sys.have_vector(nm))
157  sys.add_vector(nm);
158 
159  return sys.get_vector(nm);
160 }
161 
162 
163 const libMesh::NumericVector<Real>&
165 
166  libmesh_assert(_assembly);
167 
169 
170  std::string nm;
171  nm = sys.name();
172  nm += "imag_sol";
173  if (if_sens)
174  nm += "_sens";
175 
176  if (!sys.have_vector(nm))
177  sys.add_vector(nm);
178 
179  return sys.get_vector(nm);
180 }
181 
182 
183 
184 
185 
186 void
188 
189  libmesh_assert(_assembly);
190 
191  START_LOG("complex_solve()", "PetscFieldSplitSolver");
192 
193  // get reference to the system
194  MAST::NonlinearSystem& sys =
195  dynamic_cast<MAST::NonlinearSystem&>(_assembly->system());
196 
197  // create a petsc nested matrix of size 2x2
198  PetscErrorCode ierr;
199  Mat mat;
200  Vec res, sol, res_vec_R, res_vec_I, sol_vec_R, sol_vec_I;
201  std::vector<IS> is(2);
202  std::vector<Mat> sub_mats(4);
203  sub_mats[0] = dynamic_cast<libMesh::PetscMatrix<Real>*>(sys.matrix)->mat(); // real
204  sub_mats[1] = dynamic_cast<libMesh::PetscMatrix<Real>*>(sys.matrix_A)->mat(); // -imag
205  sub_mats[2] = dynamic_cast<libMesh::PetscMatrix<Real>*>(sys.matrix_B)->mat(); // imag
206  sub_mats[3] = dynamic_cast<libMesh::PetscMatrix<Real>*>(sys.matrix)->mat(); // real
207 
208  ierr = MatCreateNest(_assembly->system().comm().get(),
209  2, nullptr,
210  2, nullptr,
211  &sub_mats[0],
212  &mat);
213  CHKERRABORT(sys.comm().get(), ierr);
214 
215 
216  // get the IS belonging to each block
217  ierr = MatNestGetISs(mat, &is[0], nullptr); CHKERRABORT(sys.comm().get(), ierr);
218 
219 
220 
221  // setup vector for residual
222  ierr = VecCreate(sys.comm().get(), &res); CHKERRABORT(sys.comm().get(), ierr);
223  ierr = VecSetSizes(res, PETSC_DECIDE, 2*sys.solution->size()); CHKERRABORT(sys.comm().get(), ierr);
224  ierr = VecSetType(res, VECMPI); CHKERRABORT(sys.comm().get(), ierr);
225 
226 
227  // setup vector for solution
228  ierr = VecDuplicate(res, &sol); CHKERRABORT(sys.comm().get(), ierr);
229 
230 
231 
232  // get the residual subvectors
233  ierr = VecGetSubVector(res, is[0], &res_vec_R); CHKERRABORT(sys.comm().get(), ierr);
234  ierr = VecGetSubVector(res, is[1], &res_vec_I); CHKERRABORT(sys.comm().get(), ierr);
235 
236  // get the solution subvectors
237  ierr = VecGetSubVector(sol, is[0], &sol_vec_R); CHKERRABORT(sys.comm().get(), ierr);
238  ierr = VecGetSubVector(sol, is[1], &sol_vec_I); CHKERRABORT(sys.comm().get(), ierr);
239 
240 
241  // clone vectors for use as system RHS
242  std::unique_ptr<libMesh::NumericVector<Real> >
243  res_R(new libMesh::PetscVector<Real>(res_vec_R, sys.comm())),
244  res_I(new libMesh::PetscVector<Real>(res_vec_I, sys.comm())),
245  sol_R(new libMesh::PetscVector<Real>(sol_vec_R, sys.comm())),
246  sol_I(new libMesh::PetscVector<Real>(sol_vec_I, sys.comm()));
247 
248  sol_R->zero();
249  sol_R->close();
250  sol_I->zero();
251  sol_I->close();
252 
253  // assemble the matrices
255  *sol_I,
256  *res_R,
257  *res_I,
258  *sys.matrix, // J_R
259  *sys.matrix_B); // J_I
260 
261  // restore the subvectors
262  //ierr = VecRestoreSubVector(res, is[0], &res_vec_R); CHKERRABORT(sys.comm().get(), ierr);
263  //ierr = VecRestoreSubVector(res, is[1], &res_vec_I); CHKERRABORT(sys.comm().get(), ierr);
264 
265  // copy the imag Jacobian for the first row of the sys of eq.
266  sys.matrix_A->zero();
267  sys.matrix_A->close();
268  sys.matrix_A->add(-1., *sys.matrix_B);
269  sys.matrix_A->close();
270 
271  ierr = MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY); CHKERRABORT(sys.comm().get(), ierr);
272  ierr = MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY); CHKERRABORT(sys.comm().get(), ierr);
273  ierr = VecAssemblyBegin(res); CHKERRABORT(sys.comm().get(), ierr);
274  ierr = VecAssemblyEnd(res); CHKERRABORT(sys.comm().get(), ierr);
275 
276 
277 
278  // now initialize the KSP and ask for solution.
279  KSP ksp;
280  PC pc;
281 
282  // setup the KSP
283  ierr = KSPCreate(sys.comm().get(), &ksp); CHKERRABORT(sys.comm().get(), ierr);
284  ierr = KSPSetOperators(ksp, mat, mat); CHKERRABORT(sys.comm().get(), ierr);
285  ierr = KSPSetFromOptions(ksp); CHKERRABORT(sys.comm().get(), ierr);
286 
287  // setup the PC
288  ierr = KSPGetPC(ksp, &pc); CHKERRABORT(sys.comm().get(), ierr);
289  ierr = PCFieldSplitSetIS(pc, nullptr, is[0]);CHKERRABORT(sys.comm().get(), ierr);
290  ierr = PCFieldSplitSetIS(pc, nullptr, is[1]);CHKERRABORT(sys.comm().get(), ierr);
291  ierr = PCSetFromOptions(pc); CHKERRABORT(sys.comm().get(), ierr);
292 
293 
294 
295  // now solve
296  ierr = KSPSolve(ksp, res, sol);
297 
298 
299  // assemble the matrices
301  *sol_I,
302  *res_R,
303  *res_I,
304  *sys.matrix, // J_R
305  *sys.matrix_B); // J_I
306 
307  // copy solutions for output
308  /*this->real_solution() = *sol_R;
309  this->imag_solution() = *sol_I;
310  this->real_solution().close();
311  this->imag_solution().close();*/
312 
313 
314  // restore the subvectors
315  ierr = VecRestoreSubVector(res, is[0], &res_vec_R); CHKERRABORT(sys.comm().get(), ierr);
316  ierr = VecRestoreSubVector(res, is[1], &res_vec_I); CHKERRABORT(sys.comm().get(), ierr);
317  ierr = VecRestoreSubVector(sol, is[0], &sol_vec_R); CHKERRABORT(sys.comm().get(), ierr);
318  ierr = VecRestoreSubVector(sol, is[1], &sol_vec_I); CHKERRABORT(sys.comm().get(), ierr);
319 
320 
321  // destroy the objects
322  ierr = KSPDestroy(&ksp);
323  ierr = MatDestroy(&mat);
324  ierr = VecDestroy(&res);
325  ierr = VecDestroy(&sol);
326 
327  STOP_LOG("complex_solve()", "PetscFieldSplitSolver");
328 }
329 
330 
331 
332 void
334 
335  libmesh_assert(_assembly);
336 
337  START_LOG("solve_block_matrix()", "ComplexSolve");
338 
339  // get reference to the system
340  MAST::NonlinearSystem& sys =
341  dynamic_cast<MAST::NonlinearSystem&>(_assembly->system());
342 
343  libMesh::DofMap& dof_map = sys.get_dof_map();
344 
345  const PetscInt
346  my_m = dof_map.n_dofs(),
347  my_n = my_m,
348  n_l = dof_map.n_dofs_on_processor(sys.processor_id()),
349  m_l = n_l;
350 
351  const std::vector<libMesh::dof_id_type>
352  & n_nz = dof_map.get_n_nz(),
353  & n_oz = dof_map.get_n_oz();
354 
355  std::vector<libMesh::dof_id_type>
356  complex_n_nz (2*n_nz.size()),
357  complex_n_oz (2*n_oz.size());
358 
359  // create the n_nz and n_oz for the complex matrix without block format
360  for (unsigned int i=0; i<n_nz.size(); i++) {
361 
362  complex_n_nz[2*i] = 2*n_nz[i];
363  complex_n_nz[2*i+1] = 2*n_nz[i];
364  }
365 
366  for (unsigned int i=0; i<n_oz.size(); i++) {
367 
368  complex_n_oz[2*i] = 2*n_oz[i];
369  complex_n_oz[2*i+1] = 2*n_oz[i];
370  }
371 
372 
373 
374  // create the matrix
375  PetscErrorCode ierr;
376  Mat mat;
377 
378  ierr = MatCreate(sys.comm().get(), &mat); CHKERRABORT(sys.comm().get(), ierr);
379  ierr = MatSetSizes(mat, 2*m_l, 2*n_l, 2*my_m, 2*my_n); CHKERRABORT(sys.comm().get(), ierr);
380 
381  if (libMesh::on_command_line("--solver_system_names")) {
382 
383  std::string nm = _assembly->system().name() + "_complex_";
384  MatSetOptionsPrefix(mat, nm.c_str());
385  }
386  ierr = MatSetFromOptions(mat); CHKERRABORT(sys.comm().get(), ierr);
387 
388  //ierr = MatSetType(mat, MATBAIJ); CHKERRABORT(sys.comm().get(), ierr);
389  ierr = MatSetBlockSize(mat, 2); CHKERRABORT(sys.comm().get(), ierr);
390  ierr = MatSeqAIJSetPreallocation(mat,
391  2*my_m,
392  (PetscInt*)&complex_n_nz[0]); CHKERRABORT(sys.comm().get(), ierr);
393  ierr = MatMPIAIJSetPreallocation(mat,
394  0,
395  (PetscInt*)&complex_n_nz[0],
396  0,
397  (PetscInt*)&complex_n_oz[0]); CHKERRABORT(sys.comm().get(), ierr);
398  ierr = MatSeqBAIJSetPreallocation (mat, 2,
399  0, (PetscInt*)&n_nz[0]); CHKERRABORT(sys.comm().get(), ierr);
400  ierr = MatMPIBAIJSetPreallocation (mat, 2,
401  0, (PetscInt*)&n_nz[0],
402  0, (PetscInt*)&n_oz[0]); CHKERRABORT(sys.comm().get(), ierr);
403  ierr = MatSetOption(mat,
404  MAT_NEW_NONZERO_ALLOCATION_ERR,
405  PETSC_TRUE); CHKERRABORT(sys.comm().get(), ierr);
406 
407 
408  // now create the vectors
409  Vec res_vec, sol_vec;
410 
411  ierr = MatCreateVecs(mat, &res_vec, PETSC_NULL); CHKERRABORT(sys.comm().get(), ierr);
412  ierr = MatCreateVecs(mat, &sol_vec, PETSC_NULL); CHKERRABORT(sys.comm().get(), ierr);
413 
414 
415  std::unique_ptr<libMesh::SparseMatrix<Real> >
416  jac_mat(new libMesh::PetscMatrix<Real>(mat, sys.comm()));
417 
418  std::unique_ptr<libMesh::NumericVector<Real> >
419  res(new libMesh::PetscVector<Real>(res_vec, sys.comm())),
420  sol(new libMesh::PetscVector<Real>(sol_vec, sys.comm()));
421 
422 
423  // if sensitivity analysis is requested, then set the complex solution in
424  // the solution vector
425  if (p) {
426 
427  // copy the solution to separate real and imaginary vectors
428  libMesh::NumericVector<Real>
429  &sol_R = this->real_solution(),
430  &sol_I = this->imag_solution();
431 
432  unsigned int
433  first = sol_R.first_local_index(),
434  last = sol_I.last_local_index();
435 
436  for (unsigned int i=first; i<last; i++) {
437 
438  sol->set( 2*i, sol_R(i));
439  sol->set(2*i+1, sol_I(i));
440  }
441 
442  sol->close();
443  }
444 
445 
446  // assemble the matrix
448  *res,
449  *jac_mat,
450  p);
451  res->scale(-1.);
452 
453  // now initialize the KSP and ask for solution.
454  KSP ksp;
455  PC pc;
456 
457  // setup the KSP
458  ierr = KSPCreate(sys.comm().get(), &ksp); CHKERRABORT(sys.comm().get(), ierr);
459 
460  if (libMesh::on_command_line("--solver_system_names")) {
461 
462  std::string nm = _assembly->system().name() + "_complex_";
463  KSPSetOptionsPrefix(ksp, nm.c_str());
464  }
465 
466  ierr = KSPSetOperators(ksp, mat, mat); CHKERRABORT(sys.comm().get(), ierr);
467  ierr = KSPSetFromOptions(ksp); CHKERRABORT(sys.comm().get(), ierr);
468 
469  // setup the PC
470  ierr = KSPGetPC(ksp, &pc); CHKERRABORT(sys.comm().get(), ierr);
471  ierr = PCSetFromOptions(pc); CHKERRABORT(sys.comm().get(), ierr);
472 
473 
474  START_LOG("KSPSolve", "ComplexSolve");
475 
476  // now solve
477  ierr = KSPSolve(ksp, res_vec, sol_vec);
478 
479  STOP_LOG("KSPSolve", "ComplexSolve");
480 
481 
482  // copy the solution to separate real and imaginary vectors
483  libMesh::NumericVector<Real>
484  &sol_R = this->real_solution(p != nullptr),
485  &sol_I = this->imag_solution(p != nullptr);
486 
487  unsigned int
488  first = sol_R.first_local_index(),
489  last = sol_R.last_local_index();
490 
491  for (unsigned int i=first; i<last; i++) {
492  sol_R.set(i, (*sol)( 2*i));
493  sol_I.set(i, (*sol)(2*i+1));
494  }
495 
496  sol_R.close();
497  sol_I.close();
498  sol->close();
499 
500  ierr = KSPDestroy(&ksp); CHKERRABORT(sys.comm().get(), ierr);
501  ierr = MatDestroy(&mat); CHKERRABORT(sys.comm().get(), ierr);
502  ierr = VecDestroy(&res_vec); CHKERRABORT(sys.comm().get(), ierr);
503  ierr = VecDestroy(&sol_vec); CHKERRABORT(sys.comm().get(), ierr);
504 
505  STOP_LOG("solve_block_matrix()", "ComplexSolve");
506 }
507 
508 
509 
510 
void residual_and_jacobian_blocked(const libMesh::NumericVector< Real > &X, libMesh::NumericVector< Real > &R, libMesh::SparseMatrix< Real > &J, MAST::Parameter *p=nullptr)
Assembles the residual and Jacobian of the N complex system of equations split into 2N real system of...
This class implements a system for solution of nonlinear systems.
void residual_and_jacobian_field_split(const libMesh::NumericVector< Real > &X_R, const libMesh::NumericVector< Real > &X_I, libMesh::NumericVector< Real > &R_R, libMesh::NumericVector< Real > &R_I, libMesh::SparseMatrix< Real > &J_R, libMesh::SparseMatrix< Real > &J_I)
void set_assembly(MAST::ComplexAssemblyBase &assemble)
sets the assembly object for this solver
This is a scalar function whose value can be changed and one that can be used as a design variable in...
Definition: parameter.h:35
libMesh::NumericVector< Real > & imag_solution(bool if_sens=false)
void clear_assembly()
clears the assembly object from this solver
virtual void solve_pc_fieldsplit()
solves the complex system of equations using PCFieldSplit
libMesh::SparseMatrix< Real > * matrix_A
The system matrix for standard eigenvalue problems.
virtual void solve_block_matrix(MAST::Parameter *p=nullptr)
solves the complex system of equations using block matrices.
libMesh::SparseMatrix< Real > * matrix_B
A second system matrix for generalized eigenvalue problems.
virtual ~ComplexSolverBase()
destructor
ComplexSolverBase()
default constructor
libMesh::NumericVector< Real > & real_solution(bool if_sens=false)
MAST::ComplexAssemblyBase * _assembly
Associated ComplexAssembly object that provides the element level quantities.
const MAST::NonlinearSystem & system() const