source: trunk/LMDZ.COMMON/libf/evolution/outputs.F90 @ 4050

Last change on this file since 4050 was 3991, checked in by jbclement, 3 months ago

PEM:
Apply documentation template everywhere: standardized headers format with short description, separators between functions/subroutines, normalized code sections, aligned dependencies/arguments/variables declaration.
JBC

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