ROMS
Loading...
Searching...
No Matches
ana_tclima.h
Go to the documentation of this file.
1!!
2 SUBROUTINE ana_tclima (ng, tile, model)
3!
4!! git $Id$
5!!======================================================================
6!! Copyright (c) 2002-2025 The ROMS Group !
7!! Licensed under a MIT/X style license !
8!! See License_ROMS.md !
9!=======================================================================
10! !
11! This routine sets analytical tracer climatology fields. !
12! !
13!=======================================================================
14!
15 USE mod_param
16 USE mod_ncparam
17!
18! Imported variable declarations.
19!
20 integer, intent(in) :: ng, tile, model
21!
22! Local variable declarations.
23!
24 character (len=*), parameter :: MyFile = &
25 & __FILE__
26!
27#include "tile.h"
28!
29 CALL ana_tclima_tile (ng, tile, model, &
30 & lbi, ubi, lbj, ubj, &
31 & imins, imaxs, jmins, jmaxs)
32!
33! Set analytical header file name used.
34!
35#ifdef DISTRIBUTE
36 IF (lanafile) THEN
37#else
38 IF (lanafile.and.(tile.eq.0)) THEN
39#endif
40 ananame(33)=myfile
41 END IF
42!
43 RETURN
44 END SUBROUTINE ana_tclima
45!
46!***********************************************************************
47 SUBROUTINE ana_tclima_tile (ng, tile, model, &
48 & LBi, UBi, LBj, UBj, &
49 & IminS, ImaxS, JminS, JmaxS)
50!***********************************************************************
51!
52 USE mod_param
53 USE mod_clima
54 USE mod_grid
55 USE mod_scalars
56!
58#ifdef DISTRIBUTE
60#endif
61!
62! Imported variable declarations.
63!
64 integer, intent(in) :: ng, tile, model
65 integer, intent(in) :: LBi, UBi, LBj, UBj
66 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
67!
68! Local variable declarations.
69!
70 integer :: i, itrc, j, k
71!
72 real(r8) :: val1, val2, val3, val4
73
74#include "set_bounds.h"
75!
76!-----------------------------------------------------------------------
77! Set tracer climatology.
78!-----------------------------------------------------------------------
79!
80 IF (any(ltracerclm(:,ng)).or.any(lnudgetclm(:,ng))) THEN
81#if defined DOUBLE_GYRE
82 val1=(44.69_r8/39.382_r8)**2
83 val2=val1*(rho0*100.0_r8/g)* &
84 & (5.0e-5_r8/((42.689_r8/44.69_r8)**2))
85 DO k=1,n(ng)
86 DO j=jstrt,jendt
87 DO i=istrt,iendt
88 val3=t0(ng)+val2*exp(grid(ng)%z_r(i,j,k)/100.0_r8)* &
89 & (10.0_r8-0.4_r8*tanh(grid(ng)%z_r(i,j,k)/100.0_r8))
90 val4=grid(ng)%yr(i,j)/el(ng)
91 clima(ng)%tclm(i,j,k,itemp)=val3-3.0_r8*val4
92# ifdef SALINITY
93 clima(ng)%tclm(i,j,k,isalt)=34.5_r8- &
94 & 0.001_r8*grid(ng)%z_r(i,j,k)- &
95 & val4
96# endif
97 END DO
98 END DO
99 END DO
100#else
101 DO k=1,n(ng)
102 DO j=jstrt,jendt
103 DO i=istrt,iendt
104 clima(ng)%tclm(i,j,k,itemp)=???
105# ifdef SALINITY
106 clima(ng)%tclm(i,j,k,isalt)=???
107# endif
108 END DO
109 END DO
110 END DO
111#endif
112!
113! Exchange boundary data.
114!
115 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
116 DO itrc=1,ntclm(ng)
117 CALL exchange_r3d_tile (ng, tile, &
118 & lbi, ubi, lbj, ubj, 1, n(ng), &
119 & clima(ng)%tclm(:,:,:,itrc))
120 END DO
121 END IF
122
123#ifdef DISTRIBUTE
124 CALL mp_exchange4d (ng, tile, model, 1, &
125 & lbi, ubi, lbj, ubj, 1, n(ng), 1, ntclm(ng), &
126 & nghostpoints, &
127 & ewperiodic(ng), nsperiodic(ng), &
128 & clima(ng) % tclm)
129#endif
130 END IF
131!
132 RETURN
133 END SUBROUTINE ana_tclima_tile
subroutine ana_tclima(ng, tile, model)
Definition ana_tclima.h:3
subroutine ana_tclima_tile(ng, tile, model, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs)
Definition ana_tclima.h:50
subroutine exchange_r3d_tile(ng, tile, lbi, ubi, lbj, ubj, lbk, ubk, a)
type(t_clima), dimension(:), allocatable clima
Definition mod_clima.F:153
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
logical lanafile
character(len=256), dimension(39) ananame
integer, dimension(:), allocatable ntclm
Definition mod_param.F:494
integer, dimension(:), allocatable n
Definition mod_param.F:479
integer nghostpoints
Definition mod_param.F:710
real(r8), dimension(:), allocatable t0
real(r8), dimension(:), allocatable el
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
integer isalt
integer itemp
real(dp) g
logical, dimension(:,:), allocatable ltracerclm
real(dp) rho0
logical, dimension(:,:), allocatable lnudgetclm
subroutine mp_exchange4d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, lbk, ubk, lbt, ubt, nghost, ew_periodic, ns_periodic, a, b, c)