source: trunk/LMDZ.GENERIC/libf/phystd/su_gases.F90 @ 257

Last change on this file since 257 was 253, checked in by emillour, 14 years ago

Generic GCM

  • Massive update to version 0.7

EM+RW

  • Property svn:executable set to *
File size: 2.2 KB
Line 
1subroutine su_gases
2
3  implicit none
4
5#include "gases.h"
6
7       integer ngas, igas, ierr
8
9!==================================================================
10!     
11!     Purpose
12!     -------
13!     Load atmospheric composition info
14!     
15!     Authors
16!     -------
17!     R. Wordsworth (2011)
18!     
19!==================================================================
20
21       ! load gas names from file 'gases.def'
22       open(90,file='gases.def',status='old',form='formatted',iostat=ierr)
23       if (ierr.eq.0) then
24          write(*,*) "sugases.F90: reading file gases.def"
25          read(90,*)
26          read(90,*,iostat=ierr) ngas
27          if (ierr.ne.0) then
28             write(*,*) "sugases.F90: error reading number of gases"
29             write(*,*) "   (first line of gases.def) "
30             call abort
31          endif
32          ! check that the number of gases is indeed ngasmx
33          if (ngas.ne.ngasmx) then
34             write(*,*) "sugases.F90: error, wrong number of gases:"
35             write(*,*) "ngas=",ngas," whereas ngasmx=",ngasmx
36             call abort
37          endif
38
39          do igas=1,ngas
40             read(90,*,iostat=ierr) gnom(igas)
41             if (ierr.ne.0) then
42                write(*,*) 'sugases.F90: error reading gas names...'
43                call abort
44             endif
45          enddo                  !of do igas=1,ngas
46         
47          vgas=0
48          do igas=1,ngas
49             read(90,*,iostat=ierr) gfrac(igas)
50             if (ierr.ne.0) then
51                write(*,*) 'sugases.F90: error reading gas molar fractions...'
52                call abort
53             endif
54
55             ! find variable gas (if any)
56             if(gfrac(igas).eq.-1.0)then
57                if(vgas.eq.0)then
58                   vgas=igas
59                else
60                   print*,'You seem to be choosing two variable gases'
61                   print*,'Check that gases.def is correct'
62                   call abort
63                endif
64             endif
65             
66          enddo                  !of do igas=1,ngas
67
68       else
69          write(*,*) 'Cannot find required file "gases.def"'
70          call abort
71       endif
72       close(90)
73
74end subroutine su_gases
Note: See TracBrowser for help on using the repository browser.