source: trunk/LMDZ.PLUTO/libf/phypluto/writediagfi.F

Last change on this file was 3184, checked in by afalco, 10 months ago

Pluto PCM:
Added LMDZ.PLUTO, a copy of the generic model,
cleaned from some unnecessary modules (water, ...)
AF

File size: 21.2 KB
Line 
1      subroutine writediagfi(ngrid,nom,titre,unite,dim,px)
2
3!  Ecriture de variables diagnostiques au choix dans la physique
4!  dans un fichier NetCDF nomme  'diagfi'. Ces variables peuvent etre
5!  3d (ex : temperature), 2d (ex : temperature de surface), ou
6!  0d (pour un scalaire qui ne depend que du temps : ex : la longitude
7!  solaire)
8!  (ou encore 1d, dans le cas de testphys1d, pour sortir une colonne)
9!  La periode d'ecriture est donnee par
10!  "ecritphy " regle dans le fichier de controle de run :  run.def
11!
12!    writediagfi peut etre appele de n'importe quelle subroutine
13!    de la physique, plusieurs fois. L'initialisation et la creation du
14!    fichier se fait au tout premier appel.
15!
16! WARNING : les variables dynamique (u,v,t,q,ps)
17!  sauvees par writediagfi avec une
18! date donnee sont legerement differentes que dans le fichier histoire car
19! on ne leur a pas encore ajoute de la dissipation et de la physique !!!
20! IL est  RECOMMANDE d'ajouter les tendance physique a ces variables
21! avant l'ecriture dans diagfi (cf. physiq.F)
22
23! Modifs: Aug.2010 Ehouarn: enforce outputs to be real*4
24!         Oct 2011 Francois: enable having a 'diagfi.def' file to select
25!                            at runtime, which variables to put in file
26!
27!  parametres (input) :
28!  ----------
29!      ngrid : nombres de point ou est calcule la physique
30!                (ngrid = 2+(jjm-1)*iim - 1/jjm)
31!                 (= nlon ou klon dans la physique terrestre)
32!     
33!      unit : unite logique du fichier de sortie (toujours la meme)
34!      nom  : nom de la variable a sortir (chaine de caracteres)
35!      titre: titre de la variable (chaine de caracteres)
36!      unite : unite de la variable (chaine de caracteres)
37!      px : variable a sortir (real 0, 1, 2, ou 3d)
38!      dim : dimension de px : 0, 1, 2, ou 3 dimensions
39!
40!=================================================================
41      use surfdat_h, only: phisfi
42      use geometry_mod, only: cell_area
43      use time_phylmdz_mod, only: ecritphy, day_step, iphysiq, day_ini
44      USE mod_phys_lmdz_para, only : is_parallel, is_mpi_root,
45     &                               is_master, gather
46      USE mod_grid_phy_lmdz, only : klon_glo, Grid1Dto2D_glo,
47     &                              nbp_lon, nbp_lat, nbp_lev,
48     &                              grid_type, unstructured
49      implicit none
50
51! Commons
52      include "netcdf.inc"
53
54! Arguments on input:
55      integer,intent(in) :: ngrid
56      character (len=*),intent(in) :: nom,titre,unite
57      integer,intent(in) :: dim
58      real,intent(in) :: px(ngrid,nbp_lev)
59
60! Local variables:
61
62      real*4 dx3(nbp_lon+1,nbp_lat,nbp_lev) ! to store a 3D data set
63      real*4 dx2(nbp_lon+1,nbp_lat)     ! to store a 2D (surface) data set
64      real*4 dx1(nbp_lev)           ! to store a 1D (column) data set
65      real*4 dx0
66      real*4 dx3_1d(1,nbp_lev) ! to store a profile with 1D model
67      real*4 dx2_1d ! to store a surface value with 1D model
68
69      real*4,save :: date
70!$OMP THREADPRIVATE(date)
71
72      REAL phis((nbp_lon+1),nbp_lat)
73      REAL area((nbp_lon+1),nbp_lat)
74
75      integer irythme
76      integer ierr,ierr2
77      integer i,j,l, ig0
78
79      integer,save :: zitau=0
80      character(len=40),save :: firstnom='1234567890'
81!$OMP THREADPRIVATE(zitau,firstnom)
82
83! Ajouts
84      integer, save :: ntime=0
85!$OMP THREADPRIVATE(ntime)
86      integer :: idim,varid
87      integer :: nid
88      character(len=*),parameter :: fichnom="diagfi.nc"
89      integer, dimension(4) :: id
90      integer, dimension(4) :: edges,corner
91
92! Added to use diagfi.def to select output variable
93      logical,save :: diagfi_def
94      logical :: getout
95      integer,save :: n_nom_def
96      integer :: n
97      integer,parameter :: n_nom_def_max=199
98      character(len=120),save :: nom_def(n_nom_def_max)
99      logical,save :: firstcall=.true.
100!$OMP THREADPRIVATE(firstcall)  !diagfi_def,n_nom_def,nom_def read in diagfi.def
101     
102#ifdef CPP_PARA
103! Added to work in parallel mode
104      real dx3_glop(klon_glo,nbp_lev)
105      real dx3_glo(nbp_lon,nbp_lat,nbp_lev) ! to store a global 3D data set
106      real dx2_glop(klon_glo)
107      real dx2_glo(nbp_lon,nbp_lat)     ! to store a global 2D (surface) data set
108      real px2(ngrid)
109!      real dx1_glo(nbp_lev)          ! to store a 1D (column) data set
110!      real dx0_glo
111      real phisfi_glo(klon_glo) ! surface geopotential on global physics grid
112      real areafi_glo(klon_glo) ! mesh area on global physics grid
113#else
114      real phisfi_glo(ngrid) ! surface geopotential on global physics grid
115      real areafi_glo(ngrid) ! mesh area on global physics grid
116#endif
117
118      if (grid_type==unstructured) then
119           return
120      endif
121
122!***************************************************************
123!Sortie des variables au rythme voulu
124
125      irythme = int(ecritphy) ! output rate set by ecritphy
126
127!***************************************************************
128
129! At very first call, check if there is a "diagfi.def" to use and read it
130! -----------------------------------------------------------------------
131      IF (firstcall) THEN
132         firstcall=.false.
133
134!$OMP MASTER
135  !      Open diagfi.def definition file if there is one:
136         open(99,file="diagfi.def",status='old',form='formatted',
137     s   iostat=ierr2)
138
139         if(ierr2.eq.0) then
140            diagfi_def=.true.
141            write(*,*) "******************"
142            write(*,*) "Reading diagfi.def"
143            write(*,*) "******************"
144            do n=1,n_nom_def_max
145              read(99,fmt='(a)',end=88) nom_def(n)
146              write(*,*) 'Output in diagfi: ', trim(nom_def(n))
147            end do
148 88         continue
149            if (n.ge.n_nom_def_max) then
150               write(*,*)"n_nom_def_max too small in writediagfi.F:",n
151               call abort_physic("writediagfi",
152     &             "n_nom_def_max too small",1)
153            end if
154            n_nom_def=n-1
155            close(99)
156         else
157            diagfi_def=.false.
158         endif
159!$OMP END MASTER
160!$OMP BARRIER
161      END IF ! of IF (firstcall)
162
163! Get out of write_diagfi if there is diagfi.def AND variable not listed
164!  ---------------------------------------------------------------------
165      if (diagfi_def) then
166          getout=.true.
167          do n=1,n_nom_def
168             if(trim(nom_def(n)).eq.nom) getout=.false.
169          end do
170          if (getout) return
171      end if
172
173! Initialisation of 'firstnom' and create/open the "diagfi.nc" NetCDF file
174! ------------------------------------------------------------------------
175! (at very first call to the subroutine, in accordance with diagfi.def)
176
177      if (firstnom.eq.'1234567890') then ! .true. for the very first valid
178      !   call to this subroutine; now set 'firstnom'
179         firstnom = nom
180         ! just to be sure, check that firstnom is large enough to hold nom
181         if (len_trim(firstnom).lt.len_trim(nom)) then
182           write(*,*) "writediagfi: Error !!!"
183           write(*,*) "   firstnom string not long enough!!"
184           write(*,*) "   increase its size to at least ",len_trim(nom)
185           call abort_physic("writediagfi","firstnom too short",1)
186         endif
187         
188         zitau = -1 ! initialize zitau
189
190#ifdef CPP_PARA
191          ! Gather phisfi() geopotential on physics grid
192          call Gather(phisfi,phisfi_glo)
193          ! Gather cell_area() mesh area on physics grid
194          call Gather(cell_area,areafi_glo)
195#else
196         phisfi_glo(:)=phisfi(:)
197         areafi_glo(:)=cell_area(:)
198#endif
199
200         !! parallel: we cannot use the usual writediagfi method
201!!         call iophys_ini
202         if (is_master) then
203         ! only the master is required to do this
204
205         ! Create the NetCDF file
206         ierr = NF_CREATE(fichnom, NF_CLOBBER, nid)
207         ! Define the 'Time' dimension
208         ierr = nf_def_dim(nid,"Time",NF_UNLIMITED,idim)
209         ! Define the 'Time' variable
210!#ifdef NC_DOUBLE
211!         ierr = NF_DEF_VAR (nid, "Time", NF_DOUBLE, 1, idim,varid)
212!#else
213         ierr = NF_DEF_VAR (nid, "Time", NF_FLOAT, 1, idim,varid)
214!#endif
215         ! Add a long_name attribute
216         ierr = NF_PUT_ATT_TEXT (nid, varid, "long_name",
217     .          4,"Time")
218         ! Add a units attribute
219         ierr = NF_PUT_ATT_TEXT(nid, varid,'units',29,
220     .          "days since 0000-00-0 00:00:00")
221         ! Switch out of NetCDF Define mode
222         ierr = NF_ENDDEF(nid)
223
224         ! Build phis() and area()
225         IF (klon_glo>1) THEN
226          do i=1,nbp_lon+1 ! poles
227           phis(i,1)=phisfi_glo(1)
228           phis(i,nbp_lat)=phisfi_glo(klon_glo)
229           ! for area, divide at the poles by nbp_lon
230           area(i,1)=areafi_glo(1)/nbp_lon
231           area(i,nbp_lat)=areafi_glo(klon_glo)/nbp_lon
232          enddo
233          do j=2,nbp_lat-1
234           ig0= 1+(j-2)*nbp_lon
235           do i=1,nbp_lon
236              phis(i,j)=phisfi_glo(ig0+i)
237              area(i,j)=areafi_glo(ig0+i)
238           enddo
239           ! handle redundant point in longitude
240           phis(nbp_lon+1,j)=phis(1,j)
241           area(nbp_lon+1,j)=area(1,j)
242          enddo
243         ENDIF
244         
245         ! write "header" of file (longitudes, latitudes, geopotential, ...)
246         IF (klon_glo>1) THEN ! general 3D case
247           call iniwrite(nid,day_ini,phis,area,nbp_lon+1,nbp_lat)
248         ELSE ! 1D model
249           call iniwrite(nid,day_ini,phisfi_glo(1),areafi_glo(1),1,1)
250         ENDIF
251
252         endif ! of if (is_master)
253
254      else
255
256         if (is_master) then
257           ! only the master is required to do this
258
259           ! Open the NetCDF file
260           ierr = NF_OPEN(fichnom,NF_WRITE,nid)
261         endif ! of if (is_master)
262
263      endif ! if (firstnom.eq.'1234567890')
264
265! Increment time index 'zitau' if it is the "fist call" (at given time level)
266! to writediagfi
267!------------------------------------------------------------------------
268      if (nom.eq.firstnom) then
269          zitau = zitau + iphysiq
270      end if
271
272!--------------------------------------------------------
273! Write the variables to output file if it's time to do so
274!--------------------------------------------------------
275
276      if ( MOD(zitau+1,irythme) .eq.0.) then
277
278! Compute/write/extend 'Time' coordinate (date given in days)
279! (done every "first call" (at given time level) to writediagfi)
280! Note: date is incremented as 1 step ahead of physics time
281!--------------------------------------------------------
282
283        if (is_master) then
284           ! only the master is required to do this
285        if (nom.eq.firstnom) then
286        ! We have identified a "first call" (at given date)
287           ntime=ntime+1 ! increment # of stored time steps
288           ! compute corresponding date (in days and fractions thereof)
289           date=(zitau +1.)/day_step
290           ! Get NetCDF ID of 'Time' variable
291           ierr= NF_INQ_VARID(nid,"Time",varid)
292           ! Write (append) the new date to the 'Time' array
293!#ifdef NC_DOUBLE
294!           ierr= NF_PUT_VARA_DOUBLE(nid,varid,ntime,1,date)
295!#else
296           ierr= NF_PUT_VARA_REAL(nid,varid,ntime,1,date)
297!#endif
298           if (ierr.ne.NF_NOERR) then
299              write(*,*) "***** PUT_VAR matter in writediagfi_nc"
300              write(*,*) "***** with time"
301              write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr) 
302c             call abort
303           endif
304
305           write(6,*)'WRITEDIAGFI: date= ', date
306        end if ! of if (nom.eq.firstnom)
307
308        endif ! of if (is_master)
309
310!Case of a 3D variable
311!---------------------
312        if (dim.eq.3) then
313
314#ifdef CPP_PARA
315          ! Gather field on a "global" (without redundant longitude) array
316          call Gather(px,dx3_glop)
317!$OMP MASTER
318          if (is_mpi_root) then
319            call Grid1Dto2D_glo(dx3_glop,dx3_glo)
320            ! copy dx3_glo() to dx3(:) and add redundant longitude
321            dx3(1:nbp_lon,:,:)=dx3_glo(1:nbp_lon,:,:)
322            dx3(nbp_lon+1,:,:)=dx3(1,:,:)
323          endif
324!$OMP END MASTER
325!$OMP BARRIER
326#else
327!         Passage variable physique -->  variable dynamique
328!         recast (copy) variable from physics grid to dynamics grid
329          IF (klon_glo>1) THEN ! General case
330           DO l=1,nbp_lev
331             DO i=1,nbp_lon+1
332                dx3(i,1,l)=px(1,l)
333                dx3(i,nbp_lat,l)=px(ngrid,l)
334             ENDDO
335             DO j=2,nbp_lat-1
336                ig0= 1+(j-2)*nbp_lon
337                DO i=1,nbp_lon
338                   dx3(i,j,l)=px(ig0+i,l)
339                ENDDO
340                dx3(nbp_lon+1,j,l)=dx3(1,j,l)
341             ENDDO
342           ENDDO
343          ELSE ! 1D model case
344           dx3_1d(1,1:nbp_lev)=px(1,1:nbp_lev)
345          ENDIF
346#endif
347!         Ecriture du champs
348
349          if (is_master) then
350           ! only the master writes to output
351! name of the variable
352           ierr= NF_INQ_VARID(nid,nom,varid)
353           if (ierr /= NF_NOERR) then
354! corresponding dimensions
355              ierr= NF_INQ_DIMID(nid,"longitude",id(1))
356              ierr= NF_INQ_DIMID(nid,"latitude",id(2))
357              ierr= NF_INQ_DIMID(nid,"altitude",id(3))
358              ierr= NF_INQ_DIMID(nid,"Time",id(4))
359
360! Create the variable if it doesn't exist yet
361
362              write (*,*) "=========================="
363              write (*,*) "DIAGFI: creating variable ",trim(nom)
364              call def_var(nid,nom,titre,unite,4,id,varid,ierr)
365
366           else
367             if (ntime==0) then
368              write(*,*) "DIAGFI Error: failed creating variable ",
369     &                   trim(nom)
370              write(*,*) "it seems it already exists!"
371              call abort_physic("writediagfi",
372     &             trim(nom)//" already exists",1)
373             endif
374           endif
375
376           corner(1)=1
377           corner(2)=1
378           corner(3)=1
379           corner(4)=ntime
380
381           IF (klon_glo==1) THEN
382             edges(1)=1
383           ELSE
384             edges(1)=nbp_lon+1
385           ENDIF
386           edges(2)=nbp_lat
387           edges(3)=nbp_lev
388           edges(4)=1
389!#ifdef NC_DOUBLE
390!           ierr= NF_PUT_VARA_DOUBLE(nid,varid,corner,edges,dx3)
391!#else
392!           write(*,*)"test:  nid=",nid," varid=",varid
393!           write(*,*)"       corner()=",corner
394!           write(*,*)"       edges()=",edges
395!           write(*,*)"       dx3()=",dx3
396           IF (klon_glo>1) THEN ! General case
397             ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx3)
398           ELSE
399             ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx3_1d)
400           ENDIF
401!#endif
402
403           if (ierr.ne.NF_NOERR) then
404              write(*,*) "***** PUT_VAR problem in writediagfi"
405              write(*,*) "***** with dx3: ",trim(nom)
406              write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr)
407              call abort_physic("writediagfi",
408     &             "failed writing "//trim(nom),1)
409           endif
410
411          endif !of if (is_master)
412
413!Case of a 2D variable
414!---------------------
415
416        else if (dim.eq.2) then
417
418#ifdef CPP_PARA
419          ! Gather field on a "global" (without redundant longitude) array
420          px2(:)=px(:,1)
421          call Gather(px2,dx2_glop)
422!$OMP MASTER
423          if (is_mpi_root) then
424            call Grid1Dto2D_glo(dx2_glop,dx2_glo)
425            ! copy dx2_glo() to dx2(:) and add redundant longitude
426            dx2(1:nbp_lon,:)=dx2_glo(1:nbp_lon,:)
427            dx2(nbp_lon+1,:)=dx2(1,:)
428          endif
429!$OMP END MASTER
430!$OMP BARRIER
431#else
432
433!         Passage variable physique -->  physique dynamique
434!         recast (copy) variable from physics grid to dynamics grid
435          IF (klon_glo>1) THEN ! General case
436             DO i=1,nbp_lon+1
437                dx2(i,1)=px(1,1)
438                dx2(i,nbp_lat)=px(ngrid,1)
439             ENDDO
440             DO j=2,nbp_lat-1
441                ig0= 1+(j-2)*nbp_lon
442                DO i=1,nbp_lon
443                   dx2(i,j)=px(ig0+i,1)
444                ENDDO
445                dx2(nbp_lon+1,j)=dx2(1,j)
446             ENDDO
447          ELSE ! 1D model case
448            dx2_1d=px(1,1)
449          ENDIF
450#endif
451
452          if (is_master) then
453           ! only the master writes to output
454!         write (*,*) 'In  writediagfi, on sauve:  ' , nom
455!         write (*,*) 'In  writediagfi. Estimated date = ' ,date
456           ierr= NF_INQ_VARID(nid,nom,varid)
457           if (ierr /= NF_NOERR) then
458! corresponding dimensions
459              ierr= NF_INQ_DIMID(nid,"longitude",id(1))
460              ierr= NF_INQ_DIMID(nid,"latitude",id(2))
461              ierr= NF_INQ_DIMID(nid,"Time",id(3))
462
463! Create the variable if it doesn't exist yet
464
465              write (*,*) "=========================="
466              write (*,*) "DIAGFI: creating variable ",trim(nom)
467
468              call def_var(nid,nom,titre,unite,3,id,varid,ierr)
469
470           else
471             if (ntime==0) then
472              write(*,*) "DIAGFI Error: failed creating variable ",
473     &                   trim(nom)
474              write(*,*) "it seems it already exists!"
475              call abort_physic("writediagfi",
476     &             trim(nom)//" already exists",1)
477             endif
478           endif
479
480           corner(1)=1
481           corner(2)=1
482           corner(3)=ntime
483           IF (klon_glo==1) THEN
484             edges(1)=1
485           ELSE
486             edges(1)=nbp_lon+1
487           ENDIF
488           edges(2)=nbp_lat
489           edges(3)=1
490
491
492!#ifdef NC_DOUBLE
493!           ierr = NF_PUT_VARA_DOUBLE (nid,varid,corner,edges,dx2)
494!#else         
495           IF (klon_glo>1) THEN ! General case
496             ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx2)
497           ELSE
498             ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx2_1d)
499           ENDIF
500!#endif     
501
502           if (ierr.ne.NF_NOERR) then
503              write(*,*) "***** PUT_VAR matter in writediagfi"
504              write(*,*) "***** with dx2: ",trim(nom)
505              write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr)
506              call abort_physic("writediagfi",
507     &             "failed writing "//trim(nom),1)
508           endif
509
510          endif !of if (is_master)
511
512!Case of a 1D variable (ie: a column)
513!---------------------------------------------------
514
515       else if (dim.eq.1) then
516         if (is_parallel) then
517           write(*,*) "writediagfi error: dim=1 not implemented ",
518     &                 "in parallel mode. Problem for ",trim(nom)
519              call abort_physic("writediagfi",
520     &             "failed writing "//trim(nom),1)
521         endif
522!         Passage variable physique -->  physique dynamique
523!         recast (copy) variable from physics grid to dynamics grid
524          do l=1,nbp_lev
525            dx1(l)=px(1,l)
526          enddo
527         
528          ierr= NF_INQ_VARID(nid,nom,varid)
529           if (ierr /= NF_NOERR) then
530! corresponding dimensions
531              ierr= NF_INQ_DIMID(nid,"altitude",id(1))
532              ierr= NF_INQ_DIMID(nid,"Time",id(2))
533
534! Create the variable if it doesn't exist yet
535
536              write (*,*) "=========================="
537              write (*,*) "DIAGFI: creating variable ",trim(nom)
538
539              call def_var(nid,nom,titre,unite,2,id,varid,ierr)
540             
541           else
542             if (ntime==0) then
543              write(*,*) "DIAGFI Error: failed creating variable ",
544     &                   trim(nom)
545              write(*,*) "it seems it already exists!"
546              call abort_physic("writediagfi",
547     &             trim(nom)//" already exists",1)
548             endif
549           endif
550           
551           corner(1)=1
552           corner(2)=ntime
553           
554           edges(1)=nbp_lev
555           edges(2)=1
556!#ifdef NC_DOUBLE
557!           ierr= NF_PUT_VARA_DOUBLE(nid,varid,corner,edges,dx1)
558!#else
559           ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx1)
560!#endif
561
562           if (ierr.ne.NF_NOERR) then
563              write(*,*) "***** PUT_VAR problem in writediagfi"
564              write(*,*) "***** with dx1: ",trim(nom)
565              write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr)
566              call abort_physic("writediagfi",
567     &             "failed writing "//trim(nom),1)
568           endif
569
570!Case of a 0D variable (ie: a time-dependent scalar)
571!---------------------------------------------------
572
573        else if (dim.eq.0) then
574
575           dx0 = px (1,1)
576
577          if (is_master) then
578           ! only the master writes to output
579           ierr= NF_INQ_VARID(nid,nom,varid)
580           if (ierr /= NF_NOERR) then
581! corresponding dimensions
582              ierr= NF_INQ_DIMID(nid,"Time",id(1))
583
584! Create the variable if it doesn't exist yet
585
586              write (*,*) "=========================="
587              write (*,*) "DIAGFI: creating variable ",trim(nom)
588
589              call def_var(nid,nom,titre,unite,1,id,varid,ierr)
590
591           else
592             if (ntime==0) then
593              write(*,*) "DIAGFI Error: failed creating variable ",
594     &                   trim(nom)
595              write(*,*) "it seems it already exists!"
596              call abort_physic("writediagfi",
597     &             trim(nom)//" already exists",1)
598             endif
599           endif
600
601           corner(1)=ntime
602           edges(1)=1
603
604!#ifdef NC_DOUBLE
605!           ierr = NF_PUT_VARA_DOUBLE (nid,varid,corner,edges,dx0) 
606!#else
607           ierr= NF_PUT_VARA_REAL(nid,varid,corner,edges,dx0)
608!#endif
609           if (ierr.ne.NF_NOERR) then
610              write(*,*) "***** PUT_VAR matter in writediagfi"
611              write(*,*) "***** with dx0: ",trim(nom)
612              write(*,*) 'ierr=', ierr,": ",NF_STRERROR(ierr)
613              call abort_physic("writediagfi",
614     &             "failed writing "//trim(nom),1)
615           endif
616
617          endif !of if (is_master)
618
619        endif ! of if (dim.eq.3) elseif(dim.eq.2)...
620
621      endif ! of if ( MOD(zitau+1,irythme) .eq.0.)
622
623      if (is_master) then
624        ierr= NF_CLOSE(nid)
625      endif
626
627      end
Note: See TracBrowser for help on using the repository browser.