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>,
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,
GlobalMatrix::IndexType GlobalIndexType

Definition at line 95 of file ComponentGlobalIndexDict.h.

◆ GetMeshElement

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

Definition at line 16 of file ShapeMatrixCache.h.

◆ GetShapeFunctionHigherOrder

Initial value:
typename ElementTraitsLagrange::ShapeFunction

Definition at line 19 of file ShapeMatrixCache.h.

◆ GetShapeFunctionLowerOrder

Initial value:
typename ElementTraitsLagrange::LowerOrderShapeFunction

Definition at line 23 of file ShapeMatrixCache.h.

◆ GetShapeMatrix_N

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

Definition at line 34 of file ShapeMatrixCache.h.

◆ GetShapeMatrixPolicy

Initial value:
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType

Definition at line 27 of file ShapeMatrixCache.h.

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 62 of file AdvectionMatrixAssembler.h.

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

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 84 of file AdvectionMatrixAssembler.h.

87{
88 Eigen::VectorXd quasi_nodal_flux =
89 Eigen::VectorXd::Zero(laplacian_matrix.rows());
90
91 for (std::size_t ip = 0; ip < ip_flux_vector.size(); ++ip)
92 {
93 auto const& ip_data = ip_data_vector[ip];
94 auto const w = ip_data.integration_weight;
95 auto const& dNdx = ip_data.dNdx;
96
97 quasi_nodal_flux.noalias() -= ip_flux_vector[ip].transpose() * dNdx * w;
98 }
99
100 applyFullUpwind(quasi_nodal_flux, laplacian_matrix);
101}
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 22 of file AdvectionMatrixAssembler.h.

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

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 46 of file AdvectionMatrixAssembler.h.

49{
50 for (std::size_t ip = 0; ip < ip_flux_vector.size(); ++ip)
51 {
52 auto const& ip_data = ip_data_vector[ip];
53 auto const w = ip_data.integration_weight;
54 auto const& N = ip_data.N;
55 auto const& dNdx = ip_data.dNdx;
56 laplacian_matrix.noalias() +=
57 N.transpose() * ip_flux_vector[ip].transpose() * dNdx * w;
58 }
59}

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 523 of file FluxCorrectedTransport.cpp.

528{
529#ifndef USE_PETSC
530 calculateFluxCorrectedTransportSerial(t, dt, x, x_prev, process_id,
531 matrix_specification, M, K, b);
532#endif // end of ifndef USE_PETSC
533
534#ifdef USE_PETSC
535 calculateFluxCorrectedTransportPETSc(t, dt, x, x_prev, process_id,
536 matrix_specification, M, K, b);
537#endif // end of ifdef USE_PETSC
538}
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 202 of file FluxCorrectedTransport.cpp.

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

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

Referenced by calculateFluxCorrectedTransport().

◆ checkJacobianDeterminant()

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

Definition at line 85 of file NaturalCoordinatesMapping.cpp.

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

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 71 of file NaturalCoordinatesMapping.cpp.

77{
78 if constexpr (T_SHAPE_FUNC::DIM != 0)
79 {
80 double* const dNdr = shapemat.dNdr.data();
81 T_SHAPE_FUNC::computeGradShapeFunction(natural_pt, dNdr);
82 }
83}

◆ 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 60 of file NaturalCoordinatesMapping.cpp.

66{
67 T_SHAPE_FUNC::computeShapeFunction(natural_pt, shapemat.N);
68}

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 229 of file NaturalCoordinatesMapping.cpp.

235{
237 ele,
238 natural_pt,
239 ele_local_coord,
240 shapemat,
243 ele,
244 natural_pt,
245 ele_local_coord,
246 shapemat,
248}
void computeMappingMatrices(const MeshLib::Element &, const double *natural_pt, const MeshLib::ElementCoordinatesMappingLocal &, T_SHAPE_MATRICES &shapemat, FieldType< ShapeMatrixType::N >)

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 124 of file NaturalCoordinatesMapping.cpp.

130{
131 if constexpr (T_SHAPE_FUNC::DIM != 0)
132 {
134 ele,
135 natural_pt,
136 ele_local_coord,
137 shapemat,
139
140 auto constexpr nnodes = T_SHAPE_FUNC::NPOINTS;
141
142 // jacobian: J=[dx/dr dy/dr // dx/ds dy/ds]
143 for (auto k = decltype(nnodes){0}; k < nnodes; k++)
144 {
145 const MathLib::Point3d& mapped_pt =
146 ele_local_coord.getMappedCoordinates(k);
147 shapemat.J +=
148 shapemat.dNdr.col(k) *
149 mapped_pt.asEigenVector3d().head(T_SHAPE_FUNC::DIM).transpose();
150 }
151
152 shapemat.detJ = shapemat.J.determinant();
153 checkJacobianDeterminant(shapemat.detJ, ele);
154 }
155 else
156 {
157 shapemat.detJ = 1.0;
158 }
159}
Eigen::Vector3d const & asEigenVector3d() const
Definition Point3d.h:55
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 MathLib::Point3d::asEigenVector3d(), 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 184 of file NaturalCoordinatesMapping.cpp.

190{
192 ele,
193 natural_pt,
194 ele_local_coord,
195 shapemat,
197 if constexpr (T_SHAPE_FUNC::DIM != 0)
198 {
199 checkJacobianDeterminant(shapemat.detJ, ele);
200
201 // J^-1, dshape/dx
202 shapemat.invJ.noalias() = shapemat.J.inverse();
203
204 assert(shapemat.dNdr.rows() == ele.getDimension());
205 const unsigned global_dim = ele_local_coord.getGlobalDimension();
206 if (global_dim == T_SHAPE_FUNC::DIM)
207 {
208 shapemat.dNdx
209 .template topLeftCorner<T_SHAPE_FUNC::DIM,
210 T_SHAPE_FUNC::NPOINTS>()
211 .noalias() = shapemat.invJ * shapemat.dNdr;
212 }
213 else
214 {
215 auto const& matR =
216 (ele_local_coord.getRotationMatrixToGlobal().topLeftCorner(
217 global_dim, T_SHAPE_FUNC::DIM))
218 .eval();
219
220 auto const invJ_dNdr = shapemat.invJ * shapemat.dNdr;
221 auto const dshape_global = matR * invJ_dNdr;
222 shapemat.dNdx =
223 dshape_global.topLeftCorner(global_dim, T_SHAPE_FUNC::NPOINTS);
224 }
225 }
226}
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 162 of file NaturalCoordinatesMapping.cpp.

168{
170 ele,
171 natural_pt,
172 ele_local_coord,
173 shapemat,
176 ele,
177 natural_pt,
178 ele_local_coord,
179 shapemat,
181}

References computeMappingMatrices().

◆ finalize()

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

Definition at line 194 of file FluxCorrectedTransport.cpp.

195{
196 PetscCallAbort(PETSC_COMM_WORLD,
197 MatAssemblyBegin(M, MatAssemblyType::MAT_FINAL_ASSEMBLY));
198 PetscCallAbort(PETSC_COMM_WORLD,
199 MatAssemblyEnd(M, MatAssemblyType::MAT_FINAL_ASSEMBLY));
200}

Referenced by calculateFluxCorrectedTransportPETSc().

◆ 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 15 of file HydrodynamicDispersion.h.

21{
22 double const velocity_magnitude = velocity.norm();
23 if (velocity_magnitude == 0.0)
24 {
25 return porosity * pore_diffusion_coefficient;
26 }
27
28 auto const dim = velocity.size();
29 Eigen::MatrixXd const& I(Eigen::MatrixXd::Identity(dim, dim));
30 return porosity * pore_diffusion_coefficient +
31 solute_dispersivity_transverse * velocity_magnitude * I +
32 (solute_dispersivity_longitudinal - solute_dispersivity_transverse) /
33 velocity_magnitude * velocity * velocity.transpose();
34}

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 36 of file HydrodynamicDispersion.h.

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

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 32 of file NumericalDifferentiation.h.

33{
34 auto const& value = std::get<IndexInTuple>(tuple);
35
36 if constexpr (IsScalar<std::remove_cvref_t<decltype(value)>>::value)
37 {
38 return value;
39 }
40 else
41 {
42 return value[component];
43 }
44}

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 34 of file LocalDOF.h.

37{
38 static_assert(((Sizes > 0) && ...));
39 static_assert(((Offsets >= 0) && ...));
40
41 return std::tuple<Eigen::Map<typename MatrixPolicyType::VectorType<
42 NumberOfDofs<Ns_t>::value> const>...>{{x.data() + Offsets, Sizes}...};
43}
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 253 of file NaturalCoordinatesMapping.cpp.

257{
258 const MeshLib::ElementCoordinatesMappingLocal ele_local_coord(ele,
259 global_dim);
260
262 ele,
263 natural_pt,
264 ele_local_coord,
265 shapemat,
267}

References computeMappingMatrices().

Referenced by NumLib::NaturalCoordinatesMapping< ShapeFunctionType, ShapeMatrices >::computeShapeMatrices().

◆ newZeroedInstance()

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

Definition at line 24 of file FluxCorrectedTransport.cpp.

26{
28 matrix_specification);
29 result->setZero();
30 return result;
31}

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 12 of file ShapeMatrices-impl.h.

13{
14 // mat.setZero();
15 const std::size_t n = mat.rows() * mat.cols();
16 auto* v = mat.data();
17 for (std::size_t i = 0; i < n; i++)
18 {
19 v[i] = .0;
20 }
21}

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

◆ setVectorZero()

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

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

25{
26 // vec.setZero();
27 const std::size_t n = vec.size();
28 auto* v = vec.data();
29 for (std::size_t i = 0; i < n; i++)
30 {
31 v[i] = .0;
32 }
33}

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 86 of file ShapeMatrices-impl.h.

88{
91}
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 69 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 20 of file Interpolation.h.

22{
23}

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::ComponentTransport::LocalAssemblerData< ShapeFunction, GlobalDim >::calculateIntPtLiquidDensity(), 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::EmbeddedAnchor< GlobalDim >::getShapeMatricesAndGlobalIndicesAndDisplacements(), 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 28 of file Interpolation.h.

32{
33 auto const num_nodes = shape_matrix_N.size();
34 double iv = 0.0;
35
36 for (auto n = decltype(num_nodes){0}; n < num_nodes; ++n)
37 {
38 iv += nodal_values[DOFOffset * num_nodes + n] * shape_matrix_N[n];
39 }
40
41 interpolated_value = iv;
42
43 shapeFunctionInterpolate<DOFOffset + 1>(nodal_values, shape_matrix_N,
44 interpolated_values...);
45}
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)

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 51 of file NonlinearSolver.cpp.

54{
55 if (linear_solver_behaviour ==
57 linear_solver_behaviour == MathLib::LinearSolverBehaviour::REUSE)
58 {
59 WARN(
60 "The performance optimization to skip the linear solver compute() "
61 "step is not implemented for PETSc or LIS linear solvers.");
62 }
63
64 BaseLib::RunTime time_linear_solver;
65 time_linear_solver.start();
66
67 bool const iteration_succeeded = linear_solver.solve(A, rhs, x);
68
69 INFO("[time] Linear solver took {:g} s.", time_linear_solver.elapsed());
70
71 if (iteration_succeeded)
72 {
73 return true;
74 }
75
76 ERR("Picard: The linear solver failed in the solve() step.");
77 return false;
78}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:34
Count the running time.
Definition RunTime.h:18
double elapsed() const
Get the elapsed time in seconds.
Definition RunTime.h:31
void start()
Start the timer.
Definition RunTime.h:21
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().