Ignore:
Timestamp:
Feb 28, 2026, 6:36:58 PM (6 weeks ago)
Author:
emillour
Message:

Venus PCM:
Code cleanup: remove "hedin.h" and make variables included there be module
variables of compo_hedin_mod2.F90. Turn "nirdat.h" into module nirdata.F90
and incroporate initialization routine nir_leedat.F in it.
EM

File:
1 moved

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.VENUS/libf/phyvenus/nirdata.F90

    r4091 r4092  
     1module nirdata_mod
    12
    2       integer npres                ! Number of pressures in NIR correction
    3       parameter (npres=42)         ! table
     3implicit none
    44
    5       common /NIRdata/ pres1d,corgcm,oco21d,alfa,p1999
    6       real    pres1d(npres)
    7       real    corgcm(npres)
    8       real    oco21d(npres),alfa(npres),p1999(npres)
     5! Number of pressures in NIR correction table
     6integer,parameter :: npres=42   
     7
     8real,save,protected :: pres1d(npres)
     9real,save,protected :: corgcm(npres)
     10real,save,protected :: oco21d(npres)
     11real,save,protected :: alfa(npres)
     12real,save,protected :: p1999(npres)
     13!$OMP THREADPRIVATE(pres1d,corgcm,oco21d,alfa,p1999)
     14
     15contains
     16
     17  subroutine nir_leedat
     18
     19!       reads parameters for NIR NLTE calculation   
     20
     21!       nov 2011    fgg+malv    first version               
     22!***********************************************************************
     23
     24    use mod_phys_lmdz_para, only: is_master, bcast
     25
     26    implicit none
     27 
     28    integer :: ind
     29    character (len=100) :: datafile="HIGHATM"
    930   
     31    if (is_master) then ! only the master needs read the file
     32      open(43,file=trim(datafile)//'/NIRcorrection_feb2011.dat', &
     33              status='old')
     34      do ind=1,9
     35         read(43,*)
     36      enddo
     37     
     38      do ind=1,npres
     39         read(43,*)pres1d(ind),corgcm(ind),oco21d(ind),p1999(ind), &
     40                   alfa(ind)
     41         ! Tabulated pressure to Pa
     42         pres1d(ind)=pres1d(ind)*100.
     43      enddo
     44     
     45      close(43)
     46   
     47    endif ! of if (is_master)
     48 
     49    ! broadcast to all cores
     50    call bcast(pres1d)
     51    call bcast(corgcm)
     52    call bcast(oco21d)
     53    call bcast(p1999)
     54    call bcast(alfa)
     55
     56  end subroutine nir_leedat
     57
     58end module nirdata_mod
Note: See TracChangeset for help on using the changeset viewer.