source: flexpart.git/src/writeheader_txt.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 6.3 KB
Line 
1subroutine writeheader_txt
2
3  !*****************************************************************************
4  !                                                                            *
5  !  This routine produces a file header containing basic information on the   *
6  !  settings of the FLEXPART run.                                             *
7  !  The header file is essential and must be read in by any postprocessing    *
8  !  program before reading in the output data.                                *
9  !                                                                            *
10  !     Author: A. Stohl                                                       *
11  !                                                                            *
12  !     7 August 2002                                                          *
13  !     modified IP 2013 for text output                                       *
14  !*****************************************************************************
15  !                                                                            *
16  ! Variables:                                                                 *
17  !                                                                            *
18  ! xlon                   longitude                                           *
19  ! xl                     model x coordinate                                  *
20  ! ylat                   latitude                                            *
21  ! yl                     model y coordinate                                  *
22  !                                                                            *
23  !*****************************************************************************
24
25  use point_mod
26  use outg_mod
27  use par_mod
28  use com_mod
29
30  implicit none
31
32!  integer :: jjjjmmdd,ihmmss,i,ix,jy,j
33  integer :: jjjjmmdd,ihmmss,i,j
34  real :: xp1,yp1,xp2,yp2
35
36
37  !************************
38  ! Open header output file
39  !************************
40
41  open(unitheader,file=path(2)(1:length(2))//'header_txt', &
42       form='formatted',err=998)
43  open(unitheader_txt,file=path(2)(1:length(2))//'header_txt_releases', &
44       form='formatted',err=998)
45
46
47  ! Write the header information
48  !*****************************
49 
50  write(unitheader,*) '# ibdate,ibtime, iedate, ietime, flexversion'
51  write(unitheader,*) ibdate, ibtime, iedate, ietime, trim(flexversion) !  'FLEXPART V9.0'
52  !if (ldirect.eq.1) then
53  !  write(unitheader,*) ibdate,ibtime,trim(flexversion) !  'FLEXPART V9.0'
54  !else
55  !  write(unitheader,*) iedate,ietime,trim(flexversion) ! 'FLEXPART V9.0'
56  !endif
57
58  ! Write info on output interval, averaging time, sampling time
59  !*************************************************************
60 
61  write(unitheader,*) '# interval, averaging time, sampling time'
62  write(unitheader,*) loutstep,loutaver,loutsample
63
64  ! Write information on output grid setup
65  !***************************************
66 
67  write(unitheader,*) '# information on grid setup    '
68  write(unitheader,*) '#outlon0,outlat0,numxgrid,numygrid,dxout,dyout'
69  write(unitheader,*) outlon0,outlat0,numxgrid,numygrid, &
70       dxout,dyout 
71  write(unitheader,*) '# numzgrid, outheight(.) '
72  write(unitheader,*) numzgrid,(outheight(i),i=1,numzgrid)
73
74  write(unitheader,*) '# jjjjmmdd,ihmmss'
75  call caldate(bdate,jjjjmmdd,ihmmss)
76  write(unitheader,*) jjjjmmdd,ihmmss
77
78  ! Write number of species, and name for each species (+extra name for depositions)
79  ! Indicate the vertical dimension of the fields (i.e., 1 for deposition fields, numzgrid for
80  ! concentration fields
81  !*****************************************************************************
82
83  write(unitheader,*) '# information on species'
84  write(unitheader,*) '# 3*nspec,maxpointspec_act'
85  write(unitheader,*) 3*nspec,maxpointspec_act
86  write(unitheader,*) '# for nspec:'
87  write(unitheader,*) '# 1, WD_ '
88  write(unitheader,*) '# 1, DD_ '
89  write(unitheader,*) '# numzgrid,species'
90  do i=1,nspec
91    write(unitheader,*) 1,'WD_'//species(i)(1:7)
92    write(unitheader,*) 1,'DD_'//species(i)(1:7)
93    write(unitheader,*) numzgrid,species(i)
94  end do
95
96  ! Write information on release points: total number, then for each point:
97  ! start, end, coordinates, # of particles, name, mass
98  !************************************************************************
99
100
101  write(unitheader_txt,*) '# information on release points'
102  write(unitheader_txt,*) '# numpoint'
103  write(unitheader_txt,*) numpoint
104  write(unitheader_txt,*) '# for numpoint:'
105  do i=1,numpoint
106    write(unitheader_txt,*) ireleasestart(i),ireleaseend(i),kindz(i)
107    xp1=xpoint1(i)*dx+xlon0
108    yp1=ypoint1(i)*dy+ylat0
109    xp2=xpoint2(i)*dx+xlon0
110    yp2=ypoint2(i)*dy+ylat0
111    write(unitheader_txt,*) xp1,yp1,xp2,yp2,zpoint1(i),zpoint2(i)
112    write(unitheader_txt,*) npart(i),1
113    if (numpoint.le.1000) then
114      write(unitheader_txt,*) compoint(i)
115    else
116      write(unitheader_txt,*) compoint(1001)
117    endif
118    do j=1,nspec
119      write(unitheader_txt,*) xmass(i,j)
120      write(unitheader_txt,*) xmass(i,j)
121      write(unitheader_txt,*) xmass(i,j)
122    end do
123  end do
124
125  ! Write information on model switches
126  !*****************************************
127
128  write(unitheader,*) '# information on model switches'
129  write(unitheader,*) '# method,lsubgrid,lconvection, ind_source,ind_receptor'
130  write(unitheader,*) method,lsubgrid,lconvection, &
131       ind_source,ind_receptor
132
133  ! Write age class information
134  !****************************
135 
136  write(unitheader,*) '# information on age class     '
137  write(unitheader,*) nageclass,(lage(i),i=1,nageclass)
138
139
140  !Do not write topography to text output file. Keep it on the binary one
141  !********************************
142
143  !do ix=0,numxgrid-1
144  !  write(unitheader,*) (oroout(ix,jy),jy=0,numygrid-1)
145  !end do
146
147
148 
149
150
151  close(unitheader)
152  close(unitheader_txt)
153
154
155!  open(unitheader,file=path(2)(1:length(2))//'header_nml', &
156!        form='formatted',err=998)
157!  write(unitheader,NML=COMMAND)
158!  close(unitheader)
159
160  return
161
162
163998   write(*,*) ' #### FLEXPART MODEL ERROR!   THE FILE         #### '
164  write(*,*) ' #### '//path(2)(1:length(2))//'header_txt'//' #### '
165  write(*,*) ' #### CANNOT BE OPENED. IF A FILE WITH THIS    #### '
166  write(*,*) ' #### NAME ALREADY EXISTS, DELETE IT AND START #### '
167  write(*,*) ' #### THE PROGRAM AGAIN.                       #### '
168  stop
169
170end subroutine writeheader_txt
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG