1 #ifndef STAN_MATH_FWD_MAT_FUN_MDIVIDE_LEFT_TRI_LOW_HPP 2 #define STAN_MATH_FWD_MAT_FUN_MDIVIDE_LEFT_TRI_LOW_HPP 16 template <
typename T,
int R1,
int C1,
int R2,
int C2>
18 const Eigen::Matrix<
fvar<T>, R1, C1>& A,
19 const Eigen::Matrix<
fvar<T>, R2, C2>& b) {
23 Eigen::Matrix<T, R1, C2> inv_A_mult_b(A.rows(), b.cols());
24 Eigen::Matrix<T, R1, C2> inv_A_mult_deriv_b(A.rows(), b.cols());
25 Eigen::Matrix<T, R1, C1> inv_A_mult_deriv_A(A.rows(), A.cols());
26 Eigen::Matrix<T, R1, C1> val_A(A.rows(), A.cols());
27 Eigen::Matrix<T, R1, C1> deriv_A(A.rows(), A.cols());
28 Eigen::Matrix<T, R2, C2> val_b(b.rows(), b.cols());
29 Eigen::Matrix<T, R2, C2> deriv_b(b.rows(), b.cols());
33 for (
size_type j = 0; j < A.cols(); j++) {
34 for (
size_type i = j; i < A.rows(); i++) {
35 val_A(i, j) = A(i, j).val_;
36 deriv_A(i, j) = A(i, j).d_;
40 for (
size_type j = 0; j < b.cols(); j++) {
41 for (
size_type i = 0; i < b.rows(); i++) {
42 val_b(i, j) = b(i, j).val_;
43 deriv_b(i, j) = b(i, j).d_;
51 Eigen::Matrix<T, R1, C2> deriv(A.rows(), b.cols());
52 deriv = inv_A_mult_deriv_b -
multiply(inv_A_mult_deriv_A, inv_A_mult_b);
54 return to_fvar(inv_A_mult_b, deriv);
57 template <
typename T,
int R1,
int C1,
int R2,
int C2>
59 const Eigen::Matrix<double, R1, C1>& A,
60 const Eigen::Matrix<
fvar<T>, R2, C2>& b) {
64 Eigen::Matrix<T, R1, C2> inv_A_mult_b(A.rows(), b.cols());
65 Eigen::Matrix<T, R1, C2> inv_A_mult_deriv_b(A.rows(), b.cols());
66 Eigen::Matrix<T, R2, C2> val_b(b.rows(), b.cols());
67 Eigen::Matrix<T, R2, C2> deriv_b(b.rows(), b.cols());
68 Eigen::Matrix<double, R1, C1> val_A(A.rows(), A.cols());
71 for (
size_type j = 0; j < A.cols(); j++) {
72 for (
size_type i = j; i < A.rows(); i++) {
73 val_A(i, j) = A(i, j);
77 for (
size_type j = 0; j < b.cols(); j++) {
78 for (
size_type i = 0; i < b.rows(); i++) {
79 val_b(i, j) = b(i, j).val_;
80 deriv_b(i, j) = b(i, j).d_;
87 Eigen::Matrix<T, R1, C2> deriv(A.rows(), b.cols());
88 deriv = inv_A_mult_deriv_b;
90 return to_fvar(inv_A_mult_b, deriv);
93 template <
typename T,
int R1,
int C1,
int R2,
int C2>
95 const Eigen::Matrix<
fvar<T>, R1, C1>& A,
96 const Eigen::Matrix<double, R2, C2>& b) {
100 Eigen::Matrix<T, R1, C2> inv_A_mult_b(A.rows(), b.cols());
101 Eigen::Matrix<T, R1, C1> inv_A_mult_deriv_A(A.rows(), A.cols());
102 Eigen::Matrix<T, R1, C1> val_A(A.rows(), A.cols());
103 Eigen::Matrix<T, R1, C1> deriv_A(A.rows(), A.cols());
107 for (
size_type j = 0; j < A.cols(); j++) {
108 for (
size_type i = j; i < A.rows(); i++) {
109 val_A(i, j) = A(i, j).val_;
110 deriv_A(i, j) = A(i, j).d_;
117 Eigen::Matrix<T, R1, C2> deriv(A.rows(), b.cols());
118 deriv = -
multiply(inv_A_mult_deriv_A, inv_A_mult_b);
120 return to_fvar(inv_A_mult_b, deriv);
Eigen::Matrix< fvar< T >, R1, C2 > mdivide_left(const Eigen::Matrix< fvar< T >, R1, C1 > &A, const Eigen::Matrix< fvar< T >, R2, C2 > &b)
void check_square(const char *function, const char *name, const matrix_cl &y)
Check if the matrix_cl is square.
Eigen::Matrix< fvar< T >, R1, C1 > multiply(const Eigen::Matrix< fvar< T >, R1, C1 > &m, const fvar< T > &c)
std::vector< fvar< T > > to_fvar(const std::vector< T > &v)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic >::Index size_type
Type for sizes and indexes in an Eigen matrix with double e.
void check_multiplicable(const char *function, const char *name1, const T1 &y1, const char *name2, const T2 &y2)
Check if the matrices can be multiplied.
Eigen::Matrix< fvar< T >, R1, C1 > mdivide_left_tri_low(const Eigen::Matrix< fvar< T >, R1, C1 > &A, const Eigen::Matrix< fvar< T >, R2, C2 > &b)
This template class represents scalars used in forward-mode automatic differentiation, which consist of values and directional derivatives of the specified template type.