61 std::ostringstream error_message;
62 error_message <<
"The fluid and pseudo elastic mesh must be set.";
64 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
68 std::ostringstream error_message;
69 error_message <<
"The solid mesh must be set.";
71 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
75 std::ostringstream error_message;
76 error_message <<
"The Lagrange multiplier mesh must be set.";
78 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
88 unsigned nfluid_dof =
Dim + 1;
91 unsigned npseudo_elastic_dof = this->ndof_types_in_mesh(0) - nfluid_dof;
94 unsigned nsolid_dof = this->ndof_types_in_mesh(1);
97 unsigned nlagr_mult_dof = this->ndof_types_in_mesh(2);
100 unsigned ntotal_dof =
101 nfluid_dof + npseudo_elastic_dof + nsolid_dof + nlagr_mult_dof;
107 Vector<unsigned> dof_to_block_map(ntotal_dof, 0);
109 for (
unsigned i = 0; i < npseudo_elastic_dof; i++)
111 dof_to_block_map[c] = 2;
114 for (
unsigned i = 0; i < nsolid_dof; i++)
116 dof_to_block_map[c] = 1;
119 for (
unsigned i = 0; i < nlagr_mult_dof; i++)
121 dof_to_block_map[c] = 3;
126 CRDoubleMatrix* cr_matrix_pt =
dynamic_cast<CRDoubleMatrix*
>(matrix_pt());
128 if (cr_matrix_pt == 0)
130 std::ostringstream error_message;
131 error_message <<
"FSIPreconditioner only works with"
132 <<
" CRDoubleMatrix matrices" << std::endl;
134 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
139 this->block_setup(dof_to_block_map);
150 Vector<unsigned> ns_dof_list(nfluid_dof, 0);
151 for (
unsigned i = 0; i < nfluid_dof; i++)
157 ->turn_into_subsidiary_block_preconditioner(
this, ns_dof_list);
164 CRDoubleMatrix* ns_matrix_pt =
new CRDoubleMatrix;
165 this->get_block(0, 0, *ns_matrix_pt);
173 if (
dynamic_cast<BlockPreconditioner<CRDoubleMatrix>*
>(
177 GeneralPurposeBlockPreconditioner<CRDoubleMatrix>*
178 solid_block_preconditioner_pt =
179 dynamic_cast<GeneralPurposeBlockPreconditioner<CRDoubleMatrix>*
>(
182 if (solid_block_preconditioner_pt != 0)
184 unsigned offset = nfluid_dof + npseudo_elastic_dof;
185 Vector<unsigned> solid_prec_dof_list(nsolid_dof);
186 for (
unsigned i = 0; i < nsolid_dof; i++)
188 solid_prec_dof_list[i] = offset + i;
190 solid_block_preconditioner_pt
191 ->turn_into_subsidiary_block_preconditioner(
this,
192 solid_prec_dof_list);
193 solid_block_preconditioner_pt->setup(cr_matrix_pt);
197 std::ostringstream error_message;
198 error_message <<
"If the (real) solid preconditioner is a "
199 <<
"BlockPreconditioner then is must be a "
200 <<
"GeneralPurposeBlockPreconditioner";
201 throw OomphLibError(error_message.str(),
202 OOMPH_CURRENT_FUNCTION,
203 OOMPH_EXCEPTION_LOCATION);
210 CRDoubleMatrix* s_matrix_pt =
new CRDoubleMatrix;
211 this->get_block(1, 1, *s_matrix_pt);
218 unsigned ndof_for_pseudo_elastic_prec =
Dim * 3;
219 Vector<unsigned> pseudo_elastic_prec_dof_list(ndof_for_pseudo_elastic_prec,
221 for (
unsigned i = 0; i <
Dim * 2; i++)
223 pseudo_elastic_prec_dof_list[i] = nfluid_dof + i;
225 for (
unsigned i = 0; i <
Dim; i++)
227 pseudo_elastic_prec_dof_list[i +
Dim * 2] =
228 nfluid_dof + npseudo_elastic_dof + nsolid_dof + i;
231 this, pseudo_elastic_prec_dof_list);
242 CRDoubleMatrix* fp_matrix_pt =
new CRDoubleMatrix;
243 get_block(0, 2, *fp_matrix_pt);
245 this->setup_matrix_vector_product(
251 CRDoubleMatrix* sf_matrix_pt =
new CRDoubleMatrix;
252 get_block(1, 0, *sf_matrix_pt);
259 CRDoubleMatrix* sp_matrix_pt =
new CRDoubleMatrix;
260 get_block(1, 2, *sp_matrix_pt);
262 this->setup_matrix_vector_product(
268 CRDoubleMatrix* ls_matrix_pt =
new CRDoubleMatrix;
269 get_block(3, 1, *ls_matrix_pt);
271 this->setup_matrix_vector_product(
281 const DoubleVector& r, DoubleVector& z)
288 this->get_block_vector(2, z, x);
294 this->get_block_vector(0, r, x);
305 this->return_block_vector(0, x, z_copy);
315 this->return_block_vector(0, y, z);
320 this->get_block_vector(0, z, x);
323 this->get_block_vector(1, r, x);
334 DoubleVector z_copy(z);
335 this->return_block_vector(1, x, z_copy);
337 (
dynamic_cast<GeneralPurposeBlockPreconditioner<CRDoubleMatrix>*
>(
339 ->preconditioner_solve(z_copy, z);
340 this->get_block_vector(1, z, y);
346 this->return_block_vector(1, y, z);
352 this->get_block_vector(3, r, y);
356 this->return_block_vector(3, y, z_copy);
MatrixVectorProduct * Solid_pseudo_elastic_matvec_pt
solid onto pseudo solid matrix vector operatio
Mesh * Lagrange_multiplier_mesh_pt
Mesh containing the lagrange multiplier elements.
MatrixVectorProduct * Lagrange_solid_matvec_pt
void clean_up_memory()
Broken assignment operator.
Preconditioner * Solid_preconditioner_pt
pointer to the solid preconditioner
NavierStokesSchurComplementPreconditioner * Navier_stokes_schur_complement_preconditioner_pt
Navier Stokes Schur complement preconditioner.
void setup()
Setup the precoonditioner.
Preconditioner * Navier_stokes_preconditioner_pt
pointer to the navier stokes precondtioner
MatrixVectorProduct * Fluid_pseudo_elastic_matvec_pt
fluid onto pseudosolid matrix vector operator
MatrixVectorProduct * Solid_fluid_matvec_pt
solid onto fluid matrix vector operatio
Mesh * Fluid_and_pseudo_elastic_mesh_pt
Mesh containing the combined fluid and pseudo solid element.
bool Solid_preconditioner_is_block_preconditioner
boolean flag to indicate whether the Solid preconditioner is a block preconditioner
bool Use_navier_stokes_schur_complement_preconditioner
If true the Navier Stokes Schur complement preconditioner is used. Otherwise SuperLUPreconditioner is...
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply the preconditioner.
unsigned Dim
the dimension of the fluid
PseudoElasticPreconditioner * Pseudo_elastic_preconditioner_pt
pointer to the pseudo solid preconditioner
Mesh * Solid_mesh_pt
Mesh containing the solid elements.
void lagrange_multiplier_preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply the lagrange multiplier subsidiary preconditioner.
void set_lagrange_multiplier_mesh(Mesh *mesh_pt)
Access function to mesh containing the block-preconditionable lagrange multiplier elements.
void elastic_preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply the elastic subsidiary preconditioner.
void clean_up_memory()
Clears the memory.
void set_elastic_mesh(Mesh *mesh_pt)
Access function to mesh containing the block-preconditionable elastic elements.