source: codes/icosagcm/devel/src/dynamics/compute_caldyn.f90 @ 1027

Last change on this file since 1027 was 1027, checked in by dubos, 4 years ago

devel : towards conformity to F2008 standard

File size: 3.6 KB
Line 
1MODULE compute_caldyn_mod
2  USE prec, ONLY : rstd
3  IMPLICIT NONE
4  SAVE
5 
6  ! fake array dimensions, for interfaces
7  INTEGER, PARAMETER :: iim_jjm_i=1, iim_jjm_u=1, iim_jjm_v =1, llm_=1, llm1=1, nqdyn_=1
8
9  INTERFACE
10     
11    SUBROUTINE comp_pvort_only(u,rhodz,qu,qv, hv)
12      IMPORT
13      REAL(rstd), INTENT(IN)    :: u(:,:)
14      REAL(rstd), INTENT(INOUT) :: rhodz(:,:)
15      REAL(rstd), INTENT(OUT)   :: qu(:,:)
16      REAL(rstd), INTENT(OUT)   :: qv(:,:)
17      REAL(rstd), INTENT(OUT)   :: hv(:,:)
18    END SUBROUTINE comp_pvort_only
19
20    SUBROUTINE comp_theta(mass_col,theta_rhodz, rhodz,theta)
21      IMPORT
22      REAL(rstd), INTENT(IN)    :: mass_col(:)
23      REAL(rstd), INTENT(IN)    :: theta_rhodz(:,:,:)
24      REAL(rstd), INTENT(INOUT) :: rhodz(:,:)
25      REAL(rstd), INTENT(OUT)   :: theta(:,:,:)
26    END SUBROUTINE comp_theta
27
28    SUBROUTINE comp_geopot(rhodz,theta, ps,pk,geopot) 
29      IMPORT
30      REAL(rstd), INTENT(IN)    :: rhodz(:,:)    ! rho*dz = mass per unit surface in each full model level
31      REAL(rstd), INTENT(IN)    :: theta(:,:,:)  ! active scalars : theta/entropy, moisture, ...
32      REAL(rstd), INTENT(INOUT) :: ps(:)         ! surface pressure
33      REAL(rstd), INTENT(OUT)   :: pk(:,:)       ! Exner function (compressible) /Lagrange multiplier (Boussinesq)
34      REAL(rstd), INTENT(INOUT) :: geopot(:,:)   ! geopotential
35    END SUBROUTINE comp_geopot
36
37    SUBROUTINE comp_caldyn_fast(tau,theta,geopot, pk,berni,du,u)
38      IMPORT
39      REAL(rstd), INTENT(IN)    :: tau           ! "solve" u-tau*du/dt = rhs
40      REAL(rstd), INTENT(IN)    :: theta(:,:,:)
41      REAL(rstd), INTENT(IN)    :: geopot(:,:)
42      REAL(rstd), INTENT(INOUT) :: pk(:,:)
43      REAL(rstd), INTENT(INOUT) :: berni(:,:)    ! partial Bernoulli function
44      REAL(rstd), INTENT(INOUT) :: du(:,:)
45      REAL(rstd), INTENT(INOUT) :: u(:,:)        ! INOUT if tau>0
46    END SUBROUTINE comp_caldyn_fast
47
48    SUBROUTINE comp_caldyn_slow_hydro(zero, u,rhodz,hv,Kv, berni, hflux,du)
49      IMPORT
50      LOGICAL,    INTENT(IN)    :: zero
51      REAL(rstd), INTENT(IN)    :: u(:,:)      ! prognostic "velocity"
52      REAL(rstd), INTENT(IN)    :: rhodz(:,:)
53      REAL(rstd), INTENT(IN)    :: hv(:,:)     ! height/mass averaged to vertices
54      REAL(rstd), INTENT(IN)    :: Kv(:,:)     ! kinetic energy at vertices
55      REAL(rstd), INTENT(OUT)   :: berni(:,:)  ! Bernoulli function
56      REAL(rstd), INTENT(OUT)   :: hflux(:,:)  ! hflux in kg/s
57      REAL(rstd), INTENT(INOUT) :: du(:,:)
58    END SUBROUTINE comp_caldyn_slow_hydro
59   
60    SUBROUTINE comp_caldyn_coriolis(hflux,theta,qu, Ftheta, convm,dtheta_rhodz,du)
61      IMPORT
62      REAL(rstd), INTENT(IN)    :: hflux(:,:)   ! hflux in kg/s
63      REAL(rstd), INTENT(IN)    :: theta(:,:,:) ! active scalars
64      REAL(rstd), INTENT(IN)    :: qu(:,:)
65      REAL(rstd), INTENT(OUT)   :: Ftheta(:,:)  ! potential temperature flux
66      REAL(rstd), INTENT(OUT)   :: convm(:,:)   ! mass flux convergence
67      REAL(rstd), INTENT(OUT)   :: dtheta_rhodz(:,:,:)
68      REAL(rstd), INTENT(INOUT) :: du(:,:)
69    END SUBROUTINE comp_caldyn_coriolis
70
71  END INTERFACE
72
73  PROCEDURE(comp_pvort_only),        POINTER :: compute_pvort_only        => NULL()
74  PROCEDURE(comp_theta),             POINTER :: compute_theta             => NULL()
75  PROCEDURE(comp_geopot),            POINTER :: compute_geopot            => NULL()
76  PROCEDURE(comp_caldyn_fast),       POINTER :: compute_caldyn_fast       => NULL()
77  PROCEDURE(comp_caldyn_slow_hydro), POINTER :: compute_caldyn_slow_hydro => NULL()
78  PROCEDURE(comp_caldyn_coriolis),   POINTER :: compute_caldyn_coriolis   => NULL()
79
80END MODULE compute_caldyn_mod
Note: See TracBrowser for help on using the repository browser.