source: trunk/LMDZ.GENERIC/utilities/aerosols_utilities/math.f90 @ 3537

Last change on this file since 3537 was 3240, checked in by mturbet, 10 months ago

add optpropgen aerosol optical property code in the svn

File size: 55.3 KB
Line 
1module math
2    ! """
3    ! Contains useful math functions.
4    ! """
5    implicit none
6
7    integer, parameter :: &
8        dp = selected_real_kind(15,300)
9
10    real(kind=dp), parameter :: &
11        pi = 3.141592653589793238460_dp
12
13    doubleprecision, parameter :: &
14        sqrtpi = 0.5641895835477563d0  ! = 1 / sqrt(pi)
15
16    doubleprecision, parameter :: &
17        prec_high = 10d0 ** (-precision(0d0)), prec_low = 10d0 ** (-precision(0.))
18
19    save
20
21    contains
22        doubleprecision function chi2(observed_data, calculated_data)
23            ! """
24            ! Evaluate the goodeness of fit of calculated data in respect to observations.
25            !
26            ! input:
27            !   observed_data: array of observed data
28            !   calculated_data: array of calculated data
29            ! """
30            implicit None
31
32            doubleprecision, dimension(:), intent(in) :: observed_data, calculated_data
33
34            chi2 = sum((observed_data(:) - calculated_data(:))**2 / calculated_data(:))
35        end function chi2
36
37
38        doubleprecision function chi2_reduced(observed_data, calculated_data, input_deviation)
39            ! """
40            ! Evaluate the goodeness of fit of calculated data in respect to observations and observational errors.
41            !
42            ! input:
43            !   observed_data: array of observed data
44            !   calculated_data: array of calculated data
45            !   input_deviation: array of observational error
46            ! """
47            implicit None
48
49            doubleprecision, dimension(:), intent(in) :: observed_data, calculated_data, input_deviation
50
51            chi2_reduced = sum((observed_data(:) - calculated_data(:))**2 / input_deviation(:)**2)
52        end function chi2_reduced
53
54
55        doubleprecision function deg2rad(angle)
56            ! """
57            ! Return an angle in radians from an angle in degrees.
58            !
59            ! input:
60            !   angle: (degree) an angle
61            ! """
62            implicit none
63
64            doubleprecision, intent(in) :: angle
65
66            deg2rad = angle * (pi / 180.0_dp)
67        end function deg2rad
68
69
70        doubleprecision function ellipse_polar_form(semi_major_axis, semi_minor_axis, angle)
71            ! """
72            ! Return the ellipse polar form of an ellipse relative to its center.
73            !
74            ! inputs:
75            !   angle: (deg) angle between the semi-major axis, the center, and the point of the ellipse
76            !   semi_major_axis: semi-major axis of the ellipse
77            !   semi_minor_axis: semi-minor axis of the ellipse
78            ! """
79            implicit none
80
81            doubleprecision, intent(in) :: angle, semi_major_axis, semi_minor_axis
82
83            doubleprecision :: &
84                theta
85
86            theta = deg2rad(angle)
87
88            ellipse_polar_form = semi_major_axis * semi_minor_axis / &
89                sqrt((semi_minor_axis * cos(theta))**2 + (semi_major_axis * sin(theta))**2)
90
91            return
92
93        end function ellipse_polar_form
94
95
96        doubleprecision function gaussian(x, fwhm)
97            ! """
98            ! Return the gaussian of a value, for a given Full Width Half Maximum.
99            !
100            ! inputs:
101            !   x: value for which to calculate the gaussian
102            !   fwhm: Full Width Half Maximum of the gaussian function
103            ! """
104            implicit none
105
106            doubleprecision, intent(in) :: x, fwhm
107
108            doubleprecision :: &
109                sigma
110
111            sigma = fwhm / (2D0 * sqrt(2D0 * log(2D0)))
112            gaussian = 1D0 / (sigma * sqrt(2D0 * pi)) * exp(-1D0/2D0 * (x / sigma)**2D0)
113        end function gaussian
114
115
116        doubleprecision function gaussian_noise() result(n)
117            ! """
118            ! Return a random value following a standard normal distriubtion PDF.
119            !
120            ! output:
121            !   n: random gaussian noise
122            ! """
123            implicit none
124
125            doubleprecision ::&
126                r  ! random number between 0 and 1
127
128            call init_random_seed
129            call random_number(r)
130
131            n = sqrt(2D0) * erfinv(2 * r - 1)
132
133            return
134
135        end function gaussian_noise
136
137
138        doubleprecision function sec(angle)
139            ! """
140            ! Return the secant of an angle in degrees.
141            !
142            ! input:
143            !   angle: (degree) an angle
144            ! """
145            implicit none
146
147            doubleprecision, intent(in) :: angle
148
149            sec = 1D0 / cos(deg2rad(angle))
150        end function sec
151
152
153        doubleprecision function sgn(value)
154            ! """
155            ! Return the sign of a value.
156            ! Not to be confused with Fortran built-in function SIGN(A, B).
157            ! """
158            implicit none
159
160            doubleprecision, intent(in) :: value
161
162            if (value >= 0D0) then
163                sgn = 1D0
164            else
165                sgn = -1D0
166            end if
167
168            return
169        end function sgn
170
171
172        doubleprecision function sinc_fwhm(x, fwhm)
173            ! """
174            ! Return the sine cardinal of a value, for a given Full Width Half Maximum.
175            !
176            ! inputs:
177            !   x: value for which to calculate sinc
178            !   fwhm: Full Width Half Maximum of the sinc function
179            ! """
180            implicit none
181
182            doubleprecision,intent(in) :: x, fwhm
183
184            doubleprecision, parameter :: &
185                k = 1D0 / (2D0 * 1.89549D0)  ! sinc fwhm constant
186
187            if (x > 0D0 - tiny(0.) .and. x < 0D0 - tiny(0.)) then
188                sinc_fwhm = 1D0  ! avoid NaN at x = 0
189            else
190                sinc_fwhm = sin(x / (k * fwhm)) / (x / (k * fwhm))
191            end if
192        end function sinc_fwhm
193
194
195        function arange(start, stop, step) result(array)
196            ! """
197            ! Return evenly spaced values within a given interval.
198            !
199            ! inputs:
200            !   start: start of interval; the array starts with this value
201            !   stop: end of interval; the array does not inlude this value
202            !   step: spacing between values
203            !
204            ! outputs:
205            !   array: array of evenly spaced values
206            ! """
207            implicit none
208
209            doubleprecision, intent(in) :: start, stop, step
210            doubleprecision, dimension(:), allocatable :: array
211
212            integer :: &
213                i, &  ! index
214                n     ! number of elements in array
215
216            n = ceiling((stop - start) / step)
217
218            allocate(array(n))
219
220            do i = 1, n
221                array(i) = start + (i - 1) * step
222            end do
223
224            return
225        end function arange
226
227
228        function arange_include(start, stop, step) result(array)
229            ! """
230            ! Return evenly spaced values including a given interval.
231            !
232            ! inputs:
233            !   start: start of interval; the array starts with this value
234            !   stop: end of interval; the last element of array is always greater than this value by a max of 1 step
235            !   step: spacing between values
236            !
237            ! outputs:
238            !   array: array of evenly spaced values
239            ! """
240            implicit none
241
242            doubleprecision, intent(in) :: start, stop, step
243            doubleprecision, dimension(:), allocatable :: array
244
245            integer :: &
246                i, &  ! index
247                n     ! number of elements in array
248
249            n = ceiling((stop - start) / step) + 1  ! operation +1 ensures that stop is included in array
250
251            allocate(array(n))
252
253            do i = 1, n
254                array(i) = start + (i - 1) * step
255            end do
256
257            return
258        end function arange_include
259
260
261        function convolve(signal, filter) result(convolved_signal)
262            ! """
263            ! Convolve the signal by the filter using classical convolution.
264            !
265            ! inputs:
266            !   signal: the signal array
267            !   filter: the filter array
268            !
269            ! output:
270            !   convolved_signal: signal convolved by filter
271            !
272            ! source: https://fortrandev.wordpress.com/2013/04/01/fortran-convolution-algorithm/
273            ! """
274            implicit none
275
276            doubleprecision, dimension(:), intent(in) :: signal, filter
277            doubleprecision, dimension(:), allocatable :: convolved_signal, convolution
278
279            integer :: &
280                i, &
281                j, &
282                k, &
283                size_filter, &
284                size_signal, &
285                start_index
286
287            size_signal = size(signal)
288            size_filter = size(filter)
289
290            start_index = floor(size_filter / 2d0)
291
292            allocate(convolved_signal(size_signal), convolution(size_signal + size_filter))
293
294            ! Last part
295            do i = size_signal, size_signal + size_filter
296                convolution(i) = 0D0
297                j = size_signal
298
299                do k = 1, size_filter
300                    convolution(i) = convolution(i) + signal(j) * filter(k)
301                    j = j - 1
302                end do
303            end do
304
305            ! Middle part
306            do i = size_filter, size_signal
307                convolution(i) = 0D0
308                j = i
309
310                do k = 1, size_filter
311                    convolution(i) = convolution(i) + signal(j) * filter(k)
312                    j = j - 1
313                end do
314            end do
315
316            ! First part
317            do i = 1, size_filter
318                convolution(i) = 0D0
319                j = i
320                k = 1
321
322                do while (j > 0)
323                    convolution(i) = convolution(i) + signal(j) * filter(k)
324                    j = j - 1
325                    k = k + 1
326                end do
327            end do
328
329            convolved_signal(:) = convolution(start_index:start_index + size_signal - 1)
330        end function convolve
331
332
333        function slide_convolve(signal, filter) result(convolved_signal)
334            ! """
335            ! Slide convolve the signal by a sliding filter using classical convolution.
336            !
337            ! inputs:
338            !   signal: the signal array
339            !   filter: the filter 2D-array
340            !
341            ! output:
342            !   convolved_signal: signal convolved by filter
343            ! """
344            implicit none
345
346            doubleprecision, dimension(:), intent(in) :: signal
347            doubleprecision, dimension(:, :), intent(in) :: filter
348            doubleprecision, dimension(:), allocatable :: convolved_signal, convolution
349
350            integer :: &
351                i, &  ! index
352                j, &  ! index
353                k, &  ! index
354                size_filter, &  ! size of the filter
355                size_signal, &  ! size of the signal
356                start_index     ! start index of the result
357
358            size_signal = size(signal)
359            size_filter = size(filter, 1)
360
361            ! Check sizes
362            ! TODO [low] remove the size constraint
363            if (size_signal < size_filter) then
364                print '("ERROR: slide_convolve: filter size must be lower than signal size, &
365                        &but sizes are ", I10, "and ", I10)', size_filter, size_signal
366                print '("Catched error: ", ES15.8)', signal(size_filter)
367                stop  ! be sure the program stops here
368            end if
369
370            start_index = floor(size_filter / 2D0) + 1
371
372            allocate(convolved_signal(size_signal), convolution(size_signal + size_filter))
373
374            convolution(:) = 0D0
375
376            ! Last part
377            do i = size_signal + 1, start_index + size_signal - 1
378                j = size_signal
379
380                do k = 1, size_filter
381                    convolution(i) = convolution(i) + signal(j) * filter(k, j)
382                    j = j - 1
383                end do
384            end do
385
386            ! Middle part
387            do i = size_filter, size_signal
388                j = i
389
390                do k = 1, size_filter
391                    convolution(i) = convolution(i) + signal(j) * filter(k, j)
392                    j = j - 1
393                end do
394            end do
395
396            ! First part
397            do i = start_index + 1, size_filter - 1
398                j = i
399
400                do k = 1, i
401                    convolution(i) = convolution(i) + signal(j) * filter(k, j)
402                    j = j - 1
403                end do
404            end do
405
406            ! First index
407            i = start_index
408            convolution(i) = 0D0
409            j = i
410
411            do k = 1, i
412                convolution(i) = convolution(i) + signal(j) * filter(k, j)
413                j = j - 1
414            end do
415
416            k = i + 1
417            convolution(i) = convolution(i) + signal(1) * filter(k, 1)
418
419            convolved_signal(:) = convolution(start_index:start_index + size_signal - 1)
420        end function slide_convolve
421
422
423        function search_sorted(array, value) result(index)
424            ! """
425            ! Find the index into a sorted array such that the corresponding value is the closest to 'value'.
426            ! """
427            implicit none
428
429            doubleprecision, intent(in) :: value
430            doubleprecision, dimension(:), intent(in) :: array
431
432            integer :: index
433
434            integer :: &
435                i_low, &
436                i_high, &
437                i_mid
438
439            i_low = 1
440            i_mid = 0
441            i_high = size(array)
442
443            if(value < array(1)) then
444                index = 1
445                return
446            elseif(value > array(i_high)) then
447                index = i_high
448                return
449            end if
450
451            do while(i_low <= i_high)
452                i_mid = i_low + (i_high - i_low) / 2
453
454                if(value < array(i_mid)) then
455                    i_high = i_mid - 1
456                elseif(value > array(i_mid)) then
457                    i_low = i_mid + 1
458                else
459                    index = i_mid
460                    return
461                end if
462            end do
463
464            if(array(i_low) - value < value - array(i_high)) then
465                index = i_low
466            else
467                index = i_high
468            end if
469
470            return
471        end function search_sorted
472
473
474        function erfinv(x) result(r)
475            ! """
476            ! Calculate the inverse of the erf function.
477            !
478            ! input:
479            !   x: a number between -1 and 1
480            !
481            ! output:
482            !   r: x = erf(r)
483            ! """
484            implicit none
485
486            doubleprecision, parameter :: &  ! erfinv-approximation factors
487                a0 = 0.886226899, &
488                a1 = -1.645349621, &
489                a2 = 0.914624893, &
490                a3 = -0.140543331, &
491                b0 = 1, &
492                b1 = -2.118377725, &
493                b2 = 1.442710462, &
494                b3 = -0.329097515, &
495                b4 = 0.012229801, &
496                c0 = -1.970840454, &
497                c1 = -1.62490649, &
498                c2 = 3.429567803, &
499                c3 = 1.641345311, &
500                d0 = 1, &
501                d1 = 3.543889200, &
502                d2 = 1.637067800
503
504            doubleprecision :: x, r
505
506            integer ::&
507                sign_x  ! sign of elements of array x
508
509            doubleprecision :: &
510                x2, &  ! square of x
511                y  ! intermediate value
512
513            if(x < -1 .or. x > 1) then
514                print '("ERROR: erfinv(x): x must be in [-1;1]")'
515                stop
516            end if
517
518            if (x > -tiny(0.) .and. x < tiny(0.)) then
519                r = 0
520                return
521            end if
522
523            if (x > 0) then
524                sign_x = 1
525            else
526                sign_x = -1
527                x = -x
528            end if
529
530            if (x <= 0.7) then
531                x2 = x * x
532                r = x * (((a3 * x2 + a2) * x2 + a1) * x2 + a0)
533                r = r / (((b4 * x2 + b3) * x2 + b2) * x2 + b1) * x2 + b0
534            else
535                y = sqrt(-log((1 - x) / 2))
536                r = (((c3 * y + c2) * y + c1) * y + c0)
537                r = r / ((d2 * y + d1) * y + d0)
538            end if
539
540            r = r * sign_x
541            x = x * sign_x
542
543            r = r - (erf(r) - x) / (2 / sqrt(PI) * exp (-r * r))
544            r = r - (erf(r) - x) / (2 / sqrt(PI) * exp (-r * r))
545
546            return
547
548        end function erfinv
549
550
551        function interp(x_new, x, y) result(y_new)
552            ! """
553            ! Interpolate array y of abscisse x to array y_new of abscisse x_new.
554            !
555            ! inputs:
556            !   x_new: abscisses on which y will be interpolated
557            !   x: abscisses of y
558            !   y: array to interpole
559            !
560            ! output:
561            !   y_new: interpolation of y(x) on abscisses x_new
562            ! """
563
564            implicit none
565
566            doubleprecision, dimension(:), intent(in) :: x, y, x_new
567            doubleprecision, dimension(size(x_new)) :: y_new
568
569            integer ::&
570                i, &  ! index
571                j     ! index
572
573            i = 1
574            j = 1
575            y_new(:) = 0D0
576
577            if(size(x) /= size(y)) then
578                print '("ERROR: interp: x and y must have the same size")'
579                stop
580            end if
581
582            if(x(1) < x(size(x))) then  ! ascending numerical order x array
583                if(x_new(1) < x(1) .or. x_new(size(x_new)) > x(size(x))) then
584                    print '("ERROR: interp: x_new is outside x boundaries")'
585                    print *, x(1), ' < ' , x_new(1),  '--', x_new(size(x_new)), ' < ', x(size(x))
586                    stop
587                end if
588
589                do while (i <= size(x_new))
590                    if(x_new(i) <= x(j + 1)) then
591                        y_new(i) = (x_new(i) - x(j)) * (y(j + 1) - y(j)) / (x(j + 1) - x(j)) + y(j)
592
593                        i = i + 1
594                    else
595                        j = j + 1
596                    end if
597                end do
598            else  ! descending numerical order x array
599                if(x_new(1) > x(1) .or. x_new(size(x_new)) < x(size(x))) then
600                    print '("ERROR: interp: x_new is outside x boundaries")'
601                    print *,x(1), ' > ', x_new(1), '--', x_new(size(x_new)), ' > ', x(size(x))
602                    stop
603                end if
604
605                do while (i <= size(x_new))
606                    if(x_new(i) >= x(j + 1)) then
607                        y_new(i) = (x_new(i) - x(j)) * (y(j + 1) - y(j)) / (x(j + 1) - x(j)) + y(j)
608
609                        i = i + 1
610                    else
611                        j = j + 1
612                    end if
613                end do
614            end if
615
616            return
617        end function interp
618
619
620        function interp_fast(x_new, x, y) result(y_new)
621            ! """
622            ! Interpolate array y of abscisse x to array y_new of abscisse x_new.
623            ! Assumptions (no check performed):
624            !   - x_new is within x boundaries
625            !   - x and y have the same size
626            !
627            ! inputs:
628            !   x_new: abscisses on which y will be interpolated
629            !   x: abscisses of y
630            !   y: array to interpole
631            !
632            ! output:
633            !   y_new: interpolation of y(x) on abscisses x_new
634            ! """
635
636            implicit none
637
638            doubleprecision, dimension(:), intent(in) :: x, y, x_new
639            doubleprecision, dimension(size(x_new)) :: y_new
640
641            integer ::&
642                i, &  ! index
643                j     ! index
644
645            i = 1
646            j = 1
647            y_new(:) = 0D0
648
649            if(x(1) < x(size(x))) then  ! ascending numerical order x array
650                do while(i <= size(x_new))
651                    if(x_new(i) <= x(j + 1)) then
652                        y_new(i) = (x_new(i) - x(j)) * (y(j + 1) - y(j)) / (x(j + 1) - x(j)) + y(j)
653
654                        i = i + 1
655                    else
656                        j = j + 1
657                    end if
658                end do
659            else  ! descending numerical order x array
660                do while(i <= size(x_new))
661                    if(x_new(i) >= x(j + 1)) then
662                        y_new(i) = (x_new(i) - x(j)) * (y(j + 1) - y(j)) / (x(j + 1) - x(j)) + y(j)
663
664                        i = i + 1
665                    else
666                        j = j + 1
667                    end if
668                end do
669            end if
670
671            return
672        end function interp_fast
673
674
675        function interp_ex(x_new, x, y) result(y_new)
676            ! """
677            ! Interpolate array y of abscisse x to array y_new of abscisse x_new.
678            ! The points outside x range are linearly extrapolated.
679            ! Data must be in strict increasing order.
680            !
681            ! inputs:
682            !   x_new: abscisses on which y will be interpolated
683            !   x: abscisses of y
684            !   y: array to interpole
685            !
686            ! output:
687            !   y_new: interpolation of y(x) on abscisses x_new
688            ! """
689
690            implicit none
691
692            doubleprecision, dimension(:), intent(in) :: x, y, x_new
693            doubleprecision, dimension(size(x_new)) :: y_new
694
695            integer ::&
696                i, &  ! index
697                i_max, &
698                i_min, &
699                n     ! index
700
701            n = size(x)
702            i_min = 0
703            i_max = 0
704
705            ! Check order
706            if(x(1) > x(n)) then
707                if(x_new(1) < x_new(size(x_new))) then
708                    write(*, '("Error: interp_ex: x_new numerical order must be the same than x (descending)")')
709
710                    stop
711                end if
712
713                ! Find where x_new is in data range
714                do i = 1, size(x_new)
715                    if(x_new(i) <= x(1) .and. i_min == 0) then
716                        i_min = i
717                    end if
718
719                    if(x_new(i) >= x(n)) then
720                        i_max = i
721                    end if
722
723                    if(x_new(i) < x(n)) then
724                        exit
725                    end if
726                end do
727            else
728                if(x_new(1) > x_new(size(x_new))) then
729                    write(*, '("Error: interp_ex: x_new numerical order must be the same than x (ascending)")')
730
731                    stop
732                end if
733
734                ! Find where x_new is in data range
735                do i = 1, size(x_new)
736                    if(x_new(i) >= x(1) .and. i_min == 0) then
737                        i_min = i
738                    end if
739
740                    if(x_new(i) <= x(n)) then
741                        i_max = i
742                    end if
743
744                    if(x_new(i) > x(n)) then
745                        exit
746                    end if
747                end do
748            end if
749
750            ! Interpolation within data range
751            if(i_min > 0 .and. i_max > 0) then
752                y_new(i_min:i_max) = interp_fast(x_new(i_min:i_max), x, y)
753            end if
754
755            if(i_min == 0 .or. i_max == 0) then
756                if(i_min > 0) then
757                    ! Higher bound extrapolation
758                    do i = 1, size(x_new)
759                        y_new(i) = (x_new(i) - x(n - 1)) * (y(n) - y(n - 1)) / (x(n) - x(n - 1)) + y(n - 1)
760                    end do
761                else if (i_max > 0) then
762                    ! Lower bound extrapolation
763                    do i = 1, size(x_new)
764                        y_new(i) = (x_new(i) - x(1)) * (y(2) - y(1)) / (x(2) - x(1)) + y(1)
765                    end do
766                else
767                    write(*, '("Error: interp_ex: unable to determine boundary")')
768
769                    stop
770                end if
771            else
772                ! Lower bound extrapolation
773                if(i_min > 1) then
774                    do i = 1, i_min - 1
775                        y_new(i) = (x_new(i) - x(1)) * (y(2) - y(1)) / (x(2) - x(1)) + y(1)
776                    end do
777                end if
778
779                ! Higher bound extrapolation
780                if(i_max < size(x_new)) then
781                    do i = i_max + 1, size(x_new)
782                        y_new(i) = (x_new(i) - x(n - 1)) * (y(n) - y(n - 1)) / (x(n) - x(n - 1)) + y(n - 1)
783                    end do
784                end if
785            end if
786        end function interp_ex
787
788
789        function interp_ex_0d(x_new, x, y) result(y_new)
790            ! """
791            ! Interpolate array y of abscisse x to value y_new at x_new.
792            ! This function only takes a double as x_new, not an array. Return a double, not an array.
793            ! The points outside x range are linearly extrapolated.
794            ! Data must be in strict increasing order.
795            !
796            ! inputs:
797            !   x_new: value at which y will be interpolated
798            !   x: abscisses of y
799            !   y: array to interpole
800            !
801            ! output:
802            !   y_new: interpolation of y(x) at x_new
803            ! """
804            implicit none
805
806            doubleprecision, dimension(:), intent(in) :: x, y
807            doubleprecision :: x_new, y_new
808
809            doubleprecision, dimension(1) :: arr_tmp
810
811            arr_tmp = interp_ex([x_new], x, y)
812            y_new = arr_tmp(1)
813        end function interp_ex_0d
814
815
816        function mean_restep(x_new, x, y) result(y_new)
817            ! """
818            ! Change the values of y of abscisse x so that the new values y_new of abcisse x_new are the mean of y
819            ! between a step of x_new centered on x_new.
820            ! Both x and x_new must be regularly spaced.
821            ! Example:
822            ! >>> x_new = [0, 3, 6]
823            ! >>> x = [0, 1, 2, 3, 4, 5, 6]
824            ! >>> y = [1, 2, 4, 3, 6, 2, 0]
825            ! >>> mean_restep(x_new, x, y)
826            ! >>> [1.50, 4.33, 1.00]
827            !
828            ! inputs:
829            !   x_new: abscisses on which y will be interpolated
830            !   x: abscisses of y
831            !   y: array to interpole
832            !
833            ! output:
834            !   y_new: interpolation of y(x) on abscisses x_new
835            ! """
836            implicit none
837
838            doubleprecision, dimension(:), intent(in) :: x, y, x_new
839            doubleprecision, dimension(size(x_new)) :: y_new
840
841            integer :: &
842                i, &  ! index
843                j, &  ! index
844                n     ! number of elements
845
846            doubleprecision :: &
847                sum_x, &
848                min_x, &
849                max_x
850
851            call check_inputs()
852
853            j = 1
854
855            ! First value
856            i = 1
857            n = 0
858            sum_x = 0D0
859            max_x = x_new(i) + (x_new(i + 1) - x_new(i)) / 2D0
860            min_x = x_new(i)
861
862            do while (x(j) < max_x .and. j < size(x))
863                if (x(j) >= min_x) then
864                    sum_x = sum_x + y(j)
865                    n = n + 1
866                end if
867
868                j = j + 1
869            end do
870
871            y_new(i) = sum_x / n
872
873            ! Intermediate values
874            do i = 2, size(x_new) - 1
875                n = 0
876                sum_x = 0D0
877                max_x = x_new(i) + (x_new(i + 1) - x_new(i)) / 2D0
878                min_x = x_new(i) - (x_new(i) - x_new(i - 1)) / 2D0
879
880                do while (x(j) < max_x .and. j < size(x))
881                    if (x(j) >= min_x) then
882                        sum_x = sum_x + y(j)
883                        n = n + 1
884                    end if
885
886                    j = j + 1
887                end do
888
889                y_new(i) = sum_x / n
890            end do
891
892            ! Last value
893            i = size(x_new)
894            n = 0
895            sum_x = 0D0
896            max_x = x_new(i)
897            min_x = x_new(i) - (x_new(i) - x_new(i - 1)) / 2D0
898
899            do while (x(j) < max_x .and. j < size(x))
900                if (x(j) >= min_x) then
901                    sum_x = sum_x + y(j)
902                    n = n + 1
903                end if
904
905                j = j + 1
906            end do
907
908            if(n == 0) then
909                y_new(i) = y(j)
910            else
911                y_new(i) = sum_x / n
912            end if
913
914            return
915
916            contains
917                subroutine check_inputs()
918                    implicit none
919
920                    if (size(x) /= size(y)) then
921                        print '("ERROR: mean_restep: x and y must have the same size")'
922                        stop
923                    end if
924
925                    if (x(1) > x(size(x))) then
926                        print '("ERROR: mean_restep: x must be in increasing order")'
927                        stop
928                    end if
929
930                    if(x_new(1) < x(1) .or. x_new(size(x_new)) > x(size(x))) then
931                        print '("ERROR: interp: x_new is outside x boundaries")'
932                        print *,x_new(1),x(1),x_new(size(x_new)),x(size(x))
933                        stop
934                    end if
935                end subroutine check_inputs
936        end function mean_restep
937
938
939        function reverse_array(array) result(reversed_array)
940            ! """
941            ! Reverse a 1-D double precision array.
942            ! """
943            implicit none
944
945            doubleprecision, dimension(:), intent(in) :: array
946            integer :: head, tail
947
948            doubleprecision, dimension(size(array)) :: reversed_array
949
950            head = 1
951            tail = size(array)
952            reversed_array(:) = 0d0
953
954            do while(head < tail)
955                reversed_array(tail) = array(head)
956                reversed_array(head) = array(tail)
957                head = head + 1
958                tail = tail - 1
959            end do
960
961            return
962        end function reverse_array
963
964
965        function voigt(x, y)
966            ! """
967            ! Calculate the voigt function using an algorithm written by Humlicek JQSRT, 27, 437 (1982).
968            ! Calculate the complex probability function W(z) = exp(-z^2) * erfc(-z^2) in the superior complex plan
969            ! (i.e. for y >= 0). The real part of this function is the Voigt function.
970            !
971            ! The article shows that the rational function W can be written as (Eq. 11):
972            !   W_n(z) = (-iz) * sum_{k=1}^{n/2}(c_k((-iz)^2)^{k-1}) / ((-iz)^n + sum_{k=1}^{n/2}(d_k((-iz)^2)^{k-1}))
973            ! With z = x + iy, (-iz = t). The coefficent c and d are given in Table 3 of the article.
974            ! The maximal relative error on the imaginary and real parts is < 1e-4.
975            !
976            ! inputs:
977            !   x: real part coordinate
978            !   y: imaginary part coordinate
979            !
980            ! output:
981            !   voigt: value of the voigt function at z = x + iy
982            ! """
983            implicit none
984
985            doubleprecision, dimension(:), intent(in) :: x
986            doubleprecision, intent(in) :: y
987            doubleprecision, dimension(size(x)) :: voigt
988
989            integer :: &
990                i  ! index
991
992            double precision, dimension(size(x)) :: &
993                s  ! intermediate value
994
995            complex (kind=8), dimension(size(x)) :: &
996                t, &  ! intermediate value
997                u, &  ! intermediate value
998                w     ! intermediate value
999
1000            do i=1,size(x)
1001                t(i) = cmplx(y, -x(i), kind=8)
1002            end do
1003
1004            s(:) = abs(x(:)) + y
1005
1006            where (s >= 15D0)  ! n = 2 (region i)
1007                w = t * sqrtpi / (0.5D0 + t**2)
1008                voigt = dble(w * sqrtpi)
1009            elsewhere (s >= 5.5D0)  ! n = 4 (region II)
1010                u = t**2
1011                w = t * (1.410474D0 + u * sqrtpi) / (0.75D0 + u * (3D0 + u))
1012                voigt = dble(w * sqrtpi)
1013            elsewhere (0.195D0 * abs(x) - 0.176D0 <= y)  ! n = 6 (region III)
1014                w = (16.4955D0 + t * (20.20933D0 + t * (11.96482D0 + t * (3.778987D0 + t * 0.5642236D0)))) / &
1015                         (16.4955D0 + t * (38.82363D0 + t * (39.27121D0 + t * (21.69274D0 + t * (6.699398D0 + t)))))
1016                voigt = dble(w * sqrtpi)
1017            elsewhere  ! n = 8 (region IV)
1018                u = t**2
1019                w = exp(u) - &
1020                         t * (36183.31D0 - &
1021                              u * (3321.9905D0 - &
1022                                   u * (1540.787D0 - &
1023                                        u * (219.0313D0 - u *(35.76683D0 - u * (1.320522D0 - u * sqrtpi)))))) / &
1024                         (32066.6D0 - &
1025                          u * (24322.84D0 - &
1026                               u * (9022.228D0 - &
1027                                    u * (2186.181D0 - u * (364.2191D0 - u * (61.57037D0 - u * (1.841439D0 - u)))))))
1028                voigt = dble(w * sqrtpi)
1029            end where
1030
1031            return
1032        end function voigt
1033
1034
1035        doubleprecision function voigt_from_data(x, y, v) result(voigt)
1036            ! """
1037            ! """
1038            implicit none
1039
1040            doubleprecision, intent(in) :: &
1041                    x, &
1042                    y, &
1043                    v(400, 100)
1044
1045            doubleprecision, parameter :: &
1046                    dp = 0.0025d0, &
1047                    ds = 0.01d0, &
1048                    yl = 1d0 / 99d0, &
1049                    xl = 399d0
1050
1051            integer :: &
1052                    ip, &
1053                    is
1054
1055            doubleprecision :: &
1056                    a, &
1057                    b, &
1058                    pressure_space, &
1059                    x2, &
1060                    x2y2, &
1061                    s, &
1062                    y2, &
1063                    f
1064
1065            if(y > 5.4d0) then
1066                x2 = x * x
1067                y2 = y * y
1068                x2y2 = x2 + y2
1069                voigt = sqrtpi * y * (1d0 + (3d0 * x2 - y2) / (2d0 * x2y2 * x2y2) + &
1070                        0.75d0 * (5d0 * x2 * x2 + y2 * y2 - 10d0 * x2 * y2) / (x2y2 * x2y2 * x2y2 * x2y2)) / x2y2
1071            else
1072                if(y <= yl) then
1073                    b = exp(-x * x)
1074                    pressure_space = 1d0 / (1d0 + x)
1075                    ip = int(pressure_space / dp)
1076                    f = pressure_space / dp - ip
1077                    if(ip >= 1) then
1078                        if(ip < 400) then
1079                            a = (1d0 - f) * v(ip, 1) + f * v(ip + 1, 1)
1080                        else
1081                            a = v(ip, 1)
1082                        end if
1083                    else
1084                        a = f * v(ip + 1, 1)
1085                    end if
1086                    a = a * a
1087                    f = y / yl
1088                    voigt = (1d0 - f) * b + f * a
1089                else
1090                    if(x > xl) then
1091                        s = y / (1d0 + y)
1092                        is = int(s / ds)
1093                        b = v(1, is) * xl / x
1094                        b = b * b
1095                        a = v(1, is + 1) * xl / x
1096                        a = a * a
1097                        f = s / ds - is
1098                        voigt = (1d0 - f) * b + f * a
1099                    else
1100                        s = y / (1d0 + y)
1101                        pressure_space = 1d0 / (1d0 + x)
1102                        is = int(s / ds)
1103                        ip = int(pressure_space / dp)
1104                        f = pressure_space / dp - ip
1105                        if(ip < 400) then
1106                            b = (1d0 - f) * v(ip, is) + f * v(ip + 1, is)
1107                            a = (1d0 - f) * v(ip, is + 1) + f * v(ip + 1, is + 1)
1108                        else
1109                            b = v(ip, is)
1110                            a = v(ip, is + 1)
1111                        end if
1112                        b = b * b
1113                        a = a * a
1114                        f = s / ds - is
1115                        voigt = (1d0 - f) * b + f * a
1116                    end if
1117                end if
1118            end if
1119
1120            return
1121        end function voigt_from_data
1122
1123
1124        recursive subroutine fft(x)
1125            ! """
1126            ! Calculate the Cooley-Tukey FFT functions.
1127            !
1128            ! input:
1129            !   x: double complex vector of size 2^n
1130            !
1131            ! notes:
1132            !   Source: https://rosettacode.org/wiki/Fast_Fourier_transform#Fortran
1133            ! """
1134            implicit none
1135
1136            complex(kind=dp), dimension(:), intent(inout) :: x
1137            complex(kind=dp) ::&
1138                t
1139
1140            integer ::&
1141                i, &  ! index
1142                n     ! size of array x
1143
1144            complex(kind=dp), dimension(:), allocatable ::&
1145                even, &  ! even-number indexed values of x
1146                odd      ! odd-number indexed values of x
1147
1148            n = size(x)
1149
1150            if(n <= 1) return
1151
1152            allocate(odd((n+1)/2))
1153            allocate(even(n/2))
1154
1155            ! divide
1156            odd(:) = x(1:n:2)
1157            even(:) = x(2:n:2)
1158
1159            ! conquer
1160            call fft(odd)
1161            call fft(even)
1162
1163            ! combine
1164            do i = 1, n/2
1165               t = exp(cmplx(0.0_dp, -2.0_dp * pi * real(i-1, dp) / real(n, dp), kind=dp)) * even(i)
1166               x(i) = odd(i) + t
1167               x(i+n/2) = odd(i) - t
1168            end do
1169
1170            deallocate(odd)
1171            deallocate(even)
1172
1173        end subroutine fft
1174
1175
1176        recursive subroutine quicksort(array)
1177            ! """
1178            ! Sort double precision numbers into ascending numerical order.
1179            !
1180            ! input:
1181            !   array: double precision array to sort
1182            !
1183            ! output:
1184            !   array: sorted double precision array
1185            !
1186            ! notes:
1187            !   Author: t-nissie, some tweaks by 1AdAstra1
1188            !   Source: https://gist.github.com/t-nissie/479f0f16966925fa29ea
1189            ! """
1190            implicit none
1191
1192            double precision, intent(inout), dimension(:) :: array
1193
1194            double precision :: &
1195                x, &  ! pivot point
1196                tmp   ! temporary value
1197
1198            integer :: &
1199                first = 1, & ! index of the beginning of the array
1200                last, &      ! index of the end of the array
1201                insertion_size_threshold = 32
1202
1203            integer :: &
1204                i, &
1205                j
1206
1207            last = size(array, 1)
1208
1209            if(last < insertion_size_threshold) then  ! use insertion sort on small arrays
1210                do i = 2, last
1211                    tmp = array(i)
1212
1213                    do j = i - 1, 1, -1
1214                        if(array(j) < tmp) then
1215                            exit
1216                        else
1217                            array(j + 1) = array(j)
1218                        end if
1219                    end do
1220
1221                    array(j + 1) = tmp
1222                end do
1223            else
1224                x = array((first+last)/2)
1225                i = first
1226                j = last
1227
1228                do
1229                    do while (array(i) < x)
1230                        i = i + 1
1231                    end do
1232
1233                    do while (x < array(j))
1234                        j = j - 1
1235                    end do
1236
1237                    if (i >= j) then
1238                        exit
1239                    end if
1240
1241                    tmp = array(i)
1242                    array(i) = array(j)
1243                    array(j) = tmp
1244
1245                    i = i + 1
1246                    j = j - 1
1247                end do
1248
1249                if (first < i - 1) call quicksort(array(first:i-1))
1250                if (j + 1 < last)  call quicksort(array(j+1:last))
1251            end if
1252        end subroutine quicksort
1253
1254
1255        recursive subroutine quicksort_index(array, sorted_index)
1256            ! """
1257            ! Sort double precision numbers into ascending numerical order.
1258            !
1259            ! input:
1260            !   array: double precision array to sort
1261            !
1262            ! output:
1263            !   sorted_index: sorted index of the array
1264            ! """
1265            implicit none
1266
1267            double precision, intent(inout), dimension(:) :: array
1268            integer, intent(out), dimension(:) :: sorted_index
1269
1270            integer :: &
1271                i
1272
1273            do i = 1, size(array)
1274                sorted_index(i) = i
1275            end do
1276
1277            call partition(array, sorted_index)
1278
1279            contains
1280                recursive subroutine partition(array, sorted_index)
1281                    implicit none
1282
1283                    integer, intent(inout), dimension(:) :: sorted_index
1284                    double precision, intent(inout), dimension(:) :: array
1285
1286                    integer :: &
1287                        tmp_i
1288
1289                    double precision :: &
1290                        x, &  ! pivot point
1291                        tmp   ! temporary value
1292
1293                    integer :: &
1294                        first = 1, & ! index of the beginning of the array
1295                        last, &      ! index of the end of the array
1296                        insertion_size_threshold = 128
1297
1298                    integer :: &
1299                        i, &
1300                        j
1301
1302                    last = size(array, 1)
1303
1304                    if(last < insertion_size_threshold) then  ! use insertion sort on small arrays
1305                        do i = 2, last
1306                            tmp = array(i)
1307                            tmp_i = sorted_index(i)
1308
1309                            do j = i - 1, 1, -1
1310                                if(array(j) < tmp) then
1311                                    exit
1312                                else
1313                                    array(j + 1) = array(j)
1314                                    sorted_index(j + 1) = sorted_index(j)
1315                                end if
1316                            end do
1317
1318                                array(j + 1) = tmp
1319                                sorted_index(j + 1) = tmp_i
1320                        end do
1321                    else
1322                        x = array((first + last) / 2)
1323                        i = first
1324                        j = last
1325
1326                        do
1327                            do while (array(i) < x)
1328                                i = i + 1
1329                            end do
1330
1331                            do while (x < array(j))
1332                                j = j - 1
1333                            end do
1334
1335                            if (i >= j) then
1336                                exit
1337                            end if
1338
1339                            tmp = array(i)
1340                            array(i) = array(j)
1341                            array(j) = tmp
1342
1343                            tmp_i = sorted_index(i)
1344                            sorted_index(i) = sorted_index(j)
1345                            sorted_index(j) = tmp_i
1346
1347                            i = i + 1
1348                            j = j - 1
1349                        end do
1350
1351                        if (first < i - 1) call partition(array(first:i-1), sorted_index(first:i-1))
1352                        if (j + 1 < last)  call partition(array(j+1:last), sorted_index(j+1:last))
1353                    end if
1354                end subroutine partition
1355        end subroutine quicksort_index
1356
1357
1358        subroutine reallocate_1Ddouble(double_array, new_size)
1359            ! """
1360            ! Change the size of a 1D double precision array.
1361            ! :param double_array: array to change the size
1362            ! :param new_size: new size of the array
1363            ! :return double_array: resized array
1364            ! """
1365            implicit none
1366
1367            integer, intent(in) :: new_size
1368            doubleprecision, intent(inout), dimension(:), allocatable :: double_array
1369
1370            integer :: &
1371                size_array  ! size of the array
1372            doubleprecision, dimension(:), allocatable :: &
1373                tmp         ! temporary array
1374
1375            size_array = size(double_array)
1376
1377            call move_alloc(double_array, tmp)
1378            allocate(double_array(new_size))
1379
1380            double_array(:) = 0d0
1381
1382            size_array = min(new_size, size_array)
1383
1384            double_array(:size_array) = tmp(:size_array)
1385        end subroutine reallocate_1Ddouble
1386
1387
1388        subroutine reallocate_1Dinteger(integer_array, new_size)
1389            ! """
1390            ! Change the size of a 1D integer array.
1391            ! :param integer_array: array to change the size
1392            ! :param new_size: new size of the array
1393            ! :return integer_array: resized array
1394            ! """
1395            implicit none
1396
1397            integer, intent(in) :: new_size
1398            integer, intent(inout), dimension(:), allocatable :: integer_array
1399
1400            integer :: &
1401                size_array  ! size of the array
1402            integer, dimension(:), allocatable :: &
1403                tmp         ! temporary array
1404
1405            size_array = size(integer_array)
1406
1407            call move_alloc(integer_array, tmp)
1408            allocate(integer_array(new_size))
1409
1410            integer_array(:) = 0
1411
1412            size_array = min(new_size, size_array)
1413
1414            integer_array(:size_array) = tmp(:size_array)
1415        end subroutine reallocate_1Dinteger
1416
1417
1418        subroutine reallocate_2Ddouble(double_array, new_size_1, new_size_2)
1419            ! """
1420            ! Change the size of a 2D double precision array.
1421            ! :param double_array: array to change the size
1422            ! :param new_size_1: new size of dimension 1 of the array
1423            ! :param new_size_2: new size of dimension 2 the array
1424            ! :return double_array: resized array
1425            ! """
1426            implicit none
1427
1428            integer, intent(in) :: new_size_1, new_size_2
1429            doubleprecision, intent(inout), dimension(:, :), allocatable :: double_array
1430
1431            integer :: &
1432                shape_array(2)  ! shape of the array
1433
1434            doubleprecision, dimension(:, :), allocatable :: &
1435                tmp             ! temporary array
1436
1437            shape_array = shape(double_array)
1438
1439            call move_alloc(double_array, tmp)
1440            allocate(double_array(new_size_1, new_size_2))
1441
1442            double_array(:, :) = 0d0
1443
1444            shape_array(1) = min(new_size_1, shape_array(1))
1445            shape_array(2) = min(new_size_2, shape_array(2))
1446
1447            double_array(:shape_array(1), :shape_array(2)) = tmp(:shape_array(1), :shape_array(2))
1448        end subroutine reallocate_2Ddouble
1449
1450
1451        subroutine reallocate_2Dinteger(integer_array, new_size_1, new_size_2)
1452            ! """
1453            ! Change the size of a 2D integer array.
1454            ! :param integer_array: array to change the size
1455            ! :param new_size_1: new size of dimension 1 of the array
1456            ! :param new_size_2: new size of dimension 2 the array
1457            ! :return integer_array: resized array
1458            ! """
1459            implicit none
1460
1461            integer, intent(in) :: new_size_1, new_size_2
1462            integer, intent(inout), dimension(:, :), allocatable :: integer_array
1463
1464            integer :: &
1465                shape_array(2)  ! shape of the array
1466
1467            integer, dimension(:, :), allocatable :: &
1468                tmp             ! temporary array
1469
1470            shape_array = shape(integer_array)
1471
1472            call move_alloc(integer_array, tmp)
1473            allocate(integer_array(new_size_1, new_size_2))
1474
1475            integer_array(:, :) = 0
1476
1477            shape_array(1) = min(new_size_1, shape_array(1))
1478            shape_array(2) = min(new_size_2, shape_array(2))
1479
1480            integer_array(:shape_array(1), :shape_array(2)) = tmp(:shape_array(1), :shape_array(2))
1481        end subroutine reallocate_2Dinteger
1482
1483
1484        subroutine reallocate_3Ddouble(double_array, new_size_1, new_size_2, new_size_3)
1485            ! """
1486            ! Change the size of a 2D double precision array.
1487            ! :param double_array: array to change the size
1488            ! :param new_size_1: new size of dimension 1 of the array
1489            ! :param new_size_2: new size of dimension 2 of the array
1490            ! :param new_size_3: new size of dimension 3 of the array
1491            ! :return double_array: resized array
1492            ! """
1493            implicit none
1494
1495            integer, intent(in) :: new_size_1, new_size_2, new_size_3
1496            doubleprecision, intent(inout), dimension(:, :, :), allocatable :: double_array
1497
1498            integer :: &
1499                shape_array(3)  ! shape of the array
1500
1501            doubleprecision, dimension(:, :, :), allocatable :: &
1502                tmp             ! temporary array
1503
1504            shape_array = shape(double_array)
1505
1506            call move_alloc(double_array, tmp)
1507            allocate(double_array(new_size_1, new_size_2, new_size_3))
1508
1509            double_array(:, :, :) = 0d0
1510
1511            shape_array(1) = min(new_size_1, shape_array(1))
1512            shape_array(2) = min(new_size_2, shape_array(2))
1513            shape_array(3) = min(new_size_3, shape_array(3))
1514
1515            double_array(:shape_array(1), :shape_array(2), :shape_array(3)) = &
1516                tmp(:shape_array(1), :shape_array(2), :shape_array(3))
1517        end subroutine reallocate_3Ddouble
1518
1519
1520        subroutine ifft(x)
1521            ! """
1522            ! Calculate the Cooley-Tukey IFFT functions.
1523            !
1524            ! input:
1525            !   x: double complex vector of size 2^n
1526            ! """
1527            implicit none
1528
1529            complex(kind=dp), dimension(:), intent(inout)  :: x
1530
1531            call ifft_r(x)
1532
1533            x = x / real(size(x), dp)
1534
1535            contains
1536                recursive subroutine ifft_r(x)
1537                    implicit none
1538
1539                    complex(kind=dp), dimension(:), intent(inout) :: x
1540
1541                    integer ::&
1542                        i, &  ! index
1543                        n     ! size of array x
1544
1545                    complex(kind=dp) ::&
1546                        t  ! intermediate value
1547
1548                    complex(kind=dp), dimension(:), allocatable ::&
1549                        even, &  ! even-number indexed values of x
1550                        odd      ! odd-number indexed values of x
1551
1552                    n=size(x)
1553
1554                    if(n <= 1) return
1555
1556                    allocate(odd((n + 1) / 2))
1557                    allocate(even(n / 2))
1558
1559                    ! divide
1560                    odd(:) = x(1:n:2)
1561                    even(:) = x(2:n:2)
1562
1563                    ! conquer
1564                    call ifft_r(odd)
1565                    call ifft_r(even)
1566
1567                    ! combine
1568                    do i = 1, n/2
1569                       t = exp(cmplx(0.0_dp, 2.0_dp * pi * real(i - 1, dp) / real(n, dp), kind=dp)) * even(i)
1570                       x(i) = odd(i) + t
1571                       x(i+n/2) = odd(i) - t
1572                    end do
1573
1574                    deallocate(odd)
1575                    deallocate(even)
1576
1577                end subroutine ifft_r
1578
1579        end subroutine ifft
1580
1581
1582        subroutine init_random_seed()
1583            ! """
1584            ! Initialize the random seed with a varying seed in order to ensure a different random number sequence for
1585            ! each invocation of the program.
1586            !
1587            ! notes:
1588            !   Using the random number generator file takes far too long to be useful.
1589            !   Source: https://gcc.gnu.org/onlinedocs/gfortran/RANDOM_005fSEED.html
1590            ! """
1591            use iso_fortran_env, only: int64
1592
1593            implicit none
1594
1595            integer :: &
1596                i, &      ! index
1597                n, &      ! size of random seed
1598                dt(8), &  ! date and time
1599                pid       ! pid
1600
1601            integer, allocatable :: &
1602                seed(:)  ! random generato seed
1603
1604            integer(int64) :: &
1605                t  ! (s) time
1606
1607            call random_seed(size = n)
1608
1609            allocate(seed(n))
1610
1611            ! XOR:ing the current time and pid. The PID is useful in case one launches multiple instances of the same
1612            ! program in parallel.
1613            call system_clock(t)
1614
1615            if (t == 0) then
1616                call date_and_time(values=dt)
1617
1618                t = (dt(1) - 1970) * 365_int64 * 24 * 60 * 60 * 1000 &
1619                   + dt(2) * 31_int64 * 24 * 60 * 60 * 1000 &
1620                   + dt(3) * 24_int64 * 60 * 60 * 1000 &
1621                   + dt(5) * 60 * 60 * 1000 &
1622                   + dt(6) * 60 * 1000 + dt(7) * 1000 &
1623                   + dt(8)
1624            end if
1625
1626            pid = getpid()
1627            t = ieor(t, int(pid, kind(t)))
1628
1629            do i = 1, n
1630                seed(i) = lcg(t)
1631            end do
1632
1633            call random_seed(put=seed)
1634
1635            contains
1636                function lcg(s)
1637                    integer :: lcg
1638                    integer(int64) :: s
1639
1640                    if (s == 0) then
1641                        s = 104729
1642                    else
1643                        s = mod(s, 4294967296_int64)
1644                    end if
1645                        s = mod(s * 279470273_int64, 4294967291_int64)
1646                        lcg = int(mod(s, int(huge(0), int64)), kind(0))
1647                end function lcg
1648        end subroutine init_random_seed
1649
1650
1651        subroutine matinv(a, n)
1652            ! """
1653            ! Inverse a 2D matrix of dimension (n,n)
1654            ! The original matrix a(n,n) will be destroyed during the calculation
1655            ! Method: Based on Doolittle LU factorization for Ax=b
1656            ! Source: http://ww2.odu.edu/~agodunov/computing/programs/book2/Ch06/Inverse.f90
1657            ! by Alex G. December 2009
1658            !
1659            ! inputs:
1660            !   a(n,n): array of coefficients for matrix A
1661            !   n: dimension
1662            !
1663            ! output:
1664            !   c(n,n): inverse matrix of A
1665            ! """
1666            implicit none
1667
1668            integer, intent(in) :: n
1669            double precision, intent(inout) :: a(n,n)
1670
1671            double precision :: &
1672                l(n,n), &  ! lower triangular matrix
1673                U(n,n), &  ! upper triangular matrix
1674                b(n), &    ! intermediate value
1675                c(n,n), &  ! temporary matrix
1676                d(n), &    ! intermediate value
1677                x(n)       ! intermediate value
1678
1679            double precision :: &
1680                coeff  ! intermediate value
1681
1682            integer :: &
1683                i, &  ! index
1684                j, &  ! index
1685                k     ! index
1686
1687            ! Step 0: initialization for matrices l and U and b
1688            l(:, :) = 0d0
1689            U(:, :) = 0d0
1690            b = 0d0
1691            c(:, :) = 0d0
1692
1693            ! Step 1: forward elimination
1694            do k = 1, n - 1
1695                do i = k + 1, n
1696                    coeff = a(i, k) / a(k, k)
1697
1698                    l(i, k) = coeff
1699
1700                    do j = k + 1, n
1701                        a(i, j) = a(i, j) - coeff * a(k, j)
1702                    end do
1703                end do
1704            end do
1705
1706            ! Step 2: prepare l and U matrices
1707            ! l matrix is a matrix of the elimination coefficient + the diagonal elements are 1.0
1708            do i = 1, n
1709                l(i, i) = 1.0
1710            end do
1711
1712            ! U matrix is the upper triangular part of A
1713            do j = 1, n
1714                do i = 1, j
1715                    U(i, j) = a(i, j)
1716                end do
1717            end do
1718
1719            ! Step 3: compute columns of the inverse matrix C
1720            do k = 1, n
1721                b(k) = 1.0
1722                d(1) = b(1)
1723
1724                ! Step 3a: Solve Ld=b using the forward substitution
1725                do i = 2, n
1726                    d(i) = b(i)
1727
1728                    do j=1, i - 1
1729                        d(i) = d(i) - l(i, j) * d(j)
1730                    end do
1731                end do
1732
1733                ! Step 3b: Solve Ux=d using the back substitution
1734                x(n) = d(n) / U(n,n)
1735
1736                do i = n - 1, 1, -1
1737                    x(i) = d(i)
1738
1739                    do j = n, i + 1, -1
1740                        x(i) = x(i) - U(i, j) * x(j)
1741                    end do
1742
1743                    x(i) = x(i) / u(i, i)
1744                end do
1745
1746                ! Step 3c: fill the solutions x(n) into column k of C
1747                do i = 1, n
1748                    c(i, k) = x(i)
1749                end do
1750
1751                b(k) = 0d0
1752            end do
1753
1754            a(:, :) = c(:, :)
1755        end subroutine matinv
1756end module math
Note: See TracBrowser for help on using the repository browser.