SCALE-RM
Data Types | Functions/Subroutines
scale_atmos_thermodyn Module Reference

module ATMOSPHERE / Thermodynamics More...

Functions/Subroutines

subroutine, public atmos_thermodyn_setup
 Setup. More...
 
subroutine, public atmos_thermodyn_tempre (temp, pres, Ein, dens, qdry, q, CVq, Rq)
 
subroutine, public atmos_thermodyn_tempre2 (temp, pres, dens, pott, qdry, q, CPq, Rq)
 

Detailed Description

module ATMOSPHERE / Thermodynamics

Description
Thermodynamics module
Author
Team SCALE
History
  • 2011-10-24 (T.Seiki) [new] Import from NICAM
  • 2012-02-10 (H.Yashiro) [mod] Reconstruction
  • 2012-03-23 (H.Yashiro) [mod] Explicit index parameter inclusion
  • 2012-12-22 (S.Nishizawa) [mod] Use thermodyn macro set

Function/Subroutine Documentation

◆ atmos_thermodyn_setup()

subroutine, public scale_atmos_thermodyn::atmos_thermodyn_setup ( )

Setup.

Definition at line 122 of file scale_atmos_sub_thermodyn.F90.

References scale_grid_index::ieb, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, scale_grid_index::ks, and scale_tracer::qa.

Referenced by mod_rm_driver::scalerm(), and mod_rm_prep::scalerm_prep().

122  implicit none
123  !---------------------------------------------------------------------------
124 
125  if( io_l ) write(io_fid_log,*)
126  if( io_l ) write(io_fid_log,*) '++++++ Module[THERMODYN] / Categ[ATMOS SHARE] / Origin[SCALElib]'
127  if( io_l ) write(io_fid_log,*) '*** No namelists.'
128 
129  return
Here is the caller graph for this function:

◆ atmos_thermodyn_tempre()

subroutine, public scale_atmos_thermodyn::atmos_thermodyn_tempre ( real(rp), dimension(ka,ia,ja), intent(out)  temp,
real(rp), dimension(ka,ia,ja), intent(out)  pres,
real(rp), dimension (ka,ia,ja), intent(in)  Ein,
real(rp), dimension(ka,ia,ja), intent(in)  dens,
real(rp), dimension(ka,ia,ja), intent(in)  qdry,
real(rp), dimension (ka,ia,ja,qa), intent(in)  q,
real(rp), dimension (qa), intent(in)  CVq,
real(rp), dimension (qa), intent(in)  Rq 
)

Definition at line 951 of file scale_atmos_sub_thermodyn.F90.

References scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, scale_grid_index::ks, and scale_tracer::qa.

951  implicit none
952 
953  real(RP), intent(out) :: temp(KA,IA,JA) ! temperature
954  real(RP), intent(out) :: pres(KA,IA,JA) ! pressure
955  real(RP), intent(in) :: Ein (KA,IA,JA) ! internal energy
956  real(RP), intent(in) :: dens(KA,IA,JA) ! density
957  real(RP), intent(in) :: qdry(KA,IA,JA) ! dry concentration
958  real(RP), intent(in) :: q (KA,IA,JA,QA) ! water concentration
959  real(RP), intent(in) :: CVq (QA) ! specific heat
960  real(RP), intent(in) :: Rq (QA) ! gas constant
961 
962  real(RP) :: cv, Rmoist
963 
964  integer :: i, j, k, iqw
965  !---------------------------------------------------------------------------
966 
967  !$omp parallel do private(i,j,k,iqw,cv,Rmoist) OMP_SCHEDULE_ collapse(2)
968  do j = jsb, jeb
969  do i = isb, ieb
970  do k = ks, ke
971  cv = qdry(k,i,j) * cvdry
972  rmoist = qdry(k,i,j) * rdry
973  do iqw =1, qa
974  cv = cv + q(k,i,j,iqw) * cvq(iqw)
975  rmoist = rmoist + q(k,i,j,iqw) * rq(iqw)
976  enddo
977 
978  temp(k,i,j) = ein(k,i,j) / cv
979  pres(k,i,j) = dens(k,i,j) * rmoist * temp(k,i,j)
980  enddo
981  enddo
982  enddo
983 
984  return
integer, public jeb
integer, public ieb
integer, public isb
integer, public jsb

◆ atmos_thermodyn_tempre2()

subroutine, public scale_atmos_thermodyn::atmos_thermodyn_tempre2 ( real(rp), dimension(ka,ia,ja), intent(out)  temp,
real(rp), dimension(ka,ia,ja), intent(out)  pres,
real(rp), dimension(ka,ia,ja), intent(in)  dens,
real(rp), dimension(ka,ia,ja), intent(in)  pott,
real(rp), dimension(ka,ia,ja), intent(in)  qdry,
real(rp), dimension (ka,ia,ja,qa), intent(in)  q,
real(rp), dimension (qa), intent(in)  CPq,
real(rp), dimension (qa), intent(in)  Rq 
)

Definition at line 992 of file scale_atmos_sub_thermodyn.F90.

References scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, scale_grid_index::ks, and scale_tracer::qa.

992  implicit none
993 
994  real(RP), intent(out) :: temp(KA,IA,JA) ! temperature
995  real(RP), intent(out) :: pres(KA,IA,JA) ! pressure
996  real(RP), intent(in) :: dens(KA,IA,JA) ! density
997  real(RP), intent(in) :: pott(KA,IA,JA) ! potential temperature
998  real(RP), intent(in) :: qdry(KA,IA,JA) ! dry concentration
999  real(RP), intent(in) :: q (KA,IA,JA,QA) ! water concentration
1000  real(RP), intent(in) :: CPq (QA) ! specific heat
1001  real(RP), intent(in) :: Rq (QA) ! gas constant
1002 
1003  real(RP) :: Rmoist, cp
1004 
1005  integer :: i, j, k, iqw
1006  !---------------------------------------------------------------------------
1007 
1008  !$omp parallel do private(i,j,k,iqw,cp,Rmoist) OMP_SCHEDULE_ collapse(2)
1009  do j = jsb, jeb
1010  do i = isb, ieb
1011  do k = ks, ke
1012  cp = qdry(k,i,j) * cpdry
1013  rmoist = qdry(k,i,j) * rdry
1014  do iqw = 1, qa
1015  cp = cp + q(k,i,j,iqw) * cpq(iqw)
1016  rmoist = rmoist + q(k,i,j,iqw) * rq(iqw)
1017  enddo
1018  calc_pre(pres(k,i,j), dens(k,i,j), pott(k,i,j), rmoist, cp, pre00)
1019 
1020  temp(k,i,j) = pres(k,i,j) / ( dens(k,i,j) * rmoist )
1021  enddo
1022  enddo
1023  enddo
1024 
1025  return
integer, public jeb
integer, public ieb
integer, public isb
integer, public jsb