66 o2cblas_ConjTrans=113};
84 template<
class vec_t>
double dasum(
const size_t N,
const vec_t &X) {
86 for(
size_t i=0;i<N;i++) {
97 template<
class vec_t,
class vec2_t>
98 void daxpy(
const double alpha,
const size_t N,
const vec_t &X,
107 const size_t m=N % 4;
110 O2SCL_IX(Y,i)+=alpha*O2SCL_IX(X,i);
113 for (i=m;i+3<N;i+=4) {
114 O2SCL_IX(Y,i)+=alpha*O2SCL_IX(X,i);
115 O2SCL_IX(Y,i+1)+=alpha*O2SCL_IX(X,i+1);
116 O2SCL_IX(Y,i+2)+=alpha*O2SCL_IX(X,i+2);
117 O2SCL_IX(Y,i+3)+=alpha*O2SCL_IX(X,i+3);
122 template<
class vec_t,
class vec2_t>
123 double ddot(
const size_t N,
const vec_t &X,
const vec2_t &Y) {
128 const size_t m=N % 4;
131 r+=O2SCL_IX(X,i)*O2SCL_IX(Y,i);
134 for (i=m;i+3<N;i+=4) {
135 r+=O2SCL_IX(X,i)*O2SCL_IX(Y,i);
136 r+=O2SCL_IX(X,i+1)*O2SCL_IX(Y,i+1);
137 r+=O2SCL_IX(X,i+2)*O2SCL_IX(Y,i+2);
138 r+=O2SCL_IX(X,i+3)*O2SCL_IX(Y,i+3);
156 template<
class vec_t>
double dnrm2(
const size_t N,
const vec_t &X) {
165 return fabs(O2SCL_IX(X,0));
169 const double x=O2SCL_IX(X,i);
172 const double ax=fabs(x);
175 ssq=1.0+ssq*(scale/ax)*(scale/ax);
178 ssq+=(ax/scale)*(ax/scale);
184 return scale*sqrt(ssq);
189 template<
class vec_t>
190 void dscal(
const double alpha,
const size_t N, vec_t &X) {
193 const size_t m=N % 4;
196 O2SCL_IX(X,i)*=alpha;
199 for (i=m;i+3<N;i+=4) {
200 O2SCL_IX(X,i)*=alpha;
201 O2SCL_IX(X,i+1)*=alpha;
202 O2SCL_IX(X,i+2)*=alpha;
203 O2SCL_IX(X,i+3)*=alpha;
217 template<
class mat_t,
class vec_t,
class vec2_t>
220 const size_t N,
const double alpha,
const mat_t &A,
221 const vec_t &X,
const double beta, vec2_t &Y) {
227 const int Trans=(TransA != o2cblas_ConjTrans) ? TransA : o2cblas_Trans;
229 if (M == 0 || N == 0) {
233 if (alpha == 0.0 && beta == 1.0) {
237 if (Trans == o2cblas_NoTrans) {
248 for (i=0;i<lenY;i++) {
252 }
else if (beta != 1.0) {
254 for (i=0;i<lenY;i++) {
255 O2SCL_IX(Y,iy) *= beta;
264 if ((order == o2cblas_RowMajor && Trans == o2cblas_NoTrans) ||
265 (order == o2cblas_ColMajor && Trans == o2cblas_Trans)) {
269 for (i=0;i<lenY;i++) {
272 for (j=0;j<lenX;j++) {
273 temp+=O2SCL_IX(X,ix)*O2SCL_IX2(A,i,j);
276 O2SCL_IX(Y,iy)+=alpha*temp;
280 }
else if ((order == o2cblas_RowMajor && Trans == o2cblas_Trans) ||
281 (order == o2cblas_ColMajor && Trans == o2cblas_NoTrans)) {
285 for (j=0;j<lenX;j++) {
286 const double temp=alpha*O2SCL_IX(X,ix);
289 for (i=0;i<lenY;i++) {
290 O2SCL_IX(Y,iy)+=temp*O2SCL_IX2(A,j,i);
307 template<
class mat_t,
class vec_t>
312 const size_t M,
const size_t N,
const mat_t &A, vec_t &X) {
314 const int nonunit=(Diag == o2cblas_NonUnit);
317 const int Trans=(TransA != o2cblas_ConjTrans) ? TransA : o2cblas_Trans;
325 if ((order == o2cblas_RowMajor && Trans == o2cblas_NoTrans &&
326 Uplo == o2cblas_Upper) ||
327 (order == o2cblas_ColMajor && Trans == o2cblas_Trans &&
328 Uplo == o2cblas_Lower)) {
336 O2SCL_IX(X,ix)=O2SCL_IX(X,ix)/O2SCL_IX2(A,N-1,N-1);
340 for (i=(
int)(N-1);i>0 && i--;) {
341 double tmp=O2SCL_IX(X,ix);
343 for (j=i+1;j<((int)N);j++) {
344 const double Aij=O2SCL_IX2(A,i,j);
345 tmp-=Aij*O2SCL_IX(X,jx);
349 O2SCL_IX(X,ix)=tmp/O2SCL_IX2(A,i,i);
356 }
else if ((order == o2cblas_RowMajor && Trans == o2cblas_NoTrans &&
357 Uplo == o2cblas_Lower) ||
358 (order == o2cblas_ColMajor && Trans == o2cblas_Trans &&
359 Uplo == o2cblas_Upper)) {
364 O2SCL_IX(X,ix)=O2SCL_IX(X,ix)/O2SCL_IX2(A,0,0);
367 for (i=1;i<((int)N);i++) {
368 double tmp=O2SCL_IX(X,ix);
371 const double Aij=O2SCL_IX2(A,i,j);
372 tmp-=Aij*O2SCL_IX(X,jx);
376 O2SCL_IX(X,ix)=tmp/O2SCL_IX2(A,i,i);
383 }
else if ((order == o2cblas_RowMajor && Trans == o2cblas_Trans &&
384 Uplo == o2cblas_Upper) ||
385 (order == o2cblas_ColMajor && Trans == o2cblas_NoTrans &&
386 Uplo == o2cblas_Lower)) {
393 O2SCL_IX(X,ix)=O2SCL_IX(X,ix)/O2SCL_IX2(A,0,0);
396 for (i=1;i<((int)N);i++) {
397 double tmp=O2SCL_IX(X,ix);
400 const double Aji=O2SCL_IX2(A,j,i);
401 tmp-=Aji*O2SCL_IX(X,jx);
405 O2SCL_IX(X,ix)=tmp/O2SCL_IX2(A,i,i);
412 }
else if ((order == o2cblas_RowMajor && Trans == o2cblas_Trans &&
413 Uplo == o2cblas_Lower) ||
414 (order == o2cblas_ColMajor && Trans == o2cblas_NoTrans &&
415 Uplo == o2cblas_Upper)) {
422 O2SCL_IX(X,ix)=O2SCL_IX(X,ix)/O2SCL_IX2(A,N-1,N-1);
425 for (i=(
int)(N-1);i>0 && i--;) {
426 double tmp=O2SCL_IX(X,ix);
428 for (j=i+1;j<((int)N);j++) {
429 const double Aji=O2SCL_IX2(A,j,i);
430 tmp-=Aji*O2SCL_IX(X,jx);
434 O2SCL_IX(X,ix)=tmp/O2SCL_IX2(A,i,i);
449 template<
class mat_t,
class vec_t>
454 const mat_t &A, vec_t &x) {
458 const int nonunit=(Diag == o2cblas_NonUnit);
459 const int Trans=(TransA != o2cblas_ConjTrans) ? TransA : o2cblas_Trans;
461 if ((Order == o2cblas_RowMajor && Trans == o2cblas_NoTrans &&
462 Uplo == o2cblas_Upper) ||
463 (Order == o2cblas_ColMajor && Trans == o2cblas_Trans &&
464 Uplo == o2cblas_Lower)) {
470 const size_t j_min=i+1;
471 const size_t j_max=N;
473 for (j=j_min;j<j_max;j++) {
474 temp+=O2SCL_IX(x,jx)*O2SCL_IX2(A,i,j);
478 O2SCL_IX(x,i)=temp+O2SCL_IX(x,i)*O2SCL_IX2(A,i,i);
484 }
else if ((Order == o2cblas_RowMajor && Trans == o2cblas_NoTrans &&
485 Uplo == o2cblas_Lower) ||
486 (Order == o2cblas_ColMajor && Trans == o2cblas_Trans &&
487 Uplo == o2cblas_Upper)) {
492 for (i=N;i>0 && i--;) {
494 const size_t j_min=0;
495 const size_t j_max=i;
497 for (j=j_min;j<j_max;j++) {
498 temp+=O2SCL_IX(x,jx)*O2SCL_IX2(A,i,j);
502 O2SCL_IX(x,ix)=temp+O2SCL_IX(x,ix)*O2SCL_IX2(A,i,i);
504 O2SCL_IX(x,ix)+=temp;
509 }
else if ((Order == o2cblas_RowMajor && Trans == o2cblas_Trans &&
510 Uplo == o2cblas_Upper) ||
511 (Order == o2cblas_ColMajor && Trans == o2cblas_NoTrans &&
512 Uplo == o2cblas_Lower)) {
516 for (i=N;i>0 && i--;) {
518 const size_t j_min=0;
519 const size_t j_max=i;
521 for (j=j_min;j<j_max;j++) {
522 temp+=O2SCL_IX(x,jx)*O2SCL_IX2(A,j,i);
526 O2SCL_IX(x,ix)=temp+O2SCL_IX(x,ix)*O2SCL_IX2(A,i,i);
528 O2SCL_IX(x,ix)+=temp;
533 }
else if ((Order == o2cblas_RowMajor && Trans == o2cblas_Trans &&
534 Uplo == o2cblas_Lower) ||
535 (Order == o2cblas_ColMajor && Trans == o2cblas_NoTrans &&
536 Uplo == o2cblas_Upper)) {
540 const size_t j_min=i+1;
541 const size_t j_max=N;
543 for (j=j_min;j<j_max;j++) {
544 temp+=O2SCL_IX(x,jx)*O2SCL_IX2(A,j,i);
548 O2SCL_IX(x,i)=temp+O2SCL_IX(x,i)*O2SCL_IX2(A,i,i);
555 O2SCL_ERR(
"Unrecognized operation in dtrmv().",
581 template<
class mat_t>
585 const size_t N,
const size_t K,
const double alpha,
586 const mat_t &A,
const mat_t &B,
const double beta, mat_t &C) {
592 if (alpha == 0.0 && beta == 1.0) {
604 if (Order == o2cblas_RowMajor) {
613 O2SCL_IX2(C,i,j)=0.0;
616 }
else if (beta != 1.0) {
619 O2SCL_IX2(C,i,j)*=beta;
628 TransF=(TransA == o2cblas_ConjTrans) ? o2cblas_Trans : TransA;
629 TransG=(TransB == o2cblas_ConjTrans) ? o2cblas_Trans : TransB;
631 if (TransF == o2cblas_NoTrans && TransG == o2cblas_NoTrans) {
637 const double temp=alpha*O2SCL_IX2(A,i,k);
640 O2SCL_IX2(C,i,j)+=temp*O2SCL_IX2(B,k,j);
646 }
else if (TransF == o2cblas_NoTrans && TransG == o2cblas_Trans) {
654 temp+=O2SCL_IX2(A,i,k)*O2SCL_IX2(B,j,k);
656 O2SCL_IX2(C,i,j)+=alpha*temp;
660 }
else if (TransF == o2cblas_Trans && TransG == o2cblas_NoTrans) {
664 const double temp=alpha*O2SCL_IX2(A,k,i);
667 O2SCL_IX2(C,i,j)+=temp*O2SCL_IX2(B,k,j);
673 }
else if (TransF == o2cblas_Trans && TransG == o2cblas_Trans) {
679 temp+=O2SCL_IX2(A,k,i)*O2SCL_IX2(B,j,k);
681 O2SCL_IX2(C,i,j)+=alpha*temp;
700 O2SCL_IX2(C,i,j)=0.0;
703 }
else if (beta != 1.0) {
706 O2SCL_IX2(C,i,j)*=beta;
715 TransF=(TransB == o2cblas_ConjTrans) ? o2cblas_Trans : TransB;
716 TransG=(TransA == o2cblas_ConjTrans) ? o2cblas_Trans : TransA;
718 if (TransF == o2cblas_NoTrans && TransG == o2cblas_NoTrans) {
724 const double temp=alpha*O2SCL_IX2(B,i,k);
727 O2SCL_IX2(C,i,j)+=temp*O2SCL_IX2(A,k,j);
733 }
else if (TransF == o2cblas_NoTrans && TransG == o2cblas_Trans) {
741 temp+=O2SCL_IX2(B,i,k)*O2SCL_IX2(A,j,k);
743 O2SCL_IX2(C,i,j)+=alpha*temp;
747 }
else if (TransF == o2cblas_Trans && TransG == o2cblas_NoTrans) {
751 const double temp=alpha*O2SCL_IX2(B,k,i);
754 O2SCL_IX2(C,i,j)+=temp*O2SCL_IX2(A,k,j);
760 }
else if (TransF == o2cblas_Trans && TransG == o2cblas_Trans) {
766 temp+=O2SCL_IX2(B,k,i)*O2SCL_IX2(A,j,k);
768 O2SCL_IX2(C,i,j)+=alpha*temp;
796 template<
class vec_t,
class vec2_t>
798 vec2_t &Y,
const size_t ie) {
802 if (alpha == 0.0)
return;
803 #if O2SCL_NO_RANGE_CHECK 810 const size_t m=(N-ie) % 4;
812 for (i=ie;i<ie+m;i++) {
813 O2SCL_IX(Y,i)+=alpha*O2SCL_IX(X,i);
817 O2SCL_IX(Y,i)+=alpha*O2SCL_IX(X,i);
818 O2SCL_IX(Y,i+1)+=alpha*O2SCL_IX(X,i+1);
819 O2SCL_IX(Y,i+2)+=alpha*O2SCL_IX(X,i+2);
820 O2SCL_IX(Y,i+3)+=alpha*O2SCL_IX(X,i+3);
833 template<
class vec_t,
class vec2_t>
834 double ddot_subvec(
const size_t N,
const vec_t &X,
const vec2_t &Y,
839 #if O2SCL_NO_RANGE_CHECK 846 const size_t m=(N-ie) % 4;
848 for (i=ie;i<ie+m;i++) {
849 r+=O2SCL_IX(X,i)*O2SCL_IX(Y,i);
853 r+=O2SCL_IX(X,i)*O2SCL_IX(Y,i);
854 r+=O2SCL_IX(X,i+1)*O2SCL_IX(Y,i+1);
855 r+=O2SCL_IX(X,i+2)*O2SCL_IX(Y,i+2);
856 r+=O2SCL_IX(X,i+3)*O2SCL_IX(Y,i+3);
874 template<
class vec_t>
880 #if O2SCL_NO_RANGE_CHECK 888 return fabs(O2SCL_IX(X,ie));
891 for (
size_t i=ie;i<N;i++) {
892 const double x=O2SCL_IX(X,i);
895 const double ax=fabs(x);
898 ssq=1.0+ssq*(scale/ax)*(scale/ax);
901 ssq+=(ax/scale)*(ax/scale);
907 return scale*sqrt(ssq);
919 template<
class vec_t>
923 #if O2SCL_NO_RANGE_CHECK 931 const size_t m=(N-ie) % 4;
933 for (i=ie;i<ie+m;i++) {
934 O2SCL_IX(X,i)*=alpha;
938 O2SCL_IX(X,i)*=alpha;
939 O2SCL_IX(X,i+1)*=alpha;
940 O2SCL_IX(X,i+2)*=alpha;
941 O2SCL_IX(X,i+3)*=alpha;
958 template<
class mat_t,
class vec_t>
960 const size_t ir,
const size_t ic, vec_t &y) {
962 #if O2SCL_NO_RANGE_CHECK 974 const size_t m=(M-ir) % 4;
976 for (i=ir;i<m+ir;i++) {
977 O2SCL_IX(y,i)+=alpha*O2SCL_IX2(X,i,ic);
981 O2SCL_IX(y,i)+=alpha*O2SCL_IX2(X,i,ic);
982 O2SCL_IX(y,i+1)+=alpha*O2SCL_IX2(X,i+1,ic);
983 O2SCL_IX(y,i+2)+=alpha*O2SCL_IX2(X,i+2,ic);
984 O2SCL_IX(y,i+3)+=alpha*O2SCL_IX2(X,i+3,ic);
1000 template<
class mat_t,
class vec_t>
1002 const size_t ic,
const vec_t &y) {
1003 #if O2SCL_NO_RANGE_CHECK 1012 const size_t m=(M-ir) % 4;
1014 for (i=ir;i<m+ir;i++) {
1015 r+=O2SCL_IX2(X,i,ic)*O2SCL_IX(y,i);
1019 r+=O2SCL_IX2(X,i,ic)*O2SCL_IX(y,i);
1020 r+=O2SCL_IX2(X,i+1,ic)*O2SCL_IX(y,i+1);
1021 r+=O2SCL_IX2(X,i+2,ic)*O2SCL_IX(y,i+2);
1022 r+=O2SCL_IX2(X,i+3,ic)*O2SCL_IX(y,i+3);
1043 template<
class mat_t>
1051 #if O2SCL_NO_RANGE_CHECK 1060 return fabs(O2SCL_IX2(A,ir,ic));
1063 for (i=ir;i<M;i++) {
1064 const double x=O2SCL_IX2(A,i,ic);
1067 const double ax=fabs(x);
1070 ssq=1.0+ssq*(scale/ax)*(scale/ax);
1073 ssq+=(ax/scale)*(ax/scale);
1079 return scale*sqrt(ssq);
1092 template<
class mat_t>
1094 const size_t M,
const double alpha) {
1096 #if O2SCL_NO_RANGE_CHECK 1104 const size_t m=(M-ir) % 4;
1106 for (i=ir;i<m+ir;i++) {
1107 O2SCL_IX2(A,i,ic)*=alpha;
1111 O2SCL_IX2(A,i,ic)*=alpha;
1112 O2SCL_IX2(A,i+1,ic)*=alpha;
1113 O2SCL_IX2(A,i+2,ic)*=alpha;
1114 O2SCL_IX2(A,i+3,ic)*=alpha;
1130 template<
class mat_t>
1134 #if O2SCL_NO_RANGE_CHECK 1143 const size_t m=(M-ir) % 4;
1145 for (i=ir;i<m+ir;i++) {
1146 r+=fabs(O2SCL_IX2(A,i,ic));
1150 r+=fabs(O2SCL_IX2(A,i,ic));
1151 r+=fabs(O2SCL_IX2(A,i+1,ic));
1152 r+=fabs(O2SCL_IX2(A,i+2,ic));
1153 r+=fabs(O2SCL_IX2(A,i+3,ic));
1176 template<
class mat_t,
class vec_t>
1178 const size_t ir,
const size_t ic, vec_t &Y) {
1180 #if O2SCL_NO_RANGE_CHECK 1192 const size_t m=(N-ic) % 4;
1194 for (i=ic;i<m+ic;i++) {
1195 O2SCL_IX(Y,i)+=alpha*O2SCL_IX2(X,ir,i);
1199 O2SCL_IX(Y,i)+=alpha*O2SCL_IX2(X,ir,i);
1200 O2SCL_IX(Y,i+1)+=alpha*O2SCL_IX2(X,ir,i+1);
1201 O2SCL_IX(Y,i+2)+=alpha*O2SCL_IX2(X,ir,i+2);
1202 O2SCL_IX(Y,i+3)+=alpha*O2SCL_IX2(X,ir,i+3);
1222 template<
class mat_t,
class vec_t>
1224 const size_t ic,
const vec_t &Y) {
1226 #if O2SCL_NO_RANGE_CHECK 1235 const size_t m=(N-ic) % 4;
1237 for (i=ic;i<m+ic;i++) {
1238 r+=O2SCL_IX2(X,ir,i)*O2SCL_IX(Y,i);
1242 r+=O2SCL_IX2(X,ir,i)*O2SCL_IX(Y,i);
1243 r+=O2SCL_IX2(X,ir,i+1)*O2SCL_IX(Y,i+1);
1244 r+=O2SCL_IX2(X,ir,i+2)*O2SCL_IX(Y,i+2);
1245 r+=O2SCL_IX2(X,ir,i+3)*O2SCL_IX(Y,i+3);
1261 template<
class mat_t>
1270 return fabs(O2SCL_IX2(M,ir,ic));
1273 for (i=ic;i<N;i++) {
1274 const double x=O2SCL_IX2(M,ir,i);
1277 const double ax=fabs(x);
1280 ssq=1.0+ssq*(scale/ax)*(scale/ax);
1283 ssq+=(ax/scale)*(ax/scale);
1289 return scale*sqrt(ssq);
1304 template<
class mat_t>
1306 const size_t N,
const double alpha) {
1308 #if O2SCL_NO_RANGE_CHECK 1316 const size_t m=(N-ic) % 4;
1318 for (i=ic;i<m+ic;i++) {
1319 O2SCL_IX2(A,ir,i)*=alpha;
1323 O2SCL_IX2(A,ir,i)*=alpha;
1324 O2SCL_IX2(A,ir,i+1)*=alpha;
1325 O2SCL_IX2(A,ir,i+2)*=alpha;
1326 O2SCL_IX2(A,ir,i+3)*=alpha;
void daxpy_subvec(const double alpha, const size_t N, const vec_t &X, vec2_t &Y, const size_t ie)
Compute beginning with index ie and ending with index N-1.
double dasum(const size_t N, const vec_t &X)
Compute the absolute sum of vector elements.
void dtrsv(const enum o2cblas_order order, const enum o2cblas_uplo Uplo, const enum o2cblas_transpose TransA, const enum o2cblas_diag Diag, const size_t M, const size_t N, const mat_t &A, vec_t &X)
Compute .
void daxpy_subcol(const double alpha, const size_t M, const mat_t &X, const size_t ir, const size_t ic, vec_t &y)
Compute for a subcolumn of a matrix.
void dgemm(const enum o2cblas_order Order, const enum o2cblas_transpose TransA, const enum o2cblas_transpose TransB, const size_t M, const size_t N, const size_t K, const double alpha, const mat_t &A, const mat_t &B, const double beta, mat_t &C)
Compute .
invalid argument supplied by user
o2cblas_uplo
Upper- or lower-triangular.
double ddot_subrow(const size_t N, const mat_t &X, const size_t ir, const size_t ic, const vec_t &Y)
Compute for a subrow of a matrix.
o2cblas_diag
Unit or generic diagonal.
double dnrm2(const size_t N, const vec_t &X)
Compute the norm of the vector X.
void dscal_subrow(mat_t &A, const size_t ir, const size_t ic, const size_t N, const double alpha)
Compute for a subrow of a matrix.
void dgemv(const enum o2cblas_order order, const enum o2cblas_transpose TransA, const size_t M, const size_t N, const double alpha, const mat_t &A, const vec_t &X, const double beta, vec2_t &Y)
Compute .
Namespace for O2scl CBLAS function templates.
void dtrmv(const enum o2cblas_order Order, const enum o2cblas_uplo Uplo, const enum o2cblas_transpose TransA, const enum o2cblas_diag Diag, const size_t N, const mat_t &A, vec_t &x)
Compute for the triangular matrix A.
o2cblas_transpose
Transpose operations.
double dasum_subcol(mat_t &A, const size_t ir, const size_t ic, const size_t M)
Compute for a subcolumn of a matrix.
double ddot(const size_t N, const vec_t &X, const vec2_t &Y)
Compute .
double dnrm2_subcol(const mat_t &A, const size_t ir, const size_t ic, const size_t M)
Compute the norm of a subcolumn of a matrix.
double dnrm2_subvec(const size_t N, const vec_t &X, const size_t ie)
Compute the norm of the vector X beginning with index ie and ending with index N-1.
double dnrm2_subrow(const mat_t &M, const size_t ir, const size_t ic, const size_t N)
Compute the norm of a subrow of a matrix.
#define O2SCL_ERR(d, n)
Set an error with message d and code n.
double ddot_subvec(const size_t N, const vec_t &X, const vec2_t &Y, const size_t ie)
Compute beginning with index ie and ending with index N-1.
void dscal(const double alpha, const size_t N, vec_t &X)
Compute .
o2cblas_side
Left or right sided operation.
void daxpy(const double alpha, const size_t N, const vec_t &X, vec2_t &Y)
Compute .
double ddot_subcol(const size_t M, const mat_t &X, const size_t ir, const size_t ic, const vec_t &y)
Compute for a subcolumn of a matrix.
void daxpy_subrow(const double alpha, const size_t N, const mat_t &X, const size_t ir, const size_t ic, vec_t &Y)
Compute for a subrow of a matrix.
o2cblas_order
Matrix order, either column-major or row-major.
void dscal_subvec(const double alpha, const size_t N, vec_t &X, const size_t ie)
Compute beginning with index ie and ending with index N-1.
void dscal_subcol(mat_t &A, const size_t ir, const size_t ic, const size_t M, const double alpha)
Compute for a subcolumn of a matrix.