-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathharm.F
2231 lines (2111 loc) · 82.7 KB
/
harm.F
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
C******************************************************************************
C PADCIRC VERSION 45.12 03/17/2006 *
C last changes in this file prior to VERSION 44.15 *
C *
c *
c added local LOGICAL CHARMV declaration 04/21/2004 *
c******************************************************************************
c *
c PADCIRC MODULE ( HARM ) *
c *
c HA_SUBS.FOR V3.01 11/9/95 *
c *
c Least Square harmonic analysis of timeseries from ADCIRC2DDI_v27 *
c *
c Notes: *
c 1.) Both the left hand side matrix and the right hand side *
c forcing vectors are continuously updated in time. This *
c eliminates the need to store time series outputs for later *
c harmonic analysis. *
c 2.) The left hand side matrix and the right hand side forcing *
c vectors are output in the hotstart file and can be used to *
c perform harmonic analysis on an incomplete run. *
c 3.) Frequencies should be in rad/sec,times should be in sec. *
c *
c***********************************************************************
c *
c Program Written by: *
c R.A. Luettich, IMS UNC *
c J.J. Westerink, CE ND *
c *
c Program Development History: *
c 1.) lsq_stations_v004 by JJW *
c 2.) LSQEX by RL used in 2D binary extr program *
c 3.) LSQRL by RL used in 1D test codes *
c 4.) LSQ2D v1.00-v2.26 by RL real time Harmonic Analysis for ADCIRC*
c 5.) HA_SUBS v3.01 by RL real time HA for ADCIRC separate *
c subroutines for elevation station, velocity station, *
c global elevation and global velocity harmonic analysis *
c *
c***********************************************************************
c *
c SUBROUTINE LSQUPDLHS updates the LHS matrix *
c SUBROUTINE LSQUPDES updates the RHS load vector for elev stations *
c SUBROUTINE LSQUPDVS updates the RHS load vector for velocity stations*
c SUBROUTINE LSQUPDEG updates the RHS load vector for elevation global *
c SUBROUTINE LSQUPDVG updates the RHS load vector for velocity global *
c SUBROUTINE FULSOL fills out, decomposes and solves the matricies *
c SUBROUTINE LSQSOLES solves & writes output for elevation stations *
c SUBROUTINE LSQSOLVS solves & writes output for velocity stations *
c SUBROUTINE LSQSOLEG solves & writes output for elevation global *
c SUBROUTINE LSQSOLVG solves & writes output for velocity global *
c SUBROUTINE HAHOUT writes HA parameters & LHS matrix to hotstart file *
c SUBROUTINE HAHOUTES writes elev sta RHS load vector to hotstart file *
c SUBROUTINE HAHOUTVS writes vel sta RHS load vector to hotstart file *
c SUBROUTINE HAHOUTEG writes glob elev RHS load vector to hotstart file*
c SUBROUTINE HAHOUTVG writes glob vel RHS load vector to hotstart file *
c SUBROUTINE HACOLDS initializes HA param & LHS matrix for cold start *
c SUBROUTINE HACOLDSES initializes elev sta RHS load vec for cold start*
c SUBROUTINE HACOLDSVS initializes vel sta RHS load vec for cold start *
c SUBROUTINE HACOLDSEG initializes glob ele RHS load vec for cold start*
c SUBROUTINE HACOLDSVG initializes glob vel RHS load vec for cold start*
c SUBROUTINE HAHOTS initializes HA params & LHS matrix for a hot start *
c SUBROUTINE HAHOTSES initializes elev sta RHS load vec for hot start *
c SUBROUTINE HAHOTSVS initializes vel sta RHS load vec for hot start *
c SUBROUTINE HAHOTSEG initializes glob elev RHS load vec for hot start *
c SUBROUTINE HAHOTSVG initializes glob vel RHS load vec for hot start *
c *
c***********************************************************************
c *
c INPUT FILES: *
c - Frequency information is read in by ADCIRC from unit 15. *
c *
c - If the model is hot start, input is read from UNIT 67 or 68 *
c *
c OUTPUT FILES: *
C UNIT 51 : HARMONIC CONSTITUENT ELEVATION VALUES AT SPECIFIED *
C ELEVATION RECORDING STATION COORDINATES (ASCII) *
C UNIT 52 : HARMONIC CONSTITUENT VELOCITY VALUES AT SPECIFIED *
C VELOCITY RECORDING STATION COORDINATES (ASCII) *
C UNIT 53 : HARMONIC CONSTITUENT ELEVATIONS AT ALL NODES (ASCII) *
C UNIT 54 : HARMONIC CONSTITUENT VELOCITIES AT ALL NODES (ASCII) *
C UNIT 55 : COMPARISON BETWEEN THE MEAN AND VARIANCE OF THE TIME *
C SERIES GENERATED BY THE MODEL AND THE MEAN AND *
C VARIANCE OF A TIME SERIES RESYNTHESIZED FROM THE *
C COMPUTED HARMONIC CONSTITUENTS. THIS GIVES AN *
C INDICATION OF HOW COMPLETE THE HARMONIC ANALYSIS *
C WAS. (ASCII) *
C UNIT 67 or 68 : HOT START FILES (BINARY) *
c *
c***********************************************************************
C
MODULE HARM
C
USE SIZES, ONLY : SZ, MNP, MNHARF, MNSTAE, MNSTAV, LOCALDIR,
& WRITE_LOCAL_HARM_FILES
USE GLOBAL, ONLY : screenMessage, logMessage, DEBUG, ECHO, INFO,
& WARNING, ERROR, setMessageSource, pi, IFSPROTS,
& unsetMessageSource, allMessage, scratchMessage
C
C
! jgf49.44: parameters describing the harmonic constituents to be
! included in the harmonic analysis of model results
INTEGER :: IHARIND ! =1 if harmonic analysis was requested, 0 otherwise
LOGICAL :: CHARMV=.false. ! .true. for timeseries reconstruction
INTEGER :: NFREQ ! number of frequencies
CHARACTER*10,ALLOCATABLE :: NAMEFR(:) ! "M2", "S2", etc
REAL(SZ), ALLOCATABLE :: HAFREQ(:) ! frequency (rad/s)
REAL(SZ), ALLOCATABLE :: HAFF(:) ! nodal factor
REAL(SZ), ALLOCATABLE :: HAFACE(:) ! equilibrium argument (deg)
C
INTEGER :: NZ ! nz=0 if there is a steady component, nz=1 otherwise
INTEGER :: NF ! nf=0 if no steady constituent, nf=1 otherwise
C
INTEGER :: ITHAS ! time step upon which harmonic analysis starts
INTEGER :: ITHAF ! time step upon which harmonic analysis finishes
INTEGER :: ITMV ! time step upon which means and variance calcs start
INTEGER :: TIMEBEG ! model time (sec) at which means and var calcs start
INTEGER :: IHABEG ! 1 h.a. time increment past the start of h.a.
INTEGER :: ICHA ! spool counter for HA updates
INTEGER :: NTSTEPS! number of time steps into means and variance calcs
INTEGER :: ICALL ! number of times HA update sub has been called
INTEGER :: ITUD ! time step of most recent HA update
REAL(8) :: TIMEUD ! model time of most recent HA update
C
REAL(SZ), ALLOCATABLE :: HA(:,:) ! LHS matrix
INTEGER :: MM ! 2nd dim of HA matrix; 2x num freqs (+ steady)
REAL(SZ), ALLOCATABLE :: HAP(:) ! used in matrix solution
REAL(SZ), ALLOCATABLE :: HAX(:) ! used in matrix solution
C
C Load Vectors
REAL(SZ), ALLOCATABLE, TARGET :: GLOELV(:,:) ! fulldomain elevation
REAL(SZ), ALLOCATABLE, TARGET :: GLOULV(:,:) ! fulldomain u velocity
REAL(SZ), ALLOCATABLE, TARGET :: GLOVLV(:,:) ! fulldomain v velocity
REAL(SZ), ALLOCATABLE, TARGET :: STAELV(:,:) ! station elevation
REAL(SZ), ALLOCATABLE, TARGET :: STAULV(:,:) ! station u velocity
REAL(SZ), ALLOCATABLE, TARGET :: STAVLV(:,:) ! station v velocity
C
C Arrays for time series reconstruction.
REAL(SZ),ALLOCATABLE,TARGET :: XVELAV(:)
REAL(SZ),ALLOCATABLE,TARGET :: YVELAV(:)
REAL(SZ),ALLOCATABLE,TARGET :: XVELVA(:)
REAL(SZ),ALLOCATABLE,TARGET :: YVELVA(:)
REAL(SZ),ALLOCATABLE,TARGET :: ELAV(:)
REAL(SZ),ALLOCATABLE,TARGET :: ELVA(:)
C
C jgf49.44: Parameters that control the spatial locations where
C harmonic analysis is performed:
INTEGER :: NHASE ! /= 0 to perform HA at elevation recording stations
INTEGER :: NHASV ! /= 0 to perform HA at velocity recording stations
INTEGER :: NHAGE ! /= 0 to do elevation HA at all nodes in the mesh
INTEGER :: NHAGV ! /= 0 to do velocity HA at all nodes in the mesh
C
C jgf49.44: Parameters that control the calculation of harmonic
C constituents:
REAL(SZ) :: THAS ! number of days b/f harmonic analysis starts
REAL(SZ) :: THAF ! number of days at which harmonic analysis ends
INTEGER :: NHAINC ! time step increment for including data in HA
REAL(SZ) :: FMV ! fraction of HA period for means and var. calcs (0to1)
C
C jgf49.44: Harmonic analysis hotstarting:
REAL(SZ), ALLOCATABLE, TARGET :: GLOELV_g(:,:)
REAL(SZ), ALLOCATABLE, TARGET :: STAELV_g(:,:)
REAL(SZ), ALLOCATABLE, TARGET :: GLOULV_g(:,:)
REAL(SZ), ALLOCATABLE, TARGET :: GLOVLV_g(:,:)
REAL(SZ), ALLOCATABLE, TARGET :: STAULV_g(:,:)
REAL(SZ), ALLOCATABLE, TARGET :: STAVLV_g(:,:)
REAL(SZ), ALLOCATABLE, TARGET :: ELAV_g(:)
REAL(SZ), ALLOCATABLE, TARGET :: ELVA_g(:)
REAL(SZ), ALLOCATABLE, TARGET :: XVELAV_g(:)
REAL(SZ), ALLOCATABLE, TARGET :: YVELAV_g(:)
REAL(SZ), ALLOCATABLE, TARGET :: XVELVA_g(:)
REAL(SZ), ALLOCATABLE, TARGET :: YVELVA_g(:)
C
C jgf49.44: Equivalence variables for reading in hotstart data, either
C in serial or in parallel.
REAL(SZ), POINTER :: STAELV_pt(:,:),STAULV_pt(:,:),STAVLV_pt(:,:)
REAL(SZ), POINTER :: GLOELV_pt(:,:),GLOULV_pt(:,:),GLOVLV_pt(:,:)
REAL(SZ), POINTER :: ELAV_pt(:), ELVA_pt(:)
REAL(SZ), POINTER :: XVELAV_pt(:), YVELAV_pt(:)
REAL(SZ), POINTER :: XVELVA_pt(:), YVELVA_pt(:)
C
C jgf49.44: The following variables are read from the hotstart file and
C later compared with the values from the fort.15 to ensure that they
C match.
INTEGER INHASE, INHASV, INHAGE, INHAGV, IFLAG
INTEGER INFREQ, INSTAE, INSTAV, INP, INZ, INF, IMM, IICALL
REAL(SZ),ALLOCATABLE :: IFREQ(:),IFF(:),IFACE(:)
CHARACTER*10,ALLOCATABLE :: INAMEFR(:)
C
CHARACTER*16 FNAME
CHARACTER*8 FNAM8(2)
EQUIVALENCE (FNAM8(1),FNAME)
C
C jgf49.44: Raised output variables to module level and added
C full domain counterparts for globalio. These arrays represent
C magnitudes and phases for each frequency at each station or node.
C Stations:
REAL(SZ), ALLOCATABLE, TARGET :: EMAG(:,:) ! elevation magnitudes
REAL(SZ), ALLOCATABLE, TARGET :: PHASEDE(:,:) ! elevation phases
REAL(SZ), ALLOCATABLE, TARGET :: UMAG(:,:) ! u velocity magnitudes
REAL(SZ), ALLOCATABLE, TARGET :: VMAG(:,:) ! v velocity magnitudes
REAL(SZ), ALLOCATABLE, TARGET :: PHASEDU(:,:) ! u velocity phases
REAL(SZ), ALLOCATABLE, TARGET :: PHASEDV(:,:) ! v velocity phases
C Nodes:
REAL(SZ), ALLOCATABLE, TARGET :: EMAGT(:,:) ! elevation magnitudes
REAL(SZ), ALLOCATABLE, TARGET :: PHASEDEN(:,:) ! elevation phases
REAL(SZ), ALLOCATABLE, TARGET :: UMAGT(:,:) ! u velocity magnitudes
REAL(SZ), ALLOCATABLE, TARGET :: VMAGT(:,:) ! v velocity magnitudes
REAL(SZ), ALLOCATABLE, TARGET :: PHASEDUT(:,:) ! u velocity phases
REAL(SZ), ALLOCATABLE, TARGET :: PHASEDVT(:,:) ! v velocity phases
C Time series reconstruction:
REAL(SZ), ALLOCATABLE, TARGET :: EAV(:)
REAL(SZ), ALLOCATABLE, TARGET :: ESQ(:)
REAL(SZ), ALLOCATABLE, TARGET :: EAVDIF(:)
REAL(SZ), ALLOCATABLE, TARGET :: EVADIF(:)
C
REAL(SZ), ALLOCATABLE, TARGET :: UAV(:)
REAL(SZ), ALLOCATABLE, TARGET :: USQ(:)
REAL(SZ), ALLOCATABLE, TARGET :: UAVDIF(:)
REAL(SZ), ALLOCATABLE, TARGET :: UVADIF(:)
C
REAL(SZ), ALLOCATABLE, TARGET :: VAV(:)
REAL(SZ), ALLOCATABLE, TARGET :: VSQ(:)
REAL(SZ), ALLOCATABLE, TARGET :: VAVDIF(:)
REAL(SZ), ALLOCATABLE, TARGET :: VVADIF(:)
!jgf51.52.01: Variables used in adcpost.
REAL(SZ), ALLOCATABLE,TARGET :: XELAV(:)
REAL(SZ), ALLOCATABLE,TARGET :: YELAV(:)
REAL(SZ), ALLOCATABLE,TARGET :: XELVA(:)
REAL(SZ), ALLOCATABLE,TARGET :: YELVA(:)
C-----------------END OF DECLARATIONS---------------------------------------
CONTAINS
C
C***********************************************************************
C Allocate arays used by LSQ_HARM.
C
C vjp 8/99
C***********************************************************************
SUBROUTINE ALLOC_HA()
IMPLICIT NONE
C
call setMessageSource("ALLOC_HA")
#if defined(HARM_TRACE) || defined (ALL_TRACE)
call allMessage(DEBUG, "Enter.")
#endif
ALLOCATE ( HAFREQ(MNHARF),HAFF(MNHARF),HAFACE(MNHARF) )
ALLOCATE ( NAMEFR(MNHARF) )
C
ALLOCATE ( HA(2*MNHARF,2*MNHARF) )
ALLOCATE ( HAP(2*MNHARF),HAX(2*MNHARF) )
ALLOCATE ( GLOELV(2*MNHARF,MNP) )
ALLOCATE ( GLOULV(2*MNHARF,MNP),GLOVLV(2*MNHARF,MNP) )
ALLOCATE ( STAELV(2*MNHARF,MNSTAE) )
ALLOCATE ( STAULV(2*MNHARF,MNSTAV),STAVLV(2*MNHARF,MNSTAV) )
C
ALLOCATE ( IFREQ(MNHARF),IFF(MNHARF),IFACE(MNHARF) )
ALLOCATE ( INAMEFR(MNHARF) )
C
ALLOCATE ( EMAG(MNHARF,MNSTAE), PHASEDE(MNHARF,MNSTAE) )
ALLOCATE ( UMAG(MNHARF,MNSTAV), PHASEDU(MNHARF,MNSTAV) )
ALLOCATE ( VMAG(MNHARF,MNSTAV), PHASEDV(MNHARF,MNSTAV) )
ALLOCATE ( EMAGT(MNHARF,MNP), PHASEDEN(MNHARF,MNP) )
ALLOCATE ( UMAGT(MNHARF,MNP), PHASEDUT(MNHARF,MNP) )
ALLOCATE ( VMAGT(MNHARF,MNP), PHASEDVT(MNHARF,MNP) )
C
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Return.")
#endif
call unsetMessageSource()
RETURN
C----------------------------------------------------------------------
END SUBROUTINE ALLOC_HA
C----------------------------------------------------------------------
C----------------------------------------------------------------------
C...
C...Allocate space for harmonic analysis means and variance
C...calculations
C...
C----------------------------------------------------------------------
SUBROUTINE ALLOC_MAIN14()
IMPLICIT NONE
C
call setMessageSource("ALLOC_MAIN14")
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Enter.")
#endif
C
ALLOCATE ( XVELAV(MNP),YVELAV(MNP),XVELVA(MNP),YVELVA(MNP) )
ALLOCATE ( ELAV(MNP),ELVA(MNP) )
C
ALLOCATE( EAV(MNP),ESQ(MNP),EAVDIF(MNP),EVADIF(MNP) )
ALLOCATE( UAV(MNP),USQ(MNP),UAVDIF(MNP),UVADIF(MNP) )
ALLOCATE( VAV(MNP),VSQ(MNP),VAVDIF(MNP),VVADIF(MNP) )
C
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Return.")
#endif
call unsetMessageSource()
RETURN
C----------------------------------------------------------------------
END SUBROUTINE ALLOC_MAIN14
C----------------------------------------------------------------------
C--------------------------------------------------------------------
C S U B R O U T I N E
C A L L O C A T E F U L L D O M A I N H A I O A R R A Y S
C--------------------------------------------------------------------
C jgf49.44: Allocates memory for the full domain arrays for i/o
C purposes when executing in parallel.
C--------------------------------------------------------------------
SUBROUTINE allocateFullDomainHAIOArrays()
USE GLOBAL, ONLY : NSTAE_G, NSTAV_G, NP_G
IMPLICIT NONE
C
call setMessageSource("allocateFullDomainHAIOArrays")
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Enter.")
#endif
IF (IHARIND.eq.1) THEN
ALLOCATE(GLOELV_g(2*MNHARF,NP_G))
ALLOCATE(STAELV_g(2*MNHARF,NSTAE_G))
ALLOCATE(GLOULV_g(2*MNHARF,NP_G))
ALLOCATE(GLOVLV_g(2*MNHARF,NP_G))
ALLOCATE(STAULV_g(2*MNHARF,NSTAV_G))
ALLOCATE(STAVLV_g(2*MNHARF,NSTAV_G))
ENDIF
IF (CHARMV) THEN
ALLOCATE(ELAV_g(NP_G))
ALLOCATE(ELVA_g(NP_G))
ALLOCATE(XVELAV_g(NP_G))
ALLOCATE(YVELAV_g(NP_G))
ALLOCATE(XVELVA_g(NP_G))
ALLOCATE(YVELVA_g(NP_G))
ENDIF
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Return.")
#endif
call unsetMessageSource()
C--------------------------------------------------------------------
END SUBROUTINE allocateFullDomainHAIOArrays
C--------------------------------------------------------------------
C--------------------------------------------------------------------
C S U B R O U T I N E W R I T E H A
C--------------------------------------------------------------------
C jgf50.97: Writes LHS matrix to a file for use in debugging.
C--------------------------------------------------------------------
SUBROUTINE writeHA(timesec, it)
USE SIZES, ONLY : globaldir, myproc
IMPLICIT NONE
REAL(8) :: timesec ! time in seconds since cold start
INTEGER :: it ! current adcirc time step number
INTEGER :: i, j ! loop counters
C
call setMessageSource("writeHA")
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Enter.")
#endif
IF (myproc.eq.0) THEN
OPEN(56,FILE=TRIM(GLOBALDIR)//'/'//'fort.56',
& STATUS='UNKNOWN',ACCESS='SEQUENTIAL',ACTION='WRITE',
& POSITION='APPEND')
WRITE(56,2120) timesec, it
DO i=1,2*MNHARF
DO j=1,2*MNHARF
write(56,2130) i, j, ha(i,j)
END DO
END DO
CLOSE(56)
ENDIF
2120 FORMAT(2X,1pE20.10E3,5X,I10)
2130 FORMAT('ha(',I2,',',I2,')=',1pE20.10E3)
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Return.")
#endif
call unsetMessageSource()
C--------------------------------------------------------------------
END SUBROUTINE writeHA
C--------------------------------------------------------------------
C--------------------------------------------------------------------
C S U B R O U T I N E C H E C K H A R M O N I C P A R A M E T E R S
C--------------------------------------------------------------------
C jgf50.41: Checks harmonic analysis settings.
C WJP 02.20.2018: Allowing for > 1 options for netcdf output
C--------------------------------------------------------------------
SUBROUTINE checkHarmonicParameters()
USE GLOBAL, ONLY : NFOVER, useNetCDF, useNetCDFOutput
IMPLICIT NONE
C
call setMessageSource("checkHarmonicParameters")
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Enter.")
#endif
IF((NHASE.LT.0).OR.(NHASE.GT.5)) THEN
CALL allMessage(WARNING,"Input value of NHASE is not valid.")
IF (NFOVER.EQ.1) THEN
CALL allMessage(WARNING,"NHASE will be set to zero.")
ELSE
CALL harmonicTerminate()
ENDIF
ELSE IF (NHASE.EQ.1) THEN
CALL logMessage(INFO,
& 'STATION ELEVATION HARMONIC ANALYSIS'
& //' WILL BE WRITTEN TO UNIT 51 IN ASCII FORMAT.')
ELSEIF (NHASE > 1) THEN
CALL logMessage(INFO,
& 'STATION ELEVATION HARMONIC ANALYSIS'
& //' WILL BE WRITTEN TO UNIT 51 IN NETCDF FORMAT.')
useNetCDF = .true.
useNetCDFOutput = .true.
ENDIF
IF((NHASV.LT.0).OR.(NHASV.GT.5)) THEN
CALL allMessage(WARNING,"Input value of NHASV is not valid.")
IF (NFOVER.EQ.1) THEN
CALL allMessage(WARNING,"NHASV will be set to zero.")
ELSE
CALL harmonicTerminate()
ENDIF
ELSE IF (NHASV.EQ.1) THEN
CALL logMessage(INFO,
& 'STATION VELOCITY HARMONIC ANALYSIS'
& //' WILL BE WRITTEN TO UNIT 52 IN ASCII FORMAT.')
ELSE IF (NHASV > 1) THEN
CALL logMessage(INFO,
& 'STATION VELOCITY HARMONIC ANALYSIS'
& //' WILL BE WRITTEN TO UNIT 52 IN NETCDF FORMAT.')
useNetCDF = .true.
useNetCDFOutput = .true.
ENDIF
IF((NHAGE.LT.0).OR.(NHAGE.GT.5)) THEN
CALL allMessage(WARNING,"Input value of NHAGE is not valid.")
IF (NFOVER.EQ.1) THEN
CALL allMessage(WARNING,"NHAGE will be set to zero.")
ELSE
CALL harmonicTerminate()
ENDIF
ELSE IF (NHAGE.EQ.1) THEN
CALL logMessage(INFO,
& 'FULL DOMAIN ELEVATION HARMONIC ANALYSIS'
& //' WILL BE WRITTEN TO UNIT 53 IN ASCII FORMAT.')
ELSE IF (NHAGE > 1) THEN
CALL logMessage(INFO,
& 'FULL DOMAIN ELEVATION HARMONIC ANALYSIS'
& //' WILL BE WRITTEN TO UNIT 53 IN NETCDF FORMAT.')
useNetCDF = .true.
useNetCDFOutput = .true.
ENDIF
IF((NHAGV.LT.0).OR.(NHAGV.GT.5)) THEN
CALL allMessage(WARNING,"Input value of NHAGV is not valid.")
IF (NFOVER.EQ.1) THEN
CALL allMessage(WARNING,"NHAGV will be set to zero.")
ELSE
CALL harmonicTerminate()
ENDIF
ELSE IF (NHAGV.EQ.1) THEN
CALL logMessage(INFO,
& 'FULL DOMAIN VELOCITY HARMONIC ANALYSIS'
& //' WILL BE WRITTEN TO UNIT 54 ASCII FORMAT.')
ELSE IF (NHAGV > 1) THEN
CALL logMessage(INFO,
& 'FULL DOMAIN VELOCITY HARMONIC ANALYSIS'
& //' WILL BE WRITTEN TO UNIT 54 NETCDF FORMAT.')
useNetCDF = .true.
useNetCDFOutput = .true.
ENDIF
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Return.")
#endif
call unsetMessageSource()
C--------------------------------------------------------------------
END SUBROUTINE checkHarmonicParameters
C--------------------------------------------------------------------
C--------------------------------------------------------------------
C S U B R O U T I N E I N I T H A R M O N I C P A R A M E T E R S
C--------------------------------------------------------------------
C jgf49.44: Initializes harmonic analysis arrays and parameters.
C Based on HACOLDS, written by rl.
C WJP 02.20.2018: Allowing for > 1 options for netcdf output
C--------------------------------------------------------------------
SUBROUTINE initHarmonicParameters()
USE SIZES, ONLY : READ_LOCAL_HOT_START_FILES, MNPROC, myProc,
& WRITE_LOCAL_HARM_FILES,
& WRITE_LOCAL_HOT_START_FILES
USE GLOBAL, ONLY : ITHS, NSTAE, NSTAV, NHSTAR
USE MESH, ONLY : NP
IMPLICIT NONE
INTEGER :: i, j, n ! loop counters
C
call setMessageSource("initHarmonicParameters")
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Enter.")
#endif
IF (IHARIND.EQ.0) THEN
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Return.")
#endif
call unsetMessageSource()
RETURN ! EARLY RETURN if harmonic analysis is not part of this run
ENDIF
!
! jgf51.44: If subdomain harmonic analysis files will be written
! in parallel, subdomain hotstart files are also required for
! saving the load vectors and rhs of the harmonic analysis problem.
if ( (mnproc.gt.1).and.
& (write_local_harm_files.eqv..true.) ) then
call logMessage(INFO,'Subdomain harmonic analysis files '
& // 'were specified with -m on the command line.')
if ( (NHSTAR.ne.0).and.
& (write_local_hot_start_files.eqv..false.) ) then
call logMessage(INFO,'Fulldomain hotstart files will be '
& // 'created by default.')
call allMessage(ERROR,'Fulldomain hotstart output files '
& // 'cannot be used in combination with subdomain harmonic '
& // 'analysis output files. '
& // 'The -s argument can be used on the ADCIRC command '
& // 'line to specify that hotstart files should be written '
& // 'for each subdomain; this is required for the writing '
& // 'of harmonic analysis output (fort.51 etc) '
& // ' to subdomains.')
call harmonicTerminate()
endif
if ( (NHSTAR.ne.0).and.
& (write_local_hot_start_files.eqv..true.) ) then
call logMessage(INFO,'Subdomain hotstart files will '
& // 'contain subdomain harmonic analysis data.')
endif
endif
IHABEG = ITHAS + NHAINC
ICHA=0
icall=0
C
if (hafreq(1).eq.0.0) then
nz=0
nf=1
else
nz=1
nf=0
endif
c
nfreq=nfreq-nf
mm=2*nfreq+nf
C
ha(:,:)=0.0d0
C
IF (NHASE.ne.0) THEN
STAELV(:,:)=0.d0
ENDIF
IF (NHASV.ne.0) THEN
STAULV(:,:)=0.d0
STAVLV(:,:)=0.d0
ENDIF
IF (NHAGE.ne.0) THEN
GLOELV(:,:)=0.d0
ENDIF
IF (NHAGV.ne.0) THEN
GLOULV(:,:)=0.d0
GLOVLV(:,:)=0.d0
ENDIF
IF (CHARMV.eqv..true.) THEN
ELAV(:)=0.D0
XVELAV(:)=0.D0
YVELAV(:)=0.D0
ELVA(:)=0.D0
XVELVA(:)=0.D0
YVELVA(:)=0.D0
ENDIF ! charmv
C
C jgf49.44: Initialize all the full domain i/o arrays to zero, if
C they will be used (i.e., if we are in serial or if we might be
C reading or writing full domain arrays in parallel).
if ( (WRITE_LOCAL_HARM_FILES.eqv..false.).and.
& (READ_LOCAL_HOT_START_FILES.eqv..false.).and.
& (MNPROC.gt.1).and.
& (myProc.eq.0) ) then
GLOELV_g(:,:) = 0.d0
STAELV_g(:,:) = 0.d0
GLOULV_g(:,:) = 0.d0
GLOVLV_g(:,:) = 0.d0
STAULV_g(:,:) = 0.d0
STAVLV_g(:,:) = 0.d0
IF (CHARMV.eqv..true.) THEN
ELAV_g(:) = 0.d0
ELVA_g(:) = 0.d0
XVELAV_g(:) = 0.d0
YVELAV_g(:) = 0.d0
XVELVA_g(:) = 0.d0
YVELVA_g(:) = 0.d0
ENDIF
endif
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Return.")
#endif
call unsetMessageSource()
C--------------------------------------------------------------------
end subroutine initHarmonicParameters
C--------------------------------------------------------------------
C--------------------------------------------------------------------
C S U B R O U T I N E U P D A T E H A R M O N I C A N A L Y S I S
C--------------------------------------------------------------------
C jgf49.44: Removed from timestep.F and placed here.
C WJP 02.20.2018: Allowing for > 1 options for netcdf output
C...
C... IF IHARIND=1 AND THE TIME STEP FALLS WITHIN THE SPECIFIED WINDOW
C... AND ON THE SPECIFIED INCREMENT, USE MODEL RESULTS TO UPDATE
C... HARMONIC ANALYSIS MATRIX AND LOAD VECTORS. NOTE: AN 8 BYTE RECORD
C... SHOULD BE USED THROUGHOUT THE HARMONIC ANALYSIS SUBROUTINES, EVEN
C... ON 32 BIT WORKSTATIONS, SINCE IN THAT CASE THE HARMONIC ANALYSIS
C... IS DONE IN DOUBLE PRECISION.
C...
C--------------------------------------------------------------------
SUBROUTINE updateHarmonicAnalysis(it, timeh)
USE GLOBAL, ONLY : ET00, UU00, VV00, ETA2, UU2, VV2,
& STAIE1, STAIE2, STAIE3, STAIV1, STAIV2, STAIV3,
& NSTAE, NSTAV, NNE, NNV
USE MESH, ONLY : NP, NM, UVECTMP, VVECTMP, DRVMAP2DSPVEC
IMPLICIT NONE
INTEGER, intent(in) :: it ! current ADCIRC timestep
REAL(8), intent(in) :: timeh ! time (sec) incl. STATIM and REFTIM
C
INTEGER :: i, j ! loop counters
REAL(SZ) :: EE1, EE2, EE3 ! nodal elevations around an element
REAL(SZ) :: U11, U22, U33 ! nodal u velocities around an element
REAL(SZ) :: V11, V22, V33 ! nodal v velocities around an element
INTEGER :: n ! loop counter
C
call setMessageSource("updateHarmonicAnalysis")
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Enter.")
#endif
C...
IF(IHARIND.EQ.1) THEN
IF((IT.GT.ITHAS).AND.(IT.LE.ITHAF)) THEN
IF(ICHA.EQ.NHAINC) ICHA=0
ICHA=ICHA+1
IF(ICHA.EQ.NHAINC) THEN
C...
C.....UPDATE THE LHS MATRIX
C...
C.... BEG: DW 2019
IF ( (NHASV .EQ. 1) .OR. (NHAGV .EQ. 1) ) THEN
IF ( IFSPROTS .EQ. 1 ) THEN
UVECTMP = UU2 ;
VVECTMP = VV2 ;
CALL DRVMAP2DSPVEC( UU2, VV2,
& UVECTMP, VVECTMP, NP, FWD = .FALSE. )
END IF
END IF
C.... DW
CALL LSQUPDLHS(timeh,IT)
C... IF DESIRED COMPUTE ELEVATION STATION INFORMATION AND UPDATE LOAD
C.....VECTOR
C...
IF(NHASE.ne.0) THEN
DO I=1,NSTAE
EE1=ETA2(NM(NNE(I),1))
EE2=ETA2(NM(NNE(I),2))
EE3=ETA2(NM(NNE(I),3))
ET00(I)=EE1*STAIE1(I)+EE2*STAIE2(I)+EE3*STAIE3(I)
END DO
CALL LSQUPDES(ET00,NSTAE)
ENDIF
C... IF DESIRED COMPUTE VELOCITY STATION INFORMATION AND UPDATE LOAD
C.....VECTOR
C...
IF(NHASV.ne.0) THEN
DO I=1,NSTAV
U11=UU2(NM(NNV(I),1))
U22=UU2(NM(NNV(I),2))
U33=UU2(NM(NNV(I),3))
V11=VV2(NM(NNV(I),1))
V22=VV2(NM(NNV(I),2))
V33=VV2(NM(NNV(I),3))
UU00(I)=U11*STAIV1(I)+U22*STAIV2(I)+U33*STAIV3(I)
VV00(I)=V11*STAIV1(I)+V22*STAIV2(I)+V33*STAIV3(I)
END DO
CALL LSQUPDVS(UU00,VV00,NSTAV)
ENDIF
C...
C.....IF DESIRED UPDATE GLOBAL ELEVATION LOAD VECTOR
C...
IF(NHAGE.ne.0) CALL LSQUPDEG(ETA2,NP)
C...
C.....IF DESIRED UPDATE GLOBAL VELOCITY LOAD VECTOR
C...
IF(NHAGV.ne.0) CALL LSQUPDVG(UU2,VV2,NP)
C.... BEG: DW 2019
IF ( (NHASV .EQ. 1) .OR. (NHAGV .EQ. 1) ) THEN
IF ( IFSPROTS .EQ. 1 ) THEN
UU2 = UVECTMP ;
VV2 = VVECTMP ;
END IF
END IF
C.... DW
ENDIF
ENDIF
C... LINES TO COMPUTE MEANS AND VARIANCES
if (CHARMV) then
IF(IT.GT.ITMV) THEN
NTSTEPS=NTSTEPS+1
DO I=1,NP
ELAV(I)=ELAV(I)+ETA2(I)
XVELAV(I)=XVELAV(I)+UU2(I)
YVELAV(I)=YVELAV(I)+VV2(I)
ELVA(I)=ELVA(I)+ETA2(I)*ETA2(I)
XVELVA(I)=XVELVA(I)+UU2(I)*UU2(I)
YVELVA(I)=YVELVA(I)+VV2(I)*VV2(I)
END DO
ENDIF
endif ! charmv
ENDIF
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Return.")
#endif
call unsetMessageSource()
C--------------------------------------------------------------------
end subroutine updateHarmonicAnalysis
C--------------------------------------------------------------------
c***********************************************************************
c Subroutine to update the Left Hand Side Matrix *
c *
c TIMELOC - ABSOLUTE MODEL TIME (SEC) *
c IT - MODEL TIME STEP *
c icall - number of times the subroutine has been called *
c a - Left Hand Side Matrix *
c *
c RL 11/7/95 *
c***********************************************************************
c
SUBROUTINE LSQUPDLHS(TIMELOC,IT)
IMPLICIT NONE
INTEGER IT,I,J,I1,I2,J1,J2
REAL(SZ) TF1,TF2
REAL(8) TIMELOC
C
call setMessageSource("LSQUPDLHS")
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Enter.")
#endif
c
icall = icall + 1
c
c***** Update the Left Hand Side Matrix
c Note: this is a symmetric matrix and therefore only store the
c upper triangular part. The lower part will be filled out in
c SUBROUTINE FULSOL prior to the matrix's decomposition
c Take care of the steady constituent if included in the analysis
if(nf.eq.1) then
ha(1,1)=icall
do j=1,nfreq
tf1=hafreq(j+nf)*TIMELOC
ha(1,2*j) = ha(1,2*j) + cos(tf1)
ha(1,2*j+1) = ha(1,2*j+1) + sin(tf1)
end do
endif
c Take care of the other constituents
do i=1,nfreq
do j=i,nfreq
i1=2*i-(1-nf)
i2=i1+1
j1=2*j-(1-nf)
j2=j1+1
tf1=hafreq(i+nf)*TIMELOC
tf2=hafreq(j+nf)*TIMELOC
ha(i1,j1) = ha(i1,j1) + cos(tf1)*cos(tf2)
ha(i1,j2) = ha(i1,j2) + cos(tf1)*sin(tf2)
ha(i2,j2) = ha(i2,j2) + sin(tf1)*sin(tf2)
if(i2.le.j1) ha(i2,j1) = ha(i2,j1) + sin(tf1)*cos(tf2)
end do
end do
c Record update time and time step
TIMEUD = TIMELOC
ITUD = IT
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Return.")
#endif
call unsetMessageSource()
return
end subroutine LSQUPDLHS
c***********************************************************************
c Subroutine to update the Right Hand Side Load Vectors for the *
c elevation station harmonic analysis. *
c *
c STAE - STATION ELEVATION VALUES USED TO UPDATE LOAD VECTORS *
c NSTAE - NUMBER OF TIDAL ELEVATION RECORDING STATIONS *
c *
c STAELV - station elevation load vector *
c *
c RL 11/8/95 *
c***********************************************************************
c
SUBROUTINE LSQUPDES(STAE,NSTAE)
IMPLICIT NONE
INTEGER NSTAE,N,I,I1,I2
REAL(SZ) TF1,CTF1,STF1
REAL(SZ) STAE(MNSTAE)
C
call setMessageSource("LSQUPDES")
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Enter.")
#endif
c
c***** Update the Right Hand Side Load Vectors
c
c Take care of the steady constituent if included in the analysis
if(nz.eq.0) then
do n=1,NSTAE
STAELV(1,N) = STAELV(1,N) + STAE(N)
end do
endif
c Take care of the other constituents
do i=1,nfreq
i1=2*i-nz
i2=i1+1
tf1=hafreq(i+nf)*TIMEUD
ctf1 = cos(tf1)
stf1 = sin(tf1)
do n=1,NSTAE
STAELV(I1,N) = STAELV(I1,N) + STAE(N)*CTF1
STAELV(I2,N) = STAELV(I2,N) + STAE(N)*STF1
end do
end do
C
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Return.")
#endif
call unsetMessageSource()
return
end subroutine LSQUPDES
c***********************************************************************
c Subroutine to update the Right Hand Side Load Vectors for the *
c velocity station harmonic analysis. *
c *
c STAU - STATION U VELOCITY VALUES USED TO UPDATE LOAD VECTORS *
c STAV - STATION V VELOCITY VALUES USED TO UPDATE LOAD VECTORS *
c NSTAV - NUMBER OF TIDAL CURRENT RECORDING STATIONS *
c *
c STAULV - station u velocity load vector *
c STAVLV - station v velocity load vector *
c *
c RL 11/8/95 *
c***********************************************************************
c
SUBROUTINE LSQUPDVS(STAU,STAV,NSTAV)
IMPLICIT NONE
INTEGER NSTAV,N,I,I1,I2
REAL(SZ) TF1,CTF1,STF1
REAL(SZ) STAU(MNSTAV),STAV(MNSTAV)
C
call setMessageSource("LSQUPDVS")
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Enter.")
#endif
c
c***** Update the Right Hand Side Load Vectors
c
c Take care of the steady constituent if included in the analysis
if(nz.eq.0) then
do n=1,NSTAV
STAULV(1,N) = STAULV(1,N) + STAU(N)
STAVLV(1,N) = STAVLV(1,N) + STAV(N)
end do
endif
c Take care of the other constituents
do i=1,nfreq
i1=2*i-nz
i2=i1+1
tf1=hafreq(i+nf)*TIMEUD
ctf1 = cos(tf1)
stf1 = sin(tf1)
do n=1,NSTAV
STAULV(I1,N) = STAULV(I1,N) + STAU(N)*CTF1
STAVLV(I1,N) = STAVLV(I1,N) + STAV(N)*CTF1
STAULV(I2,N) = STAULV(I2,N) + STAU(N)*STF1
STAVLV(I2,N) = STAVLV(I2,N) + STAV(N)*STF1
end do
end do
C
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Return.")
#endif
call unsetMessageSource()
return
end subroutine LSQUPDVS
c***********************************************************************
c Subroutine to update the Right Hand Side Load Vectors for the *
c global elevation harmonic analysis. *
c *
c GLOE - GLOBAL ELEVATION VALUES USED TO UPDATE LOAD VECTORS *
c NP - NUMBER OF POINTS IN GLOBAL GRID *
c *
c GLOELV - global elevation load vector *
c *
c RL 11/8/95 *
c***********************************************************************
c
SUBROUTINE LSQUPDEG(GLOE,NP)
IMPLICIT NONE
INTEGER I,NP,N,I1,I2
REAL(SZ) TF1,CTF1,STF1
REAL(SZ) GLOE(MNP)
C
call setMessageSource("LSQUPDEG")
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Enter.")
#endif
c
c*****Update the Right Hand Side Load Vectors
c
c Take care of the steady constituent if included in the analysis
if(nz.eq.0) then
do n=1,np
GLOELV(1,N)=GLOELV(1,N)+GLOE(N)
end do
endif
c Take care of the other constituents
do i=1,nfreq
i1=2*i-nz
i2=i1+1
tf1=hafreq(i+nf)*TIMEUD
ctf1 = cos(tf1)
stf1 = sin(tf1)
do n=1,np
GLOELV(I1,N)=GLOELV(I1,N)+GLOE(N)*CTF1
GLOELV(I2,N)=GLOELV(I2,N)+GLOE(N)*STF1
end do
end do
C
#if defined(HARM_TRACE) || defined(ALL_TRACE)
call allMessage(DEBUG,"Return.")
#endif
call unsetMessageSource()
return
end subroutine LSQUPDEG
c***********************************************************************
c Subroutine to update the Right Hand Side Load Vectors for the *
c global velocity harmonic analysis. *
c *
c GLOU - GLOBAL U VELOCITY VALUES USED TO UPDATE LOAD VECTORS *
c GLOV - GLOBAL V VELOCITY VALUES USED TO UPDATE LOAD VECTORS *