source: LMDZ6/branches/Amaury_dev/libf/dyn3dmem/dynredem_mod.F90 @ 5153

Last change on this file since 5153 was 5128, checked in by abarral, 4 months ago

Correct bug in vlspltqs_loc.f90 from r2270 where we call SSUM with incorrect arguments.
Merge the three different versions of abort_gcm into one
Fix seq, para 3D compilation broken from r5107 onwards
(lint) usual + Remove uneeded fixed-form continuations

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