From: David Wells Date: Sat, 27 May 2017 21:56:09 +0000 (-0400) Subject: Move Trilinos index functions into a new header file. X-Git-Tag: v9.0.0-rc1~1562^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=5a5ec17249438a4b95cd30a9625aec2a57a11a40;p=dealii.git Move Trilinos index functions into a new header file. --- diff --git a/doc/news/changes/minor/20170527DavidWells b/doc/news/changes/minor/20170527DavidWells new file mode 100644 index 0000000000..5ed8c7565e --- /dev/null +++ b/doc/news/changes/minor/20170527DavidWells @@ -0,0 +1,4 @@ +New: There is a new header, lac/trilinos_index_access.h, that provides +index and size functions for some common Trilinos objects that work correctly for +both 32 and 64 bit code. +
(David Wells, 2017/05/10) diff --git a/include/deal.II/lac/trilinos_index_access.h b/include/deal.II/lac/trilinos_index_access.h new file mode 100644 index 0000000000..c57d42cad2 --- /dev/null +++ b/include/deal.II/lac/trilinos_index_access.h @@ -0,0 +1,222 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2017 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + +#ifndef dealii__trilinos_index_access_h +#define dealii__trilinos_index_access_h + +#include + +#ifdef DEAL_II_WITH_TRILINOS + +#include + +#include +#include +#include +#include + +DEAL_II_NAMESPACE_OPEN + +namespace TrilinosWrappers +{ + /** + * A helper function that queries the size of an Epetra_BlockMap object + * and calls either the 32 or 64 bit function to get the number of global + * elements in the map. + */ + inline + TrilinosWrappers::types::int_type + n_global_elements(const Epetra_BlockMap &map) + { +#ifdef DEAL_II_WITH_64BIT_INDICES + return map.NumGlobalElements64(); +#else + return map.NumGlobalElements(); +#endif + } + + /** + * A helper function that finds the minimum global index value on the + * calling processor by calling either the 32 or 64 bit function. + */ + inline + TrilinosWrappers::types::int_type + min_my_gid(const Epetra_BlockMap &map) + { +#ifdef DEAL_II_WITH_64BIT_INDICES + return map.MinMyGID64(); +#else + return map.MinMyGID(); +#endif + } + + /** + * A helper function that finds the maximum global index value on the + * calling processor by calling either the 32 or 64 bit function. + */ + inline + TrilinosWrappers::types::int_type + max_my_gid(const Epetra_BlockMap &map) + { +#ifdef DEAL_II_WITH_64BIT_INDICES + return map.MaxMyGID64(); +#else + return map.MaxMyGID(); +#endif + } + + /** + * A helper function that converts a local index to a global one calling + * either the 32 or 64 bit function. + */ + inline + TrilinosWrappers::types::int_type + global_index(const Epetra_BlockMap &map, + const dealii::types::global_dof_index i) + { +#ifdef DEAL_II_WITH_64BIT_INDICES + return map.GID64(i); +#else + return map.GID(i); +#endif + } + + /** + * A helper function that returns a pointer to the array containing the + * global indices assigned to the current process by calling either the 32 + * or 64 bit function. + */ + inline + TrilinosWrappers::types::int_type * + my_global_elements(const Epetra_BlockMap &map) + { +#ifdef DEAL_II_WITH_64BIT_INDICES + return map.MyGlobalElements64(); +#else + return map.MyGlobalElements(); +#endif + } + + /** + * A helper function that finds the global number of rows by calling + * either the 32 or 64 bit function. + */ + inline + TrilinosWrappers::types::int_type + n_global_rows(const Epetra_CrsGraph &graph) + { +#ifdef DEAL_II_WITH_64BIT_INDICES + return graph.NumGlobalRows64(); +#else + return graph.NumGlobalRows(); +#endif + } + + /** + * A helper function that finds the global number of columns by calling + * either the 32 or 64 bit function. + */ + inline + TrilinosWrappers::types::int_type + n_global_cols(const Epetra_CrsGraph &graph) + { +#ifdef DEAL_II_WITH_64BIT_INDICES + return graph.NumGlobalCols64(); +#else + return graph.NumGlobalCols(); +#endif + } + + /** + * A helper function that finds the number of global entries by calling + * either the 32 or 64 bit function. + */ + inline + TrilinosWrappers::types::int_type + n_global_entries(const Epetra_CrsGraph &graph) + { +#ifdef DEAL_II_WITH_64BIT_INDICES + return graph.NumGlobalEntries64(); +#else + return graph.NumGlobalEntries(); +#endif + } + + /** + * A helper function that finds the global row index by calling + * either the 32 or 64 bit function. + */ + inline + TrilinosWrappers::types::int_type + global_row_index(const Epetra_CrsMatrix &matrix, + const dealii::types::global_dof_index i) + { +#ifdef DEAL_II_WITH_64BIT_INDICES + return matrix.GRID64(i); +#else + return matrix.GRID(i); +#endif + } + + /** + * A helper function that finds the global column index by calling + * either the 32 or 64 bit function. + */ + inline + TrilinosWrappers::types::int_type + global_column_index(const Epetra_CrsMatrix &matrix, + const dealii::types::global_dof_index i) + { +#ifdef DEAL_II_WITH_64BIT_INDICES + return matrix.GCID64(i); +#else + return matrix.GCID(i); +#endif + } + + /** + * A helper function that finds the global length of a vector by calling + * either the 32 or 64 bit function. + */ + inline + TrilinosWrappers::types::int_type + global_length (const Epetra_MultiVector &vector) + { +#ifdef DEAL_II_WITH_64BIT_INDICES + return vector.GlobalLength64(); +#else + return vector.GlobalLength(); +#endif + } + + /** + * A helper function that finds the global number of rows by calling + * either the 32 or 64 bit function. + */ + inline + TrilinosWrappers::types::int_type + n_global_rows(const Epetra_RowMatrix &matrix) + { +#ifdef DEAL_II_WITH_64BIT_INDICES + return matrix.NumGlobalRows64(); +#else + return matrix.NumGlobalRows(); +#endif + } +} + +DEAL_II_NAMESPACE_CLOSE +#endif // DEAL_II_WITH_TRILINOS +#endif // dealii__trilinos_index_access_h diff --git a/source/lac/trilinos_block_vector.cc b/source/lac/trilinos_block_vector.cc index 3f0b2fc64a..ed1e7bcd38 100644 --- a/source/lac/trilinos_block_vector.cc +++ b/source/lac/trilinos_block_vector.cc @@ -17,33 +17,15 @@ #ifdef DEAL_II_WITH_TRILINOS -# include +#include + +#include DEAL_II_NAMESPACE_OPEN namespace TrilinosWrappers { - namespace - { - // define a helper function that queries the size of an Epetra_Map object - // by calling either the 32- or 64-bit function necessary, and returns the - // result in the correct data type so that we can use it in calling other - // Epetra member functions that are overloaded by index type -#ifndef DEAL_II_WITH_64BIT_INDICES - int n_global_elements (const Epetra_BlockMap &map) - { - return map.NumGlobalElements(); - } -#else - long long int n_global_elements (const Epetra_BlockMap &map) - { - return map.NumGlobalElements64(); - } -#endif - } - - namespace MPI { BlockVector & @@ -248,7 +230,7 @@ namespace TrilinosWrappers std::vector block_sizes (no_blocks); for (size_type i=0; iblock_indices.reinit (block_sizes); diff --git a/source/lac/trilinos_precondition_ml.cc b/source/lac/trilinos_precondition_ml.cc index cba0f83216..3912e39480 100644 --- a/source/lac/trilinos_precondition_ml.cc +++ b/source/lac/trilinos_precondition_ml.cc @@ -13,6 +13,7 @@ // // --------------------------------------------------------------------- +#include #include #ifdef DEAL_II_WITH_TRILINOS @@ -20,6 +21,7 @@ # include # include # include +# include DEAL_II_DISABLE_EXTRA_DIAGNOSTICS # include @@ -35,43 +37,6 @@ DEAL_II_NAMESPACE_OPEN namespace TrilinosWrappers { - namespace - { -#ifndef DEAL_II_WITH_64BIT_INDICES - int n_global_rows (const Epetra_RowMatrix &matrix) - { - return matrix.NumGlobalRows(); - } - - int global_length (const Epetra_MultiVector &vector) - { - return vector.GlobalLength(); - } - - int gid(const Epetra_Map &map, unsigned int i) - { - return map.GID(i); - } -#else - long long int n_global_rows (const Epetra_RowMatrix &matrix) - { - return matrix.NumGlobalRows64(); - } - - long long int global_length (const Epetra_MultiVector &vector) - { - return vector.GlobalLength64(); - } - - long long int gid(const Epetra_Map &map, dealii::types::global_dof_index i) - { - return map.GID64(i); - } -#endif - } - - - /* -------------------------- PreconditionAMG -------------------------- */ PreconditionAMG::AdditionalData:: @@ -188,12 +153,12 @@ namespace TrilinosWrappers if (constant_modes_dimension > 0) { - const size_type global_size = n_global_rows(matrix); + const size_type global_size = TrilinosWrappers::n_global_rows(matrix); (void)global_length; // work around compiler warning about unused function in release mode Assert (global_size == - static_cast(global_length(distributed_constant_modes)), + static_cast(TrilinosWrappers::global_length(distributed_constant_modes)), ExcDimensionMismatch(global_size, - global_length(distributed_constant_modes))); + TrilinosWrappers::global_length(distributed_constant_modes))); const bool constant_modes_are_global = additional_data.constant_modes[0].size() == global_size; const size_type my_size = domain_map.NumMyElements(); @@ -209,7 +174,7 @@ namespace TrilinosWrappers for (size_type row=0; row #include #ifdef DEAL_II_WITH_TRILINOS @@ -20,6 +21,7 @@ # include # include +# include # include DEAL_II_DISABLE_EXTRA_DIAGNOSTICS @@ -38,43 +40,6 @@ DEAL_II_NAMESPACE_OPEN namespace TrilinosWrappers { - namespace - { -#ifndef DEAL_II_WITH_64BIT_INDICES - int n_global_rows (const Epetra_RowMatrix &matrix) - { - return matrix.NumGlobalRows(); - } - - int global_length (const Epetra_MultiVector &vector) - { - return vector.GlobalLength(); - } - - int gid(const Epetra_Map &map, unsigned int i) - { - return map.GID(i); - } -#else - long long int n_global_rows (const Epetra_RowMatrix &matrix) - { - return matrix.NumGlobalRows64(); - } - - long long int global_length (const Epetra_MultiVector &vector) - { - return vector.GlobalLength64(); - } - - long long int gid(const Epetra_Map &map, dealii::types::global_dof_index i) - { - return map.GID64(i); - } -#endif - } - - - PreconditionAMGMueLu::AdditionalData:: AdditionalData (const bool elliptic, const unsigned int n_cycles, @@ -181,7 +146,7 @@ namespace TrilinosWrappers if (constant_modes_dimension > 0) { - const size_type n_rows = n_global_rows(matrix); + const size_type n_rows = TrilinosWrappers::n_global_rows(matrix); const bool constant_modes_are_global = additional_data.constant_modes[0].size() == n_rows; const size_type n_relevant_rows = @@ -191,9 +156,9 @@ namespace TrilinosWrappers Assert (n_relevant_rows == my_size, ExcDimensionMismatch(n_relevant_rows, my_size)); Assert (n_rows == - static_cast(global_length(distributed_constant_modes)), + static_cast(TrilinosWrappers::global_length(distributed_constant_modes)), ExcDimensionMismatch(n_rows, - global_length(distributed_constant_modes))); + TrilinosWrappers::global_length(distributed_constant_modes))); (void)n_relevant_rows; (void)global_length; @@ -204,7 +169,7 @@ namespace TrilinosWrappers for (size_type row=0; row #include #ifdef DEAL_II_WITH_TRILINOS @@ -24,6 +25,7 @@ # include # include # include +# include # include DEAL_II_DISABLE_EXTRA_DIAGNOSTICS @@ -37,79 +39,6 @@ DEAL_II_NAMESPACE_OPEN namespace TrilinosWrappers { - namespace - { -#ifndef DEAL_II_WITH_64BIT_INDICES - // define a helper function that queries the size of an Epetra_Map object - // by calling either the 32- or 64-bit function necessary, and returns the - // result in the correct data type so that we can use it in calling other - // Epetra member functions that are overloaded by index type - int n_global_elements (const Epetra_BlockMap &map) - { - return map.NumGlobalElements(); - } - - int min_my_gid(const Epetra_BlockMap &map) - { - return map.MinMyGID(); - } - - int max_my_gid(const Epetra_BlockMap &map) - { - return map.MaxMyGID(); - } - - int n_global_cols(const Epetra_CrsGraph &graph) - { - return graph.NumGlobalCols(); - } - - int global_column_index(const Epetra_CrsMatrix &matrix, int i) - { - return matrix.GCID(i); - } - - int global_row_index(const Epetra_CrsMatrix &matrix, int i) - { - return matrix.GRID(i); - } -#else - // define a helper function that queries the size of an Epetra_Map object - // by calling either the 32- or 64-bit function necessary, and returns the - // result in the correct data type so that we can use it in calling other - // Epetra member functions that are overloaded by index type - long long int n_global_elements (const Epetra_BlockMap &map) - { - return map.NumGlobalElements64(); - } - - long long int min_my_gid(const Epetra_BlockMap &map) - { - return map.MinMyGID64(); - } - - long long int max_my_gid(const Epetra_BlockMap &map) - { - return map.MaxMyGID64(); - } - - long long int n_global_cols(const Epetra_CrsGraph &graph) - { - return graph.NumGlobalCols64(); - } - - long long int global_column_index(const Epetra_CrsMatrix &matrix, int i) - { - return matrix.GCID64(i); - } - - long long int global_row_index(const Epetra_CrsMatrix &matrix, int i) - { - return matrix.GRID64(i); - } -#endif - } - namespace internal { template @@ -526,9 +455,9 @@ namespace TrilinosWrappers if (input_row_map.Comm().MyPID() == 0) { AssertDimension (sparsity_pattern.n_rows(), - static_cast(n_global_elements(input_row_map))); + static_cast(TrilinosWrappers::n_global_elements(input_row_map))); AssertDimension (sparsity_pattern.n_cols(), - static_cast(n_global_elements(input_col_map))); + static_cast(TrilinosWrappers::n_global_elements(input_col_map))); } column_space_map.reset (new Epetra_Map (input_col_map)); @@ -548,8 +477,8 @@ namespace TrilinosWrappers return; } - const size_type first_row = min_my_gid(input_row_map), - last_row = max_my_gid(input_row_map)+1; + const size_type first_row = TrilinosWrappers::min_my_gid(input_row_map), + last_row = TrilinosWrappers::max_my_gid(input_row_map)+1; std::vector n_entries_per_row(last_row-first_row); for (size_type row=first_row; row(n_global_cols(*graph))); + static_cast(TrilinosWrappers::n_global_cols(*graph))); (void)n_global_cols; // And now finally generate the matrix. @@ -659,9 +588,9 @@ namespace TrilinosWrappers nonlocal_matrix_exporter.reset(); AssertDimension (sparsity_pattern.n_rows(), - static_cast(n_global_elements(input_row_map))); + static_cast(TrilinosWrappers::n_global_elements(input_row_map))); AssertDimension (sparsity_pattern.n_cols(), - static_cast(n_global_elements(input_col_map))); + static_cast(TrilinosWrappers::n_global_elements(input_col_map))); column_space_map.reset (new Epetra_Map (input_col_map)); @@ -669,8 +598,8 @@ namespace TrilinosWrappers // serial case if (relevant_rows.size() == 0) { - relevant_rows.set_size(n_global_elements(input_row_map)); - relevant_rows.add_range(0, n_global_elements(input_row_map)); + relevant_rows.set_size(TrilinosWrappers::n_global_elements(input_row_map)); + relevant_rows.add_range(0, TrilinosWrappers::n_global_elements(input_row_map)); } relevant_rows.compress(); Assert(relevant_rows.n_elements() >= static_cast(input_row_map.NumMyElements()), @@ -786,7 +715,7 @@ namespace TrilinosWrappers graph->OptimizeStorage(); AssertDimension (sparsity_pattern.n_cols(),static_cast( - n_global_cols(*graph))); + TrilinosWrappers::n_global_cols(*graph))); matrix.reset (new Epetra_FECrsMatrix(Copy, *graph, false)); } @@ -2244,45 +2173,38 @@ namespace TrilinosWrappers inputleft.range_partitioner()); Assert (inputleft.domain_partitioner().LinearMap() == true, ExcMessage("Matrix must be partitioned contiguously between procs.")); - for (unsigned int i=0; i= 0, ExcInternalError()); -#ifndef DEAL_II_WITH_64BIT_INDICES - const size_type GID = inputleft.trilinos_matrix().RowMap().GID(i); - for (TrilinosWrappers::types::int_type j=0; j= 0, ExcInternalError()); -#ifndef DEAL_II_WITH_64BIT_INDICES - const size_type GID = inputleft.trilinos_matrix().RowMap().GID(i); - for (TrilinosWrappers::types::int_type j=0; jExtractMyRowView (i, num_entries, values, indices); for (TrilinosWrappers::types::int_type j=0; j #include #ifdef DEAL_II_WITH_TRILINOS # include # include -# include + # include +# include +# include DEAL_II_DISABLE_EXTRA_DIAGNOSTICS # include @@ -30,85 +33,6 @@ DEAL_II_NAMESPACE_OPEN namespace TrilinosWrappers { - namespace - { - // define a helper function that queries the size of an Epetra_Map object - // by calling either the 32- or 64-bit function necessary, and returns the - // result in the correct data type so that we can use it in calling other - // Epetra member functions that are overloaded by index type -#ifndef DEAL_II_WITH_64BIT_INDICES - int n_global_elements (const Epetra_BlockMap &map) - { - return map.NumGlobalElements(); - } - - int min_my_gid(const Epetra_BlockMap &map) - { - return map.MinMyGID(); - } - - int max_my_gid(const Epetra_BlockMap &map) - { - return map.MaxMyGID(); - } - - int n_global_rows(const Epetra_CrsGraph &graph) - { - return graph.NumGlobalRows(); - } - - int n_global_cols(const Epetra_CrsGraph &graph) - { - return graph.NumGlobalCols(); - } - - int n_global_entries(const Epetra_CrsGraph &graph) - { - return graph.NumGlobalEntries(); - } - - int local_to_global_index(const Epetra_BlockMap &map, int i) - { - return map.GID(i); - } -#else - long long int n_global_elements (const Epetra_BlockMap &map) - { - return map.NumGlobalElements64(); - } - - long long int min_my_gid(const Epetra_BlockMap &map) - { - return map.MinMyGID64(); - } - - long long int max_my_gid(const Epetra_BlockMap &map) - { - return map.MaxMyGID64(); - } - - long long int n_global_rows(const Epetra_CrsGraph &graph) - { - return graph.NumGlobalRows64(); - } - - long long int n_global_cols(const Epetra_CrsGraph &graph) - { - return graph.NumGlobalCols64(); - } - - long long int n_global_entries(const Epetra_CrsGraph &graph) - { - return graph.NumGlobalEntries64(); - } - - long long int local_to_global_index(const Epetra_BlockMap &map, int i) - { - return map.GID64(i); - } -#endif - } - namespace SparsityPatternIterators { void @@ -392,13 +316,13 @@ namespace TrilinosWrappers nonlocal_graph.reset(); graph.reset (); AssertDimension (n_entries_per_row.size(), - static_cast(n_global_elements(row_map))); + static_cast(TrilinosWrappers::n_global_elements(row_map))); column_space_map.reset (new Epetra_Map (col_map)); - std::vector local_entries_per_row(max_my_gid(row_map)- - min_my_gid(row_map)); + std::vector local_entries_per_row(TrilinosWrappers::max_my_gid(row_map)- + TrilinosWrappers::min_my_gid(row_map)); for (unsigned int i=0; i 1) graph.reset(new Epetra_FECrsGraph(Copy, row_map, @@ -434,17 +358,17 @@ namespace TrilinosWrappers graph.reset (); AssertDimension (sp.n_rows(), - static_cast(n_global_elements(row_map))); + static_cast(TrilinosWrappers::n_global_elements(row_map))); AssertDimension (sp.n_cols(), - static_cast(n_global_elements(col_map))); + static_cast(TrilinosWrappers::n_global_elements(col_map))); column_space_map.reset (new Epetra_Map (col_map)); Assert (row_map.LinearMap() == true, ExcMessage ("This function only works if the row map is contiguous.")); - const size_type first_row = min_my_gid(row_map), - last_row = max_my_gid(row_map)+1; + const size_type first_row = TrilinosWrappers::min_my_gid(row_map), + last_row = TrilinosWrappers::max_my_gid(row_map)+1; std::vector n_entries_per_row(last_row - first_row); // Trilinos wants the row length as an int this is hopefully never going @@ -932,7 +856,7 @@ namespace TrilinosWrappers if (graph->Filled() == true) n_cols = n_global_cols(*graph); else - n_cols = n_global_elements(*column_space_map); + n_cols = TrilinosWrappers::n_global_elements(*column_space_map); return n_cols; } @@ -953,8 +877,8 @@ namespace TrilinosWrappers SparsityPattern::local_range () const { size_type begin, end; - begin = min_my_gid(graph->RowMap()); - end = max_my_gid(graph->RowMap())+1; + begin = TrilinosWrappers::min_my_gid(graph->RowMap()); + end = TrilinosWrappers::max_my_gid(graph->RowMap())+1; return std::make_pair (begin, end); } @@ -1089,8 +1013,8 @@ namespace TrilinosWrappers { graph->ExtractMyRowView (i, num_entries, indices); for (int j=0; jRowMap(), i) - << "," << local_to_global_index(graph->ColMap(), indices[j]) << ") " + out << "(" << TrilinosWrappers::global_index(graph->RowMap(), i) + << "," << TrilinosWrappers::global_index(graph->ColMap(), indices[j]) << ") " << std::endl; } } @@ -1104,20 +1028,23 @@ namespace TrilinosWrappers SparsityPattern::print_gnuplot (std::ostream &out) const { Assert (graph->Filled() == true, ExcInternalError()); - for (unsigned int row=0; rowExtractMyRowView (row, num_entries, indices); - for (int j=0; j= 0, ExcInternalError()); + // avoid sign comparison warning + const dealii::types::global_dof_index num_entries_ = num_entries; + for (dealii::types::global_dof_index j=0; j(local_to_global_index(graph->ColMap(), indices[j])) - << " " << -static_cast(local_to_global_index(graph->RowMap(), row)) << std::endl; + out << static_cast(TrilinosWrappers::global_index(graph->ColMap(), indices[j])) + << " " << -static_cast(TrilinosWrappers::global_index(graph->RowMap(), row)) << std::endl; } AssertThrow (out, ExcIO()); diff --git a/source/lac/trilinos_vector.cc b/source/lac/trilinos_vector.cc index ac96692027..9f33b388ce 100644 --- a/source/lac/trilinos_vector.cc +++ b/source/lac/trilinos_vector.cc @@ -13,13 +13,15 @@ // // --------------------------------------------------------------------- +#include #include #ifdef DEAL_II_WITH_TRILINOS # include -# include # include +# include +# include DEAL_II_DISABLE_EXTRA_DIAGNOSTICS # include @@ -33,65 +35,8 @@ DEAL_II_NAMESPACE_OPEN namespace TrilinosWrappers { - namespace - { -#ifndef DEAL_II_WITH_64BIT_INDICES - // define a helper function that queries the size of an Epetra_BlockMap object - // by calling either the 32- or 64-bit function necessary, and returns the - // result in the correct data type so that we can use it in calling other - // Epetra member functions that are overloaded by index type - int n_global_elements (const Epetra_BlockMap &map) - { - return map.NumGlobalElements(); - } - // define a helper function that queries the pointer to internal array - // containing list of global IDs assigned to the calling processor - // by calling either the 32- or 64-bit function necessary, and returns the - // result in the correct data type so that we can use it in calling other - // Epetra member functions that are overloaded by index type - int *my_global_elements(const Epetra_BlockMap &map) - { - return map.MyGlobalElements(); - } - // define a helper function that queries the global vector length of an - // Epetra_FEVector object by calling either the 32- or 64-bit - // function necessary. - int global_length(const Epetra_FEVector &vector) - { - return vector.GlobalLength(); - } -#else - // define a helper function that queries the size of an Epetra_BlockMap object - // by calling either the 32- or 64-bit function necessary, and returns the - // result in the correct data type so that we can use it in calling other - // Epetra member functions that are overloaded by index type - long long int n_global_elements (const Epetra_BlockMap &map) - { - return map.NumGlobalElements64(); - } - // define a helper function that queries the pointer to internal array - // containing list of global IDs assigned to the calling processor - // by calling either the 32- or 64-bit function necessary, and returns the - // result in the correct data type so that we can use it in calling other - // Epetra member functions that are overloaded by index type - long long int *my_global_elements(const Epetra_BlockMap &map) - { - return map.MyGlobalElements64(); - } - // define a helper function that queries the global vector length of an - // Epetra_FEVector object by calling either the 32- or 64-bit - // function necessary. - long long int global_length(const Epetra_FEVector &vector) - { - return vector.GlobalLength64(); - } -#endif - } - namespace MPI { - - Vector::Vector () { last_action = Zero; @@ -139,9 +84,9 @@ namespace TrilinosWrappers VectorBase() { AssertThrow (parallel_partitioner.size() == - static_cast(n_global_elements(v.vector->Map())), + static_cast(TrilinosWrappers::n_global_elements(v.vector->Map())), ExcDimensionMismatch (parallel_partitioner.size(), - n_global_elements(v.vector->Map()))); + TrilinosWrappers::n_global_elements(v.vector->Map()))); last_action = Zero; @@ -307,7 +252,7 @@ namespace TrilinosWrappers for (size_type block=0; block(global_length(*actual_vec)) + AssertThrow (static_cast(TrilinosWrappers::global_length(*actual_vec)) == v.size(), - ExcDimensionMismatch (global_length(*actual_vec), + ExcDimensionMismatch (TrilinosWrappers::global_length(*actual_vec), v.size())); Epetra_Import data_exchange (vector->Map(), actual_vec->Map()); @@ -565,11 +510,11 @@ namespace TrilinosWrappers Vector::reinit (const Epetra_Map &input_map, const bool /*omit_zeroing_entries*/) { - Epetra_LocalMap map (n_global_elements(input_map), + Epetra_LocalMap map (TrilinosWrappers::n_global_elements(input_map), input_map.IndexBase(), input_map.Comm()); vector.reset (new Epetra_FEVector(map)); - owned_elements = complete_index_set(n_global_elements(input_map)); + owned_elements = complete_index_set(TrilinosWrappers::n_global_elements(input_map)); last_action = Zero; } @@ -622,7 +567,7 @@ namespace TrilinosWrappers #endif if (!same_communicators || local_range() != v.local_range()) { - Epetra_LocalMap map (global_length(*(v.vector)), + Epetra_LocalMap map (TrilinosWrappers::global_length(*(v.vector)), v.vector->Map().IndexBase(), v.vector->Comm()); vector.reset (new Epetra_FEVector(map)); @@ -679,7 +624,7 @@ namespace TrilinosWrappers { if (size() != v.size()) { - Epetra_LocalMap map (n_global_elements(v.vector->Map()), + Epetra_LocalMap map (TrilinosWrappers::n_global_elements(v.vector->Map()), v.vector->Map().IndexBase(), v.vector->Comm()); vector.reset (new Epetra_FEVector(map)); @@ -696,7 +641,7 @@ namespace TrilinosWrappers { if (size() != v.size()) { - Epetra_LocalMap map (n_global_elements(v.vector->Map()), + Epetra_LocalMap map (TrilinosWrappers::n_global_elements(v.vector->Map()), v.vector->Map().IndexBase(), v.vector->Comm()); vector.reset (new Epetra_FEVector(map));