20 data h_alc7 / 0.0, 10.0, 20.0, 35.0, 50.0, &
21 75.0, 100.0, 125.0, 150.0, 175.0, &
22 200.0, 250.0, 300.0, 350.0, 400.0, &
23 450.0, 490.0, 525.0, 560.0, 615.0, &
24 660.0, 700.0, 750.0, 800.0, 854.0, &
25 900.0, 946.0, 971.0, 1003.0, 1032.0, &
26 1065.0, 1101.0, 1143.0, 1214.0, 1299.0, &
27 1398.0, 1520.0, 1617.0, 1722.0, 1820.0, &
28 1894.0, 1946.0, 1989.0, 2024.0, 2055.0, &
29 2083.0, 2098.0, 2110.0, 2120.0, 2126.0, &
30 2130.0, 2132.0, 2134.0, 2136.0, 2138.0, &
31 2141.9, 2145.2, 2147.1, 2148.7, 2150.1, &
32 2151.2, 2152.0, 2152.9, 2153.5, 2154.1, &
33 2154.7, 2155.5, 2156.2, 2156.8, 2157.6, &
34 2158.2, 2158.8, 2159.4, 2160.0, 2160.6, &
35 2161.3, 2162.3, 2163.3, 2164.6, 2166.0, &
36 2167.7, 2170.0, 2172.3, 2175.3, 2178.7, &
37 2182.4, 2186.1, 2191.5, 2196.5, 2201.9, &
38 2207.7, 2213.7, 2220.0, 2227.7, 2235.3, &
39 2243.3, 2251.6, 2259.4, 2269.7, 2282.4, &
40 2294.4, 2308.6, 2325.8, 2346.1, 2378.2, &
41 2418.7, 2457.9, 2495.6, 2535.5, 2578.2, &
42 2628.3, 2688.1, 2761.5, 2844.4, 2970.8, &
43 3133.8, 3425.0, 3752.7, 4295.9, 4968.6/
45 data t_alc7 / 6583, 6397, 6231, 6006, 5826, &
46 5607, 5431, 5288, 5165, 5080, &
47 5010, 4907, 4805, 4700, 4590, &
48 4485, 4435, 4410, 4400, 4435, &
49 4510, 4640, 4840, 5090, 5430, &
50 5720, 5969, 6100, 6225, 6315, &
51 6400, 6474, 6531, 6576, 6598, &
52 6610, 6623, 6633, 6643, 6652, &
53 6660, 6667, 6674, 6680, 6686, &
54 6694, 6700, 6706, 6718, 6740, &
55 6768, 6800, 6870, 6992, 7248, &
56 7950, 9115, 10980, 13200, 15760, &
57 18140, 20510, 23100, 25120, 27130, &
58 29500, 32260, 34580, 36870, 39400, &
59 41450, 43400, 45140, 46800, 48490, &
60 50140, 52690, 55020, 57790, 60790, &
61 63950, 68000, 71810, 76330, 81220, &
62 86120, 90640, 96860, 102300, 107800, &
63 113200, 118700, 124000, 130200, 135800, &
64 141600, 147200, 152300, 158600, 166000, &
65 172600, 180100, 188600, 198100, 212000, &
66 227800, 241900, 254400, 266700, 279000, &
67 292500, 307300, 324000, 341300, 365000, &
68 392000, 432900, 471600, 524300, 577400/
74 data h_fontenla / -100.0, -90.0, -80.0, -70.0, -60.0, &
75 -50.0, -40.0, -30.0, -20.0, -10.0, &
76 -0.7, 9.5, 20.0, 35.0, 50.0, &
77 75.0, 100.0, 125.0, 150.0, 175.0, &
78 200.0, 250.0, 300.0, 347.0, 398.0, &
79 447.5, 485.0, 519.0, 570.0, 625.0, &
80 696.0, 772.0, 812.0, 838.0, 874.5, &
81 899.5, 926.0, 947.0, 965.0, 1020.0, &
82 1060.0, 1120.0, 1180.0, 1245.0, 1319.0, &
83 1447.0, 1571.0, 1687.0, 1801.0, 1911.0, &
84 2013.0, 2142.0, 2263.0, 2357.0, 2425.0, &
87 data t_fontenla / 9460, 9220, 8940, 8620, 8280, &
88 7940, 7600, 7270, 6960, 6690, &
89 6490, 6310, 6150, 5950, 5780, &
90 5550, 5380, 5245, 5130, 5035, &
91 4960, 4835, 4740, 4660, 4580, &
92 4510, 4450, 4390, 4300, 4200, &
93 4080, 3955, 3890, 3850, 3800, &
94 3810, 3950, 4150, 4400, 5350, &
95 5770, 6120, 6270, 6390, 6450, &
96 6490, 6530, 6570, 6620, 6670, &
97 6720, 6770, 6830, 6910, 7030, &
104 data h_valc / -75, -50, -25, 0, 50, &
105 100, 150, 250, 350, 450, &
106 515, 555, 605, 655, 705, &
107 755, 855, 905, 980, 1065, &
108 1180, 1280, 1380, 1515, 1605, &
109 1785, 1925, 1990, 2016, 2050, &
110 2070, 2080, 2090, 2104, 2107, &
111 2109, 2113, 2115, 2120, 2129, &
112 2160, 2200, 2230, 2255, 2263, &
113 2267, 2271, 2274, 2280, 2290, &
116 data t_valc / 8320, 7610, 6910, 6420, 5840, &
117 5455, 5180, 4780, 4465, 4220, &
118 4170, 4230, 4420, 4730, 5030, &
119 5280, 5650, 5755, 5925, 6040, &
120 6150, 6220, 6280, 6370, 6440, &
121 6630, 6940, 7160, 7360, 7660, &
122 7940, 8180, 8440, 9500, 10700, &
123 12300, 18500, 21000, 22500, 23000, &
124 23500, 24000, 24200, 24500, 25500, &
125 28000, 32000, 37000, 50000, 89100, &
131 data h_hong / -3.7616
d+07, -3.7384
d+07, -3.7045
d+07, -3.6561
d+07, -3.5880
d+07,&
132 -3.4950
d+07, -3.3718
d+07, -3.2153
d+07, -3.0250
d+07, -2.8050
d+07,&
133 -2.5627
d+07, -2.3088
d+07, -2.0556
d+07, -1.8151
d+07, -1.5972
d+07,&
134 -1.4004
d+07, -1.2234
d+07, -1.0625
d+07, -9.1510
d+06, -7.7773
d+06,&
135 -6.4599
d+06, -5.1584
d+06, -3.8455
d+06, -2.4946
d+06, -1.0577
d+06,&
136 4.7767
d+05, 2.1041
d+06, 3.8237
d+06, 5.6198
d+06, 7.4817
d+06,&
137 9.3930
d+06, 1.1339
d+07, 1.3304
d+07, 1.5278
d+07, 1.7250
d+07,&
138 1.9213
d+07, 2.1163
d+07, 2.3099
d+07, 2.5019
d+07, 2.6922
d+07,&
139 2.8809
d+07, 3.0678
d+07, 3.2532
d+07, 3.4371
d+07, 3.6196
d+07,&
140 3.8009
d+07, 3.9811
d+07, 4.1605
d+07, 4.3396
d+07, 4.5186
d+07,&
141 4.6977
d+07, 4.8773
d+07, 5.0576
d+07, 5.2390
d+07, 5.4218
d+07,&
142 5.6065
d+07, 5.7934
d+07, 5.9833
d+07, 6.1766
d+07, 6.3736
d+07,&
143 6.5744
d+07, 6.7791
d+07, 6.9880
d+07, 7.2011
d+07, 7.4185
d+07,&
144 7.6403
d+07, 7.8660
d+07, 8.0954
d+07, 8.3284
d+07, 8.5647
d+07,&
145 8.8038
d+07, 9.0456
d+07, 9.2896
d+07, 9.5357
d+07, 9.7838
d+07,&
146 1.0033
d+08, 1.0285
d+08, 1.0537
d+08, 1.0791
d+08, 1.1046
d+08,&
147 1.1302
d+08, 1.1559
d+08, 1.1818
d+08, 1.2077
d+08, 1.2337
d+08,&
148 1.2599
d+08, 1.2862
d+08, 1.3126
d+08, 1.3392
d+08, 1.3659
d+08,&
149 1.3927
d+08, 1.4198
d+08, 1.4470
d+08, 1.4744
d+08, 1.5019
d+08,&
150 1.5295
d+08, 1.5571
d+08, 1.5847
d+08, 1.6121
d+08, 1.6389
d+08,&
151 1.6648
d+08, 1.6895
d+08, 1.7122
d+08, 1.7325
d+08, 1.7499
d+08,&
152 1.7642
d+08, 1.7754
d+08, 1.7838
d+08, 1.7900
d+08, 1.7944
d+08,&
153 1.7974
d+08, 1.7996
d+08, 1.8010
d+08, 1.8020
d+08, 1.8027
d+08,&
154 1.8032
d+08, 1.8035
d+08, 1.8037
d+08, 1.8039
d+08, 1.8041
d+08,&
155 1.8042
d+08, 1.8044
d+08, 1.8045
d+08, 1.8046
d+08, 1.8047
d+08,&
156 1.8049
d+08, 1.8051
d+08, 1.8053
d+08, 1.8055
d+08, 1.8059
d+08,&
157 1.8063
d+08, 1.8068
d+08, 1.8074
d+08, 1.8082
d+08, 1.8092
d+08,&
158 1.8103
d+08, 1.8118
d+08, 1.8135
d+08, 1.8157
d+08, 1.8182
d+08,&
159 1.8214
d+08, 1.8252
d+08, 1.8299
d+08, 1.8357
d+08, 1.8428
d+08,&
160 1.8515
d+08, 1.8621
d+08, 1.8749
d+08, 1.8904
d+08, 1.9089
d+08,&
161 1.9307
d+08, 1.9561
d+08, 1.9852
d+08, 2.0180
d+08, 2.0544
d+08,&
162 2.0941
d+08, 2.1368
d+08, 2.1819
d+08, 2.2291
d+08, 2.2781
d+08,&
163 2.3284
d+08, 2.3798
d+08, 2.4320
d+08, 2.4849
d+08, 2.5383
d+08,&
164 2.5922
d+08, 2.6463
d+08, 2.7008
d+08, 2.7554
d+08, 2.8102
d+08,&
165 2.8651
d+08, 2.9202
d+08, 2.9753
d+08, 3.0305
d+08, 3.0858
d+08,&
166 3.1412
d+08, 3.1966
d+08, 3.2520
d+08, 3.3075
d+08, 3.3630
d+08,&
167 3.4186
d+08, 3.4741
d+08, 3.5297
d+08, 3.5854
d+08, 3.6410
d+08,&
168 3.6967
d+08, 3.7523
d+08, 3.8080
d+08, 3.8637
d+08, 3.9195
d+08,&
169 3.9752
d+08, 4.0309
d+08, 4.0867
d+08, 4.1424
d+08, 4.1982
d+08,&
170 4.2540
d+08, 4.3098
d+08, 4.3656
d+08, 4.4214
d+08, 4.4772
d+08,&
171 4.5330
d+08, 4.5888
d+08, 4.6446
d+08, 4.7004
d+08, 4.7563
d+08,&
172 4.8121
d+08, 4.8679
d+08, 4.9238
d+08, 4.9796
d+08, 5.0355
d+08,&
173 5.0913
d+08, 5.1472
d+08, 5.2030
d+08, 5.2589
d+08, 5.3148
d+08,&
174 5.3706
d+08, 5.4265
d+08, 5.4824
d+08, 5.5382
d+08, 5.5941
d+08,&
175 5.6500
d+08, 5.7059
d+08, 5.7617
d+08, 5.8176
d+08, 5.8735
d+08,&
176 5.9294
d+08, 5.9853
d+08, 6.0412
d+08, 6.0971
d+08, 6.1530
d+08,&
177 6.2089
d+08, 6.2648
d+08, 6.3207
d+08, 6.3766
d+08, 6.4325
d+08,&
178 6.4884
d+08, 6.5443
d+08, 6.6002
d+08, 6.6561
d+08, 6.7120
d+08,&
179 6.7679
d+08, 6.8238
d+08, 6.8798
d+08, 6.9357
d+08, 6.9916
d+08,&
180 7.0475
d+08, 7.1034
d+08, 7.1593
d+08, 7.2152
d+08, 7.2712
d+08,&
181 7.3271
d+08, 7.3830
d+08, 7.4389
d+08, 7.4948
d+08, 7.5508
d+08,&
182 7.6067
d+08, 7.6626
d+08, 7.7185
d+08, 7.7744
d+08, 7.8304
d+08,&
183 7.8863
d+08, 7.9422
d+08, 7.9981
d+08, 8.0540
d+08, 8.1100
d+08,&
184 8.1659
d+08, 8.2218
d+08, 8.2777
d+08, 8.3337
d+08, 8.3896
d+08,&
185 8.4455
d+08, 8.5014
d+08, 8.5574
d+08, 8.6133
d+08, 8.6692
d+08,&
186 8.7252
d+08, 8.7811
d+08, 8.8370
d+08, 8.8929
d+08, 8.9489
d+08,&
187 9.0048
d+08, 9.0607
d+08, 9.1166
d+08, 9.1726
d+08, 9.2285
d+08,&
188 9.2844
d+08, 9.3404
d+08, 9.3963
d+08, 9.4522
d+08, 9.5082
d+08,&
189 9.5641
d+08, 9.6200
d+08, 9.6759
d+08, 9.7319
d+08, 9.7878
d+08,&
190 9.8437
d+08, 9.8997
d+08, 9.9556
d+08, 1.0012
d+09, 1.0067
d+09 /
192 data t_hong / 1.0000
d+04, 1.0000
d+04, 1.0001
d+04, 1.0002
d+04, 1.0003
d+04,&
193 1.0005
d+04, 1.0008
d+04, 1.0012
d+04, 1.0018
d+04, 1.0026
d+04,&
194 1.0036
d+04, 1.0047
d+04, 1.0058
d+04, 1.0011
d+04, 9.6700
d+03,&
195 9.3468
d+03, 8.9643
d+03, 8.6116
d+03, 8.2357
d+03, 7.8466
d+03,&
196 7.4734
d+03, 7.1441
d+03, 6.8397
d+03, 6.5263
d+03, 6.2978
d+03,&
197 6.1165
d+03, 5.9274
d+03, 5.7787
d+03, 5.6233
d+03, 5.4920
d+03,&
198 5.3684
d+03, 5.2602
d+03, 5.1621
d+03, 5.0768
d+03, 4.9988
d+03,&
199 4.9202
d+03, 4.8435
d+03, 4.7743
d+03, 4.7107
d+03, 4.6492
d+03,&
200 4.5903
d+03, 4.5333
d+03, 4.4783
d+03, 4.4272
d+03, 4.3832
d+03,&
201 4.3391
d+03, 4.2948
d+03, 4.2540
d+03, 4.2277
d+03, 4.2067
d+03,&
202 4.1940
d+03, 4.1851
d+03, 4.1956
d+03, 4.2262
d+03, 4.2723
d+03,&
203 4.3381
d+03, 4.4204
d+03, 4.5278
d+03, 4.6445
d+03, 4.7643
d+03,&
204 4.8856
d+03, 5.0037
d+03, 5.1176
d+03, 5.2250
d+03, 5.3211
d+03,&
205 5.4063
d+03, 5.4888
d+03, 5.5778
d+03, 5.6443
d+03, 5.7069
d+03,&
206 5.7603
d+03, 5.8132
d+03, 5.8700
d+03, 5.9167
d+03, 5.9592
d+03,&
207 5.9938
d+03, 6.0256
d+03, 6.0561
d+03, 6.0805
d+03, 6.1051
d+03,&
208 6.1298
d+03, 6.1515
d+03, 6.1720
d+03, 6.1904
d+03, 6.2096
d+03,&
209 6.2281
d+03, 6.2458
d+03, 6.2654
d+03, 6.2877
d+03, 6.3134
d+03,&
210 6.3429
d+03, 6.3754
d+03, 6.4108
d+03, 6.4497
d+03, 6.4949
d+03,&
211 6.5456
d+03, 6.6066
d+03, 6.6837
d+03, 6.7794
d+03, 6.8983
d+03,&
212 7.0383
d+03, 7.1943
d+03, 7.3537
d+03, 7.4943
d+03, 7.6133
d+03,&
213 7.7131
d+03, 7.8049
d+03, 7.8909
d+03, 7.9677
d+03, 8.0333
d+03,&
214 8.0851
d+03, 8.1219
d+03, 8.1534
d+03, 8.2225
d+03, 8.4087
d+03,&
215 8.7744
d+03, 9.3187
d+03, 1.0097
d+04, 1.1128
d+04, 1.2398
d+04,&
216 1.3941
d+04, 1.5857
d+04, 1.8286
d+04, 2.1320
d+04, 2.4676
d+04,&
217 2.7878
d+04, 3.0947
d+04, 3.4019
d+04, 3.7182
d+04, 4.0498
d+04,&
218 4.4021
d+04, 4.7801
d+04, 5.1890
d+04, 5.6344
d+04, 6.1230
d+04,&
219 6.6619
d+04, 7.2587
d+04, 7.9209
d+04, 8.6564
d+04, 9.4709
d+04,&
220 1.0364
d+05, 1.1329
d+05, 1.2363
d+05, 1.3469
d+05, 1.4653
d+05,&
221 1.5920
d+05, 1.7275
d+05, 1.8719
d+05, 2.0250
d+05, 2.1859
d+05,&
222 2.3535
d+05, 2.5260
d+05, 2.7013
d+05, 2.8767
d+05, 3.0498
d+05,&
223 3.2184
d+05, 3.3811
d+05, 3.5372
d+05, 3.6862
d+05, 3.8283
d+05,&
224 3.9636
d+05, 4.0926
d+05, 4.2157
d+05, 4.3333
d+05, 4.4459
d+05,&
225 4.5538
d+05, 4.6575
d+05, 4.7573
d+05, 4.8535
d+05, 4.9464
d+05,&
226 5.0362
d+05, 5.1232
d+05, 5.2076
d+05, 5.2895
d+05, 5.3691
d+05,&
227 5.4467
d+05, 5.5222
d+05, 5.5958
d+05, 5.6677
d+05, 5.7379
d+05,&
228 5.8066
d+05, 5.8738
d+05, 5.9395
d+05, 6.0040
d+05, 6.0672
d+05,&
229 6.1291
d+05, 6.1900
d+05, 6.2497
d+05, 6.3084
d+05, 6.3661
d+05,&
230 6.4228
d+05, 6.4786
d+05, 6.5335
d+05, 6.5876
d+05, 6.6408
d+05,&
231 6.6933
d+05, 6.7450
d+05, 6.7960
d+05, 6.8463
d+05, 6.8959
d+05,&
232 6.9449
d+05, 6.9932
d+05, 7.0409
d+05, 7.0880
d+05, 7.1346
d+05,&
233 7.1806
d+05, 7.2260
d+05, 7.2709
d+05, 7.3154
d+05, 7.3593
d+05,&
234 7.4027
d+05, 7.4457
d+05, 7.4882
d+05, 7.5303
d+05, 7.5719
d+05,&
235 7.6131
d+05, 7.6539
d+05, 7.6943
d+05, 7.7343
d+05, 7.7739
d+05,&
236 7.8131
d+05, 7.8518
d+05, 7.8898
d+05, 7.9272
d+05, 7.9641
d+05,&
237 8.0004
d+05, 8.0362
d+05, 8.0715
d+05, 8.1063
d+05, 8.1406
d+05,&
238 8.1745
d+05, 8.2079
d+05, 8.2409
d+05, 8.2734
d+05, 8.3056
d+05,&
239 8.3373
d+05, 8.3687
d+05, 8.3997
d+05, 8.4303
d+05, 8.4606
d+05,&
240 8.4906
d+05, 8.5205
d+05, 8.5502
d+05, 8.5796
d+05, 8.6090
d+05,&
241 8.6381
d+05, 8.6670
d+05, 8.6958
d+05, 8.7244
d+05, 8.7528
d+05,&
242 8.7810
d+05, 8.8091
d+05, 8.8371
d+05, 8.8648
d+05, 8.8924
d+05,&
243 8.9199
d+05, 8.9472
d+05, 8.9743
d+05, 9.0013
d+05, 9.0281
d+05,&
244 9.0548
d+05, 9.0814
d+05, 9.1078
d+05, 9.1341
d+05, 9.1602
d+05,&
245 9.1862
d+05, 9.2121
d+05, 9.2378
d+05, 9.2634
d+05, 9.2889
d+05,&
246 9.3142
d+05, 9.3394
d+05, 9.3645
d+05, 9.3895
d+05, 9.4143
d+05,&
247 9.4391
d+05, 9.4637
d+05, 9.4882
d+05, 9.5126
d+05, 9.5368
d+05,&
248 9.5610
d+05, 9.5850
d+05, 9.6089
d+05, 9.6328
d+05, 9.6565
d+05,&
249 9.6801
d+05, 9.7036
d+05, 9.7270
d+05, 9.7503
d+05, 9.7735
d+05,&
250 9.7965
d+05, 9.8195
d+05, 9.8424
d+05, 9.8652
d+05, 9.8879
d+05,&
251 9.9105
d+05, 9.9330
d+05, 9.9554
d+05, 9.9778
d+05, 9.9778
d+05 /
264 integer,
intent(in) :: nh
265 double precision,
intent(in) :: h(nh),grav(nh)
266 double precision,
intent(out) :: rho(nh),pth(nh)
267 double precision,
intent(out),
optional :: Tem(nh)
268 double precision :: rhohc,hc
269 character(*) :: Tcurve
271 double precision :: h_cgs(nh),Te_cgs(nh),Te(nh)
272 double precision :: invT,dh,rhot,dht,ratio
280 if (
mype==0) print *,
'Temperature curve from Vernazza, Avrett & Loeser, 1981, ApJS, 45, 635'
284 if (
mype==0) print *,
'Temperature curve from Hong et al. 2017, ApJ, 845, 144'
288 if (
mype==0) print *,
'Temperature curve from Fontenla et al. 2007, ApJ, 667, 1243'
292 if (
mype==0) print *,
'Temperature curve from Avrett & Loeser 2008, ApJS, 175, 229'
295 call mpistop(
"Unknown temperature curve")
300 if(
present(tem)) tem=te
316 invt=invt+dh*(grav(j)/te(j)+grav(j-1)/te(j-1))*0.5d0
317 pth(j)=pth(1)*dexp(invt)
320 if (h(j-1)<=hc .and. h(j)>hc)
then
322 rhot=rho(j-1)+dht*(rho(j)-rho(j-1))/(h(j)-h(j-1))
336 double precision :: h(nh),Te(nh)
338 integer :: ih,j,imin,imax,n_table
339 double precision :: h_table(n_alc7),T_table(n_alc7)
340 double precision :: unit_h,unit_T,htra,Ttra,Fc,kappa,dTdh
357 if (h(ih)>h_table(n_table))
then
359 te(ih)=(3.5d0*fc*(h(ih)-htra)/kappa+ttra**3.5d0)**(2.d0/7.d0)
362 if (h(ih)<=h_table(1))
then
364 te(ih)=(h(ih)-h_table(1))*dtdh+t_table(1)
371 if (h(1)>=h_table(1) .and. h(1)<=h_table(n_table))
then
375 if (h(ih-1)<h_table(1) .and. h(ih)>=h_table(1) .and. h(ih)<=h_table(n_table)) imin=ih
378 if (h(nh)>=h_table(1) .and. h(nh)<=h_table(n_table))
then
382 if (h(ih)<=h_table(n_table) .and. h(ih+1)>h_table(n_table) .and. h(ih)>=h_table(1)) imax=ih
387 call interp_linear(h_table,t_table,n_table,h(imin:imax),te(imin:imax),imax-imin+1)
396 double precision :: h(nh),Te(nh)
398 double precision :: h_table(n_fontenla),T_table(n_fontenla)
399 double precision :: unit_h,unit_T,dTdh
400 double precision :: h1,h2,h3
401 double precision :: Tpho,Ttop,htanh,wtra
402 double precision :: htra,Ttra,Fc,kappa
403 integer :: ih,j,imin,imax,n_table
430 dtdh=(t_table(2)-t_table(1))/(h_table(2)-h_table(1))
436 te(ih)=(h(ih)-h_table(1))*dtdh+t_table(1)
440 if (h(ih)>h2 .and. h(ih)<h3)
then
441 te(ih)=tpho+0.5d0*(ttop-tpho)*(tanh((h(ih)-htanh)/wtra)+1.d0)
446 te(ih)=(3.5d0*fc*(h(ih)-htra)/kappa+ttra**3.5)**(2.d0/7.d0)
454 if (h(1)>=h_table(1) .and. h(1)<=h_table(n_table))
then
458 if (h(ih-1)<h_table(1) .and. h(ih)>=h_table(1) .and. h(ih)<=h_table(n_table)) imin=ih
461 if (h(nh)>=h_table(1) .and. h(nh)<=h_table(n_table))
then
465 if (h(ih)<=h_table(n_table) .and. h(ih+1)>h_table(n_table) .and. h(ih)>=h_table(1)) imax=ih
470 call interp_linear(h_table,t_table,n_fontenla,h(imin:imax),te(imin:imax),imax-imin+1)
479 double precision :: h(nh),Te(nh)
481 double precision :: h_table(n_valc),T_table(n_valc)
482 double precision :: unit_h,unit_T,htra,Ttra,Fc,kappa,dTdh
483 integer :: ih,j,imin,imax,n_table
494 htra=h_table(n_table)
495 ttra=t_table(n_table)
496 dtdh=(t_table(2)-t_table(1))/(h_table(2)-h_table(1))
499 if (h(ih)>h_table(n_table))
then
501 te(ih)=(3.5d0*fc*(h(ih)-htra)/kappa+ttra**3.5d0)**(2.d0/7.d0)
504 if (h(ih)<=h_table(1))
then
506 te(ih)=(h(ih)-h_table(1))*dtdh+t_table(1)
514 if (h(1)>=h_table(1) .and. h(1)<=h_table(n_table))
then
518 if (h(ih-1)<h_table(1) .and. h(ih)>=h_table(1) .and. h(ih)<=h_table(n_table)) imin=ih
521 if (h(nh)>=h_table(1) .and. h(nh)<=h_table(n_table))
then
525 if (h(ih)<=h_table(n_table) .and. h(ih+1)>h_table(n_table) .and. h(ih)>=h_table(1)) imax=ih
530 call interp_linear(h_table,t_table,n_valc,h(imin:imax),te(imin:imax),imax-imin+1)
539 double precision :: h(nh),Te(nh)
541 double precision :: h_table(n_hong),T_table(n_hong)
542 double precision :: dTdh
543 double precision :: htra,Ttra,Fc,kappa
544 integer :: ih,j,imin,imax,n_table
556 if (h(ih)>=h_table(n_table))
then
558 te(ih)=(3.5d0*fc*(h(ih)-htra)/kappa+ttra**3.5d0)**(2.d0/7.d0)
561 if (h(ih)<=h_table(1))
then
563 dtdh=(t_table(2)-t_table(1))/(h_table(2)-h_table(1))
564 te(ih)=(h(ih)-h_table(1))*dtdh+t_table(1)
572 if (h(1)>=h_table(1) .and. h(1)<=h_table(n_table))
then
576 if (h(ih-1)<h_table(1) .and. h(ih)>=h_table(1) .and. h(ih)<=h_table(n_table)) imin=ih
579 if (h(nh)>=h_table(1) .and. h(nh)<=h_table(n_table))
then
583 if (h(ih)<=h_table(n_table) .and. h(ih+1)>h_table(n_table) .and. h(ih)>=h_table(1)) imax=ih
588 call interp_linear(h_table,t_table,n_table,h(imin:imax),te(imin:imax),imax-imin+1)
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
This module contains definitions of global parameters and variables and some generic functions/subrou...
double precision unit_length
Physical scaling factor for length.
integer mype
The rank of the current MPI task.
double precision, dimension(:), allocatable, parameter d
double precision unit_temperature
Physical scaling factor for temperature.
subroutine interp_linear(x_table, y_table, n_table, x_itp, y_itp, n_itp)
module ionization degree - get ionization degree for given temperature
double precision function r_ideal_gas_law_partial_ionization(Te)
This module defines the procedures of a physics module. It contains function pointers for the various...
logical phys_partial_ionization
Solve partially ionized one-fluid plasma.
User can use subroutine get_atm_para to generate 1D solar stmosphere. User should provide heights (h)...
double precision, dimension(1:120) t_alc7
subroutine get_te_alc7(h, Te, nh)
subroutine get_te_hong(h, Te, nh)
double precision, dimension(1:300) t_hong
double precision, dimension(1:56) t_fontenla
double precision, dimension(1:52) t_valc
subroutine get_atm_para(h, rho, pth, grav, nh, Tcurve, hc, rhohc, Tem)
subroutine get_te_valc(h, Te, nh)
double precision, dimension(1:52) h_valc
double precision, dimension(1:300) h_hong
subroutine get_te_sprm(h, Te, nh)
double precision, dimension(1:56) h_fontenla
double precision, dimension(1:120) h_alc7
double precision, parameter solar_gravity