1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
use super::Alt;
use rayon::prelude::*;

/// From https://github.com/fastscape-lem/fastscapelib-fortran/blob/master/src/Diffusion.f90
///
/// See https://fastscape-lem.github.io/fastscapelib-fortran
///
/// nx = grid size x
///
/// ny = grid size y
///
/// xl = x-dimension of the model topography in meters (double precision)
///
/// yl = y-dimension of the model topography in meters (double precision)
///
/// dt = length of the time step in years (double precision)
///
/// ibc = boundary conditions for grid.  For now, we assume all four boundaries
/// are fixed height, so this parameter is ignored.
///
/// h = heights of cells at each cell in the grid.
///
/// b = basement height at each cell in the grid (see https://en.wikipedia.org/wiki/Basement_(geology)).
///
/// kd = bedrock transport coefficient (or diffusivity) for hillslope processes
/// in meter squared per year (double precision) at each cell in the grid.
///
/// kdsed = sediment transport coefficient (or diffusivity) for hillslope
/// processes in meter squared per year; (double precision;) note that when
/// kdsed < 0, its value is not used, i.e., kd for sediment and bedrock have the
/// same value, regardless of sediment thickness
/* subroutine Diffusion ()

  ! subroutine to solve the diffusion equation by ADI

  use FastScapeContext

  implicit none
*/
#[allow(clippy::needless_late_init)]
pub fn diffusion(
    nx: usize,
    ny: usize,
    xl: f64,
    yl: f64,
    dt: f64,
    _ibc: (),
    h: &mut [Alt],
    b: &mut [Alt],
    kd: impl Fn(usize) -> f64,
    kdsed: f64,
) {
    let aij = |i: usize, j: usize| j * nx + i;
    /*
      double precision, dimension(:), allocatable :: f,diag,sup,inf,res
      double precision, dimension(:,:), allocatable :: zint,kdint,zintp
      integer i,j,ij
      double precision factxp,factxm,factyp,factym,dx,dy
    */
    let mut f: Vec<f64>;
    let mut diag: Vec<f64>;
    let mut sup: Vec<f64>;
    let mut inf: Vec<f64>;
    let mut res: Vec<f64>;
    let mut zint: Vec<f64>;
    let mut kdint: Vec<f64>;
    let mut zintp: Vec<f64>;
    let mut ij: usize;
    let mut factxp: f64;
    let mut factxm: f64;
    let mut factyp: f64;
    let mut factym: f64;
    let dx: f64;
    let dy: f64;
    /*
      character cbc*4

      !print*,'Diffusion'

      write (cbc,'(i4)') ibc

      dx=xl/(nx-1)
      dy=yl/(ny-1)
    */
    // 2048*32/2048/2048
    // 1 / 64 m
    dx = xl / /*(nx - 1)*/nx as f64;
    dy = yl / /*(ny - 1)*/ny as f64;
    /*
      ! creates 2D internal arrays to store topo and kd

      allocate (zint(nx,ny),kdint(nx,ny),zintp(nx,ny))
    */
    zint = vec![Default::default(); nx * ny];
    kdint = vec![Default::default(); nx * ny];
    /*
      do j=1,ny
        do i=1,nx
          ij=(j-1)*nx+i
          zint(i,j)=h(ij)
          kdint(i,j)=kd(ij)
          if (kdsed.gt.0.d0 .and. (h(ij)-b(ij)).gt.1.d-6) kdint(i,j)=kdsed
        enddo
      enddo

      zintp = zint
    */
    for j in 0..ny {
        for i in 0..nx {
            // ij = vec2_as_uniform_idx(i, j);
            ij = aij(i, j);
            zint[ij] = h[ij];
            kdint[ij] = kd(ij);
            if kdsed > 0.0 && (h[ij] - b[ij]) > 1.0e-6 {
                kdint[ij] = kdsed;
            }
        }
    }

    zintp = zint.clone();
    /*
      ! first pass along the x-axis

      allocate (f(nx),diag(nx),sup(nx),inf(nx),res(nx))
      f=0.d0
      diag=0.d0
      sup=0.d0
      inf=0.d0
      res=0.d0
      do j=2,ny-1
    */
    f = vec![0.0; nx];
    diag = vec![0.0; nx];
    sup = vec![0.0; nx];
    inf = vec![0.0; nx];
    res = vec![0.0; nx];
    for j in 1..ny - 1 {
        /*
            do i=2,nx-1
            factxp=(kdint(i+1,j)+kdint(i,j))/2.d0*(dt/2.)/dx**2
            factxm=(kdint(i-1,j)+kdint(i,j))/2.d0*(dt/2.)/dx**2
              factyp=(kdint(i,j+1)+kdint(i,j))/2.d0*(dt/2.)/dy**2
            factym=(kdint(i,j-1)+kdint(i,j))/2.d0*(dt/2.)/dy**2
            diag(i)=1.d0+factxp+factxm
            sup(i)=-factxp
            inf(i)=-factxm
            f(i)=zintp(i,j)+factyp*zintp(i,j+1)-(factyp+factym)*zintp(i,j)+factym*zintp(i,j-1)
            enddo
        */
        for i in 1..nx - 1 {
            factxp = (kdint[aij(i + 1, j)] + kdint[aij(i, j)]) / 2.0 * (dt / 2.0) / (dx * dx);
            factxm = (kdint[aij(i - 1, j)] + kdint[aij(i, j)]) / 2.0 * (dt / 2.0) / (dx * dx);
            factyp = (kdint[aij(i, j + 1)] + kdint[aij(i, j)]) / 2.0 * (dt / 2.0) / (dy * dy);
            factym = (kdint[aij(i, j - 1)] + kdint[aij(i, j)]) / 2.0 * (dt / 2.0) / (dy * dy);
            diag[i] = 1.0 + factxp + factxm;
            sup[i] = -factxp;
            inf[i] = -factxm;
            f[i] = zintp[aij(i, j)] + factyp * zintp[aij(i, j + 1)]
                - (factyp + factym) * zintp[aij(i, j)]
                + factym * zintp[aij(i, j - 1)];
        }
        /*
        ! left bc
            if (cbc(4:4).eq.'1') then
            diag(1)=1.
            sup(1)=0.
            f(1)=zintp(1,j)
            else
            factxp=(kdint(2,j)+kdint(1,j))/2.d0*(dt/2.)/dx**2
            factyp=(kdint(1,j+1)+kdint(1,j))/2.d0*(dt/2.)/dy**2
            factym=(kdint(1,j-1)+kdint(1,j))/2.d0*(dt/2.)/dy**2
            diag(1)=1.d0+factxp
            sup(1)=-factxp
            f(1)=zintp(1,j)+factyp*zintp(1,j+1)-(factyp+factym)*zintp(1,j)+factym*zintp(1,j-1)
            endif
        */
        if true {
            diag[0] = 1.0;
            sup[0] = 0.0;
            f[0] = zintp[aij(0, j)];
        } else {
            // reflective boundary
            factxp = (kdint[aij(1, j)] + kdint[aij(0, j)]) / 2.0 * (dt / 2.0) / (dx * dx);
            factyp = (kdint[aij(0, j + 1)] + kdint[aij(0, j)]) / 2.0 * (dt / 2.0) / (dy * dy);
            factym = (kdint[aij(0, j - 1)] + kdint[aij(0, j)]) / 2.0 * (dt / 2.0) / (dy * dy);
            diag[0] = 1.0 + factxp;
            sup[0] = -factxp;
            f[0] = zintp[aij(0, j)] + factyp * zintp[aij(0, j + 1)]
                - (factyp + factym) * zintp[aij(0, j)]
                + factym * zintp[aij(0, j - 1)];
        }
        /*
        ! right bc
            if (cbc(2:2).eq.'1') then
            diag(nx)=1.
            inf(nx)=0.
            f(nx)=zintp(nx,j)
            else
            factxm=(kdint(nx-1,j)+kdint(nx,j))/2.d0*(dt/2.)/dx**2
            factyp=(kdint(nx,j+1)+kdint(nx,j))/2.d0*(dt/2.)/dy**2
            factym=(kdint(nx,j-1)+kdint(nx,j))/2.d0*(dt/2.)/dy**2
            diag(nx)=1.d0+factxm
            inf(nx)=-factxm
            f(nx)=zintp(nx,j)+factyp*zintp(nx,j+1)-(factyp+factym)*zintp(nx,j)+factym*zintp(nx,j-1)
            endif
        */
        if true {
            diag[nx - 1] = 1.0;
            inf[nx - 1] = 0.0;
            f[nx - 1] = zintp[aij(nx - 1, j)];
        } else {
            // reflective boundary
            factxm = (kdint[aij(nx - 2, j)] + kdint[aij(nx - 1, j)]) / 2.0 * (dt / 2.0) / (dx * dx);
            factyp =
                (kdint[aij(nx - 1, j + 1)] + kdint[aij(nx - 1, j)]) / 2.0 * (dt / 2.0) / (dy * dy);
            factym =
                (kdint[aij(nx - 1, j - 1)] + kdint[aij(nx - 1, j)]) / 2.0 * (dt / 2.0) / (dy * dy);
            diag[nx - 1] = 1.0 + factxm;
            inf[nx - 1] = -factxm;
            f[nx - 1] = zintp[aij(nx - 1, j)] + factyp * zintp[aij(nx - 1, j + 1)]
                - (factyp + factym) * zintp[aij(nx - 1, j)]
                + factym * zintp[aij(nx - 1, j - 1)];
        }
        /*
          call tridag (inf,diag,sup,f,res,nx)
            do i=1,nx
            zint(i,j)=res(i)
            enddo
        */
        tridag(&inf, &diag, &sup, &f, &mut res, nx);
        for i in 0..nx {
            zint[aij(i, j)] = res[i];
        }
        /*
          enddo
          deallocate (f,diag,sup,inf,res)
        */
    }

    /*
      ! second pass along y-axis

      allocate (f(ny),diag(ny),sup(ny),inf(ny),res(ny))
      f=0.d0
      diag=0.d0
      sup=0.d0
      inf=0.d0
      res=0.d0
      do i=2,nx-1
    */
    f = vec![0.0; ny];
    diag = vec![0.0; ny];
    sup = vec![0.0; ny];
    inf = vec![0.0; ny];
    res = vec![0.0; ny];
    for i in 1..nx - 1 {
        /*
            do j=2,ny-1
            factxp=(kdint(i+1,j)+kdint(i,j))/2.d0*(dt/2.)/dx**2
            factxm=(kdint(i-1,j)+kdint(i,j))/2.d0*(dt/2.)/dx**2
            factyp=(kdint(i,j+1)+kdint(i,j))/2.d0*(dt/2.)/dy**2
            factym=(kdint(i,j-1)+kdint(i,j))/2.d0*(dt/2.)/dy**2
            diag(j)=1.d0+factyp+factym
            sup(j)=-factyp
            inf(j)=-factym
            f(j)=zint(i,j)+factxp*zint(i+1,j)-(factxp+factxm)*zint(i,j)+factxm*zint(i-1,j)
            enddo
        */
        for j in 1..ny - 1 {
            factxp = (kdint[aij(i + 1, j)] + kdint[aij(i, j)]) / 2.0 * (dt / 2.0) / (dx * dx);
            factxm = (kdint[aij(i - 1, j)] + kdint[aij(i, j)]) / 2.0 * (dt / 2.0) / (dx * dx);
            factyp = (kdint[aij(i, j + 1)] + kdint[aij(i, j)]) / 2.0 * (dt / 2.0) / (dy * dy);
            factym = (kdint[aij(i, j - 1)] + kdint[aij(i, j)]) / 2.0 * (dt / 2.0) / (dy * dy);
            diag[j] = 1.0 + factyp + factym;
            sup[j] = -factyp;
            inf[j] = -factym;
            f[j] = zint[aij(i, j)] + factxp * zint[aij(i + 1, j)]
                - (factxp + factxm) * zint[aij(i, j)]
                + factxm * zint[aij(i - 1, j)];
        }
        /*
        ! bottom bc
            if (cbc(1:1).eq.'1') then
            diag(1)=1.
            sup(1)=0.
            f(1)=zint(i,1)
            else
            factxp=(kdint(i+1,1)+kdint(i,j))/2.d0*(dt/2.)/dx**2
            factxm=(kdint(i-1,1)+kdint(i,1))/2.d0*(dt/2.)/dx**2
            factyp=(kdint(i,2)+kdint(i,1))/2.d0*(dt/2.)/dy**2
            diag(1)=1.d0+factyp
            sup(1)=-factyp
            f(1)=zint(i,1)+factxp*zint(i+1,1)-(factxp+factxm)*zint(i,1)+factxm*zint(i-1,1)
            endif
        */
        if true {
            diag[0] = 1.0;
            sup[0] = 0.0;
            f[0] = zint[aij(i, 0)];
        } else {
            // reflective boundary
            // TODO: Check whether this j was actually supposed to be a 0 in the original
            // (probably).
            // factxp = (kdint[aij(i+1, 0)] + kdint[aij(i, j)]) / 2.0 * (dt / 2.0) / (dx *
            // dx);
            factxp = (kdint[aij(i + 1, 0)] + kdint[aij(i, 0)]) / 2.0 * (dt / 2.0) / (dx * dx);
            factxm = (kdint[aij(i - 1, 0)] + kdint[aij(i, 0)]) / 2.0 * (dt / 2.0) / (dx * dx);
            factyp = (kdint[aij(i, 1)] + kdint[aij(i, 0)]) / 2.0 * (dt / 2.0) / (dy * dy);
            diag[0] = 1.0 + factyp;
            sup[0] = -factyp;
            f[0] = zint[aij(i, 0)] + factxp * zint[aij(i + 1, 0)]
                - (factxp + factxm) * zint[aij(i, 0)]
                + factxm * zint[aij(i - 1, 0)];
        }
        /*
        ! top bc
            if (cbc(3:3).eq.'1') then
            diag(ny)=1.
            inf(ny)=0.
            f(ny)=zint(i,ny)
            else
            factxp=(kdint(i+1,ny)+kdint(i,ny))/2.d0*(dt/2.)/dx**2
            factxm=(kdint(i-1,ny)+kdint(i,ny))/2.d0*(dt/2.)/dx**2
            factym=(kdint(i,ny-1)+kdint(i,ny))/2.d0*(dt/2.)/dy**2
            diag(ny)=1.d0+factym
            inf(ny)=-factym
            f(ny)=zint(i,ny)+factxp*zint(i+1,ny)-(factxp+factxm)*zint(i,ny)+factxm*zint(i-1,ny)
            endif
        */
        if true {
            diag[ny - 1] = 1.0;
            inf[ny - 1] = 0.0;
            f[ny - 1] = zint[aij(i, ny - 1)];
        } else {
            // reflective boundary
            factxp =
                (kdint[aij(i + 1, ny - 1)] + kdint[aij(i, ny - 1)]) / 2.0 * (dt / 2.0) / (dx * dx);
            factxm =
                (kdint[aij(i - 1, ny - 1)] + kdint[aij(i, ny - 1)]) / 2.0 * (dt / 2.0) / (dx * dx);
            factym = (kdint[aij(i, ny - 2)] + kdint[aij(i, ny - 1)]) / 2.0 * (dt / 2.0) / (dy * dy);
            diag[ny - 1] = 1.0 + factym;
            inf[ny - 1] = -factym;
            f[ny - 1] = zint[aij(i, ny - 1)] + factxp * zint[aij(i + 1, ny - 1)]
                - (factxp + factxm) * zint[aij(i, ny - 1)]
                + factxm * zint[aij(i - 1, ny - 1)];
        }
        /*
          call tridag (inf,diag,sup,f,res,ny)
            do j=1,ny
            zintp(i,j)=res(j)
            enddo
        */
        tridag(&inf, &diag, &sup, &f, &mut res, ny);
        for j in 0..ny {
            zintp[aij(i, j)] = res[j];
        }
        /*
          enddo
          deallocate (f,diag,sup,inf,res)
        */
    }
    /*
      ! stores result in 1D array

      do j=1,ny
        do i=1,nx
          ij=(j-1)*nx+i
          etot(ij)=etot(ij)+h(ij)-zintp(i,j)
          erate(ij)=erate(ij)+(h(ij)-zintp(i,j))/dt
          h(ij)=zintp(i,j)
        enddo
      enddo

      b=min(h,b)
    */
    for j in 0..ny {
        for i in 0..nx {
            ij = aij(i, j);
            // FIXME: Track total erosion and erosion rate.
            h[ij] = zintp[ij] as Alt;
        }
    }

    b.par_iter_mut().zip(h).for_each(|(b, h)| {
        *b = h.min(*b);
    });
    /*
      deallocate (zint,kdint,zintp)

      return

    end subroutine Diffusion
    */
}

/*
!----------

! subroutine to solve a tri-diagonal system of equations (from Numerical Recipes)

      SUBROUTINE tridag(a,b,c,r,u,n)

      implicit none

      INTEGER n
      double precision a(n),b(n),c(n),r(n),u(n)
*/
#[allow(clippy::needless_late_init)]
pub fn tridag(a: &[f64], b: &[f64], c: &[f64], r: &[f64], u: &mut [f64], n: usize) {
    /*
          INTEGER j
          double precision bet
          double precision,dimension(:),allocatable::gam

          allocate (gam(n))

          if(b(1).eq.0.d0) stop 'in tridag'
    */
    let mut bet: f64;
    let mut gam: Vec<f64>;

    gam = vec![Default::default(); n];

    assert_ne!(b[0], 0.0);
    /*

    ! first pass

    bet=b(1)
    u(1)=r(1)/bet
    do 11 j=2,n
      gam(j)=c(j-1)/bet
      bet=b(j)-a(j)*gam(j)
      if(bet.eq.0.) then
        print*,'tridag failed'
        stop
      endif
      u(j)=(r(j)-a(j)*u(j-1))/bet
      11    continue
    */
    bet = b[0];
    u[0] = r[0] / bet;
    for j in 1..n {
        gam[j] = c[j - 1] / bet;
        bet = b[j] - a[j] * gam[j];
        assert_ne!(bet, 0.0);
        // Round 0: u[0] = r[0] / b[0]
        //               = r'[0] / b'[0]
        // Round j: u[j] = (r[j] - a[j] * u'[j - 1]) / b'[j]
        //               = (r[j] - a[j] * r'[j - 1] / b'[j - 1]) / b'[j]
        //               = (r[j] - (a[j] / b'[j - 1]) * r'[j - 1]) / b'[j]
        //               = (r[j] - w[j] * r'[j - 1]) / b'[j]
        //               = r'[j] / b'[j]
        u[j] = (r[j] - a[j] * u[j - 1]) / bet;
    }
    /*
      ! second pass

      do 12 j=n-1,1,-1
        u(j)=u(j)-gam(j+1)*u(j+1)
        12    continue
    */
    for j in (0..n - 1).rev() {
        u[j] -= gam[j + 1] * u[j + 1];
    }
    /*
        deallocate (gam)

        return

        END
    */
}