File indexing completed on 2023-07-14 05:10:22 UTC
view on githubraw file Latest commit de57a2ec on 2023-07-13 16:55:13 UTC
a7c06059b4 Gael*0001 #include "ECCO_OPTIONS.h"
0002
9f5240b52a Jean*0003 subroutine cost_gencost_boxmean( myThid )
a7c06059b4 Gael*0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015 implicit none
0016
0017
0018
0019 #include "EEPARAMS.h"
0020 #include "SIZE.h"
0021 #include "PARAMS.h"
0022 #include "GRID.h"
0023 #ifdef ALLOW_CAL
0024 # include "cal.h"
0025 #endif
f09238ab8f Gael*0026 #ifdef ALLOW_ECCO
13d362b8c1 Ou W*0027 # include "ECCO_SIZE.h"
0028 # include "ECCO.h"
a7c06059b4 Gael*0029 #endif
0030
0031
9f5240b52a Jean*0032 integer myThid
a7c06059b4 Gael*0033
0034 #ifdef ALLOW_GENCOST_CONTRIBUTION
0035
0036
f09238ab8f Gael*0037 integer kgen
9f5240b52a Jean*0038 _RL mybar(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
a7c06059b4 Gael*0039
0040 _RL mySumTile(nSx,nSy),myVolTile(nSx,nSy)
bcdcbe969d Gael*0041 _RL mySumGlo,myVolGlo
a7c06059b4 Gael*0042
bcdcbe969d Gael*0043 _RL tmpSumTile(nSx,nSy)
0044 _RL tmpSumGlo
a7c06059b4 Gael*0045
0046 integer bi,bj
9f5240b52a Jean*0047 integer i,j
039a16fdf5 Gael*0048 integer nrecloc,irec,il,ioUnit
de57a2ec4b Mart*0049 character*(MAX_LEN_FNAM) myfname
a7c06059b4 Gael*0050 _RL mydummy
0051 logical doglobalread
0052 logical ladinit
0053 character*(MAX_LEN_MBUF) msgbuf
0054
039a16fdf5 Gael*0055 LOGICAL exst
0056 CHARACTER*(128) tempfile
0057 _RS dummyRS(1)
0058 _RL gencost_mskTemporal
0410c4d6ff Jean*0059 _RL tmpVar(1)
039a16fdf5 Gael*0060
a7c06059b4 Gael*0061
0062 integer ilnblnk
0063 external ilnblnk
0064 LOGICAL MASTER_CPU_THREAD
0065 EXTERNAL MASTER_CPU_THREAD
0066
0067
0068
0069
bcdcbe969d Gael*0070 do kgen=1,NGENCOST
447bdc4b79 Gael*0071 if ( (gencost_flag(kgen).EQ.-3).AND.(using_gencost(kgen)) ) then
a7c06059b4 Gael*0072
0073
0074
0075
0076 doglobalread = .false.
0077 ladinit = .false.
f09238ab8f Gael*0078 mydummy=gencost_dummy(kgen)
0079 il = ilnblnk( gencost_barfile(kgen) )
de57a2ec4b Mart*0080 write(myfname,'(2a,i10.10)')
f09238ab8f Gael*0081 & gencost_barfile(kgen)(1:il),'.',eccoiter
a7c06059b4 Gael*0082
0083
0084 DO bj=myByLo(myThid),myByHi(myThid)
0085 DO bi=myBxLo(myThid),myBxHi(myThid)
0086 mySumTile(bi,bj)=0. _d 0
0087 myVolTile(bi,bj)=0. _d 0
0088 mySumGlo=0. _d 0
0089 myVolGlo=0. _d 0
0090 ENDDO
0091 ENDDO
0092
dcb4d78e77 Gael*0093 nrecloc=gencost_nrec(kgen)
0094
a7c06059b4 Gael*0095
0096
0097
dcb4d78e77 Gael*0098 do irec = 1,nrecloc
a7c06059b4 Gael*0099
0100
101f75e5cd Gael*0101 #ifdef ALLOW_AUTODIFF
bdae1843b8 Gael*0102 call active_read_xy( myfname, mybar, irec,
a7c06059b4 Gael*0103 & doglobalread, ladinit,
9f5240b52a Jean*0104 & eccoiter, myThid,
a7c06059b4 Gael*0105 & mydummy )
101f75e5cd Gael*0106 #else
bdae1843b8 Gael*0107 CALL READ_REC_XY_RL( myfname, mybar,
101f75e5cd Gael*0108 & iRec, 1, myThid )
0109 #endif
a7c06059b4 Gael*0110
0111 DO bj=myByLo(myThid),myByHi(myThid)
0112 DO bi=myBxLo(myThid),myBxHi(myThid)
0113 tmpSumTile(bi,bj)=0. _d 0
0114 tmpSumGlo=0. _d 0
0115 enddo
0116 enddo
0117
47d80787ea Gael*0118 il = ilnblnk(gencost_mask(kgen))
0119 write(tempfile(1:128),'(2A)') gencost_mask(kgen)(1:il),'T'
039a16fdf5 Gael*0120 inquire( file=tempfile(1:il+1), exist=exst )
47d80787ea Gael*0121 if ( (.NOT.exst).OR.(gencost_mask(kgen).EQ.' ') ) then
ecd3edf254 Gael*0122 gencost_mskTemporal=nrecloc
0123 gencost_mskTemporal=1. _d 0 / gencost_mskTemporal
039a16fdf5 Gael*0124 else
0125 ioUnit = 0
0126 call MDS_READVEC_LOC(tempfile,cost_iprec,ioUnit,'RL',
0410c4d6ff Jean*0127 & 1, tmpVar, dummyRS, 0, 0, iRec, myThid )
0128 gencost_mskTemporal = tmpVar(1)
039a16fdf5 Gael*0129 endif
0130
a7c06059b4 Gael*0131
039a16fdf5 Gael*0132 IF ( myProcId .EQ. 0 ) num_gencost(1,1,kgen)=
0133 & num_gencost(1,1,kgen)+gencost_mskTemporal
0134
a7c06059b4 Gael*0135 DO bj=myByLo(myThid),myByHi(myThid)
0136 DO bi=myBxLo(myThid),myBxHi(myThid)
0137 do j = 1,sNy
0138 do i = 1,sNx
0139
039a16fdf5 Gael*0140 objf_gencost(bi,bj,kgen)=objf_gencost(bi,bj,kgen)
0141 & +gencost_mskTemporal*mybar(i,j,bi,bj)
a7c06059b4 Gael*0142
0143
039a16fdf5 Gael*0144 tmpSumTile(bi,bj)=tmpSumTile(bi,bj)
0145 & +gencost_mskTemporal*mybar(i,j,bi,bj)
a7c06059b4 Gael*0146 enddo
0147 enddo
0148 enddo
0149 enddo
0150
0151
0152 CALL GLOBAL_SUM_TILE_RL( tmpSumTile, tmpSumGlo, myThid )
0153
bcdcbe969d Gael*0154 WRITE(msgBuf,'(A,I3,A,1PE21.14)')
447bdc4b79 Gael*0155 & 'boxmean/horflux :',irec,' ',tmpSumGlo
a7c06059b4 Gael*0156 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
0157 & SQUEEZE_RIGHT, myThid )
0158
0159 enddo
0160
0161
0162
0163
f09238ab8f Gael*0164 CALL GLOBAL_SUM_TILE_RL( objf_gencost(1,1,kgen),
b8bb1d0bb0 Gael*0165 & mySumGlo, myThid )
bcdcbe969d Gael*0166
447bdc4b79 Gael*0167 WRITE(msgBuf,'(A,1PE21.14)') 'boxmean/horflux fc :',mySumGlo
a7c06059b4 Gael*0168 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
0169 & SQUEEZE_RIGHT, myThid )
0170
0171
0172
447bdc4b79 Gael*0173 endif
bcdcbe969d Gael*0174 enddo
a7c06059b4 Gael*0175
0176 #endif /* ALLOW_GENCOST_CONTRIBUTION */
0177
9f5240b52a Jean*0178 return
a7c06059b4 Gael*0179 end