varNameList_mod sourceΒΆ

  1module varNameList_mod
  2  ! MODULE varNameList_mod (prefix='vnl' category='7. Low-level data objects')
  3  !
  4  !:Purpose:  Contains a list of all possible variable names that can be used
  5  !           as analysis variables along with additional information for each
  6  !           and procedures for accessing this information
  7  !
  8  use bufr_mod
  9  use midasMpi_mod
 10  use utilities_mod
 11  use MathPhysConstants_mod
 12
 13  implicit none
 14  save
 15  private
 16
 17  ! public variables (parameters)
 18  public :: vnl_numvarmax3D, vnl_numvarmax2D, vnl_numvarmaxOther, vnl_numvarmax
 19  public :: vnl_varNameList3D, vnl_varNameList2D, vnl_varNameListOther, vnl_varNameList
 20  public :: vnl_numvarmaxCloud, vnl_varNameListCloud
 21
 22  ! public procedures
 23  public :: vnl_varListIndex3d, vnl_varListIndex2d, vnl_varListIndexOther
 24  public :: vnl_varListIndex, vnl_varnameFromVarnum, vnl_varnameIsValid
 25  public :: vnl_varLevelFromVarname, vnl_varLevelFromVarnum
 26  public :: vnl_varKindFromVarname, vnl_varnumFromVarname
 27  public :: vnl_varNamesFromExistList, vnl_varMassFromVarNum, vnl_varMassFromVarName
 28  public :: vnl_isPhysicsVar, vnl_isCloudVar
 29  public :: vnl_addToVarNames
 30
 31  ! These private parameters permit side-stepping a conflict with the Sphinx documenter,
 32  ! and an infinite loop
 33  integer, parameter          :: VNLnumvarmax3D    = 52
 34  integer, parameter          :: VNLnumvarmax2D    = 37
 35  integer, parameter          :: VNLnumvarmaxOther =  6
 36  integer, parameter          :: VNLnumvarmaxCloud =  5
 37
 38  integer, parameter          :: vnl_numvarmax3D    = VNLnumvarmax3D
 39  integer, parameter          :: vnl_numvarmax2D    = VNLnumvarmax2D
 40  integer, parameter          :: vnl_numvarmaxOther = VNLnumvarmaxOther
 41  integer, parameter          :: vnl_numvarmaxCloud = VNLnumvarmaxCloud
 42
 43  character(len=4), parameter :: vnl_varNameList3D(vnl_numvarmax3D) = (/                         &
 44                                 'UU  ','VV  ','Z_T ','Z_M ','P_T ','P_M ',                      &
 45                                 'TT  ','HU  ','LQ  ','ES  ','VT  ',                             &
 46                                 'PP  ','CC  ','UC  ','UT  ','TB  ','DW  ','QR  ','DD  ',        &
 47                                 'TO3 ','O3L ','TCH4','TCO2','TCO ','TNO2','TN2O','THCH',        &
 48                                 'TSO2','TNH3','AF  ','AC  ','TNO ','ALFA','VIS ','LVIS',        &
 49                                 'HR  ','TD  ','ALFT','UV  ','LWCR','IWCR','QC  ','CH4L',        &
 50                                 'N2OL','UUW ','VVW ','TM  ','SALW','ALFO','RF  ','SF  ',        &
 51                                 'CLDR' /)
 52
 53  character(len=4), parameter :: varLevelList3D(vnl_numvarmax3D)     = (/                        &
 54                                 'MM',  'MM',  'TH',  'MM',  'TH',  'MM',                        &
 55                                 'TH',  'TH',  'TH',  'TH',  'TH',                               &
 56                                 'MM',  'MM',  'MM',  'TH',  'TH',  'TH',  'MM',  'MM',          &
 57                                 'TH',  'TH',  'TH',  'TH',  'TH',  'TH',  'TH',  'TH',          &
 58                                 'TH',  'TH',  'TH',  'TH',  'TH',  'MM',  'TH',  'TH',          &
 59                                 'TH',  'TH',  'TH',  'MM',  'TH',  'TH',  'TH',  'TH',          &
 60                                 'TH',  'DP',  'DP',  'DP',  'DP',  'DP',  'TH',  'TH',          &
 61                                 'TH' /)
 62
 63  character(len=2), parameter :: varKindList3D(vnl_numvarmax3D)     = (/                         &
 64                                 'MT',  'MT',  'MT',  'MT',  'MT',  'MT',                        &
 65                                 'MT',  'MT',  'MT',  'MT',  'MT',                               &
 66                                 'MT',  'MT',  'MT',  'MT',  'MT',  'MT',  'MT',  'MT',          &
 67                                 'CH',  'CH',  'CH',  'CH',  'CH',  'CH',  'CH',  'CH',          &
 68                                 'CH',  'CH',  'CH',  'CH',  'CH',  'MT',  'MT',  'MT',          &
 69                                 'MT',  'MT',  'MT',  'MT',  'MT',  'MT',  'MT',  'CH',          &
 70                                 'CH',  'OC',  'OC',  'OC',  'OC',  'OC',  'MT',  'MT',          &
 71                                 'MT' /)
 72
 73  character(len=4), parameter :: vnl_varNameList2D(vnl_numvarmax2D) = (/ &
 74                                 'P0  ','TG  ','UP  ','PB  ','ECO ','ENO2','EHCH','ESO2','ENH3', &
 75                                 'GL  ','WGE ','BIN ','MG  ','SSH ','QI1 ','QO1 ','STOR','ALFS', &
 76                                 'PN  ','PR  ','LPR ','I2  ','I3  ','I4  ','I5  ','I6  ','I8  ', &
 77                                 'DN  ','FB  ','FI  ','MSKC','LZS ','WT  ','LG  ','VF  ','DSLO', &
 78                                 'P0LS'/)
 79
 80  character(len=4), parameter :: varLevelList2D(vnl_numvarmax2D) = (/    &
 81                                 'SF',  'SF',  'SF',  'SF',  'SF',  'SF',  'SF',  'SF',  'SF',  &
 82                                 'SS',  'SF',  'SF',  'SF',  'SS',  'SF',  'SF',  'SF',  'SF',  &
 83                                 'SF',  'SF',  'SF',  'SF',  'SF',  'SF',  'SF',  'SF',  'SF',  &
 84                                 'SF',  'SF',  'SF',  'SF',  'SF',  'SF',  'SS',  'SS',  'SS',  &
 85                                 'SF'/)
 86
 87  character(len=2), parameter :: varKindList2D(vnl_numvarmax2D) = (/     &
 88                                 'MT',  'MT',  'MT',  'MT',  'CH',  'CH',  'CH',  'CH',  'CH', &
 89                                 'OC',  'MT',  'MT',  'MT',  'OC',  'HY',  'HY',  'HY',  'HY', &
 90                                 'MT',  'MT',  'MT',  'MT',  'MT',  'MT',  'MT',  'MT',  'MT', &
 91                                 'MT',  'MT',  'MT',  'MT',  'HY',  'MT',  'OC',  'OC',  'OC', &
 92                                 'MT'/)
 93
 94  character(len=4), parameter :: vnl_varNameListOther(vnl_numvarmaxOther) = (/ &
 95                                 'I0  ','I1  ','I7  ','I9  ','SD  ','AL  '/)
 96
 97  character(len=4), parameter :: varLevelListOther(vnl_numvarmaxOther) = (/    &
 98                                 'OT',  'OT',  'OT',  'OT',  'OT',  'OT'  /)
 99
100  character(len=2), parameter :: varKindListOther(vnl_numvarmaxOther) = (/     &
101                                 'LD',  'LD',  'LD',  'LD',  'LD',  'LD'  /) ! LD = Land
102
103  character(len=4), parameter :: vnl_varNameListCloud(vnl_numvarmaxCloud) = (/ &
104                                 'LWCR', 'IWCR', 'RF  ', 'SF  ', 'CLDR' /)                                 
105
106  integer, parameter          :: vnl_numvarmax = VNLnumvarmax3D + VNLnumvarmax2D + VNLnumvarmaxOther
107
108  character(len=4), parameter :: vnl_varNameList(vnl_numvarmax) =  &
109       (/ vnl_varNameList3D, vnl_varNameList2D, vnl_varNameListOther /)
110  character(len=4), parameter :: varLevelList   (vnl_numvarmax) =  &
111       (/ varLevelList3D   , varLevelList2D   , varLevelListOther    /)
112  character(len=2), parameter :: varKindList    (vnl_numvarmax) =  &
113       (/ varKindList3D    , varKindList2D    , varKindListOther     /)
114
115  contains
116
117   !--------------------------------------------------------------------------
118   ! vnl_varListIndex3d
119   !--------------------------------------------------------------------------
120    function vnl_varListIndex3d(varName) result(listIndex)
121      !
122      !:Purpose: To get the 3d list index from the variable name
123      !
124      implicit none
125
126      ! Arguments:
127      character(len=*), intent(in) :: varName
128      ! Result:
129      integer                      :: listIndex
130      
131      ! Locals:
132      integer                      :: jvar
133
134      listIndex=-1
135      do jvar=1,vnl_numvarmax3D
136        if( varName == vnl_varNameList3d(jvar)) then
137          listIndex=jvar
138          exit
139        end if
140      end do
141
142      if( listIndex <= 0 ) then
143        call utl_abort('vnl_varListIndex3D: Unknown variable name! ' // varName)
144      end if
145
146    end function vnl_varListIndex3d
147
148   !--------------------------------------------------------------------------
149   ! vnl_varListIndex2d
150   !--------------------------------------------------------------------------
151    function vnl_varListIndex2d(varName) result(listIndex)
152      !
153      !:Purpose: To get the 2d list index from the variable name
154      !
155      implicit none
156
157      ! Arguments:
158      character(len=*), intent(in) :: varName
159      ! Result:
160      integer                      :: listIndex
161
162      ! Locals:
163      integer                      :: jvar
164
165      listIndex=-1
166      do jvar = 1, vnl_numvarmax2D
167        if( varName == vnl_varNameList2d(jvar) ) then 
168          listIndex=jvar
169          exit
170        end if
171      end do
172
173      if( listIndex <= 0 ) then
174        call utl_abort('vnl_varListIndex2D: Unknown variable name! ' // varName)
175      end if
176
177    end function vnl_varListIndex2d
178
179   !--------------------------------------------------------------------------
180   ! vnl_varListIndexOther
181   !--------------------------------------------------------------------------
182    function vnl_varListIndexOther(varName) result(listIndex)
183      !
184      !:Purpose: To get the "Other" list index from the variable name
185      !
186      implicit none
187
188      ! Arguments:
189      character(len=*), intent(in) :: varName
190      ! Result:
191      integer                      :: listIndex
192
193      ! Locals:
194      integer                      :: jvar
195
196      listIndex=-1
197      do jvar = 1, vnl_numvarmaxOther
198        if( varName == vnl_varNameListOther(jvar) ) then 
199          listIndex=jvar
200          exit
201        end if
202      end do
203
204      if(listIndex <= 0) then
205        call utl_abort('vnl_varListIndexOther: Unknown variable name! ' // varName)
206      end if
207
208    end function vnl_varListIndexOther
209
210   !--------------------------------------------------------------------------
211   ! vnl_varListIndex
212   !--------------------------------------------------------------------------
213    function vnl_varListIndex(varName) result(listIndex)
214      !
215      !:Purpose: To get the varlist index from the variable name
216      !
217
218      implicit none
219
220      ! Arguments:
221      character(len=*), intent(in) :: varName
222      ! Result:
223      integer                      :: listIndex
224
225      ! Local:
226      integer                      :: jvar
227
228      listIndex=-1
229      do jvar=1,vnl_numvarmax
230        if(varName == vnl_varNameList(jvar)) then 
231          listIndex=jvar
232          exit
233        end if
234      end do
235
236      if(listIndex <= 0) then
237        call utl_abort('vnl_varListIndex: Unknown variable name! ' // varName)
238      end if
239
240    end function vnl_varListIndex
241
242    !--------------------------------------------------------------------------
243    ! vnl_varnameIsValid
244    !--------------------------------------------------------------------------
245    function vnl_varnameIsValid(varName) result(isValid)
246      !
247      !:Purpose: Check if the supplied variable name is known by MIDAS.
248      !
249      implicit none
250      
251      ! Arguments:
252      character(len=*), intent(in) :: varName
253      ! Result:
254      logical                      :: isValid
255      
256      ! Local::
257      integer                      :: varIndex
258
259      isValid = .false.
260      do varIndex = 1, vnl_numvarmax
261        if(varName == vnl_varNameList(varIndex)) then 
262          isValid = .true.
263          exit
264        end if
265      end do
266
267    end function vnl_varnameIsValid
268
269    !--------------------------------------------------------------------------
270    ! vnl_varnameFromVarnum
271    !--------------------------------------------------------------------------
272    function vnl_varnameFromVarnum( varNumber, varNumberChm_opt, modelName_opt ) result(varName)
273      !
274      !:Purpose: To get the variable name from the variable number
275      !
276      implicit none
277
278      ! Arguments:
279      integer,                    intent(in) :: varNumber
280      integer,          optional, intent(in) :: varNumberChm_opt
281      character(len=*), optional, intent(in) :: modelName_opt
282      ! Result:
283      character(len=4)    :: varName
284
285      varName='    '
286      select case (varNumber)
287      case ( BUFR_NEUU, BUFR_NEUS, BUFR_NEAL )
288        varName='UU'
289      case( BUFR_NEVV, BUFR_NEVS )
290        varName='VV'
291      case( BUFR_NETT, BUFR_NETS )
292        varName='TT'
293      case( BUFR_NEDZ, BUFR_NEGZ )
294        varName='Z_T'
295      case( BUFR_NEHU, BUFR_NEHS, BUFR_NEES, BUFR_NESS )
296        varName='HU'
297      case( BUFR_NEPS, BUFR_NEPN )
298        varName='P0'
299      case ( BUFR_NERF, BUFR_NEBD, BUFR_NEZD )
300        varName='TT'   ! temporarily associate refractivity and ZTD with temperature
301      case ( BUFR_NEDW )
302        varName='DW'
303      case ( BUFR_SST )
304        varname='TG'
305      case ( BUFR_ICEC, BUFR_ICEP, BUFR_ICEV, BUFR_ICES )
306        varname='GL'
307      case ( bufr_vis )
308        varname='VIS'
309      case ( bufr_logVis )
310        varname='LVIS'
311      case ( bufr_radarPrecip )
312        varname='PR'
313      case ( bufr_logRadarPrecip )
314        varname='LPR'
315      case ( bufr_gust )
316        varname='WGE'
317      case ( bufr_riverFlow )
318        varname='QO1'
319      case ( BUFR_NEFS, bufr_radvel) 
320        varname='UV'
321      case default
322        !
323        ! Search for constituents. Identification depends on value and presence of second parameter.
324        !
325        if (present(varNumberChm_opt)) then 
326           select case (varNumberChm_opt)
327              case(BUFR_NECH_O3)
328                 varname='TO3' 
329              case(BUFR_NECH_H2O)
330                 varname='HU'
331              case(BUFR_NECH_CH4)
332                 varname='TCH4'
333              case(BUFR_NECH_CO2)
334                 varname='TCO2'
335              case(BUFR_NECH_CO)
336                 varname='TCO'
337              case(BUFR_NECH_NO2)
338                 varname='TNO2'
339              case(BUFR_NECH_N2O)
340                 varname='TN2O' 
341              case(BUFR_NECH_NO)
342                 varname='TNO'
343              case(BUFR_NECH_HCHO)
344                 varname='THCH'
345              case(BUFR_NECH_SO2)
346                 varname='TSO2'
347              case(BUFR_NECH_NH3)
348                 varname='TNH3'
349              case(BUFR_NECH_PM25)
350                 varname='AF'
351              case(BUFR_NECH_PM10)
352                 varname='AC'
353              case default
354                 call utl_abort( 'vnl_varnameFromVarnum: Unknown variable number! ' // &
355                   utl_str(varNumber) // ', ' // utl_str(varNumberChm_opt) )
356           end select
357           if (present(modelName_opt)) then
358             if (trim(modelName_opt) == 'GEM') then
359               select case (varNumberChm_opt)
360                 case(BUFR_NECH_O3)
361                   varname='O3L'  
362                 case(BUFR_NECH_CH4)
363                   varname='CH4L'
364                 case(BUFR_NECH_N2O)
365                   varname='N2OL'
366                 case default
367                   call utl_abort( 'vnl_varnameFromVarnum: Unknown variable number or model! ' // &
368                      utl_str(varNumber) // ', ' // utl_str(varNumberChm_opt) // ', ' // trim(modelName_opt) )
369               end select
370             end if
371           end if      
372        else
373           write(*,*) 'vnl_varnameFromVarnum: Unknown variable number! ',varNumber
374           call utl_abort('vnl_varnameFromVarnum')
375        end if 
376      end select
377
378    end function vnl_varnameFromVarnum
379
380   !--------------------------------------------------------------------------
381   ! vnl_varnumFromVarName
382   !--------------------------------------------------------------------------
383    function vnl_varnumFromVarName(varName,varKind_opt) result(varNumber)
384      !
385      !:Purpose: Identifies varNumber from varName for use in assimilating
386      !           obs in the CH family.   
387      !           Here, for weather variables, there is a 1-1 association between
388      !           a variable name and an observation unit.
389      !           So one must provide the name directly associated to a single
390      !           BUFR code.
391      !           As such, weather variable varNames may not necessarily be a
392      !           member of the vnl_varNameList for this routine only.
393      !   
394      !           For constituents, the varNumber refers only to the field/
395      !           variable and not units. As consequence, there is a unique
396      !           pairing of varNumbers with the varNames from vnl_VarNameList.
397      !
398      implicit none
399
400      ! Arguments:
401      character(len=*),            intent(in) :: varName
402      character(len=*), optional,  intent(in) :: varKind_opt
403      ! Result:
404      integer    :: varNumber
405      
406      varNumber=0
407      select case (varName)
408      
409      ! Weather variables. Must provide name directly associated to a single
410      ! BUFR code. As such, the varName may not necessarily be a member of the
411      ! vnl_varNameList for this routine only.
412
413      case('UU')
414        varNumber=BUFR_NEUU
415      case('US')
416        varNumber=BUFR_NEUS
417      case('VV')
418        varNumber=BUFR_NEVV
419      case('VS')
420        varNumber=BUFR_NEVS
421      case('TT')
422        varNumber=BUFR_NETT
423      case('TS')
424        varNumber=BUFR_NETS
425      case('RF')            
426        varNumber=BUFR_NERF
427      case('BD')
428        varNumber=BUFR_NEBD
429      case('ZD')
430        varNumber=BUFR_NEZD
431      case('GZ')
432        varNumber=BUFR_NEGZ
433      case('DZ')
434        varNumber=BUFR_NEDZ
435      case('UV')
436        varNumber=BUFR_NEFS 
437      case('HU')
438        if (present(varKind_opt)) then
439           if (varKind_opt == 'CH') then
440              varNumber=BUFR_NECH_H2O
441           else
442              varNumber=BUFR_NEHU
443           end if
444         else
445            varNumber=BUFR_NEHU
446         end if
447      case('HS')
448        varNumber=BUFR_NEHS
449      case('SS')
450        varNumber=BUFR_NESS
451      case('P0','PS')
452        varNumber=BUFR_NEPS
453      case('PN')
454        varNumber=BUFR_NEPN
455      case('DW')
456        varNumber=BUFR_NEDW
457      case('WGE')
458        varNumber=bufr_gust
459      case('LVIS')
460        varNumber=bufr_logVis
461      case('VIS')
462        varNumber=bufr_vis
463
464      ! Atmospheric constituents other than HU
465      case('TO3','O3L')
466        varNumber=BUFR_NECH_O3
467      case('TH2O')
468        varNumber=BUFR_NECH_H2O
469      case('TCH4','CH4L')
470        varNumber=BUFR_NECH_CH4
471      case('TCO2')
472        varNumber=BUFR_NECH_CO2
473      case('TCO','ECO')
474        varNumber=BUFR_NECH_CO
475      case('TNO2','ENO2')
476        varNumber=BUFR_NECH_NO2
477      case('TN2O','N2OL')
478        varNumber=BUFR_NECH_N2O
479      case('TNO')
480        varNumber=BUFR_NECH_NO
481      case('THCH','EHCH')
482        varNumber=BUFR_NECH_HCHO
483      case('TSO2','ESO2')
484        varNumber=BUFR_NECH_SO2
485      case('TNH3','ENH3')
486        varNumber=BUFR_NECH_NH3
487      case('AF')
488        varNumber=BUFR_NECH_PM25
489      case('AC')
490        varNumber=BUFR_NECH_PM10
491        
492      case default
493         call utl_abort('vnl_varnumFromVarName: Unknown variable name ' // trim(varName) )
494      end select
495      
496    end function vnl_varnumFromVarname
497
498   !--------------------------------------------------------------------------
499   ! vnl_varLevelFromVarname
500   !--------------------------------------------------------------------------
501    function vnl_varLevelFromVarname(varName) result(varLevel)
502      !
503      !:Purpose: To get variable level list from variable name 
504      !
505      implicit none
506
507      ! Arguments:
508      character(len=*), intent(in)   :: varName
509      ! Result:
510      character(len=4)               :: varLevel
511
512      ! Locals:
513      integer                :: nulnam, ierr
514      integer, external      :: fnom, fclos
515      logical, save          :: firstTime = .true.
516
517      ! Namelist variables
518      character(len=4), save :: forceSfcOnly(vnl_numVarMax) ! List of 3D variable names only allocated at the surface
519
520      NAMELIST /namvnl/forceSfcOnly
521
522      if (firstTime) then
523        firstTime = .false.
524        ! default values (not a valid variable name)
525        forceSfcOnly(:) = 'XXXX'
526
527        if (utl_isNamelistPresent('namvnl','./flnml')) then
528          nulnam = 0
529          ierr = fnom(nulnam,'./flnml','FTN+SEQ+R/O',0)
530          read (nulnam, nml = NAMVNL, iostat = ierr)
531          if ( ierr /= 0 ) call utl_abort('vnl_varLevelFromVarname: Error reading namelist')
532          if ( mmpi_myid == 0 ) write(*,nml=namvnl)
533          ierr = fclos(nulnam)
534        else
535          write(*,*)
536          write(*,*) 'vnl_varLevelFromVarname: namvnl is missing in the namelist. The default value will be taken.'
537        end if
538      end if
539
540      varLevel = varLevelList(vnl_varListIndex(varName))
541      if (any(forceSfcOnly(:) == varName)) then
542        if (varLevel == 'TH') then
543          varLevel = 'SFTH'
544        else if (varLevel == 'MM') then
545          varLevel = 'SFMM'
546        else
547          call utl_abort('vnl_varLevelFromVarname: something is wrong')
548        end if
549      end if
550
551    end function vnl_varLevelFromVarname
552
553   !--------------------------------------------------------------------------
554   ! vnl_varLevelFromVarnum
555   !--------------------------------------------------------------------------
556    function vnl_varLevelFromVarnum(varNumber,varNumberChm_opt,modelName_opt) result(varLevel)
557      !
558      !:Purpose: To get variable level list from the variable number 
559      !
560      implicit none
561
562      ! Arguments:
563      integer,                    intent(in) :: varNumber
564      integer,          optional, intent(in) :: varNumberChm_opt
565      character(len=*), optional, intent(in) :: modelName_opt
566      ! Result:
567      character(len=4)              :: varLevel
568
569      ! Locals:
570      character(len=4)              :: varName
571
572      varName = vnl_varnameFromVarnum(varNumber,varNumberChm_opt=varNumberChm_opt,modelName_opt=modelName_opt)
573      varLevel = varLevelList(vnl_varListIndex(varName))
574
575    end function vnl_varLevelFromVarnum
576
577   !--------------------------------------------------------------------------
578   ! vnl_varKindFromVarname
579   !--------------------------------------------------------------------------
580    function vnl_varKindFromVarname(varName) result(varKind)
581      !
582      !:Purpose: To get variable kind list from the variable number 
583      !
584      implicit none
585
586      ! Arguments:
587      character(len=*), intent(in) :: varName
588      ! Result:
589      character(len=2) :: varKind
590      
591      varKind = varKindList(vnl_varListIndex(varName))
592
593    end function vnl_varKindFromVarname
594
595   !--------------------------------------------------------------------------
596   ! vnl_varNamesFromExistList
597   !--------------------------------------------------------------------------
598    subroutine vnl_varNamesFromExistList(varNames,varExistList)
599      !
600      !:Purpose: To get variable names from the variable existList 
601      !
602
603      implicit none
604
605      ! Arguments:
606      logical,                   intent(in)  :: varExistList(:)
607      character(len=4), pointer, intent(out) :: varNames(:)
608
609      ! Local:
610      integer :: varIndex, numFound
611
612      if (associated(varNames)) then
613        call utl_abort('vnl_varNamesFromExistList: varNames must be NULL pointer on input')
614      end if
615
616      numFound = 0
617      do varIndex = 1, vnl_numvarmax
618        if ( varExistList(varIndex) ) numFound = numFound + 1
619      end do
620      allocate(varNames(numFound))
621
622      numFound = 0
623      do varIndex = 1, vnl_numvarmax
624        if ( varExistList(varIndex) ) then
625          numFound = numFound + 1
626          varNames(numFound) = vnl_varNameList(varIndex)
627        end if
628      end do
629
630    end subroutine vnl_varNamesFromExistList
631 
632   !--------------------------------------------------------------------------
633   ! vnl_varMassFromVarNum
634   !--------------------------------------------------------------------------
635    function vnl_varMassFromVarNum(varNumber) result(varMass)
636      !
637      !:Purpose: Identifies constituent molar mass from varNum for use in conversions for the CH family.   
638      !
639      implicit none
640
641      ! Arguments:
642      integer, intent(in) :: varNumber
643      ! Result:
644      real(8)             :: varMass
645
646      if ( varNumber == BUFR_NECH_O3 ) then
647        varMass = MPC_MOLAR_MASS_O3_R8
648      else if ( varNumber == BUFR_NECH_H2O ) then
649        varMass = MPC_MOLAR_MASS_VAPOUR_R8
650      else if ( varNumber == BUFR_NECH_CH4 ) then
651        varMass = MPC_MOLAR_MASS_CH4_R8
652      else if ( varNumber == BUFR_NECH_CO2 ) then
653        varMass = MPC_MOLAR_MASS_CO2_R8
654      else if ( varNumber == BUFR_NECH_CO ) then
655        varMass = MPC_MOLAR_MASS_CO_R8
656      else if ( varNumber == BUFR_NECH_NO2 ) then
657        varMass = MPC_MOLAR_MASS_NO2_R8
658      else if ( varNumber == BUFR_NECH_N2O ) then
659        varMass = MPC_MOLAR_MASS_N2O_R8
660      else if ( varNumber == BUFR_NECH_HCHO ) then
661        varMass = MPC_MOLAR_MASS_HCHO_R8
662      else if ( varNumber == BUFR_NECH_SO2 ) then
663        varMass = MPC_MOLAR_MASS_SO2_R8
664      else if ( varNumber == BUFR_NECH_NH3 ) then
665        varMass = MPC_MOLAR_MASS_NH3_R8
666      else if ( varNumber == BUFR_NECH_NO ) then
667        varMass = MPC_MOLAR_MASS_NO_R8
668      else if ( varNumber == BUFR_NECH_PM25 ) then
669        varMass = 1.0d0 ! no scaling
670      else if ( varNumber == BUFR_NECH_PM10 ) then
671        varMass = 1.0d0 ! no scaling
672      else
673        call utl_abort('vnl_varMassFromVarNum: Constituent id number ' // &
674                       utl_str(varNumber) // ' not recognized' )
675      end if
676      
677    end function vnl_varMassFromVarNum
678
679   !--------------------------------------------------------------------------
680   ! vnl_varMassFromVarName
681   !--------------------------------------------------------------------------
682    function vnl_varMassFromVarName(varName) result(varMass)
683      !
684      !:Purpose: Identifies constituent molar mass from varName for use in conversions for the CH family.   
685      !
686      implicit none
687
688      ! Arguments:
689      character(len=*),  intent(in) :: varName
690      ! Result:
691      real(8)                       :: varMass
692
693      if ( varName == 'TO3' .or. varName == 'O3L'  ) then
694        varMass = MPC_MOLAR_MASS_O3_R8
695      else if ( varName == 'LQ' .or.  varName == 'HU' ) then
696        varMass = MPC_MOLAR_MASS_VAPOUR_R8
697      else if ( varName == 'TCH4' .or. varName == 'CH4L'   ) then
698        varMass = MPC_MOLAR_MASS_CH4_R8
699      else if ( varName == 'TCO2' ) then
700        varMass = MPC_MOLAR_MASS_CO2_R8
701      else if ( varName == 'TCO' ) then
702        varMass = MPC_MOLAR_MASS_CO_R8
703      else if ( varName == 'TNO2' ) then
704        varMass = MPC_MOLAR_MASS_NO2_R8
705      else if ( varName == 'TN2O' .or. varName == 'N2OL'   ) then
706        varMass = MPC_MOLAR_MASS_N2O_R8
707      else if ( varName == 'THCH' ) then
708        varMass = MPC_MOLAR_MASS_HCHO_R8
709      else if ( varName == 'TSO2' ) then
710        varMass = MPC_MOLAR_MASS_SO2_R8
711      else if ( varName == 'TNH3' ) then
712        varMass = MPC_MOLAR_MASS_NH3_R8
713      else if ( varName == 'TNO' ) then
714        varMass = MPC_MOLAR_MASS_NO_R8
715      else if ( varName == 'AF' ) then
716        varMass = 1.0d0 ! no scaling
717      else if ( varName == 'AC' ) then
718        varMass = 1.0d0 ! no scaling
719      else
720        call utl_abort('vnl_varMassFromVarName: Molar mass not found for varName ' // &
721                       trim(varName) )
722      end if
723      
724    end function vnl_varMassFromVarName
725
726    !--------------------------------------------------------------------------
727    ! vnl_isPhysicsVar
728    !--------------------------------------------------------------------------
729    function vnl_isPhysicsVar(varName) result(isPhysicsVar)
730      !
731      !:Purpose: Signals if variable is expected to be on the "physics" grid.
732      !
733      implicit none
734
735      ! Arguments:
736      character(len=*),  intent(in) :: varName
737      ! Result:
738      logical                       :: isPhysicsVar
739
740      select case (trim(varName))
741      case ( 'I0','I1','I2','I3','I4','I5','I6','I7','I8','I9', &
742             'DN','FB','FI','PR','LPR' )
743        isPhysicsVar = .true.
744      case default
745        isPhysicsVar = .false.
746      end select
747
748    end function vnl_isPhysicsVar
749
750    !-----------------------------------------------------------------------
751    ! vnl_isCloudVar
752    !----------------------------------------------------------------------
753    function vnl_isCloudVar(varName) result(isCloud)
754      !
755      !:Purpose: determine if varName is cloud variable.
756      !
757      implicit none
758  
759      ! Arguments:
760      character(len=*), intent(in) :: varName
761      ! Result:
762      logical                      :: isCloud
763  
764      ! Locals:
765      integer :: varNameIndex
766  
767      isCloud = .false.
768      do varNameIndex = 1, vnl_numvarmaxCloud
769        if (trim(varName) == trim(vnl_varNameListCloud(varNameIndex))) then
770          isCloud = .true.
771          return
772        end if
773      end do
774  
775    end function vnl_isCloudVar
776
777    !--------------------------------------------------------------------------
778    ! vnl_addToVarNames
779    !--------------------------------------------------------------------------
780    function vnl_addToVarNames(varNamesIn,varNameToAdd) result(varNamesOut)
781      !
782      !:Purpose: Add an additional varName to an existing list of varNames
783      !
784      implicit none
785
786      ! Arguments:
787      character(len=*),  intent(in) :: varNamesIn(:)
788      character(len=*),  intent(in) :: varNameToAdd
789      ! Result:
790      character(len=4), pointer     :: varNamesOut(:)
791
792      ! Locals:
793      integer :: lenVarNames
794
795      lenVarNames = size(varNamesIn)
796      allocate(varNamesOut(lenVarNames+1))
797
798      varNamesOut(1:lenVarNames) = varNamesIn(:)
799      varNamesOut(lenVarNames+1) = varNameToAdd
800
801    end function vnl_addToVarNames
802
803end module varNameList_mod