New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
flx_bulk_monthly.h90 in trunk/NEMO/OPA_SRC/SBC – NEMO

source: trunk/NEMO/OPA_SRC/SBC/flx_bulk_monthly.h90 @ 699

Last change on this file since 699 was 699, checked in by smasson, 17 years ago

insert revision Id

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 9.2 KB
Line 
1   !!----------------------------------------------------------------------
2   !!                    ***  flx_blulk_monthly.h90  ***
3   !!----------------------------------------------------------------------
4   !!   flx     : update surface thermohaline fluxes using bulk formulae
5   !!             and fields read in a NetCDF file
6   !!----------------------------------------------------------------------
7   !! * Modules used     C A U T I O N  already defined in flxmod.F90
8
9   !! * Module variables
10   
11   INTEGER ::          &
12      ji, jj,          &  ! loop indices
13      numflx,          &  ! logical unit for surface fluxes data
14      nflx1 , nflx2,   &  !  first and second record used
15      nflx11, nflx12      ! ???
16
17   INTEGER, PARAMETER :: jpf    =  7                   
18   REAL(wp), DIMENSION(jpi,jpj,2,jpf) ::   &
19      flxdta              ! 2 consecutive set of CLIO monthly fluxes
20   !!----------------------------------------------------------------------
21   !!   OPA 9.0 , LOCEAN-IPSL (2005)
22   !! $Id$
23   !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt
24   !!----------------------------------------------------------------------
25
26CONTAINS
27
28   SUBROUTINE flx( kt )
29      !!---------------------------------------------------------------------
30      !!                     ***  ROUTINE flx  ***
31      !!                   
32      !! ** Purpose :   provide the thermohaline fluxes (heat and freshwater)
33      !!      to the ocean at each time step.
34      !!
35      !! ** Method  :   Read monthly climatological fluxes in a NetCDF file
36      !!          the net downward radiative flux qsr      1 (watt/m2)
37      !!          the net downward heat flux      q        2 (watt/m2)
38      !!          the net upward water            emp      3 (mm/month)
39      !!              (evaporation - precipitation)
40      !!          the climatological ice cover    rclice   4 (0 or 1)
41      !!
42      !!     Qsr and q is obtained from Esbensen-Kushnir data (opal file) with
43      !!   some corrections :
44      !!          - Data are extended over the polar area and for the net heat
45      !!            flux, values are put at 200 w/m2 on the ice regions
46      !!          - Red sea and Mediterranean values are imposed.
47      !!
48      !!     emp is the Oberhuber climatology with a function of Levitus
49      !!   salinity
50      !!
51      !!     rclice is an handmade climalological ice cover on the polar
52      !!   regions.
53      !!
54      !!     runoff is an handmade climalological runoff.
55      !!
56      !! caution : now, in the opa global model, the net upward water flux is
57      !! -------   with mm/day unit.
58      !!
59      !! History :
60      !!        !  91-03  (O. Marti and Ph Dandin)  Original code
61      !!        !  92-07  (M. Imbard)
62      !!        !  96-11  (E. Guilyardi)  Daily AGCM input files
63      !!        !  99-11  (M. Imbard)  NetCDF FORMAT with io-ipsl
64      !!        !  00-10  (J.-P. Boulanger)  adjusted for reading any
65      !!                         daily wind stress data including a climatology
66      !!        !  01-09  (A. Lazar and C. Levy)  Daily NetCDF by default
67      !!   8.5  !  02-09  (G. Madec)  F90: Free form and module
68      !!----------------------------------------------------------------------
69      !! * modules used
70      USE iom
71      USE blk_oce         ! bulk variable
72      USE bulk            ! bulk module
73
74      !! * arguments
75      INTEGER, INTENT( in  ) ::   kt   ! ocean time step
76
77      !! * Local declarations
78      INTEGER ::   jm            ! dummy loop indices
79      INTEGER ::   &
80         imois, imois2,       &  ! temporary integers
81         i15  , iman             !    "          "
82      REAL(wp) ::   &
83         zxy   , zdtt  ,      &  !    "         "
84         zttbt , zttat ,      &  !    "         "
85         zdtts6                  !    "         "
86      !!---------------------------------------------------------------------
87
88      ! Initialization
89      ! --------------
90
91      i15 = INT( 2 * FLOAT( nday ) / ( FLOAT( nobis(nmonth) ) + 0.5 ) )
92      iman  = INT( raamo )
93      imois = nmonth + i15 - 1
94      IF( imois == 0 ) imois = iman
95      imois2 = nmonth
96
97      ! 1. first call kt=nit000
98      ! -----------------------
99
100      IF( kt == nit000 ) THEN
101         ! initializations
102         nflx1  = 0
103         nflx11 = 0
104         ! open the file
105         IF(lwp) THEN
106            WRITE(numout,*) ' '
107            WRITE(numout,*) ' **** Routine flx_bulk_monthly.h90'
108            WRITE(numout,*) ' '
109            WRITE(numout,*) ' global CLIO flx monthly fields'
110         ENDIF
111         CALL iom_open ( 'flx.nc', numflx )
112       
113         ! temperature, spline initialization, we read the first record
114         CALL iom_get ( numflx, jpdom_data, 'socliot1', flxdta(:,:,1,5), 1 )
115
116      ENDIF
117
118
119      ! Read monthly file
120      ! ----------------
121
122      IF( kt == nit000 .OR. imois /= nflx1 ) THEN
123
124         ! Calendar computation
125
126         ! nflx1 number of the first file record used in the simulation
127         ! nflx2 number of the last  file record
128
129         nflx1 = imois
130         nflx2 = nflx1+1
131         nflx1 = MOD( nflx1, iman )
132         nflx2 = MOD( nflx2, iman )
133         IF( nflx1 == 0 )   nflx1 = iman
134         IF( nflx2 == 0 )   nflx2 = iman
135         IF(lwp) WRITE(numout,*) 'first record file used nflx1 ',nflx1
136         IF(lwp) WRITE(numout,*) 'last  record file used nflx2 ',nflx2
137         
138         ! Read monthly fluxes data
139
140         ! humidity
141         CALL iom_get ( numflx, jpdom_data, 'socliohu', flxdta(:,:,1,1), nflx1 )
142         CALL iom_get ( numflx, jpdom_data, 'socliohu', flxdta(:,:,2,1), nflx2 )
143         ! 10m wind module
144         CALL iom_get ( numflx, jpdom_data, 'socliowi', flxdta(:,:,1,2), nflx1 )
145         CALL iom_get ( numflx, jpdom_data, 'socliowi', flxdta(:,:,2,2), nflx2 )
146         ! cloud cover
147         CALL iom_get ( numflx, jpdom_data, 'socliocl', flxdta(:,:,1,3), nflx1 )
148         CALL iom_get ( numflx, jpdom_data, 'socliocl', flxdta(:,:,2,3), nflx2 )
149         ! precipitations
150         CALL iom_get ( numflx, jpdom_data, 'socliopl', flxdta(:,:,1,4), nflx1 )
151         CALL iom_get ( numflx, jpdom_data, 'socliopl', flxdta(:,:,2,4), nflx2 )
152         
153         IF(lwp .AND. nitend-nit000 <= 100 ) THEN
154            WRITE(numout,*)
155            WRITE(numout,*) ' read clio flx ok'
156            WRITE(numout,*)
157            DO jm = 1, 4
158               WRITE(numout,*)
159               WRITE(numout,*) 'Clio mounth: ',nflx1,'  field: ',jm,' multiply by ',0.1
160               CALL prihre( flxdta(:,:,1,jm),jpi,jpj,1,jpi,20,1,jpj,10,.1,numout )
161            END DO
162         ENDIF
163
164      ENDIF
165
166      IF( kt == nit000 .OR. imois2 /= nflx11 ) THEN
167
168         ! calendar computation
169         
170         ! nflx1 number of the first file record used in the simulation
171         ! nflx2 number of the last  file record
172         
173         nflx11 = imois2
174         nflx12 = nflx11 + 1
175         nflx11 = MOD( nflx11, iman )
176         nflx12 = MOD( nflx12, iman )
177         IF( nflx11 == 0 )   nflx11 = iman
178         IF( nflx12 == 0 )   nflx12 = iman
179         IF(lwp) WRITE(numout,*) 'first record file used nflx11 ',nflx11
180         IF(lwp) WRITE(numout,*) 'last  record file used nflx12 ',nflx12
181         
182         ! Read monthly fluxes data Esbensen Kushnir
183         
184         ! air temperature
185         ! Utilisation d'un spline, on lit le champ a mois=nflx1 et nflx2
186         CALL iom_get (numflx,jpdom_data,'socliot1',flxdta(:,:,1,6),nflx11)
187         CALL iom_get (numflx,jpdom_data,'socliot1',flxdta(:,:,2,6),nflx12)
188         ! air temperature derivative (to reconstruct a daily field)
189         CALL iom_get (numflx,jpdom_data,'socliot2',flxdta(:,:,1,7),nflx11)
190         CALL iom_get (numflx,jpdom_data,'socliot2',flxdta(:,:,2,7),nflx12)
191         
192         IF(lwp) THEN
193            WRITE(numout,*)
194            WRITE(numout,*) ' read CLIO flx ok'
195            WRITE(numout,*)
196            DO jm = 6, jpf
197               WRITE(numout,*) 'jpf =  ', jpf !C a u t i o n : information need for SX5NEC compilo bug
198               WRITE(numout,*) 'Clio mounth: ',nflx11,'  field: ',jm,' multiply by ',0.1
199               CALL prihre(flxdta(:,:,1,jm),jpi,jpj,1,jpi,20,1,jpj,10,.1,numout)
200               WRITE(numout,*)
201            END DO
202         ENDIF
203
204      ENDIF
205
206
207      ! 3. at every time step interpolation of fluxes
208      ! ---------------------------------------------
209
210      zxy = FLOAT( nday ) / FLOAT( nobis(nflx1) ) + 0.5 - i15
211
212      zdtt = raajj / raamo
213      zttbt = (REAL(nday) - 1.)/(nobis(nmonth) - 1.)
214
215      zttat = 1. - zttbt
216      zdtts6 = zdtt/6.
217
218      hatm(:,:) = ( (1.-zxy) * flxdta(:,:,1,1) + zxy * flxdta(:,:,2,1) )
219      vatm(:,:) = ( (1.-zxy) * flxdta(:,:,1,2) + zxy * flxdta(:,:,2,2) )
220      catm(:,:) = ( (1.-zxy )* flxdta(:,:,1,3) + zxy * flxdta(:,:,2,3) )
221      watm(:,:) = ( (1.-zxy) * flxdta(:,:,1,4) + zxy * flxdta(:,:,2,4) )
222      tatm(:,:) = ( flxdta(:,:,2,6) - flxdta(:,:,1,6) )/zdtt   &
223                - ((3. * zttat * zttat - 1.) * flxdta(:,:,1,7)   &
224                - ( 3. * zttbt * zttbt - 1.) * flxdta(:,:,2,7) ) * zdtts6   &
225                + flxdta(:,:,1,5)
226 
227      CALL blk( kt )                ! bulk formulea fluxes
228
229      ! ------------------- !
230      ! Last call kt=nitend !
231      ! ------------------- !
232
233      ! Closing of the numflx file (required in mpp)
234      IF( kt == nitend ) CALL iom_close (numflx)
235
236   END SUBROUTINE flx
Note: See TracBrowser for help on using the repository browser.