-
Notifications
You must be signed in to change notification settings - Fork 571
/
Copy pathwmscrpmd.F90
1793 lines (1700 loc) · 60.6 KB
/
wmscrpmd.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
!> @brief Contains module WMSCRPMD.
!>
!> @author E. Rogers
!> @author M. Dutour
!> @author A. Roland
!> @author F. Ardhuin
!> @date 10-Dec-2014
!>
#include "w3macros.h"
!/ ------------------------------------------------------------------- /
!>
!> @brief Routines to determine and process grid dependencies in the
!> multi-grid wave model.
!>
!> @author E. Rogers
!> @author M. Dutour
!> @author A. Roland
!> @author F. Ardhuin
!> @date 10-Dec-2014
!>
MODULE WMSCRPMD
!/
!/ +-----------------------------------+
!/ | WAVEWATCH III |
!/ | E. Rogers, M. Dutour, |
!/ | A. Roland, F. Ardhuin |
!/ | FORTRAN 90 |
!/ | Last update : 10-Dec-2014 |
!/ +-----------------------------------+
!/
!/ 06-Sep-2012 : Origination, transfer from WMGRIDMD ( version 4.08 )
!/ 10-Dec-2014 : Add checks for allocate status ( version 5.04 )
!/
!/ Copyright 2012 National Weather Service (NWS),
!/ National Oceanic and Atmospheric Administration. All rights
!/ reserved. WAVEWATCH III is a trademark of the NWS.
!/ No unauthorized use without permission.
!/
! 1. Purpose :
!
! Routines to determine and process grid dependencies in the
! multi-grid wave model.
!
! 2. Variables and types :
!
! 3. Subroutines and functions :
!
! Name Type Scope Description
! ----------------------------------------------------------------
! scrip_wrapper Subr. Public as the name says ...
! get_scrip_info_structured Subr. Public as the name says ...
! get_scrip_info_unstructured Subr. Public as the name says ...
! get_scrip_info Subr. Public as the name says ...
! scrip_info_renormalization Subr. Public as the name says ...
! TRIANG_INDEXES Subr. Public as the name says ...
! get_unstructured_vertex_degree Subr. Public as the name says ...
! GET_BOUNDARY Subr. Public as the name says ...
! ----------------------------------------------------------------
!
! 4. Subroutines and functions used :
!
! Name Type Module Description
! ----------------------------------------------------------------
! get_unstructured_vertex_degree Subr. W3TRIAMD Manage data structures
! ----------------------------------------------------------------
!
! 5. Remarks :
!
! - The subroutines TRIANG_INDEXES, get_unstructured_vertex_degree, and
! GET_BOUNDARY should probably be renamed and moved to the module w3triamd
!
! 6. Switches :
!
!
! !/S Enable subroutine tracing.
! !/Tn Enable test output.
!
! 7. Source code :
!
!/ ------------------------------------------------------------------- /
!/
!/ Specify default accessibility
!/
PUBLIC
!/
!/ Module private variable for checking error returns
!/
INTEGER, PRIVATE :: ISTAT !< ISTAT Check error returns
!/
CONTAINS
!/ ------------------------------------------------------------------- /
!>
!> @brief Compute grid information required by SCRIP.
!>
!> @param[in] ID_SRC
!> @param[in] ID_DST
!> @param[in] MAPSTA_SRC
!> @param[in] MAPST2_SRC
!> @param[in] FLAGLL
!> @param[in] GRIDSHIFT
!> @param[in] L_MASTER
!> @param[in] L_READ
!> @param[in] L_TEST
!>
!> @author E. Rogers
!> @author M. Dutour
!> @author A. Roland
!> @date 10-Dec-2014
!>
SUBROUTINE SCRIP_WRAPPER (ID_SRC, ID_DST, &
MAPSTA_SRC,MAPST2_SRC,FLAGLL,GRIDSHIFT,L_MASTER, &
L_READ,L_TEST)
!/ +-----------------------------------+
!/ | WAVEWATCH III |
!/ | E. Rogers, M. Dutour, A. Roland |
!/ | FORTRAN 90 |
!/ | Last update : 10-Dec-2014 !
!/ +-----------------------------------+
!/
! 1. Original author :
!
! Erick Rogers, NRL
!
! 2. Last update :
!
! See revisions.
!
! 3. Revisions :
!
! 29-Apr-2011 : Origination ( version 4.01 )
! 20-Feb-2012 : Mathieu Dutour Sikiric, Aron Roland Z&P
! Modification is to split the code into several
! subroutines
! ---get_scrip_info_structured (the structured case)
! ---get_scrip_info (chooses according to FD/FE)
! ---scrip_info_renormalization (conv_x and all that)
! 11-Apr-2013 Kevin Lind
! Added code for reading/writing SCRIP remap files
!/ 10-Dec-2014 : Add checks for allocate status ( version 5.04 )
!
! 4. Copyright :
!
! 5. Purpose :
!
! Compute grid information required by SCRIP
!
! 6. Method :
!
! 7. Parameters, Variables and types :
!
! 8. Called by :
!
! Subroutine WMGHGH
!
! 9. Subroutines and functions used :
!
! Subroutine SCRIP
!
! 10. Error messages:
!
! 11. Remarks :
!
! 12. Structure :
!
! 13. Switches :
!
! 14. Source code :
USE SCRIP_GRIDS
USE SCRIP_REMAP_VARS
USE SCRIP_CONSTANTS
USE SCRIP_KINDSMOD
USE SCRIP_INTERFACE
USE W3SERVMD, ONLY: EXTCDE
! USE W3GDATMD, ONLY : GRIDS
IMPLICIT NONE
INTEGER(SCRIP_I4), INTENT(IN) :: ID_SRC, ID_DST
INTEGER(SCRIP_I4), INTENT(IN) :: MAPSTA_SRC(:,:)
INTEGER(SCRIP_I4), INTENT(IN) :: MAPST2_SRC(:,:)
LOGICAL(SCRIP_LOGICAL), INTENT(IN) :: FLAGLL
REAL (SCRIP_R8), INTENT(IN) :: GRIDSHIFT
LOGICAL(SCRIP_LOGICAL), INTENT(IN) :: L_MASTER ! Am I the master processor (do I/O)?
LOGICAL(SCRIP_LOGICAL), INTENT(IN) :: L_READ ! Do I read the remap file?
LOGICAL(SCRIP_LOGICAL), INTENT(IN) :: L_TEST ! Whether to include test output
! in subroutines
!/ ------------------------------------------------------------------- /
!/ local variables
!/
INTEGER(SCRIP_I4) :: IREC,I,J,NI,NJ,IDUM,NK,K, &
ILINK,IW,ICORNER,NGOODPNTS, &
NBADPNTS
INTEGER(SCRIP_I4) :: ISRC,JSRC,KSRC,IPNT,KDST, &
NI_SRC
REAL (SCRIP_R8) :: LAT_CONVERSION,OFFSET
REAL (SCRIP_R8) :: CONV_DX,CONV_DY,WEIGHT
REAL (SCRIP_R8) :: WTSUM
#ifdef W3_T38
CHARACTER (LEN=10) :: CDATE_TIME(3)
INTEGER :: DATE_TIME(8)
INTEGER :: ELAPSED_TIME, BEG_TIME, &
END_TIME
#endif
! test output for input variables
#ifdef W3_T38
if(l_master)write(*,*)'flagll = ',flagll
if(l_master)write(*,*)'gridshift = ',gridshift
#endif
!
! START: universal settings
!
! Set variables for converting to degrees
! notes: SCRIP only operates on spherical coordinates, so for the case
! where the problem is specified by the user as in a
! meters/cartesian coordinate system, it is necessary to make
! a temporary conversion to a "fake" spherical coordinate grid,
! to keep SCRIP happy. The good news here is that multi-grid
! meters-grid simulations will be very rare: we will probably only
! encounter them in the context of simple test cases. Strictly
! speaking, this conversion does not even need to be physically
! correct, e.g. we could say that 1000 km is 1 deg....as long as
! we are consistent between grids.
! Potential future improvement: make conv_dy and offset such that dst grid
! covers a specific longitude range, e.g. 1 deg east to 2 deg east
!
! START: set up src grid
!
!notes: when we work out how to interface with an unstructured grid,
! we will need to revisit this issue of how to set grid1_rank, etc.
! strategy: declare variables in grid module, but allocate them here.
!
GRID1_UNITS='degrees' ! the other option is radians...we don't use this
GRID1_NAME='src' ! this is not used, except for netcdf output
CALL GET_SCRIP_INFO(ID_SRC, &
& GRID1_CENTER_LON, GRID1_CENTER_LAT, &
& GRID1_CORNER_LON, GRID1_CORNER_LAT, GRID1_MASK, &
& GRID1_DIMS, GRID1_SIZE, GRID1_CORNERS, GRID1_RANK)
GRID2_UNITS='degrees'
GRID2_NAME='dst'
CALL GET_SCRIP_INFO(ID_DST, &
& GRID2_CENTER_LON, GRID2_CENTER_LAT, &
& GRID2_CORNER_LON, GRID2_CORNER_LAT, GRID2_MASK, &
& GRID2_DIMS, GRID2_SIZE, GRID2_CORNERS, GRID2_RANK)
IF(FLAGLL)THEN
CONV_DX=ONE
CONV_DY=ONE
OFFSET=ZERO
ELSE
LAT_CONVERSION=ZERO ! lat_conversion
! notes: this is the latitude used for conversion everywhere
! in the grid (approximation)
! (in radians)
! conv_dy=92.6*1200.0 ! physical, =92.6/(3/3600)=111000 m = 111 km
CONV_DY=1.0E+6_SCRIP_R8 ! non-physical, 1e+6=1 deg
CONV_DX=COS(LAT_CONVERSION)*CONV_DY
! notes: offset (in meters), is necessary so that our grid does
! not lie on the branch cut
OFFSET=75000.0_SCRIP_R8-MIN(MINVAL(GRID1_CENTER_LON), &
& MINVAL(GRID2_CENTER_LON))
ENDIF
!.....test output
#ifdef W3_T38
write(*,*)'l_master = ',l_master
if(l_master)then
write(*,*)'conv_dx=', conv_dx
write(*,*)'conv_dy=', conv_dy
write(*,*)'offset = ',offset
write(*,*)'grid1_size=', grid1_size
write(*,*)'grid2_size=', grid2_size
write(*,*)'l_read = ',l_read
write(*,*)'minval(grid1_center_lon) = ',minval(grid1_center_lon)
write(*,*)'maxval(grid1_center_lon) = ',maxval(grid1_center_lon)
write(*,*)'minval(grid1_center_lat) = ',minval(grid1_center_lat)
write(*,*)'maxval(grid1_center_lat) = ',maxval(grid1_center_lat)
write(*,*)'minval(grid2_center_lon) = ',minval(grid2_center_lon)
write(*,*)'maxval(grid2_center_lon) = ',maxval(grid2_center_lon)
write(*,*)'minval(grid2_center_lat) = ',minval(grid2_center_lat)
write(*,*)'maxval(grid2_center_lat) = ',maxval(grid2_center_lat)
endif
#endif
CALL SCRIP_INFO_RENORMALIZATION( &
& GRID1_CENTER_LON, GRID1_CENTER_LAT, &
& GRID1_CORNER_LON, GRID1_CORNER_LAT, GRID1_MASK, &
& GRID1_DIMS, GRID1_SIZE, GRID1_CORNERS, GRID1_RANK, &
& CONV_DX, CONV_DY, OFFSET, GRIDSHIFT)
CALL SCRIP_INFO_RENORMALIZATION( &
& GRID2_CENTER_LON, GRID2_CENTER_LAT, &
& GRID2_CORNER_LON, GRID2_CORNER_LAT, GRID2_MASK, &
& GRID2_DIMS, GRID2_SIZE, GRID2_CORNERS, GRID2_RANK, &
& CONV_DX, CONV_DY, OFFSET, ZERO)
!.....Set constants for thresholding weights:
FRAC_LOWEST =1.E-3_SCRIP_R8
FRAC_HIGHEST=ONE+1.E-3_SCRIP_R8
WT_LOWEST =ZERO
WT_HIGHEST =ONE+1.E-7_SCRIP_R8
#ifdef W3_T38
call date_and_time (CDATE_TIME(1), CDATE_TIME(2), CDATE_TIME(3), DATE_TIME)
beg_time = ((date_time(5)*60 + date_time(6))*60 +date_time(7))*1000 + date_time(8)
#endif
CALL SCRIP(ID_SRC, ID_DST, L_MASTER, L_READ, L_TEST)
#ifdef W3_T38
call date_and_time (CDATE_TIME(1), CDATE_TIME(2), CDATE_TIME(3), DATE_TIME)
end_time = ((date_time(5)*60 + date_time(6))*60 +date_time(7))*1000 + date_time(8)
elapsed_time = end_time - beg_time
write(0,*) "SCRIP: ", elapsed_time, " MSEC"
#endif
#ifdef W3_T38
if(l_master)write(*,*)'new minval(grid1_center_lon) = ',minval(grid1_center_lon)
if(l_master)write(*,*)'new maxval(grid1_center_lon) = ',maxval(grid1_center_lon)
#endif
!.....notes: at this point we have the following useful variables:
! num_wts, e.g. num_wts=3....for first order conservative remapping,
! only the first one is relevant, the other two are for second-order
! remapping.
! max_links_map1, e.g. max_links_map1=1369,
! grid2_size, e.g. grid2_size=1849,
! wts_map1(num_wts,max_links_map1), e.g. wts_map1(3,1369),
! grid1_add_map1(max_links_map1), e.g. grid1_add_map1(1369),
! grid2_add_map1(max_links_map1), e.g. grid2_add_map1(1369),
! grid2_frac(grid2_size), e.g. grid2_frac(1849),
! (see earlier versions for notes re: equivalency in netcdf/matlab)
!
!.....test output (optional)
!
!.....note re: notation: I use k for the combined i/j array, similar to isea,
! but not necessarily the same as isea since some points may
! be land etc.
#ifdef W3_T38
if(l_master)then
do k=1,grid2_size
write(403,*)grid2_frac(k)
end do
do ilink=1,max_links_map1
write(405,'(999(1x,f20.7))')(wts_map1(iw,ilink),iw=1,num_wts)
end do
do ilink=1,max_links_map1
write(406,'(i20)')grid1_add_map1(ilink) ! equivalent to
! my "src_address"
write(407,'(i20)')grid2_add_map1(ilink) ! equivalent to
! my "dst_address"
end do
endif
#endif
!.....organize results and return to wmghgh.
! what we need, for purpose of feeding back to ww3, for each dst grid node:
! a) the set of src grid nodes, in terms of isea, for which
! weights are available
! b) the corresponding set of weights
ALLOCATE ( WGTDATA(GRID2_SIZE), STAT=ISTAT )
CHECK_ALLOC_STATUS ( ISTAT )
!.....step 1: count up NR0, NR1, NR2, NRL, NLOC (NR1 and NLOC are denoted
! "n" here)
! It is especially important to determine how large npnts gets,
! so that we can allocate arrays properly
WGTDATA%NR0=0
WGTDATA%NR2=0
WGTDATA%NRL=0
WGTDATA%N=0
NI_SRC=GRID1_DIMS(1)
NGOODPNTS=0
NBADPNTS=0
DO ILINK=1,MAX_LINKS_MAP1
!........note that this pair of if-thens *must* be consistent with the
!........single if-then below
IF ((GRID2_FRAC(GRID2_ADD_MAP1(ILINK))>FRAC_LOWEST) .AND. &
(GRID2_FRAC(GRID2_ADD_MAP1(ILINK))<FRAC_HIGHEST)) THEN
! then consider this link either to include, or to print a warning
IF( (WTS_MAP1(1,ILINK)>=WT_LOWEST) .AND. &
(WTS_MAP1(1,ILINK)<=WT_HIGHEST) ) THEN
KSRC=GRID1_ADD_MAP1(ILINK)
JSRC=INT((KSRC-1)/NI_SRC)+1
ISRC=KSRC-(JSRC-1)*NI_SRC
IF (MAPSTA_SRC(JSRC,ISRC).EQ.0) THEN ! excluded point
WGTDATA(GRID2_ADD_MAP1(ILINK))%NR0 = &
WGTDATA(GRID2_ADD_MAP1(ILINK))%NR0 + 1
IF (MAPST2_SRC(JSRC,ISRC).EQ.0)THEN
WGTDATA(GRID2_ADD_MAP1(ILINK))%NRL = &
WGTDATA(GRID2_ADD_MAP1(ILINK))%NRL + 1
ENDIF
ELSE IF (ABS(MAPSTA_SRC(JSRC,ISRC)).EQ.1) THEN
! sea point
WGTDATA(GRID2_ADD_MAP1(ILINK))%N= &
WGTDATA(GRID2_ADD_MAP1(ILINK))%N+1
ELSE IF (ABS(MAPSTA_SRC(JSRC,ISRC)).EQ.2) THEN
! bnd point
WGTDATA(GRID2_ADD_MAP1(ILINK))%NR2 = &
WGTDATA(GRID2_ADD_MAP1(ILINK))%NR2 + 1
END IF
NGOODPNTS=NGOODPNTS+1
ELSEIF ( (GRID1_FRAC(GRID1_ADD_MAP1(ILINK))>FRAC_LOWEST) &
.AND. (GRID1_FRAC(GRID1_ADD_MAP1(ILINK))<FRAC_HIGHEST)&
) THEN
! note that we don't bother giving a warning for the
! cases where grid1_frac is strange.
NBADPNTS=NBADPNTS+1
! we also don't bother giving a warning if we've already
! given a lot of warnings (we keep counting, though)
IF((NBADPNTS.LT.5).AND.L_MASTER)THEN ! print warnings
WRITE(*,'(A)')'%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%'
WRITE(*,'(A)')'WARNING: SCRIP weight problem '
WRITE(*,'(4x,A,I7,A,I7)')'ilink = ',ilink,' out of ',&
max_links_map1
WRITE(*,'(4x,A,I7)')'grid1_add_map1(ilink) = ', &
grid1_add_map1(ilink)
WRITE(*,'(4x,A,I7)')'grid2_add_map1(ilink) = ', &
grid2_add_map1(ilink)
WRITE(*,'(4x,A,E12.4)')'wts_map1(1,ilink) = ', &
wts_map1(1,ilink)
WRITE(*,'(4x,A,F10.5)') &
'grid1_frac(grid1_add_map1(ilink)) = ', &
grid1_frac(grid1_add_map1(ilink))
WRITE(*,'(4x,A,F10.5)') &
'grid2_frac(grid2_add_map1(ilink)) = ', &
grid2_frac(grid2_add_map1(ilink))
WRITE(*,'(4x,A,F10.5)')'grid1_center_lat = ', &
grid1_center_lat(grid1_add_map1(ilink))
WRITE(*,'(4x,A,F10.5)')'grid1_center_lon = ', &
grid1_center_lon(grid1_add_map1(ilink))
WRITE(*,'(4x,A,F10.5)')'grid2_center_lat = ', &
grid2_center_lat(grid2_add_map1(ilink))
WRITE(*,'(4x,A,F10.5)')'grid2_center_lon = ', &
grid2_center_lon(grid2_add_map1(ilink))
WRITE(*,'(A)')'%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%'
ENDIF
ENDIF
ENDIF
END DO
IF((NBADPNTS.GT.0).AND.L_MASTER)THEN
WRITE(*,'(4x,A,I5,A)')'We had problems in ',NBADPNTS, &
' points.'
WRITE(*,'(4x,I8,A)')NGOODPNTS,' points appear to be OK.'
ENDIF
IF( (NBADPNTS.GT.(NGOODPNTS/30)) .AND.L_MASTER )THEN
WRITE(*,'(4x,A)')'Error: excessive SCRIP failure. Stopping.'
STOP 'wmscrpmd, case 1'
ENDIF
!.....step 2: allocate according to the size "n" determined above
DO KDST=1,GRID2_SIZE
ALLOCATE ( WGTDATA(KDST)%W(WGTDATA(KDST)%N), STAT=ISTAT )
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE ( WGTDATA(KDST)%K(WGTDATA(KDST)%N), STAT=ISTAT )
CHECK_ALLOC_STATUS ( ISTAT )
WGTDATA(KDST)%N=0
END DO
!.....step 3: save weights
DO ILINK=1,MAX_LINKS_MAP1
KSRC=GRID1_ADD_MAP1(ILINK)
JSRC=INT((KSRC-1)/NI_SRC)+1
ISRC=KSRC-(JSRC-1)*NI_SRC
!........note that this single if-then *must* be consistent with the
!........pair of if-thens above
IF ((GRID2_FRAC(GRID2_ADD_MAP1(ILINK))>FRAC_LOWEST) .AND. &
& (GRID2_FRAC(GRID2_ADD_MAP1(ILINK))<FRAC_HIGHEST) .AND. &
& (WTS_MAP1(1,ILINK)>=WT_LOWEST) .AND. &
& (WTS_MAP1(1,ILINK)<=WT_HIGHEST))THEN
IF (ABS(MAPSTA_SRC(JSRC,ISRC)).EQ.1) THEN ! sea point
WGTDATA(GRID2_ADD_MAP1(ILINK))%N= &
& WGTDATA(GRID2_ADD_MAP1(ILINK))%N+1
WGTDATA(GRID2_ADD_MAP1(ILINK))%W(WGTDATA( &
& GRID2_ADD_MAP1(ILINK))%N)=WTS_MAP1(1,ILINK)
WGTDATA(GRID2_ADD_MAP1(ILINK))%K(WGTDATA( &
& GRID2_ADD_MAP1(ILINK))%N)=GRID1_ADD_MAP1(ILINK)
ENDIF
ENDIF
END DO
!.....step 4: re-normalize weights. This is necessary because we called
!.....scrip without the mask. Now that we have the mask in place, we need
!.....to re-normalize the weights.
DO KDST=1,GRID2_SIZE
IF (WGTDATA(KDST)%N > 0) THEN
WTSUM=ZERO
DO IPNT=1,WGTDATA(KDST)%N
WTSUM=WTSUM+WGTDATA(KDST)%W(IPNT)
END DO
DO IPNT=1,WGTDATA(KDST)%N
WGTDATA(KDST)%W(IPNT)=WGTDATA(KDST)%W(IPNT)/WTSUM
END DO
END IF
END DO
CALL SCRIP_CLEAR
END SUBROUTINE SCRIP_WRAPPER
!/ ------------------------------------------------------------------- /
!>
!> @brief Compute grid arrays for scrip for a specific unstructured grid.
!>
!> @details For interior vertices, we select for every triangle the barycenter
!> of the triangle. So to every node contained in N triangles we associate
!> a domain with N corners. Every one of those corners is contained
!> in 3 different domains.
!>
!> For nodes that are on the boundary, we have to proceed differently.
!> For every such node, we have NEIGHBOR_PREV and NEIGHBOR_NEXT which
!> are the neighbor on each side of the boundary.
!> We put a corner on the middle of the edge. We also put a corner
!> on the node itself.
!>
!> Note that instead of taking barycenters, we could have taken
!> Voronoi vertices, but this carries danger since Voronoi vertices
!> can be outside of the triangle. And it leaves open how to treat
!> the boundary.
!>
!> @param[in] ID_GRD
!> @param[out] GRID_CENTER_LON
!> @param[out] GRID_CENTER_LAT
!> @param[out] GRID_CORNER_LON
!> @param[out] GRID_CORNER_LAT
!> @param[out] GRID_MASK
!> @param[out] GRID_DIMS
!> @param[out] GRID_SIZE
!> @param[out] GRID_CORNERS
!> @param[out] GRID_RANK
!>
!> @author M. Dutour
!> @author A. Roland
!> @date 10-Dec-2014
!>
SUBROUTINE GET_SCRIP_INFO_UNSTRUCTURED (ID_GRD, &
& GRID_CENTER_LON, GRID_CENTER_LAT, &
& GRID_CORNER_LON, GRID_CORNER_LAT, GRID_MASK, &
& GRID_DIMS, GRID_SIZE, GRID_CORNERS, GRID_RANK)
!/ +-----------------------------------+
!/ | WAVEWATCH III |
!/ | M. Dutour, A. Roland |
!/ | FORTRAN 90 |
!/ | Last update : 10-Dec-2014 !
!/ +-----------------------------------+
!/
! 1. Original author :
!
! Mathieu Dutour Sikiric, IRB & Aron Roland, Z&P
!
! 2. Last update :
!
! See revisions.
!
! 3. Revisions :
!
! 20-Feb-2012 : Origination
!/ 10-Dec-2014 : Add checks for allocate status ( version 5.04 )
!
! 4. Copyright :
!
! 5. Purpose :
!
! Compute grid arrays for scrip for a specific unstructured grid
! For interior vertices, we select for every triangle the barycenter
! of the triangle. So to every node contained in N triangles we associate
! a domain with N corners. Every one of those corners is contained
! in 3 different domains.
! For nodes that are on the boundary, we have to proceed differently.
! For every such node, we have NEIGHBOR_PREV and NEIGHBOR_NEXT which
! are the neighbor on each side of the boundary.
! We put a corner on the middle of the edge. We also put a corner
! on the node itself.
! Note that instead of taking barycenters, we could have taken
! Voronoi vertices, but this carries danger since Voronoi vertices
! can be outside of the triangle. And it leaves open how to treat
! the boundary.
!
! 6. Method :
!
! 7. Parameters, Variables and types :
!
! 8. Called by :
!
! Subroutine get_scrip_info
!
! 9. Subroutines and functions used :
!
! 10. Error messages:
!
! 11. Remarks :
!
! 12. Structure :
!
! 13. Switches :
!
! 14. Source code :
USE W3SERVMD, ONLY: EXTCDE
USE W3GDATMD, ONLY : GRIDS
IMPLICIT NONE
INTEGER, INTENT(IN) :: ID_GRD
REAL*8, INTENT(OUT), ALLOCATABLE :: GRID_CENTER_LON(:)
REAL*8, INTENT(OUT), ALLOCATABLE :: GRID_CENTER_LAT(:)
LOGICAL, INTENT(OUT), ALLOCATABLE :: GRID_MASK(:)
REAL*8, INTENT(OUT), ALLOCATABLE :: GRID_CORNER_LON(:,:)
REAL*8, INTENT(OUT), ALLOCATABLE :: GRID_CORNER_LAT(:,:)
INTEGER, INTENT(OUT), ALLOCATABLE :: GRID_DIMS(:)
INTEGER, INTENT(OUT) :: GRID_SIZE, GRID_CORNERS, GRID_RANK
INTEGER DIRAPPROACH, DUALAPPROACH, THEAPPROACH
INTEGER MNE, MNP, IE, IP, I
INTEGER NBPLUS, NBMINUS
INTEGER I1, I2, I3
REAL*8 :: ELON1, ELON2, ELON3, ELON, ELONC
REAL*8 :: ELAT1, ELAT2, ELAT3, ELAT, ELATC
REAL *8 :: DELTALON12, DELTALON13, DELTALAT12, DELTALAT13
REAL *8 :: THEDET
REAL*8 :: PT(3,2)
INTEGER, POINTER :: IOBP(:), TRIGINCD(:)
INTEGER, POINTER :: NEIGHBOR_PREV(:), NEIGHBOR_NEXT(:)
INTEGER, POINTER :: NBASSIGNEDCORNER(:), LISTNBCORNER(:)
INTEGER, POINTER :: STATUS(:), NEXTVERT(:), PREVVERT(:), FINALVERT(:)
INTEGER :: MAXCORNER, NBCORNER
INTEGER :: IDX, IPNEXT, IPPREV, NB, INEXT, IPREV
REAL*8, POINTER :: LON_CENT_TRIG(:), LAT_CENT_TRIG(:)
REAL*8 :: ELONIP, ELONNEXT, ELONPREV, ELONN, ELONP
REAL*8 :: ELATIP, ELATNEXT, ELATPREV, ELATN, ELATP
INTEGER :: ISFINISHED, ZPREV
INTEGER :: DODEBUG
GRID_RANK=2
DIRAPPROACH=1
DUALAPPROACH=2
THEAPPROACH=DUALAPPROACH
MNE=GRIDS(ID_GRD)%NTRI
MNP=GRIDS(ID_GRD)%NX
IF (THEAPPROACH .EQ. DIRAPPROACH) THEN
ALLOCATE(GRID_CENTER_LON(MNE), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(GRID_CENTER_LAT(MNE), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(GRID_CORNER_LON(3,MNE), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(GRID_CORNER_LAT(3,MNE), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(GRID_MASK(MNE), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
DO IE=1,MNE
I1=GRIDS(ID_GRD)%TRIGP(1,IE)
I2=GRIDS(ID_GRD)%TRIGP(2,IE)
I3=GRIDS(ID_GRD)%TRIGP(3,IE)
ELON1=GRIDS(ID_GRD)%XGRD(1,I1)
ELON2=GRIDS(ID_GRD)%XGRD(1,I2)
ELON3=GRIDS(ID_GRD)%XGRD(1,I3)
ELAT1=GRIDS(ID_GRD)%YGRD(1,I1)
ELAT2=GRIDS(ID_GRD)%YGRD(1,I2)
ELAT3=GRIDS(ID_GRD)%YGRD(1,I3)
ELON=(ELON1 + ELON2 + ELON3)/3
ELAT=(ELAT1 + ELAT2 + ELAT3)/3
GRID_CENTER_LON(IE)=ELON
GRID_CENTER_LAT(IE)=ELAT
GRID_CORNER_LON(1,IE)=ELON1
GRID_CORNER_LON(2,IE)=ELON2
GRID_CORNER_LON(3,IE)=ELON3
GRID_CORNER_LAT(1,IE)=ELAT1
GRID_CORNER_LAT(2,IE)=ELAT2
GRID_CORNER_LAT(3,IE)=ELAT3
GRID_MASK(IE)=.TRUE.
END DO
GRID_CORNERS=3
END IF
IF (THEAPPROACH .EQ. DUALAPPROACH) THEN
ALLOCATE(TRIGINCD(MNP), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(IOBP(MNP), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(NEIGHBOR_NEXT(MNP), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(NEIGHBOR_PREV(MNP), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(NBASSIGNEDCORNER(MNP), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(LISTNBCORNER(MNP), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(STATUS(MNP), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(NEXTVERT(MNP), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(PREVVERT(MNP), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(FINALVERT(MNP), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(LON_CENT_TRIG(MNE), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(LAT_CENT_TRIG(MNE), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
CALL GET_UNSTRUCTURED_VERTEX_DEGREE (MNP, MNE, &
GRIDS(ID_GRD)%TRIGP, TRIGINCD)
CALL GET_BOUNDARY(MNP, MNE, GRIDS(id_grd)%TRIGP, IOBP, &
NEIGHBOR_PREV, NEIGHBOR_NEXT)
! Find max number of corners
MAXCORNER=0
DO IP=1,MNP
IF (NEIGHBOR_NEXT(IP) .EQ. 0) THEN
NBCORNER=TRIGINCD(IP)
ELSE
NBCORNER=TRIGINCD(IP) + 3
END IF
LISTNBCORNER(IP)=NBCORNER
IF (NBCORNER .GT. MAXCORNER) THEN
MAXCORNER=NBCORNER
END IF
END DO
GRID_CORNERS=MAXCORNER
ALLOCATE(GRID_CENTER_LON(MNP), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(GRID_CENTER_LAT(MNP), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(GRID_CORNER_LON(MAXCORNER,MNP), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(GRID_CORNER_LAT(MAXCORNER,MNP), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
ALLOCATE(GRID_MASK(MNP), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
! Add first three corners for boundaries
NBASSIGNEDCORNER(:)=0
DO IP=1,MNP
GRID_MASK(IP)=.TRUE.
IF (NEIGHBOR_NEXT(IP) .GT. 0) THEN
IPNEXT=NEIGHBOR_NEXT(IP)
IPPREV=NEIGHBOR_PREV(IP)
ELONIP=DBLE(GRIDS(ID_GRD)%XGRD(1,IP))
ELATIP=DBLE(GRIDS(ID_GRD)%YGRD(1,IP))
ELONNEXT=DBLE(GRIDS(ID_GRD)%XGRD(1,IPNEXT))
ELATNEXT=DBLE(GRIDS(ID_GRD)%YGRD(1,IPNEXT))
ELONPREV=DBLE(GRIDS(ID_GRD)%XGRD(1,IPPREV))
ELATPREV=DBLE(GRIDS(ID_GRD)%YGRD(1,IPPREV))
! Periodicity fix for corner node
IF ( ABS(ELONIP - ELONNEXT) .GT. 180.0 ) THEN
ELONNEXT = ELONNEXT -SIGN(360.0d0,(ELONIP - ELONNEXT))
ENDIF
IF ( ABS(ELONIP - ELONPREV) .GT. 180.0 ) THEN
ELONPREV = ELONPREV -SIGN(360.0d0,(ELONIP - ELONPREV))
ENDIF
ELONN=(ELONIP+ELONNEXT)/2.0
ELATN=(ELATIP+ELATNEXT)/2.0
ELONP=(ELONIP+ELONPREV)/2.0
ELATP=(ELATIP+ELATPREV)/2.0
GRID_CORNER_LON(1,IP)=ELONN
GRID_CORNER_LAT(1,IP)=ELATN
GRID_CORNER_LON(2,IP)=ELONIP
GRID_CORNER_LAT(2,IP)=ELATIP
GRID_CORNER_LON(3,IP)=ELONP
GRID_CORNER_LAT(3,IP)=ELATP
NBASSIGNEDCORNER(IP)=3
END IF
END DO
! Compute centers
DO IP=1,MNP
GRID_CENTER_LON(IP)=DBLE(GRIDS(ID_GRD)%XGRD(1,IP))
GRID_CENTER_LAT(IP)=DBLE(GRIDS(ID_GRD)%YGRD(1,IP))
END DO
! Check triangle node orientation
! Compute triangle centers
NBPLUS=0
NBMINUS=0
DO IE=1,MNE
I1=GRIDS(ID_GRD)%TRIGP(1,IE)
I2=GRIDS(ID_GRD)%TRIGP(2,IE)
I3=GRIDS(ID_GRD)%TRIGP(3,IE)
PT(1,1)=DBLE(GRIDS(ID_GRD)%XGRD(1,I1))
PT(2,1)=DBLE(GRIDS(ID_GRD)%XGRD(1,I2))
PT(3,1)=DBLE(GRIDS(ID_GRD)%XGRD(1,I3))
PT(1,2)=DBLE(GRIDS(ID_GRD)%YGRD(1,I1))
PT(2,2)=DBLE(GRIDS(ID_GRD)%YGRD(1,I2))
PT(3,2)=DBLE(GRIDS(ID_GRD)%YGRD(1,I3))
CALL FIX_PERIODCITY(PT)
ELON1 = PT(1,1)
ELON2 = PT(2,1)
ELON3 = PT(3,1)
ELAT1 = PT(1,2)
ELAT2 = PT(2,2)
ELAT3 = PT(3,2)
DELTALON12=ELON2 - ELON1
DELTALON13=ELON3 - ELON1
DELTALAT12=ELAT2 - ELAT1
DELTALAT13=ELAT3 - ELAT1
THEDET=DELTALON12*DELTALAT13 - DELTALON13*DELTALAT12
IF (THEDET.GT.0) THEN
NBPLUS=NBPLUS+1
END IF
IF (THEDET.LT.0) THEN
NBMINUS=NBMINUS+1
END IF
ELON=(ELON1 + ELON2 + ELON3)/3.0
ELAT=(ELAT1 + ELAT2 + ELAT3)/3.0
LON_CENT_TRIG(IE)=ELON
LAT_CENT_TRIG(IE)=ELAT
END DO
DODEBUG=0
IF (DODEBUG.EQ.1) THEN
print *, 'nbplus=', nbplus, ' nbminus=', nbminus
END IF
STATUS(:) = 0
NEXTVERT(:) = 0
PREVVERT(:) = 0
DO IE=1,MNE
DO I=1,3
CALL TRIANG_INDEXES(I, INEXT, IPREV)
IP=GRIDS(ID_GRD)%TRIGP(I,IE)
IPNEXT=GRIDS(ID_GRD)%TRIGP(INEXT,IE)
IPPREV=GRIDS(ID_GRD)%TRIGP(IPREV,IE)
IF (STATUS(IP).EQ.0) THEN
IF (NEIGHBOR_NEXT(IP).EQ.0) THEN
STATUS(IP)=1
FINALVERT(IP)=IPPREV
PREVVERT(IP)=IPPREV
NEXTVERT(IP)=IPNEXT
ELSE
IF (NEIGHBOR_PREV(IP).EQ.IPPREV) THEN
STATUS(IP)=1
PREVVERT(IP)=IPPREV
NEXTVERT(IP)=IPNEXT
FINALVERT(IP)=NEIGHBOR_NEXT(IP)
END IF
END IF
END IF
END DO
END DO
STATUS(:)=0
DO
ISFINISHED=1
DO IE=1,MNE
ELON=LON_CENT_TRIG(IE)
ELAT=LAT_CENT_TRIG(IE)
DO I=1,3
CALL TRIANG_INDEXES(I, INEXT, IPREV)
IP=GRIDS(ID_GRD)%TRIGP(I,IE)
IPNEXT=GRIDS(ID_GRD)%TRIGP(INEXT,IE)
IPPREV=GRIDS(ID_GRD)%TRIGP(IPREV,IE)
IF (STATUS(IP).EQ.0) THEN
ISFINISHED=0
ZPREV=PREVVERT(IP)
IF (ZPREV.EQ.IPPREV) THEN
IDX=NBASSIGNEDCORNER(IP)
IDX=IDX+1
GRID_CORNER_LON(IDX,IP)=ELON
GRID_CORNER_LAT(IDX,IP)=ELAT
NBASSIGNEDCORNER(IP)=IDX
PREVVERT(IP)=IPNEXT
IF (IPNEXT.EQ.FINALVERT(IP)) THEN
STATUS(IP)=1
END IF
END IF
END IF
END DO
END DO
IF (ISFINISHED.EQ.1) THEN
EXIT
END IF
END DO
DO IP=1,MNP
IF (NBASSIGNEDCORNER(IP).NE.LISTNBCORNER(IP)) THEN
WRITE(*,*) 'Incoherent number at IP=', IP
WRITE(*,*) ' NbAssignedCorner(IP)=', NbAssignedCorner(IP)
WRITE(*,*) ' ListNbCorner(IP)=', ListNbCorner(IP)
WRITE(*,*) ' N_N=', NEIGHBOR_NEXT(IP), 'N_P=', NEIGHBOR_PREV(IP)
WRITE(*,*) ' TrigIncd=', TrigIncd(IP)
STOP 'wmscrpmd, case 2'
END IF
END DO
! if the number of corner is below threshold, we have to
! add some more.
DO IP=1,MNP
NB=NBASSIGNEDCORNER(IP)
IF (NB .LT. MAXCORNER) THEN
ELON=GRID_CORNER_LON(NB,IP)
ELAT=GRID_CORNER_LAT(NB,IP)
DO IDX=NB+1,MAXCORNER
GRID_CORNER_LON(IDX,IP)=ELON
GRID_CORNER_LAT(IDX,IP)=ELAT
END DO
END IF
END DO
DEALLOCATE(NBASSIGNEDCORNER, STAT=ISTAT)
CHECK_DEALLOC_STATUS ( ISTAT )
DEALLOCATE(LISTNBCORNER, STAT=ISTAT)
CHECK_DEALLOC_STATUS ( ISTAT )
DEALLOCATE(TRIGINCD, STAT=ISTAT)
CHECK_DEALLOC_STATUS ( ISTAT )
DEALLOCATE(IOBP, STAT=ISTAT)
CHECK_DEALLOC_STATUS ( ISTAT )
DEALLOCATE(NEIGHBOR_PREV, STAT=ISTAT)
CHECK_DEALLOC_STATUS ( ISTAT )
DEALLOCATE(NEIGHBOR_NEXT, STAT=ISTAT)
CHECK_DEALLOC_STATUS ( ISTAT )
DEALLOCATE(STATUS, STAT=ISTAT)
CHECK_DEALLOC_STATUS ( ISTAT )
DEALLOCATE(NEXTVERT, STAT=ISTAT)
CHECK_DEALLOC_STATUS ( ISTAT )
DEALLOCATE(PREVVERT, STAT=ISTAT)
CHECK_DEALLOC_STATUS ( ISTAT )
DEALLOCATE(FINALVERT, STAT=ISTAT)
CHECK_DEALLOC_STATUS ( ISTAT )
DEALLOCATE(LON_CENT_TRIG, STAT=ISTAT)
CHECK_DEALLOC_STATUS ( ISTAT )
DEALLOCATE(LAT_CENT_TRIG, STAT=ISTAT)
CHECK_DEALLOC_STATUS ( ISTAT )
ALLOCATE(GRID_DIMS(2), STAT=ISTAT)
CHECK_ALLOC_STATUS ( ISTAT )
GRID_DIMS(1)=MNP
GRID_DIMS(2)=1
GRID_SIZE=MNP
END IF
END SUBROUTINE GET_SCRIP_INFO_UNSTRUCTURED
!/ ------------------------------------------------------------------- /
!>
!> @brief Compute grid arrays needed for scrip for a specific
!> structured grid.
!>
!> @details This is adapted from Erick Rogers original code by
!> splitting the original scrip_wrapper function.
!>
!> @param[in] ID_GRD
!> @param[out] GRID_CENTER_LON
!> @param[out] GRID_CENTER_LAT
!> @param[out] GRID_CORNER_LON
!> @param[out] GRID_CORNER_LAT
!> @param[out] GRID_MASK
!> @param[out] GRID_DIMS
!> @param[out] GRID_SIZE
!> @param[out] GRID_CORNERS
!> @param[out] GRID_RANK
!>
!> @author M. Dutour
!> @author A. Roland
!> @date 10-Dec-2014
!>
SUBROUTINE GET_SCRIP_INFO_STRUCTURED (ID_GRD, &
& GRID_CENTER_LON, GRID_CENTER_LAT, &
& GRID_CORNER_LON, GRID_CORNER_LAT, GRID_MASK, &
& GRID_DIMS, GRID_SIZE, GRID_CORNERS, GRID_RANK)
!/ +-----------------------------------+
!/ | WAVEWATCH III |
!/ | M. Dutour, A. Roland |
!/ | FORTRAN 90 |
!/ | Last update : 10-Dec-2014 !
!/ +-----------------------------------+
!/
! 1. Original author :
!
! Mathieu Dutour Sikiric, IRB & Aron Roland, Z&P
!
! 2. Last update :
!
! See revisions.
!
! 3. Revisions :
!
! 20-Feb-2012 : Origination, this is adapted from Erick Rogers
! code by splitting the code into sections.
!/ 10-Dec-2014 : Add checks for allocate status ( version 5.04 )
!
! 4. Copyright :
!
! 5. Purpose :
!
! Compute grid arrays needed for scrip for a specific structured grid.
! This is adapted from Erick Rogers original code by splitting