source: LMDZ5/branches/LMDZ_tree_FC/libf/dyn3dmem/dynredem_mod.F90 @ 5446

Last change on this file since 5446 was 2299, checked in by dcugnet, 10 years ago

In dyn3d/:
etat0dyn_netcdf.F90: "startget_dyn3d" syntax slightly simplified.
dynredem.F90: Shortcut routines (put_var*, cre_var,
dynredem_write_*, dynredem_read_u)

modified to match dyn3dmem version and put in

module dyredem_mod.F90.
dynetat0.F90 -> *.f90: Few simplifications (no usage of NC_DOUBLE
needed => no precompilation)

Add tracers initialization in the isotope case

suppressed by accident.
dynredem_mod.F90: Created to mimic dyn3dmem equivalent.

In dyn3dmem/:
dynetat0_loc.F -> *.f90: Converted into fortran 90 to match the dyn3d
version.
dynredem_loc.F -> *.F90: Converted into fortran 90.
dynredem_mod.F90: Add some shortcut routines to match the dyn3d
version.

In phylmd/:
phyredem.F90: Bug fix: nsw instead of nsoilmx was used as
Tsoil second maximum index.

Bug fix: fevap instead of snow was saved for

"SNOW".
etat0phys_netcdf.F90: "filtreg_mod" module usage suppressed.

Local variable rugo computation removed (not

used).

In dynlonlat_phylonlat/:
grid_atob_m.F90 -> *.f90 DOUBLE PRECISION variables usage removed.

Precompilation o longer needed => .F90 extension.

  • 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
File size: 9.0 KB
Line 
1MODULE dynredem_mod
2
3  USE dimensions_mod
4  USE parallel_lmdz
5  USE mod_hallo
6  USE netcdf
7  PRIVATE
8  PUBLIC :: dynredem_write_u, dynredem_write_v, dynredem_read_u, err
9  PUBLIC :: cre_var, get_var1, put_var, fil, modname, msg
10  CHARACTER(LEN=256), SAVE :: fil, modname
11  INTEGER,            SAVE :: nvarid
12
13
14CONTAINS
15
16
17!===============================================================================
18!
19SUBROUTINE dynredem_write_u(ncid,id,var,ll)
20!
21!===============================================================================
22  IMPLICIT NONE
23!===============================================================================
24! Arguments:
25  INTEGER,          INTENT(IN) :: ncid
26  CHARACTER(LEN=*), INTENT(IN) :: id
27  REAL,             INTENT(IN) :: var(ijb_u:ije_u,ll)
28  INTEGER,          INTENT(IN) :: ll
29!===============================================================================
30! Local variables:
31  REAL, ALLOCATABLE, SAVE :: var_tmp(:,:), var_glo(:)
32  INTEGER :: start(4), count(4), l, ierr
33!===============================================================================
34  start(:)=[1,1,1,1]; count(:)=[iip1,jjp1,1,1]
35
36!$OMP MASTER
37  IF(mpi_rank==0) CALL err(NF90_INQ_VARID(ncid,id,nvarid),"inq",id)
38!$OMP END MASTER
39
40!$OMP MASTER
41  ALLOCATE(var_tmp(ijb_u:ije_u,ll),var_glo(ip1jmp1))
42!$OMP END MASTER
43!$OMP BARRIER
44
45!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
46  DO l=1,ll; var_tmp(:,l)=var(:,l); END DO
47  DO l=1,ll
48    CALL gather_field_u(var_tmp(:,l),var_glo,1)
49    IF(mpi_rank==0) THEN
50    !$OMP MASTER
51      start(3)=l
52      CALL err(NF90_PUT_VAR(ncid,nvarid,var_glo,start,count),"put",id)
53    !$OMP END MASTER
54    END IF
55  END DO
56!$OMP BARRIER
57!$OMP MASTER
58  DEALLOCATE(var_glo,var_tmp)
59!$OMP END MASTER
60!$OMP BARRIER
61 
62END SUBROUTINE dynredem_write_u
63!
64!===============================================================================
65
66
67!===============================================================================
68!
69SUBROUTINE dynredem_write_v(ncid,id,var,ll)
70!
71!===============================================================================
72  IMPLICIT NONE
73!===============================================================================
74! Arguments:
75  INTEGER,          INTENT(IN) :: ncid
76  CHARACTER(LEN=*), INTENT(IN) :: id
77  REAL,             INTENT(IN) :: var(ijb_v:ije_v,ll)
78  INTEGER,          INTENT(IN) :: ll
79!===============================================================================
80! Local variables:
81  REAL, ALLOCATABLE, SAVE :: var_tmp(:,:), var_glo(:)
82  INTEGER :: start(4), count(4), l, ierr
83!===============================================================================
84  start(:)=[1,1,1,1]; count(:)=[iip1,jjm,1,1]
85
86!$OMP MASTER
87  IF(mpi_rank==0) CALL err(NF90_INQ_VARID(ncid,id,nvarid),"inq",id)
88!$OMP END MASTER
89
90!$OMP MASTER
91  ALLOCATE(var_tmp(ijb_v:ije_v,ll),var_glo(ip1jm))
92!$OMP END MASTER
93!$OMP BARRIER
94
95!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
96  DO l=1,ll; var_tmp(:,l)=var(:,l); END DO
97  DO l=1,ll
98    CALL gather_field_v(var_tmp(:,l),var_glo,1)
99    IF(mpi_rank==0) THEN
100    !$OMP MASTER
101      start(3)=l
102      CALL err(NF90_PUT_VAR(ncid,nvarid,var_glo,start,count),"put",id)
103    !$OMP END MASTER
104    END IF
105  END DO
106!$OMP BARRIER
107!$OMP MASTER
108  DEALLOCATE(var_glo,var_tmp)
109!$OMP END MASTER
110!$OMP BARRIER
111 
112END SUBROUTINE dynredem_write_v
113!
114!===============================================================================
115
116
117!===============================================================================
118!
119SUBROUTINE dynredem_read_u(ncid,id,var,ll)
120!
121!===============================================================================
122  IMPLICIT NONE
123!===============================================================================
124! Arguments:
125  INTEGER,          INTENT(IN)  :: ncid
126  CHARACTER(LEN=*), INTENT(IN)  :: id
127  REAL,             INTENT(OUT) :: var(ijb_u:ije_u,ll)
128  INTEGER,          INTENT(IN)  :: ll
129!===============================================================================
130! Local variables:
131  REAL, ALLOCATABLE, SAVE :: var_tmp(:,:), var_glo(:)
132  INTEGER :: start(4), count(4), l, ierr
133!===============================================================================
134  start(:)=[1,1,1,1]; count(:)=[iip1,jjp1,1,1]
135
136!$OMP MASTER
137  IF(mpi_rank==0) CALL err(NF90_INQ_VARID(ncid,id,nvarid),'inq',id)
138!$OMP END MASTER
139
140!$OMP MASTER
141  ALLOCATE(var_tmp(ijb_u:ije_u,ll),var_glo(ip1jmp1))
142!$OMP END MASTER
143!$OMP BARRIER
144
145  DO l=1,ll
146    IF(mpi_rank==0) THEN
147    !$OMP MASTER
148      start(3)=l
149      CALL err(NF90_GET_VAR(ncid,nvarid,var_glo,start,count),"get",id)
150    !$OMP END MASTER
151    END IF
152    CALL scatter_field_u(var_glo,var_tmp(:,l),1)
153  END DO
154
155!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
156  DO l=1,ll; var(:,l)=var_tmp(:,l); END DO
157   
158!$OMP BARRIER
159!$OMP MASTER
160  DEALLOCATE(var_glo,var_tmp)
161!$OMP END MASTER
162!$OMP BARRIER
163 
164END SUBROUTINE dynredem_read_u   
165!
166!===============================================================================
167
168
169!===============================================================================
170!
171SUBROUTINE cre_var(ncid,var,title,did,units)
172!
173!===============================================================================
174  IMPLICIT NONE
175!===============================================================================
176! Arguments:
177  INTEGER,                    INTENT(IN) :: ncid
178  CHARACTER(LEN=*),           INTENT(IN) :: var, title
179  INTEGER,                    INTENT(IN) :: did(:)
180  CHARACTER(LEN=*), OPTIONAL, INTENT(IN) :: units
181!===============================================================================
182#ifdef NC_DOUBLE
183  CALL err(NF90_DEF_VAR(ncid,var,NF90_DOUBLE,did,nvarid),"inq",var)
184#else
185  CALL err(NF90_DEF_VAR(ncid,var,NF90_FLOAT ,did,nvarid),"inq",var)
186#endif
187  IF(title/="")      CALL err(NF90_PUT_ATT(ncid,nvarid,"title",title),var)
188  IF(PRESENT(units)) CALL err(NF90_PUT_ATT(ncid,nvarid,"units",units),var)
189
190END SUBROUTINE cre_var
191!
192!===============================================================================
193
194
195!===============================================================================
196!
197SUBROUTINE put_var(ncid,var,title,did,v,units)
198!
199!===============================================================================
200  IMPLICIT NONE
201!===============================================================================
202! Arguments:
203  INTEGER,                    INTENT(IN) :: ncid
204  CHARACTER(LEN=*),           INTENT(IN) :: var, title
205  INTEGER,                    INTENT(IN) :: did(:)
206  REAL,                       INTENT(IN) :: v(:)
207  CHARACTER(LEN=*), OPTIONAL, INTENT(IN) :: units
208!===============================================================================
209  INTEGER :: nd, k, nn(2)
210  IF(     PRESENT(units)) CALL cre_var(ncid,var,title,did,units)
211  IF(.NOT.PRESENT(units)) CALL cre_var(ncid,var,title,did)
212  CALL err(NF90_ENDDEF(ncid))
213  nd=SIZE(did)
214  DO k=1,nd; CALL err(NF90_INQUIRE_DIMENSION(ncid,did(k),len=nn(k))); END DO
215  IF(nd==1) CALL err(NF90_PUT_VAR(ncid,nvarid,RESHAPE(v,nn(1:1))),var)
216  IF(nd==2) CALL err(NF90_PUT_VAR(ncid,nvarid,RESHAPE(v,nn(1:2))),var)
217  CALL err(NF90_REDEF(ncid))
218END SUBROUTINE put_var
219!
220!===============================================================================
221
222
223!===============================================================================
224!
225FUNCTION msg(typ,nam)
226!
227!===============================================================================
228  IMPLICIT NONE
229!===============================================================================
230! Arguments:
231  CHARACTER(LEN=256)                     :: msg    !--- STANDARDIZED MESSAGE
232  CHARACTER(LEN=*),           INTENT(IN) :: typ    !--- TYPE OF OPERATION
233  CHARACTER(LEN=*), OPTIONAL, INTENT(IN) :: nam    !--- FIELD NAME
234!===============================================================================
235  SELECT CASE(typ)
236    CASE('open');  msg="Opening failed for <"//TRIM(fil)//">"
237    CASE('close'); msg="Closing failed for <"//TRIM(fil)//">"
238    CASE('get');   msg="Reading failed for <"//TRIM(nam)//">"
239    CASE('put');   msg="Writting failed for <"//TRIM(nam)//">"
240    CASE('inq');   msg="Missing field <"//TRIM(nam)//">"
241    CASE('fnd');   msg="Found field <"//TRIM(nam)//">"
242  END SELECT
243  msg=TRIM(msg)//" in file <"//TRIM(fil)//">"
244
245END FUNCTION msg
246!
247!===============================================================================
248
249
250!===============================================================================
251!
252SUBROUTINE err(ierr,typ,nam)
253!
254!===============================================================================
255  IMPLICIT NONE
256!===============================================================================
257! Arguments:
258  INTEGER,                    INTENT(IN) :: ierr   !--- NetCDF ERROR CODE
259  CHARACTER(LEN=*), OPTIONAL, INTENT(IN) :: typ    !--- TYPE OF OPERATION
260  CHARACTER(LEN=*), OPTIONAL, INTENT(IN) :: nam    !--- FIELD NAME
261!===============================================================================
262  IF(ierr==NF90_NoERR) RETURN
263  IF(.NOT.PRESENT(typ)) THEN
264    CALL ABORT_gcm(modname,NF90_STRERROR(ierr),ierr)
265  ELSE
266    CALL ABORT_gcm(modname,msg(typ,nam),ierr)
267  END IF
268
269END SUBROUTINE err
270!
271!===============================================================================
272
273END MODULE dynredem_mod   
274
275   
276   
Note: See TracBrowser for help on using the repository browser.