oomph::BlockDiagonalPreconditioner< MATRIX > Class Template Reference

#include <general_purpose_block_preconditioners.h>

+ Inheritance diagram for oomph::BlockDiagonalPreconditioner< MATRIX >:

Public Member Functions

 BlockDiagonalPreconditioner ()
 constructor - when the preconditioner is used as a master preconditioner More...
 
virtual ~BlockDiagonalPreconditioner ()
 Destructor - delete the preconditioner matrices. More...
 
virtual void clean_up_memory ()
 clean up the memory More...
 
 BlockDiagonalPreconditioner (const BlockDiagonalPreconditioner &)=delete
 Broken copy constructor. More...
 
void operator= (const BlockDiagonalPreconditioner &)=delete
 Broken assignment operator. More...
 
void preconditioner_solve (const DoubleVector &r, DoubleVector &z)
 Apply preconditioner to r. More...
 
virtual void setup ()
 Setup the preconditioner. More...
 
void enable_two_level_parallelisation ()
 Use two level parallelisation. More...
 
void disable_two_level_parallelisation ()
 Don't use two-level parallelisation. More...
 
void enable_doc_time_during_preconditioner_solve ()
 Enable Doc timings in application of block sub-preconditioners. More...
 
void disable_doc_time_during_preconditioner_solve ()
 Disable Doc timings in application of block sub-preconditioners. More...
 
void fill_in_subsidiary_preconditioners (const unsigned &nprec_needed)
 
- Public Member Functions inherited from oomph::GeneralPurposeBlockPreconditioner< MATRIX >
 GeneralPurposeBlockPreconditioner ()
 constructor More...
 
virtual ~GeneralPurposeBlockPreconditioner ()
 
 GeneralPurposeBlockPreconditioner (const GeneralPurposeBlockPreconditioner &)=delete
 Broken copy constructor. More...
 
void operator= (const GeneralPurposeBlockPreconditioner &)=delete
 Broken assignment operator. More...
 
void set_subsidiary_preconditioner_function (SubsidiaryPreconditionerFctPt sub_prec_fn)
 access function to set the subsidiary preconditioner function. More...
 
void reset_subsidiary_preconditioner_function_to_default ()
 Reset the subsidiary preconditioner function to its default. More...
 
void set_subsidiary_preconditioner_pt (Preconditioner *prec, const unsigned &i)
 
Preconditionersubsidiary_preconditioner_pt (const unsigned &i) const
 
void set_dof_to_block_map (Vector< unsigned > &dof_to_block_map)
 Specify a DOF to block map. More...
 
void add_mesh (const Mesh *mesh_pt, const bool &allow_multiple_element_type_in_mesh=false)
 
unsigned gp_nmesh ()
 
- Public Member Functions inherited from oomph::BlockPreconditioner< MATRIX >
 BlockPreconditioner ()
 Constructor. More...
 
virtual ~BlockPreconditioner ()
 Destructor. More...
 
 BlockPreconditioner (const BlockPreconditioner &)=delete
 Broken copy constructor. More...
 
void operator= (const BlockPreconditioner &)=delete
 Broken assignment operator. More...
 
MATRIX * matrix_pt () const
 
void turn_on_recursive_debug_flag ()
 
void turn_off_recursive_debug_flag ()
 
void turn_on_debug_flag ()
 Toggles on the debug flag. More...
 
void turn_off_debug_flag ()
 Toggles off the debug flag. More...
 
void turn_into_subsidiary_block_preconditioner (BlockPreconditioner< MATRIX > *master_block_prec_pt, const Vector< unsigned > &doftype_in_master_preconditioner_coarse)
 
void turn_into_subsidiary_block_preconditioner (BlockPreconditioner< MATRIX > *master_block_prec_pt, const Vector< unsigned > &doftype_in_master_preconditioner_coarse, const Vector< Vector< unsigned >> &doftype_coarsen_map_coarse)
 
virtual void block_setup ()
 
void block_setup (const Vector< unsigned > &dof_to_block_map)
 
void get_block (const unsigned &i, const unsigned &j, MATRIX &output_matrix, const bool &ignore_replacement_block=false) const
 
MATRIX get_block (const unsigned &i, const unsigned &j, const bool &ignore_replacement_block=false) const
 
void set_master_matrix_pt (MATRIX *in_matrix_pt)
 Set the matrix_pt in the upper-most master preconditioner. More...
 
void get_block_other_matrix (const unsigned &i, const unsigned &j, MATRIX *in_matrix_pt, MATRIX &output_matrix)
 
void get_blocks (DenseMatrix< bool > &required_blocks, DenseMatrix< MATRIX * > &block_matrix_pt) const
 
void get_dof_level_block (const unsigned &i, const unsigned &j, MATRIX &output_block, const bool &ignore_replacement_block=false) const
 
MATRIX get_concatenated_block (const VectorMatrix< BlockSelector > &selected_block)
 
void get_concatenated_block_vector (const Vector< unsigned > &block_vec_number, const DoubleVector &v, DoubleVector &b)
 
void return_concatenated_block_vector (const Vector< unsigned > &block_vec_number, const DoubleVector &b, DoubleVector &v) const
 Takes concatenated block ordered vector, b, and copies its. More...
 
void get_block_vectors (const Vector< unsigned > &block_vec_number, const DoubleVector &v, Vector< DoubleVector > &s) const
 
void get_block_vectors (const DoubleVector &v, Vector< DoubleVector > &s) const
 
void return_block_vectors (const Vector< unsigned > &block_vec_number, const Vector< DoubleVector > &s, DoubleVector &v) const
 
void return_block_vectors (const Vector< DoubleVector > &s, DoubleVector &v) const
 
void get_block_vector (const unsigned &n, const DoubleVector &v, DoubleVector &b) const
 
void return_block_vector (const unsigned &n, const DoubleVector &b, DoubleVector &v) const
 
void get_block_ordered_preconditioner_vector (const DoubleVector &v, DoubleVector &w)
 
void return_block_ordered_preconditioner_vector (const DoubleVector &w, DoubleVector &v) const
 
unsigned nblock_types () const
 Return the number of block types. More...
 
unsigned ndof_types () const
 Return the total number of DOF types. More...
 
const Meshmesh_pt (const unsigned &i) const
 
unsigned nmesh () const
 
int block_number (const unsigned &i_dof) const
 Return the block number corresponding to a global index i_dof. More...
 
int index_in_block (const unsigned &i_dof) const
 
const LinearAlgebraDistributionblock_distribution_pt (const unsigned &b) const
 Access function to the block distributions (const version). More...
 
LinearAlgebraDistributionblock_distribution_pt (const unsigned b)
 Access function to the block distributions (non-const version). More...
 
LinearAlgebraDistributiondof_block_distribution_pt (const unsigned &b)
 Access function to the dof-level block distributions. More...
 
const LinearAlgebraDistributionmaster_distribution_pt () const
 
unsigned ndof_types_in_mesh (const unsigned &i) const
 
bool is_subsidiary_block_preconditioner () const
 
bool is_master_block_preconditioner () const
 
void set_block_output_to_files (const std::string &basefilename)
 
void disable_block_output_to_files ()
 Turn off output of blocks (by clearing the basefilename string). More...
 
bool block_output_on () const
 Test if output of blocks is on or not. More...
 
void output_blocks_to_files (const std::string &basefilename, const unsigned &precision=8) const
 
void post_block_matrix_assembly_partial_clear ()
 
BlockPreconditioner< MATRIX > * master_block_preconditioner_pt () const
 Access function to the master block preconditioner pt. More...
 
void clear_block_preconditioner_base ()
 
void document ()
 
Vector< Vector< unsigned > > doftype_coarsen_map_fine () const
 
Vector< unsignedget_fine_grain_dof_types_in (const unsigned &i) const
 
unsigned nfine_grain_dof_types_in (const unsigned &i) const
 
MapMatrix< unsigned, CRDoubleMatrix * > replacement_dof_block_pt () const
 Access function to the replaced dof-level blocks. More...
 
void setup_matrix_vector_product (MatrixVectorProduct *matvec_prod_pt, CRDoubleMatrix *block_pt, const Vector< unsigned > &block_col_indices)
 
void setup_matrix_vector_product (MatrixVectorProduct *matvec_prod_pt, CRDoubleMatrix *block_pt, const unsigned &block_col_index)
 
void internal_get_block_ordered_preconditioner_vector (const DoubleVector &v, DoubleVector &w) const
 
void internal_return_block_ordered_preconditioner_vector (const DoubleVector &w, DoubleVector &v) const
 
unsigned internal_nblock_types () const
 
unsigned internal_ndof_types () const
 
void internal_return_block_vector (const unsigned &n, const DoubleVector &b, DoubleVector &v) const
 
void internal_get_block_vector (const unsigned &n, const DoubleVector &v, DoubleVector &b) const
 
void internal_get_block_vectors (const Vector< unsigned > &block_vec_number, const DoubleVector &v, Vector< DoubleVector > &s) const
 
void internal_get_block_vectors (const DoubleVector &v, Vector< DoubleVector > &s) const
 
void internal_return_block_vectors (const Vector< unsigned > &block_vec_number, const Vector< DoubleVector > &s, DoubleVector &v) const
 
void internal_return_block_vectors (const Vector< DoubleVector > &s, DoubleVector &v) const
 
void internal_get_block (const unsigned &i, const unsigned &j, MATRIX &output_block) const
 
int internal_block_number (const unsigned &i_dof) const
 
int internal_index_in_block (const unsigned &i_dof) const
 
const LinearAlgebraDistributioninternal_block_distribution_pt (const unsigned &b) const
 Access function to the internal block distributions. More...
 
void insert_auxiliary_block_distribution (const Vector< unsigned > &block_vec_number, LinearAlgebraDistribution *dist_pt)
 
void block_matrix_test (const unsigned &i, const unsigned &j, const MATRIX *block_matrix_pt) const
 
template<typename myType >
int get_index_of_value (const Vector< myType > &vec, const myType val, const bool sorted=false) const
 
void internal_get_block (const unsigned &block_i, const unsigned &block_j, CRDoubleMatrix &output_block) const
 
void get_dof_level_block (const unsigned &block_i, const unsigned &block_j, CRDoubleMatrix &output_block, const bool &ignore_replacement_block) const
 
- Public Member Functions inherited from oomph::Preconditioner
 Preconditioner ()
 Constructor. More...
 
 Preconditioner (const Preconditioner &)=delete
 Broken copy constructor. More...
 
void operator= (const Preconditioner &)=delete
 Broken assignment operator. More...
 
virtual ~Preconditioner ()
 Destructor (empty) More...
 
virtual void preconditioner_solve_transpose (const DoubleVector &r, DoubleVector &z)
 
void setup (DoubleMatrixBase *matrix_pt)
 
void setup (const Problem *problem_pt, DoubleMatrixBase *matrix_pt)
 
void enable_silent_preconditioner_setup ()
 Set up the block preconditioner quietly! More...
 
void disable_silent_preconditioner_setup ()
 Be verbose in the block preconditioner setup. More...
 
virtual void set_matrix_pt (DoubleMatrixBase *matrix_pt)
 Set the matrix pointer. More...
 
virtual const OomphCommunicatorcomm_pt () const
 Get function for comm pointer. More...
 
virtual void set_comm_pt (const OomphCommunicator *const comm_pt)
 Set the communicator pointer. More...
 
double setup_time () const
 Returns the time to setup the preconditioner. More...
 
virtual void turn_into_subsidiary_block_preconditioner (BlockPreconditioner< CRDoubleMatrix > *master_block_prec_pt, const Vector< unsigned > &doftype_in_master_preconditioner_coarse)
 
virtual void turn_into_subsidiary_block_preconditioner (BlockPreconditioner< CRDoubleMatrix > *master_block_prec_pt, const Vector< unsigned > &doftype_in_master_preconditioner_coarse, const Vector< Vector< unsigned >> &doftype_coarsen_map_coarse)
 
- Public Member Functions inherited from oomph::DistributableLinearAlgebraObject
 DistributableLinearAlgebraObject ()
 Default constructor - create a distribution. More...
 
 DistributableLinearAlgebraObject (const DistributableLinearAlgebraObject &matrix)=delete
 Broken copy constructor. More...
 
void operator= (const DistributableLinearAlgebraObject &)=delete
 Broken assignment operator. More...
 
virtual ~DistributableLinearAlgebraObject ()
 Destructor. More...
 
LinearAlgebraDistributiondistribution_pt () const
 access to the LinearAlgebraDistribution More...
 
unsigned nrow () const
 access function to the number of global rows. More...
 
unsigned nrow_local () const
 access function for the num of local rows on this processor. More...
 
unsigned nrow_local (const unsigned &p) const
 access function for the num of local rows on this processor. More...
 
unsigned first_row () const
 access function for the first row on this processor More...
 
unsigned first_row (const unsigned &p) const
 access function for the first row on this processor More...
 
bool distributed () const
 distribution is serial or distributed More...
 
bool distribution_built () const
 
void build_distribution (const LinearAlgebraDistribution *const dist_pt)
 
void build_distribution (const LinearAlgebraDistribution &dist)
 

Protected Member Functions

virtual unsigned get_other_diag_ds (const unsigned &i, const unsigned &nblock) const
 
- Protected Member Functions inherited from oomph::GeneralPurposeBlockPreconditioner< MATRIX >
void gp_preconditioner_set_all_meshes ()
 Set the mesh in the block preconditioning framework. More...
 
void gp_preconditioner_block_setup ()
 Modified block setup for general purpose block preconditioners. More...
 
void fill_in_subsidiary_preconditioners (const unsigned &nprec_needed)
 
- Protected Member Functions inherited from oomph::BlockPreconditioner< MATRIX >
void set_nmesh (const unsigned &n)
 
void set_mesh (const unsigned &i, const Mesh *const mesh_pt, const bool &allow_multiple_element_type_in_mesh=false)
 
void set_replacement_dof_block (const unsigned &block_i, const unsigned &block_j, CRDoubleMatrix *replacement_dof_block_pt)
 
bool any_mesh_distributed () const
 
int internal_dof_number (const unsigned &i_dof) const
 
unsigned internal_index_in_dof (const unsigned &i_dof) const
 
unsigned internal_block_dimension (const unsigned &b) const
 
unsigned internal_dof_block_dimension (const unsigned &i) const
 
unsigned master_nrow () const
 
unsigned internal_master_dof_number (const unsigned &b) const
 
const LinearAlgebraDistributioninternal_preconditioner_matrix_distribution_pt () const
 
const LinearAlgebraDistributionpreconditioner_matrix_distribution_pt () const
 
- Protected Member Functions inherited from oomph::DistributableLinearAlgebraObject
void clear_distribution ()
 

Private Attributes

PreconditionerArrayPreconditioner_array_pt
 pointer for the PreconditionerArray More...
 
bool Use_two_level_parallelisation
 Use two level parallelism using the PreconditionerArray. More...
 
bool Doc_time_during_preconditioner_solve
 Doc timings in application of block sub-preconditioners? More...
 

Additional Inherited Members

- Public Types inherited from oomph::GeneralPurposeBlockPreconditioner< MATRIX >
typedef Preconditioner *(* SubsidiaryPreconditionerFctPt) ()
 
- Protected Attributes inherited from oomph::GeneralPurposeBlockPreconditioner< MATRIX >
Vector< Preconditioner * > Subsidiary_preconditioner_pt
 List of preconditioners to use for the blocks to be solved. More...
 
SubsidiaryPreconditionerFctPt Subsidiary_preconditioner_creation_function_pt
 
- Protected Attributes inherited from oomph::BlockPreconditioner< MATRIX >
MapMatrix< unsigned, CRDoubleMatrix * > Replacement_dof_block_pt
 The replacement dof-level blocks. More...
 
Vector< LinearAlgebraDistribution * > Block_distribution_pt
 The distribution for the blocks. More...
 
Vector< Vector< unsigned > > Block_to_dof_map_coarse
 
Vector< Vector< unsigned > > Block_to_dof_map_fine
 Mapping for the block types to the most fine grain dof types. More...
 
Vector< Vector< unsigned > > Doftype_coarsen_map_coarse
 
Vector< Vector< unsigned > > Doftype_coarsen_map_fine
 
Vector< LinearAlgebraDistribution * > Internal_block_distribution_pt
 Storage for the default distribution for each internal block. More...
 
Vector< LinearAlgebraDistribution * > Dof_block_distribution_pt
 
Vector< unsignedAllow_multiple_element_type_in_mesh
 
Vector< const Mesh * > Mesh_pt
 
Vector< unsignedNdof_types_in_mesh
 
unsigned Internal_nblock_types
 
unsigned Internal_ndof_types
 
- Protected Attributes inherited from oomph::Preconditioner
bool Silent_preconditioner_setup
 Boolean to indicate whether or not the build should be done silently. More...
 
std::ostream * Stream_pt
 Pointer to the output stream – defaults to std::cout. More...
 

Detailed Description

template<typename MATRIX>
class oomph::BlockDiagonalPreconditioner< MATRIX >

Block diagonal preconditioner. By default SuperLU is used to solve the subsidiary systems, but other preconditioners can be used by setting them using passing a pointer to a function of type SubsidiaryPreconditionerFctPt to the method subsidiary_preconditioner_function_pt().

Constructor & Destructor Documentation

◆ BlockDiagonalPreconditioner() [1/2]

template<typename MATRIX >
oomph::BlockDiagonalPreconditioner< MATRIX >::BlockDiagonalPreconditioner ( )
inline

constructor - when the preconditioner is used as a master preconditioner

324  : GeneralPurposeBlockPreconditioner<MATRIX>()
325  {
326  // by default we do not use two level parallelism
328 
329  // null the Preconditioner array pt
331 
332  // Don't doc by default
334  }
bool Doc_time_during_preconditioner_solve
Doc timings in application of block sub-preconditioners?
Definition: general_purpose_block_preconditioners.h:438
PreconditionerArray * Preconditioner_array_pt
pointer for the PreconditionerArray
Definition: general_purpose_block_preconditioners.h:432
bool Use_two_level_parallelisation
Use two level parallelism using the PreconditionerArray.
Definition: general_purpose_block_preconditioners.h:435

◆ ~BlockDiagonalPreconditioner()

template<typename MATRIX >
virtual oomph::BlockDiagonalPreconditioner< MATRIX >::~BlockDiagonalPreconditioner ( )
inlinevirtual

Destructor - delete the preconditioner matrices.

338  {
339  this->clean_up_memory();
340  }
virtual void clean_up_memory()
clean up the memory
Definition: general_purpose_block_preconditioners.h:343

References oomph::GeneralPurposeBlockPreconditioner< MATRIX >::clean_up_memory().

◆ BlockDiagonalPreconditioner() [2/2]

template<typename MATRIX >
oomph::BlockDiagonalPreconditioner< MATRIX >::BlockDiagonalPreconditioner ( const BlockDiagonalPreconditioner< MATRIX > &  )
delete

Broken copy constructor.

Member Function Documentation

◆ clean_up_memory()

template<typename MATRIX >
virtual void oomph::BlockDiagonalPreconditioner< MATRIX >::clean_up_memory ( )
inlinevirtual

clean up the memory

Reimplemented from oomph::GeneralPurposeBlockPreconditioner< MATRIX >.

344  {
346  {
349  }
350 
351  // Clean up the base class too
353  }
virtual void clean_up_memory()
Definition: general_purpose_block_preconditioners.h:112

References oomph::GeneralPurposeBlockPreconditioner< MATRIX >::clean_up_memory().

◆ disable_doc_time_during_preconditioner_solve()

template<typename MATRIX >
void oomph::BlockDiagonalPreconditioner< MATRIX >::disable_doc_time_during_preconditioner_solve ( )
inline

Disable Doc timings in application of block sub-preconditioners.

392  {
394  }

◆ disable_two_level_parallelisation()

template<typename MATRIX >
void oomph::BlockDiagonalPreconditioner< MATRIX >::disable_two_level_parallelisation ( )
inline

Don't use two-level parallelisation.

380  {
382  }

◆ enable_doc_time_during_preconditioner_solve()

template<typename MATRIX >
void oomph::BlockDiagonalPreconditioner< MATRIX >::enable_doc_time_during_preconditioner_solve ( )
inline

Enable Doc timings in application of block sub-preconditioners.

386  {
388  }

◆ enable_two_level_parallelisation()

template<typename MATRIX >
void oomph::BlockDiagonalPreconditioner< MATRIX >::enable_two_level_parallelisation ( )
inline

Use two level parallelisation.

369  {
370 #ifndef OOMPH_HAS_MPI
371  throw OomphLibError("Cannot do any parallelism since we don't have MPI.",
374 #endif
376  }
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ fill_in_subsidiary_preconditioners()

template<typename MATRIX >
void oomph::BlockDiagonalPreconditioner< MATRIX >::fill_in_subsidiary_preconditioners ( const unsigned nprec_needed)
inline
397  {
398 #ifdef PARANOID
400  !this->Subsidiary_preconditioner_pt.empty())
401  {
402  std::string err_msg =
403  "Two level parallelism diagonal block preconditioners cannot have";
404  err_msg +=
405  " any preset preconditioners (due to weird memory management";
406  err_msg += " in the PreconditionerArray, you could try fixing it).";
407  throw OomphLibError(
409  }
410 #endif
411 
412  // Now call the real function
414  MATRIX>::fill_in_subsidiary_preconditioners(nprec_needed);
415  }
void fill_in_subsidiary_preconditioners(const unsigned &nprec_needed)
Definition: general_purpose_block_preconditioners.h:396
GeneralPurposeBlockPreconditioner()
constructor
Definition: general_purpose_block_preconditioners.h:87
Vector< Preconditioner * > Subsidiary_preconditioner_pt
List of preconditioners to use for the blocks to be solved.
Definition: general_purpose_block_preconditioners.h:294
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286

References oomph::GeneralPurposeBlockPreconditioner< MATRIX >::fill_in_subsidiary_preconditioners(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::Global_string_for_annotation::string(), and oomph::GeneralPurposeBlockPreconditioner< MATRIX >::Subsidiary_preconditioner_pt.

◆ get_other_diag_ds()

template<typename MATRIX >
virtual unsigned oomph::BlockDiagonalPreconditioner< MATRIX >::get_other_diag_ds ( const unsigned i,
const unsigned nblock 
) const
inlineprotectedvirtual

This is a helper function to allow us to implement AntiDiagonal preconditioner by only changing this function. Get the second index for block number i. Obviously for a diagonal preconditioner we want the blocks (i,i), (for anti diagonal we will want blocks (i, nblock - i), see that class).

Reimplemented in oomph::BlockAntiDiagonalPreconditioner< MATRIX >.

425  {
426  return i;
427  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9

References i.

◆ operator=()

template<typename MATRIX >
void oomph::BlockDiagonalPreconditioner< MATRIX >::operator= ( const BlockDiagonalPreconditioner< MATRIX > &  )
delete

Broken assignment operator.

◆ preconditioner_solve()

template<typename MATRIX >
void oomph::BlockDiagonalPreconditioner< MATRIX >::preconditioner_solve ( const DoubleVector r,
DoubleVector z 
)
virtual

Apply preconditioner to r.

Preconditioner solve for the block diagonal preconditioner.

Implements oomph::Preconditioner.

213  {
214  // Cache umber of block types
215  unsigned n_block = this->nblock_types();
216 
217  // Get the right hand side vector (b in Ax = b) in block form.
218  Vector<DoubleVector> block_r;
219  this->get_block_vectors(r, block_r);
220 
221  // Make sure z vector is built
222  if (!z.built())
223  {
224  z.build(this->distribution_pt(), 0.0);
225  }
226 
227  // Vector of vectors for storage of solution in block form.
228  Vector<DoubleVector> block_z(n_block);
229 
231  {
233  }
234  else
235  {
236  // solve each diagonal block
237  for (unsigned i = 0; i < n_block; i++)
238  {
239  double t_start = 0.0;
241  {
242  t_start = TimingHelpers::timer();
243  }
244 
245  // this->Subsidiary_preconditioner_pt[i]->preconditioner_solve(block_r[i],
246  // block_z[i]);
247 
248 
249  // See if the i-th subsidiary preconditioner is a block preconditioner
250  if (dynamic_cast<BlockPreconditioner<CRDoubleMatrix>*>(
251  this->Subsidiary_preconditioner_pt[i]) == 0)
252  {
253  this->Subsidiary_preconditioner_pt[i]->preconditioner_solve(
254  block_r[i], block_z[i]);
255  }
256  // If the subsidiary preconditioner is a block preconditioner
257  else
258  {
259  // This is pretty inefficient but there is no other way to do this
260  // with block sub preconditioners as far as I can tell because they
261  // demand to have the full r and z vectors...
262  DoubleVector block_z_with_size_of_full_z(r.distribution_pt());
263  DoubleVector r_updated(r.distribution_pt());
264 
265  // Construct the new r vector (with the update given by back subs
266  // below).
267  this->return_block_vectors(block_r, r_updated);
268 
269  // Solve (blocking is handled inside the block preconditioner).
270  this->Subsidiary_preconditioner_pt[i]->preconditioner_solve(
271  r_updated, block_z_with_size_of_full_z);
272 
273  // Extract this block's z vector because we need it for the next step
274  // (and throw away block_z_with_size_of_full_z).
275  this->get_block_vector(i, block_z_with_size_of_full_z, block_z[i]);
276  }
277 
278 
280  {
281  oomph_info << "Time for application of " << i
282  << "-th block preconditioner: "
283  << TimingHelpers::timer() - t_start << std::endl;
284  }
285  }
286  }
287 
288  // copy solution in block vectors block_z back to z
289  this->return_block_vectors(block_z, z);
290  }
std::vector< double > DoubleVector
loads clump configuration
Definition: ClumpInput.h:26
void return_block_vectors(const Vector< unsigned > &block_vec_number, const Vector< DoubleVector > &s, DoubleVector &v) const
Definition: block_preconditioner.cc:3443
void get_block_vector(const unsigned &n, const DoubleVector &v, DoubleVector &b) const
Definition: block_preconditioner.cc:4186
unsigned nblock_types() const
Return the number of block types.
Definition: block_preconditioner.h:1670
void get_block_vectors(const Vector< unsigned > &block_vec_number, const DoubleVector &v, Vector< DoubleVector > &s) const
Definition: block_preconditioner.cc:2939
LinearAlgebraDistribution * distribution_pt() const
access to the LinearAlgebraDistribution
Definition: linear_algebra_distribution.h:457
void solve_preconditioners(const Vector< DoubleVector > &r, Vector< DoubleVector > &z)
Definition: preconditioner_array.h:61
r
Definition: UniformPSDSelfTest.py:20
double timer()
returns the time in seconds after some point in past
Definition: oomph_utilities.cc:1295
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References oomph::DoubleVector::build(), oomph::DoubleVector::built(), i, oomph::oomph_info, UniformPSDSelfTest::r, and oomph::TimingHelpers::timer().

◆ setup()

template<typename MATRIX >
void oomph::BlockDiagonalPreconditioner< MATRIX >::setup
virtual

Setup the preconditioner.

setup for the block diagonal preconditioner

Implements oomph::Preconditioner.

40  {
41  // clean the memory
42  this->clean_up_memory();
43 
44  // Note: Subsidiary block preconditioners don't really use their
45  // mesh pointers (since the lookup schemes inferred from them are
46  // set up by their masters). Generally we insist that (for uniformity)
47  // a mesh should be set, but sometimes subsidiary preconditioners are
48  // set up on the fly and we can't sensibly set meshes, so we're
49  // forgiving...)
51  {
52 #ifdef PARANOID
53  if (this->gp_nmesh() == 0)
54  {
55  std::ostringstream err_msg;
56  err_msg << "There are no meshes set.\n"
57  << "Did you remember to call add_mesh(...)?";
58  throw OomphLibError(
60  }
61 #endif
62 
63  // Set all meshes if this is master block preconditioner
65  }
66 
67  // If we're meant to build silently
68  if (this->Silent_preconditioner_setup == true)
69  {
70  // Store the output stream pointer
71  this->Stream_pt = oomph_info.stream_pt();
72 
73  // Now set the oomph_info stream pointer to the null stream to
74  // disable all possible output
76  }
77 
78  // Set up the block look up schemes
80 
81  // number of types of blocks.
82  unsigned nblock_types = this->nblock_types();
83 
84  // Create any subsidiary preconditioners needed
85  this->fill_in_subsidiary_preconditioners(nblock_types);
86 
87  // The total time for extracting all the blocks from the "global" matrix
88  double t_extraction_total = 0.0;
89 
90  // The total time for setting up the subsidiary preconditioners
91  double t_subsidiary_setup_total = 0.0;
92 
93  // If using two level parallelisation then we need to use a
94  // PrecondtionerArray which requires very different setup. ??ds possibly
95  // it should have it's own class?
97  {
98  // Get the blocks. We have to use new because you can't have containers
99  // of matrices because the copy constructor is buggy (so we create a
100  // container of pointers instead). ??ds
101  Vector<CRDoubleMatrix*> block_diagonal_matrix_pt(nblock_types, 0);
102  for (unsigned i = 0; i < nblock_types; i++)
103  {
104  // Allocate space for the new matrix
105  block_diagonal_matrix_pt[i] = new CRDoubleMatrix;
106 
107  // Get the start time
108  double t_extract_start = TimingHelpers::timer();
109 
110  // Extract the i-th block
111  this->get_block(
112  i, get_other_diag_ds(i, nblock_types), *block_diagonal_matrix_pt[i]);
113 
114  // Get the end time
115  double t_extract_end = TimingHelpers::timer();
116 
117  // Update the timing total
118  t_extraction_total += (t_extract_end - t_extract_start);
119  }
120 
121  // Construct the preconditioner array
122  Preconditioner_array_pt = new PreconditionerArray;
123 
124  // Get the start time
125  double t_subsidiary_setup_start = TimingHelpers::timer();
126 
127  // Set up the preconditioners
129  block_diagonal_matrix_pt,
131  this->comm_pt());
132 
133  // Get the end time
134  double t_subsidiary_setup_end = TimingHelpers::timer();
135 
136  // Update the timing total
137  t_subsidiary_setup_total +=
138  (t_subsidiary_setup_end - t_subsidiary_setup_start);
139 
140  // and delete the blocks
141  for (unsigned i = 0; i < nblock_types; i++)
142  {
143  delete block_diagonal_matrix_pt[i];
144  block_diagonal_matrix_pt[i] = 0;
145  }
146 
147  // Preconditioner array is weird, it calls delete on all the
148  // preconditioners you give it and requires new ones each time!
149  this->Subsidiary_preconditioner_pt.clear();
150  }
151  // Otherwise just set up each block's preconditioner in order
152  else
153  {
154  for (unsigned i = 0; i < nblock_types; i++)
155  {
156  // Get the block
157  unsigned j = get_other_diag_ds(i, nblock_types);
158 
159  // Get the start time
160  double t_extract_start = TimingHelpers::timer();
161 
162  // Get the block
163  CRDoubleMatrix block = this->get_block(i, j);
164 
165  // Get the end time
166  double t_extract_end = TimingHelpers::timer();
167 
168  // Update the timing total
169  t_extraction_total += (t_extract_end - t_extract_start);
170 
171  // Get the start time
172  double t_subsidiary_setup_start = TimingHelpers::timer();
173 
174  // Set the (subsidiary) preconditioner up for this block
175  this->Subsidiary_preconditioner_pt[i]->setup(&block);
176 
177  // Get the end time
178  double t_subsidiary_setup_end = TimingHelpers::timer();
179 
180  // Tell the user
181  oomph_info << "Took "
182  << t_subsidiary_setup_end - t_subsidiary_setup_start
183  << "s to setup." << std::endl;
184 
185  // Update the timing total
186  t_subsidiary_setup_total +=
187  (t_subsidiary_setup_end - t_subsidiary_setup_start);
188  }
189  }
190 
191  // If we're meant to build silently, reassign the oomph stream pointer
192  if (this->Silent_preconditioner_setup == true)
193  {
194  // Store the output stream pointer
195  oomph_info.stream_pt() = this->Stream_pt;
196 
197  // Reset our own stream pointer
198  this->Stream_pt = 0;
199  }
200 
201  // Tell the user
202  oomph_info << "Total block extraction time [sec]: " << t_extraction_total
203  << "\nTotal subsidiary preconditioner setup time [sec]: "
204  << t_subsidiary_setup_total << std::endl;
205  }
m m block(1, 0, 2, 2)<< 4
virtual unsigned get_other_diag_ds(const unsigned &i, const unsigned &nblock) const
Definition: general_purpose_block_preconditioners.h:423
void get_block(const unsigned &i, const unsigned &j, MATRIX &output_matrix, const bool &ignore_replacement_block=false) const
Definition: block_preconditioner.h:673
bool is_master_block_preconditioner() const
Definition: block_preconditioner.h:2068
void gp_preconditioner_set_all_meshes()
Set the mesh in the block preconditioning framework.
Definition: general_purpose_block_preconditioners.h:218
void gp_preconditioner_block_setup()
Modified block setup for general purpose block preconditioners.
Definition: general_purpose_block_preconditioners.h:241
unsigned gp_nmesh()
Definition: general_purpose_block_preconditioners.h:211
std::ostream *& stream_pt()
Access function for the stream pointer.
Definition: oomph_definitions.h:464
void setup_preconditioners(Vector< CRDoubleMatrix * > matrix_pt, Vector< Preconditioner * > prec_pt, const OomphCommunicator *comm_pt)
Definition: preconditioner_array.h:52
std::ostream * Stream_pt
Pointer to the output stream – defaults to std::cout.
Definition: preconditioner.h:232
bool Silent_preconditioner_setup
Boolean to indicate whether or not the build should be done silently.
Definition: preconditioner.h:229
virtual const OomphCommunicator * comm_pt() const
Get function for comm pointer.
Definition: preconditioner.h:171
Nullstream oomph_nullstream
Single (global) instantiation of the Nullstream.
Definition: oomph_definitions.cc:313
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References block(), oomph::TerminateHelper::clean_up_memory(), i, j, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::oomph_info, oomph::oomph_nullstream, GlobalParameters::Stream_pt, oomph::OomphInfo::stream_pt(), and oomph::TimingHelpers::timer().

Member Data Documentation

◆ Doc_time_during_preconditioner_solve

template<typename MATRIX >
bool oomph::BlockDiagonalPreconditioner< MATRIX >::Doc_time_during_preconditioner_solve
private

Doc timings in application of block sub-preconditioners?

◆ Preconditioner_array_pt

template<typename MATRIX >
PreconditionerArray* oomph::BlockDiagonalPreconditioner< MATRIX >::Preconditioner_array_pt
private

pointer for the PreconditionerArray

◆ Use_two_level_parallelisation

template<typename MATRIX >
bool oomph::BlockDiagonalPreconditioner< MATRIX >::Use_two_level_parallelisation
private

Use two level parallelism using the PreconditionerArray.


The documentation for this class was generated from the following files: