A Modular Arbitrary-Order Ocean-Atmosphere Model -- Stochastic implementation
inprod_analytic.f90
Go to the documentation of this file.
1 
2 ! inprod_analytic.f90
3 !
4 !> Inner products between the truncated set of basis functions for the
5 !> ocean and atmosphere streamfunction fields.
6 !> These are partly calculated using the analytical expressions from
7 !> Cehelsky, P., & Tung, K. K. : Theories of multiple equilibria and
8 !> weather regimes-A critical reexamination. Part II: Baroclinic two-layer
9 !> models. Journal of the atmospheric sciences, 44(21), 3282-3303, 1987.
10 !
11 !> @copyright
12 !> 2015 Lesley De Cruz & Jonathan Demaeyer.
13 !> See LICENSE.txt for license information.
14 !
15 !---------------------------------------------------------------------------!
16 !
17 !> @remark
18 !> Generated Fortran90/95 code
19 !> from inprod_analytic.lua
20 !
21 !---------------------------------------------------------------------------!
22 
24 
25  !-----------------------------------------------------!
26  ! !
27  ! Preamble and variables declaration !
28  ! !
29  !-----------------------------------------------------!
30 
31  USE params, only: nbatm, nboc, natm, noc, n, oms, ams, pi
32  IMPLICIT NONE
33 
34  PRIVATE
35 
37 
38  !> Atmospheric bloc specification type
39  TYPE :: atm_wavenum
40  CHARACTER :: typ
41  INTEGER :: m=0,p=0,h=0
42  REAL(KIND=8) :: nx=0.,ny=0.
43  END TYPE atm_wavenum
44 
45  !> Oceanic bloc specification type
46  TYPE :: ocean_wavenum
47  INTEGER :: p,h
48  REAL(KIND=8) :: nx,ny
49  END TYPE ocean_wavenum
50 
51  !> Type holding the atmospheric inner products tensors
52  TYPE :: atm_tensors
53  REAL(KIND=8), DIMENSION(:,:), ALLOCATABLE :: a,c,d,s
54  REAL(KIND=8), DIMENSION(:,:,:), ALLOCATABLE :: b,g
55  END TYPE atm_tensors
56 
57  !> Type holding the oceanic inner products tensors
58  TYPE :: ocean_tensors
59  REAL(KIND=8), DIMENSION(:,:), ALLOCATABLE :: k,m,n,w
60  REAL(KIND=8), DIMENSION(:,:,:), ALLOCATABLE :: o,c
61  END TYPE ocean_tensors
62 
63  !> Atmospheric blocs specification
64  TYPE(atm_wavenum), DIMENSION(:), ALLOCATABLE, PUBLIC :: awavenum
65  !> Oceanic blocs specification
66  TYPE(ocean_wavenum), DIMENSION(:), ALLOCATABLE, PUBLIC :: owavenum
67 
68  !> Atmospheric tensors
69  TYPE(atm_tensors), PUBLIC :: atmos
70  !> Oceanic tensors
71  TYPE(ocean_tensors), PUBLIC :: ocean
72 
73 
74  !-----------------------------------------------------!
75  ! !
76  ! End of preamble !
77  ! !
78  !-----------------------------------------------------!
79 
80 CONTAINS
81 
82  !-----------------------------------------------------!
83  ! !
84  ! Definition of the Helper functions from Cehelsky !
85  ! & Tung !
86  ! !
87  !-----------------------------------------------------!
88 
89  !> Cehelsky & Tung Helper functions
90  REAL(KIND=8) FUNCTION b1(Pi, Pj, Pk)
91  INTEGER :: Pi,Pj,Pk
92  b1 = (pk + pj) / REAL(pi)
93  END FUNCTION b1
94 
95  !> Cehelsky & Tung Helper functions
96  REAL(KIND=8) FUNCTION b2(Pi, Pj, Pk)
97  INTEGER :: Pi,Pj,Pk
98  b2 = (pk - pj) / REAL(pi)
99  END FUNCTION b2
100 
101  !> Integer Dirac delta function
102  REAL(KIND=8) FUNCTION delta(r)
103  INTEGER :: r
104  IF (r==0) THEN
105  delta = 1.d0
106  ELSE
107  delta = 0.d0
108  ENDIF
109  END FUNCTION delta
110 
111  !> "Odd or even" function
112  REAL(KIND=8) FUNCTION flambda(r)
113  INTEGER :: r
114  IF (mod(r,2)==0) THEN
115  flambda = 0.d0
116  ELSE
117  flambda = 1.d0
118  ENDIF
119  END FUNCTION flambda
120 
121  !> Cehelsky & Tung Helper functions
122  REAL(KIND=8) FUNCTION s1(Pj, Pk, Mj, Hk)
123  INTEGER :: Pk,Pj,Mj,Hk
124  s1 = -((pk * mj + pj * hk)) / 2.d0
125  END FUNCTION s1
126 
127  !> Cehelsky & Tung Helper functions
128  REAL(KIND=8) FUNCTION s2(Pj, Pk, Mj, Hk)
129  INTEGER :: Pk,Pj,Mj,Hk
130  s2 = (pk * mj - pj * hk) / 2.d0
131  END FUNCTION s2
132 
133  !> Cehelsky & Tung Helper functions
134  REAL(KIND=8) FUNCTION s3(Pj, Pk, Hj, Hk)
135  INTEGER :: Pj,Pk,Hj,Hk
136  s3 = (pk * hj + pj * hk) / 2.d0
137  END FUNCTION s3
138 
139  !> Cehelsky & Tung Helper functions
140  REAL(KIND=8) FUNCTION s4(Pj, Pk, Hj, Hk)
141  INTEGER :: Pj,Pk,Hj,Hk
142  s4 = (pk * hj - pj * hk) / 2.d0
143  END FUNCTION s4
144 
145  !-----------------------------------------------------!
146  ! Inner products definition routines !
147  !--------------------------------------------------------!
148  ! 1. Inner products in the equations for the atmosphere !
149  !--------------------------------------------------------!
150 
151  !> Eigenvalues of the Laplacian (atmospheric)
152  !>
153  !> \f$ a_{i,j} = (F_i, \nabla^2 F_j)\f$ .
154  SUBROUTINE calculate_a
155  INTEGER :: i
156  TYPE(atm_wavenum) :: Ti
157  INTEGER :: AllocStat
158  IF (natm == 0 ) THEN
159  stop "*** Problem with calculate_a : natm==0 ! ***"
160  ELSE
161  IF (.NOT. ALLOCATED(atmos%a)) THEN
162  ALLOCATE(atmos%a(natm,natm), stat=allocstat)
163  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
164  END IF
165  END IF
166  atmos%a=0.d0
167 
168  DO i=1,natm
169  ti = awavenum(i)
170  atmos%a(i,i) = -(n**2) * ti%Nx**2 - ti%Ny**2
171  ENDDO
172  END SUBROUTINE calculate_a
173 
174  !> Streamfunction advection terms (atmospheric)
175  !>
176  !> \f$ b_{i,j,k} = (F_i, J(F_j, \nabla^2 F_k))\f$ .
177  !
178  !> @remark
179  !> Atmospheric g and a tensors must be computed before calling
180  !> this routine
181  SUBROUTINE calculate_b
182  INTEGER :: i,j,k
183  INTEGER :: AllocStat
184 
185  IF ((.NOT. ALLOCATED(atmos%a)) .OR. (.NOT. ALLOCATED(atmos%g))) THEN
186  stop "*** atmos%a and atmos%g must be defined before calling calculate_b ! ***"
187  END IF
188 
189  IF (natm == 0 ) THEN
190  stop "*** Problem with calculate_b : natm==0 ! ***"
191  ELSE
192  IF (.NOT. ALLOCATED(atmos%b)) THEN
193  ALLOCATE(atmos%b(natm,natm,natm), stat=allocstat)
194  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
195  END IF
196  END IF
197  atmos%b=0.d0
198 
199  DO i=1,natm
200  DO j=1,natm
201  DO k=1,natm
202  atmos%b(i,j,k)= atmos%a(k,k) * atmos%g(i,j,k)
203  END DO
204  END DO
205  END DO
206  END SUBROUTINE calculate_b
207 
208  !> Beta term for the atmosphere
209  !>
210  !> \f$ c_{i,j} = (F_i, \partial_x F_j)\f$ .
211  !
212  !> @remark
213  !> Strict function !! Only accepts KL type.
214  !> For any other combination, it will not calculate anything
215  SUBROUTINE calculate_c_atm
216  INTEGER :: i,j
217  TYPE(atm_wavenum) :: Ti, Tj
218  REAL(KIND=8) :: val
219  INTEGER :: AllocStat
220 
221  IF (natm == 0 ) THEN
222  stop "*** Problem with calculate_c_atm : natm==0 ! ***"
223  ELSE
224  IF (.NOT. ALLOCATED(atmos%c)) THEN
225  ALLOCATE(atmos%c(natm,natm), stat=allocstat)
226  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
227  END IF
228  END IF
229  atmos%c=0.d0
230 
231  DO i=1,natm
232  DO j=1,natm
233  ti = awavenum(i)
234  tj = awavenum(j)
235  val = 0.d0
236  IF ((ti%typ == "K") .AND. (tj%typ == "L")) THEN
237  val = n * ti%M * delta(ti%M - tj%H) * delta(ti%P - tj%P)
238  END IF
239  IF (val /= 0.d0) THEN
240  atmos%c(i,j)=val
241  atmos%c(j,i)= - val
242  ENDIF
243  END DO
244  END DO
245  END SUBROUTINE calculate_c_atm
246 
247  !> Forcing of the ocean on the atmosphere.
248  !>
249  !> \f$ d_{i,j} = (F_i, \nabla^2 \eta_j)\f$ .
250  !
251  !> @remark
252  !> Atmospheric s tensor and oceanic M tensor must be computed before
253  !> calling this routine !
254  SUBROUTINE calculate_d
255  INTEGER :: i,j
256  INTEGER :: AllocStat
257 
258  IF ((.NOT. ALLOCATED(atmos%s)) .OR. (.NOT. ALLOCATED(ocean%M))) THEN
259  stop "*** atmos%s and ocean%M must be defined before calling calculate_d ! ***"
260  END IF
261 
262 
263  IF (natm == 0 ) THEN
264  stop "*** Problem with calculate_d : natm==0 ! ***"
265  ELSE
266  IF (.NOT. ALLOCATED(atmos%d)) THEN
267  ALLOCATE(atmos%d(natm,noc), stat=allocstat)
268  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
269  END IF
270  END IF
271  atmos%d=0.d0
272 
273  DO i=1,natm
274  DO j=1,noc
275  atmos%d(i,j)=atmos%s(i,j) * ocean%M(j,j)
276  END DO
277  END DO
278  END SUBROUTINE calculate_d
279 
280  !> Temperature advection terms (atmospheric)
281  !>
282  !> \f$ g_{i,j,k} = (F_i, J(F_j, F_k))\f$ .
283  !
284  ! @remark
285  ! This is a strict function: it only accepts AKL KKL and LLL types.
286  ! For any other combination, it will not calculate anything.
287  SUBROUTINE calculate_g
288  INTEGER :: i,j,k
289  TYPE(atm_wavenum) :: Ti, Tj, Tk
290  REAL(KIND=8) :: val,vb1, vb2, vs1, vs2, vs3, vs4
291  INTEGER :: AllocStat
292 
293  IF (natm == 0 ) THEN
294  stop "*** Problem with calculate_g : natm==0 ! ***"
295  ELSE
296  IF (.NOT. ALLOCATED(atmos%g)) THEN
297  ALLOCATE(atmos%g(natm,natm,natm), stat=allocstat)
298  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
299  END IF
300  END IF
301  atmos%g=0.d0
302 
303  DO i=1,natm
304  DO j=1,natm
305  DO k=1,natm
306  ti = awavenum(i)
307  tj = awavenum(j)
308  tk = awavenum(k)
309  val=0.d0
310  IF ((ti%typ == "A") .AND. (tj%typ == "K") .AND. (tk%typ == "L")) THEN
311  vb1 = b1(ti%P,tj%P,tk%P)
312  vb2 = b2(ti%P,tj%P,tk%P)
313  val = -2 * sqrt(2.) / pi * tj%M * delta(tj%M - tk%H) * flambda(ti%P + tj%P + tk%P)
314  IF (val /= 0.d0) val = val * (vb1**2 / (vb1**2 - 1) - vb2**2 / (vb2**2 - 1))
315  ELSEIF ((ti%typ == "K") .AND. (tj%typ == "K") .AND. (tk%typ == "L")) THEN
316  vs1 = s1(tj%P,tk%P,tj%M,tk%H)
317  vs2 = s2(tj%P,tk%P,tj%M,tk%H)
318  val = vs1 * (delta(ti%M - tk%H - tj%M) * delta(ti%P -&
319  & tk%P + tj%P) - delta(ti%M- tk%H - tj%M) *&
320  & delta(ti%P + tk%P - tj%P) + (delta(tk%H - tj%M&
321  & + ti%M) + delta(tk%H - tj%M - ti%M)) *&
322  & delta(tk%P + tj%P - ti%P)) + vs2 * (delta(ti%M&
323  & - tk%H - tj%M) * delta(ti%P - tk%P - tj%P) +&
324  & (delta(tk%H - tj%M - ti%M) + delta(ti%M + tk%H&
325  & - tj%M)) * (delta(ti%P - tk%P + tj%P) -&
326  & delta(tk%P - tj%P + ti%P)))
327  END IF
328  val=val*n
329  IF (val /= 0.d0) THEN
330  atmos%g(i,j,k) = val
331  atmos%g(j,k,i) = val
332  atmos%g(k,i,j) = val
333  atmos%g(i,k,j) = -val
334  atmos%g(j,i,k) = -val
335  atmos%g(k,j,i) = -val
336  ENDIF
337  END DO
338  END DO
339  END DO
340 
341  DO i=1,natm
342  DO j=i,natm
343  DO k=j,natm
344  ti = awavenum(i)
345  tj = awavenum(j)
346  tk = awavenum(k)
347  val=0.d0
348 
349  IF ((ti%typ == "L") .AND. (tj%typ == "L") .AND. (tk%typ == "L")) THEN
350  vs3 = s3(tj%P,tk%P,tj%H,tk%H)
351  vs4 = s4(tj%P,tk%P,tj%H,tk%H)
352  val = vs3 * ((delta(tk%H - tj%H - ti%H) - delta(tk%H &
353  &- tj%H + ti%H)) * delta(tk%P + tj%P - ti%P) +&
354  & delta(tk%H + tj%H - ti%H) * (delta(tk%P - tj%P&
355  & + ti%P) - delta(tk%P - tj%P - ti%P))) + vs4 *&
356  & ((delta(tk%H + tj%H - ti%H) * delta(tk%P - tj&
357  &%P - ti%P)) + (delta(tk%H - tj%H + ti%H) -&
358  & delta(tk%H - tj%H - ti%H)) * (delta(tk%P - tj&
359  &%P - ti%P) - delta(tk%P - tj%P + ti%P)))
360  ENDIF
361  val=val*n
362  IF (val /= 0.d0) THEN
363  atmos%g(i,j,k) = val
364  atmos%g(j,k,i) = val
365  atmos%g(k,i,j) = val
366  atmos%g(i,k,j) = -val
367  atmos%g(j,i,k) = -val
368  atmos%g(k,j,i) = -val
369  ENDIF
370  ENDDO
371  ENDDO
372  ENDDO
373 
374  END SUBROUTINE calculate_g
375 
376  !> Forcing (thermal) of the ocean on the atmosphere.
377  !>
378  !> \f$ s_{i,j} = (F_i, \eta_j)\f$ .
379  SUBROUTINE calculate_s
380  INTEGER :: i,j
381  TYPE(atm_wavenum) :: Ti
382  TYPE(ocean_wavenum) :: Dj
383  REAL(KIND=8) :: val
384  INTEGER :: AllocStat
385  IF (natm == 0 ) THEN
386  stop "*** Problem with calculate_s : natm==0 ! ***"
387  ELSEIF (noc == 0) then
388  stop "*** Problem with calculate_s : noc==0 ! ***"
389  ELSE
390  IF (.NOT. ALLOCATED(atmos%s)) THEN
391  ALLOCATE(atmos%s(natm,noc), stat=allocstat)
392  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
393  END IF
394  END IF
395  atmos%s=0.d0
396 
397  DO i=1,natm
398  DO j=1,noc
399  ti = awavenum(i)
400  dj = owavenum(j)
401  val=0.d0
402  IF (ti%typ == "A") THEN
403  val = flambda(dj%H) * flambda(dj%P + ti%P)
404  IF (val /= 0.d0) THEN
405  val = val*8*sqrt(2.)*dj%P/(pi**2 * (dj%P**2 - ti%P**2) * dj%H)
406  END IF
407  ELSEIF (ti%typ == "K") THEN
408  val = flambda(2 * ti%M + dj%H) * delta(dj%P - ti%P)
409  IF (val /= 0.d0) THEN
410  val = val*4*dj%H/(pi * (-4 * ti%M**2 + dj%H**2))
411  END IF
412  ELSEIF (ti%typ == "L") THEN
413  val = delta(dj%P - ti%P) * delta(2 * ti%H - dj%H)
414  END IF
415  IF (val /= 0.d0) THEN
416  atmos%s(i,j)=val
417  ENDIF
418  END DO
419  END DO
420  END SUBROUTINE calculate_s
421 
422  !--------------------------------------------------------!
423  ! 2. Inner products in the equations for the ocean !
424  !--------------------------------------------------------!
425 
426  !> Forcing of the atmosphere on the ocean.
427  !>
428  !> \f$ K_{i,j} = (\eta_i, \nabla^2 F_j)\f$ .
429  !
430  !> @remark
431  !> atmospheric a and s tensors must be computed before calling
432  !> this function !
433  SUBROUTINE calculate_k
434  INTEGER :: i,j
435  INTEGER :: AllocStat
436 
437  IF ((.NOT. ALLOCATED(atmos%a)) .OR. (.NOT. ALLOCATED(atmos%s))) THEN
438  stop "*** atmos%a and atmos%s must be defined before calling calculate_K ! ***"
439  END IF
440 
441  IF (noc == 0 ) THEN
442  stop "*** Problem with calculate_K : noc==0 ! ***"
443  ELSEIF (natm == 0 ) THEN
444  stop "*** Problem with calculate_K : natm==0 ! ***"
445  ELSE
446  IF (.NOT. ALLOCATED(ocean%K)) THEN
447  ALLOCATE(ocean%K(noc,natm), stat=allocstat)
448  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
449  END IF
450  END IF
451  ocean%K=0.d0
452 
453  DO i=1,noc
454  DO j=1,natm
455  ocean%K(i,j) = atmos%s(j,i) * atmos%a(j,j)
456  END DO
457  END DO
458  END SUBROUTINE calculate_k
459 
460  !> Forcing of the ocean fields on the ocean.
461  !>
462  !> \f$ M_{i,j} = (eta_i, \nabla^2 \eta_j)\f$ .
463  SUBROUTINE calculate_m
464  INTEGER :: i
465  TYPE(ocean_wavenum) :: Di
466  INTEGER :: AllocStat
467  IF (noc == 0 ) THEN
468  stop "*** Problem with calculate_M : noc==0 ! ***"
469  ELSE
470  IF (.NOT. ALLOCATED(ocean%M)) THEN
471  ALLOCATE(ocean%M(noc,noc), stat=allocstat)
472  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
473  END IF
474  END IF
475  ocean%M=0.d0
476 
477  DO i=1,noc
478  di = owavenum(i)
479  ocean%M(i,i) = -(n**2) * di%Nx**2 - di%Ny**2
480  END DO
481  END SUBROUTINE calculate_m
482 
483  !> Beta term for the ocean
484  !>
485  !> \f$ N_{i,j} = (\eta_i, \partial_x \eta_j) \f$.
486  SUBROUTINE calculate_n
487  INTEGER :: i,j
488  TYPE(ocean_wavenum) :: Di,Dj
489  REAL(KIND=8) :: val
490  INTEGER :: AllocStat
491  IF (noc == 0 ) THEN
492  stop "*** Problem with calculate_N : noc==0 ! ***"
493  ELSE
494  IF (.NOT. ALLOCATED(ocean%N)) THEN
495  ALLOCATE(ocean%N(noc,noc), stat=allocstat)
496  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
497  END IF
498  END IF
499  ocean%N=0.d0
500  val=0.d0
501 
502  DO i=1,noc
503  DO j=1,noc
504  di = owavenum(i)
505  dj = owavenum(j)
506  val = delta(di%P - dj%P) * flambda(di%H + dj%H)
507  IF (val /= 0.d0) ocean%N(i,j) = val * (-2) * dj%H * di%H * n / ((dj%H**2 - di%H**2) * pi)
508  END DO
509  END DO
510  END SUBROUTINE calculate_n
511 
512  !> Temperature advection term (passive scalar)
513  !>
514  !> \f$ O_{i,j,k} = (\eta_i, J(\eta_j, \eta_k))\f$ .
515  SUBROUTINE calculate_o
516  INTEGER :: i,j,k
517  REAL(KIND=8) :: vs3,vs4,val
518  TYPE(ocean_wavenum) :: Di,Dj,Dk
519  INTEGER :: AllocStat
520  IF (noc == 0 ) THEN
521  stop "*** Problem with calculate_O : noc==0 ! ***"
522  ELSE
523  IF (.NOT. ALLOCATED(ocean%O)) THEN
524  ALLOCATE(ocean%O(noc,noc,noc), stat=allocstat)
525  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
526  END IF
527  END IF
528  ocean%O=0.d0
529  val=0.d0
530 
531  DO i=1,noc
532  DO j=i,noc
533  DO k=j,noc
534  di = owavenum(i)
535  dj = owavenum(j)
536  dk = owavenum(k)
537  vs3 = s3(dj%P,dk%P,dj%H,dk%H)
538  vs4 = s4(dj%P,dk%P,dj%H,dk%H)
539  val = vs3*((delta(dk%H - dj%H - di%H) - delta(dk%H - dj&
540  &%H + di%H)) * delta(dk%P + dj%P - di%P) + delta(dk&
541  &%H + dj%H - di%H) * (delta(dk%P - dj%P + di%P) -&
542  & delta(dk%P - dj%P - di%P))) + vs4 * ((delta(dk%H &
543  &+ dj%H - di%H) * delta(dk%P - dj%P - di%P)) +&
544  & (delta(dk%H - dj%H + di%H) - delta(dk%H - dj%H -&
545  & di%H)) * (delta(dk%P - dj%P - di%P) - delta(dk%P &
546  &- dj%P + di%P)))
547  val = val * n / 2
548  IF (val /= 0.d0) THEN
549  ocean%O(i,j,k) = val
550  ocean%O(j,k,i) = val
551  ocean%O(k,i,j) = val
552  ocean%O(i,k,j) = -val
553  ocean%O(j,i,k) = -val
554  ocean%O(k,j,i) = -val
555  END IF
556  END DO
557  END DO
558  END DO
559  END SUBROUTINE calculate_o
560 
561  !> Streamfunction advection terms (oceanic)
562  !>
563  !> \f$ C_{i,j,k} = (\eta_i, J(\eta_j,\nabla^2 \eta_k))\f$ .
564  !
565  !> @remark
566  !> Requires O_{i,j,k} and M_{i,j} to be calculated beforehand.
567  SUBROUTINE calculate_c_oc
568  INTEGER :: i,j,k
569  REAL(KIND=8) :: val
570  INTEGER :: AllocStat
571 
572  IF ((.NOT. ALLOCATED(ocean%O)) .OR. (.NOT. ALLOCATED(ocean%M))) THEN
573  stop "*** ocean%O and ocean%M must be defined before calling calculate_C ! ***"
574  END IF
575 
576  IF (noc == 0 ) THEN
577  stop "*** Problem with calculate_C : noc==0 ! ***"
578  ELSE
579  IF (.NOT. ALLOCATED(ocean%C)) THEN
580  ALLOCATE(ocean%C(noc,noc,noc), stat=allocstat)
581  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
582  END IF
583  END IF
584  ocean%C=0.d0
585  val=0.d0
586 
587  DO i=1,noc
588  DO j=1,noc
589  DO k=1,noc
590  val = ocean%M(k,k) * ocean%O(i,j,k)
591  IF (val /= 0.d0) ocean%C(i,j,k) = val
592  END DO
593  END DO
594  END DO
595  END SUBROUTINE calculate_c_oc
596 
597  !> Short-wave radiative forcing of the ocean
598  !>
599  !> \f$ W_{i,j} = (\eta_i, F_j)\f$ .
600  !
601  !> @remark
602  !> atmospheric s tensor must be computed before calling
603  !> this function !
604  SUBROUTINE calculate_w
605  INTEGER :: i,j
606  INTEGER :: AllocStat
607 
608  IF (.NOT. ALLOCATED(atmos%s)) THEN
609  stop "*** atmos%s must be defined before calling calculate_W ! ***"
610  END IF
611 
612  IF (noc == 0 ) THEN
613  stop "*** Problem with calculate_W : noc==0 ! ***"
614  ELSEIF (natm == 0 ) THEN
615  stop "*** Problem with calculate_W : natm==0 ! ***"
616  ELSE
617  IF (.NOT. ALLOCATED(ocean%W)) THEN
618  ALLOCATE(ocean%W(noc,natm), stat=allocstat)
619  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
620  END IF
621  END IF
622  ocean%W=0.d0
623 
624  DO i=1,noc
625  DO j=1,natm
626  ocean%W(i,j) = atmos%s(j,i)
627  END DO
628  END DO
629  END SUBROUTINE calculate_w
630 
631  !-----------------------------------------------------!
632  ! !
633  ! Initialisation routine !
634  ! !
635  !-----------------------------------------------------!
636 
637  !> Initialisation of the inner product
638  SUBROUTINE init_inprod
639  INTEGER :: i,j
640  INTEGER :: AllocStat
641 
642  ! Definition of the types and wave numbers tables
643 
644  ALLOCATE(owavenum(noc),awavenum(natm), stat=allocstat)
645  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
646 
647  j=0
648  DO i=1,nbatm
649  IF (ams(i,1)==1) THEN
650  awavenum(j+1)%typ='A'
651  awavenum(j+2)%typ='K'
652  awavenum(j+3)%typ='L'
653 
654  awavenum(j+1)%P=ams(i,2)
655  awavenum(j+2)%M=ams(i,1)
656  awavenum(j+2)%P=ams(i,2)
657  awavenum(j+3)%H=ams(i,1)
658  awavenum(j+3)%P=ams(i,2)
659 
660  awavenum(j+1)%Ny=REAL(ams(i,2))
661  awavenum(j+2)%Nx=REAL(ams(i,1))
662  awavenum(j+2)%Ny=REAL(ams(i,2))
663  awavenum(j+3)%Nx=REAL(ams(i,1))
664  awavenum(j+3)%Ny=REAL(ams(i,2))
665 
666  j=j+3
667  ELSE
668  awavenum(j+1)%typ='K'
669  awavenum(j+2)%typ='L'
670 
671  awavenum(j+1)%M=ams(i,1)
672  awavenum(j+1)%P=ams(i,2)
673  awavenum(j+2)%H=ams(i,1)
674  awavenum(j+2)%P=ams(i,2)
675 
676  awavenum(j+1)%Nx=REAL(ams(i,1))
677  awavenum(j+1)%Ny=REAL(ams(i,2))
678  awavenum(j+2)%Nx=REAL(ams(i,1))
679  awavenum(j+2)%Ny=REAL(ams(i,2))
680 
681  j=j+2
682 
683  ENDIF
684  ENDDO
685 
686  DO i=1,noc
687  owavenum(i)%H=oms(i,1)
688  owavenum(i)%P=oms(i,2)
689 
690  owavenum(i)%Nx=oms(i,1)/2.d0
691  owavenum(i)%Ny=oms(i,2)
692 
693  ENDDO
694 
695  ! Computation of the atmospheric inner products tensors
696 
697  CALL calculate_a
698  CALL calculate_g
699  CALL calculate_s
700  CALL calculate_b
701  CALL calculate_c_atm
702 
703  ! Computation of the oceanic inner products tensors
704 
705  CALL calculate_m
706  CALL calculate_n
707  CALL calculate_o
708  CALL calculate_c_oc
709  CALL calculate_w
710  CALL calculate_k
711 
712  ! A last atmospheric one that needs ocean%M
713 
714  CALL calculate_d
715 
716 
717 
718  END SUBROUTINE init_inprod
719 
720  !> Deallocation of the inner products
721  SUBROUTINE deallocate_inprod
722  INTEGER :: AllocStat
723 
724  ! Deallocation of atmospheric inprod
725  allocstat=0
726  IF (ALLOCATED(atmos%a)) DEALLOCATE(atmos%a, stat=allocstat)
727  IF (allocstat /= 0) stop "*** Problem to deallocate ! ***"
728 
729  allocstat=0
730  IF (ALLOCATED(atmos%c)) DEALLOCATE(atmos%c, stat=allocstat)
731  IF (allocstat /= 0) stop "*** Problem to deallocate ! ***"
732 
733  allocstat=0
734  IF (ALLOCATED(atmos%d)) DEALLOCATE(atmos%d, stat=allocstat)
735  IF (allocstat /= 0) stop "*** Problem to deallocate ! ***"
736 
737  allocstat=0
738  IF (ALLOCATED(atmos%s)) DEALLOCATE(atmos%s, stat=allocstat)
739  IF (allocstat /= 0) stop "*** Problem to deallocate ! ***"
740 
741  allocstat=0
742  IF (ALLOCATED(atmos%g)) DEALLOCATE(atmos%g, stat=allocstat)
743  IF (allocstat /= 0) stop "*** Problem to deallocate ! ***"
744 
745  allocstat=0
746  IF (ALLOCATED(atmos%b)) DEALLOCATE(atmos%b, stat=allocstat)
747  IF (allocstat /= 0) stop "*** Problem to deallocate ! ***"
748 
749  ! Deallocation of oceanic inprod
750  allocstat=0
751  IF (ALLOCATED(ocean%K)) DEALLOCATE(ocean%K, stat=allocstat)
752  IF (allocstat /= 0) stop "*** Problem to deallocate ! ***"
753 
754  allocstat=0
755  IF (ALLOCATED(ocean%M)) DEALLOCATE(ocean%M, stat=allocstat)
756  IF (allocstat /= 0) stop "*** Problem to deallocate ! ***"
757 
758  allocstat=0
759  IF (ALLOCATED(ocean%N)) DEALLOCATE(ocean%N, stat=allocstat)
760  IF (allocstat /= 0) stop "*** Problem to deallocate ! ***"
761 
762  allocstat=0
763  IF (ALLOCATED(ocean%W)) DEALLOCATE(ocean%W, stat=allocstat)
764  IF (allocstat /= 0) stop "*** Problem to deallocate ! ***"
765 
766  allocstat=0
767  IF (ALLOCATED(ocean%O)) DEALLOCATE(ocean%O, stat=allocstat)
768  IF (allocstat /= 0) stop "*** Problem to deallocate ! ***"
769 
770  allocstat=0
771  IF (ALLOCATED(ocean%C)) DEALLOCATE(ocean%C, stat=allocstat)
772  IF (allocstat /= 0) stop "*** Problem to deallocate ! ***"
773  END SUBROUTINE deallocate_inprod
774 
775 END MODULE inprod_analytic
776 
integer noc
Number of oceanic basis functions.
Definition: params.f90:78
Atmospheric bloc specification type.
subroutine, public init_inprod
Initialisation of the inner product.
integer nbatm
Number of oceanic blocks.
Definition: params.f90:76
Inner products between the truncated set of basis functions for the ocean and atmosphere streamfuncti...
Statistics accumulators.
Definition: stat.f90:14
real(kind=8) function b1(Pi, Pj, Pk)
Cehelsky & Tung Helper functions.
real(kind=8) function s4(Pj, Pk, Hj, Hk)
Cehelsky & Tung Helper functions.
subroutine calculate_g
Temperature advection terms (atmospheric)
real(kind=8) pi
Definition: params.f90:48
subroutine calculate_m
Forcing of the ocean fields on the ocean.
real(kind=8) function s2(Pj, Pk, Mj, Hk)
Cehelsky & Tung Helper functions.
real(kind=8) function b2(Pi, Pj, Pk)
Cehelsky & Tung Helper functions.
subroutine calculate_a
Eigenvalues of the Laplacian (atmospheric)
real(kind=8) function delta(r)
Integer Dirac delta function.
subroutine calculate_k
Forcing of the atmosphere on the ocean.
subroutine, public deallocate_inprod
Deallocation of the inner products.
subroutine calculate_c_atm
Beta term for the atmosphere.
type(ocean_wavenum), dimension(:), allocatable, public owavenum
Oceanic blocs specification.
subroutine calculate_s
Forcing (thermal) of the ocean on the atmosphere.
Type holding the oceanic inner products tensors.
type(atm_wavenum), dimension(:), allocatable, public awavenum
Atmospheric blocs specification.
Oceanic bloc specification type.
real(kind=8) function s1(Pj, Pk, Mj, Hk)
Cehelsky & Tung Helper functions.
real(kind=8) function s3(Pj, Pk, Hj, Hk)
Cehelsky & Tung Helper functions.
integer natm
Number of atmospheric basis functions.
Definition: params.f90:77
type(atm_tensors), public atmos
Atmospheric tensors.
subroutine calculate_o
Temperature advection term (passive scalar)
subroutine calculate_n
Beta term for the ocean.
integer nboc
Number of atmospheric blocks.
Definition: params.f90:75
real(kind=8) function flambda(r)
"Odd or even" function
The model parameters module.
Definition: params.f90:18
subroutine calculate_w
Short-wave radiative forcing of the ocean.
Type holding the atmospheric inner products tensors.
integer, dimension(:,:), allocatable ams
Atmospheric mode selection array.
Definition: params.f90:81
subroutine calculate_b
Streamfunction advection terms (atmospheric)
type(ocean_tensors), public ocean
Oceanic tensors.
integer, dimension(:,:), allocatable oms
Ocean mode selection array.
Definition: params.f90:80
subroutine calculate_c_oc
Streamfunction advection terms (oceanic)
subroutine calculate_d
Forcing of the ocean on the atmosphere.
real(kind=8) n
- Aspect ratio
Definition: params.f90:24