Ignore:
Timestamp:
Feb 24, 2026, 9:59:11 AM (5 days ago)
Author:
emillour
Message:

Generic PCM:
OpenMP bug fix in "rad_correlatedk_stellar_spectrum", all intermediate computations
should be done by all OpenMp? threads.
While at it did some cleanup:

  • added some OMP threadprivate statements in radcommon_h.F90 (not alwyas necessary, but best practice is that all saved variables be threadprivate)
  • turned rad_blackbody.F into a module and modernized routines.
  • use clearer stategy (wrt OpenMP) in "rad_correlatedk_init_thermal.F90" and "rad_correlatedk_init_stellar.F90": have the master read in file and data and then broadcast to all cores.

EM

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/libf/phygeneric/rad_blackbody.F

    r4077 r4081  
     1      module rad_blackbody_mod
     2     
     3      implicit none
     4     
     5      contains
     6     
    17      subroutine rad_blackbody_planck_law_wavelength(blalong,blat,blae)
    28
    3       implicit double precision (a-h,o-z)
     9      implicit none
     10     
     11      double precision, intent(in) :: blalong
     12      double precision, intent(in) :: blat
     13      double precision, intent(out) :: blae
     14     
     15!      double precision :: sigma,pi,c0,h,cbol,rind,c,c1,c2
    416
    517      ! physical constants
    6       sigma=5.670374D-8
    7       pi=datan(1.d0)*4.d0
    8       c0=2.997925d+08
    9       h=6.62607d-34
    10       cbol=1.380649d-23
    11       rind=1.d0
    12       c=c0/rind
    13       c1=h*(c**2)
    14       c2=h*c/cbol
     18      double precision, parameter :: sigma=5.670374D-8
     19      double precision, parameter :: pi=atan(1.d0)*4.d0
     20      double precision, parameter :: c0=2.997925d+08
     21      double precision, parameter :: h=6.62607d-34
     22      double precision, parameter :: cbol=1.380649d-23
     23      double precision, parameter :: rind=1.d0
     24      double precision, parameter :: c=c0/rind
     25      double precision, parameter :: c1=h*(c**2)
     26      double precision, parameter :: c2=h*c/cbol
    1527
    1628
    17       blae=2.d0*pi*c1/blalong**5/(dexp(c2/blalong/blat)-1.d0)
     29      blae=2.d0*pi*c1/blalong**5/(exp(c2/blalong/blat)-1.d0)
    1830
    19 
    20       return
    21       end
     31      end subroutine rad_blackbody_planck_law_wavelength
    2232
    2333      subroutine rad_blackbody_planck_law_wavenumber(blalong,blat,blae)
    2434
    25       implicit double precision (a-h,o-z)
     35      implicit none
     36
     37      double precision, intent(in) :: blalong
     38      double precision, intent(in) :: blat
     39      double precision, intent(out) :: blae
     40     
     41!      double precision :: sigma,pi,c0,h,cbol,rind,c,c1,c2
    2642
    2743      ! physical constants
    28       sigma=5.670374D-8
    29       pi=datan(1.d0)*4.d0
    30       c0=2.997925d+08
    31       h=6.62607d-34
    32       cbol=1.380649d-23
    33       rind=1.d0
    34       c=c0/rind
    35       c1=h*(c**2)
    36       c2=h*c/cbol
     44      double precision, parameter :: sigma=5.670374D-8
     45      double precision, parameter :: pi=atan(1.d0)*4.d0
     46      double precision, parameter :: c0=2.997925d+08
     47      double precision, parameter :: h=6.62607d-34
     48      double precision, parameter :: cbol=1.380649d-23
     49      double precision, parameter :: rind=1.d0
     50      double precision, parameter :: c=c0/rind
     51      double precision, parameter :: c1=h*(c**2)
     52      double precision, parameter :: c2=h*c/cbol
    3753
    3854
    39       blae=2.d0*pi*c1*blalong**3/(dexp(c2*blalong/blat)-1.d0)
     55      blae=2.d0*pi*c1*blalong**3/(exp(c2*blalong/blat)-1.d0)
    4056
    41 
    42       return
    43       end
     57      end subroutine rad_blackbody_planck_law_wavenumber
     58     
     59      end module rad_blackbody_mod
Note: See TracChangeset for help on using the changeset viewer.