SCALE-RM
scale_atmos_solarins.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
12 !-------------------------------------------------------------------------------
13 #include "scalelib.h"
15  !-----------------------------------------------------------------------------
16  !
17  !++ used modules
18  !
19  use scale_precision
20  use scale_io
21  use scale_prof
22  !-----------------------------------------------------------------------------
23  implicit none
24  private
25  !-----------------------------------------------------------------------------
26  !
27  !++ Public procedure
28  !
29  public :: atmos_solarins_setup
30  public :: atmos_solarins_orbit
31  public :: atmos_solarins_insolation
32 
33  interface atmos_solarins_insolation
34  module procedure atmos_solarins_insolation_0d
35  module procedure atmos_solarins_insolation_2d
36  end interface atmos_solarins_insolation
37 
38  !-----------------------------------------------------------------------------
39  !
40  !++ Public parameters & variables
41  !
42  real(RP), public :: atmos_solarins_constant = 1360.250117_rp ! Solar constant [W/m2]
43 
44  logical, public :: atmos_solarins_set_ve = .false. ! Set vernal equinox condition?
45 
46  logical, public :: atmos_solarins_set_ideal = .false. ! Set obliquity and eccentricity?
47  real(RP), public :: atmos_solarins_obliquity = 0.0_rp ! Obliquity [deg]
48  real(RP), public :: atmos_solarins_eccentricity = 0.0_rp ! Eccentricity
49 
50  logical, public :: atmos_solarins_fixedlatlon = .false. ! Latitude/Longitude is fixed?
51  real(RP), public :: atmos_solarins_lon ! Longitude for radiation [rad]
52  real(RP), public :: atmos_solarins_lat ! Latitude for radiation [rad]
53 
54  logical, public :: atmos_solarins_fixeddate = .false. ! Date is fixed?
55  integer, public :: atmos_solarins_date(6) ! Date for radiation [Y,M,D,H,M,S]
56 
57  !-----------------------------------------------------------------------------
58  !
59  !++ Private procedure
60  !
61  !-----------------------------------------------------------------------------
62  !
63  !++ Private parameters & variables
64  !
65  real(RP), private :: obliquity ! obliquity [rad]
66  real(RP), private :: e ! eccentricity
67  real(RP), private :: omega ! longitude of perigee [rad]
68  real(RP), private :: lambda_m0 ! longitude at the vernal equinox [rad]
69 
70  integer, private, parameter :: year_ref = 1950 ! reference year [year]
71  integer, private :: ve_date(6) ! reference date of vernal equinox
72  data ve_date / 1950, 3, 21, 0, 0, 0 /
73 
74  real(RP), private, parameter :: obliquity_ref = 23.320556_rp ! initial condition of obliquity (epsilon_star)
75  real(RP), private, parameter :: psi_bar = 50.439273_rp ! parameter for general precession [second of arc]
76  real(RP), private, parameter :: zeta = 3.392506_rp ! parameter for general precession [degree]
77 
78  !-----< Parameter tables from Berger(1978b) >-----
79  integer, private, parameter :: nobliq = 47 ! # of terms of the series expansion of epsilon
80  real(RP), private :: obliq_amp (nobliq) ! amplitude [second of arc]
81  real(RP), private :: obliq_rate (nobliq) ! mean rate [second of arc/year]
82  real(RP), private :: obliq_phase(nobliq) ! phase [degree of arc]
83 
84  data obliq_amp / &
85  -2462.2214466_rp, & ! 1
86  -857.3232075_rp, & ! 2
87  -629.3231835_rp, & ! 3
88  -414.2804924_rp, & ! 4
89  -311.7632587_rp, & ! 5
90  308.9408604_rp, & ! 6
91  -162.5533601_rp, & ! 7
92  -116.1077911_rp, & ! 8
93  101.1189923_rp, & ! 9
94  -67.6856209_rp, & ! 10
95  24.9079067_rp, & ! 11
96  22.5811241_rp, & ! 12
97  -21.1648355_rp, & ! 13
98  -15.6549876_rp, & ! 14
99  15.3936813_rp, & ! 15
100  14.6660938_rp, & ! 16
101  -11.7273029_rp, & ! 17
102  10.2742696_rp, & ! 18
103  6.4914588_rp, & ! 19
104  5.8539148_rp, & ! 20
105  -5.4872205_rp, & ! 21
106  -5.4290191_rp, & ! 22
107  5.1609570_rp, & ! 23
108  5.0786314_rp, & ! 24
109  -4.0735782_rp, & ! 25
110  3.7227167_rp, & ! 26
111  3.3971932_rp, & ! 27
112  -2.8347004_rp, & ! 28
113  -2.6550721_rp, & ! 29
114  -2.5717867_rp, & ! 30
115  -2.4712188_rp, & ! 31
116  2.4625410_rp, & ! 32
117  2.2464112_rp, & ! 33
118  -2.0755511_rp, & ! 34
119  -1.9713669_rp, & ! 35
120  -1.8813061_rp, & ! 36
121  -1.8468785_rp, & ! 37
122  1.8186742_rp, & ! 38
123  1.7601888_rp, & ! 39
124  -1.5428851_rp, & ! 40
125  1.4738838_rp, & ! 41
126  -1.4593669_rp, & ! 42
127  1.4192259_rp, & ! 43
128  -1.1818980_rp, & ! 44
129  1.1756474_rp, & ! 45
130  -1.1316126_rp, & ! 46
131  1.0896928_rp / ! 47
132 
133  data obliq_rate / &
134  31.609974_rp, & ! 1
135  32.620504_rp, & ! 2
136  24.172203_rp, & ! 3
137  31.983787_rp, & ! 4
138  44.828336_rp, & ! 5
139  30.973257_rp, & ! 6
140  43.668246_rp, & ! 7
141  32.246691_rp, & ! 8
142  30.599444_rp, & ! 9
143  42.681324_rp, & ! 10
144  43.836462_rp, & ! 11!
145  47.439436_rp, & ! 12
146  63.219948_rp, & ! 13
147  64.230478_rp, & ! 14
148  1.010530_rp, & ! 15
149  7.437771_rp, & ! 16
150  55.782177_rp, & ! 17
151  0.373813_rp, & ! 18
152  13.218362_rp, & ! 19
153  62.583231_rp, & ! 20
154  63.593761_rp, & ! 21
155  76.438310_rp, & ! 22
156  45.815258_rp, & ! 23
157  8.448301_rp, & ! 24
158  56.792707_rp, & ! 25
159  49.747842_rp, & ! 26
160  12.058272_rp, & ! 27
161  75.278220_rp, & ! 28
162  65.241008_rp, & ! 29
163  64.604291_rp, & ! 30
164  1.647247_rp, & ! 31
165  7.811584_rp, & ! 32
166  12.207832_rp, & ! 33
167  63.856665_rp, & ! 34
168  56.155990_rp, & ! 35
169  77.448840_rp, & ! 36
170  6.801054_rp, & ! 37
171  62.209418_rp, & ! 38
172  20.656133_rp, & ! 39
173  48.344406_rp, & ! 40
174  55.145460_rp, & ! 41
175  69.000539_rp, & ! 42
176  11.071350_rp, & ! 43
177  74.291298_rp, & ! 44
178  11.047742_rp, & ! 45
179  0.636717_rp, & ! 46
180  12.844549_rp / ! 47
181 
182  data obliq_phase / &
183  251.9025_rp, & ! 1
184  280.8325_rp, & ! 2
185  128.3057_rp, & ! 3
186  292.7252_rp, & ! 4
187  15.3747_rp, & ! 5
188  263.7951_rp, & ! 6
189  308.4258_rp, & ! 7
190  240.0099_rp, & ! 8
191  222.9725_rp, & ! 9
192  268.7809_rp, & ! 10
193  316.7998_rp, & ! 11
194  319.6024_rp, & ! 12
195  143.8050_rp, & ! 13
196  172.7351_rp, & ! 14
197  28.9300_rp, & ! 15
198  123.5968_rp, & ! 16
199  20.2082_rp, & ! 17
200  40.8226_rp, & ! 18
201  123.4722_rp, & ! 19
202  155.6977_rp, & ! 20
203  184.6277_rp, & ! 21
204  267.2772_rp, & ! 22
205  55.0196_rp, & ! 23
206  152.5268_rp, & ! 24
207  49.1382_rp, & ! 25
208  204.6609_rp, & ! 26
209  56.5233_rp, & ! 27
210  200.3284_rp, & ! 28
211  201.6651_rp, & ! 29
212  213.5577_rp, & ! 30
213  17.0374_rp, & ! 31
214  164.4194_rp, & ! 32
215  94.5422_rp, & ! 33
216  131.9124_rp, & ! 34
217  61.0309_rp, & ! 35
218  296.2073_rp, & ! 36
219  135.4894_rp, & ! 37
220  114.8750_rp, & ! 38
221  247.0691_rp, & ! 39
222  256.6114_rp, & ! 40
223  32.1008_rp, & ! 41
224  143.6804_rp, & ! 42
225  16.8784_rp, & ! 43
226  160.6835_rp, & ! 44
227  27.5932_rp, & ! 45
228  348.1074_rp, & ! 46
229  82.6496_rp / ! 47
230 
231  integer, private, parameter :: neclip = 19 ! # of terms of the series expansion of ecliptic
232  real(RP), private :: eclip_amp (neclip) ! amplitude
233  real(RP), private :: eclip_rate (neclip) ! mean rate [second of arc/year]
234  real(RP), private :: eclip_phase(neclip) ! phase [degree of arc]
235 
236  data eclip_amp / &
237  0.01860798_rp, & ! 1
238  0.01627522_rp, & ! 2
239  -0.01300660_rp, & ! 3
240  0.00988829_rp, & ! 4
241  -0.00336700_rp, & ! 5
242  0.00333077_rp, & ! 6
243  -0.00235400_rp, & ! 7
244  0.00140015_rp, & ! 8
245  0.00100700_rp, & ! 9
246  0.00085700_rp, & ! 10
247  0.00064990_rp, & ! 11
248  0.00059900_rp, & ! 12
249  0.00037800_rp, & ! 13
250  -0.00033700_rp, & ! 14
251  0.00027600_rp, & ! 15
252  0.00018200_rp, & ! 16
253  -0.00017400_rp, & ! 17
254  -0.00012400_rp, & ! 18
255  0.00001250_rp / ! 19
256 
257  data eclip_rate / &
258  4.2072050_rp, & ! 1
259  7.3460910_rp, & ! 2
260  17.8572630_rp, & ! 3
261  17.2205460_rp, & ! 4
262  16.8467330_rp, & ! 5
263  5.1990790_rp, & ! 6
264  18.2310760_rp, & ! 7
265  26.2167580_rp, & ! 8
266  6.3591690_rp, & ! 9
267  16.2100160_rp, & ! 10
268  3.0651810_rp, & ! 11
269  16.5838290_rp, & ! 12
270  18.4939800_rp, & ! 13
271  6.1909530_rp, & ! 14
272  18.8677930_rp, & ! 15
273  17.4255670_rp, & ! 16
274  6.1860010_rp, & ! 17
275  18.4174410_rp, & ! 18
276  0.6678630_rp / ! 19
277 
278  data eclip_phase / &
279  28.620089_rp, & ! 1
280  193.788772_rp, & ! 2
281  308.307024_rp, & ! 3
282  320.199637_rp, & ! 4
283  279.376984_rp, & ! 5
284  87.195000_rp, & ! 6
285  349.129677_rp, & ! 7
286  128.443387_rp, & ! 8
287  154.143880_rp, & ! 9
288  291.269597_rp, & ! 10
289  114.860583_rp, & ! 11
290  332.092251_rp, & ! 12
291  296.414411_rp, & ! 13
292  145.769910_rp, & ! 14
293  337.237063_rp, & ! 15
294  152.092288_rp, & ! 16
295  126.839891_rp, & ! 17
296  210.667199_rp, & ! 18
297  72.108838_rp / ! 19
298 
299  integer, private, parameter :: nprece = 78 ! # of terms of the series expansion of general precession
300  real(RP), private :: prece_amp (nprece) ! amplitude
301  real(RP), private :: prece_rate (nprece) ! mean rate [second of arc/year]
302  real(RP), private :: prece_phase(nprece) ! phase [degree of arc]
303 
304  data prece_amp / &
305  7391.0225890_rp, & ! 1
306  2555.1526947_rp, & ! 2
307  2022.7629188_rp, & ! 3
308  -1973.6517951_rp, & ! 4
309  1240.2321818_rp, & ! 5
310  953.8679112_rp, & ! 6
311  -931.7537108_rp, & ! 7
312  872.3795383_rp, & ! 8
313  606.3544732_rp, & ! 9
314  -496.0274038_rp, & ! 10
315  456.9608039_rp, & ! 11
316  346.9462320_rp, & ! 12
317  -305.8412902_rp, & ! 13
318  249.6173246_rp, & ! 14
319  -199.1027200_rp, & ! 15
320  191.0560889_rp, & ! 16
321  -175.2936572_rp, & ! 17
322  165.9068833_rp, & ! 18
323  161.1285917_rp, & ! 19
324  139.7878093_rp, & ! 20
325  -133.5228399_rp, & ! 21
326  117.0673811_rp, & ! 22
327  104.6907281_rp, & ! 23
328  95.3227476_rp, & ! 24
329  86.7824524_rp, & ! 25
330  86.0857729_rp, & ! 26
331  70.5893698_rp, & ! 27
332  -69.9719343_rp, & ! 28
333  -62.5817473_rp, & ! 29
334  61.5450059_rp, & ! 30
335  -57.9364011_rp, & ! 31
336  57.1899832_rp, & ! 32
337  -57.0236109_rp, & ! 33
338  -54.2119253_rp, & ! 34
339  53.2834147_rp, & ! 35
340  52.1223575_rp, & ! 36
341  -49.0059908_rp, & ! 37
342  -48.3118757_rp, & ! 38
343  -45.4191685_rp, & ! 39
344  -42.2357920_rp, & ! 40
345  -34.7971099_rp, & ! 41
346  34.4623613_rp, & ! 42
347  -33.8356643_rp, & ! 43
348  33.6689362_rp, & ! 44
349  -31.2521586_rp, & ! 45
350  -30.8798701_rp, & ! 46
351  28.4640769_rp, & ! 47
352  -27.1960802_rp, & ! 48
353  27.0860736_rp, & ! 49
354  -26.3437456_rp, & ! 50
355  24.7253740_rp, & ! 51
356  24.6732126_rp, & ! 52
357  24.4272733_rp, & ! 53
358  24.0127327_rp, & ! 54
359  21.7150294_rp, & ! 55
360  -21.5375347_rp, & ! 56
361  18.1148363_rp, & ! 57
362  -16.9603104_rp, & ! 58
363  -16.1765215_rp, & ! 59
364  15.5567653_rp, & ! 60
365  15.4846529_rp, & ! 61
366  15.2150632_rp, & ! 62
367  14.5047426_rp, & ! 63
368  -14.3873316_rp, & ! 64
369  13.1351419_rp, & ! 65
370  12.8776311_rp, & ! 66
371  11.9867234_rp, & ! 67
372  11.9385578_rp, & ! 68
373  11.7030822_rp, & ! 69
374  11.6018181_rp, & ! 70
375  -11.2617293_rp, & ! 71
376  -10.4664199_rp, & ! 72
377  10.4333970_rp, & ! 73
378  -10.2377466_rp, & ! 74
379  10.1934446_rp, & ! 75
380  -10.1280191_rp, & ! 76
381  10.0289441_rp, & ! 77
382  -10.0034259_rp / ! 78
383 
384  data prece_rate / &
385  31.609974_rp, & ! 1
386  32.620504_rp, & ! 2
387  24.172203_rp, & ! 3
388  0.636717_rp, & ! 4
389  31.983787_rp, & ! 5
390  3.138886_rp, & ! 6
391  30.973257_rp, & ! 7
392  44.828336_rp, & ! 8
393  0.991874_rp, & ! 9
394  0.373813_rp, & ! 10
395  43.668246_rp, & ! 11
396  32.246691_rp, & ! 12
397  30.599444_rp, & ! 13
398  2.147012_rp, & ! 14
399  10.511172_rp, & ! 15
400  42.681324_rp, & ! 16
401  13.650058_rp, & ! 17
402  0.986922_rp, & ! 18
403  9.874455_rp, & ! 19
404  13.013341_rp, & ! 20
405  0.262904_rp, & ! 21
406  0.004952_rp, & ! 22
407  1.142024_rp, & ! 23
408  63.219948_rp, & ! 24
409  0.205021_rp, & ! 25
410  2.151964_rp, & ! 26
411  64.230478_rp, & ! 27
412  43.836462_rp, & ! 28
413  47.439436_rp, & ! 29
414  1.384343_rp, & ! 30
415  7.437771_rp, & ! 31
416  18.829299_rp, & ! 32
417  9.500642_rp, & ! 33
418  0.431696_rp, & ! 34
419  1.160090_rp, & ! 35
420  55.782177_rp, & ! 36
421  12.639528_rp, & ! 37
422  1.155138_rp, & ! 38
423  0.168216_rp, & ! 39
424  1.647247_rp, & ! 40
425  10.884985_rp, & ! 41
426  5.610937_rp, & ! 42
427  12.658184_rp, & ! 43
428  1.010530_rp, & ! 44
429  1.983748_rp, & ! 45
430  14.023871_rp, & ! 46
431  0.560178_rp, & ! 47
432  1.273434_rp, & ! 48
433  12.021467_rp, & ! 49
434  62.583231_rp, & ! 50
435  63.593761_rp, & ! 51
436  76.438310_rp, & ! 52
437  4.280910_rp, & ! 53
438  13.218362_rp, & ! 54
439  17.818769_rp, & ! 55
440  8.359495_rp, & ! 56
441  56.792707_rp, & ! 57
442  8.448301_rp, & ! 58
443  1.978796_rp, & ! 59
444  8.863925_rp, & ! 60
445  0.186365_rp, & ! 61
446  8.996212_rp, & ! 62
447  6.771027_rp, & ! 63
448  45.815258_rp, & ! 64
449  12.002811_rp, & ! 65
450  75.278220_rp, & ! 66
451  65.241008_rp, & ! 67
452  18.870667_rp, & ! 68
453  22.009553_rp, & ! 69
454  64.604291_rp, & ! 70
455  11.498094_rp, & ! 71
456  0.578834_rp, & ! 72
457  9.237738_rp, & ! 73
458  49.747842_rp, & ! 74
459  2.147012_rp, & ! 75
460  1.196895_rp, & ! 76
461  2.133898_rp, & ! 77
462  0.173168_rp / ! 78
463 
464  data prece_phase / &
465  251.9025_rp, & ! 1
466  280.8325_rp, & ! 2
467  128.3057_rp, & ! 3
468  348.1074_rp, & ! 4
469  292.7252_rp, & ! 5
470  165.1686_rp, & ! 6
471  263.7951_rp, & ! 7
472  15.3747_rp, & ! 8
473  58.5749_rp, & ! 9
474  40.8226_rp, & ! 10
475  308.4258_rp, & ! 11
476  240.0099_rp, & ! 12
477  222.9725_rp, & ! 13
478  106.5937_rp, & ! 14
479  114.5182_rp, & ! 15
480  268.7809_rp, & ! 16
481  279.6869_rp, & ! 17
482  39.6448_rp, & ! 18
483  126.4108_rp, & ! 19
484  291.5795_rp, & ! 20
485  307.2848_rp, & ! 21
486  18.9300_rp, & ! 22
487  273.7596_rp, & ! 23
488  143.8050_rp, & ! 24
489  191.8927_rp, & ! 25
490  125.5237_rp, & ! 26
491  172.7351_rp, & ! 27
492  316.7998_rp, & ! 28
493  319.6024_rp, & ! 29
494  69.7526_rp, & ! 30
495  123.5968_rp, & ! 31
496  217.6432_rp, & ! 32
497  85.5882_rp, & ! 33
498  156.2147_rp, & ! 34
499  66.9489_rp, & ! 35
500  20.2082_rp, & ! 36
501  250.7568_rp, & ! 37
502  48.0188_rp, & ! 38
503  8.3739_rp, & ! 39
504  17.0374_rp, & ! 40
505  155.3409_rp, & ! 41
506  94.1709_rp, & ! 42
507  221.1120_rp, & ! 43
508  28.9300_rp, & ! 44
509  117.1498_rp, & ! 45
510  320.5095_rp, & ! 46
511  262.3602_rp, & ! 47
512  336.2148_rp, & ! 48
513  233.0046_rp, & ! 49
514  155.6977_rp, & ! 50
515  184.6277_rp, & ! 51
516  267.2772_rp, & ! 52
517  78.9281_rp, & ! 53
518  123.4722_rp, & ! 54
519  188.7132_rp, & ! 55
520  180.1364_rp, & ! 56
521  49.1382_rp, & ! 57
522  152.5268_rp, & ! 58
523  98.2198_rp, & ! 59
524  97.4808_rp, & ! 60
525  221.5376_rp, & ! 61
526  168.2438_rp, & ! 62
527  161.1199_rp, & ! 63
528  55.0196_rp, & ! 64
529  262.6495_rp, & ! 65
530  200.3284_rp, & ! 66
531  201.6651_rp, & ! 67
532  294.6547_rp, & ! 68
533  99.8233_rp, & ! 69
534  213.5577_rp, & ! 70
535  154.1631_rp, & ! 71
536  232.7153_rp, & ! 72
537  138.3034_rp, & ! 73
538  204.6609_rp, & ! 74
539  106.5938_rp, & ! 75
540  250.4676_rp, & ! 76
541  332.3345_rp, & ! 77
542  27.3039_rp / ! 78
543 
544  real(RP), private :: arcsec2d, arcsec2r ! unit converter
545 
546  logical, private :: debug = .false.
547 
548  !-----------------------------------------------------------------------------
549 contains
550  !-----------------------------------------------------------------------------
552  !-----------------------------------------------------------------------------
553  subroutine atmos_solarins_setup( &
554  basepoint_lon, basepoint_lat, &
555  iyear )
556  use scale_prc, only: &
557  prc_abort
558  use scale_const, only: &
559  const_d2r
560  implicit none
561  integer, intent(in) :: iyear ! year at setup
562  real(RP), intent(in) :: basepoint_lon
563  real(RP), intent(in) :: basepoint_lat
564 
565  namelist / param_atmos_solarins / &
576  debug
577 
578  real(RP) :: dyear ! delta t [year]
579  integer :: year ! used year
580 
581  integer :: ierr
582  !---------------------------------------------------------------------------
583 
584  log_newline
585  log_info("ATMOS_SOLARINS_setup",*) 'Setup'
586 
587  atmos_solarins_lon = basepoint_lon
588  atmos_solarins_lat = basepoint_lat
589  atmos_solarins_date(:) = -1
590 
591  !--- read namelist
592  rewind(io_fid_conf)
593  read(io_fid_conf,nml=param_atmos_solarins,iostat=ierr)
594  if( ierr < 0 ) then !--- missing
595  log_info("ATMOS_SOLARINS_setup",*) 'Not found namelist. Default used.'
596  elseif( ierr > 0 ) then !--- fatal error
597  log_error("ATMOS_SOLARINS_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_SOLARINS. Check!'
598  call prc_abort
599  endif
600  log_nml(param_atmos_solarins)
601 
602  if ( atmos_solarins_set_ve ) then
603  if ( .NOT. atmos_solarins_set_ideal ) then
604  atmos_solarins_set_ideal = .true.
605  atmos_solarins_obliquity = 0.0_rp
607  endif
608  if ( .NOT. atmos_solarins_fixedlatlon ) then
610  atmos_solarins_lon = 0.0_rp
611  atmos_solarins_lat = 0.0_rp
612  endif
613  if ( .NOT. atmos_solarins_fixeddate ) then
614  atmos_solarins_fixeddate = .true.
615  atmos_solarins_date(1) = ve_date(1)
616  atmos_solarins_date(2) = ve_date(2)
617  atmos_solarins_date(3) = ve_date(3)
618  atmos_solarins_date(4) = 12
619  atmos_solarins_date(5) = 0
620  atmos_solarins_date(6) = 0
621  endif
622  endif
623 
624  arcsec2d = 1.0_rp / (60.0_rp*60.0_rp)
625  arcsec2r = arcsec2d * const_d2r
626 
627  year = iyear
628 
629  if ( atmos_solarins_fixeddate ) then
630  if( atmos_solarins_date(1) >= 0 ) year = atmos_solarins_date(1)
631  endif
632 
633  call atmos_solarins_orbit( year )
634 
635  dyear = real( year-year_ref, kind=rp )
636 
637  !----- report data -----
638  log_newline
639  log_info("ATMOS_SOLARINS_setup",'(1x,A)') 'Solar insolation parameters '
640  log_info("ATMOS_SOLARINS_setup",'(1x,A,I7)') 'Reference year : ', year_ref
641  log_info("ATMOS_SOLARINS_setup",'(1x,A,I7)') 'Current year : ', year
642  log_info("ATMOS_SOLARINS_setup",'(1x,A,I7)') 'Difference from ref. : ', int(dyear)
643  log_newline
644  log_info("ATMOS_SOLARINS_setup",'(1x,A,F12.7)') 'Solar constant [W/m2] : ', atmos_solarins_constant
645  log_info("ATMOS_SOLARINS_setup",'(1x,A,F12.7)') 'Obliquity [deg] : ', obliquity / const_d2r
646  log_info("ATMOS_SOLARINS_setup",'(1x,A,F12.7)') 'Eccentricity : ', e
647  log_info("ATMOS_SOLARINS_setup",'(1x,A,F12.7)') 'Longitude of perihelion [deg] : ', omega / const_d2r
648  log_info("ATMOS_SOLARINS_setup",'(1x,A,F12.7)') 'Longitude at the vernal equinox [deg] : ', lambda_m0 / const_d2r
649  log_newline
650 
651  log_info("ATMOS_SOLARINS_setup",*) 'Latitude/Longitude is fixed? : ', atmos_solarins_fixedlatlon
652  if ( atmos_solarins_fixedlatlon ) then
653  log_info("ATMOS_SOLARINS_setup",*) 'Longitude [deg] : ', atmos_solarins_lon
654  log_info("ATMOS_SOLARINS_setup",*) 'Latitude [deg] : ', atmos_solarins_lat
657  endif
658 
659  log_info("ATMOS_SOLARINS_setup",*) 'Date is fixed? : ', atmos_solarins_fixeddate
660  if ( atmos_solarins_fixeddate ) then
661  log_info("ATMOS_SOLARINS_setup",*) 'Date : ', atmos_solarins_date
662  endif
663 
664  return
665  end subroutine atmos_solarins_setup
666 
667  !-----------------------------------------------------------------------------
669  !-----------------------------------------------------------------------------
670  subroutine atmos_solarins_orbit( &
671  iyear )
672  use scale_const, only: &
673  pi => const_pi, &
674  d2r => const_d2r
675  implicit none
676 
677  integer, intent(in) :: iyear ! year at setup
678 
679  real(RP) :: dyear ! delta t [year]
680  real(RP) :: EsinOMG ! e * sin(w) : ecliptic parameter
681  real(RP) :: EcosOMG ! e * cos(w) : ecliptic parameter
682  real(RP) :: Perih_pi ! pi : longitude of fixed perihelion [rad]
683  real(RP) :: Perih_psi ! psi : general precession [degree]
684  real(RP) :: Perih_omega ! omega_bar : longitude of perihelion [degree]
685 
686  real(RP) :: temp, beta, EcosOMG_mod
687  integer :: i
688  !---------------------------------------------------------------------------
689 
690  ! time from reference year(1950.0 AD)
691  dyear = real( iyear-year_ref, kind=rp )
692 
693  ! obliquity
694  temp = 0.0_rp
695  do i = 1, nobliq
696  temp = temp + obliq_amp(i)*arcsec2d * cos( obliq_rate(i)*arcsec2r*dyear + obliq_phase(i)*d2r )
697  enddo
698  obliquity = ( obliquity_ref + temp ) * d2r
699 
700  ! eccentricity
701  esinomg = 0.0_rp
702  ecosomg = 0.0_rp
703  do i = 1, neclip
704  esinomg = esinomg + eclip_amp(i) * sin( eclip_rate(i)*arcsec2r*dyear + eclip_phase(i)*d2r )
705  ecosomg = ecosomg + eclip_amp(i) * cos( eclip_rate(i)*arcsec2r*dyear + eclip_phase(i)*d2r )
706  enddo
707  e = sqrt( esinomg*esinomg + ecosomg*ecosomg )
708 
709  if ( atmos_solarins_set_ideal ) then
710  obliquity = atmos_solarins_obliquity * d2r
712  endif
713 
714  ! longitude of fixed perihelion
715  if ( ecosomg == 0.0_rp ) then
716  ecosomg_mod = 1.e-30_rp
717  else
718  ecosomg_mod = ecosomg
719  endif
720  perih_pi = atan(esinomg/ecosomg_mod) + pi * ( 0.5_rp - sign(0.5_rp, ecosomg) )
721 
722 ! if( abs(EcosOMG) < 1.E-8_RP ) then
723 ! if ( EsinOMG == 0.0_RP ) then
724 ! Perih_pi = 0.0_RP
725 ! elseif( EsinOMG < 0.0_RP ) then
726 ! Perih_pi = 1.5_RP * PI
727 ! elseif( EsinOMG > 0.0_RP ) then
728 ! Perih_pi = 0.5_RP * PI
729 ! endif
730 !
731 ! elseif( EcosOMG < 0.0_RP ) then
732 !
733 ! Perih_pi = atan(EsinOMG/EcosOMG) + PI
734 !
735 ! elseif( EcosOMG > 0.0_RP ) then
736 !
737 ! if ( EsinOMG < 0.0_RP ) then
738 ! Perih_pi = atan(EsinOMG/EcosOMG) + 2.0_RP * PI
739 ! else
740 ! Perih_pi = atan(EsinOMG/EcosOMG)
741 ! endif
742 !
743 ! endif
744  perih_pi = perih_pi / d2r ! [rad]->[degree]
745 
746  ! general precession
747  temp = 0.0_rp
748  do i = 1, nprece
749  temp = temp + prece_amp(i)*arcsec2d * sin( prece_rate(i)*arcsec2r*dyear + prece_phase(i)*d2r )
750  enddo
751  perih_psi = psi_bar * arcsec2d * dyear + zeta + temp
752 
753  ! longitude of perihelion
754  perih_omega = perih_pi + perih_psi
755 
756  do i = 1, 1000
757  if ( perih_omega + 180.0_rp < 0.0_rp ) then
758  perih_omega = perih_omega + 360.0_rp
759  elseif( perih_omega + 180.0_rp >= 360.0_rp ) then
760  perih_omega = perih_omega - 360.0_rp
761  else
762  exit
763  endif
764  enddo
765 
766  ! longitude of perigee (see berger et al.(1993))
767  omega = ( perih_omega + 180.0_rp ) * d2r
768 
769  beta = sqrt( 1.0_rp - e*e )
770 
771  ! longitude at the vernal equinox
772  lambda_m0 = 2.0_rp * ( ( 1.0_rp/2.0_rp*e + 1.0_rp/8.0_rp*e*e*e ) * ( 1.0_rp + beta ) * sin( omega) &
773  - ( 1.0_rp/4.0_rp*e*e ) * ( 1.0_rp/2.0_rp + beta ) * sin(2.0_rp*omega) &
774  + ( 1.0_rp/8.0_rp*e*e*e ) * ( 1.0_rp/3.0_rp + beta ) * sin(3.0_rp*omega) )
775 
776  return
777  end subroutine atmos_solarins_orbit
778 
779  !-----------------------------------------------------------------------------
782  Re_factor, &
783  sinDEC, &
784  cosDEC, &
785  hourangle, &
786  now_date, &
787  offset_year )
788  use scale_const, only: &
789  pi => const_pi
790  use scale_calendar, only: &
794  i_year, i_month, i_day, &
795  i_hour, i_min, i_sec
796  implicit none
797 
798  real(RP), intent(out) :: re_factor ! factor of the distance of Earth from the sun (1/rho2)
799  real(RP), intent(out) :: sinDEC ! sin/cos(solar declination)
800  real(RP), intent(out) :: cosDEC ! sin/cos(solar declination)
801  real(RP), intent(out) :: hourangle ! hour angle: relative longitude of subsolar point
802  integer, intent(in) :: now_date(6) ! date(yyyy,mm,dd,hh,mm,ss)
803  integer, intent(in) :: offset_year ! year offset
804 
805  integer :: date(6)
806  integer :: oyear
807  integer :: absday, absday_ve
808  real(DP) :: DayOfYear, abssec
809 
810  real(RP) :: lambda_m ! mean longitude from vernal equinox
811  real(RP) :: lambda ! actual longitude from vernal equinox
812  real(RP) :: nu
813  !---------------------------------------------------------------------------
814 
815  date(:) = now_date(:)
816  oyear = offset_year
817 
818  if ( atmos_solarins_fixeddate ) then
819  if( atmos_solarins_date(1) >= 0 ) oyear = 0
820  if( atmos_solarins_date(1) >= 0 ) date(1) = atmos_solarins_date(1)
821  if( atmos_solarins_date(2) >= 1 ) date(2) = atmos_solarins_date(2)
822  if( atmos_solarins_date(3) >= 1 ) date(3) = atmos_solarins_date(3)
823  if( atmos_solarins_date(4) >= 0 ) date(4) = atmos_solarins_date(4)
824  if( atmos_solarins_date(5) >= 0 ) date(5) = atmos_solarins_date(5)
825  if( atmos_solarins_date(6) >= 0 ) date(6) = atmos_solarins_date(6)
826  endif
827 
828  call calendar_getdayofyear( dayofyear, date(i_year) )
829 
830  call calendar_ymd2absday( absday, & ! [OUT]
831  date(i_year), & ! [IN]
832  date(i_month), & ! [IN]
833  date(i_day), & ! [IN]
834  oyear ) ! [IN]
835 
836  call calendar_ymd2absday( absday_ve, & ! [OUT]
837  date(i_year), & ! [IN]
838  ve_date(i_month), & ! [IN]
839  ve_date(i_day), & ! [IN]
840  oyear ) ! [IN]
841 
842  call calendar_hms2abssec( abssec, & ! [OUT]
843  date(i_hour), & ! [IN]
844  date(i_min), & ! [IN]
845  date(i_sec), & ! [IN]
846  0.0_dp ) ! [IN]
847 
848  lambda_m = lambda_m0 + 2.0_rp * pi * real(absday-absday_ve,kind=RP) / DayOfYear
849 
850  nu = lambda_m - omega
851 
852  ! 1 / (rho*rho)
853  re_factor = 1.0_rp &
854  + 2.0_rp*e * cos( nu ) &
855  + 0.5_rp*e*e * cos( 2.0_rp*nu ) &
856  + 2.5_rp*e*e &
857  + 4.0_rp*e*e*e * cos( 3.0_rp*nu )
858 
859  ! actual longitude from vernal equinox
860  lambda = lambda_m &
861  + ( 2.0_rp*e - 1.0_rp/ 4.0_rp*e*e*e ) * sin( nu ) &
862  + ( 5.0_rp/ 4.0_rp*e*e ) * sin( 2.0_rp*nu ) &
863  + ( 13.0_rp/12.0_rp*e*e*e ) * sin( 3.0_rp*nu )
864 
865  ! solar declination
866  sindec = sin(lambda) * sin(obliquity)
867  cosdec = sqrt( 1.0_rp - sindec*sindec )
868 
869  ! hour angle
870  hourangle = 2.0_rp * pi * abssec / (24.0_rp*60.0_rp*60.0_rp)
871 
872  if ( debug ) then
873  log_info("ATMOS_SOLARINS_ecliptic_longitude",*) lambda_m, nu, lambda
874  endif
875 
876  return
877  end subroutine atmos_solarins_ecliptic_longitude
878 
879  !-----------------------------------------------------------------------------
881  subroutine atmos_solarins_insolation_0d( &
882  real_lon, real_lat, &
883  now_date, offset_year, &
884  solins, cosSZA )
885  use scale_const, only: &
886  eps => const_eps
887  implicit none
888  real(RP), intent(in) :: real_lon ! longitude
889  real(RP), intent(in) :: real_lat ! latitude
890  integer, intent(in) :: now_date(6) ! date(yyyy,mm,dd,hh,mm,ss)
891  integer, intent(in) :: offset_year ! year offset
892 
893  real(RP), intent(out) :: solins ! solar insolation
894  real(RP), intent(out) :: cosSZA ! cos(Solar Zenith Angle)
895 
896  real(RP) :: Re_factor ! factor of the distance of Earth from the sun (1/rho2)
897  real(RP) :: sinDEC, cosDEC ! sin/cos(solar declination)
898  real(RP) :: hourangle ! hour angle: relative longitude of subsolar point
899 
900  real(RP) :: lon
901  real(RP) :: lat
902  !---------------------------------------------------------------------------
903 
904  call atmos_solarins_ecliptic_longitude( re_factor, & ! [OUT]
905  sindec, & ! [OUT]
906  cosdec, & ! [OUT]
907  hourangle, & ! [OUT]
908  now_date(:), & ! [IN]
909  offset_year ) ! [IN]
910 
911  if ( atmos_solarins_fixedlatlon ) then
912  lon = atmos_solarins_lon
913  lat = atmos_solarins_lat
914  else
915  lon = real_lon
916  lat = real_lat
917  endif
918 
919  cossza = sin(lat)*sindec - cos(lat)*cosdec*cos(lon+hourangle)
920  solins = atmos_solarins_constant * re_factor * ( 0.5_rp + sign(0.5_rp,cossza-eps) )
921 
922  return
923  end subroutine atmos_solarins_insolation_0d
924 
925  !-----------------------------------------------------------------------------
927  subroutine atmos_solarins_insolation_2d( &
928  IA, IS, IE, JA, JS, JE, &
929  real_lon, real_lat, &
930  now_date, offset_year, &
931  solins, cosSZA )
932  use scale_const, only: &
933  eps => const_eps
934  implicit none
935  integer, intent(in) :: IA, IS, IE
936  integer, intent(in) :: JA, JS, JE
937 
938  real(RP), intent(in) :: real_lon(ia,ja) ! longitude [rad]
939  real(RP), intent(in) :: real_lat(ia,ja) ! latitude [rad]
940  integer, intent(in) :: now_date(6) ! date(yyyy,mm,dd,hh,mm,ss)
941  integer, intent(in) :: offset_year ! year offset
942 
943  real(RP), intent(out) :: solins (ia,ja) ! solar insolation
944  real(RP), intent(out) :: cosSZA (ia,ja) ! cos(Solar Zenith Angle)
945 
946  real(RP) :: Re_factor ! factor of the distance of Earth from the sun (1/rho2)
947  real(RP) :: sinDEC, cosDEC ! sin/cos(solar declination)
948  real(RP) :: hourangle ! hour angle: relative longitude of subsolar point
949 
950  real(RP) :: lon(ia,ja)
951  real(RP) :: lat(ia,ja)
952 
953  integer :: i, j
954  !---------------------------------------------------------------------------
955 
956  call atmos_solarins_ecliptic_longitude( re_factor, & ! [OUT]
957  sindec, & ! [OUT]
958  cosdec, & ! [OUT]
959  hourangle, & ! [OUT]
960  now_date(:), & ! [IN]
961  offset_year ) ! [IN]
962 
963  if ( atmos_solarins_fixedlatlon ) then
964  !$omp parallel do OMP_SCHEDULE_ collapse(2)
965  do j = js, je
966  do i = is, ie
967  lon(i,j) = atmos_solarins_lon
968  lat(i,j) = atmos_solarins_lat
969  enddo
970  enddo
971  else
972  !$omp parallel do OMP_SCHEDULE_ collapse(2)
973  do j = js, je
974  do i = is, ie
975  lon(i,j) = real_lon(i,j)
976  lat(i,j) = real_lat(i,j)
977  enddo
978  enddo
979  endif
980 
981  !$omp parallel do OMP_SCHEDULE_ collapse(2)
982  do j = js, je
983  do i = is, ie
984  cossza(i,j) = sin(lat(i,j))*sindec - cos(lat(i,j))*cosdec*cos(lon(i,j)+hourangle)
985  solins(i,j) = atmos_solarins_constant * re_factor * ( 0.5_rp + sign(0.5_rp,cossza(i,j)-eps) )
986  enddo
987  enddo
988 
989  return
990  end subroutine atmos_solarins_insolation_2d
991 
992 end module scale_atmos_solarins
subroutine atmos_solarins_insolation_2d(IA, IS, IE, JA, JS, JE, real_lon, real_lat, now_date, offset_year, solins, cosSZA)
calc factor of Earths solar insolation
logical, public atmos_solarins_set_ve
integer, parameter, public i_month
[index] month
integer, parameter, public i_year
[index] year
real(rp), public atmos_solarins_constant
logical, public atmos_solarins_set_ideal
subroutine, public calendar_ymd2absday(absday, gyear, gmonth, gday, oyear)
Convert from gregorian date to absolute day, DAY 0 is AD1/1/1.
subroutine, public atmos_solarins_setup(basepoint_lon, basepoint_lat, iyear)
setup solar incidence module
logical, public atmos_solarins_fixeddate
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
real(rp), public atmos_solarins_lat
subroutine, public atmos_solarins_orbit(iyear)
setup solar incidence module
logical, public atmos_solarins_fixedlatlon
subroutine atmos_solarins_insolation_0d(real_lon, real_lat, now_date, offset_year, solins, cosSZA)
calc factor of Earths solar insolation
subroutine atmos_solarins_ecliptic_longitude(Re_factor, sinDEC, cosDEC, hourangle, now_date, offset_year)
calc factor of Earths solar insolation
integer, parameter, public i_min
[index] minute
module PROCESS
Definition: scale_prc.F90:11
integer, parameter, public i_sec
[index] second
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
integer, parameter, public i_hour
[index] hour
integer, parameter, public i_day
[index] day
real(rp), public atmos_solarins_lon
module profiler
Definition: scale_prof.F90:11
real(rp), public const_eps
small number
Definition: scale_const.F90:33
module PRECISION
real(rp), public const_pi
pi
Definition: scale_const.F90:31
module atmosphere / SOLARINS
module CALENDAR
integer, dimension(6), public atmos_solarins_date
subroutine, public calendar_getdayofyear(DayOfYear, iyear)
Get day of year.
subroutine, public calendar_hms2abssec(abssec, hour, minute, second, subsec)
Hour, minute, second, subsecond -> absolute second.
module STDIO
Definition: scale_io.F90:10
integer, parameter, public rp
real(rp), public atmos_solarins_eccentricity
real(rp), public atmos_solarins_obliquity