mirror of
https://github.com/OPM/opm-simulators.git
synced 2024-11-22 09:16:27 -06:00
Correct a few embarrasing matrix dimension errors (in the "linear
pressure" part).
This commit is contained in:
parent
da327ba410
commit
ac97aa8c23
@ -78,16 +78,16 @@ mim_ip_linpress_exact(int nf, int d, double vol, double *N, double *K,
|
||||
double a1, a2, t;
|
||||
|
||||
t = 0.0;
|
||||
for (i = 0; i < nf; i++) {
|
||||
t += K[i + i*nf];
|
||||
for (i = 0; i < d; i++) {
|
||||
t += K[i + i*d];
|
||||
}
|
||||
|
||||
/* Step 4) T <- N*K */
|
||||
m = nf ; n = nf ; k = d;
|
||||
m = nf ; n = d ; k = d;
|
||||
ld1 = nf ; ld2 = d ;
|
||||
a1 = 1.0; a2 = 0.0;
|
||||
dgemm_("No Transpose", "No Transpose", &m, &n, &k,
|
||||
&a1, N, &ld1, K, &ld1, &a2, T, &ld1);
|
||||
&a1, N, &ld1, K, &ld2, &a2, T, &ld1);
|
||||
|
||||
/* Step 5) Binv <- (N*K*N' + t*X) / vol */
|
||||
a1 = 1.0 / vol ;
|
||||
|
Loading…
Reference in New Issue
Block a user