updated top-level README and version_decl for V4.4.2 (#1795)
[WRF.git] / external / io_grib2 / g2lib / gribend.F
blobf2b6c3a5562fdcae76aff3811f8ec0fbf7123af5
1       subroutine gribend(cgrib,lcgrib,lengrib,ierr)
2 !$$$  SUBPROGRAM DOCUMENTATION BLOCK
3 !                .      .    .                                       .
4 ! SUBPROGRAM:    gribend 
5 !   PRGMMR: Gilbert         ORG: W/NP11    DATE: 2000-05-02
7 ! ABSTRACT: This subroutine finalizes a GRIB message after all grids
8 !   and fields have been added.  It adds the End Section ( "7777" )
9 !   to the end of the GRIB message and calculates the length and stores
10 !   it in the appropriate place in Section 0.
11 !   This routine is used with routines "gribcreate", "addlocal", "addgrid",
12 !   and "addfield" to create a complete GRIB2 message.  Subroutine
13 !   gribcreate must be called first to initialize a new GRIB2 message.
15 ! PROGRAM HISTORY LOG:
16 ! 2000-05-02  Gilbert
18 ! USAGE:    CALL gribend(cgrib,lcgrib,lengrib,ierr)
19 !   INPUT ARGUMENT LIST:
20 !     cgrib    - Character array to contain the GRIB2 message
21 !     lcgrib   - Maximum length (bytes) of array cgrib.
23 !   OUTPUT ARGUMENT LIST:      
24 !     cgrib    - Character array to contain the GRIB2 message
25 !     lengrib  - Length of the final GRIB2 message in octets (bytes)
26 !     ierr     - Error return code.
27 !                0 = no error
28 !                1 = GRIB message was not initialized.  Need to call
29 !                    routine gribcreate first.
30 !                2 = GRIB message already complete.  
31 !                3 = Sum of Section byte counts doesn't add to total byte count.
32 !                4 = Previous Section was not 7.
34 ! REMARKS: This routine is intended for use with routines "gribcreate", 
35 !          "addlocal", "addgrid", and "addfield" to create a complete 
36 !          GRIB2 message.
38 ! ATTRIBUTES:
39 !   LANGUAGE: Fortran 90
40 !   MACHINE:  IBM SP
42 !$$$
44       character(len=1),intent(inout) :: cgrib(lcgrib)
45       integer,intent(in) :: lcgrib
46       integer,intent(out) :: lengrib,ierr
47       
48       character(len=4),parameter :: grib='GRIB',c7777='7777'
49       character(len=4):: ctemp
50       integer iofst,ibeg,lencurr,len
52       ierr=0
54 !  Check to see if beginning of GRIB message exists
56       ctemp=cgrib(1)//cgrib(2)//cgrib(3)//cgrib(4)
57       if ( ctemp.ne.grib ) then
58         print *,'gribend: GRIB not found in given message.'
59         ierr=1
60         return
61       endif
63 !  Get current length of GRIB message
64 !  
65       call g2lib_gbyte(cgrib,lencurr,96,32)
67 !  Check to see if GRIB message is already complete
69 !      ctemp=cgrib(lencurr-3)//cgrib(lencurr-2)//cgrib(lencurr-1)
70 !     &      //cgrib(lencurr)
71 !      if ( ctemp.eq.c7777 ) then
72 !        print *,'gribend: GRIB message already complete.'
73 !        ierr=2
74 !        return
75 !      endif
77 !  Loop through all current sections of the GRIB message to
78 !  find the last section number.
80       len=16    ! Length of Section 0
81       do 
82       !    Get number and length of next section
83         iofst=len*8
84         call g2lib_gbyte(cgrib,ilen,iofst,32)
85         iofst=iofst+32
86         call g2lib_gbyte(cgrib,isecnum,iofst,8)
87         len=len+ilen
88       !    Exit loop if last section reached
89         if ( len.eq.lencurr ) exit
90       !    If byte count for each section doesn't match current
91       !    total length, then there is a problem.
92         if ( len.gt.lencurr ) then
93           print *,'gribend: Section byte counts don''t add to total.'
94           print *,'gribend: Sum of section byte counts = ',len
95           print *,'gribend: Total byte count in Section 0 = ',lencurr
96           ierr=3
97           return
98         endif
99       enddo
101 !  Can only add End Section (Section 8) after Section 7.
103       if ( isecnum.ne.7 ) then
104         print *,'gribend: Section 8 can only be added after Section 7.'
105         print *,'gribend: Section ',isecnum,' was the last found in',
106      &          ' given GRIB message.'
107         ierr=4
108         return
109       endif
111 !  Add Section 8  - End Section
113       cgrib(lencurr+1:lencurr+4)=c7777
116 !  Update current byte total of message in Section 0
118       lengrib=lencurr+4
119       call g2lib_sbyte(cgrib,lengrib,96,32)
121       return
122       end