source: trunk/LMDZ.COMMON/libf/evolution/writediagpem.F90 @ 3383

Last change on this file since 3383 was 3345, checked in by jbclement, 6 months ago

PEM:
Making "writediagsoilpem.F90" work in 3D with the correct PEM soil depth dimension.
JBC

File size: 34.4 KB
Line 
1MODULE writediagpem_mod
2
3implicit none
4
5!=======================================================================
6contains
7!=======================================================================
8
9SUBROUTINE writediagpem(ngrid,nom,titre,unite,dim,px)
10
11!  Ecriture de variables diagnostiques au choix dans la physique
12!  dans un fichier NetCDF nomme  'diagpem'. Ces variables peuvent etre
13!  3d (ex : temperature), 2d (ex : temperature de surface), ou
14!  0d (pour un scalaire qui ne depend que du temps : ex : la longitude
15!  solaire)
16!  (ou encore 1d, dans le cas de testphys1d, pour sortir une colonne)
17!  La periode d'ecriture est donnee par
18!  "ecritpem " regle dans le fichier de controle de run :  run.def
19!
20!    writediagpem peut etre appele de n'importe quelle subroutine
21!    de la physique, plusieurs fois. L'initialisation et la creation du
22!    fichier se fait au tout premier appel.
23!
24! WARNING : les variables dynamique (u,v,t,q,ps)
25!  sauvees par writediagpem avec une
26! date donnee sont legerement differentes que dans le fichier histoire car
27! on ne leur a pas encore ajoute de la dissipation et de la physique !!!
28! IL est  RECOMMANDE d'ajouter les tendance physique a ces variables
29! avant l'ecriture dans diagpem (cf. physiq.F)
30!
31! Modifs: Aug.2010 Ehouarn: enforce outputs to be real*4
32!         Oct 2011 Francois: enable having a 'diagpem.def' file to select
33!                            at runtime, which variables to put in file
34!         Oct 2023 JBC: conversion into Fortran 90 with module for the PEM
35!
36!  parametres (input) :
37!  ----------
38!      ngrid : nombres de point ou est calcule la physique
39!                (ngrid = 2+(jjm-1)*iim - 1/jjm)
40!                (= nlon ou klon dans la physique terrestre)
41!      unit : unite logique du fichier de sortie (toujours la meme)
42!      nom  : nom de la variable a sortir (chaine de caracteres)
43!      titre: titre de la variable (chaine de caracteres)
44!      unite : unite de la variable (chaine de caracteres)
45!      px : variable a sortir (real 0, 1, 2, ou 3d)
46!      dim : dimension de px : 0, 1, 2, ou 3 dimensions
47!
48!=================================================================
49use surfdat_h,          only: phisfi
50use geometry_mod,       only: cell_area
51use mod_phys_lmdz_para, only: is_parallel, is_mpi_root, is_master, gather
52use mod_grid_phy_lmdz,  only: klon_glo, Grid1Dto2D_glo, nbp_lon, nbp_lat, nbp_lev, grid_type, unstructured
53use time_evol_mod,      only: ecritpem, dt_pem
54
55implicit none
56
57! Commons
58include "netcdf.inc"
59
60! Arguments on input:
61integer,                        intent(in) :: ngrid
62character(len=*),               intent(in) :: nom, titre, unite
63integer,                        intent(in) :: dim
64real, dimension(ngrid,nbp_lev), intent(in) :: px
65
66! Local variables:
67real*4, dimension(nbp_lon + 1,nbp_lat,nbp_lev) :: dx3    ! to store a 3D data set
68real*4, dimension(nbp_lon + 1,nbp_lat)         :: dx2    ! to store a 2D (surface) data set
69real*4, dimension(nbp_lev)                     :: dx1    ! to store a 1D (column) data set
70real*4                                         :: dx0
71real*4, dimension(1,nbp_lev)                   :: dx3_1d ! to store a profile with 1D model
72real*4                                         :: dx2_1d ! to store a surface value with 1D model
73real*4, save                                   :: date
74
75!$OMP THREADPRIVATE(date)
76real, dimension((nbp_lon + 1),nbp_lat) :: phis
77real, dimension((nbp_lon + 1),nbp_lat) :: area
78integer                                :: irythme, ierr, ierr2, i, j, l, ig0
79integer,                          save :: zitau = 0
80character(27),                    save :: firstnom='1234567890'
81!$OMP THREADPRIVATE(zitau,firstnom)
82
83! Ajouts
84integer,           save :: ntime = 0
85!$OMP THREADPRIVATE(ntime)
86integer                 :: idim, varid
87integer                 :: nid
88character(*), parameter :: fichnom = "diagpem.nc"
89integer, dimension(4)   :: id
90integer, dimension(4)   :: edges, corner
91
92! Added to use diagpem.def to select output variable
93logical,                                  save :: diagpem_def
94logical                                        :: getout
95integer,                                  save :: n_nom_def
96integer                                        :: n
97integer, parameter                             :: n_nom_def_max = 199
98character(120), dimension(n_nom_def_max), save :: nom_def
99logical,                                  save :: firstcall = .true.
100!$OMP THREADPRIVATE(firstcall) !diagpem_def,n_nom_def,nom_def read in diagpem.def
101
102#ifdef CPP_PARA
103! Added to work in parallel mode
104    real, dimension(klon_glo,nbp_lev)        :: dx3_glop
105    real, dimension(nbp_lon,nbp_lat,nbp_lev) :: dx3_glo  ! to store a global 3D data set
106    real, dimension(klon_glo)                :: dx2_glop
107    real, dimension(nbp_lon,nbp_lat)         :: dx2_glo  ! to store a global 2D (surface) data set
108    real, dimension(ngrid)                   :: px2
109!    real, dimension(nbp_lev)                :: dx1_glo   ! to store a 1D (column) data set
110!    real :: dx0_glo
111    real, dimension(klon_glo) :: phisfi_glo ! surface geopotential on global physics grid
112    real, dimension(klon_glo) :: areafi_glo ! mesh area on global physics grid
113#else
114    real, dimension(ngrid) :: phisfi_glo ! surface geopotential on global physics grid
115    real, dimension(ngrid) :: areafi_glo ! mesh area on global physics grid
116#endif
117
118if (grid_type == unstructured) return
119
120!***************************************************************
121!Sortie des variables au rythme voulu
122irythme = int(ecritpem) ! output rate set by ecritpem
123
124!***************************************************************
125! At very first call, check if there is a "diagpem.def" to use and read it
126! ------------------------------------------------------------------------
127IF (firstcall) THEN
128    firstcall=.false.
129
130!$OMP MASTER
131    ! Open diagpem.def definition file if there is one:
132        open(99,file="diagpem.def",status='old',form='formatted',iostat=ierr2)
133
134        if (ierr2 == 0) then
135            diagpem_def=.true.
136            write(*,*) "*******************"
137            write(*,*) "Reading diagpem.def"
138            write(*,*) "*******************"
139            do n=1,n_nom_def_max
140                read(99,fmt='(a)',end=88) nom_def(n)
141                write(*,*) 'Output in diagpem: ', trim(nom_def(n))
142            enddo
143 88         continue
144            if (n.ge.n_nom_def_max) then
145                write(*,*)"n_nom_def_max too small in writediagpem.F:",n
146                call abort_physic("writediagpem","n_nom_def_max too small",1)
147            end if
148            n_nom_def=n-1
149            close(99)
150        else
151            diagpem_def=.false.
152        endif
153!$OMP END MASTER
154!$OMP BARRIER
155ENDIF ! of IF (firstcall)
156
157! Get out of write_diagpem if there is diagpem.def AND variable not listed
158!  -----------------------------------------------------------------------
159if (diagpem_def) then
160    getout=.true.
161    do n=1,n_nom_def
162        if (trim(nom_def(n)).eq.nom) getout=.false.
163    enddo
164    if (getout) return
165endif
166
167! Initialisation of 'firstnom' and create/open the "diagpem.nc" NetCDF file
168! -------------------------------------------------------------------------
169! (at very first call to the subroutine, in accordance with diagpem.def)
170if (firstnom.eq.'1234567890') then ! .true. for the very first valid
171    ! call to this subroutine; now set 'firstnom'
172    firstnom = nom
173    ! just to be sure, check that firstnom is large enough to hold nom
174    if (len_trim(firstnom).lt.len_trim(nom)) then
175        write(*,*) "writediagpem: Error !!!"
176        write(*,*) "   firstnom string not long enough!!"
177        write(*,*) "   increase its size to at least ",len_trim(nom)
178        call abort_physic("writediagpem","firstnom too short",1)
179    endif
180
181    zitau = -1 ! initialize zitau
182
183#ifdef CPP_PARA
184        ! Gather phisfi() geopotential on physics grid
185        call Gather(phisfi,phisfi_glo)
186        ! Gather cell_area() mesh area on physics grid
187        call Gather(cell_area,areafi_glo)
188#else
189        phisfi_glo(:)=phisfi(:)
190        areafi_glo(:)=cell_area(:)
191#endif
192
193    !! parallel: we cannot use the usual writediagpem method
194    !! call iophys_ini
195    if (is_master) then
196        ! only the master is required to do this
197        ! Create the NetCDF file
198        ierr = NF_CREATE(fichnom,NF_CLOBBER,nid)
199        ! Define the 'Time' dimension
200        ierr = nf_def_dim(nid,"Time",NF_UNLIMITED,idim)
201        ! Define the 'Time' variable
202    !#ifdef NC_DOUBLE
203    !    ierr = NF_DEF_VAR (nid, "Time", NF_DOUBLE, 1, idim,varid)
204    !#else
205        ierr = NF_DEF_VAR (nid, "Time", NF_FLOAT, 1, idim,varid)
206    !#endif
207        ! Add a long_name attribute
208        ierr = NF_PUT_ATT_TEXT (nid, varid, "long_name",4,"Time")
209        ! Add a units attribute
210        ierr = NF_PUT_ATT_TEXT(nid, varid,'units',29,"days since 0000-00-0 00:00:00")
211        ! Switch out of NetCDF Define mode
212        ierr = NF_ENDDEF(nid)
213
214        ! Build phis() and area()
215        IF (klon_glo>1) THEN
216            do i=1,nbp_lon+1 ! poles
217            phis(i,1)=phisfi_glo(1)
218            phis(i,nbp_lat)=phisfi_glo(klon_glo)
219            ! for area, divide at the poles by nbp_lon
220            area(i,1)=areafi_glo(1)/nbp_lon
221            area(i,nbp_lat)=areafi_glo(klon_glo)/nbp_lon
222            enddo
223            do j=2,nbp_lat-1
224            ig0= 1+(j-2)*nbp_lon
225            do i=1,nbp_lon
226                phis(i,j)=phisfi_glo(ig0+i)
227                area(i,j)=areafi_glo(ig0+i)
228            enddo
229            ! handle redundant point in longitude
230            phis(nbp_lon+1,j)=phis(1,j)
231            area(nbp_lon+1,j)=area(1,j)
232            enddo
233        ENDIF
234
235        ! write "header" of file (longitudes, latitudes, geopotential, ...)
236        IF (klon_glo>1) THEN ! general 3D case
237            call iniwrite(nid,0,phis,area,nbp_lon+1,nbp_lat)
238        ELSE ! 1D model
239            call iniwrite(nid,0,phisfi_glo(1),areafi_glo(1),1,1)
240        ENDIF
241
242    endif ! of if (is_master)
243
244    else
245
246    if (is_master) then
247        ! only the master is required to do this
248        ! Open the NetCDF file
249        ierr = NF_OPEN(fichnom,NF_WRITE,nid)
250    endif ! of if (is_master)
251
252endif ! if (firstnom.eq.'1234567890')
253
254! Increment time index 'zitau' if it is the "fist call" (at given time level)
255! to writediagpem
256!------------------------------------------------------------------------
257if (nom == firstnom) zitau = zitau + dt_pem
258
259!--------------------------------------------------------
260! Write the variables to output file if it's time to do so
261!--------------------------------------------------------
262if (MOD(zitau+1,irythme) == 0.) then
263! Compute/write/extend 'Time' coordinate (date given in days)
264! (done every "first call" (at given time level) to writediagpem)
265! Note: date is incremented as 1 step ahead of physics time
266!--------------------------------------------------------
267    if (is_master) then
268        ! only the master is required to do this
269        if (nom.eq.firstnom) then
270        ! We have identified a "first call" (at given date)
271            ntime=ntime+1 ! increment # of stored time steps
272            ! compute corresponding date (in days and fractions thereof)
273            date = float(zitau + 1)
274            ! Get NetCDF ID of 'Time' variable
275            ierr= NF_INQ_VARID(nid,"Time",varid)
276            ! Write (append) the new date to the 'Time' array
277!#ifdef NC_DOUBLE
278!    ierr= NF_PUT_VARA_DOUBLE(nid,varid,[ntime],[1],[date])
279!#else
280            ierr= NF_PUT_VARA_REAL(nid,varid,[ntime],[1],[date])
281!#endif
282            if (ierr.ne.NF_NOERR) then
283                write(*,*) "***** PUT_VAR matter in writediagpem_nc"
284                write(*,*) "***** with time"
285                write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr)
286!                call abort
287            endif
288
289            write(6,*)'WRITEDIAGPEM: date= ', date
290        endif ! of if (nom.eq.firstnom)
291    endif ! of if (is_master)
292
293!Case of a 3D variable
294!---------------------
295    if (dim == 3) then
296        IF (klon_glo>1) THEN ! General case
297#ifdef CPP_PARA
298        ! Gather field on a "global" (without redundant longitude) array
299        call Gather(px,dx3_glop)
300!$OMP MASTER
301        if (is_mpi_root) then
302            call Grid1Dto2D_glo(dx3_glop,dx3_glo)
303            ! copy dx3_glo() to dx3(:) and add redundant longitude
304            dx3(1:nbp_lon,:,:)=dx3_glo(1:nbp_lon,:,:)
305            dx3(nbp_lon+1,:,:)=dx3(1,:,:)
306        endif
307!$OMP END MASTER
308!$OMP BARRIER
309#else
310!         Passage variable physique -->  variable dynamique
311!         recast (copy) variable from physics grid to dynamics grid
312            DO l=1,nbp_lev
313                DO i=1,nbp_lon+1
314                    dx3(i,1,l)=px(1,l)
315                    dx3(i,nbp_lat,l)=px(ngrid,l)
316                ENDDO
317                DO j=2,nbp_lat-1
318                    ig0= 1+(j-2)*nbp_lon
319                    DO i=1,nbp_lon
320                        dx3(i,j,l)=px(ig0+i,l)
321                    ENDDO
322                    dx3(nbp_lon+1,j,l)=dx3(1,j,l)
323                ENDDO
324            ENDDO
325#endif
326        ELSE ! 1D model case
327            dx3_1d(1,1:nbp_lev)=px(1,1:nbp_lev)
328        ENDIF
329!       Ecriture du champs
330        if (is_master) then
331! only the master writes to output
332! name of the variable
333            ierr= NF_INQ_VARID(nid,nom,varid)
334            if (ierr /= NF_NOERR) then
335! corresponding dimensions
336                ierr = NF_INQ_DIMID(nid,"longitude",id(1))
337                ierr = NF_INQ_DIMID(nid,"latitude",id(2))
338                ierr = NF_INQ_DIMID(nid,"altitude",id(3))
339                ierr = NF_INQ_DIMID(nid,"Time",id(4))
340
341! Create the variable if it doesn't exist yet
342                write (*,*) "==========================="
343                write (*,*) "DIAGPEM: creating variable ",trim(nom)
344                call def_var(nid,nom,titre,unite,4,id,varid,ierr)
345
346            else
347                if (ntime==0) then
348                    write(*,*) "DIAGPEM Error: failed creating variable ",trim(nom)
349                    write(*,*) "it seems it already exists!"
350                    call abort_physic("writediagpem",trim(nom)//" already exists",1)
351                endif
352            endif
353
354            corner(1)=1
355            corner(2)=1
356            corner(3)=1
357            corner(4)=ntime
358
359            IF (klon_glo==1) THEN
360                edges(1)=1
361            ELSE
362                edges(1)=nbp_lon+1
363            ENDIF
364            edges(2)=nbp_lat
365            edges(3)=nbp_lev
366            edges(4)=1
367!#ifdef NC_DOUBLE
368!            ierr= NF_PUT_VARA_DOUBLE(nid,varid,corner,edges,dx3)
369!#else
370!            write(*,*)"test:  nid=",nid," varid=",varid
371!            write(*,*)"       corner()=",corner
372!            write(*,*)"       edges()=",edges
373!            write(*,*)"       dx3()=",dx3
374            IF (klon_glo>1) THEN ! General case
375                ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx3)
376            ELSE
377                ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx3_1d)
378            ENDIF
379!#endif
380
381            if (ierr.ne.NF_NOERR) then
382                write(*,*) "***** PUT_VAR problem in writediagpem"
383                write(*,*) "***** with dx3: ",trim(nom)
384                write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr)
385                call abort_physic("writediagpem","failed writing "//trim(nom),1)
386            endif
387
388        endif !of if (is_master)
389
390!Case of a 2D variable
391!---------------------
392    else if (dim == 2) then
393
394        IF (klon_glo>1) THEN ! General case
395#ifdef CPP_PARA
396        ! Gather field on a "global" (without redundant longitude) array
397        px2(:)=px(:,1)
398        call Gather(px2,dx2_glop)
399!$OMP MASTER
400        if (is_mpi_root) then
401            call Grid1Dto2D_glo(dx2_glop,dx2_glo)
402            ! copy dx2_glo() to dx2(:) and add redundant longitude
403            dx2(1:nbp_lon,:)=dx2_glo(1:nbp_lon,:)
404            dx2(nbp_lon+1,:)=dx2(1,:)
405        endif
406!$OMP END MASTER
407!$OMP BARRIER
408#else
409!         Passage variable physique -->  physique dynamique
410!         recast (copy) variable from physics grid to dynamics grid
411            DO i=1,nbp_lon+1
412                dx2(i,1)=px(1,1)
413                dx2(i,nbp_lat)=px(ngrid,1)
414            ENDDO
415            DO j=2,nbp_lat-1
416                ig0= 1+(j-2)*nbp_lon
417                DO i=1,nbp_lon
418                    dx2(i,j)=px(ig0+i,1)
419                ENDDO
420                dx2(nbp_lon+1,j)=dx2(1,j)
421            ENDDO
422#endif
423        ELSE ! 1D model case
424            dx2_1d=px(1,1)
425        ENDIF
426
427        if (is_master) then
428        ! only the master writes to output
429!            write (*,*) 'In  writediagpem, on sauve:  ' , nom
430!            write (*,*) 'In  writediagpem. Estimated date = ' ,date
431            ierr= NF_INQ_VARID(nid,nom,varid)
432            if (ierr /= NF_NOERR) then
433! corresponding dimensions
434                ierr= NF_INQ_DIMID(nid,"longitude",id(1))
435                ierr= NF_INQ_DIMID(nid,"latitude",id(2))
436                ierr= NF_INQ_DIMID(nid,"Time",id(3))
437
438! Create the variable if it doesn't exist yet
439                write (*,*) "==========================="
440                write (*,*) "DIAGPEM: creating variable ",trim(nom)
441
442                call def_var(nid,nom,titre,unite,3,id,varid,ierr)
443
444            else
445                if (ntime==0) then
446                    write(*,*) "DIAGPEM Error: failed creating variable ",trim(nom)
447                    write(*,*) "it seems it already exists!"
448                    call abort_physic("writediagpem",trim(nom)//" already exists",1)
449                endif
450            endif
451
452            corner(1)=1
453            corner(2)=1
454            corner(3)=ntime
455            IF (klon_glo==1) THEN
456                edges(1)=1
457            ELSE
458                edges(1)=nbp_lon+1
459            ENDIF
460            edges(2)=nbp_lat
461            edges(3)=1
462
463
464!#ifdef NC_DOUBLE
465!           ierr = NF_PUT_VARA_DOUBLE (nid,varid,corner,edges,dx2)
466!#else
467            IF (klon_glo>1) THEN ! General case
468                ierr = NF_PUT_VARA_REAL(nid,varid,corner,edges,dx2)
469            ELSE
470                ierr = NF_PUT_VARA_REAL(nid,varid,corner,edges,[dx2_1d])
471            ENDIF
472!#endif
473
474            if (ierr.ne.NF_NOERR) then
475                write(*,*) "***** PUT_VAR matter in writediagpem"
476                write(*,*) "***** with dx2: ",trim(nom)
477                write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr)
478                call abort_physic("writediagpem","failed writing "//trim(nom),1)
479            endif
480
481        endif !of if (is_master)
482
483!Case of a 1D variable (ie: a column)
484!---------------------------------------------------
485    else if (dim.eq.1) then
486        if (is_parallel) then
487            write(*,*) "writediagpem error: dim=1 not implemented ","in parallel mode. Problem for ",trim(nom)
488                call abort_physic("writediagpem","failed writing "//trim(nom),1)
489        endif
490!       Passage variable physique -->  physique dynamique
491!       recast (copy) variable from physics grid to dynamics grid
492        do l=1,nbp_lev
493            dx1(l)=px(1,l)
494        enddo
495
496        ierr= NF_INQ_VARID(nid,nom,varid)
497        if (ierr /= NF_NOERR) then
498! corresponding dimensions
499            ierr= NF_INQ_DIMID(nid,"altitude",id(1))
500            ierr= NF_INQ_DIMID(nid,"Time",id(2))
501
502! Create the variable if it doesn't exist yet
503
504            write (*,*) "==========================="
505            write (*,*) "DIAGPEM: creating variable ",trim(nom)
506
507            call def_var(nid,nom,titre,unite,2,id,varid,ierr)
508
509        else
510            if (ntime==0) then
511                write(*,*) "DIAGPEM Error: failed creating variable ",trim(nom)
512                write(*,*) "it seems it already exists!"
513                call abort_physic("writediagpem",trim(nom)//" already exists",1)
514            endif
515        endif
516
517        corner(1)=1
518        corner(2)=ntime
519        edges(1)=nbp_lev
520        edges(2)=1
521!#ifdef NC_DOUBLE
522!           ierr= NF_PUT_VARA_DOUBLE(nid,varid,corner,edges,dx1)
523!#else
524        ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx1)
525!#endif
526
527        if (ierr /= NF_NOERR) then
528            write(*,*) "***** PUT_VAR problem in writediagpem"
529            write(*,*) "***** with dx1: ",trim(nom)
530            write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr)
531            call abort_physic("writediagpem","failed writing "//trim(nom),1)
532        endif
533
534!Case of a 0D variable (ie: a time-dependent scalar)
535!---------------------------------------------------
536    else if (dim == 0) then
537
538        dx0 = px (1,1)
539
540        if (is_master) then
541            ! only the master writes to output
542            ierr= NF_INQ_VARID(nid,nom,varid)
543            if (ierr /= NF_NOERR) then
544! corresponding dimensions
545                ierr= NF_INQ_DIMID(nid,"Time",id(1))
546
547! Create the variable if it doesn't exist yet
548                write (*,*) "==========================="
549                write (*,*) "DIAGPEM: creating variable ",trim(nom)
550                call def_var(nid,nom,titre,unite,1,id,varid,ierr)
551
552            else
553                if (ntime==0) then
554                    write(*,*) "DIAGPEM Error: failed creating variable ",trim(nom)
555                    write(*,*) "it seems it already exists!"
556                    call abort_physic("writediagpem",trim(nom)//" already exists",1)
557                endif
558            endif
559
560            corner(1)=ntime
561            edges(1)=1
562
563!#ifdef NC_DOUBLE
564!            ierr = NF_PUT_VARA_DOUBLE (nid,varid,[corner(1)],[edges(1)],[dx0])
565!#else
566            ierr= NF_PUT_VARA_REAL(nid,varid,[corner(1)],[edges(1)],[dx0])
567!#endif
568            if (ierr.ne.NF_NOERR) then
569                write(*,*) "***** PUT_VAR matter in writediagpem"
570                write(*,*) "***** with dx0: ",trim(nom)
571                write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr)
572                call abort_physic("writediagpem","failed writing "//trim(nom),1)
573            endif
574        endif !of if (is_master)
575    endif ! of if (dim.eq.3) elseif(dim.eq.2)...
576endif ! of if ( MOD(zitau+1,irythme) .eq.0.)
577
578if (is_master) ierr= NF_CLOSE(nid)
579
580END SUBROUTINE writediagpem
581
582!=================================================================
583
584SUBROUTINE writediagsoilpem(ngrid,name,title,units,dimpx,px)
585
586! Write variable 'name' to NetCDF file 'diagsoilpem.nc'.
587! The variable may be 3D (lon,lat,depth) subterranean field,
588! a 2D (lon,lat) surface field, or a simple scalar (0D variable).
589!
590! Calls to 'writediagsoilpem' can originate from anywhere in the program;
591! An initialisation of variable 'name' is done if it is the first time
592! that this routine is called with given 'name'; otherwise data is appended
593! (yielding the sought time series of the variable)
594
595! Modifs: Aug.2010 Ehouarn: enforce outputs to be real*4
596
597use comsoil_h_PEM,      only: mlayer_PEM, nsoilmx_PEM, inertiedat_PEM
598use geometry_mod,       only: cell_area
599use mod_phys_lmdz_para, only: is_mpi_root, is_master, gather
600use mod_grid_phy_lmdz,  only: klon_glo, Grid1Dto2D_glo, nbp_lon, nbp_lat
601use mod_grid_phy_lmdz,  only: grid_type, unstructured
602use time_evol_mod,      only: ecritpem, dt_pem
603use iniwritesoil_mod,   only: iniwritesoil
604
605implicit none
606
607include"netcdf.inc"
608
609! Arguments:
610integer,intent(in) :: ngrid ! number of (horizontal) points of physics grid
611! i.e. ngrid = 2+(jjm-1)*iim - 1/jjm
612character(len=*),intent(in) :: name ! 'name' of the variable
613character(len=*),intent(in) :: title ! 'long_name' attribute of the variable
614character(len=*),intent(in) :: units ! 'units' attribute of the variable
615integer,intent(in) :: dimpx ! dimension of the variable (3,2 or 0)
616real,dimension(ngrid,nsoilmx_PEM),intent(in) :: px ! variable
617
618! Local variables:
619real*4,dimension(nbp_lon+1,nbp_lat,nsoilmx_PEM) :: data3 ! to store 3D data
620real*4,dimension(nbp_lon+1,nbp_lat) :: data2 ! to store 2D data
621real*4 :: data0 ! to store 0D data
622real*4 :: data3_1d(1,nsoilmx_PEM) ! to store a profile in 1D model
623real*4 :: data2_1d ! to store surface value with 1D model
624integer :: i,j,l ! for loops
625integer :: ig0
626
627real*4,save :: date ! time counter (in elapsed days)
628
629real :: inertia((nbp_lon+1),nbp_lat,nsoilmx_PEM)
630real :: area((nbp_lon+1),nbp_lat)
631
632real :: inertiafi_glo(klon_glo,nsoilmx_PEM)
633real :: areafi_glo(klon_glo)
634
635integer,save :: isample ! sample rate at which data is to be written to output
636integer,save :: ntime=0 ! counter to internally store time steps
637character(len=20),save :: firstname="1234567890"
638integer,save :: zitau=0
639!$OMP THREADPRIVATE(date,isample,ntime,firstname,zitau)
640
641character(len=30) :: filename="diagsoilpem.nc"
642
643! NetCDF stuff:
644integer :: nid ! NetCDF output file ID
645integer :: varid ! NetCDF ID of a variable
646integer :: ierr ! NetCDF routines return code
647integer,dimension(4) :: id ! NetCDF IDs of the dimensions of the variable
648integer,dimension(4) :: edges,corners
649
650#ifdef CPP_PARA
651! Added to work in parallel mode
652real dx3_glop(klon_glo,nsoilmx_PEM)
653real dx3_glo(nbp_lon,nbp_lat,nsoilmx_PEM) ! to store a global 3D data set
654real dx2_glop(klon_glo)
655real dx2_glo(nbp_lon,nbp_lat)     ! to store a global 2D (surface) data set
656real px2(ngrid)
657#endif
658
659! 0. Do we ouput a diagsoilpem.nc file? If not just bail out now.
660
661! additional check: one can only output diagsoilpem.nc files
662! in lon-lat case (or 1D)
663if (grid_type==unstructured) then
664  write(*,*) "writediagsoilpem: Error !!!"
665  write(*,*) "diagsoilpem.nc outputs not possible on unstructured grids!!"
666  call abort_physic("writediagsoilpem","impossible on unstructured grid",1)
667endif
668
669! 1. Initialization step
670if (firstname.eq."1234567890") then
671  ! Store 'name' as 'firstname'
672  firstname=name
673  ! From now on, if 'name'.eq.'firstname', then it is a new time cycle
674
675  ! just to be sure, check that firstnom is large enough to hold nom
676  if (len_trim(firstname).lt.len_trim(name)) then
677    write(*,*) "writediagsoilpem: Error !!!"
678    write(*,*) "   firstname string not long enough!!"
679    write(*,*) "   increase its size to at least ",len_trim(name)
680    call abort_physic("writediagsoilpem","firstname too short",1)
681  endif
682 
683  ! Set output sample rate
684  isample = int(ecritpem) ! same as for diagpem outputs
685 
686  ! Create output NetCDF file
687  if (is_master) then
688   ierr=NF_CREATE(filename,NF_CLOBBER,nid)
689   if (ierr.ne.NF_NOERR) then
690    write(*,*)'writediagsoilpem: Error, failed creating file '//trim(filename)
691    call abort_physic("writediagsoilpem","failed creating"//trim(filename),1)
692   endif
693  endif
694
695#ifdef CPP_PARA
696   ! Gather inertiedat_PEM() soil thermal inertia on physics grid
697   call Gather(inertiedat_PEM,inertiafi_glo)
698   ! Gather cell_area() mesh area on physics grid
699   call Gather(cell_area,areafi_glo)
700#else
701         inertiafi_glo(:,:)=inertiedat_PEM(:,:)
702         areafi_glo(:)=cell_area(:)
703#endif
704
705  if (is_master) then
706   ! build inertia() and area()
707   if (klon_glo>1) then
708    do i=1,nbp_lon+1 ! poles
709     inertia(i,1,1:nsoilmx_PEM)=inertiafi_glo(1,1:nsoilmx_PEM)
710     inertia(i,nbp_lat,1:nsoilmx_PEM)=inertiafi_glo(klon_glo,1:nsoilmx_PEM)
711     ! for area, divide at the poles by nbp_lon
712     area(i,1)=areafi_glo(1)/nbp_lon
713     area(i,nbp_lat)=areafi_glo(klon_glo)/nbp_lon
714    enddo
715    do j=2,nbp_lat-1
716     ig0= 1+(j-2)*nbp_lon
717     do i=1,nbp_lon
718        inertia(i,j,1:nsoilmx_PEM)=inertiafi_glo(ig0+i,1:nsoilmx_PEM)
719        area(i,j)=areafi_glo(ig0+i)
720     enddo
721     ! handle redundant point in longitude
722     inertia(nbp_lon+1,j,1:nsoilmx_PEM)=inertia(1,j,1:nsoilmx_PEM)
723     area(nbp_lon+1,j)=area(1,j)
724    enddo
725   endif
726   
727   ! write "header" of file (longitudes, latitudes, geopotential, ...)
728   if (klon_glo>1) then ! general 3D case
729     call iniwritesoil(nid,ngrid,inertia,area,nbp_lon+1,nbp_lat,nsoilmx_PEM,mlayer_PEM)
730   else ! 1D model
731     call iniwritesoil(nid,ngrid,inertiafi_glo(1,:),areafi_glo(1),1,1,nsoilmx_PEM,mlayer_PEM)
732   endif
733
734  endif ! of if (is_master)
735 
736  ! set zitau to -1 to be compatible with zitau incrementation step below
737  zitau=-1
738 
739else
740  ! If not an initialization call, simply open the NetCDF file
741  if (is_master) then
742   ierr=NF_OPEN(filename,NF_WRITE,nid)
743  endif
744endif ! of if (firstname.eq."1234567890")
745
746! 2. Increment local time counter, if necessary
747if (name.eq.firstname) then
748  ! if we run across 'firstname', then it is a new time step
749  zitau = zitau + dt_pem
750endif
751
752! 3. Write data, if the time index matches the sample rate
753if (mod(zitau+1,isample).eq.0) then
754
755! 3.1 If first call at this date, update 'time' variable
756  if (name.eq.firstname) then
757    ntime=ntime+1
758    date = float(zitau + 1)
759   
760    if (is_master) then
761     ! Get NetCDF ID for "time"
762     ierr=NF_INQ_VARID(nid,"time",varid)
763     ! Add the current value of date to the "time" array
764!#ifdef NC_DOUBLE
765!     ierr=NF_PUT_VARA_DOUBLE(nid,varid,ntime,1,date)
766!#else
767     ierr=NF_PUT_VARA_REAL(nid,varid,[ntime],[1],[date])
768!#endif
769     if (ierr.ne.NF_NOERR) then
770      write(*,*)"writediagsoilpem: Failed writing date to time variable"
771      call abort_physic("writediagsoilpem","failed writing time",1)
772     endif
773    endif ! of if (is_master)
774  endif ! of if (name.eq.firstname)
775
776! 3.2 Write the variable to the NetCDF file
777if (dimpx.eq.3) then ! Case of a 3D variable
778  ! A. Recast data along 'dynamics' grid
779  if (klon_glo>1) then ! General case
780#ifdef CPP_PARA
781  ! gather field on a "global" (without redundant longitude) array
782  call Gather(px,dx3_glop)
783!$OMP MASTER
784  if (is_mpi_root) then
785    call Grid1Dto2D_glo(dx3_glop,dx3_glo)
786    ! copy dx3_glo() to dx3(:) and add redundant longitude
787    data3(1:nbp_lon,:,:)=dx3_glo(1:nbp_lon,:,:)
788    data3(nbp_lon+1,:,:)=data3(1,:,:)
789  endif
790!$OMP END MASTER
791!$OMP BARRIER
792#else
793   do l=1,nsoilmx_PEM
794    ! handle the poles
795    do i=1,nbp_lon+1
796      data3(i,1,l)=px(1,l)
797      data3(i,nbp_lat,l)=px(ngrid,l)
798    enddo
799    ! rest of the grid
800    do j=2,nbp_lat-1
801      ig0=1+(j-2)*nbp_lon
802      do i=1,nbp_lon
803        data3(i,j,l)=px(ig0+i,l)
804      enddo
805      data3(nbp_lon+1,j,l)=data3(1,j,l) ! extra (modulo) longitude
806    enddo
807   enddo
808#endif
809  else ! 1D model case
810   data3_1d(1,1:nsoilmx_PEM)=px(1,1:nsoilmx_PEM)
811  endif
812 
813  ! B. Write (append) the variable to the NetCDF file
814  if (is_master) then
815  ! B.1. Get the ID of the variable
816  ierr=NF_INQ_VARID(nid,name,varid)
817  if (ierr.ne.NF_NOERR) then
818    ! If we failed geting the variable's ID, we assume it is because
819    ! the variable doesn't exist yet and must be created.
820    ! Start by obtaining corresponding dimensions IDs
821    ierr=NF_INQ_DIMID(nid,"longitude",id(1))
822    ierr=NF_INQ_DIMID(nid,"latitude",id(2))
823    ierr=NF_INQ_DIMID(nid,"depth",id(3))
824    ierr=NF_INQ_DIMID(nid,"time",id(4))
825    ! Tell the world about it
826    write(*,*) "====================="
827    write(*,*) "writediagsoilpem: creating variable "//trim(name)
828    call def_var(nid,name,title,units,4,id,varid,ierr)
829  endif ! of if (ierr.ne.NF_NOERR)
830 
831  ! B.2. Prepare things to be able to write/append the variable
832  corners(1)=1
833  corners(2)=1
834  corners(3)=1
835  corners(4)=ntime
836 
837  if (klon_glo==1) then
838    edges(1)=1
839  else
840    edges(1)=nbp_lon+1
841  endif
842  edges(2)=nbp_lat
843  edges(3)=nsoilmx_PEM
844  edges(4)=1
845 
846  ! B.3. Write the slab of data
847!#ifdef NC_DOUBLE
848!  ierr=NF_PUT_VARA_DOUBLE(nid,varid,corners,edges,data3)
849!#else
850  if (klon_glo>1) then
851    ierr=NF_PUT_VARA_REAL(nid,varid,corners,edges,data3)
852  else
853    ierr=NF_PUT_VARA_REAL(nid,varid,corners,edges,data3_1d)
854  endif
855!#endif
856  if (ierr.ne.NF_NOERR) then
857    write(*,*) "writediagsoilpem: Error: Failed writing "//trim(name)//&
858               " to file "//trim(filename)//" at time",date
859  endif
860  endif ! of if (is_master)
861
862elseif (dimpx.eq.2) then ! Case of a 2D variable
863
864  ! A. Recast data along 'dynamics' grid
865  if (klon_glo>1) then ! General case
866#ifdef CPP_PARA
867  ! gather field on a "global" (without redundant longitude) array
868  px2(:)=px(:,1)
869  call Gather(px2,dx2_glop)
870!$OMP MASTER
871  if (is_mpi_root) then
872    call Grid1Dto2D_glo(dx2_glop,dx2_glo)
873    ! copy dx3_glo() to dx3(:) and add redundant longitude
874    data2(1:nbp_lon,:)=dx2_glo(1:nbp_lon,:)
875    data2(nbp_lon+1,:)=data2(1,:)
876  endif
877!$OMP END MASTER
878!$OMP BARRIER
879#else
880    ! handle the poles
881    do i=1,nbp_lon+1
882      data2(i,1)=px(1,1)
883      data2(i,nbp_lat)=px(ngrid,1)
884    enddo
885    ! rest of the grid
886    do j=2,nbp_lat-1
887      ig0=1+(j-2)*nbp_lon
888      do i=1,nbp_lon
889        data2(i,j)=px(ig0+i,1)
890      enddo
891      data2(nbp_lon+1,j)=data2(1,j) ! extra (modulo) longitude
892    enddo
893#endif
894  else ! 1D model case
895    data2_1d=px(1,1)
896  endif
897
898  ! B. Write (append) the variable to the NetCDF file
899  if (is_master) then
900  ! B.1. Get the ID of the variable
901  ierr=NF_INQ_VARID(nid,name,varid)
902  if (ierr.ne.NF_NOERR) then
903    ! If we failed geting the variable's ID, we assume it is because
904    ! the variable doesn't exist yet and must be created.
905    ! Start by obtaining corresponding dimensions IDs
906    ierr=NF_INQ_DIMID(nid,"longitude",id(1))
907    ierr=NF_INQ_DIMID(nid,"latitude",id(2))
908    ierr=NF_INQ_DIMID(nid,"time",id(3))
909    ! Tell the world about it
910    write(*,*) "====================="
911    write(*,*) "writediagsoilpem: creating variable "//trim(name)
912    call def_var(nid,name,title,units,3,id,varid,ierr)
913  endif ! of if (ierr.ne.NF_NOERR)
914
915  ! B.2. Prepare things to be able to write/append the variable
916  corners(1)=1
917  corners(2)=1
918  corners(3)=ntime
919 
920  if (klon_glo==1) then
921    edges(1)=1
922  else
923    edges(1)=nbp_lon+1
924  endif
925  edges(2)=nbp_lat
926  edges(3)=1
927 
928  ! B.3. Write the slab of data
929!#ifdef NC_DOUBLE
930!  ierr=NF_PUT_VARA_DOUBLE(nid,varid,corners,edges,data2)
931!#else
932  if (klon_glo>1) then ! General case
933    ierr=NF_PUT_VARA_REAL(nid,varid,corners,edges,data2)
934  else
935    ierr=NF_PUT_VARA_REAL(nid,varid,corners,edges,[data2_1d])
936  endif
937!#endif
938  if (ierr.ne.NF_NOERR) then
939    write(*,*) "writediagsoilpem: Error: Failed writing "//trim(name)//&
940               " to file "//trim(filename)//" at time",date
941  endif
942  endif ! of if (is_master)
943
944elseif (dimpx.eq.0) then ! Case of a 0D variable
945#ifdef CPP_PARA
946  write(*,*) "writediagsoilpem: dimps==0 case not implemented in // mode!!"
947  call abort_physic("writediagsoilpem","dimps==0 not implemented",1)
948#endif
949  ! A. Copy data value
950  data0=px(1,1)
951
952  ! B. Write (append) the variable to the NetCDF file
953  ! B.1. Get the ID of the variable
954  ierr=NF_INQ_VARID(nid,name,varid)
955  if (ierr.ne.NF_NOERR) then
956    ! If we failed geting the variable's ID, we assume it is because
957    ! the variable doesn't exist yet and must be created.
958    ! Start by obtaining corresponding dimensions IDs
959    ierr=NF_INQ_DIMID(nid,"time",id(1))
960    ! Tell the world about it
961    write(*,*) "====================="
962    write(*,*) "writediagsoilpem: creating variable "//trim(name)
963    call def_var(nid,name,title,units,1,id,varid,ierr)
964  endif ! of if (ierr.ne.NF_NOERR)
965
966  ! B.2. Prepare things to be able to write/append the variable
967  corners(1)=ntime
968 
969  edges(1)=1
970
971  ! B.3. Write the data
972!#ifdef NC_DOUBLE
973!  ierr=NF_PUT_VARA_DOUBLE(nid,varid,corners,edges,data0)
974!#else
975  ierr=NF_PUT_VARA_REAL(nid,varid,corners,edges,[data0])
976!#endif
977  if (ierr.ne.NF_NOERR) then
978    write(*,*) "writediagsoilpem: Error: Failed writing "//trim(name)//&
979               " to file "//trim(filename)//" at time",date
980  endif
981
982endif ! of if (dimpx.eq.3) elseif (dimpx.eq.2) ...
983endif ! of if (mod(zitau+1,isample).eq.0)
984
985! 4. Close the NetCDF file
986if (is_master) then
987  ierr=NF_CLOSE(nid)
988endif
989
990END SUBROUTINE writediagsoilpem
991
992END MODULE writediagpem_mod
Note: See TracBrowser for help on using the repository browser.