Ignore:
Timestamp:
Dec 12, 2025, 5:36:11 PM (3 months ago)
Author:
emillour
Message:

Generic PCM:
OpenMP bug fix in ave_stelspec.F90 (line counting strategy should only be done
by the master thread). While at it turned ave_stelspec and setspv into modules.
EM

Location:
trunk/LMDZ.GENERIC/libf/phystd
Files:
3 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/libf/phystd/ave_stelspec.F90

    r3942 r3990  
     1      module ave_stelspec_mod
     2     
     3      implicit none
     4     
     5      contains
     6     
    17      subroutine ave_stelspec(STELLAR)
    28
     
    3844      integer ifine,band
    3945
    40       real,allocatable,save :: lam(:),stel_f(:)         !read by master
     46      real,allocatable,save :: lam(:),stel_f(:) ! read by master thread
     47                                                ! but used by all threads
    4148      real lamm,lamp
    4249      real dl
     
    5057     
    5158      integer :: ios ! file opening/reading status
     59      logical :: file_exists
    5260
    5361      STELLAR(:)=0.0
     
    93101         write(*,*) 'Please use ',1,' and only ',1,' header line in ',trim(stelspec_file)
    94102
    95          ! Opening file
     103         ! Check the target file is there
    96104         file_path = trim(datadir)//'/stellar_spectra/'//stelspec_file
    97105         print*, 'stellar flux : ', file_path
    98          OPEN(UNIT=110,FILE=file_path,STATUS='old',iostat=ios)
     106         inquire(FILE=file_path,EXIST=file_exists)         
    99107   
    100          if (ios /= 0) THEN
     108         if (.not.file_exists) THEN
    101109           write(*,*)'Beware that startype is now deprecated, you should use '
    102110           write(*,*)'stelspec_file and tstellar to define the input stellar spectrum.'
     
    116124         end if
    117125
     126!$OMP MASTER
     127         ! Open the file
     128         OPEN(UNIT=110,FILE=file_path,STATUS='old',iostat=ios)
    118129         ! Get number of line in the file
    119130         READ(110,*) ! skip first line header just in case
     
    127138         READ(110,*) ! skip first line header just in case
    128139
    129 !$OMP MASTER
    130140         allocate(lam(Nfine),stel_f(Nfine))
    131141
     
    158168!$OMP BARRIER
    159169!$OMP MASTER
    160          if (allocated(lam)) deallocate(lam)
    161          if (allocated(stel_f)) deallocate(stel_f)
     170         deallocate(lam)
     171         deallocate(stel_f)
    162172!$OMP END MASTER
    163173!$OMP BARRIER         
     
    165175
    166176      end subroutine ave_stelspec
     177     
     178      end module ave_stelspec_mod
  • trunk/LMDZ.GENERIC/libf/phystd/physiq_mod.F90

    r3928 r3990  
    2424      use radcommon_h, only: sigma, glat, grav, BWNV, WNOI, DWNI, DWNV, WNOV
    2525      use suaer_corrk_mod, only: suaer_corrk
     26      use setspv_mod, only: setspv
    2627      use radii_mod, only: h2o_reffrad, co2_reffrad
    2728      use aerosol_mod, only: iniaerosol, iaero_co2, iaero_h2o
  • trunk/LMDZ.GENERIC/libf/phystd/setspv.F90

    r3893 r3990  
     1      module setspv_mod
     2     
     3      implicit none
     4     
     5      contains
     6     
    17      subroutine setspv
    28
     
    2733      use datafile_mod, only: datadir
    2834      use callkeys_mod, only: Fat1AU
     35      use ave_stelspec_mod, only: ave_stelspec
    2936
    3037      implicit none
     
    132139      print*,' '
    133140
    134       RETURN
    135141    END subroutine setspv
     142   
     143    end module setspv_mod
     144   
Note: See TracChangeset for help on using the changeset viewer.