source: LMDZ5/branches/IPSLCM5A2.1_ISO/libf/phylmd/limit_read_mod.F90 @ 5360

Last change on this file since 5360 was 2311, checked in by Ehouarn Millour, 9 years ago

Further modifications to enforce physics/dynamics separation:

  • moved iniprint.h and misc_mod back to dyn3d_common, as these should only be used by dynamics.
  • created print_control_mod in the physics to store flags prt_level, lunout, debug to be local to physics (should be used rather than iniprint.h)
  • created abort_physic.F90 , which does the same job as abort_gcm() did, but should be used instead when in physics.
  • reactivated inifis (turned it into a module, inifis_mod.F90) to initialize physical constants and print_control_mod flags.

EM

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 13.1 KB
Line 
1!
2! $Header$
3!
4MODULE limit_read_mod
5!
6! This module reads the fichier "limit.nc" containing fields for surface forcing.
7!
8! Module subroutines :
9!  limit_read_frac    : call limit_read_tot and return the fractions
10!  limit_read_rug_alb : return rugosity and albedo, if coupled ocean call limit_read_tot first
11!  limit_read_sst     : return sea ice temperature   
12!  limit_read_tot     : read limit.nc and store the fields in local modules variables
13!
14  IMPLICIT NONE
15
16  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE, PRIVATE :: pctsrf
17!$OMP THREADPRIVATE(pctsrf)
18  REAL, ALLOCATABLE, DIMENSION(:),   SAVE, PRIVATE :: rugos
19!$OMP THREADPRIVATE(rugos)
20  REAL, ALLOCATABLE, DIMENSION(:),   SAVE, PRIVATE :: albedo
21!$OMP THREADPRIVATE(albedo) 
22  REAL, ALLOCATABLE, DIMENSION(:),   SAVE, PRIVATE :: sst
23!$OMP THREADPRIVATE(sst) 
24  LOGICAL,SAVE :: read_continents=.FALSE.
25!$OMP THREADPRIVATE(read_continents)
26
27CONTAINS
28!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
29!!
30!! Public subroutines :
31!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
32
33  SUBROUTINE limit_read_frac(itime, dtime, jour, pctsrf_new, is_modified)
34!
35! This subroutine is called from "change_srf_frac" for case of
36! ocean=force or from ocean_slab_frac for ocean=slab.
37! The fraction for all sub-surfaces at actual time step is returned.
38
39    USE dimphy
40    USE indice_sol_mod
41
42! Input arguments
43!****************************************************************************************
44    INTEGER, INTENT(IN) :: itime   ! time step
45    INTEGER, INTENT(IN) :: jour    ! current day
46    REAL   , INTENT(IN) :: dtime   ! length of time step
47 
48! Output arguments
49!****************************************************************************************
50    REAL, DIMENSION(klon,nbsrf), INTENT(OUT) :: pctsrf_new  ! sub surface fractions
51    LOGICAL, INTENT(OUT)                     :: is_modified ! true if pctsrf is modified at this time step
52
53! End declaration
54!****************************************************************************************
55
56! 1) Read file limit.nc
57    CALL limit_read_tot(itime, dtime, jour, is_modified)
58
59! 2) Return the fraction read in limit_read_tot
60    pctsrf_new(:,:) = pctsrf(:,:)
61   
62  END SUBROUTINE limit_read_frac
63
64!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
65
66  SUBROUTINE limit_read_rug_alb(itime, dtime, jour, &
67       knon, knindex, &
68       rugos_out, alb_out)
69!
70! This subroutine is called from surf_land_bucket.
71! The flag "ok_veget" must can not be true. If coupled run, "ocean=couple"
72! then this routine will call limit_read_tot.
73!
74    USE dimphy
75    USE surface_data
76
77! Input arguments
78!****************************************************************************************
79    INTEGER, INTENT(IN) :: itime                     ! numero du pas de temps courant
80    INTEGER, INTENT(IN) :: jour                      ! jour a lire dans l'annee
81    REAL   , INTENT(IN) :: dtime                     ! pas de temps de la physique (en s)
82    INTEGER, INTENT(IN) :: knon                      ! nomber of points on compressed grid
83    INTEGER, DIMENSION(klon), INTENT(IN) :: knindex  ! grid point number for compressed grid
84! Output arguments
85!****************************************************************************************
86    REAL, DIMENSION(klon), INTENT(OUT) :: rugos_out
87    REAL, DIMENSION(klon), INTENT(OUT) :: alb_out
88   
89! Local variables
90!****************************************************************************************
91    INTEGER :: i
92    LOGICAL :: is_modified
93!****************************************************************************************
94
95IF (type_ocean == 'couple'.OR. &
96         (type_ocean == 'slab' .AND. version_ocean == 'sicINT')) THEN
97       ! limit.nc has not yet been read. Do it now!
98       CALL limit_read_tot(itime, dtime, jour, is_modified)
99    END IF
100
101    DO i=1,knon
102       rugos_out(i) = rugos(knindex(i))
103       alb_out(i)  = albedo(knindex(i))
104    END DO
105
106  END SUBROUTINE limit_read_rug_alb
107
108!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
109
110  SUBROUTINE limit_read_sst(knon, knindex, sst_out)
111!
112! This subroutine returns the sea surface temperature already read from limit.nc.
113!
114    USE dimphy, ONLY : klon
115
116    INTEGER, INTENT(IN)                  :: knon     ! nomber of points on compressed grid
117    INTEGER, DIMENSION(klon), INTENT(IN) :: knindex  ! grid point number for compressed grid
118    REAL, DIMENSION(klon), INTENT(OUT)   :: sst_out
119
120    INTEGER :: i
121
122    DO i = 1, knon
123       sst_out(i) = sst(knindex(i))
124    END DO
125
126  END SUBROUTINE limit_read_sst
127
128!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
129!!
130!! Private subroutine :
131!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
132
133  SUBROUTINE limit_read_tot(itime, dtime, jour, is_modified)
134!
135! Read everything needed from limit.nc
136!
137! 0) Initialize
138! 1) Open the file limit.nc, if it is time
139! 2) Read fraction, if not type_ocean=couple
140! 3) Read sea surface temperature, if not type_ocean=couple
141! 4) Read albedo and rugosity for land surface, only in case of no vegetation model
142! 5) Close file and distribuate variables to all processus
143
144    USE dimphy
145    USE mod_grid_phy_lmdz
146    USE mod_phys_lmdz_para
147    USE surface_data, ONLY : type_ocean, ok_veget
148    USE netcdf
149    USE indice_sol_mod
150
151    IMPLICIT NONE
152   
153! In- and ouput arguments
154!****************************************************************************************
155    INTEGER, INTENT(IN) :: itime   ! numero du pas de temps courant
156    INTEGER, INTENT(IN) :: jour    ! jour a lire dans l'annee
157    REAL   , INTENT(IN) :: dtime   ! pas de temps de la physique (en s)
158
159    LOGICAL, INTENT(OUT) :: is_modified  ! true if pctsrf is modified at this time step
160
161! Locals variables with attribute SAVE
162!****************************************************************************************
163! frequence de lecture des conditions limites (en pas de physique)
164    INTEGER,SAVE                              :: lmt_pas
165!$OMP THREADPRIVATE(lmt_pas)
166    LOGICAL, SAVE                             :: first_call=.TRUE.
167!$OMP THREADPRIVATE(first_call) 
168    INTEGER, SAVE                             :: jour_lu = -1
169!$OMP THREADPRIVATE(jour_lu) 
170! Locals variables
171!****************************************************************************************
172    INTEGER                                   :: nid, nvarid
173    INTEGER                                   :: ii, ierr
174    INTEGER, DIMENSION(2)                     :: start, epais
175    REAL, DIMENSION(klon_glo,nbsrf)           :: pct_glo  ! fraction at global grid
176    REAL, DIMENSION(klon_glo)                 :: sst_glo  ! sea-surface temperature at global grid
177    REAL, DIMENSION(klon_glo)                 :: rug_glo  ! rugosity at global grid
178    REAL, DIMENSION(klon_glo)                 :: alb_glo  ! albedo at global grid
179    CHARACTER(len=20)                         :: modname='limit_read_mod'     
180
181! End declaration
182!****************************************************************************************
183
184!****************************************************************************************
185! 0) Initialization
186!
187!****************************************************************************************
188    IF (first_call) THEN
189       ! calculate number of time steps for one day
190       lmt_pas = NINT(86400./dtime * 1.0)
191       
192       ! Allocate module save variables
193       IF ( type_ocean /= 'couple' ) THEN
194          ALLOCATE(pctsrf(klon,nbsrf), sst(klon), stat=ierr)
195          IF (ierr /= 0) CALL abort_physic(modname, 'PB in allocating pctsrf and sst',1)
196       END IF
197
198       IF ( .NOT. ok_veget ) THEN
199          ALLOCATE(rugos(klon), albedo(klon), stat=ierr)
200          IF (ierr /= 0) CALL abort_physic(modname, 'PB in allocating rugos and albedo',1)
201       END IF
202
203       first_call=.FALSE.
204    ENDIF
205 
206!****************************************************************************************
207! 1) Open the file limit.nc if it is the right moment to read, once a day.
208!    The file is read only by the master thread of the master mpi process(is_mpi_root)
209!
210!****************************************************************************************
211
212    is_modified = .FALSE.
213    IF (MOD(itime-1, lmt_pas) == 0 .OR. jour_lu /= jour ) THEN   ! time to read
214       jour_lu = jour
215       is_modified = .TRUE.
216!$OMP MASTER  ! Only master thread
217       IF (is_mpi_root) THEN ! Only master processus
218
219          ierr = NF90_OPEN ('limit.nc', NF90_NOWRITE, nid)
220          IF (ierr /= NF90_NOERR) CALL abort_physic(modname,&
221               'Pb d''ouverture du fichier de conditions aux limites',1)
222         
223          ! La tranche de donnees a lire:
224          start(1) = 1
225          start(2) = jour
226          epais(1) = klon_glo
227          epais(2) = 1
228
229
230!****************************************************************************************
231! 2) Read fraction if not type_ocean=couple
232!
233!****************************************************************************************
234
235          IF ( type_ocean /= 'couple') THEN
236!
237! Ocean fraction
238             ierr = NF90_INQ_VARID(nid, 'FOCE', nvarid)
239             IF (ierr /= NF90_NOERR) CALL abort_physic(modname, 'Le champ <FOCE> est absent',1)
240             
241             ierr = NF90_GET_VAR(nid,nvarid,pct_glo(:,is_oce),start,epais)
242             IF (ierr /= NF90_NOERR) CALL abort_physic(modname,'Lecture echouee pour <FOCE>' ,1)
243!
244! Sea-ice fraction
245             ierr = NF90_INQ_VARID(nid, 'FSIC', nvarid)
246             IF (ierr /= NF90_NOERR) CALL abort_physic(modname,'Le champ <FSIC> est absent',1)
247
248             ierr = NF90_GET_VAR(nid,nvarid,pct_glo(:,is_sic),start,epais)
249             IF (ierr /= NF90_NOERR) CALL abort_physic(modname,'Lecture echouee pour <FSIC>' ,1)
250
251
252! Read land and continentals fraction only if asked for
253             IF (read_continents .OR. itime == 1) THEN
254!
255! Land fraction
256                ierr = NF90_INQ_VARID(nid, 'FTER', nvarid)
257                IF (ierr /= NF90_NOERR) CALL abort_physic(modname,'Le champ <FTER> est absent',1)
258               
259                ierr = NF90_GET_VAR(nid,nvarid,pct_glo(:,is_ter),start,epais)
260                IF (ierr /= NF90_NOERR) CALL abort_physic(modname,'Lecture echouee pour <FTER>',1)
261!
262! Continentale ice fraction
263                ierr = NF90_INQ_VARID(nid, 'FLIC', nvarid)
264                IF (ierr /= NF90_NOERR) CALL abort_physic(modname,'Le champ <FLIC> est absent',1)
265
266                ierr = NF90_GET_VAR(nid,nvarid,pct_glo(:,is_lic),start,epais)
267                IF (ierr /= NF90_NOERR) CALL abort_physic(modname,'Lecture echouee pour <FLIC>',1)
268             END IF
269
270          END IF ! type_ocean /= couple
271
272!****************************************************************************************
273! 3) Read sea-surface temperature, if not coupled ocean
274!
275!****************************************************************************************
276          IF ( type_ocean /= 'couple') THEN
277
278             ierr = NF90_INQ_VARID(nid, 'SST', nvarid)
279             IF (ierr /= NF90_NOERR) CALL abort_physic(modname,'Le champ <SST> est absent',1)
280
281             ierr = NF90_GET_VAR(nid,nvarid,sst_glo,start,epais)
282             IF (ierr /= NF90_NOERR) CALL abort_physic(modname,'Lecture echouee pour <SST>',1)
283         
284          END IF
285
286!****************************************************************************************
287! 4) Read albedo and rugosity for land surface, only in case of no vegetation model
288!
289!****************************************************************************************
290
291          IF (.NOT. ok_veget) THEN
292!
293! Read albedo
294             ierr = NF90_INQ_VARID(nid, 'ALB', nvarid)
295             IF (ierr /= NF90_NOERR) CALL abort_physic(modname,'Le champ <ALB> est absent',1)
296
297             ierr = NF90_GET_VAR(nid,nvarid,alb_glo,start,epais)
298             IF (ierr /= NF90_NOERR) CALL abort_physic(modname,'Lecture echouee pour <ALB>',1)
299!
300! Read rugosity
301             ierr = NF90_INQ_VARID(nid, 'RUG', nvarid)
302             IF (ierr /= NF90_NOERR) CALL abort_physic(modname,'Le champ <RUG> est absent',1)
303
304             ierr = NF90_GET_VAR(nid,nvarid,rug_glo,start,epais)
305             IF (ierr /= NF90_NOERR) CALL abort_physic(modname,'Lecture echouee pour <RUG>',1)
306
307          END IF
308
309!****************************************************************************************
310! 5) Close file and distribuate variables to all processus
311!
312!****************************************************************************************
313          ierr = NF90_CLOSE(nid)
314          IF (ierr /= NF90_NOERR) CALL abort_physic(modname,'Pb when closing file', 1)
315       ENDIF ! is_mpi_root
316
317!$OMP END MASTER
318!$OMP BARRIER
319
320       IF ( type_ocean /= 'couple') THEN
321          CALL Scatter(sst_glo,sst)
322          CALL Scatter(pct_glo(:,is_oce),pctsrf(:,is_oce))
323          CALL Scatter(pct_glo(:,is_sic),pctsrf(:,is_sic))
324          IF (read_continents .OR. itime == 1) THEN
325             CALL Scatter(pct_glo(:,is_ter),pctsrf(:,is_ter))
326             CALL Scatter(pct_glo(:,is_lic),pctsrf(:,is_lic))
327          END IF
328       END IF
329
330       IF (.NOT. ok_veget) THEN
331          CALL Scatter(alb_glo, albedo)
332          CALL Scatter(rug_glo, rugos)
333       END IF
334
335    ENDIF ! time to read
336
337  END SUBROUTINE limit_read_tot
338
339
340END MODULE limit_read_mod
Note: See TracBrowser for help on using the repository browser.