source: branches/petra/src/writeheader_nest.f90 @ 36

Last change on this file since 36 was 36, checked in by pesei, 9 years ago

Implement switch for incremental deposition, see ticket:113 and many small changes, see changelog.txt

File size: 6.4 KB
Line 
1!**********************************************************************
2! Copyright 1998-2015                                                 *
3! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa,             *
4! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann   *
5!                                                                     *
6! This file is part of FLEXPART.                                      *
7!                                                                     *
8! FLEXPART is free software: you can redistribute it and/or modify    *
9! it under the terms of the GNU General Public License as published by*
10! the Free Software Foundation, either version 3 of the License, or   *
11! (at your option) any later version.                                 *
12!                                                                     *
13! FLEXPART is distributed in the hope that it will be useful,         *
14! but WITHOUT ANY WARRANTY; without even the implied warranty of      *
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
16! GNU General Public License for more details.                        *
17!                                                                     *
18! You should have received a copy of the GNU General Public License   *
19! along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
20!**********************************************************************
21
22subroutine writeheader_nest
23
24  !*****************************************************************************
25  !                                                                            *
26  !  This routine produces a file header containing basic information on the   *
27  !  settings of the FLEXPART run.                                             *
28  !  The header file is essential and must be read in by any postprocessing    *
29  !  program before reading in the output data.                                *
30  !                                                                            *
31  !     Author: A. Stohl                                                       *
32  !                                                                            *
33  !     7 August 2002                                                          *
34  !                                                                            *
35  !     2/2015, AP+PS: version string length written 
36  !                                                                            *
37  !*****************************************************************************
38  !                                                                            *
39  ! Variables:                                                                 *
40  !                                                                            *
41  ! xlon                   longitude                                           *
42  ! xl                     model x coordinate                                  *
43  ! ylat                   latitude                                            *
44  ! yl                     model y coordinate                                  *
45  !                                                                            *
46  !*****************************************************************************
47
48  use point_mod
49  use outg_mod
50  use par_mod
51  use com_mod
52
53  implicit none
54
55  integer :: jjjjmmdd,ihmmss,i,ix,jy,j
56  real :: xp1,yp1,xp2,yp2
57
58
59  !************************
60  ! Open header output file
61  !************************
62
63  open(unitheader,file=path(2)(1:length(2))//'header_nest', &
64       form='unformatted',err=998)
65
66
67  ! Write the header information
68  !*****************************
69
70  if (ldirect.eq.1) then
71    write(unitheader) ibdate,ibtime, len_trim(flexversion), trim(flexversion)
72  else
73    write(unitheader) iedate,ietime, len_trim(flexversion), trim(flexversion)
74  endif
75
76  ! Write info on output interval, averaging time, sampling time
77  !*************************************************************
78
79  write(unitheader) loutstep,loutaver,loutsample
80
81  ! Write information on output grid setup
82  !***************************************
83
84  write(unitheader) outlon0n,outlat0n,numxgridn,numygridn, &
85       dxoutn,dyoutn
86  write(unitheader) numzgrid,(outheight(i),i=1,numzgrid)
87
88  call caldate(bdate,jjjjmmdd,ihmmss)
89  write(unitheader) jjjjmmdd,ihmmss
90
91  ! Write number of species, and name for each species (+extra name for depositions)
92  ! Indicate the dimension of the fields (i.e., 1 for deposition fields, numzgrid for
93  ! concentration fields
94  !*****************************************************************************
95
96  write(unitheader) 3*nspec,maxpointspec_act
97  do i=1,nspec
98    write(unitheader) 1,'WD_'//species(i)(1:7)
99    write(unitheader) 1,'DD_'//species(i)(1:7)
100    write(unitheader) numzgrid,species(i)
101  end do
102
103  ! Write information on release points: total number, then for each point:
104  ! start, end, coordinates, # of particles, name, mass
105  !************************************************************************
106
107  write(unitheader) numpoint
108  do i=1,numpoint
109    write(unitheader) ireleasestart(i),ireleaseend(i),kindz(i)
110    xp1=xpoint1(i)*dx+xlon0
111    yp1=ypoint1(i)*dy+ylat0
112    xp2=xpoint2(i)*dx+xlon0
113    yp2=ypoint2(i)*dy+ylat0
114    write(unitheader) xp1,yp1,xp2,yp2,zpoint1(i),zpoint2(i)
115    write(unitheader) npart(i),1
116    if (numpoint.le.1000) then
117      write(unitheader) compoint(i)
118    else
119      write(unitheader) compoint(1001)
120   endif
121    do j=1,nspec
122      write(unitheader) xmass(i,j)
123      write(unitheader) xmass(i,j)
124      write(unitheader) xmass(i,j)
125    end do
126  end do
127
128  ! Write information on some model switches
129  !*****************************************
130
131  write(unitheader) method,lsubgrid,lconvection, &
132       ind_source,ind_receptor
133
134  ! Write age class information
135  !****************************
136
137  write(unitheader) nageclass,(lage(i),i=1,nageclass)
138
139
140  ! Write topography to output file
141  !********************************
142
143  do ix=0,numxgridn-1
144    write(unitheader) (orooutn(ix,jy),jy=0,numygridn-1)
145  end do
146  close(unitheader)
147
148  return
149
150
151998   write(*,*) ' #### FLEXPART MODEL ERROR!   THE FILE         #### '
152  write(*,*) ' #### '//path(2)(1:length(2))//'header'//' #### '
153  write(*,*) ' #### CANNOT BE OPENED. IF A FILE WITH THIS    #### '
154  write(*,*) ' #### NAME ALREADY EXISTS, DELETE IT AND START #### '
155  write(*,*) ' #### THE PROGRAM AGAIN.                       #### '
156  stop
157
158end subroutine writeheader_nest
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG