-
Notifications
You must be signed in to change notification settings - Fork 1
/
interface.f90
2077 lines (1616 loc) · 54 KB
/
interface.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
module pop_interface
! !DESCRIPTION:
! This is the interface to the main driver for the Parallel Ocean Program (POP).
!
! !USES:
use POP_KindsMod
use POP_ErrorMod
use POP_FieldMod
use POP_GridHorzMod
use POP_CommMod
use POP_InitMod
use POP_HaloMod
use kinds_mod, only: int_kind, r8
use blocks
use communicate, only: my_task, master_task, MPI_COMM_OCN
use gather_scatter
use grid
use constants, only: grav
use domain, only: distrb_clinic, nprocs_clinic, nprocs_tropic, clinic_distribution_type, &
tropic_distribution_type, ew_boundary_type, ns_boundary_type
use forcing_fields, only: SMF, SMFT, lsmft_avail, STF
use forcing_shf, only: set_shf, SHF_QSW, shf_filename, shf_data_type, &
shf_formulation, shf_interp_freq, shf_interp_type, shf_restore_tau
use forcing_ws, only: set_ws, ws_filename, ws_data_type, ws_data_next, ws_data_update, ws_interp_freq, ws_interp_type, &
ws_interp_next, ws_interp_last, ws_interp_inc
use forcing_sfwf, only: sfwf_filename, sfwf_data_type, sfwf_formulation, sfwf_interp_freq, sfwf_interp_type, fwf_imposed
use forcing_tools, only: never
use initial, only: init_ts_option, init_ts_file, init_ts_file_fmt
use io_types, only: nml_filename
use operators, only: wcalc
use prognostic, only: TRACER, PSURF, UVEL, VVEL, UBTROP, VBTROP, RHO, curtime
use restart, only: restart_freq_opt, restart_freq, restart_outfile
use tavg, only: tavg_freq_opt, tavg_freq, tavg_outfile
use movie, only: movie_freq_opt, movie_freq, movie_outfile
use timers, only: timer_print_all, get_timer, timer_start, timer_stop
use time_management
! use time_management, only: get_time_flag_id, check_time_flag, &
! nsteps_run, stdout, exit_pop, override_time_flag
use step_mod, only: step
use diagnostics, only: check_KE
use output, only: output_driver
use exit_mod, only: sigAbort, sigExit
! use io, only:
implicit none
!interface
! subroutine start_timer ( ) bind (c)
! use iso_c_binding
! end subroutine start_timer
! subroutine stop_timer ( time ) bind (c)
! use iso_c_binding
! real (c_float) :: time
! end subroutine stop_timer
!end interface
!-----------------------------------------------------------------------
!
! module private variables
!
!-----------------------------------------------------------------------
integer (POP_i4) :: &
errorCode ! error flag
integer (int_kind) :: &
timer_step, &! timer number for step
fstop_now ! flag id for stop_now flag
logical :: &
fupdate_coriolis, &! flag for update coriolis force
fupdate_wind_stress ! flag for update wind stress
logical :: initialized = .false.
real (r8), dimension (nx_block, ny_block, max_blocks_clinic) :: &
tau_x, tau_y ! wind stress on t grid in N/m**2
real (r8), dimension (nx_block, ny_block, km, max_blocks_clinic) :: &
WVEL ! 3D vertical velocity field, derived from UVEL and VVEL using wcalc
real (r8), dimension (:,:), allocatable :: &
WORK_G ! temporary 2D array for gathered data
contains
!-----------------------------------------------------------------------
!
! pre-initialize the model run
!
!-----------------------------------------------------------------------
function initialize_code() result(ret)
integer :: ret
errorCode = POP_Success
call POP_Initialize0(errorCode)
! always allocate
!do same for horizontal grid(allocate first, deallocate if not amuse)
allocate(KMT_G(nx_global,ny_global))
if (errorCode /= POP_Success) then
ret=-1
else
ret=0
endif
end function
!-----------------------------------------------------------------------
!
! fully initialize the model run
!
!-----------------------------------------------------------------------
function commit_parameters() result(ret)
integer :: ret
!call POP_CommInitMessageEnvironment
ret=0
errorCode = POP_Success
if(topography_opt.NE.'amuse') then
deallocate(KMT_G) ! will be reallocated later
endif
!if(horiz_grid_opt.EQ.'amuse') then
! allocate (ULAT_G(nx_global, ny_global), &
! ULON_G(nx_global, ny_global))
call POP_Initialize(errorCode)
fstop_now = get_time_flag_id('stop_now')
fupdate_coriolis = .false.
fupdate_wind_stress = .false.
!allocate space for gather of global field
if (my_task == master_task) then
allocate(WORK_G(nx_global,ny_global))
endif
!-----------------------------------------------------------------------
!
! start up the main timer
!
!-----------------------------------------------------------------------
call get_timer(timer_total,'TOTAL',1,distrb_clinic%nprocs)
call timer_start(timer_total)
call get_timer(timer_step,'STEP',1,distrb_clinic%nprocs)
!-----------------------------------------------------------------------
!
! set windstress to that it can be read before the first step is run
!
!-----------------------------------------------------------------------
if (lsmft_avail) then
call set_ws(SMF, SMFT=SMFT)
else
call set_ws(SMF)
endif
!-----------------------------------------------------------------------
!
! set surface heat flux so that it can be read before the first step
! this also sets SHF_QSW, although not passed as an argument
!
!-----------------------------------------------------------------------
call set_shf(STF)
initialized = .true.
! ensure the forcings can be read
ret = prepare_parameters()
! ensure full grid info is present on master_task
call initialize_global_grid
if (errorCode /= POP_Success) then
ret=-1
endif
end function
!-----------------------------------------------------------------------
!
! advance the model in time until the model has run for 'tend' seconds
!
!-----------------------------------------------------------------------
function evolve_model(tend) result(ret)
integer :: ret
real(8), intent(in) :: tend
! stepsize_next is the time in seconds for next timestep
! it is set and updated by POP's time management, this way
! (instead of using dtt) half steps are correctly supported
do while ( tsecond < tend - stepsize_next*0.5 .and. &
errorCode == POP_Success )
call timer_start(timer_step)
call step(errorCode)
call timer_stop(timer_step)
!***
!*** exit if energy is blowing
!***
if (check_KE(100.0_r8)) then
call override_time_flag(fstop_now,value=.true.)
call output_driver(errorCode)
call exit_POP(sigAbort,'ERROR: k.e. > 100 ')
endif
!-----------------------------------------------------------------------
!
! write restart dumps and archiving
!
!-----------------------------------------------------------------------
call output_driver(errorCode)
enddo
if (errorCode == POP_Success) then
ret=0
else
ret=-1
endif
end function
!-----------------------------------------------------------------------
!
! print timing information and clean up various environments if
! they have been used
!
!-----------------------------------------------------------------------
function cleanup_code() result(ret)
integer :: ret
!override time flag fstop_now and produce final output
call override_time_flag(fstop_now,value=.true.)
call output_driver(errorCode)
!produce final output
call timer_stop(timer_total)
call timer_print_all(stats=.true.)
call POP_ErrorPrint(errorCode, printTask=POP_masterTask)
!call exit_POP(sigExit,'Successful completion of POP run')
if (my_task == master_task) then
write(6,*) 'Successful completion of POP run'
endif
ret=0
end function
!-----------------------------------------------------------------------
!
! Returns the elapsed time in the model in seconds
!
!-----------------------------------------------------------------------
function get_model_time(tout) result(ret)
integer :: ret
real(8), intent(out) :: tout
tout = tsecond
ret=0
end function
!-----------------------------------------------------------------------
!
! Returns the default number of seconds in a full time step
!
!-----------------------------------------------------------------------
function get_timestep(dtout) result(ret)
integer :: ret
real(8), intent(out) :: dtout
dtout = dtt
ret=0
end function
!-----------------------------------------------------------------------
!
! Returns number of seconds in the next step, can be less than default for half steps
!
!-----------------------------------------------------------------------
function get_timestep_next(dtout) result(ret)
integer :: ret
real(8), intent(out) :: dtout
dtout = stepsize_next
ret=0
end function
!-----------------------------------------------------------------------
!
! Getter and setter for wind stress per grid point
! Be sure to call prepare_wind_stress() before calling this
!
!-----------------------------------------------------------------------
function get_node_wind_stress(g_i, g_j, tau_x_, tau_y_, n) result(ret)
integer :: ret,n
integer, dimension(n), intent(in) :: g_i, g_j
real*8, dimension(n), intent(out) :: tau_x_, tau_y_
call get_gridded_variable_vector(g_i, g_j, tau_x, tau_x_, n)
call get_gridded_variable_vector(g_i, g_j, tau_y, tau_y_, n)
ret=0
end function
!----------------------------------------------------------------------
!
! This function converts SMF or SMFT to wind stress
! The result is stored in the internal variables tau_x and tau_y
! which can in turn be accessed by get_node_wind_stress()
!
!----------------------------------------------------------------------
function prepare_wind_stress() result(ret)
integer :: ret
tau_x = SMF(:,:,1,:)/momentum_factor * cos(ANGLE) - SMF(:,:,2,:)/momentum_factor * sin(ANGLE)
tau_y = ( SMF(:,:,2,:)/momentum_factor + tau_x * sin(ANGLE) ) / cos (ANGLE)
if (errorCode == POP_Success) then
ret=0
else
ret=-1
endif
end function
function set_node_wind_stress (g_i, g_j, tau_x_, tau_y_, n) result(ret)
integer :: ret, n
real*8, dimension(n), intent(in) :: tau_x_, tau_y_
integer, dimension(n), intent(in) :: g_i, g_j ! global i and j indexes for the gridpoint
call set_gridded_variable_vector(g_i, g_j, tau_x, tau_x_, n)
call set_gridded_variable_vector(g_i, g_j, tau_y, tau_y_, n)
fupdate_wind_stress = .true.
ret=0
end function
function set_wind_stress() result(ret)
integer :: ret
fupdate_wind_stress = .false.
!if we want to overwrite the wind stress, then we have to make sure that SMF is
!not overwritten by some later interpolation routine, therefore, we have to
!enforce the following settings
ws_data_type = 'none'
ws_data_next = never
ws_data_update = never
ws_interp_freq = 'never'
ws_interp_next = never
ws_interp_last = never
ws_interp_inc = c0
! rotate_wind_stress applies tau_x and tau_y to the U-grid
! variable SMF (surface momentum fluxes (wind stress))
! wind stress in N/m^2 is converted to vel flux (cm^2/s^2)
! the method expects tau_x and tau_y to be supplied for the U grid
! set_gridded_variable only sets grid points in the physical domain
! still need to update halo values
call POP_HaloUpdate(tau_x(:,:,:), POP_haloClinic, &
POP_gridHorzLocCenter, &
POP_fieldKindVector, errorCode, &
fillValue = 0.0_r8)
call POP_HaloUpdate(tau_y(:,:,:), POP_haloClinic, &
POP_gridHorzLocCenter, &
POP_fieldKindVector, errorCode, &
fillValue = 0.0_r8)
call rotate_wind_stress(tau_x, tau_y)
if (errorCode == POP_Success) then
ret=0
else
ret=-1
endif
end function
! rotate_wind_stress()
! rotates true zonal/meridional wind stress into local
! coordinates, converts to dyne/cm**2
!
! adapted from forcing_coupled.F90, but modified to:
! 1. not apply the land mask RCALC
! 2. to have tau_x and tau_y be supplied on the U-grid instead of the T-grid
subroutine rotate_wind_stress(WORK1, WORK2)
real (r8), dimension(nx_block,ny_block,max_blocks_clinic), intent(in) :: &
WORK1, WORK2 ! contains tau_x and tau_y from coupler on U-grid
integer :: errorCode
!-----------------------------------------------------------------------
!
! rotate and convert
!
!-----------------------------------------------------------------------
SMF(:,:,1,:) = (WORK1(:,:,:)*cos(ANGLE(:,:,:)) + &
WORK2(:,:,:)*sin(ANGLE(:,:,:)))* &
momentum_factor
SMF(:,:,2,:) = (WORK2(:,:,:)*cos(ANGLE(:,:,:)) - &
WORK1(:,:,:)*sin(ANGLE(:,:,:)))* &
momentum_factor
!-----------------------------------------------------------------------
!
! perform halo updates following the vector rotation
!
!-----------------------------------------------------------------------
call POP_HaloUpdate(SMF(:,:,1,:),POP_haloClinic, &
POP_gridHorzLocCenter, &
POP_fieldKindVector, errorCode, &
fillValue = 0.0_POP_r8)
call POP_HaloUpdate(SMF(:,:,2,:),POP_haloClinic, &
POP_gridHorzLocCenter, &
POP_fieldKindVector, errorCode, &
fillValue = 0.0_POP_r8)
end subroutine rotate_wind_stress
!-----------------------------------------------------------------------
!
! Generic getter and setters for gridded variables
!
!-----------------------------------------------------------------------
subroutine get_gridded_variable(g_i, g_j, grid, value)
integer, intent(in) :: g_i, g_j
real*8, intent(in), dimension(:,:,:) :: grid
real*8, intent(out) :: value
integer :: i,j,iblock
type (block) :: this_block
value = 0.0 !important for MPI_SUM used later
! search rather than iterate over all points
do iblock=1, nblocks_clinic
this_block = get_block(blocks_clinic(iblock),iblock)
! i_glob and j_glob denote the global index of gridpoint within a block
! the ghost cells are also given indices but they may have very different values on cyclic borders etc
! for now, only look at the grid points that are not ghost cells
!
! nx_block and ny_block are the dimensions of the block
! nghost+1 is the first index within the block that is part of the physical domain
! nx_block-nghost is the last index within the block that is part of the physical domain
if (this_block%i_glob(nghost+1) <= g_i .and. this_block%i_glob(nx_block-nghost) >= g_i .and. &
this_block%j_glob(nghost+1) <= g_j .and. this_block%j_glob(ny_block-nghost) >= g_j ) then
!find index in this block from global index (1+ because indexes start at 1 in Fortran)
i = 1+nghost + g_i - this_block%i_glob(1+nghost)
j = 1+nghost + g_j - this_block%j_glob(1+nghost)
!debugging: check if this went ok
!if (this_block%i_glob(i) /= g_i .or. this_block%j_glob(j) /= g_j) then
! write (*,*) 'Error: g_i =', g_i, ', g_j=', g_j, ' and found i=', this_block%i_glob(i), ' j=', this_block%j_glob(j)
!endif
value = grid(i,j,iblock)
endif
enddo
end subroutine get_gridded_variable
subroutine set_gridded_variable(g_i, g_j, grid, value)
integer, intent(in) :: g_i, g_j
real*8, intent(inout), dimension(:,:,:) :: grid
real*8, intent(in) :: value
integer :: i,j,iblock
type (block) :: this_block
! search rather than iterate over all points
do iblock=1, nblocks_clinic
this_block = get_block(blocks_clinic(iblock),iblock)
if (this_block%i_glob(nghost+1) <= g_i .and. this_block%i_glob(nx_block-nghost) >= g_i .and. &
this_block%j_glob(nghost+1) <= g_j .and. this_block%j_glob(ny_block-nghost) >= g_j ) then
i = 1+nghost + g_i - this_block%i_glob(1+nghost)
j = 1+nghost + g_j - this_block%j_glob(1+nghost)
!debugging: check if this went ok
!if (this_block%i_glob(i) /= g_i .or. this_block%j_glob(j) /= g_j) then
! write (*,*) 'Error: g_i =', g_i, ', g_j=', g_j, ' and found i=', this_block%i_glob(i), ' j=', this_block%j_glob(j)
!endif
grid(i,j,iblock) = value
endif
enddo
end subroutine set_gridded_variable
subroutine get_gridded_variable_3D(g_i, g_j, k, grid, value)
integer, intent(in) :: g_i, g_j, k
real*8, intent(in), dimension(:,:,:,:) :: grid
real*8, intent(out) :: value
integer :: i,j,iblock
type (block) :: this_block
do iblock=1, nblocks_clinic
this_block = get_block(blocks_clinic(iblock),iblock)
if (this_block%i_glob(nghost+1) <= g_i .and. this_block%i_glob(nx_block-nghost) >= g_i .and. &
this_block%j_glob(nghost+1) <= g_j .and. this_block%j_glob(ny_block-nghost) >= g_j ) then
i = 1+nghost + g_i - this_block%i_glob(1+nghost)
j = 1+nghost + g_j - this_block%j_glob(1+nghost)
value = grid(i,j,k,iblock)
endif
enddo
end subroutine get_gridded_variable_3D
subroutine set_gridded_variable_3D(g_i, g_j, k, grid, value)
integer, intent(in) :: g_i, g_j, k
real*8, intent(inout), dimension(:,:,:,:) :: grid
real*8, intent(in) :: value
integer :: i,j,iblock
type (block) :: this_block
do iblock=1, nblocks_clinic
this_block = get_block(blocks_clinic(iblock),iblock)
if (this_block%i_glob(nghost+1) <= g_i .and. this_block%i_glob(nx_block-nghost) >= g_i .and. &
this_block%j_glob(nghost+1) <= g_j .and. this_block%j_glob(ny_block-nghost) >= g_j ) then
i = 1+nghost + g_i - this_block%i_glob(1+nghost)
j = 1+nghost + g_j - this_block%j_glob(1+nghost)
grid(i,j,k,iblock) = value
endif
enddo
end subroutine set_gridded_variable_3D
subroutine get_gather(g_i, g_j, grid, value_, n)
integer, intent(in) :: n
integer, dimension(n), intent(in) :: g_i, g_j
real*8, dimension(n), intent(out) :: value_
real*8, dimension(:,:,:), intent(inout) :: grid
integer :: ii
value_ = 0.0 !first zero the array
call gather_global(WORK_G, grid, master_task, distrb_clinic)
if (my_task == master_task) then
do ii=1,n
value_(ii) = WORK_G(g_i(ii),g_j(ii))
enddo
endif
end subroutine get_gather
subroutine get_gather_3D(g_i, g_j, k, grid, value_, n)
integer, intent(in) :: n
integer, dimension(n), intent(in) :: g_i, g_j, k
real*8, dimension(n), intent(out) :: value_
real*8, dimension(:,:,:,:), intent(inout) :: grid
integer :: ii, ki
do ki=1,km
if (ANY(k(:) == ki)) then
call gather_global(WORK_G, grid(:,:,ki,:), master_task, distrb_clinic)
if (my_task == master_task) then
do ii=1,n
if (k(ii) == ki) then
value_(ii) = WORK_G(g_i(ii),g_j(ii))
endif
enddo
endif
endif
enddo
end subroutine get_gather_3D
subroutine get_gridded_variable_vector(g_i, g_j, grid, value, n)
integer, intent(in) :: n
integer, dimension(n), intent(in) :: g_i, g_j
real*8, intent(in), dimension(:,:,:) :: grid
real*8, dimension(n), intent(out) :: value
integer :: ii,ierr
include 'mpif.h'
value = 0.0 !important for MPI_SUM used later
do ii=1,n
call get_gridded_variable(g_i(ii), g_j(ii), grid, value(ii))
enddo
!only the output of the node with MPI rank 0 is significant in AMUSE, reduce to single value per element
if(my_task == master_task) then
call MPI_REDUCE(MPI_IN_PLACE, value, n, MPI_DBL, MPI_SUM, master_task, MPI_COMM_OCN, ierr)
else
call MPI_REDUCE(value, 0, n, MPI_DBL, MPI_SUM, master_task, MPI_COMM_OCN, ierr)
endif
end subroutine get_gridded_variable_vector
subroutine set_gridded_variable_vector(g_i, g_j, grid, value, n)
integer, intent(in) :: n
integer, dimension(n), intent(in) :: g_i, g_j
real*8, intent(inout), dimension(:,:,:) :: grid
real*8, dimension(n), intent(in) :: value
integer :: ii
do ii=1,n
call set_gridded_variable(g_i(ii), g_j(ii), grid, value(ii))
enddo
end subroutine set_gridded_variable_vector
subroutine get_gridded_variable_vector_3D(g_i, g_j, k, grid, value, n)
integer, intent(in) :: n
integer, dimension(n), intent(in) :: g_i, g_j, k
real*8, intent(in), dimension(:,:,:,:) :: grid
real*8, dimension(n), intent(out) :: value
integer :: ii,ierr
include 'mpif.h'
value = 0.0 !important for MPI_SUM used later
do ii=1,n
call get_gridded_variable_3D(g_i(ii), g_j(ii), k(ii), grid, value(ii))
enddo
!only the output of the node with MPI rank 0 is significant in AMUSE, reduce to single value per element
if(my_task == master_task) then
call MPI_REDUCE(MPI_IN_PLACE, value, n, MPI_DBL, MPI_SUM, master_task, MPI_COMM_OCN, ierr)
else
call MPI_REDUCE(value, 0, n, MPI_DBL, MPI_SUM, master_task, MPI_COMM_OCN, ierr)
endif
end subroutine get_gridded_variable_vector_3D
subroutine set_gridded_variable_vector_3D(g_i, g_j, k, grid, value, n)
integer, intent(in) :: n
integer, dimension(n), intent(in) :: g_i, g_j, k
real*8, intent(inout), dimension(:,:,:,:) :: grid
real*8, dimension(n), intent(in) :: value
integer :: ii
do ii=1,n
call set_gridded_variable_3D(g_i(ii), g_j(ii), k(ii), grid, value(ii))
enddo
end subroutine set_gridded_variable_vector_3D
!-----------------------------------------------------------------------
!
! Getter and setter for Coriolis Force
!
! these are currently implemented using the FCOR on the U grid
! set_coriolis_f() can be used to update the FCORT on the T grid based
! on the modified FCOR for the U grid
!
!-----------------------------------------------------------------------
function get_node_coriolis_f(g_i, g_j, corif_, n) result(ret)
integer :: ret, n
integer, dimension(n), intent(in) :: g_i, g_j
real*8, dimension(n), intent(out) :: corif_
if (initialized .eqv. .false.) then
call exit_POP(sigAbort, 'Error: get_node_coriolis_f() called before initialize_code()')
endif
call get_gridded_variable_vector(g_i, g_j, FCOR, corif_, n)
ret=0
end function
function set_node_coriolis_f(g_i, g_j, corif_, n) result(ret)
integer :: ret, n
integer, dimension(n), intent(in) :: g_i, g_j
real*8, dimension(n), intent(in) :: corif_
call set_gridded_variable_vector(g_i, g_j, FCOR, corif_, n)
fupdate_coriolis = .true.
ret=0
end function
! set_coriolis_f updates FCORT (T-grid) based on FCOR (U-grid)
function set_coriolis_f() result(ret)
integer :: ret
integer :: iblock
fupdate_coriolis = .false.
!first update ghost cells for FCOR on U-grid
call POP_HaloUpdate(FCOR, POP_haloClinic, &
POP_gridHorzLocCenter, &
POP_fieldKindVector, errorCode, &
fillValue = 0.0_r8)
!convert to T-grid
do iblock=1, nblocks_clinic
call ugrid_to_tgrid(FCORT, FCOR, iblock)
enddo
!update halo cells on T-grid
call POP_HaloUpdate(FCORT, POP_haloClinic, &
POP_gridHorzLocCenter, &
POP_fieldKindVector, errorCode, &
fillValue = 0.0_r8)
if (errorCode == POP_Success) then
ret=0
else
ret=-1
endif
end function
!-----------------------------------------------------------------------
!
! Recommit forcings will call all functions that had forcings set to
! temporary locations and will ensure these are properly converted to
! internal variables used by POP
!
!-----------------------------------------------------------------------
function recommit_forcings() result(ret)
integer :: ret
ret = 0
if (fupdate_wind_stress) then
ret = set_wind_stress()
endif
if (ret /= 0) return
if (fupdate_coriolis) then
ret = set_coriolis_f()
endif
end function
!-----------------------------------------------------------------------
!
! Getters for node and element surface state
!
!-----------------------------------------------------------------------
function get_node_barotropic_vel(g_i, g_j, uvel_, vvel_, n) result (ret)
integer :: ret
integer, intent(in) :: n
integer, dimension(n), intent(in) :: g_i, g_j
real*8, dimension(n), intent(out) :: uvel_, vvel_
if (n < nx_global*ny_global) then
call get_gridded_variable_vector(g_i, g_j, UBTROP(:,:,curtime,:), uvel_, n)
call get_gridded_variable_vector(g_i, g_j, VBTROP(:,:,curtime,:), vvel_, n)
else
call get_gather(g_i, g_j, UBTROP(:,:,curtime,:), uvel_, n)
call get_gather(g_i, g_j, VBTROP(:,:,curtime,:), vvel_, n)
endif
ret=0
end function
function get_node_surface_state(g_i, g_j, ssh_, uvel_, vvel_, n) result (ret)
integer :: ret
integer, intent(in) :: n
integer, dimension(n), intent(in) :: g_i, g_j
real*8, dimension(n), intent(out) :: ssh_, uvel_, vvel_
if (n < nx_global*ny_global) then
call get_gridded_variable_vector(g_i, g_j, PSURF(:,:,curtime,:), ssh_, n)
ssh_ = ssh_ / grav
call get_gridded_variable_vector(g_i, g_j, UVEL(:,:,1,curtime,:), uvel_, n)
call get_gridded_variable_vector(g_i, g_j, VVEL(:,:,1,curtime,:), vvel_, n)
else
call get_gather(g_i, g_j, PSURF(:,:,curtime,:), ssh_, n)
ssh_ = ssh_ / grav
call get_gather(g_i, g_j, UVEL(:,:,1,curtime,:), uvel_, n)
call get_gather(g_i, g_j, VVEL(:,:,1,curtime,:), vvel_, n)
endif
ret=0
end function
function get_element_surface_state(g_i, g_j, temp_, salt_, n) result (ret)
integer :: ret
integer, intent(in) :: n
integer, dimension(n), intent(in) :: g_i, g_j
real*8, dimension(n), intent(out) :: temp_, salt_
if (my_task == master_task) then
write (*,*) 'get_element_surface_state() called with n=', n
endif
if (n < nx_global*ny_global) then
call get_gridded_variable_vector(g_i, g_j, TRACER(:,:,1,1,curtime,:), temp_, n)
call get_gridded_variable_vector(g_i, g_j, TRACER(:,:,1,2,curtime,:), salt_, n)
else
call get_gather(g_i, g_j, TRACER(:,:,1,1,curtime,:), temp_, n)
call get_gather(g_i, g_j, TRACER(:,:,1,2,curtime,:), salt_, n)
endif
ret=0
end function
function get_element_surface_heat_flux(g_i, g_j, shf_, n) result (ret)
integer :: ret
integer, intent(in) :: n
integer, dimension(n), intent(in) :: g_i, g_j
real*8, dimension(n), intent(out) :: shf_
real*8, dimension(nx_block, ny_block, max_blocks_clinic) :: WORK1
integer :: iblock
do iblock=1, nblocks_clinic
where (KMT(:,:,iblock) > 0)
WORK1(:,:,iblock) = (STF(:,:,1,iblock)+SHF_QSW(:,:,iblock))/hflux_factor ! W/m^2
elsewhere
WORK1(:,:,iblock) = c0
end where
end do
call get_gridded_variable_vector(g_i, g_j, WORK1, shf_, n)
ret=0
end function
!-----------------------------------------------------------------------
!
! Getters for grid info
!
! Nodes are on the U grid
! Elements on the T grid
!
!-----------------------------------------------------------------------
function get_node_position(g_i, g_j, lat_, lon_, n) result(ret)
integer :: ret,n,ii
integer, dimension(n), intent(in) :: g_i, g_j
real*8, dimension(n), intent(out) :: lat_, lon_
! real :: time
integer :: ierr
! if (my_task == master_task) then
! write (*,*) 'get_node_position() called n=', n
! endif
! if (n < nx_global*ny_global) then
! call get_gridded_variable_vector(g_i, g_j, ULAT, lat_, n)
! call get_gridded_variable_vector(g_i, g_j, ULON, lon_, n)
! else
! time = 0.0
! call MPI_Barrier(MPI_COMM_OCN, ierr)
! call start_timer
! call get_gather(g_i, g_j, ULAT, lat_, n)
! call MPI_Barrier(MPI_COMM_OCN, ierr)
! call stop_timer(time)
! if (my_task == master_task) then
! write(*,*) 'get_gather took: ', time, 'ms .'
! endif
! time = 0.0
! call MPI_Barrier(MPI_COMM_OCN, ierr)
! call start_timer
! call get_gridded_variable_vector(g_i, g_j, ULAT, lat_, n)
! call MPI_Barrier(MPI_COMM_OCN, ierr)
! call stop_timer(time)
! if (my_task == master_task) then
! write(*,*) 'get_gridded_variable_vector took: ', time, ' ms.'
! endif
! call get_gather(g_i, g_j, ULON, lon_, n)
if (my_task == master_task) then
do ii=1,n
lon_(ii) = ULON_G(g_i(ii),g_j(ii))
lat_(ii) = ULAT_G(g_i(ii),g_j(ii))
enddo
endif
! endif
ret=0
end function
function get_element_position(g_i, g_j, lat_, lon_, n) result(ret)
integer :: ret,n
integer :: ii
integer, dimension(n), intent(in) :: g_i, g_j
real*8, dimension(n), intent(out) :: lat_, lon_
!if (n < nx_global*ny_global) then
! call get_gridded_variable_vector(g_i, g_j, TLAT, lat_, n)
! call get_gridded_variable_vector(g_i, g_j, TLON, lon_, n)
!else
! call get_gather(g_i, g_j, TLAT, lat_, n)
! call get_gather(g_i, g_j, TLON, lon_, n)
!endif
if (my_task == master_task) then
do ii=1,n
lon_(ii) = TLON_G(g_i(ii),g_j(ii))
lat_(ii) = TLAT_G(g_i(ii),g_j(ii))
enddo
endif
ret=0
end function
function get_node_vposition(g_i, g_j, k, depth_, n) result(ret)
integer :: ret,n
integer, dimension(n), intent(in) :: g_i, g_j, k
real*8, dimension(n), intent(out) :: depth_
ret=0
if (partial_bottom_cells) then
call get_gridded_variable_vector_3D(g_i, g_j, k, DZU, depth_, n)
else
ret = get_zt(k, depth_, n)
endif
end function
function get_element_vposition(g_i, g_j, k, depth_, n) result(ret)
integer :: ret,n
integer, dimension(n), intent(in) :: g_i, g_j, k
real*8, dimension(n), intent(out) :: depth_
ret=0
if (partial_bottom_cells) then
call get_gridded_variable_vector_3D(g_i, g_j, k, DZT, depth_, n)
else
ret = get_zt(k, depth_, n)
endif
end function
function get_number_of_nodes(num) result(ret)
integer :: ret
integer, intent(out) :: num
num = nx_global * ny_global
ret=0
end function
function get_domain_size(x_, y_) result(ret)
integer :: ret
integer, intent(out) :: x_, y_
x_ = nx_global