source: flexpart.git/src/writeheader_nest.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: 5.4 KB
Line 
1subroutine writeheader_nest
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  !                                                                            *
14  !*****************************************************************************
15  !                                                                            *
16  !  Modified to remove TRIM around the output of flexversion so that          *
17  !  it will be a constant length (defined in com_mod.f90) in output header    *
18  !                                                                            *
19  !     Don Morton, Boreal Scientific Computing                                *
20  !     07 May 2017                                                            *
21  !                                                                            *
22  !*****************************************************************************
23  !                                                                            *
24  ! Variables:                                                                 *
25  !                                                                            *
26  ! xlon                   longitude                                           *
27  ! xl                     model x coordinate                                  *
28  ! ylat                   latitude                                            *
29  ! yl                     model y coordinate                                  *
30  !                                                                            *
31  !*****************************************************************************
32
33  use point_mod
34  use outg_mod
35  use par_mod
36  use com_mod
37
38  implicit none
39
40  integer :: jjjjmmdd,ihmmss,i,ix,jy,j
41  real :: xp1,yp1,xp2,yp2
42
43
44  !************************
45  ! Open header output file
46  !************************
47
48  open(unitheader,file=path(2)(1:length(2))//'header_nest', &
49       form='unformatted',err=998)
50
51
52  ! Write the header information
53  !*****************************
54
55  if (ldirect.eq.1) then
56    write(unitheader) ibdate,ibtime, flexversion
57  else
58    write(unitheader) iedate,ietime, flexversion
59  endif
60
61  ! Write info on output interval, averaging time, sampling time
62  !*************************************************************
63
64  write(unitheader) loutstep,loutaver,loutsample
65
66  ! Write information on output grid setup
67  !***************************************
68
69  write(unitheader) outlon0n,outlat0n,numxgridn,numygridn, &
70       dxoutn,dyoutn
71  write(unitheader) numzgrid,(outheight(i),i=1,numzgrid)
72
73  call caldate(bdate,jjjjmmdd,ihmmss)
74  write(unitheader) jjjjmmdd,ihmmss
75
76  ! Write number of species, and name for each species (+extra name for depositions)
77  ! Indicate the dimension of the fields (i.e., 1 for deposition fields, numzgrid for
78  ! concentration fields
79  !*****************************************************************************
80
81  write(unitheader) 3*nspec,maxpointspec_act
82  do i=1,nspec
83    write(unitheader) 1,'WD_'//species(i)(1:7)
84    write(unitheader) 1,'DD_'//species(i)(1:7)
85    write(unitheader) numzgrid,species(i)
86  end do
87
88  ! Write information on release points: total number, then for each point:
89  ! start, end, coordinates, # of particles, name, mass
90  !************************************************************************
91
92  write(unitheader) numpoint
93  do i=1,numpoint
94    write(unitheader) ireleasestart(i),ireleaseend(i),kindz(i)
95    xp1=xpoint1(i)*dx+xlon0
96    yp1=ypoint1(i)*dy+ylat0
97    xp2=xpoint2(i)*dx+xlon0
98    yp2=ypoint2(i)*dy+ylat0
99    write(unitheader) xp1,yp1,xp2,yp2,zpoint1(i),zpoint2(i)
100    write(unitheader) npart(i),1
101    if (numpoint.le.1000) then
102      write(unitheader) compoint(i)
103    else
104      write(unitheader) compoint(1001)
105   endif
106    do j=1,nspec
107      write(unitheader) xmass(i,j)
108      write(unitheader) xmass(i,j)
109      write(unitheader) xmass(i,j)
110    end do
111  end do
112
113  ! Write information on some model switches
114  !*****************************************
115
116  write(unitheader) method,lsubgrid,lconvection, &
117       ind_source,ind_receptor
118
119  ! Write age class information
120  !****************************
121
122  write(unitheader) nageclass,(lage(i),i=1,nageclass)
123
124
125  ! Write topography to output file
126  !********************************
127
128  do ix=0,numxgridn-1
129    write(unitheader) (orooutn(ix,jy),jy=0,numygridn-1)
130  end do
131  close(unitheader)
132
133  return
134
135
136998   write(*,*) ' #### FLEXPART MODEL ERROR!   THE FILE         #### '
137  write(*,*) ' #### '//path(2)(1:length(2))//'header'//' #### '
138  write(*,*) ' #### CANNOT BE OPENED. IF A FILE WITH THIS    #### '
139  write(*,*) ' #### NAME ALREADY EXISTS, DELETE IT AND START #### '
140  write(*,*) ' #### THE PROGRAM AGAIN.                       #### '
141  stop
142
143end subroutine writeheader_nest
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG