source: trunk/LMDZ.COMMON/libf/dyn3dpar/dynredem_p.F90 @ 3567

Last change on this file since 3567 was 3510, checked in by jbclement, 6 weeks ago

Dynamic:
Following of r3509, the description of the 'controle' array in the "start.nc" file is adapted to the planet type (earth, mars or titan).
JBC

File size: 18.7 KB
Line 
1!
2! $Id: dynredem_p.F 1635 2012-07-12 11:37:16Z lguez $
3!
4SUBROUTINE dynredem0_p(fichnom,iday_end,phis)
5#ifdef CPP_IOIPSL
6  USE IOIPSL
7#endif
8  USE parallel_lmdz, ONLY: mpi_rank
9  USE infotrac, ONLY: nqtot, tname, ttext
10  USE netcdf, ONLY: NF90_CREATE, NF90_DEF_DIM, NF90_INQ_VARID, NF90_GLOBAL,    &
11                    NF90_CLOSE,  NF90_PUT_ATT, NF90_UNLIMITED, NF90_CLOBBER
12  USE dynredem_mod, ONLY: cre_var, put_var1, put_var2, put_char1, err, modname, fil, int2fmtstr
13  use netcdf95, only: NF95_PUT_VAR
14  use control_mod, only : planet_type
15  USE comvert_mod, ONLY: ap,bp,aps,bps,presnivs,pseudoalt,pa,preff, &
16                        nivsig,nivsigs
17  USE comconst_mod, ONLY: daysec,dtvr,rad,omeg,g,cpp,kappa,pi
18  USE logic_mod, ONLY: fxyhypb,ysinus
19  USE serre_mod, ONLY: clon,clat,grossismx,grossismy,dzoomx,dzoomy, &
20                        taux,tauy
21  USE temps_mod, ONLY: annee_ref,day_ref,itau_dyn,itaufin, &
22                        start_time,hour_ini
23  USE ener_mod, ONLY: etot0,ptot0,ztot0,stot0,ang0
24
25  IMPLICIT NONE
26!=======================================================================
27! Writting the NetCDF restart file (initialisation)
28!=======================================================================
29!   Declarations:
30!   -------------
31  include "dimensions.h"
32  include "paramet.h"
33  include "comgeom2.h"
34  include "netcdf.inc"
35  include "iniprint.h"
36
37!===============================================================================
38! Arguments:
39  CHARACTER(LEN=*), INTENT(IN) :: fichnom          !--- FILE NAME
40  INTEGER,          INTENT(IN) :: iday_end         !---
41  REAL,             INTENT(IN) :: phis(iip1, jjp1) !--- GROUND GEOPOTENTIAL
42!===============================================================================
43! Local variables:
44  INTEGER :: iq,l
45  INTEGER, PARAMETER :: length=100
46  INTEGER, PARAMETER :: ndscrpt = 42
47  INTEGER :: ldscrpt
48  REAL :: tab_cntrl(length) ! run parameters
49  character(ndscrpt), dimension(:), allocatable :: dscrpt_tab_cntrl
50  INTEGER :: ierr
51  CHARACTER(LEN=80) :: abort_message
52
53!   For NetCDF:
54  CHARACTER(LEN=30) :: unites
55  INTEGER :: indexID, descrptID, dscrpt_sID
56  INTEGER :: rlonuID, rlonvID, rlatuID, rlatvID
57  INTEGER :: sID, sigID, nID, vID, timID
58  INTEGER :: yyears0, jjour0, mmois0
59  REAL    :: zan0, zjulian, hours
60
61  CHARACTER(len=12) :: start_file_type="earth" ! default start file type
62  INTEGER :: idecal
63
64!===============================================================================
65  if (mpi_rank==0) then ! only the master reads input file
66  ! fill dynredem_mod module variables
67  modname='dynredem0_p'; fil=fichnom
68
69#ifdef CPP_IOIPSL
70  call ymds2ju(annee_ref, 1, iday_end, 0.0, zjulian)
71  call ju2ymds(zjulian, yyears0, mmois0, jjour0, hours)
72#else
73! set yyears0, mmois0, jjour0 to 0,1,1 (hours is not used)
74  yyears0=0
75  mmois0=1
76  jjour0=1
77#endif       
78
79  !!! AS: idecal is a hack to be able to read planeto starts...
80  !!!     .... while keeping everything OK for LMDZ EARTH
81  if ((planet_type.eq."generic").or.(planet_type.eq."mars")) then
82    write(lunout,*) trim(modname),' : Planeto-like start file'
83    start_file_type="planeto"
84    idecal = 4
85    ldscrpt = 28
86  else
87    if ( planet_type.eq."titan" ) then
88      ! Titan inherited Earth-like start files with idecal=5
89      write(lunout,*) trim(modname),' : Titan start file'
90      ldscrpt = 29
91    else
92      write(lunout,*) trim(modname),' : Earth-like start file'
93      ldscrpt = 32
94    endif
95    idecal = 5
96  endif
97
98  tab_cntrl(:)  = 0.
99  allocate(dscrpt_tab_cntrl(ldscrpt))
100  dscrpt_tab_cntrl(1:5) = (/                        &
101      "(1)  Number of nodes along longitude      ", &
102      "(2)  Number of nodes along latitude       ", &
103      "(3)  Number of atmospheric layers         ", &
104      "(4)  Reference day                        ", &
105      "(5)  Reference year                       "/)
106  tab_cntrl(1)  = REAL(iim)
107  tab_cntrl(2)  = REAL(jjm)
108  tab_cntrl(3)  = REAL(llm)
109  if (start_file_type.eq."earth") then
110    tab_cntrl(4)=REAL(day_ref)
111  else
112    !tab_cntrl(4)=REAL(day_end)
113    tab_cntrl(4)=REAL(iday_end)
114    dscrpt_tab_cntrl(4) = "(4)  Initial day                          "
115  endif
116  tab_cntrl(5)  = REAL(annee_ref)
117  tab_cntrl(idecal+1)  = rad
118  tab_cntrl(idecal+2)  = omeg
119  tab_cntrl(idecal+3)  = g
120  tab_cntrl(idecal+4)  = cpp
121  tab_cntrl(idecal+5) = kappa
122  tab_cntrl(idecal+6) = daysec
123  tab_cntrl(idecal+7) = dtvr
124  tab_cntrl(idecal+8) = etot0
125  tab_cntrl(idecal+9) = ptot0
126  tab_cntrl(idecal+10) = ztot0
127  tab_cntrl(idecal+11) = stot0
128  tab_cntrl(idecal+12) = ang0
129  tab_cntrl(idecal+13) = pa
130  tab_cntrl(idecal+14) = preff
131
132!    .....    parameters for the zoom      ......   
133  tab_cntrl(idecal+15)  = clon
134  tab_cntrl(idecal+16)  = clat
135  tab_cntrl(idecal+17)  = grossismx
136  tab_cntrl(idecal+18)  = grossismy
137  dscrpt_tab_cntrl(idecal + 1:idecal + 18) = (/     &
138      int2fmtstr(idecal + 1)//" Radius of the planet                 ", &
139      int2fmtstr(idecal + 2)//" Rotation of the planet (rad/s)       ", &
140      int2fmtstr(idecal + 3)//" Gravity (m/s2)                       ", &
141      int2fmtstr(idecal + 4)//" Specific heat Cp (J.kg-1.K-1)        ", &
142      int2fmtstr(idecal + 5)//" = r/Cp (=kappa)                      ", &
143      int2fmtstr(idecal + 6)//" Lenght of a sol (s)                  ", &
144      int2fmtstr(idecal + 7)//" Dynamical time step (s)              ", &
145      int2fmtstr(idecal + 8)//" Total energy                         ", &
146      int2fmtstr(idecal + 9)//" Total pressure                       ", &
147      int2fmtstr(idecal + 10)//" Total enstrophy                      ", &
148      int2fmtstr(idecal + 11)//" Total enthalpy                       ", &
149      int2fmtstr(idecal + 12)//" Total angular momentum               ", &
150      int2fmtstr(idecal + 13)//" Reference pressure (Pa)              ", &
151      int2fmtstr(idecal + 14)//" Reference surface pressure (Pa)      ", &
152      int2fmtstr(idecal + 15)//" Longitude of center of zoom          ", &
153      int2fmtstr(idecal + 16)//" Latitude of center of zoom           ", &
154      int2fmtstr(idecal + 17)//" Zooming factor, along longitude      ", &
155      int2fmtstr(idecal + 18)//" Zooming factor, along latitude       "/)
156!
157  IF ( fxyhypb )   THEN
158    tab_cntrl(idecal+19) = 1.
159    tab_cntrl(idecal+20) = dzoomx
160    tab_cntrl(idecal+21) = dzoomy
161    tab_cntrl(idecal+22) = 0.
162    tab_cntrl(idecal+23) = taux
163    tab_cntrl(idecal+24) = tauy
164    dscrpt_tab_cntrl(idecal + 19:idecal + 24) = (/  &
165      int2fmtstr(idecal + 19)//" 1.                                   ", &
166      int2fmtstr(idecal + 20)//" Extention (in longitude) of zoom     ", &
167      int2fmtstr(idecal + 21)//" Extention (in latitude) of zoom      ", &
168      int2fmtstr(idecal + 22)//" -                                    ", &
169      int2fmtstr(idecal + 23)//" Stiffness factor of zoom in longitude", &
170      int2fmtstr(idecal + 24)//" Stiffness factor of zoom in latitude "/)
171  ELSE
172    tab_cntrl(idecal+19) = 0.
173    tab_cntrl(idecal+20) = dzoomx
174    tab_cntrl(idecal+21) = dzoomy
175    tab_cntrl(idecal+22) = 0.
176    tab_cntrl(idecal+23) = 0.
177    tab_cntrl(idecal+24) = 0.
178    IF( ysinus )  tab_cntrl(idecal+22) = 1.
179    dscrpt_tab_cntrl(idecal + 19:idecal + 24) = (/  &
180      int2fmtstr(idecal + 19)//" -                                    ", &
181      int2fmtstr(idecal + 20)//" Extention (in longitude) of zoom     ", &
182      int2fmtstr(idecal + 21)//" Extention (in latitude) of zoom      ", &
183      int2fmtstr(idecal + 22)//" 1.                                   ", &
184      int2fmtstr(idecal + 23)//" Stiffness factor of zoom in longitude", &
185      int2fmtstr(idecal + 24)//" Stiffness factor of zoom in latitude "/)
186  ENDIF
187
188  if (start_file_type.eq."earth") then
189    tab_cntrl(idecal+25) = REAL(iday_end)
190    tab_cntrl(idecal+26) = REAL(itau_dyn + itaufin)
191! start_time: start_time of simulation (not necessarily 0.)
192    tab_cntrl(idecal+27) = start_time
193    dscrpt_tab_cntrl(idecal + 25:idecal + 27) = (/  &
194      "(30) Initial day                          ", &
195      "(31) Extention (in longitude) of zoom     ", &
196      "(32) Start time of simulation             "/)
197  endif
198
199  if (planet_type=="mars") then ! For Mars only
200!    tab_cntrl(29)=hour_ini
201    tab_cntrl(29)=0
202  endif
203
204!--- File creation
205  CALL err(NF90_CREATE(fichnom,NF90_CLOBBER,nid))
206
207!--- Some global attributes
208  CALL err(NF90_PUT_ATT(nid,NF90_GLOBAL,"title","Fichier demarrage dynamique"))
209
210!--- Dimensions
211  if (start_file_type.eq."earth") then
212    CALL err(NF90_DEF_DIM(nid,"index", length, indexID))
213    CALL err(NF90_DEF_DIM(nid,"rlonu", iip1,   rlonuID))
214    CALL err(NF90_DEF_DIM(nid,"rlatu", jjp1,   rlatuID))
215    CALL err(NF90_DEF_DIM(nid,"rlonv", iip1,   rlonvID))
216    CALL err(NF90_DEF_DIM(nid,"rlatv", jjm,    rlatvID))
217    CALL err(NF90_DEF_DIM(nid,"sigs",  llm,        sID))
218    CALL err(NF90_DEF_DIM(nid,"sig",   llmp1,    sigID))
219    CALL err(NF90_DEF_DIM(nid,"temps", NF90_UNLIMITED, timID))
220    CALL err(NF90_DEF_DIM(nid,"descriptor", ldscrpt, descrptID))
221    CALL err(NF90_DEF_DIM(nid,"description_size", ndscrpt, dscrpt_sID))
222  else
223    CALL err(NF90_DEF_DIM(nid,"index", length, indexID))
224    CALL err(NF90_DEF_DIM(nid,"rlonu", iip1,   rlonuID))
225    CALL err(NF90_DEF_DIM(nid,"latitude", jjp1,   rlatuID))
226    CALL err(NF90_DEF_DIM(nid,"longitude", iip1,   rlonvID))
227    CALL err(NF90_DEF_DIM(nid,"rlatv", jjm,    rlatvID))
228    CALL err(NF90_DEF_DIM(nid,"altitude",  llm,        sID))
229    CALL err(NF90_DEF_DIM(nid,"interlayer",   llmp1,    sigID))
230    CALL err(NF90_DEF_DIM(nid,"Time", NF90_UNLIMITED, timID))
231    CALL err(NF90_DEF_DIM(nid,"descriptor", ldscrpt, descrptID))
232    CALL err(NF90_DEF_DIM(nid,"description_size", ndscrpt, dscrpt_sID))
233  endif
234
235!--- Define and save invariant fields
236  CALL put_var1(nid,"controle","Parametres de controle" ,[indexID],tab_cntrl)
237  CALL put_char1(nid,"controle_descriptor","Description of control parameters",[dscrpt_sID,descrptID],dscrpt_tab_cntrl)
238  CALL put_var1(nid,"rlonu"   ,"Longitudes des points U",[rlonuID],rlonu)
239  CALL put_var1(nid,"rlatu"   ,"Latitudes des points U" ,[rlatuID],rlatu)
240  CALL put_var1(nid,"rlonv"   ,"Longitudes des points V",[rlonvID],rlonv)
241  CALL put_var1(nid,"rlatv"   ,"Latitudes des points V" ,[rlatvID],rlatv)
242  if (start_file_type.eq."earth") then
243    CALL put_var1(nid,"nivsigs" ,"Numero naturel des couches s"    ,[sID]  ,nivsigs)
244    CALL put_var1(nid,"nivsig"  ,"Numero naturel des couches sigma",[sigID],nivsig)
245  endif ! of if (start_file_type.eq."earth")
246  CALL put_var1(nid,"ap"      ,"Coefficient A pour hybride"      ,[sigID],ap)
247  CALL put_var1(nid,"bp"      ,"Coefficient B pour hybride"      ,[sigID],bp)
248  if (start_file_type.ne."earth") then
249    CALL put_var1(nid,"aps","Coef AS: hybrid pressure at midlayers",[sID],aps)
250    CALL put_var1(nid,"bps","Coef BS: hybrid sigma at midlayers",[sID],bps)
251  endif ! of if (start_file_type.eq."earth")
252  CALL put_var1(nid,"presnivs",""                                ,[sID]  ,presnivs)
253  if (start_file_type.ne."earth") then
254        ierr = NF_REDEF (nid)
255#ifdef NC_DOUBLE
256        ierr = NF_DEF_VAR(nid,"latitude",NF_DOUBLE,1,rlatuID,vID)
257#else
258        ierr = NF_DEF_VAR(nid,"latitude",NF_FLOAT,1,rlatuID,vID)
259#endif
260        ierr =NF_PUT_ATT_TEXT(nid,vID,'units',13,"degrees_north")
261        ierr = NF_PUT_ATT_TEXT (nid,vID,"long_name", 14, &
262              "North latitude")
263        ierr = NF_ENDDEF(nid)
264        call NF95_PUT_VAR(nid,vID,rlatu*180/pi)
265!
266        ierr = NF_REDEF (nid)
267#ifdef NC_DOUBLE
268        ierr=NF_DEF_VAR(nid,"longitude",NF_DOUBLE,1,rlonvID,vID)
269#else
270        ierr=NF_DEF_VAR(nid,"longitude",NF_FLOAT,1,rlonvID,vID)
271#endif
272        ierr = NF_PUT_ATT_TEXT (nid,vID,"long_name", 14, &
273              "East longitude")
274        ierr = NF_PUT_ATT_TEXT(nid,vID,'units',12,"degrees_east")
275        ierr = NF_ENDDEF(nid)
276        call NF95_PUT_VAR(nid,vID,rlonv*180/pi)
277!
278        ierr = NF_REDEF (nid)
279#ifdef NC_DOUBLE
280        ierr = NF_DEF_VAR (nid, "altitude", NF_DOUBLE, 1, &
281             sID,vID)
282#else
283        ierr = NF_DEF_VAR (nid, "altitude", NF_FLOAT, 1, &
284             sID,vID)
285#endif
286        ierr = NF_PUT_ATT_TEXT(nid,vID,"long_name",10,"pseudo-alt")
287        ierr = NF_PUT_ATT_TEXT (nid,vID,'units',2,"km")
288        ierr = NF_PUT_ATT_TEXT (nid,vID,'positive',2,"up")
289        ierr = NF_ENDDEF(nid)
290        call NF95_PUT_VAR(nid,vID,pseudoalt)
291        CALL err(NF_REDEF(nid))
292  endif ! of if (start_file_type.ne."earth")
293
294! covariant <-> contravariant <-> natural conversion coefficients
295  CALL put_var2(nid,"cu","Coefficient de passage pour U",[rlonuID,rlatuID],cu)
296  CALL put_var2(nid,"cv","Coefficient de passage pour V",[rlonvID,rlatvID],cv)
297  CALL put_var2(nid,"aire","Aires de chaque maille"     ,[rlonvID,rlatuID],aire)
298  CALL put_var2(nid,"phisinit","Geopotentiel au sol"    ,[rlonvID,rlatuID],phis)
299
300
301! Define variables that will be stored later:
302  WRITE(unites,"('days since ',i4,'-',i2.2,'-',i2.2,' 00:00:00')") &
303               yyears0,mmois0,jjour0
304  IF (start_file_type.eq."earth") THEN
305    CALL cre_var(nid,"temps","Temps de simulation",[timID],unites)
306  ELSE
307    CALL cre_var(nid,"Time","Temps de simulation",[timID],unites)
308  ENDIF
309
310  CALL cre_var(nid,"ucov" ,"Vitesse U"  ,[rlonuID,rlatuID,sID,timID])
311  CALL cre_var(nid,"vcov" ,"Vitesse V"  ,[rlonvID,rlatvID,sID,timID])
312  CALL cre_var(nid,"teta" ,"Temperature",[rlonvID,rlatuID,sID,timID])
313
314  IF(nqtot.GE.1) THEN
315    DO iq=1,nqtot
316      CALL cre_var(nid,tname(iq),ttext(iq),[rlonvID,rlatuID,sID,timID])
317    END DO
318  ENDIF
319
320  CALL cre_var(nid,"masse","Masse d air"    ,[rlonvID,rlatuID,sID,timID])
321  CALL cre_var(nid,"ps"   ,"Pression au sol",[rlonvID,rlatuID    ,timID])
322
323  CALL err(NF90_CLOSE (nid)) ! close file
324
325  WRITE(lunout,*)TRIM(modname)//': iim,jjm,llm,iday_end',iim,jjm,llm,iday_end
326  WRITE(lunout,*)TRIM(modname)//': rad,omeg,g,cpp,kappa',rad,omeg,g,cpp,kappa
327
328  deallocate(dscrpt_tab_cntrl)
329
330  endif ! of if (mpi_rank==0)
331
332END SUBROUTINE dynredem0_p
333
334!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
335
336SUBROUTINE dynredem1_p(fichnom,time,vcov,ucov,teta,q,masse,ps)
337!
338!-------------------------------------------------------------------------------
339! Purpose: Write the NetCDF restart file (append).
340!-------------------------------------------------------------------------------
341  USE parallel_lmdz, ONLY: mpi_rank, gather_field
342  USE infotrac, ONLY: nqtot, tname, type_trac
343  USE control_mod, only : planet_type
344  USE netcdf,   ONLY: NF90_OPEN,  NF90_NOWRITE, NF90_GET_VAR, NF90_INQ_VARID,  &
345                      NF90_CLOSE, NF90_WRITE,   NF90_PUT_VAR, NF90_NoErr
346  use netcdf95, only: NF95_PUT_VAR
347  USE temps_mod, ONLY: itaufin,itau_dyn
348  USE dynredem_mod, ONLY: dynredem_write_u, dynredem_write_v, dynredem_read_u, &
349                          err, modname, fil, msg
350 
351  IMPLICIT NONE
352  include "dimensions.h"
353  include "paramet.h"
354  include "netcdf.inc"
355  include "comgeom.h"
356  include "iniprint.h"
357!===============================================================================
358! Arguments:
359  CHARACTER(LEN=*), INTENT(IN) :: fichnom              !-- FILE NAME
360  REAL, INTENT(IN)    ::  time                         !-- TIME
361  REAL, INTENT(IN)    ::  vcov(iip1,jjm, llm)          !-- V COVARIANT WIND
362  REAL, INTENT(IN)    ::  ucov(iip1,jjp1,llm)          !-- U COVARIANT WIND
363  REAL, INTENT(IN)    ::  teta(iip1,jjp1,llm)          !-- POTENTIAL TEMPERATURE
364  REAL, INTENT(INOUT) ::     q(iip1,jjp1,llm,nqtot)    !-- TRACERS
365  REAL, INTENT(IN)    :: masse(iip1,jjp1,llm)          !-- MASS PER CELL
366  REAL, INTENT(IN)    ::    ps(iip1,jjp1)              !-- GROUND PRESSURE
367!===============================================================================
368! Local variables:
369  INTEGER :: l, iq, nid, vID, ierr, nid_trac, vID_trac
370  INTEGER,SAVE :: nb=0
371  INTEGER, PARAMETER :: length=100
372  REAL               :: tab_cntrl(length) ! tableau des parametres du run
373  CHARACTER(LEN=256) :: var, dum
374  LOGICAL            :: lread_inca
375  CHARACTER(LEN=80) :: abort_message
376  CHARACTER(len=12) :: start_file_type="earth" ! default start file type
377
378  ! fill dynredem_mod module variables
379  modname='dynredem1_p'; fil=fichnom
380
381  ! Gather datasets
382  call Gather_Field(ucov,ip1jmp1,llm,0)
383  call Gather_Field(vcov,ip1jm,llm,0)
384  call Gather_Field(teta,ip1jmp1,llm,0)
385  call Gather_Field(masse,ip1jmp1,llm,0)
386  call Gather_Field(ps,ip1jmp1,1,0)
387     
388  do iq=1,nqtot
389    call Gather_Field(q(:,:,:,iq),ip1jmp1,llm,0)
390  enddo
391
392  IF (mpi_rank==0) THEN ! only the master writes restart file
393
394  if ((planet_type.eq."generic").or.(planet_type.eq."mars")) then
395      write(lunout,*) trim(modname),' : Planeto-like start file'
396      start_file_type="planeto"
397  else
398      write(lunout,*) trim(modname),' : Earth-like start file'
399  endif
400
401  CALL err(NF90_OPEN(fil,NF90_WRITE,nid),"open",fil)
402
403!--- Write/extend time coordinate
404  nb = nb + 1
405  if (start_file_type.eq."earth") then
406        ierr = NF_INQ_VARID(nid, "temps", vID)
407        IF (ierr .NE. NF_NOERR) THEN
408          write(lunout,*) NF_STRERROR(ierr)
409          abort_message='Variable temps n est pas definie'
410          CALL abort_gcm(modname,abort_message,ierr)
411        ENDIF
412 else
413        ierr = NF_INQ_VARID(nid,"Time", vID)
414        IF (ierr .NE. NF_NOERR) THEN
415          write(lunout,*) NF_STRERROR(ierr)
416          abort_message='Variable Time not found'
417          CALL abort_gcm(modname,abort_message,ierr)
418        ENDIF
419  endif ! of if (start_file_type.eq."earth")
420  call NF95_PUT_VAR(nid,vID,time,start=(/nb/))
421  WRITE(lunout,*)TRIM(modname)//": Saving for ", nb, time
422
423!--- Rewrite control table (itaufin undefined in dynredem0)
424  var="controle"
425  CALL err(NF90_INQ_VARID(nid,var,vID),"inq",var)
426  CALL err(NF90_GET_VAR(nid,vID,tab_cntrl),"get",var)
427  if (start_file_type=="earth") then
428    tab_cntrl(31) = REAL(itau_dyn + itaufin)
429  else
430    tab_cntrl(31) = 0
431  endif
432  CALL err(NF90_INQ_VARID(nid,var,vID),"inq",var)
433  CALL err(NF90_PUT_VAR(nid,vID,tab_cntrl),"put",var)
434
435!--- Save fields
436  CALL dynredem_write_u(nid,"ucov" ,ucov ,llm, nb)
437  CALL dynredem_write_v(nid,"vcov" ,vcov ,llm, nb)
438  CALL dynredem_write_u(nid,"teta" ,teta ,llm, nb)
439  CALL dynredem_write_u(nid,"masse",masse,llm, nb)
440  CALL dynredem_write_u(nid,"ps"   ,ps   ,1, nb)
441
442!--- Tracers in file "start_trac.nc" (added by Anne)
443  lread_inca=.FALSE.; fil="start_trac.nc"
444  IF(type_trac=='inca') INQUIRE(FILE=fil,EXIST=lread_inca)
445  IF(lread_inca) CALL err(NF90_OPEN(fil,NF90_NOWRITE,nid_trac),"open")
446
447!--- Save tracers
448  IF(nqtot.GE.1) THEN
449    DO iq=1,nqtot
450      var=tname(iq); ierr=-1
451      IF(lread_inca) THEN                  !--- Possibly read from "start_trac.nc"
452        fil="start_trac.nc"
453        ierr=NF90_INQ_VARID(nid_trac,var,vID_trac)
454        dum='inq'; IF(ierr==NF90_NoErr) dum='fnd'
455        WRITE(lunout,*)msg(dum,var)
456
457
458        IF(ierr==NF90_NoErr) CALL dynredem_read_u(nid_trac,var,q(:,:,:,iq),llm)
459      END IF ! of IF(lread_inca)
460      fil=fichnom
461      CALL dynredem_write_u(nid,var,q(:,:,:,iq),llm,nb)
462    END DO ! of DO iq=1,nqtot
463  ENDIF ! of IF(nqtot.GE.1)
464
465  CALL err(NF90_CLOSE(nid),"close")
466  fil="start_trac.nc"
467  IF(lread_inca) CALL err(NF90_CLOSE(nid_trac),"close")
468
469  ENDIF ! of IF (mpi_rank==0)
470
471END SUBROUTINE dynredem1_p
472
Note: See TracBrowser for help on using the repository browser.