-
Notifications
You must be signed in to change notification settings - Fork 0
/
phys.f90
executable file
·1121 lines (999 loc) · 42.2 KB
/
phys.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
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! regional model module - phys
! usage: include all the subroutines
! and functions to calculate physical equations
! Yun Zhang 04/29/2015
! @stanford
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
module phys
use constant_parameter
use boundary
use turbulence
use source
use initialization
implicit none
contains
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! calculate_air_density
! usage: calculate air density for each cell
! each layer
! Yun Zhang 04/30/2015
! @stanford
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine calculate_air_density(rhoa_c,Pa_c,qv_c,Temp_c)
real(dp),dimension(NLAT,NLONG,NVERT),intent(in)::Pa_c,qv_c,Temp_c
real(dp),dimension(NLAT,NLONG,NVERT),intent(inout)::rhoa_c
rhoa_c=Pa_c/R_prime/(1+0.608*qv_c)/Temp_c
end subroutine calculate_air_density
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! calculate_mass_weighed_P
! usage: calculate P_c P_bot
! each layer
! Yun Zhang 04/30/2015
! @stanford
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine calculate_mass_weighed_P(Pa_bot,P_bot,P_c)
real(dp),dimension(NLAT,NLONG,0:NVERT),intent(in)::Pa_bot
real(dp),dimension(NLAT,NLONG,0:NVERT),intent(inout)::P_bot
real(dp),dimension(NLAT,NLONG,NVERT),intent(inout)::P_c
integer::i,j,k
do i=1,NLAT
do j=1,NLONG
P_bot(i,j,:)=(Pa_bot(i,j,:)/1000)**k_therm
P_c(i,j,:)=1.0_dp/(1+k_therm)*(P_bot(i,j,1:NVERT)*Pa_bot(i,j,1:NVERT)-Pa_bot(i,j,0:(NVERT-1))*Pa_bot(i,j,0:(NVERT-1)))
P_c(i,j,:)=P_c(i,j,:)/(Pa_bot(i,j,1:NVERT)-Pa_bot(i,j,0:NVERT-1))
enddo
enddo
end subroutine calculate_mass_weighed_P
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! calculate_fluxface_column_pressure
! usage: calculate column pressure at each edge
! of each cell. The values will be used to calculate
! cell-centered column pressure
! Yun Zhang 04/30/2015
! @stanford
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine calculate_fluxface_column_pressure(pi_c_old,pi_uface,pi_vface)
integer::i,j
real(dp), dimension(NLAT,NLONG), intent(in)::pi_c_old
real(dp), dimension(NLAT,NLONG+1), intent(inout):: pi_uface
real(dp), dimension(NLAT+1,NLONG), intent(inout):: pi_vface
! calculate uface
! use central differencing
do i=1,NLAT
if(periodicBC==1) then
pi_uface(i,1)=0.5*(pi_c_old(i,1)+pi_c_old(i,NLONG))
pi_uface(i,NLONG+1)=0.5*(pi_c_old(i,1)+pi_c_old(i,NLONG))
else
pi_uface(i,1)=pi_c_old(i,1)
pi_uface(i,NLONG+1)=pi_c_old(i,NLONG)
endif
pi_uface(i,2:NLONG)=0.5*(pi_c_old(i,1:(NLONG-1))+pi_c_old(i,2:NLONG))
enddo
! calculate vface
! use central differencing
do i=1,NLONG
if(periodicBC==1) then
pi_vface(1,i)=0.5*(pi_c_old(1,i)+pi_c_old(NLAT,i))
pi_vface(NLAT+1,i)=0.5*(pi_c_old(1,i)+pi_c_old(NLAT,i))
else
pi_vface(1,i)=pi_c_old(1,i)
pi_vface(NLAT+1,i)=pi_c_old(NLAT,i)
endif
pi_vface(2:NLAT,i)=0.5*(pi_c_old(1:(NLAT-1),NLONG)+pi_c_old(2:NLAT,NLONG))
enddo
end subroutine calculate_fluxface_column_pressure
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! calculate_fluxsum
! usage:calculate fluxsum for each cell and each layer
! fluxsum(i,j,k) means the sum of flux of 4 edges
! from layer 1 to layer k at cell i,j
! and also calculate flux_uface and flux_vface
! Yun Zhang 04/30/2015
! @stanford
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine calculate_fluxsum(fluxsum,flux_uface,flux_vface,u,v,pi_uface,pi_vface,lat_c,lat_vface,dsigma)
integer:: k
real(dp), dimension(NVERT),intent(in):: dsigma
real(dp), dimension(NLAT,NLONG,NVERT),intent(inout):: fluxsum
real(dp), dimension(NLAT,NLONG+1,NVERT), intent(in):: u
real(dp), dimension(NLAT+1,NLONG,NVERT), intent(in):: v
real(dp), dimension(NLAT,NLONG+1), intent(in):: pi_uface
real(dp), dimension(NLAT+1,NLONG), intent(in):: pi_vface,lat_vface
real(dp), dimension(NLAT,NLONG),intent(in):: lat_c
real(dp), dimension(NLAT,NLONG+1,NVERT),intent(inout):: flux_uface
real(dp), dimension(NLAT+1,NLONG,NVERT),intent(inout):: flux_vface
! reset all flux as zero
fluxsum=0.0_dp
flux_vface=0.0_dp
flux_uface=0.0_dp
! set for the first layer 7.15 7.16
do k=1,NVERT
flux_uface(:,:,k)=u(:,:,k)*pi_uface*Re*dphi
flux_vface(:,:,k)=v(:,:,k)*pi_vface*Re*dlamda_e*cos(lat_vface)
enddo
k=1
fluxsum(:,:,k)=(flux_uface(:,2:(NLONG+1),k)-flux_uface(:,1:NLONG,k))*dsigma(k)
fluxsum(:,:,k)=fluxsum(:,:,k)+(flux_vface(2:(NLAT+1),:,k)-flux_vface(1:NLAT,:,k))*dsigma(k)
! for other layers
do k=2,NVERT
fluxsum(:,:,k)=fluxsum(:,:,k-1)+(flux_uface(:,2:(NLONG+1),k)-flux_uface(:,1:NLONG,k))*dsigma(k)
fluxsum(:,:,k)=fluxsum(:,:,k)+(flux_vface(2:(NLAT+1),:,k)-flux_vface(1:NLAT,:,k))*dsigma(k)
enddo
!print *, fluxsum(1,1,NVERT),flux_uface(1,1,1),flux_uface(1,2,1),flux_vface(1,1,1),flux_vface(2,1,1)
end subroutine calculate_fluxsum
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! calculate_column_pressure
! usage: update column pressure for each cell
! and call by regional model every time step
! Yun Zhang 04/30/2015
! @stanford
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine calculate_column_pressure(pi_c_old,pi_c_new,lat_c,fluxsum)
integer:: i,j,k
real(dp), dimension(NLAT,NLONG,NVERT),intent(in):: fluxsum
real(dp), dimension(NLAT,NLONG),intent(inout):: pi_c_new
real(dp), dimension(NLAT,NLONG),intent(in):: pi_c_old
real(dp), dimension(NLAT,NLONG),intent(in):: lat_c
! calculate new column pressure
pi_c_new=pi_c_old-dt/(Re**2)/dphi/dlamda_e/cos(lat_c)*fluxsum(:,:,NVERT)
end subroutine calculate_column_pressure
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! calculate_w_sigma
! usage: update w_sigma after each time step
! for each layer using continuity equation for
! each layer
! Yun Zhang 04/30/2015
! @stanford
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine calculate_w_sigma(w_sigma, pi_c_old,pi_c_new,lat_c,fluxsum,sigma_bot,dsigma)
integer:: i,j,k
real(dp),dimension(NLAT,NLONG,0:NVERT),intent(inout)::w_sigma
real(dp),dimension(NLAT,NLONG),intent(in):: pi_c_old,pi_c_new,lat_c
real(dp),dimension(NLAT,NLONG,NVERT),intent(in):: fluxsum
real(dp),dimension(NVERT),intent(in):: dsigma
real(dp),dimension(0:NVERT),intent(in):: sigma_bot
! assume the top and bottom w_sigma=0
w_sigma(:,:,0)=0.0_dp
w_sigma(:,:,NVERT)=0.0_dp
! calculate interior layers boundaries 7.21
do i=1,NLAT
do j=1,NLONG
w_sigma(i,j,1:NVERT)=-1.0_dp/(pi_c_new(i,j)*Re*Re*cos(lat_c(i,j))*dlamda_e*dphi)*fluxsum(i,j,1:NVERT)
w_sigma(i,j,1:NVERT)=w_sigma(i,j,1:NVERT)-sigma_bot(1:NVERT)*(pi_c_new(i,j)&
-pi_c_old(i,j))/dt/pi_c_new(i,j)
enddo
enddo
end subroutine calculate_w_sigma
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! calculate_fluxface_PVT
! usage: calculate the potential virtual temperature
! at flux face
! Yun Zhang 05/01/2015
! @STANFORD
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine calculate_fluxface_PVT(PVT_c_old,PVT_uface,PVT_vface,PVT_bot,&
P_c,P_bot,lat_uface,long_uface,lat_vface,long_vface,nt)
integer::i,j,k
integer,intent(in):: nt
real(dp), dimension(NLAT,NLONG,NVERT), intent(in):: PVT_c_old
real(dp), dimension(NLAT,NLONG+1,NVERT), intent(inout):: PVT_uface
real(dp), dimension(NLAT+1,NLONG,NVERT), intent(inout):: PVT_vface
real(dp),dimension(NLAT,NLONG,0:NVERT),intent(inout)::PVT_bot
real(dp),dimension(NLAT+1,NLONG),intent(in):: lat_vface, long_vface
real(dp),dimension(NLAT,NLONG+1),intent(in):: lat_uface, long_uface
real(dp), dimension(NLAT,NLONG,NVERT),intent(in)::P_c
real(dp), dimension(NLAT,NLONG,0:NVERT),intent(in)::P_bot
! calculate uface
! use central differencing
do i=1,NLAT
do k=1,NVERT
if(periodicBC==1) then
PVT_uface(i,1,k)=0.5*(PVT_c_old(i,1,k)+PVT_c_old(i,NLONG,k))
PVT_uface(i,NLONG+1,k)=0.5*(PVT_c_old(i,1,k)+PVT_c_old(i,NLONG,k))
else
PVT_uface(i,1,k)=PVT_c_old(i,1,k)
PVT_uface(i,NLONG+1,k)=PVT_c_old(i,NLONG,k)
endif
PVT_uface(i,2:NLONG,k)=0.5*(PVT_c_old(i,1:(NLONG-1),k)+PVT_c_old(i,2:NLONG,k))
enddo
enddo
! calculate vface
! use central differencing
do i=1,NLONG
do k=1,NVERT
if(periodicBC==1) then
PVT_vface(1,i,k)=0.5*(PVT_c_old(1,i,k)+PVT_c_old(NLAT,i,k))
PVT_vface(NLAT+1,i,k)=0.5*(PVT_c_old(1,i,k)+PVT_c_old(NLAT,i,k))
else
PVT_vface(1,i,k)=PVT_c_old(1,i,k)
PVT_vface(NLAT+1,i,k)=PVT_c_old(NLAT,i,k)
endif
PVT_vface(2:NLAT,i,k)=0.5*(PVT_c_old(1:(NLAT-1),i,k)+PVT_c_old(2:NLAT,i,k))
enddo
enddo
! set boundary PTV value
if(PTVbound==1) then
! vface value
do j=1,NLONG
do k=1,NVERT
PVT_vface(1,j,k)=boundary_potential_virtual_temp(lat_vface(1,j),long_vface(1,j),k,nt*dt)
PVT_vface(NLAT+1,j,k)=boundary_potential_virtual_temp(lat_vface(NLAT+1,j),long_vface(NLAT+1,j),k,nt*dt)
enddo
enddo
! uface value
do i=1,NLAT
do k=1,NVERT
PVT_uface(i,1,k)=boundary_potential_virtual_temp(lat_vface(i,1),long_uface(i,1),k,nt*dt)
PVT_uface(i,NLONG+1,k)=boundary_potential_virtual_temp(lat_uface(i,NLONG+1),long_uface(i,NLONG+1),k,nt*dt)
enddo
enddo
endif
! layer bottom eqn 7.11
do i=1,NLAT
do j=1,NLONG
! for zero 0, P_bot(i,j,0)=P_c(i,j,0)
PVT_bot(i,j,0)=PVT_c_old(i,j,1)
PVT_bot(i,j,1:(NVERT-1))=(P_bot(i,j,1:(NVERT-1))-P_c(i,j,1:(NVERT-1)))*PVT_c_old(i,j,1:(NVERT-1))
PVT_bot(i,j,1:(NVERT-1))=PVT_bot(i,j,1:(NVERT-1))+(P_c(i,j,2:NVERT)-P_bot(i,j,1:(NVERT-1)))*PVT_c_old(i,j,2:NVERT)
PVT_bot(i,j,1:(NVERT-1))=PVT_bot(i,j,1:(NVERT-1))/(P_c(i,j,2:NVERT)-P_c(i,j,1:(NVERT-1)))
PVT_bot(i,j,NVERT)=PVT_c_old(i,j,NVERT)
enddo
enddo
end subroutine calculate_fluxface_PVT
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! calculate_fluxface_qv
! usage: calculate specific humidity
! at flux face
! Yun Zhang 05/01/2015
! @STANFORD
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine calculate_fluxface_qv(qv_c_old,qv_uface,qv_vface,qv_bot,P_c,P_bot,&
lat_uface,long_uface,lat_vface,long_vface,nt)
integer::i,j,k
integer,intent(in):: nt
real(dp), dimension(NLAT,NLONG,NVERT), intent(in):: qv_c_old
real(dp), dimension(NLAT,NLONG+1,NVERT), intent(inout):: qv_uface
real(dp), dimension(NLAT+1,NLONG,NVERT), intent(inout):: qv_vface
real(dp), dimension(NLAT,NLONG,0:NVERT),intent(inout)::qv_bot
real(dp), dimension(NLAT+1,NLONG),intent(in):: lat_vface, long_vface
real(dp), dimension(NLAT,NLONG+1),intent(in):: lat_uface, long_uface
real(dp), dimension(NLAT,NLONG,NVERT),intent(in)::P_c
real(dp), dimension(NLAT,NLONG,0:NVERT),intent(in)::P_bot
! calculate uface
! use central differencing
do i=1,NLAT
do k=1,NVERT
if(periodicBC==1) then
qv_uface(i,1,k)=0.5*(qv_c_old(i,1,k)+qv_c_old(i,NLONG,k))
qv_uface(i,NLONG+1,k)=0.5*(qv_c_old(i,1,k)+qv_c_old(i,NLONG,k))
else
qv_uface(i,1,k)=qv_c_old(i,1,k)
qv_uface(i,NLONG+1,k)=qv_c_old(i,NLONG,k)
endif
qv_uface(i,2:NLONG,k)=0.5*(qv_c_old(i,1:(NLONG-1),k)+qv_c_old(i,2:NLONG,k))
enddo
enddo
! calculate vface
! use central differencing
do i=1,NLONG
do k=1,NVERT
if(periodicBC==1) then
qv_vface(1,i,k)=0.5*(qv_c_old(1,i,k)+qv_c_old(NLAT,i,k))
qv_vface(NLAT+1,i,k)=0.5*(qv_c_old(1,i,k)+qv_c_old(NLAT,i,k))
else
qv_vface(1,i,k)=qv_c_old(1,i,k)
qv_vface(NLAT+1,i,k)=qv_c_old(NLAT,i,k)
endif
qv_vface(2:NLAT,i,k)=0.5*(qv_c_old(1:(NLAT-1),i,k)+qv_c_old(2:NLAT,i,k))
enddo
enddo
! set boundary PTV value
if(qvbound==1) then
! vface value
do j=1,NLONG
do k=1,NVERT
qv_vface(1,j,k)=boundary_specific_humidity(lat_vface(1,j),long_vface(1,j),k,nt*dt)
qv_vface(NLAT+1,j,k)=boundary_specific_humidity(lat_vface(NLAT+1,j),long_vface(NLAT+1,j),k,nt*dt)
enddo
enddo
! uface value
do i=1,NLAT
do k=1,NVERT
qv_uface(i,1,k)=boundary_specific_humidity(lat_vface(i,1),long_uface(i,1),k,nt*dt)
qv_uface(i,NLONG+1,k)=boundary_specific_humidity(lat_uface(i,NLONG+1),long_uface(i,NLONG+1),k,nt*dt)
enddo
enddo
endif
! layer bottom 7.25
do i=1,NLAT
do j=1,NLONG
! for zero 0, P_bot(i,j,0)=P_c(i,j,0)
qv_bot(i,j,0)=qv_c_old(i,j,1)
do k=1,(NVERT-1)
if(qv_c_old(i,j,k)==qv_c_old(i,j,k+1)) then
qv_bot(i,j,k)=qv_c_old(i,j,k)
else if(qv_c_old(i,j,k)==0 .or. qv_c_old(i,j,k+1)==0) then
qv_bot(i,j,k)=0.5*(qv_c_old(i,j,k)+qv_c_old(i,j,k+1))
else
qv_bot(i,j,k)=(log(qv_c_old(i,j,k))-log(qv_c_old(i,j,k+1)))&
/(1/qv_c_old(i,j,k+1)-1/qv_c_old(i,j,k))
endif
enddo
qv_bot(i,j,NVERT)=qv_c_old(i,j,NVERT)
enddo
enddo
end subroutine calculate_fluxface_qv
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! calculate_fluxface_gas
! usage: calculate gas concentration
! at flux face
! Yun Zhang 05/01/2015
! @STANFORD
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine calculate_fluxface_gas(gas_c_old,gas_uface,gas_vface,gas_bot,P_c,P_bot,&
lat_uface,long_uface,lat_vface,long_vface,nt)
integer::i,j,k
integer,intent(in):: nt
real(dp), dimension(NLAT,NLONG,NVERT), intent(in):: gas_c_old
real(dp), dimension(NLAT,NLONG+1,NVERT), intent(inout):: gas_uface
real(dp), dimension(NLAT+1,NLONG,NVERT), intent(inout):: gas_vface
real(dp), dimension(NLAT,NLONG,0:NVERT),intent(inout)::gas_bot
real(dp), dimension(NLAT+1,NLONG),intent(in):: lat_vface, long_vface
real(dp), dimension(NLAT,NLONG+1),intent(in):: lat_uface, long_uface
real(dp), dimension(NLAT,NLONG,NVERT),intent(in)::P_c
real(dp), dimension(NLAT,NLONG,0:NVERT),intent(in)::P_bot
! calculate uface
! use central differencing
do i=1,NLAT
do k=1,NVERT
if(periodicBC==1) then
gas_uface(i,1,k)=0.5*(gas_c_old(i,1,k)+gas_c_old(i,NLONG,k))
gas_uface(i,NLONG+1,k)=0.5*(gas_c_old(i,1,k)+gas_c_old(i,NLONG,k))
else
gas_uface(i,1,k)=gas_c_old(i,1,k)
gas_uface(i,NLONG+1,k)=gas_c_old(i,NLONG,k)
endif
gas_uface(i,2:NLONG,k)=0.5*(gas_c_old(i,1:(NLONG-1),k)+gas_c_old(i,2:NLONG,k))
enddo
enddo
! calculate vface
! use central differencing
do i=1,NLONG
do k=1,NVERT
if(periodicBC==1) then
gas_vface(1,i,k)=0.5*(gas_c_old(1,i,k)+gas_c_old(NLAT,i,k))
gas_vface(NLAT+1,i,k)=0.5*(gas_c_old(1,i,k)+gas_c_old(NLAT,i,k))
else
gas_vface(1,i,k)=gas_c_old(1,i,k)
gas_vface(NLAT+1,i,k)=gas_c_old(NLAT,i,k)
endif
gas_vface(2:NLAT,i,k)=0.5*(gas_c_old(1:(NLAT-1),i,k)+gas_c_old(2:NLAT,i,k))
enddo
enddo
! set boundary PTV value
if(gasbound==1) then
! vface value
do j=1,NLONG
do k=1,NVERT
gas_vface(1,j,k)=boundary_gas(lat_vface(1,j),long_vface(1,j),k,nt*dt)
gas_vface(NLAT+1,j,k)=boundary_gas(lat_vface(NLAT+1,j),long_vface(NLAT+1,j),k,nt*dt)
enddo
enddo
! uface value
do i=1,NLAT
do k=1,NVERT
gas_uface(i,1,k)=boundary_gas(lat_vface(i,1),long_uface(i,1),k,nt*dt)
gas_uface(i,NLONG+1,k)=boundary_gas(lat_uface(i,NLONG+1),long_uface(i,NLONG+1),k,nt*dt)
enddo
enddo
endif
! layer bottom 7.25
do i=1,NLAT
do j=1,NLONG
gas_bot(i,j,0)=gas_c_old(i,j,1)
gas_bot(i,j,NVERT)=gas_c_old(i,j,NVERT)
do k=1,(NVERT-1)
gas_bot(i,j,k)=0.5*(gas_c_old(i,j,k)+gas_c_old(i,j,k+1))
enddo
enddo
enddo
end subroutine calculate_fluxface_gas
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! calculate_heat_source
! usage: calculate heat source value for each time
! step for each cell at each layer
! the heat source is treated as explicit
! and defined in the heat_source function
! in source.f90
! Yun Zhang 05/01/2015
! @stanford
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine calculate_heat_source(Pa_c,Q,lat_c,long_c,nt)
real(dp),dimension(NLAT,NLONG),intent(in)::lat_c,long_c
real(dp),dimension(NLAT,NLONG,NVERT),intent(inout):: Q
real(dp),dimension(NLAT,NLONG,NVERT),intent(in):: Pa_c
integer, intent(in):: nt
integer:: i,j,k
do i=1,NLAT
do j=1,NLONG
do k=1,NVERT
Q(i,j,k)=heat_source(lat_c(i,j),long_c(i,j),k,nt*dt)
Q(i,j,k)=((1000/Pa_c(i,j,k))**k_therm)/Cp_d*Q(i,j,k)
enddo
enddo
enddo
end subroutine calculate_heat_source
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! calculate_gas_source
! usage: calculate gas source value for each time
! step for each cell at each layer
! the heat source is treated as explicit
! and defined in the heat_source function
! in source.f90
! Yun Zhang 05/01/2015
! @stanford
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine calculate_gas_source(Pa_c,Q,lat_c,long_c,nt)
real(dp),dimension(NLAT,NLONG),intent(in)::lat_c,long_c
real(dp),dimension(NLAT,NLONG,NVERT),intent(inout):: Q
real(dp),dimension(NLAT,NLONG,NVERT),intent(in):: Pa_c
integer, intent(in):: nt
integer:: i,j,k
do i=1,NLAT
do j=1,NLONG
do k=1,NVERT
Q(i,j,k)=gas_source(lat_c(i,j),long_c(i,j),k,nt*dt)
enddo
enddo
enddo
end subroutine calculate_gas_source
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! calculate_gas_source
! usage: calculate gas source value for each time
! step for each cell at each layer
! the heat source is treated as explicit
! and defined in the heat_source function
! in source.f90
! Yun Zhang 05/01/2015
! @stanford
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine calculate_qv_source(Pa_c,Q,lat_c,long_c,nt)
real(dp),dimension(NLAT,NLONG),intent(in)::lat_c,long_c
real(dp),dimension(NLAT,NLONG,NVERT),intent(inout):: Q
real(dp),dimension(NLAT,NLONG,NVERT),intent(in):: Pa_c
integer, intent(in):: nt
integer:: i,j,k
do i=1,NLAT
do j=1,NLONG
do k=1,NVERT
Q(i,j,k)=qv_source(lat_c(i,j),long_c(i,j),k,nt*dt)
enddo
enddo
enddo
end subroutine calculate_qv_source
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! calculate_scalar_field
! usage: update scalar transport for each cell
! at every time step
! eqn 7.24 7.27
! use for PTV, qv and other components
! code use PTV as example
! Yun Zhang 04/30/2015
! @stanford
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine calculate_scalar_field(PVT_c_old, PVT_c_new, PVT_bot, &
pi_c_old,pi_c_new,PVT_uface, PVT_vface,flux_uface,flux_vface,lat_c, rhoa_c,&
dsigma, w_sigma_new,K_t,Pa_c,Q)
real(dp),dimension(NLAT,NLONG,NVERT),intent(in):: PVT_c_old, K_t, Q
real(dp),dimension(NLAT,NLONG,NVERT),intent(inout):: PVT_c_new
real(dp),dimension(NLAT,NLONG),intent(in)::pi_c_old, pi_c_new, lat_c
real(dp), dimension(NLAT+1,NLONG,NVERT),intent(in):: PVT_vface, flux_vface
real(dp), dimension(NLAT,NLONG+1,NVERT),intent(in):: PVT_uface, flux_uface
real(dp), dimension(NLAT,NLONG,0:NVERT),intent(in):: w_sigma_new
real(dp), dimension(NVERT),intent(in):: dsigma
real(dp), dimension(NLAT,NLONG,0:NVERT),intent(in):: PVT_bot
real(dp), dimension(NLAT,NLONG,NVERT),intent(in):: Pa_c, rhoa_c
integer:: i,j,k
! old value
do k=1,NVERT
PVT_c_new(:,:,k)=pi_c_old*PVT_c_old(:,:,k)/pi_c_new
enddo
! horizontal flux at boundary
do k=1,NVERT
PVT_c_new(:,:,k)=PVT_c_new(:,:,k)+dt/(pi_c_new*Re*Re*cos(lat_c)*dlamda_e*dphi)&
*(flux_uface(:,1:NLONG,k)*PVT_uface(:,1:NLONG,k)-flux_uface(:,2:(NLONG+1),k)*PVT_uface(:,2:(NLONG+1),k)&
+flux_vface(1:NLAT,:,k)*PVT_vface(1:NLAT,:,k)-flux_vface(2:(NLAT+1),:,k)*PVT_vface(2:(NLAT+1),:,k))
enddo
! vertical flux
do k=1,NVERT
PVT_c_new(:,:,k)=PVT_c_new(:,:,k)+dt/dsigma(k)*(w_sigma_new(:,:,k-1)*PVT_bot(:,:,k-1)&
-w_sigma_new(:,:,k)*PVT_bot(:,:,k))
enddo
! turbulence vertical flux
! not finished
if(turbmodel==1) then
do k=1,NVERT
PVT_c_new(:,:,k)=PVT_c_new(:,:,k)
enddo
endif
! heat source
if(heatsource==1) then
do k=1,NVERT
PVT_c_new(:,:,k)=PVT_c_new(:,:,k)+dt*pi_c_old/pi_c_new*Q(:,:,k)
enddo
endif
end subroutine calculate_scalar_field
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! calculate_geopotential
! usage: calculate geopotential for each cell center
! and layer bottom
! need to specify the boundary geopotential at bottom
! use boundary_surf_geopotential
! Yun 05/02/2015
! @stanford
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine calculate_geopotential(geopot_bot,geopot_c,geopot_c_tminus1,P_c,P_bot,PVT_c_old,lat_c,long_c,nt)
real(dp), dimension(NLAT,NLONG,NVERT),intent(in)::P_c,PVT_c_old
real(dp), dimension(NLAT,NLONG,NVERT),intent(inout)::geopot_c,geopot_c_tminus1
real(dp), dimension(NLAT,NLONG,0:NVERT),intent(inout):: geopot_bot
real(dp), dimension(NLAT,NLONG,0:NVERT),intent(in):: P_bot
real(dp), dimension(NLAT,NLONG), intent(in):: lat_c,long_c
integer,intent(in):: nt
integer::i,j,k
if(nt>0) geopot_c_tminus1=geopot_c
do i=1,NLAT
do j=1,NLONG
geopot_bot(i,j,NVERT)=boundary_surf_geopotential(lat_c(i,j),long_c(i,j),nt*dt)
geopot_c(i,j,NVERT)=geopot_bot(i,j,NVERT)-Cp_d*(PVT_c_old(i,j,NVERT)*(P_c(i,j,NVERT)-P_bot(i,j,NVERT)))
do k=NVERT-1,0,-1
geopot_bot(i,j,k)=geopot_c(i,j,k+1)-Cp_d*(PVT_c_old(i,j,k+1)*(P_bot(i,j,k)-P_c(i,j,k+1)))
if(k/=0) then
geopot_c(i,j,k)=geopot_bot(i,j,k)-Cp_d*(PVT_c_old(i,j,k)*(P_c(i,j,k)-P_bot(i,j,k)))
endif
enddo
enddo
enddo
if(nt==0) geopot_c_tminus1=geopot_c
end subroutine calculate_geopotential
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! calculate_velocity_field
! usage: update velocity field for each flux face
! Yun Zhang 05/02/2015
! @stanford
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine calculate_velocity_field(pi_c_old,pi_c_new,pi_c_tmp,pi_c_tminus1,lat_c,lat_vface,&
u_old,u_new,u_tmp,v_old,v_new,v_tmp,flux_uface,flux_vface,dsigma,&
sigma_bot,w_sigma_new,f_c,geopot_c,geopot_c_tminus1,&
PVT_c_old,P_bot,P_c,nu_t,rhoa_c)
real(dp),dimension(NLAT,NLONG),intent(in)::pi_c_old,pi_c_new,pi_c_tmp,pi_c_tminus1,lat_c,f_c
real(dp),dimension(NLAT,NLONG+1,NVERT),intent(in)::u_old,flux_uface,u_tmp
real(dp),dimension(NLAT,NLONG+1,NVERT),intent(inout)::u_new
real(dp),dimension(NLAT+1,NLONG,NVERT),intent(in)::v_old,flux_vface,v_tmp
real(dp),dimension(NLAT+1,NLONG,NVERT),intent(inout)::v_new
real(dp),dimension(0:NVERT),intent(in)::sigma_bot
real(dp),dimension(NVERT),intent(in)::dsigma
real(dp),dimension(NLAT,NLONG,0:NVERT),intent(in)::w_sigma_new,P_bot
real(dp),dimension(NLAT,NLONG,NVERT),intent(in)::geopot_c,geopot_c_tminus1,PVT_c_old,P_c,nu_t,rhoa_c
real(dp),dimension(NLAT+1,NLONG),intent(in)::lat_vface
real(dp),dimension(NLAT,NLONG+1,0:NVERT)::u_bot_old
real(dp),dimension(NLAT+1,NLONG,0:NVERT)::v_bot_old
! u field
real(dp),dimension(NLAT,NLONG+1)::A_old,A_new
real(dp),dimension(NLAT,NLONG+1,0:NVERT)::F
real(dp),dimension(NLAT,0:NLONG+1,NVERT)::B
real(dp),dimension(NLAT+1,NLONG+1,NVERT)::C
real(dp),dimension(NLAT+1,0:NLONG+1,NVERT)::D,E
! v field
real(dp),dimension(NLAT+1,NLONG)::P_old,P_new
real(dp),dimension(NLAT+1,NLONG,0:NVERT)::O
real(dp),dimension(0:NLAT+1,NLONG,NVERT)::R
real(dp),dimension(NLAT+1,NLONG+1,NVERT)::Q
real(dp),dimension(0:NLAT+1,NLONG+1,NVERT)::S,T
integer::i,j,k
! prepare for the boundary values
! use new virtual matrix to calculate velocity to simplify equation
real(dp),dimension(0:NLAT+1,0:NLONG+1)::v_f_c,v_pi_c_new,v_pi_c_old,v_lat_c
real(dp),dimension(0:NLAT+1,0:NLONG+1,0:NVERT)::v_F_new
real(dp),dimension(0:NLAT+1,0:NLONG+1,0:NVERT)::v_w_sigma_new
!real(dp),dimension(0:NLAT+1,0:NLONG+1,NVERT):: v_PVT_c_old, v_P_c, v_geopot_c
!real(dp),dimension(0:NLAT+1,0:NLONG+1,0:NVERT):: v_P_bot
real(dp),dimension(0:NLAT+1,0:NLONG+2,NVERT)::v_u_old, v_flux_uface
real(dp),dimension(0:NLAT+2,0:NLONG+1,NVERT)::v_v_old, v_flux_vface
v_f_c(1:NLAT,1:NLONG)=f_c
v_pi_c_old(1:NLAT,1:NLONG)=pi_c_old
v_pi_c_new(1:NLAT,1:NLONG)=pi_c_new
v_lat_c(1:NLAT,1:NLONG)=lat_c
v_w_sigma_new(1:NLAT,1:NLONG,0:NVERT)=w_sigma_new
if(periodicBC==1) then
v_pi_c_old(0,1:NLONG)=pi_c_old(NLAT,:)
v_pi_c_old(NLAT+1,1:NLONG)=pi_c_old(1,:)
v_pi_c_old(1:NLAT,0)=pi_c_old(:,NLONG)
v_pi_c_old(1:NLAT,NLONG+1)=pi_c_old(:,1)
v_pi_c_old(0,0)=0.5*(pi_c_old(1,NLONG)+pi_c_old(NLAT,1))
v_pi_c_old(NLAT,0)=0.5*(pi_c_old(1,NLONG)+pi_c_old(1,1))
v_pi_c_old(0,NLONG+1)=0.5*(pi_c_old(1,1)+pi_c_old(NLAT,NLONG))
v_pi_c_old(NLAT,NLONG+1)=0.5*(pi_c_old(1,NLONG)+pi_c_old(NLAT,1))
else
v_pi_c_old(0,1:NLONG)=pi_c_old(1,:)
v_pi_c_old(NLAT+1,1:NLONG)=pi_c_old(NLAT,:)
v_pi_c_old(:,0)=v_pi_c_old(:,1)
v_pi_c_old(:,NLONG+1)=v_pi_c_old(:,NLONG)
endif
if(periodicBC==1) then
v_pi_c_new(0,1:NLONG)=pi_c_new(NLAT,:)
v_pi_c_new(NLAT+1,1:NLONG)=pi_c_new(1,:)
v_pi_c_new(1:NLAT,0)=pi_c_new(:,NLONG)
v_pi_c_new(1:NLAT,NLONG+1)=pi_c_new(:,1)
v_pi_c_new(0,0)=0.5*(pi_c_new(1,NLONG)+pi_c_new(NLAT,1))
v_pi_c_new(NLAT,0)=0.5*(pi_c_new(1,NLONG)+pi_c_new(1,1))
v_pi_c_new(0,NLONG+1)=0.5*(pi_c_new(1,1)+pi_c_new(NLAT,NLONG))
v_pi_c_new(NLAT,NLONG+1)=0.5*(pi_c_new(1,NLONG)+pi_c_new(NLAT,1))
else
v_pi_c_new(0,1:NLONG)=pi_c_new(1,:)
v_pi_c_new(NLAT+1,1:NLONG)=pi_c_new(NLAT,:)
v_pi_c_new(:,0)=v_pi_c_new(:,1)
v_pi_c_new(:,NLONG+1)=v_pi_c_new(:,NLONG)
endif
v_lat_c(0,:)=lat_c(1,1)-dphi
v_lat_c(NLAT+1,:)=lat_c(NLAT,1)+dphi
v_lat_c(1:NLAT,0)=lat_c(1:NLAT,1)
v_lat_c(1:NLAT,NLONG+1)=lat_c(1:NLAT,NLONG)
if(coriolis==1) then
v_f_c=2*Omega*sin(v_lat_c)
else
v_f_c=0.0_dp
endif
if(periodicBC==1) then
v_w_sigma_new(0,1:NLONG,:)=w_sigma_new(NLAT,:,:)
v_w_sigma_new(NLAT+1,1:NLONG,:)=w_sigma_new(1,:,:)
v_w_sigma_new(1:NLAT,0,:)=w_sigma_new(:,NLONG,:)
v_w_sigma_new(1:NLAT,NLONG+1,:)=w_sigma_new(:,1,:)
v_w_sigma_new(0,0,:)=0.5*(w_sigma_new(1,NLONG,:)+w_sigma_new(NLAT,1,:))
v_w_sigma_new(NLAT,0,:)=0.5*(w_sigma_new(1,NLONG,:)+w_sigma_new(1,1,:))
v_w_sigma_new(0,NLONG+1,:)=0.5*(w_sigma_new(1,1,:)+w_sigma_new(NLAT,NLONG,:))
v_w_sigma_new(NLAT,NLONG+1,:)=0.5*(w_sigma_new(1,NLONG,:)+w_sigma_new(NLAT,1,:))
else
v_w_sigma_new(0,1:NLONG,:)=w_sigma_new(1,:,:)
v_w_sigma_new(NLAT+1,1:NLONG,:)=w_sigma_new(NLAT,:,:)
v_w_sigma_new(:,0,:)=v_w_sigma_new(:,1,:)
v_w_sigma_new(:,NLONG+1,:)=v_w_sigma_new(:,NLONG,:)
endif
v_u_old(1:NLAT,1:NLONG+1,:)=u_tmp
v_flux_uface(1:NLAT,1:NLONG+1,:)=flux_uface
v_v_old(1:NLAT+1,1:NLONG,:)=v_tmp
v_flux_vface(1:NLAT+1,1:NLONG,:)=flux_vface
v_u_old(0,1:NLONG+1,:)=u_tmp(1,:,:)
v_u_old(NLAT+1,1:NLONG+1,:)=u_tmp(NLAT,:,:)
v_u_old(:,0,:)=v_u_old(:,1,:)
v_u_old(:,NLONG+2,:)=v_u_old(:,NLONG+1,:)
v_flux_uface(0,1:NLONG+1,:)=flux_uface(1,:,:)
v_flux_uface(NLAT+1,1:NLONG+1,:)=flux_uface(NLAT,:,:)
v_flux_uface(:,0,:)=v_flux_uface(:,1,:)
v_flux_uface(:,NLONG+2,:)=v_flux_uface(:,NLONG+1,:)
v_v_old(1:NLAT+1,0,:)=v_tmp(:,1,:)
v_v_old(1:NLAT+1,NLONG+1,:)=v_tmp(:,NLONG,:)
v_v_old(0,:,:)=v_v_old(1,:,:)
v_v_old(NLAT+2,:,:)=v_v_old(NLAT+1,:,:)
v_flux_vface(1:NLAT+1,0,:)=flux_vface(:,1,:)
v_flux_vface(1:NLAT+1,NLONG+1,:)=flux_vface(:,NLONG,:)
v_flux_vface(0,:,:)=v_flux_vface(1,:,:)
v_flux_vface(NLAT+2,:,:)=v_flux_vface(NLAT+1,:,:)
! Update u field
! calculate column pressure multiplied by grid-cell area 7.38
! interior points
A_old=1.0_dp/8.0_dp*Re*Re*dphi*dlamda_e*(v_pi_c_old(2:NLAT+1,0:NLONG)*cos(v_lat_c(2:NLAT+1,0:NLONG))&
+v_pi_c_old(2:NLAT+1,1:NLONG+1)*cos(v_lat_c(2:NLAT+1,1:NLONG+1))+2*v_pi_c_old(1:NLAT,0:NLONG)*cos(v_lat_c(1:NLAT,0:NLONG))&
+2*v_pi_c_old(1:NLAT,1:NLONG+1)*cos(v_lat_c(1:NLAT,1:NLONG+1))+v_pi_c_old(0:NLAT-1,0:NLONG)*cos(v_lat_c(0:NLAT-1,0:NLONG))&
+v_pi_c_old(0:NLAT-1,1:NLONG+1)*cos(v_lat_c(0:NLAT-1,1:NLONG+1)))
A_new=1.0_dp/8.0_dp*Re*Re*dphi*dlamda_e*(v_pi_c_new(2:NLAT+1,0:NLONG)*cos(v_lat_c(2:NLAT+1,0:NLONG))&
+v_pi_c_new(2:NLAT+1,1:NLONG+1)*cos(v_lat_c(2:NLAT+1,1:NLONG+1))&
+2*v_pi_c_new(1:NLAT,0:NLONG)*cos(v_lat_c(1:NLAT,0:NLONG))&
+2*v_pi_c_new(1:NLAT,1:NLONG+1)*cos(v_lat_c(1:NLAT,1:NLONG+1))&
+v_pi_c_new(0:NLAT-1,0:NLONG)*cos(v_lat_c(0:NLAT-1,0:NLONG))&
+v_pi_c_new(0:NLAT-1,1:NLONG+1)*cos(v_lat_c(0:NLAT-1,1:NLONG+1)))
! equation 7.41
B=1.0_dp/12.0_dp*(v_flux_uface(0:NLAT-1,0:NLONG+1,:)+v_flux_uface(0:NLAT-1,1:NLONG+2,:)&
+2*v_flux_uface(1:NLAT,0:NLONG+1,:)+2*v_flux_uface(1:NLAT,1:NLONG+2,:)&
+v_flux_uface(2:NLAT+1,0:NLONG+1,:)+v_flux_uface(2:NLAT+1,1:NLONG+2,:))
! equation 7.42
C=1.0_dp/12.0_dp*(v_flux_vface(0:NLAT,0:NLONG,:)+v_flux_vface(0:NLAT,1:NLONG+1,:)&
+2*v_flux_vface(1:NLAT+1,0:NLONG,:)+2*v_flux_vface(1:NLAT+1,1:NLONG+1,:)&
+v_flux_vface(2:NLAT+2,0:NLONG,:)+v_flux_vface(2:NLAT+2,1:NLONG+1,:))
! equation 7.43
D=1.0_dp/24.0_dp*(v_flux_vface(0:NLAT,0:NLONG+1,:)+2*v_flux_vface(1:NLAT+1,0:NLONG+1,:)+v_flux_vface(2:NLAT+2,0:NLONG+1,:)&
+v_flux_uface(0:NLAT,0:NLONG+1,:)+v_flux_uface(1:NLAT+1,0:NLONG+1,:)+v_flux_uface(0:NLAT,1:NLONG+2,:)&
+v_flux_uface(1:NLAT+1,1:NLONG+2,:))
! equation 7.44
E=1.0_dp/24.0_dp*(v_flux_vface(0:NLAT,0:NLONG+1,:)+2*v_flux_vface(1:NLAT+1,0:NLONG+1,:)+v_flux_vface(2:NLAT+2,0:NLONG+1,:)&
-v_flux_uface(0:NLAT,0:NLONG+1,:)-v_flux_uface(1:NLAT+1,0:NLONG+1,:)-v_flux_uface(0:NLAT,1:NLONG+2,:)&
-v_flux_uface(1:NLAT+1,1:NLONG+2,:))
! equation 7.45 get velocity at layer bottom
! the top most and bottom most u_bot are not important because w_sigma is zero
do i=1,NLAT
do j=1,NLONG+1
u_bot_old(i,j,1:NVERT-1)=(u_tmp(i,j,1:NVERT-1)*dsigma(1:NVERT-1)+u_tmp(i,j,2:NVERT)*dsigma(2:NVERT))&
/(dsigma(1:NVERT-1)+dsigma(2:NVERT))
u_bot_old(i,j,0)=u_tmp(i,j,1)
u_bot_old(i,j,NVERT)=u_tmp(i,j,NVERT)
enddo
enddo
! equation 7.40
do i=0,NLAT+1
do j=0,NLONG+1
v_F_new(i,j,:)=v_pi_c_new(i,j)*Re*Re*dphi*dlamda_e*cos(v_lat_c(i,j))*v_w_sigma_new(i,j,:)
enddo
enddo
F=1.0_dp/8.0_dp*(v_F_new(2:NLAT+1,0:NLONG,:)+v_F_new(2:NLAT+1,1:NLONG+1,:)&
+2*v_F_new(1:NLAT,0:NLONG,:)+2*v_F_new(1:NLAT,1:NLONG+1,:)&
+v_F_new(0:NLAT-1,0:NLONG,:)+v_F_new(0:NLAT-1,1:NLONG+1,:))
! equation 7.32 time difference term
do k=1,NVERT
u_new(:,:,k)=A_old/A_new*u_old(:,:,k)
enddo
! equation 7.33 horizontal advection B
do k=1,NVERT
u_new(:,:,k)=u_new(:,:,k)+dt/A_new*(B(:,0:NLONG,k)&
*0.5*(v_u_old(1:NLAT,0:NLONG,k)+v_u_old(1:NLAT,1:NLONG+1,k))&
-B(:,1:NLONG+1,k)*0.5*(v_u_old(1:NLAT,1:NLONG+1,k)+v_u_old(1:NLAT,2:NLONG+2,k)))
enddo
! equation 7.33 horizontal advection C
do k=1,NVERT
u_new(:,:,k)=u_new(:,:,k)+dt/A_new*(C(1:NLAT,1:NLONG+1,k)&
*0.5*(v_u_old(0:NLAT-1,1:NLONG+1,k)+v_u_old(1:NLAT,1:NLONG+1,k))&
-C(2:NLAT+1,1:NLONG+1,k)*0.5*(v_u_old(1:NLAT,1:NLONG+1,k)+v_u_old(2:NLAT+1,1:NLONG+1,k)))
enddo
! equation 7.33 horizontal advection D
do k=1,NVERT
u_new(:,:,k)=u_new(:,:,k)+dt/A_new*(D(1:NLAT,0:NLONG,k)&
*0.5*(v_u_old(0:NLAT-1,0:NLONG,k)+v_u_old(1:NLAT,1:NLONG+1,k))&
-D(2:NLAT+1,1:NLONG+1,k)*0.5*(v_u_old(1:NLAT,1:NLONG+1,k)+v_u_old(2:NLAT+1,2:NLONG+2,k)))
enddo
! equation 7.33 horizontal advection E
do k=1,NVERT
u_new(:,:,k)=u_new(:,:,k)+dt/A_new*(E(1:NLAT,1:NLONG+1,k)&
*0.5*(v_u_old(0:NLAT-1,2:NLONG+2,k)+v_u_old(1:NLAT,1:NLONG+1,k))&
-E(2:NLAT+1,0:NLONG,k)*0.5*(v_u_old(1:NLAT,1:NLONG+1,k)+v_u_old(2:NLAT+1,0:NLONG,k)))
enddo
! equation 7.34 vertical transport F
do k=1,NVERT
u_new(:,:,k)=u_new(:,:,k)+dt/A_new/dsigma(k)*(F(:,:,k-1)*u_bot_old(:,:,k-1)&
-F(:,:,k)*u_bot_old(:,:,k))
enddo
! equation 7.35 coriolis and sperical grid conversion
do k=1,NVERT
u_new(:,:,k)=u_new(:,:,k)+0.5*dt/A_new*Re*dlamda_e*dphi&
*(v_pi_c_old(1:NLAT,0:NLONG)*0.5*(v_v_old(1:NLAT,0:NLONG,k)+v_v_old(2:NLAT+1,0:NLONG,k))&
*(v_f_c(1:NLAT,0:NLONG)*Re*cos(v_lat_c(1:NLAT,0:NLONG))&
+0.5*(v_u_old(1:NLAT,0:NLONG,k)+v_u_old(1:NLAT,1:NLONG+1,k))*sin(v_lat_c(1:NLAT,0:NLONG)))&
+v_pi_c_old(1:NLAT,1:NLONG+1)*0.5*(v_v_old(1:NLAT,1:NLONG+1,k)+v_v_old(2:NLAT+1,1:NLONG+1,k))&
*(v_f_c(1:NLAT,1:NLONG+1)*Re*cos(v_lat_c(1:NLAT,1:NLONG+1))&
+0.5*(v_u_old(1:NLAT,1:NLONG+1,k)+v_u_old(1:NLAT,2:NLONG+2,k))*sin(v_lat_c(1:NLAT,1:NLONG+1))))
enddo
! equation 7.36 pressure gradient
do k=1,NVERT
u_new(:,2:NLONG,k)=u_new(:,2:NLONG,k)-dt/A_new(:,2:NLONG)*Re*dphi*((geopot_c(:,2:NLONG,k)&
-geopot_c(:,1:NLONG-1,k))*0.5*(pi_c_tmp(:,1:NLONG-1)+pi_c_tmp(:,2:NLONG))&
+0.5*(-pi_c_tmp(:,1:NLONG-1)+pi_c_tmp(:,2:NLONG))*Cp_d&
*(PVT_c_old(:,1:NLONG-1,k)/dsigma(k)*(sigma_bot(k)&
*(P_bot(:,1:NLONG-1,k)-P_c(:,1:NLONG-1,k))+sigma_bot(k-1)&
*(P_c(:,1:NLONG-1,k)-P_bot(:,1:NLONG-1,k-1)))+PVT_c_old(:,2:NLONG,k)/dsigma(k)&
*(sigma_bot(k)*(P_bot(:,2:NLONG,k)-P_c(:,2:NLONG,k))&
+sigma_bot(k-1)*(P_c(:,2:NLONG,k)-P_bot(:,2:NLONG,k-1)))))
! west boundary
u_new(:,1,k)=u_new(:,1,k)-dt/A_new(:,1)*Re*dphi*((-geopot_c_tminus1(:,1,k)+geopot_c(:,1,k))*pi_c_tmp(:,1)&
+(-pi_c_tminus1(:,1)+pi_c_tmp(:,1))*Cp_d*(PVT_c_old(:,1,k)/dsigma(k)*(sigma_bot(k)*(P_bot(:,1,k)&
-P_c(:,1,k))+sigma_bot(k-1)*(P_c(:,1,k)-P_bot(:,1,k-1)))))
! east boundary
if(periodicBC==1) then
u_new(:,NLONG+1,k)=u_new(:,1,k);
else
u_new(:,NLONG+1,k)=u_new(:,NLONG+1,k)-dt/A_new(:,NLONG+1)*Re*dphi*((geopot_c_tminus1(:,NLONG,k)&
-geopot_c(:,NLONG,k))*pi_c_tmp(:,NLONG)&
+(pi_c_tminus1(:,NLONG)-pi_c_tmp(:,NLONG))*Cp_d*(PVT_c_old(:,NLONG,k)/dsigma(k)*(sigma_bot(k)*(P_bot(:,NLONG,k)&
-P_c(:,NLONG,k))+sigma_bot(k-1)*(P_c(:,NLONG,k)-P_bot(:,NLONG,k-1)))))
endif
enddo
! equation 7.37 eddy vicousity not finished!
if(turbmodel==1) then
do k=1,NVERT
u_new(:,:,k)=u_new(:,:,k)
enddo
endif
! Update v field
! column pressure multiplied by the area eqn 7.53
P_old=1.0_dp/8.0_dp*Re*Re*dphi*dlamda_e*(v_pi_c_old(0:NLAT,2:NLONG+1)*cos(v_lat_c(0:NLAT,2:NLONG+1))&
+v_pi_c_old(1:NLAT+1,2:NLONG+1)*cos(v_lat_c(1:NLAT+1,2:NLONG+1))&
+2*v_pi_c_old(0:NLAT,1:NLONG)*cos(v_lat_c(0:NLAT,1:NLONG))&
+2*v_pi_c_old(1:NLAT+1,1:NLONG)*cos(v_lat_c(1:NLAT+1,1:NLONG))&
+v_pi_c_old(0:NLAT,0:NLONG-1)*cos(v_lat_c(0:NLAT,0:NLONG-1))&
+v_pi_c_old(1:NLAT+1,0:NLONG-1)*cos(v_lat_c(1:NLAT+1,0:NLONG-1)))
P_new=1.0_dp/8.0_dp*Re*Re*dphi*dlamda_e*(v_pi_c_new(0:NLAT,2:NLONG+1)*cos(v_lat_c(0:NLAT,2:NLONG+1))&
+v_pi_c_new(1:NLAT+1,2:NLONG+1)*cos(v_lat_c(1:NLAT+1,2:NLONG+1))&
+2*v_pi_c_new(0:NLAT,1:NLONG)*cos(v_lat_c(0:NLAT,1:NLONG))&
+2*v_pi_c_new(1:NLAT+1,1:NLONG)*cos(v_lat_c(1:NLAT+1,1:NLONG))&
+v_pi_c_new(0:NLAT,0:NLONG-1)*cos(v_lat_c(0:NLAT,0:NLONG-1))&
+v_pi_c_new(1:NLAT+1,0:NLONG-1)*cos(v_lat_c(1:NLAT+1,0:NLONG-1)))
! equation 7.56
R=1.0_dp/12.0_dp*(v_flux_vface(0:NLAT+1,0:NLONG-1,:)+v_flux_vface(1:NLAT+2,0:NLONG-1,:)&
+2*v_flux_vface(0:NLAT+1,1:NLONG,:)+2*v_flux_vface(1:NLAT+2,1:NLONG,:)&
+v_flux_vface(0:NLAT+1,2:NLONG+1,:)+v_flux_vface(1:NLAT+2,2:NLONG+1,:))
! equation 7.55
Q=1.0_dp/12.0_dp*(v_flux_uface(0:NLAT,0:NLONG,:)+v_flux_uface(1:NLAT+1,0:NLONG,:)&
+2*v_flux_uface(0:NLAT,1:NLONG+1,:)+2*v_flux_uface(1:NLAT+1,1:NLONG+1,:)&
+v_flux_uface(0:NLAT,2:NLONG+2,:)+v_flux_uface(1:NLAT+1,2:NLONG+2,:))
! equation 7.57
S=1.0_dp/24.0_dp*(v_flux_uface(0:NLAT+1,0:NLONG,:)+2*v_flux_uface(0:NLAT+1,1:NLONG+1,:)+v_flux_uface(0:NLAT+1,2:NLONG+2,:)&
+v_flux_vface(0:NLAT+1,0:NLONG,:)+v_flux_vface(0:NLAT+1,1:NLONG+1,:)+v_flux_vface(1:NLAT+2,0:NLONG,:)&
+v_flux_vface(1:NLAT+2,1:NLONG+1,:))
! equation 7.58
T=1.0_dp/24.0_dp*(-v_flux_uface(0:NLAT+1,0:NLONG,:)-2*v_flux_uface(0:NLAT+1,1:NLONG+1,:)-v_flux_uface(0:NLAT+1,2:NLONG+2,:)&
+v_flux_vface(0:NLAT+1,0:NLONG,:)+v_flux_vface(0:NLAT+1,1:NLONG+1,:)+v_flux_vface(1:NLAT+2,0:NLONG,:)&
+v_flux_vface(1:NLAT+2,1:NLONG+1,:))
! equation 7.54
O=1.0_dp/8.0_dp*(v_F_new(0:NLAT,2:NLONG+1,:)+v_F_new(1:NLAT+1,2:NLONG+1,:)&
+2*v_F_new(0:NLAT,1:NLONG,:)+2*v_F_new(1:NLAT+1,1:NLONG,:)&
+v_F_new(0:NLAT,0:NLONG-1,:)+v_F_new(1:NLAT+1,0:NLONG-1,:))
! equation 7.45 get velocity at layer bottom
! the top most and bottom most u_bot are not important because w_sigma is zero
do i=1,NLAT+1
do j=1,NLONG
v_bot_old(i,j,1:NVERT-1)=(v_tmp(i,j,1:NVERT-1)*dsigma(1:NVERT-1)+v_tmp(i,j,2:NVERT)*dsigma(2:NVERT))&
/(dsigma(1:NVERT-1)+dsigma(2:NVERT))
v_bot_old(i,j,0)=v_tmp(i,j,1)
v_bot_old(i,j,NVERT)=v_tmp(i,j,NVERT)
enddo
enddo
! equation 7.47 time difference term
do k=1,NVERT
v_new(:,:,k)=P_old/P_new*v_old(:,:,k)
enddo
! equation 7.48 horizontal advection R
do k=1,NVERT
v_new(:,:,k)=v_new(:,:,k)+dt/P_new*(R(0:NLAT,:,k)&
*0.5*(v_v_old(0:NLAT,1:NLONG,k)+v_v_old(1:NLAT+1,1:NLONG,k))&
-R(1:NLAT+1,:,k)*0.5*(v_v_old(1:NLAT+1,1:NLONG,k)+v_v_old(2:NLAT+2,1:NLONG,k)))
enddo
! equation 7.48 horizontal advection Q
do k=1,NVERT
v_new(:,:,k)=v_new(:,:,k)+dt/P_new*(Q(1:NLAT+1,1:NLONG,k)&
*0.5*(v_v_old(1:NLAT+1,0:NLONG-1,k)+v_v_old(1:NLAT+1,1:NLONG,k))&
-Q(1:NLAT+1,2:NLONG+1,k)*0.5*(v_v_old(1:NLAT+1,1:NLONG,k)+v_v_old(1:NLAT+1,2:NLONG+1,k)))
enddo
! equation 7.48 horizontal advection S
do k=1,NVERT
v_new(:,:,k)=v_new(:,:,k)+dt/P_new*(S(0:NLAT,1:NLONG,k)&
*0.5*(v_v_old(0:NLAT,0:NLONG-1,k)+v_v_old(1:NLAT+1,1:NLONG,k))&
-S(1:NLAT+1,2:NLONG+1,k)*0.5*(v_v_old(1:NLAT+1,1:NLONG,k)+v_v_old(2:NLAT+2,2:NLONG+1,k)))
enddo
! equation 7.48 horizontal advection T
do k=1,NVERT
v_new(:,:,k)=v_new(:,:,k)+dt/P_new*(-T(1:NLAT+1,1:NLONG,k)&
*0.5*(v_v_old(2:NLAT+2,0:NLONG-1,k)+v_v_old(1:NLAT+1,1:NLONG,k))&
+T(0:NLAT,2:NLONG+1,k)*0.5*(v_v_old(1:NLAT+1,1:NLONG,k)+v_v_old(0:NLAT,2:NLONG+1,k)))
enddo
! equation 7.49 vertical transport O
do k=1,NVERT
v_new(:,:,k)=v_new(:,:,k)+dt/P_new/dsigma(k)*(O(:,:,k-1)*v_bot_old(:,:,k-1)&
-O(:,:,k)*v_bot_old(:,:,k))
enddo
! equation 7.50 coriolis and sperical grid conversion
do k=1,NVERT
v_new(:,:,k)=v_new(:,:,k)-0.5*dt/P_new*Re*dlamda_e*dphi&
*(v_pi_c_old(0:NLAT,1:NLONG)*0.5*(v_u_old(0:NLAT,1:NLONG,k)+v_u_old(0:NLAT,2:NLONG+1,k))&
*(v_f_c(0:NLAT,1:NLONG)*Re*cos(v_lat_c(0:NLAT,1:NLONG))&
+0.5*(v_u_old(0:NLAT,1:NLONG,k)+v_u_old(0:NLAT,2:NLONG+1,k))*sin(v_lat_c(0:NLAT,1:NLONG)))&
+v_pi_c_old(1:NLAT+1,1:NLONG)*0.5*(v_u_old(1:NLAT+1,1:NLONG,k)+v_u_old(1:NLAT+1,2:NLONG+1,k))&
*(v_f_c(1:NLAT+1,1:NLONG)*Re*cos(v_lat_c(1:NLAT+1,1:NLONG))&
+0.5*(v_u_old(1:NLAT+1,1:NLONG,k)+v_u_old(1:NLAT+1,2:NLONG+1,k))*sin(v_lat_c(1:NLAT+1,1:NLONG))))
enddo
! equation 7.51 pressure gradient
do k=1,NVERT
v_new(2:NLAT,:,k)=v_new(2:NLAT,:,k)-dt/P_new(2:NLAT,:)*Re*dlamda_e*cos(lat_vface(2:NLAT,:))&
*((geopot_c(2:NLAT,:,k)-geopot_c(1:NLAT-1,:,k))*0.5*(pi_c_old(1:NLAT-1,:)+pi_c_old(2:NLAT,:))&
+0.5*(-pi_c_old(1:NLAT-1,:)+pi_c_old(2:NLAT,:))*Cp_d&
*(PVT_c_old(1:NLAT-1,:,k)/dsigma(k)*(sigma_bot(k)&
*(P_bot(1:NLAT-1,:,k)-P_c(1:NLAT-1,:,k))+sigma_bot(k-1)&
*(P_c(1:NLAT-1,:,k)-P_bot(1:NLAT-1,:,k-1)))+PVT_c_old(2:NLAT,:,k)/dsigma(k)&
*(sigma_bot(k)*(P_bot(2:NLAT,:,k)-P_c(2:NLAT,:,k))&
+sigma_bot(k-1)*(P_c(2:NLAT,:,k)-P_bot(2:NLAT,:,k-1)))))