aboutsummaryrefslogtreecommitdiff
path: root/src/EOS_Omni_MultiVarCalls.F90
blob: 664a31f34b30e3fdb7fadbb62ab1739373ef9322 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
#include "cctk_Functions.h"

! eoskey:
! 1 --- polytropic EOS
! 2 --- gamma-law EOS
! 3 --- hybrid EOS
! 4 --- finite-T microphysical NSE EOS


subroutine EOS_Omni_EOS_short(eoskey,keytemp,rf_precision,npoints,&
     rho,eps,temp,ye,press,entropy,cs2,dedt,dpderho,dpdrhoe,munu,&
     keyerr,anyerr)

  use EOS_Omni_Module
  implicit none
  DECLARE_CCTK_PARAMETERS

  CCTK_INT, intent(in)     :: eoskey,keytemp,npoints
  CCTK_INT, intent(out)    :: keyerr(npoints)
  CCTK_INT, intent(out)    :: anyerr
  CCTK_REAL, intent(in)    :: rf_precision
  CCTK_REAL, intent(in)    :: rho(npoints),ye(npoints)
  CCTK_REAL, intent(inout) :: eps(npoints), temp(npoints)
  CCTK_REAL, intent(out)   :: press(npoints)
  CCTK_REAL, intent(inout) :: entropy(npoints)
  CCTK_REAL, intent(out)   :: cs2(npoints)
  CCTK_REAL, intent(out)   :: dedt(npoints)
  CCTK_REAL, intent(out)   :: dpderho(npoints)
  CCTK_REAL, intent(out)   :: dpdrhoe(npoints)
  CCTK_REAL, intent(out)   :: munu(npoints)

  ! local vars
  integer          :: i
  character(256)   :: warnstring
  ! temporary vars for nuc_eos
  real*8           :: xrho,xye,xtemp,xenr,xent
  real*8           :: xprs,xmunu,xcs2
  real*8           :: xdedt,xdpderho,xdpdrhoe

  if(eoskey.ne.4) then
     call CCTK_WARN(0,"EOS_Omni_EOS_short currently does not work for this eoskey")
  endif

  anyerr    = 0
  keyerr(:) = 0
  
  do i=1,npoints

     xrho = rho(i) * inv_rho_gf
     xtemp = temp(i)
     xye = ye(i)
     xenr = eps(i) * inv_eps_gf
     xent = entropy(i)
     call nuc_eos_short(xrho,xtemp,xye,xenr,xprs,& 
          xent,xcs2,xdedt,xdpderho,xdpdrhoe,xmunu,&
          keytemp,keyerr(i),rf_precision)
     
     if(keyerr(i).ne.0) then
        anyerr = 1
     endif
     
     if(keytemp.eq.1) then
        eps(i) = xenr * eps_gf
     else if(keytemp.eq.2) then
        eps(i) = xenr * eps_gf
        temp(i) = xtemp
     else
        temp(i) = xtemp
     endif
     
     press(i)     = xprs * press_gf
     entropy(i)   = xent
     cs2(i)       = xcs2
     dedt(i)      = xdedt * eps_gf
     dpderho(i)   = xdpderho * press_gf * inv_eps_gf
     dpdrhoe(i)   = xdpdrhoe * press_gf * inv_rho_gf
     munu(i)      = xmunu
     
  enddo

end subroutine EOS_Omni_EOS_short


subroutine EOS_Omni_EOS_dpderho_dpdrhoe(eoskey,keytemp,rf_precision,npoints,&
     rho,eps,temp,ye,dpderho,dpdrhoe,keyerr,anyerr)

  use EOS_Omni_Module
  implicit none
  DECLARE_CCTK_PARAMETERS

  CCTK_INT, intent(in)     :: eoskey,keytemp,npoints
  CCTK_INT, intent(out)    :: keyerr(npoints)
  CCTK_INT, intent(out)    :: anyerr
  CCTK_REAL, intent(in)    :: rf_precision
  CCTK_REAL, intent(in)    :: rho(npoints),ye(npoints)
  CCTK_REAL, intent(inout) :: eps(npoints), temp(npoints)
  CCTK_REAL, intent(out)   :: dpderho(npoints)
  CCTK_REAL, intent(out)   :: dpdrhoe(npoints)

  ! local vars
  integer          :: i
  character(256)   :: warnstring
  real*8           :: hybrid_local_gamma
  real*8           :: hybrid_local_k_cgs
  real*8           :: hybrid_dp_poly,hybrid_dp_th1,hybrid_dp_th2
  ! temporary vars for nuc_eos
  real*8           :: xrho,xye,xtemp,xenr,xent
  real*8           :: xprs,xmunu,xcs2
  real*8           :: xdedt,xdpderho,xdpdrhoe


  anyerr    = 0
  keyerr(:) = 0

  select case (eoskey)
  case (1)
     ! polytropic EOS                                                                                      
        if(keytemp.eq.1) then
           do i=1,npoints
              eps(i) = press_gf * poly_k_cgs * &
                   (rho(i)*inv_rho_gf)**(poly_gamma) / &
                   (poly_gamma - 1.0d0) / rho(i)
           enddo
        endif
        do i=1,npoints
           dpdrhoe(i) = press_gf * poly_k_cgs *  &
                poly_gamma * inv_rho_gf *        &
                (rho(i)*inv_rho_gf) ** (poly_gamma - 1.d0)
           dpderho(i) = 0.0d0
        enddo
  case (2)
     ! gamma-law EOS                                                                                       
        if(keytemp.eq.1) then
           do i=1,npoints
              eps(i) = press_gf * gl_k_cgs * &
                   (rho(i)*inv_rho_gf)**(gl_gamma) / &
                   (gl_gamma - 1.0d0) / rho(i)
           enddo
        endif
        do i=1,npoints
           dpdrhoe(i) = (gl_gamma-1.0d0) * &
                eps(i)
           dpderho(i) = (gl_gamma - 1.0d0) * &
                rho(i)
        enddo
  case (3)
     ! hybrid EOS                                                                                          
        do i=1,npoints
           if(rho(i).gt.hybrid_rho_nuc) then
              hybrid_local_gamma = hybrid_gamma2
              hybrid_local_k_cgs = hybrid_k2_cgs
           else
              hybrid_local_gamma = hybrid_gamma1
              hybrid_local_k_cgs = hybrid_k1_cgs
           endif
           hybrid_dp_poly = hybrid_local_gamma * press_gf *                 &
                hybrid_local_k_cgs * rho(i)**(hybrid_local_gamma - 1.0d0) * &
                inv_rho_gf**hybrid_local_gamma

           hybrid_dp_th1 = - hybrid_local_gamma * press_gf * hybrid_local_k_cgs * &
                (hybrid_gamma_th - 1.d0) / (hybrid_local_gamma - 1.d0) *      &
                rho(i)**(hybrid_local_gamma - 1.d0) * inv_rho_gf**hybrid_local_gamma

           hybrid_dp_th2  = (hybrid_gamma_th - 1.d0) * eps(i)                        &
                - (hybrid_gamma_th - 1.d0) * (hybrid_local_gamma - hybrid_gamma1) /  &
                (hybrid_gamma1 - 1.d0) / (hybrid_gamma2 - 1.d0) *                    &
                press_gf * hybrid_k1_cgs * inv_rho_gf**hybrid_gamma1 *               &
                hybrid_rho_nuc**(hybrid_gamma1 - 1.d0)

           dpdrhoe(i) = hybrid_dp_poly + hybrid_dp_th1 + hybrid_dp_th2
           dpderho(i) = (hybrid_gamma_th - 1.0d0) * rho(i)
        enddo
    case (4)
        do i=1,npoints
           xrho = rho(i) * inv_rho_gf
           xtemp = temp(i)
           xye = ye(i)
           xenr = eps(i) * inv_eps_gf
           call nuc_eos_dpdr_dpde(xrho,xtemp,xye,xenr,xdpdrhoe,&
                xdpderho,&
                keytemp,keyerr(i),rf_precision)

           if(keyerr(i).ne.0) then
              anyerr = 1
           endif

           if(keytemp.eq.1) then
              eps(i) = xenr * eps_gf
           else
              temp(i) = xtemp
           endif

           dpdrhoe(i) = xdpdrhoe * press_gf * inv_rho_gf
           dpderho(i) = xdpderho * press_gf * inv_eps_gf
        enddo
   case DEFAULT
      write(warnstring,*) "eoskey ",eoskey," not implemented!"
      call CCTK_WARN(0,warnstring)
   end select

 end subroutine EOS_Omni_EOS_dpderho_dpdrhoe