28 #include "libmesh/libmesh_logging.h" 29 #include "libmesh/numeric_vector.h" 30 #include "libmesh/system.h" 31 #include "libmesh/petsc_matrix.h" 32 #include "libmesh/petsc_vector.h" 33 #include "libmesh/libmesh_common.h" 34 #include "libmesh/dof_map.h" 75 if (sys.have_vector(nm))
76 sys.remove_vector(nm);
81 if (sys.have_vector(nm))
82 sys.remove_vector(nm);
88 if (sys.have_vector(nm))
89 sys.remove_vector(nm);
94 if (sys.have_vector(nm))
95 sys.remove_vector(nm);
103 libMesh::NumericVector<Real>&
116 if (!sys.have_vector(nm))
119 return sys.get_vector(nm);
124 const libMesh::NumericVector<Real>&
137 if (!sys.have_vector(nm))
140 return sys.get_vector(nm);
144 libMesh::NumericVector<Real>&
157 if (!sys.have_vector(nm))
160 return sys.get_vector(nm);
164 const libMesh::NumericVector<Real>&
177 if (!sys.have_vector(nm))
180 return sys.get_vector(nm);
192 START_LOG(
"complex_solve()",
"PetscFieldSplitSolver");
201 Vec res, sol, res_vec_R, res_vec_I, sol_vec_R, sol_vec_I;
202 std::vector<IS> is(2);
203 std::vector<Mat> sub_mats(4);
204 sub_mats[0] =
dynamic_cast<libMesh::PetscMatrix<Real>*
>(sys.matrix)->mat();
205 sub_mats[1] =
dynamic_cast<libMesh::PetscMatrix<Real>*
>(sys.
matrix_A)->mat();
206 sub_mats[2] =
dynamic_cast<libMesh::PetscMatrix<Real>*
>(sys.
matrix_B)->mat();
207 sub_mats[3] =
dynamic_cast<libMesh::PetscMatrix<Real>*
>(sys.matrix)->mat();
214 CHKERRABORT(sys.comm().get(), ierr);
218 ierr = MatNestGetISs(mat, &is[0],
nullptr); CHKERRABORT(sys.comm().get(), ierr);
223 ierr = VecCreate(sys.comm().get(), &res); CHKERRABORT(sys.comm().get(), ierr);
224 ierr = VecSetSizes(res, PETSC_DECIDE, 2*sys.solution->size()); CHKERRABORT(sys.comm().get(), ierr);
225 ierr = VecSetType(res, VECMPI); CHKERRABORT(sys.comm().get(), ierr);
229 ierr = VecDuplicate(res, &sol); CHKERRABORT(sys.comm().get(), ierr);
234 ierr = VecGetSubVector(res, is[0], &res_vec_R); CHKERRABORT(sys.comm().get(), ierr);
235 ierr = VecGetSubVector(res, is[1], &res_vec_I); CHKERRABORT(sys.comm().get(), ierr);
238 ierr = VecGetSubVector(sol, is[0], &sol_vec_R); CHKERRABORT(sys.comm().get(), ierr);
239 ierr = VecGetSubVector(sol, is[1], &sol_vec_I); CHKERRABORT(sys.comm().get(), ierr);
243 std::unique_ptr<libMesh::NumericVector<Real> >
244 res_R(
new libMesh::PetscVector<Real>(res_vec_R, sys.comm())),
245 res_I(
new libMesh::PetscVector<Real>(res_vec_I, sys.comm())),
246 sol_R(
new libMesh::PetscVector<Real>(sol_vec_R, sys.comm())),
247 sol_I(
new libMesh::PetscVector<Real>(sol_vec_I, sys.comm()));
272 ierr = MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY); CHKERRABORT(sys.comm().get(), ierr);
273 ierr = MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY); CHKERRABORT(sys.comm().get(), ierr);
274 ierr = VecAssemblyBegin(res); CHKERRABORT(sys.comm().get(), ierr);
275 ierr = VecAssemblyEnd(res); CHKERRABORT(sys.comm().get(), ierr);
284 ierr = KSPCreate(sys.comm().get(), &ksp); CHKERRABORT(sys.comm().get(), ierr);
285 ierr = KSPSetOperators(ksp, mat, mat); CHKERRABORT(sys.comm().get(), ierr);
286 ierr = KSPSetFromOptions(ksp); CHKERRABORT(sys.comm().get(), ierr);
289 ierr = KSPGetPC(ksp, &pc); CHKERRABORT(sys.comm().get(), ierr);
290 ierr = PCFieldSplitSetIS(pc,
nullptr, is[0]);CHKERRABORT(sys.comm().get(), ierr);
291 ierr = PCFieldSplitSetIS(pc,
nullptr, is[1]);CHKERRABORT(sys.comm().get(), ierr);
292 ierr = PCSetFromOptions(pc); CHKERRABORT(sys.comm().get(), ierr);
297 ierr = KSPSolve(ksp, res, sol);
316 ierr = VecRestoreSubVector(res, is[0], &res_vec_R); CHKERRABORT(sys.comm().get(), ierr);
317 ierr = VecRestoreSubVector(res, is[1], &res_vec_I); CHKERRABORT(sys.comm().get(), ierr);
318 ierr = VecRestoreSubVector(sol, is[0], &sol_vec_R); CHKERRABORT(sys.comm().get(), ierr);
319 ierr = VecRestoreSubVector(sol, is[1], &sol_vec_I); CHKERRABORT(sys.comm().get(), ierr);
323 ierr = KSPDestroy(&ksp);
324 ierr = MatDestroy(&mat);
325 ierr = VecDestroy(&res);
326 ierr = VecDestroy(&sol);
328 STOP_LOG(
"complex_solve()",
"PetscFieldSplitSolver");
338 START_LOG(
"solve_block_matrix()",
"ComplexSolve");
344 libMesh::DofMap& dof_map = sys.get_dof_map();
347 my_m = dof_map.n_dofs(),
349 n_l = dof_map.n_dofs_on_processor(sys.processor_id()),
352 const std::vector<libMesh::dof_id_type>
353 & n_nz = dof_map.get_n_nz(),
354 & n_oz = dof_map.get_n_oz();
356 std::vector<libMesh::dof_id_type>
357 complex_n_nz (2*n_nz.size()),
358 complex_n_oz (2*n_oz.size());
361 for (
unsigned int i=0; i<n_nz.size(); i++) {
363 complex_n_nz[2*i] = 2*n_nz[i];
364 complex_n_nz[2*i+1] = 2*n_nz[i];
367 for (
unsigned int i=0; i<n_oz.size(); i++) {
369 complex_n_oz[2*i] = 2*n_oz[i];
370 complex_n_oz[2*i+1] = 2*n_oz[i];
379 ierr = MatCreate(sys.comm().get(), &mat); CHKERRABORT(sys.comm().get(), ierr);
380 ierr = MatSetSizes(mat, 2*m_l, 2*n_l, 2*my_m, 2*my_n); CHKERRABORT(sys.comm().get(), ierr);
382 if (libMesh::on_command_line(
"--solver_system_names")) {
385 MatSetOptionsPrefix(mat, nm.c_str());
387 ierr = MatSetFromOptions(mat); CHKERRABORT(sys.comm().get(), ierr);
390 ierr = MatSetBlockSize(mat, 2); CHKERRABORT(sys.comm().get(), ierr);
391 ierr = MatSeqAIJSetPreallocation(mat,
393 (PetscInt*)&complex_n_nz[0]); CHKERRABORT(sys.comm().get(), ierr);
394 ierr = MatMPIAIJSetPreallocation(mat,
396 (PetscInt*)&complex_n_nz[0],
398 (PetscInt*)&complex_n_oz[0]); CHKERRABORT(sys.comm().get(), ierr);
399 ierr = MatSeqBAIJSetPreallocation (mat, 2,
400 0, (PetscInt*)&n_nz[0]); CHKERRABORT(sys.comm().get(), ierr);
401 ierr = MatMPIBAIJSetPreallocation (mat, 2,
402 0, (PetscInt*)&n_nz[0],
403 0, (PetscInt*)&n_oz[0]); CHKERRABORT(sys.comm().get(), ierr);
404 ierr = MatSetOption(mat,
405 MAT_NEW_NONZERO_ALLOCATION_ERR,
406 PETSC_TRUE); CHKERRABORT(sys.comm().get(), ierr);
410 Vec res_vec, sol_vec;
412 ierr = MatCreateVecs(mat, &res_vec, PETSC_NULL); CHKERRABORT(sys.comm().get(), ierr);
413 ierr = MatCreateVecs(mat, &sol_vec, PETSC_NULL); CHKERRABORT(sys.comm().get(), ierr);
416 std::unique_ptr<libMesh::SparseMatrix<Real> >
417 jac_mat(
new libMesh::PetscMatrix<Real>(mat, sys.comm()));
419 std::unique_ptr<libMesh::NumericVector<Real> >
420 res(
new libMesh::PetscVector<Real>(res_vec, sys.comm())),
421 sol(
new libMesh::PetscVector<Real>(sol_vec, sys.comm()));
429 libMesh::NumericVector<Real>
434 first = sol_R.first_local_index(),
435 last = sol_I.last_local_index();
437 for (
unsigned int i=first; i<last; i++) {
439 sol->set( 2*i, sol_R(i));
440 sol->set(2*i+1, sol_I(i));
459 ierr = KSPCreate(sys.comm().get(), &ksp); CHKERRABORT(sys.comm().get(), ierr);
461 if (libMesh::on_command_line(
"--solver_system_names")) {
464 KSPSetOptionsPrefix(ksp, nm.c_str());
467 ierr = KSPSetOperators(ksp, mat, mat); CHKERRABORT(sys.comm().get(), ierr);
468 ierr = KSPSetFromOptions(ksp); CHKERRABORT(sys.comm().get(), ierr);
471 ierr = KSPGetPC(ksp, &pc); CHKERRABORT(sys.comm().get(), ierr);
472 ierr = PCSetFromOptions(pc); CHKERRABORT(sys.comm().get(), ierr);
475 START_LOG(
"KSPSolve",
"ComplexSolve");
478 ierr = KSPSolve(ksp, res_vec, sol_vec);
480 STOP_LOG(
"KSPSolve",
"ComplexSolve");
484 libMesh::NumericVector<Real>
489 first = sol_R.first_local_index(),
490 last = sol_R.last_local_index();
492 for (
unsigned int i=first; i<last; i++) {
493 sol_R.set(i, (*sol)( 2*i));
494 sol_I.set(i, (*sol)(2*i+1));
501 ierr = KSPDestroy(&ksp); CHKERRABORT(sys.comm().get(), ierr);
502 ierr = MatDestroy(&mat); CHKERRABORT(sys.comm().get(), ierr);
503 ierr = VecDestroy(&res_vec); CHKERRABORT(sys.comm().get(), ierr);
504 ierr = VecDestroy(&sol_vec); CHKERRABORT(sys.comm().get(), ierr);
506 STOP_LOG(
"solve_block_matrix()",
"ComplexSolve");
const MAST::NonlinearSystem & system() const
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...
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.