source: LMDZ5/trunk/libf/phylmd/iotd_ecrit.F90 @ 4735

Last change on this file since 4735 was 2239, checked in by Ehouarn Millour, 10 years ago

Reorganizing physics/dynamics interface:

  • what is related to dynamics-physics interface is now in a seperate directory: dynlmdz_phy* for physics in phy*
  • 1d model and related dependencies (including a couple from "dynamics", set up as symbolic links) is now in subdirectory "dyn1d" of phy*.
  • "bibio" directory is now "misc" and should only contain autonomous utilities.
  • "cosp" is now a subdirectory of phylmd.

EM

  • Property svn:executable set to *
File size: 4.7 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      implicit none
25
26! Commons
27
28#include "netcdf.inc"
29#include "iotd.h"
30
31
32! Arguments on input:
33      integer llm
34      character (len=*) :: nom,titre,unite
35      integer imjmax
36      parameter (imjmax=100000)
37      real px(imjmax*llm)
38
39! Local variables:
40
41      real*4 date
42      real*4 zx(imjmax*llm)
43
44
45      integer ierr,ndim,dim_cc(4)
46      integer iq
47      integer i,j,l
48
49      integer zitau
50      character firstnom*20
51      SAVE firstnom
52      SAVE zitau
53      SAVE date
54      data firstnom /'1234567890'/
55      data zitau /0/
56
57! Ajouts
58      integer, save :: ntime=0
59      integer :: idim,varid
60      character (len =50):: fichnom
61      integer, dimension(4) :: id
62      integer, dimension(4) :: edges,corner
63     
64
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=ntime
96!          date= float (zitau +1)/float (day_step)
97
98           ! compute corresponding date (in days and fractions thereof)
99           ! Get NetCDF ID of 'Time' variable
100
101           ierr=NF_SYNC(nid)
102
103           ierr= NF_INQ_VARID(nid,"Time",varid)
104           ! Write (append) the new date to the 'Time' array
105
106
107           ierr= NF_PUT_VARA_REAL(nid,varid,ntime,1,date)
108
109!          print*,'date ',date,ierr,nid
110!        print*,'IOTD Date ,varid,nid,ntime,date',varid,nid,ntime,date
111
112           if (ierr.ne.NF_NOERR) then
113              write(*,*) "***** PUT_VAR matter in writediagfi_nc"
114              write(*,*) "***** with time"
115              write(*,*) 'ierr=', ierr   
116           endif
117
118!          write(6,*)'WRITEDIAGFI: date= ', date
119        end if ! of if (nom.eq.firstnom)
120
121
122!Case of a 3D variable
123!---------------------
124        if (llm==lmax) then
125           ndim=4
126           corner(1)=1
127           corner(2)=1
128           corner(3)=1
129           corner(4)=ntime
130           edges(1)=imax
131           edges(2)=jmax
132           edges(3)=llm
133           edges(4)=1
134           dim_cc=dim_coord
135
136
137!Case of a 2D variable
138!---------------------
139
140        else if (llm==1) then
141           ndim=3
142           corner(1)=1
143           corner(2)=1
144           corner(3)=ntime
145           corner(4)=1
146           edges(1)=imax
147           edges(2)=jmax
148           edges(3)=1
149           edges(4)=1
150           dim_cc(1:2)=dim_coord(1:2)
151           dim_cc(3)=dim_coord(4)
152
153        endif ! of if llm=1 ou llm
154
155! AU premier pas de temps, on crée les variables
156!-----------------------------------------------
157
158      if (ntime==1) then
159          ierr = NF_REDEF (nid)
160          ierr = NF_DEF_VAR(nid,nom,NF_FLOAT,ndim,dim_cc,varid)
161          print*,'DEF ',nom,nid,varid
162          ierr = NF_ENDDEF(nid)
163      else
164         ierr= NF_INQ_VARID(nid,nom,varid)
165          print*,'INQ ',nom,nid,varid
166! Commandes pour recuperer automatiquement les coordonnees
167!             ierr= NF_INQ_DIMID(nid,"longitude",id(1))
168      endif
169
170
171      ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,zx)
172
173      if (ierr.ne.NF_NOERR) then
174           write(*,*) "***** PUT_VAR problem in writediagfi"
175           write(*,*) "***** with ",nom
176           write(*,*) 'ierr=', ierr
177      endif
178
179
180      end
Note: See TracBrowser for help on using the repository browser.