source: LMDZ6/trunk/libf/phylmd/iotd_ecrit.f90 @ 5291

Last change on this file since 5291 was 5291, checked in by abarral, 5 hours ago

Move thermcell_old.h iotd.h to module

File size: 5.0 KB
Line 
1      subroutine iotd_ecrit(nom,llm,titre,unite,px)
2
3
4!=======================================================================
5!
6!   Auteur:  F. Hourdin
7!   -------
8!
9!   Objet:
10!   ------
11!   Light interface for netcdf outputs. can be used outside LMDZ
12!
13!=======================================================================
14!-----------------------------------------------------------------------
15!  ----------
16!      nom  : nom de la variable a sortir (chaine de caracteres)
17!      llm  : nombre de couches
18!      titre: titre de la variable (chaine de caracteres)
19!      unite : unite de la variable (chaine de caracteres)
20!      px : variable a sortir
21!
22!=================================================================
23 
24      USE netcdf, ONLY: nf90_put_var, nf90_inq_varid, nf90_enddef, nf90_redef, nf90_sync, nf90_noerr, &
25            nf90_float, nf90_def_var
26      USE iotd_mod_h
27      implicit none
28
29
30! Arguments on input:
31      integer llm
32      character (len=*) :: nom,titre,unite
33      integer imjmax
34      parameter (imjmax=100000)
35      real px(imjmax*llm)
36
37! Local variables:
38
39      real*4 date
40      real*4 zx(imjmax*llm)
41
42
43      integer ierr,ndim,dim_cc(4)
44      integer iq
45      integer i,j,l
46
47      integer zitau
48      character firstnom*20
49      SAVE firstnom
50      SAVE zitau
51      SAVE date
52      data firstnom /'1234567890'/
53      data zitau /0/
54
55! Ajouts
56      integer, save :: ntime=0
57      integer :: idim,varid
58      character (len =50):: fichnom
59      integer, dimension(4) :: id
60      integer, dimension(4) :: edges,corner
61     
62
63
64       if (n_names_iotd_def>0 .and..not.any(names_iotd_def==nom)) return
65!***************************************************************
66! Initialisation of 'firstnom' and create/open the "diagfi.nc" NetCDF file
67! ------------------------------------------------------------------------
68! (Au tout premier appel de la subroutine durant le run.)
69
70
71!--------------------------------------------------------
72! Write the variables to output file if it's time to do so
73!--------------------------------------------------------
74
75
76! Compute/write/extend 'time' coordinate (date given in days)
77! (done every "first call" (at given time level) to writediagfi)
78! Note: date is incremented as 1 step ahead of physics time
79!--------------------------------------------------------
80
81        zx(1:imax*jmax*llm)=px(1:imax*jmax*llm)
82        if (firstnom =='1234567890') then
83            firstnom=nom
84        endif
85
86       !print*,'nom ',nom,firstnom
87
88!! Quand on tombe sur la premiere variable on ajoute un pas de temps
89        if (nom.eq.firstnom) then
90        ! We have identified a "first call" (at given date)
91
92           ntime=ntime+1 ! increment # of stored time steps
93
94!!          print*,'ntime ',ntime
95           date=iotd_t0+ntime*iotd_ts
96           !print*,'iotd_ecrit ',iotd_ts,ntime, date
97!          date= float (zitau +1)/float (day_step)
98
99           ! compute corresponding date (in days and fractions thereof)
100           ! Get NetCDF ID of 'time' variable
101
102           ierr=nf90_sync(nid)
103
104           ierr= nf90_inq_varid(nid,"time",varid)
105           ! Write (append) the new date to the 'time' array
106
107
108           ierr= NF90_PUT_VAR(nid,varid,date,[ntime])
109
110!          print*,'date ',date,ierr,nid
111!        print*,'IOTD Date ,varid,nid,ntime,date',varid,nid,ntime,date
112
113           if (ierr.ne.nf90_noerr) then
114              write(*,*) "***** PUT_VAR matter in writediagfi_nc"
115              write(*,*) "***** with time"
116              write(*,*) 'ierr=', ierr   
117           endif
118
119!          write(6,*)'WRITEDIAGFI: date= ', date
120        end if ! of if (nom.eq.firstnom)
121
122
123!Case of a 3D variable
124!---------------------
125        if (llm==lmax) then
126           ndim=4
127           corner(1)=1
128           corner(2)=1
129           corner(3)=1
130           corner(4)=ntime
131           edges(1)=imax
132           edges(2)=jmax
133           edges(3)=llm
134           edges(4)=1
135           dim_cc=dim_coord
136
137
138!Case of a 2D variable
139!---------------------
140
141        else if (llm==1) then
142           ndim=3
143           corner(1)=1
144           corner(2)=1
145           corner(3)=ntime
146           corner(4)=1
147           edges(1)=imax
148           edges(2)=jmax
149           edges(3)=1
150           edges(4)=1
151           dim_cc(1:2)=dim_coord(1:2)
152           dim_cc(3)=dim_coord(4)
153
154        endif ! of if llm=1 ou llm
155
156! AU premier pas de temps, on crée les variables
157!-----------------------------------------------
158
159      if (ntime==1) then
160          ierr = nf90_redef (nid)
161          ierr = nf90_def_var(nid,nom,nf90_float,dim_cc,varid)
162          !print*,'DEF ',nom,nid,varid
163          ierr = nf90_enddef(nid)
164      else
165         ierr= nf90_inq_varid(nid,nom,varid)
166          !print*,'INQ ',nom,nid,varid
167! Commandes pour recuperer automatiquement les coordonnees
168!             ierr= nf90_inq_dimid(nid,"longitude",id(1))
169      endif
170
171
172      ierr= NF90_PUT_VAR(nid,varid,zx,corner,edges)
173
174      if (ierr.ne.nf90_noerr) then
175           write(*,*) "***** PUT_VAR problem in writediagfi"
176           write(*,*) "***** with ",nom
177           write(*,*) 'ierr=', ierr
178      endif
179
180
181      end
Note: See TracBrowser for help on using the repository browser.