-
Notifications
You must be signed in to change notification settings - Fork 6
/
Copy pathmo_utils.f90
2426 lines (1905 loc) · 61.6 KB
/
mo_utils.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!> \file mo_utils.f90
!> \brief General utilities for the JAMS library
!> \details This module provides general utilities such as comparisons of two reals.
!> \authors Matthias Cuntz, Juliane Mai
!> \date Feb 2014
MODULE mo_utils
! Written Matthias Cuntz, Juliane Mai, Feb 2014
! Modified Matthias Cuntz, Juliane Mai, Feb 2014 - equal, notequal
! Matthias Cuntz, May 2014 - swap
! Matthias Cuntz, May 2014 - is_finite, is_nan, is_normal, special_value
! Matthias Cuntz, Jun 2016 - special_value as elemental function
! Matthias Cuntz, Jun 2016 - cumsum, arange, linspace, imaxloc/iminloc
! Matthias Cuntz, Jun 2016 - copy toupper of mo_string_utils into module
! Matthias Cuntz, Jan 2017 - isin, isinloc
! Matthias Cuntz, Jan 2025 - remove isnan. Use version of user hkvzjal at
! https://fortran-lang.discourse.group/t/challenge-testing-inf-and-nan-with-gfortran-13-ofast/6481/30
! License
! -------
! This file is part of the JAMS Fortran package, distributed under the MIT License.
!
! Copyright (c) 2014-2017 Matthias Cuntz, Juliane Mai - mc (at) macu (dot) de
!
! Permission is hereby granted, free of charge, to any person obtaining a copy
! of this software and associated documentation files (the "Software"), to deal
! in the Software without restriction, including without limitation the rights
! to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
! copies of the Software, and to permit persons to whom the Software is
! furnished to do so, subject to the following conditions:
!
! The above copyright notice and this permission notice shall be included in all
! copies or substantial portions of the Software.
!
! THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
! IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
! FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
! AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
! LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
! OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
! SOFTWARE.
USE mo_kind, only: sp, dp, i4, i8, spc, dpc
IMPLICIT NONE
PUBLIC :: arange ! Natural numbers within interval
PUBLIC :: cumsum ! Cumulative sum
PUBLIC :: eq ! a == b, a .eq. b
#ifndef __PYTHON__
PUBLIC :: equal ! a == b, a .eq. b
#endif
PUBLIC :: ge ! a >= b, a .ge. b
#ifndef __PYTHON__
PUBLIC :: greaterequal ! a >= b, a .ge. b
#endif
PUBLIC :: imaxloc ! maxloc(arr)(1)
PUBLIC :: iminloc ! maxloc(arr)(1)
PUBLIC :: isin ! .true. if scalar present in array
PUBLIC :: isinloc ! first index of scalar in an array
PUBLIC :: is_finite ! .true. if not IEEE Inf and not IEEE NaN
PUBLIC :: is_nan ! .true. if IEEE NaN
PUBLIC :: is_normal ! .true. if not IEEE Inf and not IEEE NaN
PUBLIC :: le ! a <= b, a .le. b
#ifndef __PYTHON__
PUBLIC :: lesserequal ! a <= b, a .le. b
#endif
PUBLIC :: linspace ! Evenly spaced numbers in interval
PUBLIC :: locate ! Find closest values in a monotonic series
PUBLIC :: ne ! a /= b, a .ne. b
#ifndef __PYTHON__
PUBLIC :: notequal ! a /= b, a .ne. b
#endif
PUBLIC :: special_value ! Special IEEE values
PUBLIC :: swap ! Swaps arrays or elements of an array
! ------------------------------------------------------------------
!
! NAME
! cumsum
!
! PURPOSE
! Calculate the cumulative sum
!
!> \brief Cumulative sum.
!
!> \details The cumulative sum of the elements of an array
!> \f[ cumsum(i) = \sum_{j=1}^i array(j) \f]
!
! INTENT(IN)
!> \param[in] "integer(i4/i8)/real(sp/dp)/complex(spc/dpc) :: arr(:)" 1D array
!
! INTENT(INOUT)
! None
!
! INTENT(OUT)
! None
!
! INTENT(IN), OPTIONAL
! None
!
! INTENT(INOUT), OPTIONAL
! None
!
! INTENT(OUT), OPTIONAL
! None
!
! RETURN
!> \return kind(arr) :: cumsum(size(arr)) — Cumulative sum
!
! RESTRICTIONS
! None
!
! EXAMPLE
! vec = (/ 1., 2., 3., 4., 5., 6. /)
! cum = cumsum(vec)
! -> see also example in test directory
!
! LITERATURE
! None
!
! HISTORY
!> \authors Matthias Cuntz
!> \date Jun 2016
INTERFACE cumsum
MODULE PROCEDURE cumsum_i4, cumsum_i8, cumsum_dp, cumsum_sp, cumsum_dpc, cumsum_spc
END INTERFACE cumsum
! ------------------------------------------------------------------
!
! NAME
! equal / notequal / greaterequal / lesserequal
!
! PURPOSE
! Elemental function returning .true. or .false. depending if the reals are equal or not.
!
!> \brief Comparison of real values.
!
!> \details Compares two reals if they are numerically equal or not, i.e.
!> equal: \f[ |\frac{a-b}{b}| < \epsilon \f]
!
! INTENT(IN)
!> \param[in] "real(sp/dp) :: a" First number to compare
!> \param[in] "real(sp/dp) :: b" Second number to compare
!
! INTENT(INOUT)
! None
!
! INTENT(OUT)
! None
!
! INTENT(IN), OPTIONAL
! None
!
! INTENT(INOUT), OPTIONAL
! None
!
! INTENT(OUT), OPTIONAL
! None
!
! RETURN
!> \return real(sp/dp) :: equal — \f$ a == b \f$ logically true or false
!
! RESTRICTIONS
! None
!
! EXAMPLE
! vec1 = (/ 1., 2., 3., -999., 5., 6. /)
! vec2 = (/ 1., 1., 3., -999., 10., 6. /)
! isequal = equal(vec1, vec2)
! -> see also example in test directory
!
! LITERATURE
! None
!
! HISTORY
!> \authors Matthias Cuntz, Juliane Mai
!> \date Feb 2014
! Modified, Matthias Cuntz, Juliane Mai, Feb 2014 - sp, dp
INTERFACE eq
MODULE PROCEDURE equal_sp, equal_dp
END INTERFACE eq
INTERFACE ge
MODULE PROCEDURE greaterequal_sp, greaterequal_dp
END INTERFACE ge
INTERFACE le
MODULE PROCEDURE lesserequal_sp, lesserequal_dp
END INTERFACE le
INTERFACE ne
MODULE PROCEDURE notequal_sp, notequal_dp
END INTERFACE ne
#ifndef __PYTHON__
INTERFACE equal
MODULE PROCEDURE equal_sp, equal_dp
END INTERFACE equal
INTERFACE greaterequal
MODULE PROCEDURE greaterequal_sp, greaterequal_dp
END INTERFACE greaterequal
INTERFACE lesserequal
MODULE PROCEDURE lesserequal_sp, lesserequal_dp
END INTERFACE lesserequal
INTERFACE notequal
MODULE PROCEDURE notequal_sp, notequal_dp
END INTERFACE notequal
#endif
! ------------------------------------------------------------------
!
! NAME
! imaxloc / iminloc
!
! PURPOSE
! First index location in array of element with the maximum/minimum value.
!
!> \brief First index location in an array of the element with the maximum/minimum value.
!
!> \details Fortran intrinsics maxloc and minloc return arrays with all indexes
!> corresponding to the maximum/minimum value in an array.\n
!> This routine returns only the first entry as scalar integer.
!
! INTENT(IN)
!> \param[in] "integer(i4/i8)/real(sp/dp)/complex(spc/dpc) :: array(:)" Input array
!
! INTENT(INOUT)
! None
!
! INTENT(OUT)
! None
!
! INTENT(IN), OPTIONAL
!> \param[in] "logical :: mask(:)" If present, only those locations in array corresponding to
!> the true values in mask are searched for the maximum value.
!
! INTENT(INOUT), OPTIONAL
! None
!
! INTENT(OUT), OPTIONAL
! None
!
! RETURN
!> \return integer(i4) :: imaxloc/iminloc — First index location of maximum/minimum
!
! RESTRICTIONS
! None
!
! EXAMPLE
! integer(i4) :: imin
! imin = iminloc(vec, mask=mask)
! -> see also example in test directory
!
! LITERATURE
! None
!
! HISTORY
!> \authors Matthias Cuntz, Juliane Mai
!> \date Feb 2014
! Modified, Matthias Cuntz, Juliane Mai, Feb 2014 - sp, dp
INTERFACE imaxloc
MODULE PROCEDURE imaxloc_i4, imaxloc_i8, imaxloc_sp, imaxloc_dp
END INTERFACE imaxloc
INTERFACE iminloc
MODULE PROCEDURE iminloc_i4, iminloc_i8, iminloc_sp, iminloc_dp
END INTERFACE iminloc
! ------------------------------------------------------------------
!
! NAME
! isin
!
! PURPOSE
! Return true if one element of an array corresponds to a scalar,
! false otherwise.
!
!> \brief True if scalar is present in array.
!
!> \details Ask 'Is this scalar present in the array?'
!> Returns .true. if one element of the array corresponds
!> to the scalar value.\n
!> This is basically any(array == scalar) but works also with
!> floating point variables and character strings.\n
!> Leading and trailing blank characters are removed from strings.
!
! INTENT(IN)
!> \param[in] "integer(i4/i8)/real(sp/dp)/character(len=*) :: scalar" Single scalar value
!> \param[in] "integer(i4/i8)/real(sp/dp)/character(len=*) :: array(:)" Input vector
!
! INTENT(INOUT)
! None
!
! INTENT(OUT)
! None
!
! INTENT(IN), OPTIONAL
!> \param[in] "logical :: mask(:)" If present, only those locations in array corresponding to
!> the true values in mask are searched for the scalar value.
!> \param[in] "logical :: ignore_case" If .true., ignore case in comparison of strings;
!> if .false., comparison case sensitive (default: .true.)
!
! INTENT(INOUT), OPTIONAL
! None
!
! INTENT(OUT), OPTIONAL
! None
!
! RETURN
!> \return logical :: isin — .true. if scalar present in array, .false. otherwise
!
! RESTRICTIONS
! Only 1D-arrays.
!
! EXAMPLE
! sca = 1.1
! vec = (/ 0.0, 1.1, 2.2, 3.3 /)
! if (isin(sca, vec)) write(*,*) 'It is in.'
! -> see also example in test directory
!
! LITERATURE
! None
!
! HISTORY
!> \authors Matthias Cuntz
!> \date Jan 2017
!> Modified Matthias Cuntz, Mar 2018 - ignore_case for isin_char
INTERFACE isin
MODULE PROCEDURE isin_i4, isin_i8, isin_sp, isin_dp, isin_char
END INTERFACE isin
! ------------------------------------------------------------------
!
! NAME
! isinloc
!
! PURPOSE
! Returns the first index location of the element of an array that corresponds
! to a given scalar, 0 otherwise.
!
!> \brief First index location of scalar in array.
!
!> \details Returns the first index location in an array where an array element
!> matches a given scalar. Returns 0 if the element is not present in
!> the array.\n
!> Leading and trailing blank characters are removed from strings.
!
! INTENT(IN)
!> \param[in] "integer(i4/i8)/real(sp/dp)/character(len=*) :: scalar" Single scalar value
!> \param[in] "integer(i4/i8)/real(sp/dp)/character(len=*) :: array(:)" Input vector
!
! INTENT(INOUT)
! None
!
! INTENT(OUT)
! None
!
! INTENT(IN), OPTIONAL
!> \param[in] "logical :: mask(:)" If present, only those locations in array corresponding to
!> the true values in mask are searched for the scalar value.
!> \param[in] "logical :: ignore_case" If .true., ignore case in comparison of strings;
!> if .false., comparison case sensitive (default: .true.)
!
! INTENT(INOUT), OPTIONAL
! None
!
! INTENT(OUT), OPTIONAL
! None
!
! RETURN
!> \return integer(i4) :: isinloc — First index location of scalar in array, 0 otherwise
!
! RESTRICTIONS
! Only 1D-arrays.
!
! EXAMPLE
! sca = 1.1
! vec = (/ 0.0, 1.1, 2.2, 3.3 /)
! ii = isinloc(sca, vec)
! -> see also example in test directory
!
! LITERATURE
! None
!
! HISTORY
!> \authors Matthias Cuntz
!> \date Jan 2017
INTERFACE isinloc
MODULE PROCEDURE isinloc_i4, isinloc_i8, isinloc_sp, isinloc_dp, isinloc_char
END INTERFACE isinloc
! ------------------------------------------------------------------
!
! NAME
! is_finite / is_nan / is_normal
!
! PURPOSE
! Elemental inquiry functions returning .true. if the argument has a value
! implied by the name of the function.
!
!> \brief .true. if not IEEE Inf, IEEE NaN, nor IEEE Inf nor IEEE NaN, respectively.
!
!> \details Checks for IEEE Inf and IEEE NaN, i.e. Infinity and Not-a-Number.\n
!> Wraps to functions of the intrinsic module ieee_arithmetic
!> but gives alternatives for gfortran, which does not provide ieee_arithmetic.
!
! INTENT(IN)
!> \param[in] "real(sp/dp) :: x" Number to check
!
! INTENT(INOUT)
! None
!
! INTENT(OUT)
! None
!
! INTENT(IN), OPTIONAL
! None
!
! INTENT(INOUT), OPTIONAL
! None
!
! INTENT(OUT), OPTIONAL
! None
!
! RETURN
!> \return logical :: is_finite/is_nan/is_normal — \f$ a /= Inf, a == NaN, a /= Inf and a == NaN \f$,
!> logically true or false
!
! RESTRICTIONS
! None
!
! EXAMPLE
! vec1 = (/ NaN, 2., 3., Inf, 5., 6. /)
! is_finite = equal(vec1)
! is_nan = equal(vec1)
! is_normal = equal(vec1)
! -> see also example in test directory
!
! LITERATURE
! None
!
! HISTORY
!> \authors Matthias Cuntz
!> \date Mar 2015
INTERFACE is_finite
MODULE PROCEDURE is_finite_sp, is_finite_dp
END INTERFACE is_finite
INTERFACE is_nan
MODULE PROCEDURE is_nan_sp, is_nan_dp
END INTERFACE is_nan
INTERFACE is_normal
MODULE PROCEDURE is_normal_sp, is_normal_dp
END INTERFACE is_normal
! ------------------------------------------------------------------
!
! NAME
! linspace
!
! PURPOSE
! Return evenly spaced numbers over a specified interval.
!
!> \brief Evenly spaced numbers in interval.
!
!> \details Return N evenly spaced numbers over a specified interval [lower,upper].
!> \f[ linspace(lower,upper,N) = lower + arange(0,N-1)/(N-1) * (upper-lower) \f]
!
!> Output array has kind of lower.
!
! INTENT(IN)
!> \param[in] "integer(i4/i8)/real(sp/dp) :: lower" Start of interval.
!> \param[in] "integer(i4/i8)/real(sp/dp) :: upper" End of interval.
!> \param[in] "integer(i4) :: nstep" Number of steps.
!
! INTENT(INOUT)
! None
!
! INTENT(OUT)
! None
!
! INTENT(IN), OPTIONAL
! None
!
! INTENT(INOUT), OPTIONAL
! None
!
! INTENT(OUT), OPTIONAL
! None
!
! RETURN
!> \return kind(lower) :: linspace(N) — 1D array with evenly spaced numbers between lower and upper.
!
! RESTRICTIONS
! None
!
! EXAMPLE
! rr = linspace(1.0_dp,11._dp,101)
! -> see also example in test directory
!
! LITERATURE
! None
!
! HISTORY
!> \authors Matthias Cuntz
!> \date Jun 2016
INTERFACE linspace
MODULE PROCEDURE linspace_i4, linspace_i8, linspace_dp, linspace_sp
END INTERFACE linspace
! ------------------------------------------------------------------
!
! NAME
! locate
!
! PURPOSE
! Find closest values in a monotonic series
!
!> \brief Find closest values in a monotonic series, returns the indexes.
!
!> \details Given an array x(1:n), and given a value y,
!> returns a value j such that y is between
!> x(j) and x(j+1).\n
!
!> x must be monotonically increasing.\n
!> j=0 or j=N is returned to indicate that x is out of range.
!
! INTENT(IN)
!> \param[in] "real(dp/sp) :: x(:)" Sorted array
!> \param[in] "real(dp/sp) :: y[(:)]" Value(s) of which the closest match in x(:) is wanted
!
! INTENT(INOUT)
! None
!
! INTENT(OUT)
! None
!
! INTENT(IN), OPTIONAL
! None
!
! INTENT(INOUT), OPTIONAL
! None
!
! INTENT(OUT), OPTIONAL
! None
!
! RETURN
!> \return integer(i4) :: index[(:)] — index(es) of x so that y is between x(index) and x(index+1)
!
! RESTRICTIONS
!> \note x must be monotonically increasing.\n
!
! EXAMPLE
! x = (/ 1., 2., 3., -999., 5., 6. /)
! y = (/ 1.1, 5.6 /)
! ii = locate(x, y)
! -> ii == (/ 1, 5 /)
! y = 1.1
! ii = locate(x, y)
! -> ii == 1
! -> see also example in test directory
!
! LITERATURE
! None
!
! HISTORY
!> \author Matthias Cuntz
!> \date May 2014
INTERFACE locate
MODULE PROCEDURE locate_0d_dp, locate_0d_sp, locate_1d_dp, locate_1d_sp
END INTERFACE locate
! ------------------------------------------------------------------
!
! NAME
! arange
!
! PURPOSE
! Gives natural numbers within a given interval.
!
!> \brief Numbers within a given range.
!
!> \details Gives array with numbers in a given interval, i.e.
!> \f[ arange(1) = lower \f]
!> \f[ arange(2) = lower+1 \f]
!> ...
!> \f[ arange(n) = upper \f]
!
!> Default is lower=1.
!
!> Output array has kind of lower.
!
! INTENT(IN)
!> \param[in] "integer(i4/i8)/real(sp/dp) :: lower" Start of interval if upper is given,
!> Otherwise end of interval and start of interval is 1.
!
! INTENT(INOUT)
! None
!
! INTENT(OUT)
! None
!
! INTENT(IN), OPTIONAL
!> \param[in] "integer(i4/i8)/real(sp/dp) :: upper End of interval"
!
! INTENT(INOUT), OPTIONAL
! None
!
! INTENT(OUT), OPTIONAL
! None
!
! RETURN
!> \return kind(arr) :: arange(upper-lower+1) — 1D array with values within given interval.
!
! RESTRICTIONS
! None
!
! EXAMPLE
! rr = arange(100._dp)
! -> see also example in test directory
!
! LITERATURE
! None
!
! HISTORY
!> \authors Matthias Cuntz
!> \date Jun 2016
INTERFACE arange
MODULE PROCEDURE arange_i4, arange_i8, arange_dp, arange_sp
END INTERFACE arange
! ------------------------------------------------------------------
!
! NAME
! swap
!
! PURPOSE
! Swap two values/arrays or two elements in 1D-array.
!
!> \brief Swap two values or exchange two elements in array.
!
!> \details Swaps either two entities, i.e. scalars, vectors, matrices,
!> or exchanges two elements in a vector.
!> If an optinal mask is given, the only elements with mask==.true. will be exchanged.\n
!> The call is either \n
!> call swap(x, y, mask=mask) \n
!> or \n
!> call swap(vec, i, j, mask=mask)
!
! INTENT(IN)
!> \param[in] "integer(i4) :: i" Index of first element to be swapped with second [case swap(vec,i,j)]
!> \param[in] "integer(i4) :: j" Index of second element to be swapped with first [case swap(vec,i,j)]
!
! INTENT(INOUT)
!> \param[inout] "real(sp/dp)/integer(i4)/complex(spc/dpc) :: x[(:,...)]"
!> First scalar or array to swap with second [case swap(x,y)]
!> \param[inout] "real(sp/dp)/integer(i4)/complex(spc/dpc) :: y[(:[,:])]"
!> Second scalar or array to swap with first [case swap(x,y)]
!>
!> \param[inout] "real(sp/dp)/integer(i4)/complex(spc/dpc) :: x(:)"
!> Vector of which to elements are swapped [case swap(vec,i,j)]
!
! INTENT(OUT)
! None
!
! INTENT(IN), OPTIONAL
!> \param[in] "logical, optional :: mask[(:,...)]" scalar or array logical mask\n
!> If present, only those elements will be swapped
!> where mask==.true.
!
! INTENT(INOUT), OPTIONAL
! None
!
! INTENT(OUT), OPTIONAL
! None
!
! RESTRICTIONS
! None
!
! EXAMPLE
! vec1 = (/ 1., 2., 3., -999., 5., 6. /)
! vec2 = (/ 1., 1., 3., -999., 10., 6. /)
! call swap(vec1, vec2, mask=(vec==-999.))
! call swap(vec1, 1, 3)
! -> see also example in test directory
!
! LITERATURE
! None
!
! HISTORY
!> \author Matthias Cuntz
!> \date May 2014
INTERFACE swap
MODULE PROCEDURE &
swap_xy_dp, swap_xy_sp, swap_xy_i4, swap_xy_dpc, swap_xy_spc, &
swap_xy_mask_dp, swap_xy_mask_sp, swap_xy_mask_i4, swap_xy_mask_dpc, swap_xy_mask_spc, &
swap_vec_dp, swap_vec_sp, swap_vec_i4, swap_vec_dpc, swap_vec_spc,&
swap_vec_mask_dp, swap_vec_mask_sp, swap_vec_mask_i4, swap_vec_mask_dpc, swap_vec_mask_spc
END INTERFACE swap
! ------------------------------------------------------------------
!
! NAME
! special_value
!
! PURPOSE
! Mimics the function ieee_value of the intrinsic module ieee_arithmetic.
!
!> \brief Special IEEE values.
!
!> \details Returns special IEEE values such as Infinity or Not-a-Number.\n
!> Wraps to function ieee_value of the intrinsic module ieee_arithmetic
!> but gives alternatives for gfortran, which does not provide ieee_arithmetic.\n
!> Quiet and signaling NaN are the same in case of gfortran;\n
!> also denormal values are the same as inf.
!>
!> Current special values are:\n
!> IEEE_SIGNALING_NAN\n
!> IEEE_QUIET_NAN\n
!> IEEE_NEGATIVE_INF\n
!> IEEE_POSITIVE_INF\n
!> IEEE_NEGATIVE_DENORMAL\n
!> IEEE_POSITIVE_DENORMAL\n
!> IEEE_NEGATIVE_NORMAL\n
!> IEEE_POSITIVE_NORMAL\n
!> IEEE_NEGATIVE_ZERO\n
!> IEEE_POSITIVE_ZERO
!
! INTENT(IN)
!> \param[in] "real(sp/dp) :: x" dummy for kind of output
!> \param[in] "character(le=*) :: name ieee signal nanme
!
! INTENT(INOUT)
! None
!
! INTENT(OUT)
!
! INTENT(IN), OPTIONAL
! None
!
! INTENT(INOUT), OPTIONAL
! None
!
! INTENT(OUT), OPTIONAL
! None
!
! RETURN
!> \return real(sp/dp) :: special_value — IEEE special value\n
!> IEEE_SIGNALING_NAN\n
!> IEEE_QUIET_NAN (==IEEE_SIGNALING_NAN for gfortran)\n
!> IEEE_NEGATIVE_INF\n
!> IEEE_POSITIVE_INF\n
!> IEEE_NEGATIVE_DENORMAL (==-0.0 for gfortran)\n
!> IEEE_POSITIVE_DENORMAL (==0.0 for gfortran)\n
!> IEEE_NEGATIVE_NORMAL (==-1.0 for gfortran)\n
!> IEEE_POSITIVE_NORMAL (==1.0 for gfortran)\n
!> IEEE_NEGATIVE_ZERO\n
!> IEEE_POSITIVE_ZERO\n
!
! RESTRICTIONS
! None
!
! EXAMPLE
! NaN = special_value(1.0, 'IEEE_QUIET_NAN')
! nan = special_value(1.0_dp, 'ieee_quiet_nan')
! -> see also example in test directory
!
! LITERATURE
! None
!
! HISTORY
!> \authors Matthias Cuntz
!> \date Mar 2015
INTERFACE special_value
MODULE PROCEDURE special_value_sp, special_value_dp
END INTERFACE special_value
! ------------------------------------------------------------------
PRIVATE
! ------------------------------------------------------------------
CONTAINS
! ------------------------------------------------------------------
function arange_i4(lower, upper)
implicit none
integer(i4), intent(in) :: lower
integer(i4), intent(in), optional :: upper
integer(i4), dimension(:), allocatable :: arange_i4
integer(i4) :: istart, istop
integer(i4) :: i
if (present(upper)) then
istart = lower
istop = upper
else
istart = 1_i4
istop = lower
endif
allocate(arange_i4(istop-istart+1_i4))
forall(i=istart:istop) arange_i4(i-istart+1) = i
end function arange_i4
function arange_i8(lower, upper)
implicit none
integer(i8), intent(in) :: lower
integer(i8), intent(in), optional :: upper
integer(i8), dimension(:), allocatable :: arange_i8
integer(i8) :: istart, istop
integer(i8) :: i
if (present(upper)) then
istart = lower
istop = upper
else
istart = 1_i8
istop = lower
endif
allocate(arange_i8(istop-istart+1_i8))
forall(i=istart:istop) arange_i8(i-istart+1) = i
end function arange_i8
function arange_dp(lower, upper)
implicit none
real(dp), intent(in) :: lower
real(dp), intent(in), optional :: upper
real(dp), dimension(:), allocatable :: arange_dp
integer(i8) :: istart, istop
integer(i8) :: i
if (present(upper)) then
istart = int(lower,i8)
istop = int(upper,i8)
else
istart = 1_i8
istop = int(lower,i8)
endif
allocate(arange_dp(istop-istart+1_i8))
forall(i=istart:istop) arange_dp(i-istart+1) = real(i,dp)
end function arange_dp
function arange_sp(lower, upper)
implicit none
real(sp), intent(in) :: lower
real(sp), intent(in), optional :: upper
real(sp), dimension(:), allocatable :: arange_sp
integer(i8) :: istart, istop
integer(i8) :: i
if (present(upper)) then
istart = int(lower,i8)
istop = int(upper,i8)
else
istart = 1_i8
istop = int(lower,i8)
endif
allocate(arange_sp(istop-istart+1_i8))
forall(i=istart:istop) arange_sp(i-istart+1) = real(i,sp)
end function arange_sp
! ------------------------------------------------------------------
function cumsum_i4(arr)
implicit none
integer(i4), dimension(:), intent(in) :: arr
integer(i4), dimension(size(arr,1)) :: cumsum_i4
integer(i4) :: i
cumsum_i4(1) = arr(1)
do i=2, size(arr)
cumsum_i4(i) = cumsum_i4(i-1) + arr(i)
end do
end function cumsum_i4
function cumsum_i8(arr)
implicit none
integer(i8), dimension(:), intent(in) :: arr
integer(i8), dimension(size(arr,1)) :: cumsum_i8
integer(i4) :: i
cumsum_i8(1) = arr(1)
do i=2, size(arr)
cumsum_i8(i) = cumsum_i8(i-1) + arr(i)
end do
end function cumsum_i8
function cumsum_dp(arr)
implicit none
real(dp), dimension(:), intent(in) :: arr
real(dp), dimension(size(arr,1)) :: cumsum_dp
integer(i4) :: i
cumsum_dp(1) = arr(1)
do i=2, size(arr)
cumsum_dp(i) = cumsum_dp(i-1) + arr(i)
end do
end function cumsum_dp
function cumsum_dpc(arr)
implicit none
complex(dpc), dimension(:), intent(in) :: arr
complex(dpc), dimension(size(arr,1)) :: cumsum_dpc
integer(i4) :: i
cumsum_dpc(1) = arr(1)
do i=2, size(arr)
cumsum_dpc(i) = cumsum_dpc(i-1) + arr(i)
end do
end function cumsum_dpc
function cumsum_sp(arr)
implicit none
real(sp), dimension(:), intent(in) :: arr
real(sp), dimension(size(arr,1)) :: cumsum_sp
integer(i4) :: i
cumsum_sp(1) = arr(1)
do i=2, size(arr)
cumsum_sp(i) = cumsum_sp(i-1) + arr(i)
end do
end function cumsum_sp
function cumsum_spc(arr)
implicit none
complex(spc), dimension(:), intent(in) :: arr
complex(spc), dimension(size(arr,1)) :: cumsum_spc
integer(i4) :: i
cumsum_spc(1) = arr(1)
do i=2, size(arr)
cumsum_spc(i) = cumsum_spc(i-1) + arr(i)
end do
end function cumsum_spc
! ------------------------------------------------------------------
ELEMENTAL PURE FUNCTION equal_dp(a, b)