Newer
Older
MODULE diadct
!!======================================================================
!! *** MODULE diadct ***
!! Ocean diagnostics: Compute the transport trough a sec.
!!======================================================================
!! History : OPA ! 02/1999 (Y Drillet) original code
!! ! 10/2001 (Y Drillet, R Bourdalle Badie)
!! NEMO 1.0 ! 10/2005 (M Laborie) F90
!! 3.0 ! 04/2007 (G Garric) Ice sections
!! - ! 04/2007 (C Bricaud) test on sec%nb_point, initialisation of ztransp1,ztransp2,...
!! 3.4 ! 09/2011 (C Bricaud)
!!----------------------------------------------------------------------
#if ! defined key_agrif
!! ==>> CAUTION: does not work with agrif
!!----------------------------------------------------------------------
!! dia_dct : Compute the transport through a sec.
!! dia_dct_init : Read namelist.
!! readsec : Read sections description and pathway
!! removepoints : Remove points which are common to 2 procs
!! transport : Compute transport for each sections
!! dia_dct_wri : Write tranports results in ascii files
!! interp : Compute temperature/salinity/density at U-point or V-point
!!
!!----------------------------------------------------------------------
USE oce ! ocean dynamics and tracers
USE dom_oce ! ocean space and time domain
USE phycst ! physical constants
USE in_out_manager ! I/O manager
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
USE daymod ! calendar
USE dianam ! build name of file
USE lib_mpp ! distributed memory computing library
#if defined key_si3
USE ice
#endif
USE domvvl
USE timing ! preformance summary
IMPLICIT NONE
PRIVATE
PUBLIC dia_dct ! routine called by step.F90
PUBLIC dia_dct_init ! routine called by nemogcm.F90
! !!** namelist variables **
INTEGER :: nn_dct ! Frequency of computation
INTEGER :: nn_secdebug ! Number of the section to debug
INTEGER, PARAMETER :: nb_class_max = 10
INTEGER, PARAMETER :: nb_sec_max = 150
INTEGER, PARAMETER :: nb_point_max = 2000
INTEGER, PARAMETER :: nb_type_class = 10
INTEGER, PARAMETER :: nb_3d_vars = 3
INTEGER, PARAMETER :: nb_2d_vars = 2
TYPE POINT_SECTION
INTEGER :: I,J
END TYPE POINT_SECTION
TYPE COORD_SECTION
REAL(wp) :: lon,lat
END TYPE COORD_SECTION
TYPE SECTION
CHARACTER(len=60) :: name ! name of the sec
LOGICAL :: llstrpond ! true if you want the computation of salt and heat transports
LOGICAL :: ll_ice_section ! ice surface and ice volume computation
LOGICAL :: ll_date_line ! = T if the section crosses the date-line
TYPE(COORD_SECTION), DIMENSION(2) :: coordSec ! longitude and latitude of the extremities of the sec
INTEGER :: nb_class ! number of boundaries for density classes
INTEGER, DIMENSION(nb_point_max) :: direction ! vector direction of the point in the section
CHARACTER(len=40),DIMENSION(nb_class_max) :: classname ! characteristics of the class
REAL(wp), DIMENSION(nb_class_max) :: zsigi ! in-situ density classes (99 if you don't want)
REAL(wp), DIMENSION(nb_class_max) :: zsigp ! potential density classes (99 if you don't want)
REAL(wp), DIMENSION(nb_class_max) :: zsal ! salinity classes (99 if you don't want)
REAL(wp), DIMENSION(nb_class_max) :: ztem ! temperature classes(99 if you don't want)
REAL(wp), DIMENSION(nb_class_max) :: zlay ! level classes (99 if you don't want)
REAL(wp), DIMENSION(nb_type_class,nb_class_max) :: transport ! transport output
REAL(wp) :: slopeSection ! slope of the section
INTEGER :: nb_point ! number of points in the section
TYPE(POINT_SECTION),DIMENSION(nb_point_max) :: listPoint ! list of points in the sections
END TYPE SECTION
TYPE(SECTION),DIMENSION(nb_sec_max) :: secs ! Array of sections
REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) :: transports_3d
REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: transports_2d
#if defined key_xios
REAL(wp), ALLOCATABLE, DIMENSION(:) :: heat_transport
REAL(wp), ALLOCATABLE, DIMENSION(:) :: salt_transport
REAL(wp), ALLOCATABLE, DIMENSION(:) :: vol_transport
#endif
# include "single_precision_substitute.h90"
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
# include "domzgr_substitute.h90"
!!----------------------------------------------------------------------
!! NEMO/OCE 4.0 , NEMO Consortium (2018)
!! $Id: diadct.F90 13286 2020-07-09 15:48:29Z smasson $
!! Software governed by the CeCILL license (see ./LICENSE)
!!----------------------------------------------------------------------
CONTAINS
INTEGER FUNCTION diadct_alloc()
!!----------------------------------------------------------------------
!! *** FUNCTION diadct_alloc ***
!!----------------------------------------------------------------------
ALLOCATE( transports_3d(nb_3d_vars,nb_sec_max,nb_point_max,jpk), &
& transports_2d(nb_2d_vars,nb_sec_max,nb_point_max) , STAT=diadct_alloc )
CALL mpp_sum( 'diadct', diadct_alloc )
IF( diadct_alloc /= 0 ) CALL ctl_stop( 'STOP', 'diadct_alloc: failed to allocate arrays' )
END FUNCTION diadct_alloc
SUBROUTINE dia_dct_init
!!---------------------------------------------------------------------
!! *** ROUTINE diadct ***
!!
!! ** Purpose: Read the namelist parameters
!! Open output files
!!
!!---------------------------------------------------------------------
INTEGER :: ios ! Local integer output status for namelist read
!!
NAMELIST/nam_diadct/ln_diadct, nn_dct, nn_dctwri, nn_secdebug
!!---------------------------------------------------------------------
READ ( numnam_ref, nam_diadct, IOSTAT = ios, ERR = 901)
901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'nam_diadct in reference namelist' )
READ ( numnam_cfg, nam_diadct, IOSTAT = ios, ERR = 902 )
902 IF( ios > 0 ) CALL ctl_nam ( ios , 'nam_diadct in configuration namelist' )
IF(lwm) WRITE ( numond, nam_diadct )
IF( lwp ) THEN
WRITE(numout,*) " "
WRITE(numout,*) "diadct_init: compute transports through sections "
WRITE(numout,*) "~~~~~~~~~~~~~~~~~~~~~"
WRITE(numout,*) " Calculate transport thru sections: ln_diadct = ", ln_diadct
WRITE(numout,*) " Frequency of computation: nn_dct = ", nn_dct
WRITE(numout,*) " Frequency of write: nn_dctwri = ", nn_dctwri
IF ( nn_secdebug .GE. 1 .AND. nn_secdebug .LE. nb_sec_max )THEN
WRITE(numout,*)" Debug section number: ", nn_secdebug
ELSE IF ( nn_secdebug == 0 )THEN ; WRITE(numout,*)" No section to debug"
ELSE IF ( nn_secdebug == -1 )THEN ; WRITE(numout,*)" Debug all sections"
ELSE ; WRITE(numout,*)" Wrong value for nn_secdebug : ",nn_secdebug
ENDIF
ENDIF
IF( ln_diadct ) THEN
! control
IF(nn_dct .GE. nn_dctwri .AND. MOD(nn_dct,nn_dctwri) .NE. 0) &
& CALL ctl_stop( 'diadct: nn_dct should be smaller and a multiple of nn_dctwri' )
! allocate dia_dct arrays
IF( diadct_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'diadct_alloc: failed to allocate arrays' )
!Read section_ijglobal.diadct
CALL readsec
#if ! defined key_xios
!open output file
IF( lwm ) THEN
CALL ctl_opn( numdct_vol, 'volume_transport', 'NEW', 'FORMATTED', 'SEQUENTIAL', -1, numout, .FALSE. )
CALL ctl_opn( numdct_heat, 'heat_transport' , 'NEW', 'FORMATTED', 'SEQUENTIAL', -1, numout, .FALSE. )
CALL ctl_opn( numdct_salt, 'salt_transport' , 'NEW', 'FORMATTED', 'SEQUENTIAL', -1, numout, .FALSE. )
ENDIF
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
! Initialise arrays to zero
transports_3d(:,:,:,:)=0.0
transports_2d(:,:,:) =0.0
!
ENDIF
!
END SUBROUTINE dia_dct_init
SUBROUTINE dia_dct( kt, Kmm )
!!---------------------------------------------------------------------
!! *** ROUTINE diadct ***
!!
!! Purpose :: Compute section transports and write it in numdct files
!!
!! Method :: All arrays initialised to zero in dct_init
!! Each nn_dct time step call subroutine 'transports' for
!! each section to sum the transports over each grid cell.
!! Each nn_dctwri time step:
!! Divide the arrays by the number of summations to gain
!! an average value
!! Call dia_dct_sum to sum relevant grid boxes to obtain
!! totals for each class (density, depth, temp or sal)
!! Call dia_dct_wri to write the transports into file
!! Reinitialise all relevant arrays to zero
!!---------------------------------------------------------------------
INTEGER, INTENT(in) :: kt ! ocean time step
INTEGER, INTENT(in) :: Kmm ! time level index
!
INTEGER :: jsec ! loop on sections
INTEGER :: itotal ! nb_sec_max*nb_type_class*nb_class_max
LOGICAL :: lldebug =.FALSE. ! debug a section
INTEGER , DIMENSION(1) :: ish ! work array for mpp_sum
INTEGER , DIMENSION(3) :: ish2 ! "
REAL(wp), ALLOCATABLE, DIMENSION(:) :: zwork ! "
REAL(wp), ALLOCATABLE, DIMENSION(:,:,:):: zsum ! "
!!---------------------------------------------------------------------
!
IF( ln_timing ) CALL timing_start('dia_dct')
IF( lk_mpp )THEN
itotal = nb_sec_max*nb_type_class*nb_class_max
ALLOCATE( zwork(itotal) , zsum(nb_sec_max,nb_type_class,nb_class_max) )
ENDIF
! Initialise arrays
zwork(:) = 0.0
zsum(:,:,:) = 0.0
IF( lwp .AND. kt==nit000+nn_dct-1 ) THEN
WRITE(numout,*) " "
WRITE(numout,*) "diadct: compute transport"
WRITE(numout,*) "~~~~~~~~~~~~~~~~~~~~~~~~~"
WRITE(numout,*) "nb_sec = ",nb_sec
ENDIF
! Compute transport and write only at nn_dctwri
IF( MOD(kt,nn_dct)==0 ) THEN
DO jsec=1,nb_sec
!debug this section computing ?
lldebug=.FALSE.
IF( (jsec==nn_secdebug .OR. nn_secdebug==-1) .AND. kt==nit000+nn_dct-1 ) lldebug=.TRUE.
!Compute transport through section
CALL transport(Kmm,secs(jsec),lldebug,jsec)
ENDDO
IF( MOD(kt,nn_dctwri)==0 )THEN
IF( kt==nit000+nn_dctwri-1 )WRITE(numout,*)" diadct: average transports and write at kt = ",kt
!! divide arrays by nn_dctwri/nn_dct to obtain average
transports_3d(:,:,:,:)=transports_3d(:,:,:,:)/(nn_dctwri/nn_dct)
transports_2d(:,:,:) =transports_2d(:,:,:) /(nn_dctwri/nn_dct)
! Sum over each class
DO jsec=1,nb_sec
CALL dia_dct_sum(Kmm,secs(jsec),jsec)
ENDDO
!Sum on all procs
! IF( lk_mpp )THEN
! removed by AV do not know the significance of this test, following test taken from diaprt
#if ! defined key_mpi_off
IF( .NOT. l_istiled .OR. ntile == nijtile ) THEN
ish(1) = nb_sec_max*nb_type_class*nb_class_max
ish2 = (/nb_sec_max,nb_type_class,nb_class_max/)
DO jsec=1,nb_sec ; zsum(jsec,:,:) = secs(jsec)%transport(:,:) ; ENDDO
zwork(:)= RESHAPE(zsum(:,:,:), ish )
CALL mpp_sum('diadct', zwork, ish(1))
zsum(:,:,:)= RESHAPE(zwork,ish2)
DO jsec=1,nb_sec ; secs(jsec)%transport(:,:) = zsum(jsec,:,:) ; ENDDO
ENDIF
#if defined key_xios
! IF( lwm ) sec_transport(:) = 0
ALLOCATE(heat_transport(nb_sec),salt_transport(nb_sec),vol_transport(nb_sec))
#endif
#if defined key_xios
! xios waits for a send from all procs
CALL dia_dct_wri(kt,jsec,secs(jsec))
#else
!nullify transports values after writing
transports_3d(:,jsec,:,:)=0.
transports_2d(:,jsec,: )=0.
secs(jsec)%transport(:,:)=0.
ENDDO
#if defined key_xios
IF ( .NOT. l_istiled .OR. ntile == nijtile ) THEN
CALL iom_put('mfo' , vol_transport )
CALL iom_put('sfo' , salt_transport )
CALL iom_put('hfo' , heat_transport )
ENDIF
DEALLOCATE(heat_transport, salt_transport, vol_transport)
#endif
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
ENDIF
ENDIF
IF( lk_mpp )THEN
itotal = nb_sec_max*nb_type_class*nb_class_max
DEALLOCATE( zwork , zsum )
ENDIF
IF( ln_timing ) CALL timing_stop('dia_dct')
!
END SUBROUTINE dia_dct
SUBROUTINE readsec
!!---------------------------------------------------------------------
!! *** ROUTINE readsec ***
!!
!! ** Purpose:
!! Read a binary file(section_ijglobal.diadct)
!! generated by the tools "NEMOGCM/TOOLS/SECTIONS_DIADCT"
!!
!!
!!---------------------------------------------------------------------
INTEGER :: iptglo , iptloc ! Global and local number of points for a section
INTEGER :: isec, iiglo, ijglo, iiloc, ijloc,iost,i1 ,i2 ! temporary integer
INTEGER :: jsec, jpt ! dummy loop indices
INTEGER, DIMENSION(2) :: icoord
LOGICAL :: llbon, lldebug ! local logical
CHARACTER(len=160) :: clname ! filename
CHARACTER(len=200) :: cltmp
CHARACTER(len=200) :: clformat !automatic format
TYPE(POINT_SECTION),DIMENSION(nb_point_max) ::coordtemp !contains listpoints coordinates read in the file
INTEGER, DIMENSION(nb_point_max) :: directemp !contains listpoints directions read in the files
!!-------------------------------------------------------------------------------------
!open input file
!---------------
CALL ctl_opn( numdct_in, 'section_ijglobal.diadct', 'OLD', 'UNFORMATTED', 'SEQUENTIAL', -1, numout, .FALSE. )
!---------------
!Read input file
!---------------
DO jsec=1,nb_sec_max !loop on the nb_sec sections
IF ( jsec==nn_secdebug .OR. nn_secdebug==-1 ) &
& WRITE(numout,*)'debuging for section number: ',jsec
!initialization
!---------------
secs(jsec)%name=''
secs(jsec)%llstrpond = .FALSE. ; secs(jsec)%ll_ice_section = .FALSE.
secs(jsec)%ll_date_line = .FALSE. ; secs(jsec)%nb_class = 0
secs(jsec)%zsigi = 99._wp ; secs(jsec)%zsigp = 99._wp
secs(jsec)%zsal = 99._wp ; secs(jsec)%ztem = 99._wp
secs(jsec)%zlay = 99._wp
secs(jsec)%transport = 0._wp ; secs(jsec)%nb_point = 0
!read section's number / name / computing choices / classes / slopeSection / points number
!-----------------------------------------------------------------------------------------
READ(numdct_in,iostat=iost)isec
IF (iost .NE. 0 )EXIT !end of file
WRITE(cltmp,'(a,i4.4,a,i4.4)')'diadct: read sections : Problem of section number: isec= ',isec,' and jsec= ',jsec
IF( jsec .NE. isec ) CALL ctl_stop( cltmp )
IF( jsec==nn_secdebug .OR. nn_secdebug==-1 )WRITE(numout,*)"isec ",isec
READ(numdct_in)secs(jsec)%name
READ(numdct_in)secs(jsec)%llstrpond
READ(numdct_in)secs(jsec)%ll_ice_section
READ(numdct_in)secs(jsec)%ll_date_line
READ(numdct_in)secs(jsec)%coordSec
READ(numdct_in)secs(jsec)%nb_class
READ(numdct_in)secs(jsec)%zsigi
READ(numdct_in)secs(jsec)%zsigp
READ(numdct_in)secs(jsec)%zsal
READ(numdct_in)secs(jsec)%ztem
READ(numdct_in)secs(jsec)%zlay
READ(numdct_in)secs(jsec)%slopeSection
READ(numdct_in)iptglo
!debug
!-----
IF( jsec==nn_secdebug .OR. nn_secdebug==-1 )THEN
WRITE(clformat,'(a,i2,a)') '(A40,', nb_class_max,'(f8.3,1X))'
WRITE(numout,*) " Section name : ",TRIM(secs(jsec)%name)
WRITE(numout,*) " Compute heat and salt transport ? ",secs(jsec)%llstrpond
WRITE(numout,*) " Compute ice transport ? ",secs(jsec)%ll_ice_section
WRITE(numout,*) " Section crosses date-line ? ",secs(jsec)%ll_date_line
WRITE(numout,*) " Slope section : ",secs(jsec)%slopeSection
WRITE(numout,*) " Number of points in the section: ",iptglo
WRITE(numout,*) " Number of classes ",secs(jsec)%nb_class
WRITE(numout,clformat)" Insitu density classes : ",secs(jsec)%zsigi
WRITE(numout,clformat)" Potential density classes : ",secs(jsec)%zsigp
WRITE(numout,clformat)" Salinity classes : ",secs(jsec)%zsal
WRITE(numout,clformat)" Temperature classes : ",secs(jsec)%ztem
WRITE(numout,clformat)" Depth classes : ",secs(jsec)%zlay
ENDIF
IF( iptglo /= 0 )THEN
!read points'coordinates and directions
!--------------------------------------
coordtemp(:) = POINT_SECTION(0,0) !list of points read
directemp(:) = 0 !value of directions of each points
DO jpt=1,iptglo
READ(numdct_in) i1, i2
coordtemp(jpt)%I = i1
coordtemp(jpt)%J = i2
ENDDO
READ(numdct_in) directemp(1:iptglo)
!debug
!-----
IF( jsec==nn_secdebug .OR. nn_secdebug==-1 )THEN
WRITE(numout,*)" List of points in global domain:"
DO jpt=1,iptglo
WRITE(numout,*)' # I J ',jpt,coordtemp(jpt),directemp(jpt)
ENDDO
ENDIF
!Now each proc selects only points that are in its domain:
!--------------------------------------------------------
iptloc = 0 ! initialize number of points selected
DO jpt = 1, iptglo ! loop on listpoint read in the file
!
iiglo=coordtemp(jpt)%I ! global coordinates of the point
ijglo=coordtemp(jpt)%J ! "
IF( iiglo==jpiglo .AND. nimpp==1 ) iiglo = 2 !!gm BUG: Hard coded periodicity !
iiloc=iiglo-nimpp+1 ! local coordinates of the point
ijloc=ijglo-njmpp+1 ! "
!verify if the point is on the local domain:(1,Nie0)*(1,Nje0)
IF( iiloc >= 1 .AND. iiloc <= Nie0 .AND. &
ijloc >= 1 .AND. ijloc <= Nje0 )THEN
iptloc = iptloc + 1 ! count local points
secs(jsec)%listPoint(iptloc) = POINT_SECTION(mi0(iiglo),mj0(ijglo)) ! store local coordinates
secs(jsec)%direction(iptloc) = directemp(jpt) ! store local direction
ENDIF
!
END DO
secs(jsec)%nb_point=iptloc !store number of section's points
!debug
!-----
IF( jsec==nn_secdebug .OR. nn_secdebug==-1 )THEN
WRITE(numout,*)" List of points selected by the proc:"
DO jpt = 1,iptloc
iiglo = secs(jsec)%listPoint(jpt)%I + nimpp - 1
ijglo = secs(jsec)%listPoint(jpt)%J + njmpp - 1
WRITE(numout,*)' # I J : ',iiglo,ijglo
ENDDO
ENDIF
IF(jsec==nn_secdebug .AND. secs(jsec)%nb_point .NE. 0)THEN
DO jpt = 1,iptloc
iiglo = secs(jsec)%listPoint(jpt)%I + nimpp - 1
ijglo = secs(jsec)%listPoint(jpt)%J + njmpp - 1
ENDDO
ENDIF
!remove redundant points between processors
!------------------------------------------
lldebug = .FALSE. ; IF ( jsec==nn_secdebug .OR. nn_secdebug==-1 ) lldebug = .TRUE.
IF( iptloc .NE. 0 )THEN
CALL removepoints(secs(jsec),'I','top_list',lldebug)
CALL removepoints(secs(jsec),'I','bot_list',lldebug)
CALL removepoints(secs(jsec),'J','top_list',lldebug)
CALL removepoints(secs(jsec),'J','bot_list',lldebug)
ENDIF
IF(jsec==nn_secdebug .AND. secs(jsec)%nb_point .NE. 0)THEN
DO jpt = 1,secs(jsec)%nb_point
iiglo = secs(jsec)%listPoint(jpt)%I + nimpp - 1
ijglo = secs(jsec)%listPoint(jpt)%J + njmpp - 1
ENDDO
ENDIF
!debug
!-----
IF( jsec==nn_secdebug .OR. nn_secdebug==-1 )THEN
WRITE(numout,*)" List of points after removepoints:"
iptloc = secs(jsec)%nb_point
DO jpt = 1,iptloc
iiglo = secs(jsec)%listPoint(jpt)%I + nimpp - 1
ijglo = secs(jsec)%listPoint(jpt)%J + njmpp - 1
WRITE(numout,*)' # I J : ',iiglo,ijglo
CALL FLUSH(numout)
ENDDO
ENDIF
ELSE ! iptglo = 0
IF( jsec==nn_secdebug .OR. nn_secdebug==-1 )&
WRITE(numout,*)' No points for this section.'
ENDIF
ENDDO !end of the loop on jsec
nb_sec = jsec-1 !number of section read in the file
!
END SUBROUTINE readsec
SUBROUTINE removepoints(sec,cdind,cdextr,ld_debug)
!!---------------------------------------------------------------------------
!! *** function removepoints
!!
!! ** Purpose :: Remove points which are common to 2 procs
!!
!----------------------------------------------------------------------------
!! * arguments
TYPE(SECTION),INTENT(INOUT) :: sec
CHARACTER(len=1),INTENT(IN) :: cdind ! = 'I'/'J'
CHARACTER(len=8),INTENT(IN) :: cdextr ! = 'top_list'/'bot_list'
LOGICAL,INTENT(IN) :: ld_debug
!! * Local variables
INTEGER :: iextr ,& !extremity of listpoint that we verify
iind ,& !coord of listpoint that we verify
itest ,& !indice value of the side of the domain
!where points could be redundant
isgn ,& ! isgn= 1 : scan listpoint from start to end
! isgn=-1 : scan listpoint from end to start
istart,iend !first and last points selected in listpoint
INTEGER :: jpoint !loop on list points
INTEGER, DIMENSION(nb_point_max) :: idirec !contains temporary sec%direction
INTEGER, DIMENSION(2,nb_point_max) :: icoord !contains temporary sec%listpoint
!----------------------------------------------------------------------------
!
IF( ld_debug )WRITE(numout,*)' -------------------------'
IF( ld_debug )WRITE(numout,*)' removepoints in listpoint'
!iextr=extremity of list_point that we verify
IF ( cdextr=='bot_list' )THEN ; iextr=1 ; isgn=1
ELSE IF ( cdextr=='top_list' )THEN ; iextr=sec%nb_point ; isgn=-1
ELSE ; CALL ctl_stop("removepoints :Wrong value for cdextr")
ENDIF
!which coordinate shall we verify ?
IF ( cdind=='I' )THEN ; itest=Nie0 ; iind=1
ELSE IF ( cdind=='J' )THEN ; itest=Nje0 ; iind=2
ELSE ; CALL ctl_stop("removepoints :Wrong value for cdind")
ENDIF
IF( ld_debug )THEN
WRITE(numout,*)' case: coord/list extr/domain side'
WRITE(numout,*)' ', cdind,' ',cdextr,' ',itest
WRITE(numout,*)' Actual number of points: ',sec%nb_point
ENDIF
icoord(1,1:nb_point_max) = sec%listPoint%I
icoord(2,1:nb_point_max) = sec%listPoint%J
idirec = sec%direction
sec%listPoint = POINT_SECTION(0,0)
sec%direction = 0
jpoint=iextr+isgn
DO WHILE( jpoint .GE. 1 .AND. jpoint .LE. sec%nb_point )
IF( icoord( iind,jpoint-isgn ) == itest .AND. icoord( iind,jpoint ) == itest )THEN ; jpoint=jpoint+isgn
ELSE ; EXIT
ENDIF
ENDDO
IF( cdextr=='bot_list')THEN ; istart=jpoint-1 ; iend=sec%nb_point
ELSE ; istart=1 ; iend=jpoint+1
ENDIF
sec%listPoint(1:1+iend-istart)%I = icoord(1,istart:iend)
sec%listPoint(1:1+iend-istart)%J = icoord(2,istart:iend)
sec%direction(1:1+iend-istart) = idirec(istart:iend)
sec%nb_point = iend-istart+1
IF( ld_debug )THEN
WRITE(numout,*)' Number of points after removepoints :',sec%nb_point
WRITE(numout,*)' sec%direction after removepoints :',sec%direction(1:sec%nb_point)
ENDIF
!
END SUBROUTINE removepoints
SUBROUTINE transport(Kmm,sec,ld_debug,jsec)
!!-------------------------------------------------------------------------------------------
!! *** ROUTINE transport ***
!!
!! Purpose :: Compute the transport for each point in a section
!!
!! Method :: Loop over each segment, and each vertical level and add the transport
!! Be aware :
!! One section is a sum of segments
!! One segment is defined by 2 consecutive points in sec%listPoint
!! All points of sec%listPoint are positioned on the F-point of the cell
!!
!! There are two loops:
!! loop on the segment between 2 nodes
!! loop on the level jk !!
!!
!! Output :: Arrays containing the volume,density,heat,salt transports for each i
!! point in a section, summed over each nn_dct.
!!
!!-------------------------------------------------------------------------------------------
INTEGER ,INTENT(IN) :: Kmm ! time level index
TYPE(SECTION),INTENT(INOUT) :: sec
LOGICAL ,INTENT(IN) :: ld_debug
INTEGER ,INTENT(IN) :: jsec ! numeric identifier of section
!
INTEGER :: jk, jseg, jclass,jl, isgnu, isgnv ! loop on level/segment/classes/ice categories
REAL(wp):: zumid, zvmid, zumid_ice, zvmid_ice ! U/V ocean & ice velocity on a cell segment
REAL(wp):: zTnorm ! transport of velocity through one cell's sides
REAL(wp):: ztn, zsn, zrhoi, zrhop, zsshn, zdep ! temperature/salinity/potential density/ssh/depth at u/v point
TYPE(POINT_SECTION) :: k
!!--------------------------------------------------------
!
IF( ld_debug )WRITE(numout,*)' Compute transport'
!---------------------------!
! COMPUTE TRANSPORT !
!---------------------------!
IF(sec%nb_point .NE. 0)THEN
!----------------------------------------------------------------------------------------------------
!Compute sign for velocities:
!
!convention:
! non horizontal section: direction + is toward left hand of section
! horizontal section: direction + is toward north of section
!
!
! slopeSection < 0 slopeSection > 0 slopeSection=inf slopeSection=0
! ---------------- ----------------- --------------- --------------
!
! isgnv=1 direction +
! ______ _____ ______
! | //| | | direction +
! | isgnu=1 // | |isgnu=1 |isgnu=1 /|\
! |_______ // ______| \\ | ---\ |
! | | isgnv=-1 \\ | | ---/ direction + ____________
! | | __\\| |
! | | direction + | isgnv=1
!
!----------------------------------------------------------------------------------------------------
isgnu = 1
IF( sec%slopeSection .GT. 0 ) THEN ; isgnv = -1
ELSE ; isgnv = 1
ENDIF
IF( sec%slopeSection .GE. 9999. ) isgnv = 1
IF( ld_debug )write(numout,*)"sec%slopeSection isgnu isgnv ",sec%slopeSection,isgnu,isgnv
!--------------------------------------!
! LOOP ON THE SEGMENT BETWEEN 2 NODES !
!--------------------------------------!
DO jseg=1,MAX(sec%nb_point-1,0)
!-------------------------------------------------------------------------------------------
! Select the appropriate coordinate for computing the velocity of the segment
!
! CASE(0) Case (2)
! ------- --------
! listPoint(jseg) listPoint(jseg+1) listPoint(jseg) F(i,j)
! F(i,j)----------V(i+1,j)-------F(i+1,j) |
! |
! |
! |
! Case (3) U(i,j)
! -------- |
! |
! listPoint(jseg+1) F(i,j+1) |
! | |
! | |
! | listPoint(jseg+1) F(i,j-1)
! |
! |
! U(i,j+1)
! | Case(1)
! | ------
! |
! | listPoint(jseg+1) listPoint(jseg)
! | F(i-1,j)-----------V(i,j) -------f(jseg)
! listPoint(jseg) F(i,j)
!
!-------------------------------------------------------------------------------------------
SELECT CASE( sec%direction(jseg) )
CASE(0) ; k = sec%listPoint(jseg)
CASE(1) ; k = POINT_SECTION(sec%listPoint(jseg)%I+1,sec%listPoint(jseg)%J)
CASE(2) ; k = sec%listPoint(jseg)
CASE(3) ; k = POINT_SECTION(sec%listPoint(jseg)%I,sec%listPoint(jseg)%J+1)
END SELECT
!---------------------------|
! LOOP ON THE LEVEL |
!---------------------------|
DO jk = 1, mbkt(k%I,k%J) !Sum of the transport on the vertical
! ! compute temperature, salinity, insitu & potential density, ssh and depth at U/V point
SELECT CASE( sec%direction(jseg) )
CASE(0,1)
ztn = interp(Kmm,k%I,k%J,jk,'V',ts(:,:,:,jp_tem,Kmm) )
zsn = interp(Kmm,k%I,k%J,jk,'V',ts(:,:,:,jp_sal,Kmm) )
zrhop = interp(Kmm,k%I,k%J,jk,'V',rhop)

sparonuz
committed
zrhoi =interp(Kmm,k%I,k%J,jk,'V',CASTDP(rhd*rho0+rho0))
zsshn = 0.5*( ssh(k%I,k%J,Kmm) + ssh(k%I,k%J+1,Kmm) ) * vmask(k%I,k%J,1)
CASE(2,3)
ztn = interp(Kmm,k%I,k%J,jk,'U',ts(:,:,:,jp_tem,Kmm) )
zsn = interp(Kmm,k%I,k%J,jk,'U',ts(:,:,:,jp_sal,Kmm) )
zrhop = interp(Kmm,k%I,k%J,jk,'U',rhop)

sparonuz
committed
zrhoi =interp(Kmm,k%I,k%J,jk,'U',CASTDP(rhd*rho0+rho0))
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
zsshn = 0.5*( ssh(k%I,k%J,Kmm) + ssh(k%I+1,k%J,Kmm) ) * umask(k%I,k%J,1)
END SELECT
!
zdep= gdept(k%I,k%J,jk,Kmm)
SELECT CASE( sec%direction(jseg) ) !compute velocity with the correct direction
CASE(0,1)
zumid=0._wp
zvmid=isgnv*vv(k%I,k%J,jk,Kmm)*vmask(k%I,k%J,jk)
CASE(2,3)
zumid=isgnu*uu(k%I,k%J,jk,Kmm)*umask(k%I,k%J,jk)
zvmid=0._wp
END SELECT
!zTnorm=transport through one cell;
!velocity* cell's length * cell's thickness
zTnorm = zumid*e2u(k%I,k%J) * e3u(k%I,k%J,jk,Kmm) &
& + zvmid*e1v(k%I,k%J) * e3v(k%I,k%J,jk,Kmm)
!!gm THIS is WRONG no transport due to ssh in linear free surface case !!!!!
IF( ln_linssh ) THEN !add transport due to free surface
IF( jk==1 ) THEN
zTnorm = zTnorm + zumid* e2u(k%I,k%J) * zsshn * umask(k%I,k%J,jk) &
& + zvmid* e1v(k%I,k%J) * zsshn * vmask(k%I,k%J,jk)
ENDIF
ENDIF
!!gm end
!COMPUTE TRANSPORT
transports_3d(1,jsec,jseg,jk) = transports_3d(1,jsec,jseg,jk) + zTnorm
IF( sec%llstrpond ) THEN
transports_3d(2,jsec,jseg,jk) = transports_3d(2,jsec,jseg,jk) + zTnorm * ztn * zrhop * rcp
transports_3d(3,jsec,jseg,jk) = transports_3d(3,jsec,jseg,jk) + zTnorm * zsn * zrhop * 0.001
ENDIF
END DO !end of loop on the level
#if defined key_si3
!ICE CASE
!------------
IF( sec%ll_ice_section )THEN
SELECT CASE (sec%direction(jseg))
CASE(0)
zumid_ice = 0
zvmid_ice = isgnv*0.5*(v_ice(k%I,k%J+1)+v_ice(k%I+1,k%J+1))
CASE(1)
zumid_ice = 0
zvmid_ice = isgnv*0.5*(v_ice(k%I,k%J+1)+v_ice(k%I+1,k%J+1))
CASE(2)
zvmid_ice = 0
zumid_ice = isgnu*0.5*(u_ice(k%I+1,k%J)+u_ice(k%I+1,k%J+1))
CASE(3)
zvmid_ice = 0
zumid_ice = isgnu*0.5*(u_ice(k%I+1,k%J)+u_ice(k%I+1,k%J+1))
END SELECT
zTnorm=zumid_ice*e2u(k%I,k%J)+zvmid_ice*e1v(k%I,k%J)
#if defined key_si3
DO jl=1,jpl
transports_2d(1,jsec,jseg) = transports_2d(1,jsec,jseg) + (zTnorm)* &
a_i(sec%listPoint(jseg)%I,sec%listPoint(jseg)%J,jl) * &
( h_i(sec%listPoint(jseg)%I,sec%listPoint(jseg)%J,jl) + &
h_s(sec%listPoint(jseg)%I,sec%listPoint(jseg)%J,jl) )
transports_2d(2,jsec,jseg) = transports_2d(2,jsec,jseg) + (zTnorm)* &
a_i(sec%listPoint(jseg)%I,sec%listPoint(jseg)%J,jl)
END DO
#endif
ENDIF !end of ice case
#endif
END DO !end of loop on the segment
ENDIF !end of sec%nb_point =0 case
!
END SUBROUTINE transport
SUBROUTINE dia_dct_sum(Kmm,sec,jsec)
!!-------------------------------------------------------------
!! Purpose: Average the transport over nn_dctwri time steps
!! and sum over the density/salinity/temperature/depth classes
!!
!! Method: Sum over relevant grid cells to obtain values
!! for each class
!! There are several loops:
!! loop on the segment between 2 nodes
!! loop on the level jk
!! loop on the density/temperature/salinity/level classes
!! test on the density/temperature/salinity/level
!!
!! Note: Transport through a given section is equal to the sum of transports
!! computed on each proc.
!! On each proc,transport is equal to the sum of transport computed through
!! segments linking each point of sec%listPoint with the next one.
!!
!!-------------------------------------------------------------
INTEGER ,INTENT(IN) :: Kmm ! time level index
TYPE(SECTION),INTENT(INOUT) :: sec
INTEGER ,INTENT(IN) :: jsec ! numeric identifier of section
TYPE(POINT_SECTION) :: k
INTEGER :: jk,jseg,jclass ! dummy variables for looping on level/segment/classes
REAL(wp) :: ztn, zsn, zrhoi, zrhop, zsshn, zdep ! temperature/salinity/ssh/potential density /depth at u/v point
!!-------------------------------------------------------------
!! Sum the relevant segments to obtain values for each class
IF(sec%nb_point .NE. 0)THEN
!--------------------------------------!
! LOOP ON THE SEGMENT BETWEEN 2 NODES !
!--------------------------------------!
DO jseg=1,MAX(sec%nb_point-1,0)
!-------------------------------------------------------------------------------------------
! Select the appropriate coordinate for computing the velocity of the segment
!
! CASE(0) Case (2)
! ------- --------
! listPoint(jseg) listPoint(jseg+1) listPoint(jseg) F(i,j)
! F(i,j)----------V(i+1,j)-------F(i+1,j) |
! |
! |
! |
! Case (3) U(i,j)
! -------- |
! |
! listPoint(jseg+1) F(i,j+1) |
! | |
! | |
! | listPoint(jseg+1) F(i,j-1)
! |
! |
! U(i,j+1)
! | Case(1)
! | ------
! |
! | listPoint(jseg+1) listPoint(jseg)
! | F(i-1,j)-----------V(i,j) -------f(jseg)
! listPoint(jseg) F(i,j)
!
!-------------------------------------------------------------------------------------------
SELECT CASE( sec%direction(jseg) )
CASE(0) ; k = sec%listPoint(jseg)
CASE(1) ; k = POINT_SECTION(sec%listPoint(jseg)%I+1,sec%listPoint(jseg)%J)
CASE(2) ; k = sec%listPoint(jseg)
CASE(3) ; k = POINT_SECTION(sec%listPoint(jseg)%I,sec%listPoint(jseg)%J+1)
END SELECT
!---------------------------|
! LOOP ON THE LEVEL |
!---------------------------|
!Sum of the transport on the vertical
DO jk=1,mbkt(k%I,k%J)
! compute temperature, salinity, insitu & potential density, ssh and depth at U/V point
SELECT CASE( sec%direction(jseg) )
CASE(0,1)
ztn = interp(Kmm,k%I,k%J,jk,'V',ts(:,:,:,jp_tem,Kmm) )
zsn = interp(Kmm,k%I,k%J,jk,'V',ts(:,:,:,jp_sal,Kmm) )
zrhop = interp(Kmm,k%I,k%J,jk,'V',rhop)

sparonuz
committed
zrhoi =interp(Kmm,k%I,k%J,jk,'V',CASTDP(rhd*rho0+rho0))
CASE(2,3)
ztn = interp(Kmm,k%I,k%J,jk,'U',ts(:,:,:,jp_tem,Kmm) )
zsn = interp(Kmm,k%I,k%J,jk,'U',ts(:,:,:,jp_sal,Kmm) )
zrhop = interp(Kmm,k%I,k%J,jk,'U',rhop)

sparonuz
committed
zrhoi =interp(Kmm,k%I,k%J,jk,'U',CASTDP(rhd*rho0+rho0))
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
zsshn = 0.5*( ssh(k%I,k%J,Kmm) + ssh(k%I+1,k%J,Kmm) ) * umask(k%I,k%J,1)
END SELECT
zdep= gdept(k%I,k%J,jk,Kmm)
!-------------------------------
! LOOP ON THE DENSITY CLASSES |
!-------------------------------
!The computation is made for each density/temperature/salinity/depth class
DO jclass=1,MAX(1,sec%nb_class-1)
!----------------------------------------------!
!TEST ON THE DENSITY/SALINITY/TEMPERATURE/LEVEL!
!----------------------------------------------!
IF ( ( &
((( zrhop .GE. (sec%zsigp(jclass)+1000. )) .AND. &
( zrhop .LE. (sec%zsigp(jclass+1)+1000. ))) .OR. &
( sec%zsigp(jclass) .EQ. 99.)) .AND. &
((( zrhoi .GE. (sec%zsigi(jclass) + 1000. )) .AND. &
( zrhoi .LE. (sec%zsigi(jclass+1)+1000. ))) .OR. &
( sec%zsigi(jclass) .EQ. 99.)) .AND. &
((( zsn .GT. sec%zsal(jclass)) .AND. &
( zsn .LE. sec%zsal(jclass+1))) .OR. &
( sec%zsal(jclass) .EQ. 99.)) .AND. &
((( ztn .GE. sec%ztem(jclass)) .AND. &
( ztn .LE. sec%ztem(jclass+1))) .OR. &
( sec%ztem(jclass) .EQ.99.)) .AND. &
((( zdep .GE. sec%zlay(jclass)) .AND. &
( zdep .LE. sec%zlay(jclass+1))) .OR. &
( sec%zlay(jclass) .EQ. 99. )) &
)) THEN
!SUM THE TRANSPORTS FOR EACH CLASSES FOR THE POSITIVE AND NEGATIVE DIRECTIONS
!----------------------------------------------------------------------------
IF (transports_3d(1,jsec,jseg,jk) .GE. 0.0) THEN
sec%transport(1,jclass) = sec%transport(1,jclass)+transports_3d(1,jsec,jseg,jk)*1.E-6
ELSE
sec%transport(2,jclass) = sec%transport(2,jclass)+transports_3d(1,jsec,jseg,jk)*1.E-6
ENDIF
IF( sec%llstrpond )THEN
IF ( transports_3d(2,jsec,jseg,jk) .GE. 0.0 ) THEN
sec%transport(3,jclass) = sec%transport(3,jclass)+transports_3d(2,jsec,jseg,jk)
ELSE
sec%transport(4,jclass) = sec%transport(4,jclass)+transports_3d(2,jsec,jseg,jk)
ENDIF
IF ( transports_3d(3,jsec,jseg,jk) .GE. 0.0 ) THEN
sec%transport(5,jclass) = sec%transport(5,jclass)+transports_3d(3,jsec,jseg,jk)
ELSE
sec%transport(6,jclass) = sec%transport(6,jclass)+transports_3d(3,jsec,jseg,jk)
ENDIF
ELSE
sec%transport( 3,jclass) = 0._wp
sec%transport( 4,jclass) = 0._wp
sec%transport( 5,jclass) = 0._wp
sec%transport( 6,jclass) = 0._wp
ENDIF
ENDIF ! end of test if point is in class
END DO ! end of loop on the classes
END DO ! loop over jk
#if defined key_si3
!ICE CASE
IF( sec%ll_ice_section )THEN
IF ( transports_2d(1,jsec,jseg) .GE. 0.0 ) THEN
sec%transport( 7,1) = sec%transport( 7,1)+transports_2d(1,jsec,jseg)*1.E-6
ELSE
sec%transport( 8,1) = sec%transport( 8,1)+transports_2d(1,jsec,jseg)*1.E-6
ENDIF
IF ( transports_2d(3,jsec,jseg) .GE. 0.0 ) THEN
sec%transport( 9,1) = sec%transport( 9,1)+transports_2d(2,jsec,jseg)*1.E-6
ELSE
sec%transport(10,1) = sec%transport(10,1)+transports_2d(2,jsec,jseg)*1.E-6
ENDIF
ENDIF !end of ice case
#endif
END DO !end of loop on the segment
ELSE !if sec%nb_point =0
sec%transport(1:2,:)=0.
IF (sec%llstrpond) sec%transport(3:6,:)=0.
IF (sec%ll_ice_section) sec%transport(7:10,:)=0.
ENDIF !end of sec%nb_point =0 case
END SUBROUTINE dia_dct_sum
SUBROUTINE dia_dct_wri(kt,ksec,sec)
!!-------------------------------------------------------------
!! Write transport output in numdct
!!
!! Purpose: Write transports in ascii files
!!
!! Method:
!! 1. Write volume transports in "volume_transport"
!! Unit: Sv : area * Velocity / 1.e6
!!
!! 2. Write heat transports in "heat_transport"
!! Unit: Peta W : area * Velocity * T * rhop * Cp * 1.e-15