forked from joernc/QGModel
-
Notifications
You must be signed in to change notification settings - Fork 0
/
model.py
556 lines (478 loc) · 20.9 KB
/
model.py
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
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
# Copyright (C) 2013,2014,2015 Joern Callies
#
# This file is part of QGModel.
#
# QGModel is free software: you can redistribute it and/or modify
# it under the terms of the GNU General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# QGModel is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU General Public License for more details.
#
# You should have received a copy of the GNU General Public License
# along with QGModel. If not, see <http://www.gnu.org/licenses/>.
import sys
import os
import pickle
import numpy as np
import pyfftw as fftw
import matplotlib.pyplot as plt
import matplotlib.colors as mcl
# Define blue colormap.
cdict = {
'red': ((0.0, 0.0, 0.0), (0.5, 0.216, 0.216), (1.0, 1.0, 1.0)),
'green': ((0.0, 0.0, 0.0), (0.5, 0.494, 0.494), (1.0, 1.0, 1.0)),
'blue': ((0.0, 0.0, 0.0), (0.5, 0.722, 0.722), (1.0, 1.0, 1.0))}
cm_blues = mcl.LinearSegmentedColormap('cm_blues', cdict, 256)
class Model:
"""
General QG model
This is the skeleton of a QG model that consists of a number of conserved
quantities that are advected horizontally. Implementations of this model
need to specify the number of conserved quantities (nz) and supply an
inversion relation that yields the streamfuncion given the conserved quan-
tities. The model geometry is doubly periodic in the perturbations; mean
flow and gradients in the conserved quantities can be prescribed.
"""
def __init__(self, nz):
self.nz = nz # number of levels
self.u = np.zeros(nz) # mean zonal velocities
self.v = np.zeros(nz) # mean meridional velocities
self.qx = np.zeros(nz) # mean zonal PV gradient
self.qy = np.zeros(nz) # mean meridional PV gradient
def linstab(self, k, l):
"""
Perform linear stability analysis for wavenumbers k and l.
Returns complex eigen-frequenies omega of the eigenvalue problem
[k U + l V + (k Qy - l Qx) L^-1] q = omega q.
The phase speeds are Re omega/k and the growth rates Im omega.
"""
# Set up inversion matrix with specified wavenumbers.
L = self.invmatrix(k[np.newaxis,:,np.newaxis],
l[:,np.newaxis,np.newaxis])
# Allow proper broadcasting over matrices.
kk = k[np.newaxis,:,np.newaxis,np.newaxis]
ll = l[:,np.newaxis,np.newaxis,np.newaxis]
# Compute (k Gy - l Gx) L^-1.
GL = np.einsum('...ij,...jk->...ik', kk*np.diag(self.qy)
- ll*np.diag(self.qx), np.linalg.inv(L))
# Solve eigenvalue problem.
w, v = np.linalg.eig(kk*np.diag(self.u) + ll*np.diag(self.v) + GL)
# Sort eigenvalues.
w.sort()
return w
def initnum(self, a, n, dt):
"""Initialize numerics."""
self.a = a # domain size
self.n = n # number of Fourier modes per direction
self.dt = dt # time step
self.diffexp = 2 # exponent of diffusion operator
self.hypodiff = 0. # hypodiffusion coefficient
self.threads = 1 # number of threads for FFT
self.clock = 0. # initial simulation time
# Set up grid.
self.grid()
# Set up inversion matrix.
self.L = self.invmatrix(self.k, self.l)
def grid(self):
"""Set up spectral and physical grid."""
# Set up spectral grid.
k = abs(np.fft.fftfreq(self.n, d=self.a/(2*np.pi*self.n))[:self.n/2+1])
l = np.fft.fftfreq(self.n, d=self.a/(2*np.pi*self.n))
self.k = k[np.newaxis,:,np.newaxis]
self.l = l[:,np.newaxis,np.newaxis]
# Set up physical grid.
x = np.arange(self.n) * self.a / self.n
y = np.arange(self.n) * self.a / self.n
self.x = x[np.newaxis,:,np.newaxis]
self.y = y[:,np.newaxis,np.newaxis]
def initq(self, qp):
"""Transform qp to spectral space and initialize q."""
self.q = fftw.interfaces.numpy_fft.rfft2(qp, axes=(0, 1),
threads=self.threads)
self.q[0,0,:] = 0. # ensuring zero mean
def timestep(self):
"""Perform time step."""
self.advection()
self.diffusion()
self.clock += self.dt
def advection(self):
"""Perform RK4 step for advective terms (linear and nonlinear)."""
q1 = self.advrhs(self.q)
q2 = self.advrhs(self.q + self.dt*q1/2)
q3 = self.advrhs(self.q + self.dt*q2/2)
q4 = self.advrhs(self.q + self.dt*q3)
self.q += self.dt*(q1 + 2*q2 + 2*q3 + q4)/6
def diffusion(self):
"""Perform implicit (hyper- and hypo-) diffusion step."""
k2 = self.k**2 + self.l**2
k2[0,0,:] = 1. # preventing div. by zero for wavenumber 0
self.q *= np.exp(-self.nu * k2**(self.diffexp/2.) * self.dt)
if self.hypodiff > 0:
self.q *= np.exp(-self.hypodiff / k2 * self.dt)
def advrhs(self, q):
"""
Calculate advective terms on RHS of PV equation.
Calculate mean-eddy and eddy-eddy advection terms:
u q'_x + v q'_y + u' q_x + v' q_y + J(p', q')
"""
# Perform inversion.
p = np.linalg.solve(self.L, q)
# Calculate RHS.
rhs = \
- 1j * (self.k*self.u + self.l*self.v) * q \
- 1j * (self.k*self.qy - self.l*self.qx) * p \
- self.jacobian(p, q)
return rhs
def jacobian(self, A, B):
"""
Calculate Jacobian A_x B_y - A_y B_x.
Transform Ax, Ay, Bx, By to physical space, perform multi-
plication and subtraction, and transform back to spectral space.
To avoid aliasing, apply 3/2 padding.
"""
Axp = self.ifft_pad(1j * self.k * A)
Ayp = self.ifft_pad(1j * self.l * A)
Bxp = self.ifft_pad(1j * self.k * B)
Byp = self.ifft_pad(1j * self.l * B)
return self.fft_truncate(Axp * Byp - Ayp * Bxp)
def fft_truncate(self, up):
"""Perform forward FFT on physical field up and truncate (3/2 rule)."""
us = fftw.interfaces.numpy_fft.rfft2(up, axes=(0, 1),
threads=self.threads)
u = np.zeros((self.n, self.n/2 + 1, self.nz), dtype=complex)
u[: self.n/2, :, :] = us[: self.n/2, : self.n/2 + 1, :]
u[self.n/2 :, :, :] = us[self.n : 3*self.n/2, : self.n/2 + 1, :]
return u/2.25 # accounting for normalization
def ifft_pad(self, u):
"""Pad spectral field u (3/2 rule) and perform inverse FFT."""
us = np.zeros((3*self.n/2, 3*self.n/4 + 1, self.nz), dtype=complex)
us[: self.n/2, : self.n/2 + 1, :] = u[: self.n/2, :, :]
us[self.n : 3*self.n/2, : self.n/2 + 1, :] = u[self.n/2 :, :, :]
return fftw.interfaces.numpy_fft.irfft2(2.25*us, axes=(0, 1),
threads=self.threads)
def doubleres(self):
"""Double the resolution, interpolate fields."""
self.n *= 2
# Pad spectral field.
qs = np.zeros((self.n, self.n/2 + 1, self.nz), dtype=complex)
qs[: self.n/4, : self.n/4 + 1, :] = self.q[: self.n/4, :, :]
qs[3*self.n/4 : self.n, : self.n/4 + 1, :] = self.q[self.n/4 :, :, :]
# Account for normalization.
self.q = 4*qs
# Update grid.
self.grid()
# Update inversion matrix.
self.L = self.invmatrix(self.k, self.l)
def screenlog(self):
"""Print model state info on screen."""
# Write time (in seconds).
sys.stdout.write(' {:15.0f}'.format(self.clock))
# Write mean enstrophy for each layer.
for i in range(self.nz):
sys.stdout.write(' {:5e}'.format(np.mean(np.abs(self.q[:,:,i])**2)
/self.n**2))
sys.stdout.write('\n')
def snapshot(self, name):
"""Save snapshots of total q (mean added) for each level."""
# Check whether directory exists.
if not os.path.isdir(name + '/snapshots'):
os.makedirs(name + '/snapshots')
# Transform to physical space.
qp = fftw.interfaces.numpy_fft.irfft2(self.q, axes=(0, 1),
threads=self.threads)
# Add mean gradients
qp += self.qx * (self.x - self.a / 2)
qp += self.qy * (self.y - self.a / 2)
# Determine range of colorbars.
m = np.max([np.abs(self.qx * self.a), np.abs(self.qy * self.a)], axis=0)
# Save image for each layer.
for i in range(self.nz):
plt.imsave(
name + '/snapshots/{:03d}_{:015.0f}.png'.format(i, self.clock),
qp[:,:,i], origin='lower', vmin=-m[i], vmax=m[i], cmap=cm_blues)
def save(self, name):
"""Save model state."""
# Check whether directory exists.
if not os.path.isdir(name + '/data'):
os.makedirs(name + '/data')
# Save.
with open(name + '/data/{:015.0f}'.format(self.clock), 'w') as f:
pickle.dump(self, f)
def load(name, time):
"""Load model state."""
with open(name + '/data/{:015.0f}'.format(time), 'r') as f:
m = pickle.load(f)
return m
class TwoDim(Model):
"""
Two dimensional model
This implements a two-dimensional model that consists of vorticity conser-
vation. The inversion relation is simply a Poisson equation.
"""
def __init__(self):
Model.__init__(self, 1)
def initmean(self, qx, qy):
"""Set up the mean state."""
# Set up PV gradients.
self.qx[0] = qx
self.qy[0] = qy
def invmatrix(self, k, l):
"""Set up the inversion matrix L."""
k2 = (k**2 + l**2)[:,:,0]
k2[k2 == 0.] = 1. # preventing div. by zero for wavenumber 0
L = np.empty((l.size, k.size, 1, 1))
L[:,:,0,0] = - k2
return L
class Surface(Model):
"""
Surface QG model
This implements an SQG model that consists of surface buoyancy conserva-
tion and implicit dynamics in an infinitely deep interior determined by
zero PV there. The conserved quantity here is PV-like and implements the
"oceanographic" case, where the surface is an upper surface. The dynamics
can also be used in an "atmospheric" case or a case with an interface
between two semi-infinite layers (see Held et al., 1995). The conserved
quantity is
q[0] = - f b(0) / N^2.
"""
def __init__(self):
Model.__init__(self, 1)
def initmean(self, f, N, Sx, Sy):
"""Set up the mean state."""
self.f = f # Coriolis parameter
self.N = N # buoyancy frequency
# Set up mean flow.
self.u = np.array([0])
self.v = np.array([0])
# Set up mean PV gradients.
self.qx = np.array([- f**2 * Sy / N**2])
self.qy = np.array([+ f**2 * Sx / N**2])
def invmatrix(self, k, l):
"""Set up the inversion matrix L."""
kh = np.hypot(k, l)[:,:,0]
kh[kh == 0.] = 1. # preventing div. by zero for wavenumber 0
L = np.empty((l.size, k.size, 1, 1))
L[:,:,0,0] = - self.f * kh / self.N
return L
class Layered(Model):
"""
Multi-layer model
This implements a multi-layer model with a rigid lid. See Vallis (2006)
for the formulation and notation. The layer thicknisses h, buoyancy
jumps g, mean flows (u, v), and beta can be passed to initmean.
"""
def initmean(self, f, h, g, u, v, beta):
"""Set up the mean state."""
self.f = f # Coriolis parameter
self.h = np.array(h) # layer thicknesses
self.g = np.array(g) # buoyancy jumps at interfaces
self.u = np.array(u) # mean zonal flow
self.v = np.array(v) # mean meridional flow
# Set up mean PV gradients.
self.qx = np.zeros(self.nz)
self.qx[:-1] -= f**2 * (self.v[:-1] - self.v[1:]) / (self.h[:-1] * g)
self.qx[1:] += f**2 * (self.v[:-1] - self.v[1:]) / (self.h[1::] * g)
self.qy = beta * np.ones(self.nz)
self.qy[:-1] += f**2 * (self.u[:-1] - self.u[1:]) / (self.h[:-1] * g)
self.qy[1:] -= f**2 * (self.u[:-1] - self.u[1:]) / (self.h[1::] * g)
def invmatrix(self, k, l):
"""Set up the inversion matrix L."""
k2 = (k**2 + l**2)[:,:,0]
k2[k2 == 0.] = 1. # preventing div. by zero for wavenumber 0
L = np.zeros((l.size, k.size, self.nz, self.nz))
for i in range(self.nz):
L[:,:,i,i] -= k2
if i > 0:
L[:,:,i,i-1] += self.f**2 / (self.h[i] * self.g[i-1])
L[:,:,i,i] -= self.f**2 / (self.h[i] * self.g[i-1])
if i < self.nz - 1:
L[:,:,i,i] -= self.f**2 / (self.h[i] * self.g[i])
L[:,:,i,i+1] += self.f**2 / (self.h[i] * self.g[i])
return L
class Eady(Model):
"""
Eady model
This implements an Eady model that consists of surface and bottom buoyancy
conservation and implicit interior dynamics determined by zero PV there.
The conserved quantities here are PV-like:
q[0] = - f b(0) / N^2,
q[1] = + f b(-H) / N^2.
"""
def __init__(self):
Model.__init__(self, 2)
def initmean(self, f, N, H, Sx, Sy):
"""Set up the mean state."""
self.f = f # Coriolis parameter
self.N = N # buoyancy frequency
self.H = H # depth
# Set up mean flow.
self.u = np.array([0, - Sx * H])
self.v = np.array([0, - Sy * H])
# Set up mean PV gradients.
self.qx = np.array([- f**2 * Sy / N**2, + f**2 * Sy / N**2])
self.qy = np.array([+ f**2 * Sx / N**2, - f**2 * Sx / N**2])
def invmatrix(self, k, l):
"""Set up the inversion matrix L."""
kh = np.hypot(k, l)[:,:,0]
kh[kh == 0.] = 1. # preventing div. by zero for wavenumber 0
mu = self.N * kh * self.H / self.f
L = np.empty((l.size, k.size, 2, 2))
L[:,:,0,0] = - self.f * kh / (self.N * np.tanh(mu))
L[:,:,0,1] = + self.f * kh / (self.N * np.sinh(mu))
L[:,:,1,0] = + self.f * kh / (self.N * np.sinh(mu))
L[:,:,1,1] = - self.f * kh / (self.N * np.tanh(mu))
return L
class FloatingEady(Model):
"""
Floating Eady model
This implements a "floating" Eady model that consists of a layer of
constant PV coupled to an infinitely deep layer below that also has
constant PV (see Callies, Flierl, Ferrari, Fox-Kemper, 2015). The model
consists of two conserved PV-like quantities at the surface and the inter-
face between the layers:
q[0] = - f b(0) / N[0]^2,
q[1] = + f [b^+(-H) / N[0]^2 - b^-(-H) / N[1]^2,
where N[0] and N[1] are the buoyancy frequencies of the Eady and deep
layers, respectively.
"""
def __init__(self):
Model.__init__(self, 2)
def initmean(self, f, N, H, Sx, Sy):
"""Set up the mean state."""
self.f = f # Coriolis parameter
self.N = np.array(N) # buoyancy frequencies of the two layers
self.H = H # depth of upper layer
# Set up mean flow.
self.u = np.array([0, - Sx[0] * H])
self.v = np.array([0, - Sy[0] * H])
# Set up mean PV gradients.
self.qx = np.array([
- f**2 * Sy[0] / N[0]**2,
+ f**2 * (Sy[0] / N[0]**2 - Sy[1] / N[1]**2)])
self.qy = np.array([
+ f**2 * Sx[0] / N[0]**2,
- f**2 * (Sx[0] / N[0]**2 - Sx[1] / N[1]**2)])
def invmatrix(self, k, l):
"""Set up the inversion matrix L."""
kh = np.hypot(k, l)[:,:,0]
kh[kh == 0.] = 1. # preventing div. by zero for wavenumber 0
mu = self.N[0] * kh * self.H / self.f
L = np.empty((l.size, k.size, 2, 2))
L[:,:,0,0] = - self.f * kh / (self.N[0] * np.tanh(mu))
L[:,:,0,1] = + self.f * kh / (self.N[0] * np.sinh(mu))
L[:,:,1,0] = + self.f * kh / (self.N[0] * np.sinh(mu))
L[:,:,1,1] = - self.f * kh / (self.N[0] * np.tanh(mu)) \
- self.f * kh / self.N[1]
return L
class TwoEady(Model):
"""
Two-Eady model
This implements a two-Eady model that consists of two layers of
constant PV coupled at a deformable interface (see Callies, Flierl,
Ferrari, Fox-Kemper, 2015). The model consists of three conserved PV-like
quantities at the surface, the interface between the layers, and the
bottom:
q[0] = - f b(0) / N[0]^2,
q[1] = + f [b^+(-H[0]) / N[0]^2 - b^-(-H[0]) / N[1]^2,
q[0] = + f b(-H[0]-H[1]) / N[1]^2,
where N[0] and N[1] are the buoyancy frequencies of the two layers and H[0]
and H[1] are their depths.
"""
def __init__(self):
Model.__init__(self, 3)
def initmean(self, f, N, H, Sx, Sy):
"""Set up the mean state."""
self.f = f # Coriolis parameter
self.N = np.array(N) # buoyancy frequencies of the two layers
self.H = np.array(H) # depths of the two layers
# Set up mean flow.
self.u = np.array([0, - Sx[0] * H[0], - Sx[0] * H[0] - Sx[1] * H[1]])
self.v = np.array([0, - Sy[0] * H[0], - Sy[0] * H[0] - Sy[1] * H[1]])
# Set up mean PV gradients.
self.qx = np.array([
- f**2 * Sy[0] / N[0]**2,
+ f**2 * (Sy[0] / N[0]**2 - Sy[1] / N[1]**2),
+ f**2 * Sy[1] / N[1]**2])
self.qy = np.array([
+ f**2 * Sx[0] / N[0]**2,
- f**2 * (Sx[0] / N[0]**2 - Sx[1] / N[1]**2),
- f**2 * Sx[1] / N[1]**2])
def invmatrix(self, k, l):
"""Set up the inversion matrix L."""
kh = np.hypot(k, l)
kh[kh == 0.] = 1. # preventing div. by zero for wavenumber 0
mu = self.N * kh * self.H / self.f
kh = kh[:,:,0]
L = np.zeros((l.size, k.size, 3, 3))
L[:,:,0,0] = - self.f * kh / (self.N[0] * np.tanh(mu[:,:,0]))
L[:,:,0,1] = + self.f * kh / (self.N[0] * np.sinh(mu[:,:,0]))
L[:,:,1,0] = + self.f * kh / (self.N[0] * np.sinh(mu[:,:,0]))
L[:,:,1,1] = - self.f * kh / (self.N[0] * np.tanh(mu[:,:,0])) \
- self.f * kh / (self.N[1] * np.tanh(mu[:,:,1]))
L[:,:,1,2] = + self.f * kh / (self.N[1] * np.sinh(mu[:,:,1]))
L[:,:,2,1] = + self.f * kh / (self.N[1] * np.sinh(mu[:,:,1]))
L[:,:,2,2] = - self.f * kh / (self.N[1] * np.tanh(mu[:,:,1]))
return L
class TwoEadyJump(Model):
"""
Two-Eady model with jump
This is a two-Eady model with an additional buoyancy and velocity jump
at the interface. The model then has an additional conserved quantity,
as if there was an additional layer at the interface. (The jump can be
thought of as an additional, infinitesimally thin layer with infinite
stratification and shear.) The conserved quantities are now
q[0] = - f b(0) / N[0]^2,
q[1] = + f b^+(-H[0]) / N[0]^2 + f \eta,
q[1] = - f b^-(-H[0]) / N[1]^2 - f \eta,
q[0] = + f b(-H[0]-H[1]) / N[1]^2,
where \eta is the interface displacement. The buoyancy jump g' and the
velocity jump (U, V) can be passed to initmean.
"""
def __init__(self):
Model.__init__(self, 4)
def initmean(self, f, N, H, g, Sx, Sy, U, V):
"""Set up the mean state."""
self.f = f # Coriolis parameter
self.N = np.array(N) # buoyancy frequencies of the two layers
self.H = np.array(H) # depths of the two layers
self.g = g # buoyancy jump at interface
# Set up mean flow.
self.u = np.array([0, - Sx[0] * H[0], - Sx[0] * H[0] - U,
- Sx[0] * H[0] - U - Sx[1] * H[1]])
self.v = np.array([0, - Sy[0] * H[0], - Sy[0] * H[0] - V,
- Sy[0] * H[0] - V - Sy[1] * H[1]])
# Set up mean PV gradients.
self.qx = np.array([
- f**2 * Sy[0] / N[0]**2,
+ f**2 * Sy[0] / N[0]**2 - f**2 * V / g,
- f**2 * Sy[1] / N[1]**2 + f**2 * V / g,
+ f**2 * Sy[1] / N[1]**2])
self.qy = np.array([
+ f**2 * Sx[0] / N[0]**2,
- f**2 * Sx[0] / N[0]**2 + f**2 * U / g,
+ f**2 * Sx[1] / N[1]**2 - f**2 * U / g,
- f**2 * Sx[1] / N[1]**2])
def invmatrix(self, k, l):
"""Set up the inversion matrix L."""
kh = np.hypot(k, l)
kh[kh == 0.] = 1. # preventing div. by zero for wavenumber 0
mu = self.N * kh * self.H / self.f
kh = kh[:,:,0]
L = np.zeros((l.size, k.size, 4, 4))
L[:,:,0,0] = - self.f * kh / (self.N[0] * np.tanh(mu[:,:,0]))
L[:,:,0,1] = + self.f * kh / (self.N[0] * np.sinh(mu[:,:,0]))
L[:,:,1,0] = + self.f * kh / (self.N[0] * np.sinh(mu[:,:,0]))
L[:,:,1,1] = - self.f * kh / (self.N[0] * np.tanh(mu[:,:,0])) \
- self.f**2 / self.g
L[:,:,1,2] = self.f**2 / self.g
L[:,:,2,1] = self.f**2 / self.g
L[:,:,2,2] = - self.f * kh / (self.N[1] * np.tanh(mu[:,:,1])) \
- self.f**2 / self.g
L[:,:,2,3] = + self.f * kh / (self.N[1] * np.sinh(mu[:,:,1]))
L[:,:,3,2] = + self.f * kh / (self.N[1] * np.sinh(mu[:,:,1]))
L[:,:,3,3] = - self.f * kh / (self.N[1] * np.tanh(mu[:,:,1]))
return L