Loading [MathJax]/extensions/tex2jax.js
OGS
NumLib::detail Namespace Reference

Detailed Description

Classes

struct  ByComponent
 
struct  ByGlobalIndex
 
struct  ByLocation
 
struct  ByLocationAndComponent
 
struct  ComputeDerivativeWrtOneScalar_CD
 
struct  ComputeDerivativeWrtOneScalar_FD
 
struct  DefaultPerturbationStrategy
 
struct  FieldType
 
struct  IsScalar
 
struct  IsScalar< Eigen::Matrix< double, N, 1, Eigen::ColMajor, N, 1 > >
 
struct  Line
 
struct  LineByLocationAndComponentComparator
 
struct  LineByLocationComparator
 
struct  LowerOrderShapeFunctionOrSame
 
struct  LowerOrderShapeFunctionOrSame< ShapeFunction, std::void_t< typename NumLib::LowerDim< ShapeFunction >::type > >
 
struct  NumberOfDofs
 
struct  NumberOfDofs< Vectorial< ShapeFunction, Dim > >
 
struct  ShapeDataFieldType
 
struct  ShapeFunctionTraits
 

Typedefs

using ComponentGlobalIndexDict
 
template<typename ElementTraitsLagrange >
using GetMeshElement = typename ElementTraitsLagrange::Element
 
template<typename ElementTraitsLagrange >
using GetShapeFunctionHigherOrder
 
template<typename ElementTraitsLagrange >
using GetShapeFunctionLowerOrder
 
template<typename ShapeFunction >
using GetShapeMatrixPolicy
 
template<typename ShapeMatrixPolicy >
using GetShapeMatrix_N = typename ShapeMatrixPolicy::ShapeMatrices::ShapeType
 

Functions

template<typename... Ns_t, typename ElementDOFVector , int... Offsets, int... Sizes>
auto localDOFImpl (ElementDOFVector const &x, boost::mp11::mp_list_c< int, Offsets... >, boost::mp11::mp_list_c< int, Sizes... >)
 
template<class T_SHAPE_FUNC , class T_SHAPE_MATRICES >
void computeMappingMatrices (const MeshLib::Element &, const double *natural_pt, const MeshLib::ElementCoordinatesMappingLocal &, T_SHAPE_MATRICES &shapemat, FieldType< ShapeMatrixType::N >)
 
template<class T_SHAPE_FUNC , class T_SHAPE_MATRICES >
void computeMappingMatrices (const MeshLib::Element &, const double *natural_pt, const MeshLib::ElementCoordinatesMappingLocal &, T_SHAPE_MATRICES &shapemat, FieldType< ShapeMatrixType::DNDR >)
 
static void checkJacobianDeterminant (const double detJ, MeshLib::Element const &element)
 
template<class T_SHAPE_FUNC , class T_SHAPE_MATRICES >
void computeMappingMatrices (const MeshLib::Element &ele, const double *natural_pt, const MeshLib::ElementCoordinatesMappingLocal &ele_local_coord, T_SHAPE_MATRICES &shapemat, FieldType< ShapeMatrixType::DNDR_J >)
 
template<class T_SHAPE_FUNC , class T_SHAPE_MATRICES >
void computeMappingMatrices (const MeshLib::Element &ele, const double *natural_pt, const MeshLib::ElementCoordinatesMappingLocal &ele_local_coord, T_SHAPE_MATRICES &shapemat, FieldType< ShapeMatrixType::N_J >)
 
template<class T_SHAPE_FUNC , class T_SHAPE_MATRICES >
void computeMappingMatrices (const MeshLib::Element &ele, const double *natural_pt, const MeshLib::ElementCoordinatesMappingLocal &ele_local_coord, T_SHAPE_MATRICES &shapemat, FieldType< ShapeMatrixType::DNDX >)
 
template<class T_SHAPE_FUNC , class T_SHAPE_MATRICES >
void computeMappingMatrices (const MeshLib::Element &ele, const double *natural_pt, const MeshLib::ElementCoordinatesMappingLocal &ele_local_coord, T_SHAPE_MATRICES &shapemat, FieldType< ShapeMatrixType::ALL >)
 
template<class T_SHAPE_FUNC , class T_SHAPE_MATRICES , ShapeMatrixType T_SHAPE_MATRIX_TYPE>
void naturalCoordinatesMappingComputeShapeMatrices (const MeshLib::Element &ele, const double *natural_pt, T_SHAPE_MATRICES &shapemat, const unsigned global_dim)
 Used to explicitly instantiate the NaturalCoordinatesMapping class template.
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN (ShapeHex20)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN (ShapeHex8)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN (ShapeLine2)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN (ShapeLine3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN (ShapePoint1)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN (ShapePrism15)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN (ShapePrism6)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN (ShapePyra13)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN (ShapePyra5)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN (ShapeQuad4)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN (ShapeQuad8)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN (ShapeQuad9)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN (ShapeTet10)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN (ShapeTet4)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN (ShapeTri3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN (ShapeTri6)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapePoint1, 0)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeLine2, 1)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeLine3, 1)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapePoint1, 1)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeLine2, 2)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeLine3, 2)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapePoint1, 2)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeQuad4, 2)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeQuad8, 2)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeQuad9, 2)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeTri3, 2)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeTri6, 2)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeHex20, 3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeHex8, 3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeLine2, 3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeLine3, 3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapePoint1, 3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapePrism15, 3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapePrism6, 3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapePyra13, 3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapePyra5, 3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeQuad4, 3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeQuad8, 3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeQuad9, 3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeTet10, 3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeTet4, 3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeTri3, 3)
 
 OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX (ShapeTri6, 3)
 
template<class T >
void setMatrixZero (T &mat)
 
template<class T >
void setVectorZero (T &vec)
 
template<class T_N , class T_DNDR , class T_J , class T_DNDX >
void setZero (ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX > &shape, ShapeDataFieldType< ShapeMatrixType::N >)
 
template<class T_N , class T_DNDR , class T_J , class T_DNDX >
void setZero (ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX > &shape, ShapeDataFieldType< ShapeMatrixType::DNDR >)
 
template<class T_N , class T_DNDR , class T_J , class T_DNDX >
void setZero (ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX > &shape, ShapeDataFieldType< ShapeMatrixType::DNDR_J >)
 
template<class T_N , class T_DNDR , class T_J , class T_DNDX >
void setZero (ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX > &shape, ShapeDataFieldType< ShapeMatrixType::N_J >)
 
template<class T_N , class T_DNDR , class T_J , class T_DNDX >
void setZero (ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX > &shape, ShapeDataFieldType< ShapeMatrixType::DNDX >)
 
template<class T_N , class T_DNDR , class T_J , class T_DNDX >
void setZero (ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX > &shape, ShapeDataFieldType< ShapeMatrixType::ALL >)
 
template<unsigned DOFOffset, typename NodalValues , typename ShapeMatrix >
void shapeFunctionInterpolate (const NodalValues &, const ShapeMatrix &)
 
template<unsigned DOFOffset, typename NodalValues , typename ShapeMatrix , typename... ScalarTypes>
void shapeFunctionInterpolate (const NodalValues &nodal_values, const ShapeMatrix &shape_matrix_N, double &interpolated_value, ScalarTypes &... interpolated_values)
 
template<std::size_t IndexInTuple, typename Tuple >
double getScalarOrVectorComponent (Tuple const &tuple, Eigen::Index component)
 
template<typename MeshElementType , typename IPData , typename FluxVectorType , typename Derived >
void assembleAdvectionMatrix (IPData const &ip_data_vector, NumLib::ShapeMatrixCache const &shape_matrix_cache, std::vector< FluxVectorType > const &ip_flux_vector, Eigen::MatrixBase< Derived > &laplacian_matrix)
 
template<typename IPData , typename FluxVectorType , typename Derived >
void assembleAdvectionMatrix (IPData const &ip_data_vector, std::vector< FluxVectorType > const &ip_flux_vector, Eigen::MatrixBase< Derived > &laplacian_matrix)
 
template<typename Derived >
void applyFullUpwind (Eigen::VectorXd const &quasi_nodal_flux, Eigen::MatrixBase< Derived > &laplacian_matrix)
 
template<typename IPData , typename FluxVectorType , typename Derived >
void applyFullUpwind (IPData const &ip_data_vector, std::vector< FluxVectorType > const &ip_flux_vector, Eigen::MatrixBase< Derived > &laplacian_matrix)
 
template<typename MatrixVectorType >
std::unique_ptr< MatrixVectorType > newZeroedInstance (MathLib::MatrixSpecifications const &matrix_specification)
 
void finalize (Mat &M)
 
void calculateFluxCorrectedTransportPETSc (const double t, const double dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, const MathLib::MatrixSpecifications &matrix_specification, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)
 
void calculateFluxCorrectedTransport (const double t, const double dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, const MathLib::MatrixSpecifications &matrix_specification, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)
 
Eigen::MatrixXd getHydrodynamicDispersion (Eigen::MatrixXd const &pore_diffusion_coefficient, Eigen::VectorXd const &velocity, double const porosity, double const solute_dispersivity_transverse, double const solute_dispersivity_longitudinal)
 
Eigen::MatrixXd getHydrodynamicDispersionWithArtificialDiffusion (IsotropicDiffusionStabilization const &stabilizer, std::size_t const element_id, Eigen::MatrixXd const &pore_diffusion_coefficient, Eigen::VectorXd const &velocity, double const porosity, double const solute_dispersivity_transverse, double const solute_dispersivity_longitudinal)
 
bool solvePicard (GlobalLinearSolver &linear_solver, GlobalMatrix &A, GlobalVector &rhs, GlobalVector &x, MathLib::LinearSolverBehaviour const linear_solver_behaviour)
 

Typedef Documentation

◆ ComponentGlobalIndexDict

Initial value:
boost::multi_index::multi_index_container<
Line, boost::multi_index::indexed_by<
boost::multi_index::ordered_unique<
boost::multi_index::tag<ByLocationAndComponent>,
boost::multi_index::identity<Line>,
LineByLocationAndComponentComparator>,
boost::multi_index::ordered_non_unique<
boost::multi_index::tag<ByLocation>,
boost::multi_index::identity<Line>, LineByLocationComparator>,
boost::multi_index::ordered_non_unique<
boost::multi_index::tag<ByComponent>,
boost::multi_index::member<Line, int, &Line::comp_id>>,
boost::multi_index::ordered_non_unique<
boost::multi_index::tag<ByGlobalIndex>,
boost::multi_index::member<Line, GlobalIndexType,
&Line::global_index>>>>
GlobalMatrix::IndexType GlobalIndexType

Definition at line 104 of file ComponentGlobalIndexDict.h.

◆ GetMeshElement

template<typename ElementTraitsLagrange >
using NumLib::detail::GetMeshElement = typename ElementTraitsLagrange::Element

Definition at line 23 of file ShapeMatrixCache.h.

◆ GetShapeFunctionHigherOrder

Initial value:
typename ElementTraitsLagrange::ShapeFunction

Definition at line 26 of file ShapeMatrixCache.h.

◆ GetShapeFunctionLowerOrder

Initial value:
typename ElementTraitsLagrange::LowerOrderShapeFunction

Definition at line 30 of file ShapeMatrixCache.h.

◆ GetShapeMatrix_N

template<typename ShapeMatrixPolicy >
using NumLib::detail::GetShapeMatrix_N = typename ShapeMatrixPolicy::ShapeMatrices::ShapeType

Definition at line 41 of file ShapeMatrixCache.h.

◆ GetShapeMatrixPolicy

template<typename ShapeFunction >
using NumLib::detail::GetShapeMatrixPolicy

Function Documentation

◆ applyFullUpwind() [1/2]

template<typename Derived >
void NumLib::detail::applyFullUpwind ( Eigen::VectorXd const & quasi_nodal_flux,
Eigen::MatrixBase< Derived > & laplacian_matrix )

Definition at line 68 of file AdvectionMatrixAssembler.h.

70{
71 Eigen::VectorXd const down_mask =
72 (quasi_nodal_flux.array() < 0).cast<double>();
73 Eigen::VectorXd const down = quasi_nodal_flux.cwiseProduct(down_mask);
74
75 double const q_in = -down.sum();
76 if (q_in < std::numeric_limits<double>::epsilon())
77 {
78 return;
79 }
80
81 Eigen::VectorXd const up_mask =
82 (quasi_nodal_flux.array() >= 0).cast<double>();
83 Eigen::VectorXd const up = quasi_nodal_flux.cwiseProduct(up_mask);
84
85 laplacian_matrix.diagonal().noalias() += up;
86 laplacian_matrix.noalias() += down * up.transpose() / q_in;
87}

Referenced by applyFullUpwind(), NumLib::assembleAdvectionMatrix(), and NumLib::assembleAdvectionMatrix().

◆ applyFullUpwind() [2/2]

template<typename IPData , typename FluxVectorType , typename Derived >
void NumLib::detail::applyFullUpwind ( IPData const & ip_data_vector,
std::vector< FluxVectorType > const & ip_flux_vector,
Eigen::MatrixBase< Derived > & laplacian_matrix )

Definition at line 90 of file AdvectionMatrixAssembler.h.

93{
94 Eigen::VectorXd quasi_nodal_flux =
95 Eigen::VectorXd::Zero(laplacian_matrix.rows());
96
97 for (std::size_t ip = 0; ip < ip_flux_vector.size(); ++ip)
98 {
99 auto const& ip_data = ip_data_vector[ip];
100 auto const w = ip_data.integration_weight;
101 auto const& dNdx = ip_data.dNdx;
102
103 quasi_nodal_flux.noalias() -= ip_flux_vector[ip].transpose() * dNdx * w;
104 }
105
106 applyFullUpwind(quasi_nodal_flux, laplacian_matrix);
107}
void applyFullUpwind(Eigen::VectorXd const &quasi_nodal_flux, Eigen::MatrixBase< Derived > &laplacian_matrix)

References applyFullUpwind().

◆ assembleAdvectionMatrix() [1/2]

template<typename MeshElementType , typename IPData , typename FluxVectorType , typename Derived >
void NumLib::detail::assembleAdvectionMatrix ( IPData const & ip_data_vector,
NumLib::ShapeMatrixCache const & shape_matrix_cache,
std::vector< FluxVectorType > const & ip_flux_vector,
Eigen::MatrixBase< Derived > & laplacian_matrix )

Definition at line 28 of file AdvectionMatrixAssembler.h.

32{
33 auto const& Ns = shape_matrix_cache.NsHigherOrder<MeshElementType>();
34
35 for (std::size_t ip = 0; ip < ip_flux_vector.size(); ++ip)
36 {
37 auto const& ip_data = ip_data_vector[ip];
38 auto const w = ip_data.integration_weight;
39 auto const& dNdx = ip_data.dNdx;
40 auto const& N = Ns[ip];
41 laplacian_matrix.noalias() +=
42 (N.transpose() * ip_flux_vector[ip].transpose() * dNdx * w)
43#if defined(__GNUC__) && (__GNUC__ == 14) && (__GNUC_MINOR__ == 2) && \
44 (__GNUC_PATCHLEVEL__ == 1)
45 .eval()
46#endif
47 ;
48 }
49}

References NumLib::N, and NumLib::ShapeMatrixCache::NsHigherOrder().

Referenced by ProcessLib::HT::MonolithicHTFEM< ShapeFunction, GlobalDim >::assemble(), NumLib::assembleAdvectionMatrix(), NumLib::assembleAdvectionMatrix(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::assembleBlockMatrices(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::assembleComponentTransportEquation(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::assembleHeatTransportEquation(), ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::assembleHeatTransportEquation(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobian(), and ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::assembleWithJacobianComponentTransportEquation().

◆ assembleAdvectionMatrix() [2/2]

template<typename IPData , typename FluxVectorType , typename Derived >
void NumLib::detail::assembleAdvectionMatrix ( IPData const & ip_data_vector,
std::vector< FluxVectorType > const & ip_flux_vector,
Eigen::MatrixBase< Derived > & laplacian_matrix )

Definition at line 52 of file AdvectionMatrixAssembler.h.

55{
56 for (std::size_t ip = 0; ip < ip_flux_vector.size(); ++ip)
57 {
58 auto const& ip_data = ip_data_vector[ip];
59 auto const w = ip_data.integration_weight;
60 auto const& N = ip_data.N;
61 auto const& dNdx = ip_data.dNdx;
62 laplacian_matrix.noalias() +=
63 N.transpose() * ip_flux_vector[ip].transpose() * dNdx * w;
64 }
65}

References NumLib::N.

◆ calculateFluxCorrectedTransport()

void NumLib::detail::calculateFluxCorrectedTransport ( const double t,
const double dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id,
const MathLib::MatrixSpecifications & matrix_specification,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b )

Definition at line 529 of file FluxCorrectedTransport.cpp.

534{
535#ifndef USE_PETSC
536 calculateFluxCorrectedTransportSerial(t, dt, x, x_prev, process_id,
537 matrix_specification, M, K, b);
538#endif // end of ifndef USE_PETSC
539
540#ifdef USE_PETSC
541 calculateFluxCorrectedTransportPETSc(t, dt, x, x_prev, process_id,
542 matrix_specification, M, K, b);
543#endif // end of ifdef USE_PETSC
544}
void calculateFluxCorrectedTransportPETSc(const double t, const double dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, const MathLib::MatrixSpecifications &matrix_specification, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)

References calculateFluxCorrectedTransportPETSc().

Referenced by NumLib::computeFluxCorrectedTransport().

◆ calculateFluxCorrectedTransportPETSc()

void NumLib::detail::calculateFluxCorrectedTransportPETSc ( const double t,
const double dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id,
const MathLib::MatrixSpecifications & matrix_specification,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b )

Definition at line 208 of file FluxCorrectedTransport.cpp.

214{
215 BaseLib::MPI::Mpi mpi(PETSC_COMM_WORLD);
216 if (mpi.size > 1)
217 {
218 OGS_FATAL(
219 "The current parallel implementation of flux corrected transport "
220 "supports only one MPI rank, but the simulation uses {} ranks",
221 mpi.size);
222 }
223 Mat K_raw = K.getRawMatrix();
224 Mat D;
225 PetscCallAbort(PETSC_COMM_WORLD,
226 MatDuplicate(K_raw, MAT_DO_NOT_COPY_VALUES, &D));
227 PetscCallAbort(PETSC_COMM_WORLD, MatZeroEntries(D));
228 finalize(D);
229 PetscInt start, end;
230 PetscCallAbort(PETSC_COMM_WORLD,
231 MatGetOwnershipRange(K.getRawMatrix(), &start, &end));
232
233 PetscCallAbort(PETSC_COMM_WORLD, MatScale(K_raw, -1.0));
234 K.finalizeAssembly();
235 // compute artificial diffusion operator D
236 for (PetscInt i = start; i < end; ++i)
237 {
238 PetscInt number_of_columns;
239 PetscInt const* columns;
240 PetscScalar const* values;
241
242 PetscCallAbort(PETSC_COMM_WORLD, MatGetRow(K_raw, i, &number_of_columns,
243 &columns, &values));
244
245 auto const columns_span =
246 std::span<PetscInt const>(columns, number_of_columns);
247 auto const values_span =
248 std::span<PetscScalar const>(values, number_of_columns);
249 for (auto const [j, kij] : std::views::zip(columns_span, values_span))
250 {
251 if (i == j) // skip diagonal entries
252 {
253 continue;
254 }
255 PetscScalar kji;
256 // MatGetValue is not supported for all PETSc matrix types
257 PetscCallAbort(PETSC_COMM_WORLD, MatGetValue(K_raw, j, i, &kji));
258 PetscScalar dij = std::max({-kij, 0.0, -kji});
259 PetscCallAbort(PETSC_COMM_WORLD,
260 MatSetValue(D, i, j, dij, INSERT_VALUES));
261 }
262 // MatRestoreRow() is required after MatGetRow()
263 PetscCallAbort(
264 PETSC_COMM_WORLD,
265 MatRestoreRow(K_raw, i, &number_of_columns, &columns, &values));
266 }
267 finalize(D);
268 auto row_sums = newZeroedInstance<GlobalVector>(matrix_specification);
269 Vec row_sums_raw = row_sums->getRawVector();
270 PetscCallAbort(PETSC_COMM_WORLD, MatGetRowSum(D, row_sums_raw));
271 PetscCallAbort(PETSC_COMM_WORLD, VecScale(row_sums_raw, -1.0));
272 PetscCallAbort(PETSC_COMM_WORLD, VecAssemblyBegin(row_sums_raw));
273 PetscCallAbort(PETSC_COMM_WORLD, VecAssemblyEnd(row_sums_raw));
274 PetscCallAbort(PETSC_COMM_WORLD,
275 MatDiagonalSet(D, row_sums_raw, ADD_VALUES));
276 finalize(D);
277 // need to access x and x_prev in the following loops
278 x[process_id]->setLocalAccessibleVector();
279 x_prev[process_id]->setLocalAccessibleVector();
280
281 // compute F
282 Mat F;
283 PetscCallAbort(PETSC_COMM_WORLD,
284 MatDuplicate(K.getRawMatrix(), MAT_DO_NOT_COPY_VALUES, &F));
285 PetscCallAbort(PETSC_COMM_WORLD, MatZeroEntries(F));
286 for (PetscInt i = start; i < end; ++i)
287 {
288 PetscInt cols;
289 const PetscInt* col_indices;
290 const PetscScalar* values;
291
292 PetscCallAbort(PETSC_COMM_WORLD,
293 MatGetRow(D, i, &cols, &col_indices, &values));
294
295 auto const columns = std::span<PetscInt const>(col_indices, cols);
296 auto const values_D = std::span<PetscScalar const>(values, cols);
297
298 for (auto const [j, dij] : std::views::zip(columns, values_D))
299 {
300 PetscScalar const xi = x[process_id]->get(i);
301 PetscScalar const xj = x[process_id]->get(j);
302 PetscScalar const fij = -dij * (xj - xi);
303 PetscCallAbort(PETSC_COMM_WORLD,
304 MatSetValue(F, i, j, fij, INSERT_VALUES));
305 }
306 PetscCallAbort(PETSC_COMM_WORLD,
307 MatRestoreRow(D, i, &cols, &col_indices, &values));
308 }
309 finalize(F);
310 auto& M_raw = M.getRawMatrix();
311
312 auto xdot = newZeroedInstance<GlobalVector>(matrix_specification);
313 xdot->setLocalAccessibleVector();
314 for (PetscInt i = start; i < end; ++i)
315 {
316 xdot->add(i, (x[process_id]->get(i) - x_prev[process_id]->get(i)) / dt);
317 }
318 xdot->finalizeAssembly();
319 xdot->setLocalAccessibleVector();
320
321 for (PetscInt i = start; i < end; ++i)
322 {
323 PetscInt cols;
324 PetscInt const* col_indices;
325 PetscScalar const* values;
326
327 PetscCallAbort(PETSC_COMM_WORLD,
328 MatGetRow(M_raw, i, &cols, &col_indices, &values));
329
330 auto const columns = std::span<PetscInt const>(col_indices, cols);
331 auto const values_M = std::span<PetscScalar const>(values, cols);
332
333 for (auto const [j, mij] : std::views::zip(columns, values_M))
334 {
335 PetscScalar const fij = -mij * (xdot->get(j) - xdot->get(i));
336 PetscCallAbort(PETSC_COMM_WORLD,
337 MatSetValue(F, i, j, fij, ADD_VALUES));
338 }
339 PetscCallAbort(PETSC_COMM_WORLD,
340 MatRestoreRow(M_raw, i, &cols, &col_indices, &values));
341 }
342 finalize(F);
343 auto P_plus = newZeroedInstance<GlobalVector>(matrix_specification);
344 auto P_minus = newZeroedInstance<GlobalVector>(matrix_specification);
345
346 P_plus->setLocalAccessibleVector();
347 P_minus->setLocalAccessibleVector();
348 // Fill P_plus, P_minus, Q_plus, Q_minus
349 for (PetscInt i = start; i < end; ++i)
350 {
351 PetscInt cols;
352 PetscInt const* col_indices;
353 PetscScalar const* values;
354
355 PetscCallAbort(PETSC_COMM_WORLD,
356 MatGetRow(F, i, &cols, &col_indices, &values));
357
358 auto const column_indices =
359 std::span<PetscInt const>(col_indices, cols);
360 auto const values_span = std::span<PetscScalar const>(values, cols);
361 for (auto const [j, fij] : std::views::zip(column_indices, values_span))
362 {
363 // skip diagonal
364 if (i == j)
365 {
366 continue;
367 }
368 P_plus->add(i, std::max(0., fij));
369 P_minus->add(i, std::min(0., fij));
370 }
371 PetscCallAbort(PETSC_COMM_WORLD,
372 MatRestoreRow(F, i, &cols, &col_indices, &values));
373 }
374 P_plus->finalizeAssembly();
375 P_plus->setLocalAccessibleVector();
376 P_minus->finalizeAssembly();
377 P_minus->setLocalAccessibleVector();
378 auto Q_plus = newZeroedInstance<GlobalVector>(matrix_specification);
379 auto Q_minus = newZeroedInstance<GlobalVector>(matrix_specification);
380 Q_plus->setLocalAccessibleVector();
381 Q_minus->setLocalAccessibleVector();
382
383 for (PetscInt i = start; i < end; ++i)
384 {
385 PetscInt cols;
386 PetscInt const* col_indices;
387
388 PetscScalar Q_plus_i = 0.0;
389 PetscScalar Q_minus_i = 0.0;
390
391 PetscCallAbort(PETSC_COMM_WORLD,
392 MatGetRow(F, i, &cols, &col_indices, nullptr));
393 auto const column_indices =
394 std::span<PetscInt const>(col_indices, cols);
395 for (auto const j : column_indices)
396 {
397 // skip diagonal
398 if (i == j)
399 {
400 continue;
401 }
402 PetscScalar const x_prev_i = x_prev[process_id]->get(i);
403 PetscScalar const x_prev_j = x_prev[process_id]->get(j);
404
405 Q_plus_i = std::max({Q_plus_i, 0., x_prev_j - x_prev_i});
406 Q_minus_i = std::min({Q_minus_i, 0., x_prev_j - x_prev_i});
407 }
408 Q_plus->set(i, Q_plus_i);
409 Q_minus->set(i, Q_minus_i);
410 PetscCallAbort(PETSC_COMM_WORLD,
411 MatRestoreRow(F, i, &cols, &col_indices, nullptr));
412 }
413 Q_plus->finalizeAssembly();
414 Q_plus->setLocalAccessibleVector();
415 Q_minus->finalizeAssembly();
416 Q_minus->setLocalAccessibleVector();
417
418 auto R_plus = newZeroedInstance<GlobalVector>(matrix_specification);
419 R_plus->setLocalAccessibleVector();
420
421 auto row_sums_M = newZeroedInstance<GlobalVector>(matrix_specification);
422 Vec M_L = row_sums_M->getRawVector();
423 PetscCallAbort(PETSC_COMM_WORLD, MatGetRowSum(M_raw, M_L));
424 row_sums_M->finalizeAssembly();
425 row_sums_M->setLocalAccessibleVector();
426
427 for (auto k = R_plus->getRangeBegin(); k < R_plus->getRangeEnd(); ++k)
428 {
429 PetscScalar const P_plus_i = P_plus->get(k);
430 if (P_plus_i == 0.0)
431 {
432 continue;
433 }
434 PetscScalar const mi = row_sums_M->get(k);
435 PetscScalar const Q_plus_i = Q_plus->get(k);
436 R_plus->set(k, std::min(1.0, mi * Q_plus_i / dt / P_plus_i));
437 }
438 R_plus->finalizeAssembly();
439 R_plus->setLocalAccessibleVector();
440 auto R_minus = newZeroedInstance<GlobalVector>(matrix_specification);
441 R_minus->setLocalAccessibleVector();
442
443 for (auto k = R_minus->getRangeBegin(); k < R_minus->getRangeEnd(); ++k)
444 {
445 PetscScalar const P_minus_i = P_minus->get(k);
446 if (P_minus_i == 0.0)
447 {
448 continue;
449 }
450 PetscScalar const mi = row_sums_M->get(k);
451 PetscScalar const Q_minus_i = Q_minus->get(k);
452 R_minus->set(k, std::min(1.0, mi * Q_minus_i / dt / P_minus_i));
453 }
454 R_minus->finalizeAssembly();
455 R_minus->setLocalAccessibleVector();
456 // walk over F, R_plus, and R_minus and compute alpha values; set entries in
457 // the rhs b that limit the antidiffusive flux
458 for (PetscInt i = start; i < end; ++i)
459 {
460 PetscInt cols;
461 PetscInt const* col_indices;
462 PetscScalar const* values;
463
464 PetscCallAbort(PETSC_COMM_WORLD,
465 MatGetRow(F, i, &cols, &col_indices, &values));
466
467 auto col_indices_span = std::span<PetscInt const>(col_indices, cols);
468 auto values_span = std::span<PetscScalar const>(values, cols);
469
470 double alpha_ij;
471
472 for (auto [j, fij] : std::views::zip(col_indices_span, values_span))
473 {
474 if (i == j)
475 {
476 continue;
477 }
478 if (fij > 0.)
479 {
480 alpha_ij = std::min(R_plus->get(i), R_minus->get(j));
481 }
482 else
483 {
484 alpha_ij = std::min(R_minus->get(i), R_plus->get(j));
485 }
486 b.add(i, alpha_ij * fij);
487 }
488 PetscCallAbort(PETSC_COMM_WORLD,
489 MatRestoreRow(F, i, &cols, &col_indices, &values));
490 }
491
492 // compute low-order operator
493 // update K_raw with D: K_raw += D;
494 for (PetscInt i = start; i < end; ++i)
495 {
496 PetscInt cols_D;
497 PetscInt const* col_indices_D;
498 PetscScalar const* values_D;
499 PetscCallAbort(PETSC_COMM_WORLD,
500 MatGetRow(D, i, &cols_D, &col_indices_D, &values_D));
501
502 auto const column_indices =
503 std::span<PetscInt const>(col_indices_D, cols_D);
504 auto const values = std::span<PetscScalar const>(values_D, cols_D);
505 for (auto [j, Dij] : std::views::zip(column_indices, values))
506 {
507 PetscCallAbort(PETSC_COMM_WORLD,
508 MatSetValues(K_raw, 1, &i, 1, &j, &Dij, ADD_VALUES));
509 }
510 PetscCallAbort(PETSC_COMM_WORLD,
511 MatRestoreRow(D, i, &cols_D, &col_indices_D, &values_D));
512 }
513 finalize(K_raw);
514 PetscCallAbort(PETSC_COMM_WORLD, MatScale(K_raw, -1.0));
515 K.finalizeAssembly();
516 M.setZero();
517 for (PetscInt i = start; i < end; ++i)
518 {
519 PetscScalar const row_sum_M_i = row_sums_M->get(i);
520 MatSetValues(M.getRawMatrix(), 1, &i, 1, &i, &row_sum_M_i,
521 INSERT_VALUES);
522 }
523 M.finalizeAssembly();
524 MatDestroy(&D);
525 MatDestroy(&F);
526}
#define OGS_FATAL(...)
Definition Error.h:26
void setZero()
reset data entries to zero.
Definition EigenMatrix.h:65
RawMatrixType & getRawMatrix()
void add(IndexType rowId, double v)
add entry
Definition EigenVector.h:76
auto & get(Tuples &... ts)
Definition Get.h:59
void finalize()

References MathLib::EigenVector::add(), finalize(), MathLib::EigenMatrix::getRawMatrix(), newZeroedInstance(), OGS_FATAL, MathLib::EigenMatrix::setZero(), and BaseLib::MPI::Mpi::size.

Referenced by calculateFluxCorrectedTransport().

◆ checkJacobianDeterminant()

static void NumLib::detail::checkJacobianDeterminant ( const double detJ,
MeshLib::Element const & element )
static

Definition at line 92 of file NaturalCoordinatesMapping.cpp.

94{
95 if (detJ > 0)
96 { // The usual case
97 return;
98 }
99
100 if (detJ < 0)
101 {
102 ERR("det J = {:g} is negative for element {:d}.",
103 detJ,
104 element.getID());
105#ifndef NDEBUG
106 std::cerr << element << "\n";
107#endif // NDEBUG
108 OGS_FATAL(
109 "Please check whether the node numbering of the element is correct,"
110 "or additional elements (like boundary elements) are still present "
111 "in the mesh.");
112 }
113
114 if (detJ == 0)
115 {
116 ERR("det J is zero for element {:d}.", element.getID());
117#ifndef NDEBUG
118 std::cerr << element << "\n";
119#endif // NDEBUG
120 OGS_FATAL(
121 "Please check whether:\n"
122 "\t the element nodes may have the same coordinates,\n"
123 "\t or the coordinates of all nodes are not given on the x-axis "
124 "for a 1D problem or in the xy-plane in the 2D case.\n"
125 "The first case can occur, if not all boundary elements"
126 "were removed from the bulk mesh.");
127 }
128}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45

References ERR(), MeshLib::Element::getID(), and OGS_FATAL.

Referenced by computeMappingMatrices(), and computeMappingMatrices().

◆ computeMappingMatrices() [1/6]

template<class T_SHAPE_FUNC , class T_SHAPE_MATRICES >
void NumLib::detail::computeMappingMatrices ( const MeshLib::Element & ,
const double * natural_pt,
const MeshLib::ElementCoordinatesMappingLocal & ,
T_SHAPE_MATRICES & shapemat,
FieldType< ShapeMatrixType::DNDR >  )
inline

Definition at line 78 of file NaturalCoordinatesMapping.cpp.

84{
85 if constexpr (T_SHAPE_FUNC::DIM != 0)
86 {
87 double* const dNdr = shapemat.dNdr.data();
88 T_SHAPE_FUNC::computeGradShapeFunction(natural_pt, dNdr);
89 }
90}

◆ computeMappingMatrices() [2/6]

template<class T_SHAPE_FUNC , class T_SHAPE_MATRICES >
void NumLib::detail::computeMappingMatrices ( const MeshLib::Element & ,
const double * natural_pt,
const MeshLib::ElementCoordinatesMappingLocal & ,
T_SHAPE_MATRICES & shapemat,
FieldType< ShapeMatrixType::N >  )
inline

Definition at line 67 of file NaturalCoordinatesMapping.cpp.

73{
74 T_SHAPE_FUNC::computeShapeFunction(natural_pt, shapemat.N);
75}

Referenced by computeMappingMatrices(), computeMappingMatrices(), computeMappingMatrices(), computeMappingMatrices(), and naturalCoordinatesMappingComputeShapeMatrices().

◆ computeMappingMatrices() [3/6]

template<class T_SHAPE_FUNC , class T_SHAPE_MATRICES >
void NumLib::detail::computeMappingMatrices ( const MeshLib::Element & ele,
const double * natural_pt,
const MeshLib::ElementCoordinatesMappingLocal & ele_local_coord,
T_SHAPE_MATRICES & shapemat,
FieldType< ShapeMatrixType::ALL >  )
inline

Definition at line 243 of file NaturalCoordinatesMapping.cpp.

249{
250 computeMappingMatrices<T_SHAPE_FUNC, T_SHAPE_MATRICES>(
251 ele,
252 natural_pt,
253 ele_local_coord,
254 shapemat,
256 computeMappingMatrices<T_SHAPE_FUNC, T_SHAPE_MATRICES>(
257 ele,
258 natural_pt,
259 ele_local_coord,
260 shapemat,
262}

References computeMappingMatrices().

◆ computeMappingMatrices() [4/6]

template<class T_SHAPE_FUNC , class T_SHAPE_MATRICES >
void NumLib::detail::computeMappingMatrices ( const MeshLib::Element & ele,
const double * natural_pt,
const MeshLib::ElementCoordinatesMappingLocal & ele_local_coord,
T_SHAPE_MATRICES & shapemat,
FieldType< ShapeMatrixType::DNDR_J >  )
inline

Definition at line 131 of file NaturalCoordinatesMapping.cpp.

137{
138 if constexpr (T_SHAPE_FUNC::DIM != 0)
139 {
140 computeMappingMatrices<T_SHAPE_FUNC, T_SHAPE_MATRICES>(
141 ele,
142 natural_pt,
143 ele_local_coord,
144 shapemat,
146
147 auto const dim = T_SHAPE_FUNC::DIM;
148 auto const nnodes = T_SHAPE_FUNC::NPOINTS;
149
150 // jacobian: J=[dx/dr dy/dr // dx/ds dy/ds]
151 for (auto k = decltype(nnodes){0}; k < nnodes; k++)
152 {
153 const MathLib::Point3d& mapped_pt =
154 ele_local_coord.getMappedCoordinates(k);
155 // outer product of dNdr and mapped_pt for a particular node
156 for (auto i_r = decltype(dim){0}; i_r < dim; i_r++)
157 {
158 for (auto j_x = decltype(dim){0}; j_x < dim; j_x++)
159 {
160 shapemat.J(i_r, j_x) +=
161 shapemat.dNdr(i_r, k) * mapped_pt[j_x];
162 }
163 }
164 }
165
166 shapemat.detJ = shapemat.J.determinant();
167 checkJacobianDeterminant(shapemat.detJ, ele);
168 }
169 else
170 {
171 shapemat.detJ = 1.0;
172 }
173}
MathLib::Point3d const & getMappedCoordinates(std::size_t node_id) const
return mapped coordinates of the node
static void checkJacobianDeterminant(const double detJ, MeshLib::Element const &element)

References checkJacobianDeterminant(), computeMappingMatrices(), and MeshLib::ElementCoordinatesMappingLocal::getMappedCoordinates().

◆ computeMappingMatrices() [5/6]

template<class T_SHAPE_FUNC , class T_SHAPE_MATRICES >
void NumLib::detail::computeMappingMatrices ( const MeshLib::Element & ele,
const double * natural_pt,
const MeshLib::ElementCoordinatesMappingLocal & ele_local_coord,
T_SHAPE_MATRICES & shapemat,
FieldType< ShapeMatrixType::DNDX >  )
inline

Definition at line 198 of file NaturalCoordinatesMapping.cpp.

204{
205 computeMappingMatrices<T_SHAPE_FUNC, T_SHAPE_MATRICES>(
206 ele,
207 natural_pt,
208 ele_local_coord,
209 shapemat,
211 if constexpr (T_SHAPE_FUNC::DIM != 0)
212 {
213 checkJacobianDeterminant(shapemat.detJ, ele);
214
215 // J^-1, dshape/dx
216 shapemat.invJ.noalias() = shapemat.J.inverse();
217
218 assert(shapemat.dNdr.rows() == ele.getDimension());
219 const unsigned global_dim = ele_local_coord.getGlobalDimension();
220 if (global_dim == T_SHAPE_FUNC::DIM)
221 {
222 shapemat.dNdx
223 .template topLeftCorner<T_SHAPE_FUNC::DIM,
224 T_SHAPE_FUNC::NPOINTS>()
225 .noalias() = shapemat.invJ * shapemat.dNdr;
226 }
227 else
228 {
229 auto const& matR =
230 (ele_local_coord.getRotationMatrixToGlobal().topLeftCorner(
231 global_dim, T_SHAPE_FUNC::DIM))
232 .eval();
233
234 auto const invJ_dNdr = shapemat.invJ * shapemat.dNdr;
235 auto const dshape_global = matR * invJ_dNdr;
236 shapemat.dNdx =
237 dshape_global.topLeftCorner(global_dim, T_SHAPE_FUNC::NPOINTS);
238 }
239 }
240}
const RotationMatrix & getRotationMatrixToGlobal() const
return a rotation matrix converting to global coordinates
unsigned getGlobalDimension() const
return the global dimension
virtual constexpr unsigned getDimension() const =0
Get dimension of the mesh element.

References checkJacobianDeterminant(), computeMappingMatrices(), MeshLib::Element::getDimension(), MeshLib::ElementCoordinatesMappingLocal::getGlobalDimension(), and MeshLib::ElementCoordinatesMappingLocal::getRotationMatrixToGlobal().

◆ computeMappingMatrices() [6/6]

template<class T_SHAPE_FUNC , class T_SHAPE_MATRICES >
void NumLib::detail::computeMappingMatrices ( const MeshLib::Element & ele,
const double * natural_pt,
const MeshLib::ElementCoordinatesMappingLocal & ele_local_coord,
T_SHAPE_MATRICES & shapemat,
FieldType< ShapeMatrixType::N_J >  )
inline

Definition at line 176 of file NaturalCoordinatesMapping.cpp.

182{
183 computeMappingMatrices<T_SHAPE_FUNC, T_SHAPE_MATRICES>(
184 ele,
185 natural_pt,
186 ele_local_coord,
187 shapemat,
189 computeMappingMatrices<T_SHAPE_FUNC, T_SHAPE_MATRICES>(
190 ele,
191 natural_pt,
192 ele_local_coord,
193 shapemat,
195}

References computeMappingMatrices().

◆ finalize()

void NumLib::detail::finalize ( Mat & M)

Definition at line 200 of file FluxCorrectedTransport.cpp.

201{
202 PetscCallAbort(PETSC_COMM_WORLD,
203 MatAssemblyBegin(M, MatAssemblyType::MAT_FINAL_ASSEMBLY));
204 PetscCallAbort(PETSC_COMM_WORLD,
205 MatAssemblyEnd(M, MatAssemblyType::MAT_FINAL_ASSEMBLY));
206}

◆ getHydrodynamicDispersion()

Eigen::MatrixXd NumLib::detail::getHydrodynamicDispersion ( Eigen::MatrixXd const & pore_diffusion_coefficient,
Eigen::VectorXd const & velocity,
double const porosity,
double const solute_dispersivity_transverse,
double const solute_dispersivity_longitudinal )
inline

Definition at line 21 of file HydrodynamicDispersion.h.

27{
28 double const velocity_magnitude = velocity.norm();
29 if (velocity_magnitude == 0.0)
30 {
31 return porosity * pore_diffusion_coefficient;
32 }
33
34 auto const dim = velocity.size();
35 Eigen::MatrixXd const& I(Eigen::MatrixXd::Identity(dim, dim));
36 return porosity * pore_diffusion_coefficient +
37 solute_dispersivity_transverse * velocity_magnitude * I +
38 (solute_dispersivity_longitudinal - solute_dispersivity_transverse) /
39 velocity_magnitude * velocity * velocity.transpose();
40}

Referenced by NumLib::computeHydrodynamicDispersion().

◆ getHydrodynamicDispersionWithArtificialDiffusion()

Eigen::MatrixXd NumLib::detail::getHydrodynamicDispersionWithArtificialDiffusion ( IsotropicDiffusionStabilization const & stabilizer,
std::size_t const element_id,
Eigen::MatrixXd const & pore_diffusion_coefficient,
Eigen::VectorXd const & velocity,
double const porosity,
double const solute_dispersivity_transverse,
double const solute_dispersivity_longitudinal )
inline

Definition at line 42 of file HydrodynamicDispersion.h.

50{
51 double const velocity_magnitude = velocity.norm();
52 if (velocity_magnitude == 0.0)
53 {
54 return porosity * pore_diffusion_coefficient;
55 }
56
57 double const artificial_diffusion =
58 stabilizer.computeArtificialDiffusion(element_id, velocity_magnitude);
59
60 auto const dim = velocity.size();
61 Eigen::MatrixXd const& I(Eigen::MatrixXd::Identity(dim, dim));
62 return porosity * pore_diffusion_coefficient +
63 (solute_dispersivity_transverse * velocity_magnitude +
64 artificial_diffusion) *
65 I +
66 (solute_dispersivity_longitudinal - solute_dispersivity_transverse) /
67 velocity_magnitude * velocity * velocity.transpose();
68}

References NumLib::IsotropicDiffusionStabilization::computeArtificialDiffusion().

Referenced by NumLib::computeHydrodynamicDispersion().

◆ getScalarOrVectorComponent()

template<std::size_t IndexInTuple, typename Tuple >
double NumLib::detail::getScalarOrVectorComponent ( Tuple const & tuple,
Eigen::Index component )

Definition at line 39 of file NumericalDifferentiation.h.

40{
41 auto const& value = std::get<IndexInTuple>(tuple);
42
43 if constexpr (IsScalar<std::remove_cvref_t<decltype(value)>>::value)
44 {
45 return value;
46 }
47 else
48 {
49 return value[component];
50 }
51}

Referenced by NumLib::detail::ComputeDerivativeWrtOneScalar_CD::operator()(), and NumLib::detail::ComputeDerivativeWrtOneScalar_FD< Value >::operator()().

◆ localDOFImpl()

template<typename... Ns_t, typename ElementDOFVector , int... Offsets, int... Sizes>
auto NumLib::detail::localDOFImpl ( ElementDOFVector const & x,
boost::mp11::mp_list_c< int, Offsets... > ,
boost::mp11::mp_list_c< int, Sizes... >  )

Definition at line 42 of file LocalDOF.h.

45{
46 static_assert(((Sizes > 0) && ...));
47 static_assert(((Offsets >= 0) && ...));
48
49 return std::tuple<Eigen::Map<typename MatrixPolicyType::VectorType<
50 NumberOfDofs<Ns_t>::value> const>...>{{x.data() + Offsets, Sizes}...};
51}
typename ::detail::EigenMatrixType< N, 1 >::type VectorType

Referenced by NumLib::localDOF().

◆ naturalCoordinatesMappingComputeShapeMatrices()

template<class T_SHAPE_FUNC , class T_SHAPE_MATRICES , ShapeMatrixType T_SHAPE_MATRIX_TYPE>
void NumLib::detail::naturalCoordinatesMappingComputeShapeMatrices ( const MeshLib::Element & ele,
const double * natural_pt,
T_SHAPE_MATRICES & shapemat,
const unsigned global_dim )

Used to explicitly instantiate the NaturalCoordinatesMapping class template.

Definition at line 267 of file NaturalCoordinatesMapping.cpp.

271{
272 const MeshLib::ElementCoordinatesMappingLocal ele_local_coord(ele,
273 global_dim);
274
275 detail::computeMappingMatrices<T_SHAPE_FUNC, T_SHAPE_MATRICES>(
276 ele,
277 natural_pt,
278 ele_local_coord,
279 shapemat,
281}

References computeMappingMatrices().

Referenced by NumLib::NaturalCoordinatesMapping< T_SHAPE_FUNC, T_SHAPE_MATRICES >::computeShapeMatrices().

◆ newZeroedInstance()

template<typename MatrixVectorType >
std::unique_ptr< MatrixVectorType > NumLib::detail::newZeroedInstance ( MathLib::MatrixSpecifications const & matrix_specification)

Definition at line 30 of file FluxCorrectedTransport.cpp.

32{
34 matrix_specification);
35 result->setZero();
36 return result;
37}

Referenced by calculateFluxCorrectedTransportPETSc().

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN() [1/16]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN ( ShapeHex20 )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN() [2/16]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN ( ShapeHex8 )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN() [3/16]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN ( ShapeLine2 )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN() [4/16]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN ( ShapeLine3 )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN() [5/16]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN ( ShapePoint1 )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN() [6/16]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN ( ShapePrism15 )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN() [7/16]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN ( ShapePrism6 )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN() [8/16]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN ( ShapePyra13 )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN() [9/16]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN ( ShapePyra5 )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN() [10/16]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN ( ShapeQuad4 )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN() [11/16]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN ( ShapeQuad8 )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN() [12/16]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN ( ShapeQuad9 )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN() [13/16]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN ( ShapeTet10 )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN() [14/16]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN ( ShapeTet4 )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN() [15/16]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN ( ShapeTri3 )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN() [16/16]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_DYN ( ShapeTri6 )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [1/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeHex20 ,
3  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [2/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeHex8 ,
3  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [3/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeLine2 ,
1  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [4/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeLine2 ,
2  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [5/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeLine2 ,
3  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [6/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeLine3 ,
1  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [7/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeLine3 ,
2  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [8/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeLine3 ,
3  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [9/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapePoint1 ,
0  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [10/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapePoint1 ,
1  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [11/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapePoint1 ,
2  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [12/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapePoint1 ,
3  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [13/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapePrism15 ,
3  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [14/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapePrism6 ,
3  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [15/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapePyra13 ,
3  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [16/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapePyra5 ,
3  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [17/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeQuad4 ,
2  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [18/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeQuad4 ,
3  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [19/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeQuad8 ,
2  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [20/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeQuad8 ,
3  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [21/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeQuad9 ,
2  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [22/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeQuad9 ,
3  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [23/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeTet10 ,
3  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [24/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeTet4 ,
3  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [25/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeTri3 ,
2  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [26/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeTri3 ,
3  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [27/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeTri6 ,
2  )

◆ OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX() [28/28]

NumLib::detail::OGS_INSTANTIATE_NATURAL_COORDINATES_MAPPING_FIX ( ShapeTri6 ,
3  )

◆ setMatrixZero()

template<class T >
void NumLib::detail::setMatrixZero ( T & mat)

Definition at line 21 of file ShapeMatrices-impl.h.

22{
23 // mat.setZero();
24 const std::size_t n = mat.rows() * mat.cols();
25 auto* v = mat.data();
26 for (std::size_t i = 0; i < n; i++)
27 {
28 v[i] = .0;
29 }
30}

Referenced by setZero(), setZero(), and setZero().

◆ setVectorZero()

template<class T >
void NumLib::detail::setVectorZero ( T & vec)

Definition at line 33 of file ShapeMatrices-impl.h.

34{
35 // vec.setZero();
36 const std::size_t n = vec.size();
37 auto* v = vec.data();
38 for (std::size_t i = 0; i < n; i++)
39 {
40 v[i] = .0;
41 }
42}

Referenced by setZero().

◆ setZero() [1/6]

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
void NumLib::detail::setZero ( ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX > & shape,
ShapeDataFieldType< ShapeMatrixType::ALL >  )
inline

Definition at line 95 of file ShapeMatrices-impl.h.

97{
100}
void setZero(ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX > &shape, ShapeDataFieldType< ShapeMatrixType::N >)

References setZero().

◆ setZero() [2/6]

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
void NumLib::detail::setZero ( ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX > & shape,
ShapeDataFieldType< ShapeMatrixType::DNDR >  )
inline

◆ setZero() [3/6]

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
void NumLib::detail::setZero ( ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX > & shape,
ShapeDataFieldType< ShapeMatrixType::DNDR_J >  )
inline

◆ setZero() [4/6]

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
void NumLib::detail::setZero ( ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX > & shape,
ShapeDataFieldType< ShapeMatrixType::DNDX >  )
inline

◆ setZero() [5/6]

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
void NumLib::detail::setZero ( ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX > & shape,
ShapeDataFieldType< ShapeMatrixType::N >  )
inline

◆ setZero() [6/6]

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
void NumLib::detail::setZero ( ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX > & shape,
ShapeDataFieldType< ShapeMatrixType::N_J >  )
inline

Definition at line 78 of file ShapeMatrices-impl.h.

References setZero().

◆ shapeFunctionInterpolate() [1/2]

template<unsigned DOFOffset, typename NodalValues , typename ShapeMatrix >
void NumLib::detail::shapeFunctionInterpolate ( const NodalValues & ,
const ShapeMatrix &  )
See also
NumLib::shapeFunctionInterpolate()

Definition at line 27 of file Interpolation.h.

29{
30}

Referenced by ProcessLib::HCNonAdvectiveFreeComponentFlowBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::assemble(), ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::assemble(), ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::assemble(), ProcessLib::HT::MonolithicHTFEM< ShapeFunction, GlobalDim >::assemble(), ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::assemble(), ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::assemble(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assemble(), ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::assemble(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::assemble(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowLocalAssembler< ShapeFunction, GlobalDim >::assemble(), ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPLocalAssembler< ShapeFunction, GlobalDim >::assemble(), ProcessLib::VariableDependentNeumannBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::assemble(), ProcessLib::WellboreCompensateNeumannBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::assemble(), ProcessLib::WellboreSimulator::WellboreSimulatorFEM< ShapeFunction, GlobalDim >::assemble(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::assembleBlockMatrices(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::assembleHeatTransportEquation(), ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::assembleHeatTransportEquation(), ProcessLib::HT::StaggeredHTFEM< ShapeFunction, GlobalDim >::assembleHydraulicEquation(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::assembleReactionEquationConcrete(), ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::assembleWithJacobian(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobian(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowLocalAssembler< ShapeFunction, GlobalDim >::assembleWithJacobian(), ProcessLib::HydroMechanics::HydroMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobianForDeformationEquations(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::calculateIntPtDarcyVelocity(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowLocalAssembler< ShapeFunction, GlobalDim >::computeSecondaryVariableConcrete(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::getFlux(), ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::getFlux(), ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::getFlux(), ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::getFlux(), ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::getIntPtDarcyVelocity(), ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, GlobalDim >::getIntPtDarcyVelocity(), ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::getIntPtDarcyVelocity(), ProcessLib::HT::HTFEM< ShapeFunction, GlobalDim >::getIntPtDarcyVelocityLocal(), ProcessLib::RichardsComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::getIntPtSaturation(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::initializeChemicalSystemConcrete(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::setChemicalSystemConcrete(), ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::setInitialConditionsConcrete(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowLocalAssembler< ShapeFunction, GlobalDim >::setInitialConditionsConcrete(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunction, DisplacementDim, ConstitutiveTraits >::setInitialConditionsConcrete(), shapeFunctionInterpolate(), and NumLib::shapeFunctionInterpolate().

◆ shapeFunctionInterpolate() [2/2]

template<unsigned DOFOffset, typename NodalValues , typename ShapeMatrix , typename... ScalarTypes>
void NumLib::detail::shapeFunctionInterpolate ( const NodalValues & nodal_values,
const ShapeMatrix & shape_matrix_N,
double & interpolated_value,
ScalarTypes &... interpolated_values )
See also
NumLib::shapeFunctionInterpolate()

Definition at line 35 of file Interpolation.h.

39{
40 auto const num_nodes = shape_matrix_N.size();
41 double iv = 0.0;
42
43 for (auto n = decltype(num_nodes){0}; n < num_nodes; ++n)
44 {
45 iv += nodal_values[DOFOffset * num_nodes + n] * shape_matrix_N[n];
46 }
47
48 interpolated_value = iv;
49
50 shapeFunctionInterpolate<DOFOffset + 1>(nodal_values, shape_matrix_N,
51 interpolated_values...);
52}

References shapeFunctionInterpolate().

◆ solvePicard()

bool NumLib::detail::solvePicard ( GlobalLinearSolver & linear_solver,
GlobalMatrix & A,
GlobalVector & rhs,
GlobalVector & x,
MathLib::LinearSolverBehaviour const linear_solver_behaviour )

Definition at line 58 of file NonlinearSolver.cpp.

61{
62 if (linear_solver_behaviour ==
64 linear_solver_behaviour == MathLib::LinearSolverBehaviour::REUSE)
65 {
66 WARN(
67 "The performance optimization to skip the linear solver compute() "
68 "step is not implemented for PETSc or LIS linear solvers.");
69 }
70
71 BaseLib::RunTime time_linear_solver;
72 time_linear_solver.start();
73
74 bool const iteration_succeeded = linear_solver.solve(A, rhs, x);
75
76 INFO("[time] Linear solver took {:g} s.", time_linear_solver.elapsed());
77
78 if (iteration_succeeded)
79 {
80 return true;
81 }
82
83 ERR("Picard: The linear solver failed in the solve() step.");
84 return false;
85}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
Count the running time.
Definition RunTime.h:29
double elapsed() const
Get the elapsed time in seconds.
Definition RunTime.h:42
void start()
Start the timer.
Definition RunTime.h:32
bool solve(EigenMatrix &A, EigenVector &b, EigenVector &x)

References BaseLib::RunTime::elapsed(), ERR(), INFO(), MathLib::RECOMPUTE_AND_STORE, MathLib::REUSE, MathLib::EigenLisLinearSolver::solve(), BaseLib::RunTime::start(), and WARN().

Referenced by NumLib::NonlinearSolver< NonlinearSolverTag::Picard >::solve().