A Modular Arbitrary-Order Ocean-Atmosphere Model -- Stochastic implementation
Data Types | Functions/Subroutines | Variables
inprod_analytic Module Reference

Inner products between the truncated set of basis functions for the ocean and atmosphere streamfunction fields. These are partly calculated using the analytical expressions from Cehelsky, P., & Tung, K. K. : Theories of multiple equilibria and weather regimes-A critical reexamination. Part II: Baroclinic two-layer models. Journal of the atmospheric sciences, 44(21), 3282-3303, 1987. More...

Data Types

type  atm_tensors
 Type holding the atmospheric inner products tensors. More...
 
type  atm_wavenum
 Atmospheric bloc specification type. More...
 
type  ocean_tensors
 Type holding the oceanic inner products tensors. More...
 
type  ocean_wavenum
 Oceanic bloc specification type. More...
 

Functions/Subroutines

real(kind=8) function b1 (Pi, Pj, Pk)
 Cehelsky & Tung Helper functions. More...
 
real(kind=8) function b2 (Pi, Pj, Pk)
 Cehelsky & Tung Helper functions. More...
 
real(kind=8) function delta (r)
 Integer Dirac delta function. More...
 
real(kind=8) function flambda (r)
 "Odd or even" function More...
 
real(kind=8) function s1 (Pj, Pk, Mj, Hk)
 Cehelsky & Tung Helper functions. More...
 
real(kind=8) function s2 (Pj, Pk, Mj, Hk)
 Cehelsky & Tung Helper functions. More...
 
real(kind=8) function s3 (Pj, Pk, Hj, Hk)
 Cehelsky & Tung Helper functions. More...
 
real(kind=8) function s4 (Pj, Pk, Hj, Hk)
 Cehelsky & Tung Helper functions. More...
 
subroutine calculate_a
 Eigenvalues of the Laplacian (atmospheric) More...
 
subroutine calculate_b
 Streamfunction advection terms (atmospheric) More...
 
subroutine calculate_c_atm
 Beta term for the atmosphere. More...
 
subroutine calculate_d
 Forcing of the ocean on the atmosphere. More...
 
subroutine calculate_g
 Temperature advection terms (atmospheric) More...
 
subroutine calculate_s
 Forcing (thermal) of the ocean on the atmosphere. More...
 
subroutine calculate_k
 Forcing of the atmosphere on the ocean. More...
 
subroutine calculate_m
 Forcing of the ocean fields on the ocean. More...
 
subroutine calculate_n
 Beta term for the ocean. More...
 
subroutine calculate_o
 Temperature advection term (passive scalar) More...
 
subroutine calculate_c_oc
 Streamfunction advection terms (oceanic) More...
 
subroutine calculate_w
 Short-wave radiative forcing of the ocean. More...
 
subroutine, public init_inprod
 Initialisation of the inner product. More...
 
subroutine, public deallocate_inprod
 Deallocation of the inner products. More...
 

Variables

type(atm_wavenum), dimension(:), allocatable, public awavenum
 Atmospheric blocs specification. More...
 
type(ocean_wavenum), dimension(:), allocatable, public owavenum
 Oceanic blocs specification. More...
 
type(atm_tensors), public atmos
 Atmospheric tensors. More...
 
type(ocean_tensors), public ocean
 Oceanic tensors. More...
 

Detailed Description

Inner products between the truncated set of basis functions for the ocean and atmosphere streamfunction fields. These are partly calculated using the analytical expressions from Cehelsky, P., & Tung, K. K. : Theories of multiple equilibria and weather regimes-A critical reexamination. Part II: Baroclinic two-layer models. Journal of the atmospheric sciences, 44(21), 3282-3303, 1987.

Remarks
Generated Fortran90/95 code from inprod_analytic.lua

Function/Subroutine Documentation

real(kind=8) function inprod_analytic::b1 ( integer  Pi,
integer  Pj,
integer  Pk 
)
private

Cehelsky & Tung Helper functions.

Definition at line 91 of file inprod_analytic.f90.

91  INTEGER :: pi,pj,pk
92  b1 = (pk + pj) / REAL(pi)
real(kind=8) function inprod_analytic::b2 ( integer  Pi,
integer  Pj,
integer  Pk 
)
private

Cehelsky & Tung Helper functions.

Definition at line 97 of file inprod_analytic.f90.

97  INTEGER :: pi,pj,pk
98  b2 = (pk - pj) / REAL(pi)
subroutine inprod_analytic::calculate_a ( )
private

Eigenvalues of the Laplacian (atmospheric)

\( a_{i,j} = (F_i, \nabla^2 F_j)\) .

Definition at line 155 of file inprod_analytic.f90.

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
Statistics accumulators.
Definition: stat.f90:14
subroutine inprod_analytic::calculate_b ( )
private

Streamfunction advection terms (atmospheric)

\( b_{i,j,k} = (F_i, J(F_j, \nabla^2 F_k))\) .

Remarks
Atmospheric g and a tensors must be computed before calling this routine

Definition at line 182 of file inprod_analytic.f90.

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
Statistics accumulators.
Definition: stat.f90:14
subroutine inprod_analytic::calculate_c_atm ( )
private

Beta term for the atmosphere.

\( c_{i,j} = (F_i, \partial_x F_j)\) .

Remarks
Strict function !! Only accepts KL type. For any other combination, it will not calculate anything

Definition at line 216 of file inprod_analytic.f90.

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
Statistics accumulators.
Definition: stat.f90:14
subroutine inprod_analytic::calculate_c_oc ( )
private

Streamfunction advection terms (oceanic)

\( C_{i,j,k} = (\eta_i, J(\eta_j,\nabla^2 \eta_k))\) .

Remarks
Requires O_{i,j,k} and M_{i,j} to be calculated beforehand.

Definition at line 568 of file inprod_analytic.f90.

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
Statistics accumulators.
Definition: stat.f90:14
subroutine inprod_analytic::calculate_d ( )
private

Forcing of the ocean on the atmosphere.

\( d_{i,j} = (F_i, \nabla^2 \eta_j)\) .

Remarks
Atmospheric s tensor and oceanic M tensor must be computed before calling this routine !

Definition at line 255 of file inprod_analytic.f90.

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
Statistics accumulators.
Definition: stat.f90:14
subroutine inprod_analytic::calculate_g ( )
private

Temperature advection terms (atmospheric)

\( g_{i,j,k} = (F_i, J(F_j, F_k))\) .

Definition at line 288 of file inprod_analytic.f90.

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 
Statistics accumulators.
Definition: stat.f90:14
subroutine inprod_analytic::calculate_k ( )
private

Forcing of the atmosphere on the ocean.

\( K_{i,j} = (\eta_i, \nabla^2 F_j)\) .

Remarks
atmospheric a and s tensors must be computed before calling this function !

Definition at line 434 of file inprod_analytic.f90.

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
Statistics accumulators.
Definition: stat.f90:14
subroutine inprod_analytic::calculate_m ( )
private

Forcing of the ocean fields on the ocean.

\( M_{i,j} = (eta_i, \nabla^2 \eta_j)\) .

Definition at line 464 of file inprod_analytic.f90.

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
Statistics accumulators.
Definition: stat.f90:14
subroutine inprod_analytic::calculate_n ( )
private

Beta term for the ocean.

\( N_{i,j} = (\eta_i, \partial_x \eta_j) \).

Definition at line 487 of file inprod_analytic.f90.

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
Statistics accumulators.
Definition: stat.f90:14
subroutine inprod_analytic::calculate_o ( )
private

Temperature advection term (passive scalar)

\( O_{i,j,k} = (\eta_i, J(\eta_j, \eta_k))\) .

Definition at line 516 of file inprod_analytic.f90.

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
Statistics accumulators.
Definition: stat.f90:14
subroutine inprod_analytic::calculate_s ( )
private

Forcing (thermal) of the ocean on the atmosphere.

\( s_{i,j} = (F_i, \eta_j)\) .

Definition at line 380 of file inprod_analytic.f90.

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
Statistics accumulators.
Definition: stat.f90:14
subroutine inprod_analytic::calculate_w ( )
private

Short-wave radiative forcing of the ocean.

\( W_{i,j} = (\eta_i, F_j)\) .

Remarks
atmospheric s tensor must be computed before calling this function !

Definition at line 605 of file inprod_analytic.f90.

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
Statistics accumulators.
Definition: stat.f90:14
subroutine, public inprod_analytic::deallocate_inprod ( )

Deallocation of the inner products.

Definition at line 722 of file inprod_analytic.f90.

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 ! ***"
Statistics accumulators.
Definition: stat.f90:14
real(kind=8) function inprod_analytic::delta ( integer  r)
private

Integer Dirac delta function.

Definition at line 103 of file inprod_analytic.f90.

103  INTEGER :: r
104  IF (r==0) THEN
105  delta = 1.d0
106  ELSE
107  delta = 0.d0
108  ENDIF
real(kind=8) function inprod_analytic::flambda ( integer  r)
private

"Odd or even" function

Definition at line 113 of file inprod_analytic.f90.

113  INTEGER :: r
114  IF (mod(r,2)==0) THEN
115  flambda = 0.d0
116  ELSE
117  flambda = 1.d0
118  ENDIF
subroutine, public inprod_analytic::init_inprod ( )

Initialisation of the inner product.

Definition at line 639 of file inprod_analytic.f90.

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 
Statistics accumulators.
Definition: stat.f90:14
real(kind=8) function inprod_analytic::s1 ( integer  Pj,
integer  Pk,
integer  Mj,
integer  Hk 
)
private

Cehelsky & Tung Helper functions.

Definition at line 123 of file inprod_analytic.f90.

123  INTEGER :: pk,pj,mj,hk
124  s1 = -((pk * mj + pj * hk)) / 2.d0
real(kind=8) function inprod_analytic::s2 ( integer  Pj,
integer  Pk,
integer  Mj,
integer  Hk 
)
private

Cehelsky & Tung Helper functions.

Definition at line 129 of file inprod_analytic.f90.

129  INTEGER :: pk,pj,mj,hk
130  s2 = (pk * mj - pj * hk) / 2.d0
real(kind=8) function inprod_analytic::s3 ( integer  Pj,
integer  Pk,
integer  Hj,
integer  Hk 
)
private

Cehelsky & Tung Helper functions.

Definition at line 135 of file inprod_analytic.f90.

135  INTEGER :: pj,pk,hj,hk
136  s3 = (pk * hj + pj * hk) / 2.d0
real(kind=8) function inprod_analytic::s4 ( integer  Pj,
integer  Pk,
integer  Hj,
integer  Hk 
)
private

Cehelsky & Tung Helper functions.

Definition at line 141 of file inprod_analytic.f90.

141  INTEGER :: pj,pk,hj,hk
142  s4 = (pk * hj - pj * hk) / 2.d0

Variable Documentation

type(atm_tensors), public inprod_analytic::atmos

Atmospheric tensors.

Definition at line 69 of file inprod_analytic.f90.

69  TYPE(atm_tensors), PUBLIC :: atmos
type(atm_wavenum), dimension(:), allocatable, public inprod_analytic::awavenum

Atmospheric blocs specification.

Definition at line 64 of file inprod_analytic.f90.

64  TYPE(atm_wavenum), DIMENSION(:), ALLOCATABLE, PUBLIC :: awavenum
type(ocean_tensors), public inprod_analytic::ocean

Oceanic tensors.

Definition at line 71 of file inprod_analytic.f90.

71  TYPE(ocean_tensors), PUBLIC :: ocean
type(ocean_wavenum), dimension(:), allocatable, public inprod_analytic::owavenum

Oceanic blocs specification.

Definition at line 66 of file inprod_analytic.f90.

66  TYPE(ocean_wavenum), DIMENSION(:), ALLOCATABLE, PUBLIC :: owavenum