const Epetra_Comm &
comm_self();
+ /**
+ * Return a Teuchos::Comm object needed for creation of Tpetra::Maps.
+ *
+ * If deal.II has been configured to use a compiler that does not support
+ * MPI then the resulting communicator will be a serial one. Otherwise,
+ * the communicator will correspond to MPI_COMM_SELF, i.e. a communicator
+ * that comprises only this one processor.
+ */
const Teuchos::RCP<const Teuchos::Comm<int>> &
tpetra_comm_self();
Assert(size == 0 || values != nullptr, ExcInternalError("Export failed."));
AssertDimension(size, stored_elements.n_elements());
- if (operation == VectorOperation::insert)
- {
- for (size_type i = 0; i < size; ++i)
- values[i] = new_values[i];
- }
- else if (operation == VectorOperation::add)
- {
- for (size_type i = 0; i < size; ++i)
- values[i] += new_values[i];
- }
- else if (operation == VectorOperation::min)
+ switch (operation)
{
- for (size_type i = 0; i < size; ++i)
- if (std::real(new_values[i]) - std::real(values[i]) < 0.0)
- values[i] = new_values[i];
- }
- else if (operation == VectorOperation::max)
- {
- for (size_type i = 0; i < size; ++i)
- if (std::real(new_values[i]) - std::real(values[i]) > 0.0)
+ case VectorOperation::insert:
+ for (size_type i = 0; i < size; ++i)
values[i] = new_values[i];
+ break;
+
+ case VectorOperation::add:
+ for (size_type i = 0; i < size; ++i)
+ values[i] += new_values[i];
+ break;
+
+ case VectorOperation::min:
+ for (size_type i = 0; i < size; ++i)
+ if (std::real(new_values[i]) - std::real(values[i]) < 0.0)
+ values[i] = new_values[i];
+ break;
+
+ case VectorOperation::max:
+ for (size_type i = 0; i < size; ++i)
+ if (std::real(new_values[i]) - std::real(values[i]) > 0.0)
+ values[i] = new_values[i];
+ break;
+
+ default:
+ AssertThrow(false, ExcNotImplemented());
}
- else
- AssertThrow(false, ExcNotImplemented());
}
* initialized with the same IndexSet that was used for the column indices
* of the matrix.
*
+ * This function will be called when the underlying number type for the
+ * matrix object and the one one for the vector object are the same.
+ *
* In case of a serial vector, this function will only work when
* running on one processor, since the matrix object is inherently
* distributed. Otherwise, an exception will be thrown.
TrilinosScalar>::value>::type
vmult(VectorType &dst, const VectorType &src) const;
+ /**
+ * Same as the function above for the case that the underlying number type
+ * for the matrix object and the one for the vector object do not coincide.
+ * This case is not implemented. Calling it will result in an runtime error.
+ */
template <typename VectorType>
typename std::enable_if<!std::is_same<typename VectorType::value_type,
TrilinosScalar>::value>::type
*
* This function can be called with several types of vector objects,
* see the discussion about @p VectorType in vmult().
+ *
+ * This function will be called when the underlying number type for the
+ * matrix object and the one one for the vector object are the same.
*/
template <typename VectorType>
typename std::enable_if<std::is_same<typename VectorType::value_type,
TrilinosScalar>::value>::type
Tvmult(VectorType &dst, const VectorType &src) const;
+ /**
+ * Same as the function above for the case that the underlying number type
+ * for the matrix object and the one for the vector object do not coincide.
+ * This case is not implemented. Calling it will result in an runtime error.
+ */
template <typename VectorType>
typename std::enable_if<!std::is_same<typename VectorType::value_type,
TrilinosScalar>::value>::type
return communicator;
}
+
+
const Epetra_Comm &
comm_self()
{