-
Notifications
You must be signed in to change notification settings - Fork 14
/
module_map_utils.F90
2214 lines (1707 loc) · 74 KB
/
module_map_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
MODULE map_utils_mod
! Module that defines constants, data structures, and
! subroutines used to convert grid indices to lat/lon
! and vice versa.
!
! SUPPORTED PROJECTIONS
! ---------------------
! Cylindrical Lat/Lon (code = PROJ_LATLON)
! Mercator (code = PROJ_MERC)
! Lambert Conformal (code = PROJ_LC)
! Gaussian (code = PROJ_GAUSS)
! Polar Stereographic (code = PROJ_PS)
! Rotated Lat/Lon (code = PROJ_ROTLL)
!
! REMARKS
! -------
! The routines contained within were adapted from routines
! obtained from NCEP's w3 library. The original NCEP routines were less
! flexible (e.g., polar-stereo routines only supported truelat of 60N/60S)
! than what we needed, so modifications based on equations in Hoke, Hayes, and
! Renninger (AFGWC/TN/79-003) were added to improve the flexibility.
! Additionally, coding was improved to F90 standards and the routines were
! combined into this module.
!
! ASSUMPTIONS
! -----------
! Grid Definition:
! For mercator, lambert conformal, and polar-stereographic projections,
! the routines within assume the following:
!
! 1. Grid is dimensioned (i,j) where i is the East-West direction,
! positive toward the east, and j is the north-south direction,
! positive toward the north.
! 2. Origin is at (1,1) and is located at the southwest corner,
! regardless of hemispere.
! 3. Grid spacing (dx) is always positive.
! 4. Values of true latitudes must be positive for NH domains
! and negative for SH domains.
!
! For the latlon and Gaussian projection, the grid origin may be at any
! of the corners, and the deltalat and deltalon values can be signed to
! account for this using the following convention:
! Origin Location Deltalat Sign Deltalon Sign
! --------------- ------------- -------------
! SW Corner + +
! NE Corner - -
! NW Corner - +
! SE Corner + -
!
! Data Definitions:
! 1. Any arguments that are a latitude value are expressed in
! degrees north with a valid range of -90 -> 90
! 2. Any arguments that are a longitude value are expressed in
! degrees east with a valid range of -180 -> 180.
! 3. Distances are in meters and are always positive.
! 4. The standard longitude (stdlon) is defined as the longitude
! line which is parallel to the grid's y-axis (j-direction), along
! which latitude increases (NOT the absolute value of latitude, but
! the actual latitude, such that latitude increases continuously
! from the south pole to the north pole) as j increases.
! 5. One true latitude value is required for polar-stereographic and
! mercator projections, and defines at which latitude the
! grid spacing is true. For lambert conformal, two true latitude
! values must be specified, but may be set equal to each other to
! specify a tangent projection instead of a secant projection.
!
! USAGE
! -----
! To use the routines in this module, the calling routines must have the
! following statement at the beginning of its declaration block:
! USE map_utils
!
! The use of the module not only provides access to the necessary routines,
! but also defines a structure of TYPE (proj_info) that can be used
! to declare a variable of the same type to hold your map projection
! information. It also defines some integer parameters that contain
! the projection codes so one only has to use those variable names rather
! than remembering the acutal code when using them. The basic steps are
! as follows:
!
! 1. Ensure the "USE map_utils" is in your declarations.
! 2. Declare the projection information structure as type(proj_info):
! TYPE(proj_info) :: proj
! 3. Populate your structure by calling the map_set routine:
! CALL map_set(code,lat1,lon1,knowni,knownj,dx,stdlon,truelat1,truelat2,proj)
! where:
! code (input) = one of PROJ_LATLON, PROJ_MERC, PROJ_LC, PROJ_PS,
! PROJ_GAUSS, or PROJ_ROTLL
! lat1 (input) = Latitude of grid origin point (i,j)=(1,1)
! (see assumptions!)
! lon1 (input) = Longitude of grid origin
! knowni (input) = origin point, x-location
! knownj (input) = origin point, y-location
! dx (input) = grid spacing in meters (ignored for LATLON projections)
! stdlon (input) = Standard longitude for PROJ_PS and PROJ_LC,
! deltalon (see assumptions) for PROJ_LATLON,
! ignored for PROJ_MERC
! truelat1 (input) = 1st true latitude for PROJ_PS, PROJ_LC, and
! PROJ_MERC, deltalat (see assumptions) for PROJ_LATLON
! truelat2 (input) = 2nd true latitude for PROJ_LC,
! ignored for all others.
! proj (output) = The structure of type (proj_info) that will be fully
! populated after this call
!
! 4. Now that the proj structure is populated, you may call either
! of the following routines:
!
! latlon_to_ij(proj, lat, lon, i, j)
! ij_to_latlon(proj, i, j, lat, lon)
!
! It is incumbent upon the calling routine to determine whether or
! not the values returned are within your domain's bounds. All values
! of i, j, lat, and lon are REAL values.
!
!
! REFERENCES
! ----------
! Hoke, Hayes, and Renninger, "Map Preojections and Grid Systems for
! Meteorological Applications." AFGWC/TN-79/003(Rev), Air Weather
! Service, 1985.
!
! NCAR MM5v3 Modeling System, REGRIDDER program, module_first_guess_map.F
! NCEP routines w3fb06, w3fb07, w3fb08, w3fb09, w3fb11, w3fb12
!
! HISTORY
! -------
! 27 Mar 2001 - Original Version
! Brent L. Shaw, NOAA/FSL (CSU/CIRA)
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
use constants_module
use misc_definitions_module
use utils_mod
! Define some private constants
INTEGER, PRIVATE, PARAMETER :: HIGH = 8
TYPE proj_info
INTEGER :: code ! Integer code for projection TYPE
INTEGER :: nlat ! For Gaussian -- number of latitude points
! north of the equator
INTEGER :: nlon !
!
INTEGER :: nxmin ! Starting x-coordinate of periodic, regular lat/lon dataset
INTEGER :: nxmax ! Ending x-coordinate of periodic, regular lat/lon dataset
INTEGER :: ixdim ! For Rotated Lat/Lon -- number of mass points
! in an odd row
INTEGER :: jydim ! For Rotated Lat/Lon -- number of rows
INTEGER :: stagger ! For Rotated Lat/Lon -- mass or velocity grid
REAL :: phi ! For Rotated Lat/Lon -- domain half-extent in
! degrees latitude
REAL :: lambda ! For Rotated Lat/Lon -- domain half-extend in
! degrees longitude
REAL :: lat1 ! SW latitude (1,1) in degrees (-90->90N)
REAL :: lon1 ! SW longitude (1,1) in degrees (-180->180E)
REAL :: lat0 ! For Cassini, latitude of projection pole
REAL :: lon0 ! For Cassini, longitude of projection pole
REAL :: dx ! Grid spacing in meters at truelats, used
! only for ps, lc, and merc projections
REAL :: dy ! Grid spacing in meters at truelats, used
! only for ps, lc, and merc projections
REAL :: latinc ! Latitude increment for cylindrical lat/lon
REAL :: loninc ! Longitude increment for cylindrical lat/lon
! also the lon increment for Gaussian grid
REAL :: dlat ! Lat increment for lat/lon grids
REAL :: dlon ! Lon increment for lat/lon grids
REAL :: stdlon ! Longitude parallel to y-axis (-180->180E)
REAL :: truelat1 ! First true latitude (all projections)
REAL :: truelat2 ! Second true lat (LC only)
REAL :: hemi ! 1 for NH, -1 for SH
REAL :: cone ! Cone factor for LC projections
REAL :: polei ! Computed i-location of pole point
REAL :: polej ! Computed j-location of pole point
REAL :: rsw ! Computed radius to SW corner
REAL :: rebydx ! Earth radius divided by dx
REAL :: knowni ! X-location of known lat/lon
REAL :: knownj ! Y-location of known lat/lon
REAL :: re_m ! Radius of spherical earth, meters
REAL :: rho0 ! For Albers equal area
REAL :: nc ! For Albers equal area
REAL :: bigc ! For Albers equal area
LOGICAL :: init ! Flag to indicate if this struct is
! ready for use
LOGICAL :: wrap ! For Gaussian -- flag to indicate wrapping
! around globe?
LOGICAL :: comp_ll ! Work in computational lat/lon space for Cassini
REAL, POINTER, DIMENSION(:) :: gauss_lat ! Latitude array for Gaussian grid
END TYPE proj_info
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
CONTAINS
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
SUBROUTINE map_init(proj)
! Initializes the map projection structure to missing values
IMPLICIT NONE
TYPE(proj_info), INTENT(INOUT) :: proj
proj%lat1 = -999.9
proj%lon1 = -999.9
proj%lat0 = -999.9
proj%lon0 = -999.9
proj%dx = -999.9
proj%dy = -999.9
proj%latinc = -999.9
proj%loninc = -999.9
proj%stdlon = -999.9
proj%truelat1 = -999.9
proj%truelat2 = -999.9
proj%phi = -999.9
proj%lambda = -999.9
proj%ixdim = -999
proj%jydim = -999
proj%stagger = HH
proj%nlat = 0
proj%nlon = 0
proj%nxmin = 1
proj%nxmax = 43200
proj%hemi = 0.0
proj%cone = -999.9
proj%polei = -999.9
proj%polej = -999.9
proj%rsw = -999.9
proj%knowni = -999.9
proj%knownj = -999.9
proj%re_m = EARTH_RADIUS_M
proj%init = .FALSE.
proj%wrap = .FALSE.
proj%rho0 = 0.
proj%nc = 0.
proj%bigc = 0.
proj%comp_ll = .FALSE.
nullify(proj%gauss_lat)
END SUBROUTINE map_init
SUBROUTINE map_set(proj_code, proj, lat1, lon1, lat0, lon0, knowni, knownj, dx, dy, latinc, &
loninc, stdlon, truelat1, truelat2, nlat, nlon, ixdim, jydim, nxmin, nxmax, &
stagger, phi, lambda, r_earth)
! Given a partially filled proj_info structure, this routine computes
! polei, polej, rsw, and cone (if LC projection) to complete the
! structure. This allows us to eliminate redundant calculations when
! calling the coordinate conversion routines multiple times for the
! same map.
! This will generally be the first routine called when a user wants
! to be able to use the coordinate conversion routines, and it
! will call the appropriate subroutines based on the
! proj%code which indicates which projection type this is.
IMPLICIT NONE
! Declare arguments
INTEGER, INTENT(IN) :: proj_code
INTEGER, INTENT(IN), OPTIONAL :: nlat
INTEGER, INTENT(IN), OPTIONAL :: nlon
INTEGER, INTENT(IN), OPTIONAL :: ixdim
INTEGER, INTENT(IN), OPTIONAL :: jydim
INTEGER, INTENT(IN), OPTIONAL :: nxmin
INTEGER, INTENT(IN), OPTIONAL :: nxmax
INTEGER, INTENT(IN), OPTIONAL :: stagger
REAL, INTENT(IN), OPTIONAL :: latinc
REAL, INTENT(IN), OPTIONAL :: loninc
REAL, INTENT(IN), OPTIONAL :: lat1
REAL, INTENT(IN), OPTIONAL :: lon1
REAL, INTENT(IN), OPTIONAL :: lat0
REAL, INTENT(IN), OPTIONAL :: lon0
REAL, INTENT(IN), OPTIONAL :: dx
REAL, INTENT(IN), OPTIONAL :: dy
REAL, INTENT(IN), OPTIONAL :: stdlon
REAL, INTENT(IN), OPTIONAL :: truelat1
REAL, INTENT(IN), OPTIONAL :: truelat2
REAL, INTENT(IN), OPTIONAL :: knowni
REAL, INTENT(IN), OPTIONAL :: knownj
REAL, INTENT(IN), OPTIONAL :: phi
REAL, INTENT(IN), OPTIONAL :: lambda
REAL, INTENT(IN), OPTIONAL :: r_earth
TYPE(proj_info), INTENT(OUT) :: proj
INTEGER :: iter
REAL :: dummy_lon1
REAL :: dummy_lon0
REAL :: dummy_stdlon
! First, verify that mandatory parameters are present for the specified proj_code
IF ( proj_code == PROJ_LC ) THEN
IF ( .NOT.PRESENT(truelat1) .OR. &
.NOT.PRESENT(truelat2) .OR. &
.NOT.PRESENT(lat1) .OR. &
.NOT.PRESENT(lon1) .OR. &
.NOT.PRESENT(knowni) .OR. &
.NOT.PRESENT(knownj) .OR. &
.NOT.PRESENT(stdlon) .OR. &
.NOT.PRESENT(dx) ) THEN
PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
PRINT '(A)', ' truelat1, truelat2, lat1, lon1, knowni, knownj, stdlon, dx'
call error_handler('MAP_INIT',ERROR_CODE)
END IF
ELSE IF ( proj_code == PROJ_PS ) THEN
IF ( .NOT.PRESENT(truelat1) .OR. &
.NOT.PRESENT(lat1) .OR. &
.NOT.PRESENT(lon1) .OR. &
.NOT.PRESENT(knowni) .OR. &
.NOT.PRESENT(knownj) .OR. &
.NOT.PRESENT(stdlon) .OR. &
.NOT.PRESENT(dx) ) THEN
PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
PRINT '(A)', ' truelat1, lat1, lon1, knonwi, knownj, stdlon, dx'
call error_handler('MAP_INIT',ERROR_CODE)
END IF
ELSE IF ( proj_code == PROJ_PS_WGS84 ) THEN
IF ( .NOT.PRESENT(truelat1) .OR. &
.NOT.PRESENT(lat1) .OR. &
.NOT.PRESENT(lon1) .OR. &
.NOT.PRESENT(knowni) .OR. &
.NOT.PRESENT(knownj) .OR. &
.NOT.PRESENT(stdlon) .OR. &
.NOT.PRESENT(dx) ) THEN
PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
PRINT '(A)', ' truelat1, lat1, lon1, knonwi, knownj, stdlon, dx'
call error_handler('MAP_INIT',ERROR_CODE)
END IF
ELSE IF ( proj_code == PROJ_ALBERS_NAD83 ) THEN
IF ( .NOT.PRESENT(truelat1) .OR. &
.NOT.PRESENT(truelat2) .OR. &
.NOT.PRESENT(lat1) .OR. &
.NOT.PRESENT(lon1) .OR. &
.NOT.PRESENT(knowni) .OR. &
.NOT.PRESENT(knownj) .OR. &
.NOT.PRESENT(stdlon) .OR. &
.NOT.PRESENT(dx) ) THEN
PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
PRINT '(A)', ' truelat1, truelat2, lat1, lon1, knonwi, knownj, stdlon, dx'
call error_handler('MAP_INIT',ERROR_CODE)
END IF
ELSE IF ( proj_code == PROJ_MERC ) THEN
IF ( .NOT.PRESENT(truelat1) .OR. &
.NOT.PRESENT(lat1) .OR. &
.NOT.PRESENT(lon1) .OR. &
.NOT.PRESENT(knowni) .OR. &
.NOT.PRESENT(knownj) .OR. &
.NOT.PRESENT(dx) ) THEN
PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
PRINT '(A)', ' truelat1, lat1, lon1, knowni, knownj, dx'
call error_handler('MAP_INIT',ERROR_CODE)
END IF
ELSE IF ( proj_code == PROJ_LATLON ) THEN
IF ( .NOT.PRESENT(latinc) .OR. &
.NOT.PRESENT(loninc) .OR. &
.NOT.PRESENT(knowni) .OR. &
.NOT.PRESENT(knownj) .OR. &
.NOT.PRESENT(lat1) .OR. &
.NOT.PRESENT(lon1) ) THEN
PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
PRINT '(A)', ' latinc, loninc, knowni, knownj, lat1, lon1'
call error_handler('MAP_INIT',ERROR_CODE)
END IF
ELSE IF ( proj_code == PROJ_CYL ) THEN
IF ( .NOT.PRESENT(latinc) .OR. &
.NOT.PRESENT(loninc) .OR. &
.NOT.PRESENT(stdlon) ) THEN
PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
PRINT '(A)', ' latinc, loninc, stdlon'
call error_handler('MAP_INIT',ERROR_CODE)
END IF
ELSE IF ( proj_code == PROJ_CASSINI ) THEN
IF ( .NOT.PRESENT(latinc) .OR. &
.NOT.PRESENT(loninc) .OR. &
.NOT.PRESENT(lat1) .OR. &
.NOT.PRESENT(lon1) .OR. &
.NOT.PRESENT(lat0) .OR. &
.NOT.PRESENT(lon0) .OR. &
.NOT.PRESENT(knowni) .OR. &
.NOT.PRESENT(knownj) .OR. &
.NOT.PRESENT(stdlon) ) THEN
PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
PRINT '(A)', ' latinc, loninc, lat1, lon1, knowni, knownj, lat0, lon0, stdlon'
call error_handler('MAP_INIT',ERROR_CODE)
END IF
ELSE IF ( proj_code == PROJ_GAUSS ) THEN
IF ( .NOT.PRESENT(nlat) .OR. &
.NOT.PRESENT(lat1) .OR. &
.NOT.PRESENT(lon1) .OR. &
.NOT.PRESENT(loninc) ) THEN
PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
PRINT '(A)', ' nlat, lat1, lon1, loninc'
call error_handler('MAP_INIT',ERROR_CODE)
END IF
ELSE IF ( proj_code == PROJ_ROTLL ) THEN
IF ( .NOT.PRESENT(ixdim) .OR. &
.NOT.PRESENT(jydim) .OR. &
.NOT.PRESENT(phi) .OR. &
.NOT.PRESENT(lambda) .OR. &
.NOT.PRESENT(lat1) .OR. &
.NOT.PRESENT(lon1) .OR. &
.NOT.PRESENT(stagger) ) THEN
PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
PRINT '(A)', ' ixdim, jydim, phi, lambda, lat1, lon1, stagger'
call error_handler('MAP_INIT',ERROR_CODE)
END IF
ELSE
PRINT '(A,I2)', 'Unknown projection code: ', proj_code
call error_handler('MAP_INIT',ERROR_CODE)
END IF
! Check for validity of mandatory variables in proj
IF ( PRESENT(lat1) ) THEN
IF ( ABS(lat1) .GT. 90. ) THEN
PRINT '(A)', 'Latitude of origin corner required as follows:'
PRINT '(A)', ' -90N <= lat1 < = 90.N'
call error_handler('MAP_INIT',ERROR_CODE)
ENDIF
ENDIF
IF ( PRESENT(lon1) ) THEN
dummy_lon1 = lon1
IF ( ABS(dummy_lon1) .GT. 180.) THEN
iter = 0
DO WHILE (ABS(dummy_lon1) > 180. .AND. iter < 10)
IF (dummy_lon1 < -180.) dummy_lon1 = dummy_lon1 + 360.
IF (dummy_lon1 > 180.) dummy_lon1 = dummy_lon1 - 360.
iter = iter + 1
END DO
IF (abs(dummy_lon1) > 180.) THEN
PRINT '(A)', 'Longitude of origin required as follows:'
PRINT '(A)', ' -180E <= lon1 <= 180W'
call error_handler('MAP_INIT',ERROR_CODE)
ENDIF
ENDIF
ENDIF
IF ( PRESENT(lon0) ) THEN
dummy_lon0 = lon0
IF ( ABS(dummy_lon0) .GT. 180.) THEN
iter = 0
DO WHILE (ABS(dummy_lon0) > 180. .AND. iter < 10)
IF (dummy_lon0 < -180.) dummy_lon0 = dummy_lon0 + 360.
IF (dummy_lon0 > 180.) dummy_lon0 = dummy_lon0 - 360.
iter = iter + 1
END DO
IF (abs(dummy_lon0) > 180.) THEN
PRINT '(A)', 'Longitude of pole required as follows:'
PRINT '(A)', ' -180E <= lon0 <= 180W'
call error_handler('MAP_INIT',ERROR_CODE)
ENDIF
ENDIF
ENDIF
IF ( PRESENT(dx) ) THEN
IF ((dx .LE. 0.).AND.(proj_code .NE. PROJ_LATLON)) THEN
PRINT '(A)', 'Require grid spacing (dx) in meters be positive!'
call error_handler('MAP_INIT',ERROR_CODE)
ENDIF
ENDIF
IF ( PRESENT(stdlon) ) THEN
dummy_stdlon = stdlon
IF ((ABS(dummy_stdlon) > 180.).AND.(proj_code /= PROJ_MERC)) THEN
iter = 0
DO WHILE (ABS(dummy_stdlon) > 180. .AND. iter < 10)
IF (dummy_stdlon < -180.) dummy_stdlon = dummy_stdlon + 360.
IF (dummy_stdlon > 180.) dummy_stdlon = dummy_stdlon - 360.
iter = iter + 1
END DO
IF (abs(dummy_stdlon) > 180.) THEN
PRINT '(A)', 'Need orientation longitude (stdlon) as: '
PRINT '(A)', ' -180E <= stdlon <= 180W'
call error_handler('MAP_INIT',ERROR_CODE)
ENDIF
ENDIF
ENDIF
IF ( PRESENT(truelat1) ) THEN
IF (ABS(truelat1).GT.90.) THEN
PRINT '(A)', 'Set true latitude 1 for all projections!'
call error_handler('MAP_INIT',ERROR_CODE)
ENDIF
ENDIF
CALL map_init(proj)
proj%code = proj_code
IF ( PRESENT(lat1) ) proj%lat1 = lat1
IF ( PRESENT(lon1) ) proj%lon1 = dummy_lon1
IF ( PRESENT(lat0) ) proj%lat0 = lat0
IF ( PRESENT(lon0) ) proj%lon0 = dummy_lon0
IF ( PRESENT(latinc) ) proj%latinc = latinc
IF ( PRESENT(loninc) ) proj%loninc = loninc
IF ( PRESENT(knowni) ) proj%knowni = knowni
IF ( PRESENT(knownj) ) proj%knownj = knownj
IF ( PRESENT(nxmin) ) proj%nxmin = nxmin
IF ( PRESENT(nxmax) ) proj%nxmax = nxmax
IF ( PRESENT(dx) ) proj%dx = dx
IF ( PRESENT(dy) ) THEN
proj%dy = dy
ELSE IF ( PRESENT(dx) ) THEN
proj%dy = dx
END IF
IF ( PRESENT(stdlon) ) proj%stdlon = dummy_stdlon
IF ( PRESENT(truelat1) ) proj%truelat1 = truelat1
IF ( PRESENT(truelat2) ) proj%truelat2 = truelat2
IF ( PRESENT(nlat) ) proj%nlat = nlat
IF ( PRESENT(nlon) ) proj%nlon = nlon
IF ( PRESENT(ixdim) ) proj%ixdim = ixdim
IF ( PRESENT(jydim) ) proj%jydim = jydim
IF ( PRESENT(stagger) ) proj%stagger = stagger
IF ( PRESENT(phi) ) proj%phi = phi
IF ( PRESENT(lambda) ) proj%lambda = lambda
IF ( PRESENT(r_earth) ) proj%re_m = r_earth
IF ( PRESENT(dx) ) THEN
IF ( (proj_code == PROJ_LC) .OR. (proj_code == PROJ_PS) .OR. &
(proj_code == PROJ_PS_WGS84) .OR. (proj_code == PROJ_ALBERS_NAD83) .OR. &
(proj_code == PROJ_MERC) ) THEN
proj%dx = dx
IF (truelat1 .LT. 0.) THEN
proj%hemi = -1.0
ELSE
proj%hemi = 1.0
ENDIF
proj%rebydx = proj%re_m / dx
ENDIF
ENDIF
pick_proj: SELECT CASE(proj%code)
CASE(PROJ_PS)
CALL set_ps(proj)
CASE(PROJ_PS_WGS84)
CALL set_ps_wgs84(proj)
CASE(PROJ_ALBERS_NAD83)
CALL set_albers_nad83(proj)
CASE(PROJ_LC)
IF (ABS(proj%truelat2) .GT. 90.) THEN
proj%truelat2=proj%truelat1
ENDIF
CALL set_lc(proj)
CASE (PROJ_MERC)
CALL set_merc(proj)
CASE (PROJ_LATLON)
CASE (PROJ_GAUSS)
CALL set_gauss(proj)
CASE (PROJ_CYL)
CALL set_cyl(proj)
CASE (PROJ_CASSINI)
CALL set_cassini(proj)
CASE (PROJ_ROTLL)
END SELECT pick_proj
proj%init = .TRUE.
RETURN
END SUBROUTINE map_set
SUBROUTINE latlon_to_ij(proj, lat, lon, i, j)
! Converts input lat/lon values to the cartesian (i,j) value
! for the given projection.
IMPLICIT NONE
TYPE(proj_info), INTENT(IN) :: proj
REAL, INTENT(IN) :: lat
REAL, INTENT(IN) :: lon
REAL, INTENT(OUT) :: i
REAL, INTENT(OUT) :: j
IF (.NOT.proj%init) THEN
PRINT '(A)', 'You have not called map_set for this projection!'
call error_handler('LATLON_TO_IJ',ERROR_CODE)
ENDIF
SELECT CASE(proj%code)
CASE(PROJ_LATLON)
CALL llij_latlon(lat,lon,proj,i,j)
CASE(PROJ_MERC)
CALL llij_merc(lat,lon,proj,i,j)
CASE(PROJ_PS)
CALL llij_ps(lat,lon,proj,i,j)
CASE(PROJ_PS_WGS84)
CALL llij_ps_wgs84(lat,lon,proj,i,j)
CASE(PROJ_ALBERS_NAD83)
CALL llij_albers_nad83(lat,lon,proj,i,j)
CASE(PROJ_LC)
CALL llij_lc(lat,lon,proj,i,j)
CASE(PROJ_GAUSS)
CALL llij_gauss(lat,lon,proj,i,j)
CASE(PROJ_CYL)
CALL llij_cyl(lat,lon,proj,i,j)
CASE(PROJ_CASSINI)
CALL llij_cassini(lat,lon,proj,i,j)
CASE(PROJ_ROTLL)
CALL llij_rotlatlon(lat,lon,proj,i,j)
CASE DEFAULT
PRINT '(A,I2)', 'Unrecognized map projection code: ', proj%code
call error_handler('LATLON_TO_IJ',ERROR_CODE)
END SELECT
RETURN
END SUBROUTINE latlon_to_ij
SUBROUTINE ij_to_latlon(proj, i, j, lat, lon)
! Computes geographical latitude and longitude for a given (i,j) point
! in a grid with a projection of proj
IMPLICIT NONE
TYPE(proj_info),INTENT(IN) :: proj
REAL, INTENT(IN) :: i
REAL, INTENT(IN) :: j
REAL, INTENT(OUT) :: lat
REAL, INTENT(OUT) :: lon
IF (.NOT.proj%init) THEN
PRINT '(A)', 'You have not called map_set for this projection!'
call error_handler('IJ_TO_LATLON',ERROR_CODE)
ENDIF
SELECT CASE (proj%code)
CASE (PROJ_LATLON)
CALL ijll_latlon(i, j, proj, lat, lon)
CASE (PROJ_MERC)
CALL ijll_merc(i, j, proj, lat, lon)
CASE (PROJ_PS)
CALL ijll_ps(i, j, proj, lat, lon)
CASE (PROJ_PS_WGS84)
CALL ijll_ps_wgs84(i, j, proj, lat, lon)
CASE (PROJ_ALBERS_NAD83)
CALL ijll_albers_nad83(i, j, proj, lat, lon)
CASE (PROJ_LC)
CALL ijll_lc(i, j, proj, lat, lon)
CASE (PROJ_CYL)
CALL ijll_cyl(i, j, proj, lat, lon)
CASE (PROJ_CASSINI)
CALL ijll_cassini(i, j, proj, lat, lon)
CASE (PROJ_ROTLL)
CALL ijll_rotlatlon(i, j, proj, lat, lon)
CASE DEFAULT
PRINT '(A,I2)', 'Unrecognized map projection code: ', proj%code
call error_handler('IJ_TO_LATLON',ERROR_CODE)
END SELECT
RETURN
END SUBROUTINE ij_to_latlon
SUBROUTINE set_ps(proj)
! Initializes a polar-stereographic map projection from the partially
! filled proj structure. This routine computes the radius to the
! southwest corner and computes the i/j location of the pole for use
! in llij_ps and ijll_ps.
IMPLICIT NONE
! Declare args
TYPE(proj_info), INTENT(INOUT) :: proj
! Local vars
REAL :: ala1
REAL :: alo1
REAL :: reflon
REAL :: scale_top
! Executable code
reflon = proj%stdlon + 90.
! Compute numerator term of map scale factor
scale_top = 1. + proj%hemi * SIN(proj%truelat1 * rad_per_deg)
! Compute radius to lower-left (SW) corner
ala1 = proj%lat1 * rad_per_deg
proj%rsw = proj%rebydx*COS(ala1)*scale_top/(1.+proj%hemi*SIN(ala1))
! Find the pole point
alo1 = (proj%lon1 - reflon) * rad_per_deg
proj%polei = proj%knowni - proj%rsw * COS(alo1)
proj%polej = proj%knownj - proj%hemi * proj%rsw * SIN(alo1)
RETURN
END SUBROUTINE set_ps
SUBROUTINE llij_ps(lat,lon,proj,i,j)
! Given latitude (-90 to 90), longitude (-180 to 180), and the
! standard polar-stereographic projection information via the
! public proj structure, this routine returns the i/j indices which
! if within the domain range from 1->nx and 1->ny, respectively.
IMPLICIT NONE
! Delcare input arguments
REAL, INTENT(IN) :: lat
REAL, INTENT(IN) :: lon
TYPE(proj_info),INTENT(IN) :: proj
! Declare output arguments
REAL, INTENT(OUT) :: i !(x-index)
REAL, INTENT(OUT) :: j !(y-index)
! Declare local variables
REAL :: reflon
REAL :: scale_top
REAL :: ala
REAL :: alo
REAL :: rm
! BEGIN CODE
reflon = proj%stdlon + 90.
! Compute numerator term of map scale factor
scale_top = 1. + proj%hemi * SIN(proj%truelat1 * rad_per_deg)
! Find radius to desired point
ala = lat * rad_per_deg
rm = proj%rebydx * COS(ala) * scale_top/(1. + proj%hemi *SIN(ala))
alo = (lon - reflon) * rad_per_deg
i = proj%polei + rm * COS(alo)
j = proj%polej + proj%hemi * rm * SIN(alo)
RETURN
END SUBROUTINE llij_ps
SUBROUTINE ijll_ps(i, j, proj, lat, lon)
! This is the inverse subroutine of llij_ps. It returns the
! latitude and longitude of an i/j point given the projection info
! structure.
IMPLICIT NONE
! Declare input arguments
REAL, INTENT(IN) :: i ! Column
REAL, INTENT(IN) :: j ! Row
TYPE (proj_info), INTENT(IN) :: proj
! Declare output arguments
REAL, INTENT(OUT) :: lat ! -90 -> 90 north
REAL, INTENT(OUT) :: lon ! -180 -> 180 East
! Local variables
REAL :: reflon
REAL :: scale_top
REAL :: xx,yy
REAL(KIND=HIGH) :: gi2, r2
REAL :: arccos
! Begin Code
! Compute the reference longitude by rotating 90 degrees to the east
! to find the longitude line parallel to the positive x-axis.
reflon = proj%stdlon + 90.
! Compute numerator term of map scale factor
scale_top = 1. + proj%hemi * SIN(proj%truelat1 * rad_per_deg)
! Compute radius to point of interest
xx = i - proj%polei
yy = (j - proj%polej) * proj%hemi
r2 = xx**2 + yy**2
! Now the magic code
IF (r2 .EQ. 0.) THEN
lat = proj%hemi * 90.
lon = reflon
ELSE
gi2 = (proj%rebydx * scale_top)**2.
lat = deg_per_rad * proj%hemi * ASIN((gi2-r2)/(gi2+r2))
arccos = ACOS(MIN(MAX(xx/SQRT(r2),-1.0_HIGH),1.0_HIGH))
IF (yy .GT. 0) THEN
lon = reflon + deg_per_rad * arccos
ELSE
lon = reflon - deg_per_rad * arccos
ENDIF
ENDIF
! Convert to a -180 -> 180 East convention
IF (lon .GT. 180.) lon = lon - 360.
IF (lon .LT. -180.) lon = lon + 360.
RETURN
END SUBROUTINE ijll_ps
SUBROUTINE set_ps_wgs84(proj)
! Initializes a polar-stereographic map projection (WGS84 ellipsoid)
! from the partially filled proj structure. This routine computes the
! radius to the southwest corner and computes the i/j location of the
! pole for use in llij_ps and ijll_ps.
IMPLICIT NONE
! Arguments
TYPE(proj_info), INTENT(INOUT) :: proj
! Local variables
real :: h, mc, tc, t, rho
h = proj%hemi
mc = cos(h*proj%truelat1*rad_per_deg)/sqrt(1.0-(E_WGS84*sin(h*proj%truelat1*rad_per_deg))**2.0)
tc = sqrt(((1.0-sin(h*proj%truelat1*rad_per_deg))/(1.0+sin(h*proj%truelat1*rad_per_deg)))* &
(((1.0+E_WGS84*sin(h*proj%truelat1*rad_per_deg))/(1.0-E_WGS84*sin(h*proj%truelat1*rad_per_deg)))**E_WGS84 ))
! Find the i/j location of reference lat/lon with respect to the pole of the projection
t = sqrt(((1.0-sin(h*proj%lat1*rad_per_deg))/(1.0+sin(h*proj%lat1*rad_per_deg)))* &
(((1.0+E_WGS84*sin(h*proj%lat1*rad_per_deg))/(1.0-E_WGS84*sin(h*proj%lat1*rad_per_deg)) )**E_WGS84 ) )
rho = h * (A_WGS84 / proj%dx) * mc * t / tc
proj%polei = rho * sin((h*proj%lon1 - h*proj%stdlon)*rad_per_deg)
proj%polej = -rho * cos((h*proj%lon1 - h*proj%stdlon)*rad_per_deg)
RETURN
END SUBROUTINE set_ps_wgs84
SUBROUTINE llij_ps_wgs84(lat,lon,proj,i,j)
! Given latitude (-90 to 90), longitude (-180 to 180), and the
! standard polar-stereographic projection information via the
! public proj structure, this routine returns the i/j indices which
! if within the domain range from 1->nx and 1->ny, respectively.
IMPLICIT NONE
! Arguments
REAL, INTENT(IN) :: lat
REAL, INTENT(IN) :: lon
REAL, INTENT(OUT) :: i !(x-index)
REAL, INTENT(OUT) :: j !(y-index)
TYPE(proj_info),INTENT(IN) :: proj
! Local variables
real :: h, mc, tc, t, rho
h = proj%hemi
mc = cos(h*proj%truelat1*rad_per_deg)/sqrt(1.0-(E_WGS84*sin(h*proj%truelat1*rad_per_deg))**2.0)
tc = sqrt(((1.0-sin(h*proj%truelat1*rad_per_deg))/(1.0+sin(h*proj%truelat1*rad_per_deg)))* &
(((1.0+E_WGS84*sin(h*proj%truelat1*rad_per_deg))/(1.0-E_WGS84*sin(h*proj%truelat1*rad_per_deg)))**E_WGS84 ))
t = sqrt(((1.0-sin(h*lat*rad_per_deg))/(1.0+sin(h*lat*rad_per_deg))) * &
(((1.0+E_WGS84*sin(h*lat*rad_per_deg))/(1.0-E_WGS84*sin(h*lat*rad_per_deg)))**E_WGS84))
! Find the x/y location of the requested lat/lon with respect to the pole of the projection
rho = (A_WGS84 / proj%dx) * mc * t / tc
i = h * rho * sin((h*lon - h*proj%stdlon)*rad_per_deg)
j = h *(-rho)* cos((h*lon - h*proj%stdlon)*rad_per_deg)
! Get i/j relative to reference i/j
i = proj%knowni + (i - proj%polei)
j = proj%knownj + (j - proj%polej)
RETURN
END SUBROUTINE llij_ps_wgs84
SUBROUTINE ijll_ps_wgs84(i, j, proj, lat, lon)
! This is the inverse subroutine of llij_ps. It returns the
! latitude and longitude of an i/j point given the projection info
! structure.
implicit none
! Arguments
REAL, INTENT(IN) :: i ! Column
REAL, INTENT(IN) :: j ! Row
REAL, INTENT(OUT) :: lat ! -90 -> 90 north
REAL, INTENT(OUT) :: lon ! -180 -> 180 East
TYPE (proj_info), INTENT(IN) :: proj
! Local variables
real :: h, mc, tc, t, rho, x, y
real :: chi, a, b, c, d
h = proj%hemi
x = (i - proj%knowni + proj%polei)
y = (j - proj%knownj + proj%polej)
mc = cos(h*proj%truelat1*rad_per_deg)/sqrt(1.0-(E_WGS84*sin(h*proj%truelat1*rad_per_deg))**2.0)
tc = sqrt(((1.0-sin(h*proj%truelat1*rad_per_deg))/(1.0+sin(h*proj%truelat1*rad_per_deg))) * &
(((1.0+E_WGS84*sin(h*proj%truelat1*rad_per_deg))/(1.0-E_WGS84*sin(h*proj%truelat1*rad_per_deg)))**E_WGS84 ))
rho = sqrt((x*proj%dx)**2.0 + (y*proj%dx)**2.0)
t = rho * tc / (A_WGS84 * mc)
lon = h*proj%stdlon*rad_per_deg + h*atan2(h*x,h*(-y))
chi = PI/2.0-2.0*atan(t)
a = 1./2.*E_WGS84**2. + 5./24.*E_WGS84**4. + 1./40.*E_WGS84**6. + 73./2016.*E_WGS84**8.
b = 7./24.*E_WGS84**4. + 29./120.*E_WGS84**6. + 54113./40320.*E_WGS84**8.
c = 7./30.*E_WGS84**6. + 81./280.*E_WGS84**8.
d = 4279./20160.*E_WGS84**8.
lat = chi + sin(2.*chi)*(a + cos(2.*chi)*(b + cos(2.*chi)*(c + d*cos(2.*chi))))
lat = h * lat
lat = lat*deg_per_rad
lon = lon*deg_per_rad
RETURN
END SUBROUTINE ijll_ps_wgs84
SUBROUTINE set_albers_nad83(proj)
! Initializes an Albers equal area map projection (NAD83 ellipsoid)
! from the partially filled proj structure. This routine computes the
! radius to the southwest corner and computes the i/j location of the
! pole for use in llij_albers_nad83 and ijll_albers_nad83.
IMPLICIT NONE
! Arguments
TYPE(proj_info), INTENT(INOUT) :: proj
! Local variables
real :: h, m1, m2, q1, q2, theta, q, sinphi
h = proj%hemi
m1 = cos(h*proj%truelat1*rad_per_deg)/sqrt(1.0-(E_NAD83*sin(h*proj%truelat1*rad_per_deg))**2.0)
m2 = cos(h*proj%truelat2*rad_per_deg)/sqrt(1.0-(E_NAD83*sin(h*proj%truelat2*rad_per_deg))**2.0)
sinphi = sin(proj%truelat1*rad_per_deg)
q1 = (1.0-E_NAD83**2.0) * &
((sinphi/(1.0-(E_NAD83*sinphi)**2.0)) - 1.0/(2.0*E_NAD83) * log((1.0-E_NAD83*sinphi)/(1.0+E_NAD83*sinphi)))
sinphi = sin(proj%truelat2*rad_per_deg)
q2 = (1.0-E_NAD83**2.0) * &
((sinphi/(1.0-(E_NAD83*sinphi)**2.0)) - 1.0/(2.0*E_NAD83) * log((1.0-E_NAD83*sinphi)/(1.0+E_NAD83*sinphi)))
if (proj%truelat1 == proj%truelat2) then
proj%nc = sin(proj%truelat1*rad_per_deg)
else
proj%nc = (m1**2.0 - m2**2.0) / (q2 - q1)
end if
proj%bigc = m1**2.0 + proj%nc*q1
! Find the i/j location of reference lat/lon with respect to the pole of the projection
sinphi = sin(proj%lat1*rad_per_deg)
q = (1.0-E_NAD83**2.0) * &
((sinphi/(1.0-(E_NAD83*sinphi)**2.0)) - 1.0/(2.0*E_NAD83) * log((1.0-E_NAD83*sinphi)/(1.0+E_NAD83*sinphi)))
proj%rho0 = h * (A_NAD83 / proj%dx) * sqrt(proj%bigc - proj%nc * q) / proj%nc
theta = proj%nc*(proj%lon1 - proj%stdlon)*rad_per_deg
proj%polei = proj%rho0 * sin(h*theta)
proj%polej = proj%rho0 - proj%rho0 * cos(h*theta)
RETURN
END SUBROUTINE set_albers_nad83
SUBROUTINE llij_albers_nad83(lat,lon,proj,i,j)
! Given latitude (-90 to 90), longitude (-180 to 180), and the
! standard projection information via the