1 | MODULE mod_oasis_advance |
---|
2 | |
---|
3 | USE mod_oasis_kinds |
---|
4 | USE mod_oasis_data |
---|
5 | USE mod_oasis_parameters |
---|
6 | USE mod_oasis_namcouple, ONLY: namscrmet |
---|
7 | USE mod_oasis_coupler |
---|
8 | USE mod_oasis_part |
---|
9 | USE mod_oasis_timer |
---|
10 | USE mod_oasis_var |
---|
11 | USE mod_oasis_sys |
---|
12 | USE mod_oasis_io |
---|
13 | USE mod_oasis_mpi |
---|
14 | USE mct_mod |
---|
15 | |
---|
16 | IMPLICIT NONE |
---|
17 | |
---|
18 | private |
---|
19 | |
---|
20 | public oasis_advance_init |
---|
21 | public oasis_advance_run |
---|
22 | |
---|
23 | |
---|
24 | contains |
---|
25 | |
---|
26 | !--------------------------------------------------------------------- |
---|
27 | SUBROUTINE oasis_advance_init(kinfo) |
---|
28 | |
---|
29 | ! ---------------------------------------------------------------- |
---|
30 | ! This routine handles initial restart and communication |
---|
31 | ! of data for fields with positive lags |
---|
32 | ! ---------------------------------------------------------------- |
---|
33 | |
---|
34 | IMPLICIT none |
---|
35 | ! ---------------------------------------------------------------- |
---|
36 | INTEGER(kind=ip_i4_p), intent(inout) :: kinfo ! status |
---|
37 | ! ---------------------------------------------------------------- |
---|
38 | integer(kind=ip_i4_p) :: cplid,partid,varid |
---|
39 | INTEGER(kind=ip_i4_p) :: nf,lsize,nflds,icount |
---|
40 | integer(kind=ip_i4_p) :: dt,ltime,lag,getput |
---|
41 | integer(kind=ip_i4_p) :: msec |
---|
42 | real (kind=ip_r8_p), allocatable :: array(:) ! data |
---|
43 | real (kind=ip_r8_p), allocatable :: array2(:) ! data |
---|
44 | real (kind=ip_r8_p), allocatable :: array3(:) ! data |
---|
45 | real (kind=ip_r8_p), allocatable :: array4(:) ! data |
---|
46 | real (kind=ip_r8_p), allocatable :: array5(:) ! data |
---|
47 | logical :: a2on,a3on,a4on,a5on ! data 2-5 logicals |
---|
48 | integer(kind=ip_i4_p) :: mseclag ! model time + lag |
---|
49 | character(len=ic_xl) :: rstfile ! restart filename |
---|
50 | character(len=ic_med) :: vstring ! temporary string |
---|
51 | character(len=*),parameter :: subname = 'oasis_advance_init' |
---|
52 | |
---|
53 | call oasis_debug_enter(subname) |
---|
54 | |
---|
55 | if (mpi_comm_local == MPI_COMM_NULL) then |
---|
56 | write(nulprt,*) subname,' ERROR called on non coupling task' |
---|
57 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
58 | CALL oasis_flush(nulprt) |
---|
59 | call oasis_abort_noarg() |
---|
60 | endif |
---|
61 | |
---|
62 | kinfo = OASIS_OK |
---|
63 | |
---|
64 | call oasis_timer_start ('advance_init') |
---|
65 | |
---|
66 | if (OASIS_debug >= 2) then |
---|
67 | write(nulprt,*) ' subname at time time+lag act: field ' |
---|
68 | write(nulprt,*) ' diags : fldname min max sum ' |
---|
69 | endif |
---|
70 | |
---|
71 | call oasis_debug_note(subname//' loop over cplid') |
---|
72 | cplid=1 |
---|
73 | DO WHILE (cplid <= prism_ncoupler) |
---|
74 | dt = prism_coupler(cplid)%dt |
---|
75 | lag = prism_coupler(cplid)%lag |
---|
76 | ltime = prism_coupler(cplid)%ltime |
---|
77 | getput= prism_coupler(cplid)%getput |
---|
78 | rstfile=TRIM(prism_coupler(cplid)%rstfile) |
---|
79 | partid= prism_coupler(cplid)%partID |
---|
80 | msec = 0 ! reasonable default to start with |
---|
81 | mseclag = msec |
---|
82 | |
---|
83 | IF (OASIS_Debug >= 2) THEN |
---|
84 | WRITE(nulprt,*) '----------------------------------------------------------------' |
---|
85 | WRITE(nulprt,*) subname,' Field cplid :',cplid,TRIM(prism_coupler(cplid)%fldlist) |
---|
86 | WRITE(nulprt,*) subname,' lag prism_ncoupler :',lag,prism_ncoupler |
---|
87 | WRITE(nulprt,*) '----------------------------------------------------------------' |
---|
88 | CALL oasis_flush(nulprt) |
---|
89 | ENDIF |
---|
90 | |
---|
91 | !------------------------------------------------ |
---|
92 | ! check that lag is reasonable |
---|
93 | !------------------------------------------------ |
---|
94 | |
---|
95 | IF (lag > dt .OR. lag <= -dt) THEN |
---|
96 | WRITE(nulprt,*) subname,' ERROR lag out of dt range cplid/dt/lag=',cplid,dt,lag |
---|
97 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
98 | CALL oasis_flush(nulprt) |
---|
99 | CALL oasis_abort_noarg() |
---|
100 | ENDIF |
---|
101 | |
---|
102 | !------------------------------------------------ |
---|
103 | ! read restart and call advance for the current fields |
---|
104 | ! right now, do not know whether any hot map terms are there |
---|
105 | ! assume they are if something is read, otherwise not |
---|
106 | !------------------------------------------------ |
---|
107 | IF ( (getput == OASIS3_PUT .AND. lag > 0) ) THEN |
---|
108 | msec=0-lag |
---|
109 | lsize = mct_aVect_lsize(prism_coupler(cplid)%aVect1) |
---|
110 | nflds = mct_aVect_nRAttr(prism_coupler(cplid)%aVect1) |
---|
111 | |
---|
112 | ALLOCATE(array(lsize)) |
---|
113 | ALLOCATE(array2(lsize)) |
---|
114 | ALLOCATE(array3(lsize)) |
---|
115 | ALLOCATE(array4(lsize)) |
---|
116 | ALLOCATE(array5(lsize)) |
---|
117 | |
---|
118 | DO nf = 1,nflds |
---|
119 | varid = prism_coupler(cplid)%varid(nf) |
---|
120 | CALL oasis_advance_run(OASIS_Out,varid,msec,kinfo,icount=icount,nff=nf,array1din=array,& |
---|
121 | readrest=.TRUE., array2=array2,array3=array3, & |
---|
122 | array4=array4,array5=array5) |
---|
123 | ENDDO |
---|
124 | cplid=cplid+icount |
---|
125 | |
---|
126 | DEALLOCATE(array) |
---|
127 | DEALLOCATE(array2) |
---|
128 | DEALLOCATE(array3) |
---|
129 | DEALLOCATE(array4) |
---|
130 | DEALLOCATE(array5) |
---|
131 | ELSE |
---|
132 | cplid=cplid+1 |
---|
133 | ENDIF |
---|
134 | ENDDO ! cplid |
---|
135 | ! |
---|
136 | DO cplid=1, prism_ncoupler |
---|
137 | dt = prism_coupler(cplid)%dt |
---|
138 | lag = prism_coupler(cplid)%lag |
---|
139 | ltime = prism_coupler(cplid)%ltime |
---|
140 | getput= prism_coupler(cplid)%getput |
---|
141 | rstfile=TRIM(prism_coupler(cplid)%rstfile) |
---|
142 | partid= prism_coupler(cplid)%partID |
---|
143 | msec = 0 ! reasonable default to start with |
---|
144 | mseclag = msec |
---|
145 | |
---|
146 | IF (OASIS_Debug >= 2) THEN |
---|
147 | WRITE(nulprt,*) '----------------------------------------------------------------' |
---|
148 | WRITE(nulprt,*) subname,' Field cplid :',cplid,TRIM(prism_coupler(cplid)%fldlist) |
---|
149 | WRITE(nulprt,*) '----------------------------------------------------------------' |
---|
150 | CALL oasis_flush(nulprt) |
---|
151 | ENDIF |
---|
152 | !------------------------------------------------ |
---|
153 | ! read restart for LOCTRANS fields |
---|
154 | ! do after restart and advance above because prism_advance_run |
---|
155 | ! fills in the avect with the array info |
---|
156 | !------------------------------------------------ |
---|
157 | |
---|
158 | call oasis_debug_note(subname//' check for loctrans restart') |
---|
159 | IF (getput == OASIS3_PUT .AND. prism_coupler(cplid)%trans /= ip_instant) THEN |
---|
160 | if (len_trim(rstfile) < 1) then |
---|
161 | write(nulprt,*) subname,' ERROR restart undefined' |
---|
162 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
163 | CALL oasis_flush(nulprt) |
---|
164 | call oasis_abort_noarg() |
---|
165 | endif |
---|
166 | if (OASIS_debug >= 2) then |
---|
167 | write(nulprt,*) subname,' at ',msec,mseclag,' RTRN: ',& |
---|
168 | trim(prism_coupler(cplid)%fldlist),' ',trim(rstfile) |
---|
169 | endif |
---|
170 | lsize = mct_aVect_lsize(prism_coupler(cplid)%aVect1) |
---|
171 | |
---|
172 | write(vstring,'(a,i2.2,a)') 'loc',prism_coupler(cplid)%trans,'_cnt' |
---|
173 | call oasis_io_read_array(rstfile,iarray=prism_coupler(cplid)%avcnt,& |
---|
174 | ivarname=trim(vstring),abort=.false.) |
---|
175 | |
---|
176 | write(vstring,'(a,i2.2,a)') 'loc',prism_coupler(cplid)%trans,'_' |
---|
177 | call oasis_io_read_avfile(rstfile,prism_coupler(cplid)%avect1,& |
---|
178 | prism_part(partid)%gsmap,abort=.false.,nampre=trim(vstring)) |
---|
179 | |
---|
180 | call mct_aVect_init(prism_coupler(cplid)%aVect2,prism_coupler(cplid)%aVect1,lsize) |
---|
181 | call mct_aVect_zero(prism_coupler(cplid)%aVect2) |
---|
182 | write(vstring,'(a,i2.2,a)') 'av2loc',prism_coupler(cplid)%trans,'_' |
---|
183 | call oasis_io_read_avfile(rstfile,prism_coupler(cplid)%avect2,& |
---|
184 | prism_part(partid)%gsmap,abort=.false.,nampre=trim(vstring),& |
---|
185 | didread=prism_coupler(cplid)%aVon(2)) |
---|
186 | if (.not. prism_coupler(cplid)%aVon(2)) then |
---|
187 | call mct_aVect_clean(prism_coupler(cplid)%avect2) |
---|
188 | endif |
---|
189 | |
---|
190 | call mct_aVect_init(prism_coupler(cplid)%aVect3,prism_coupler(cplid)%aVect1,lsize) |
---|
191 | call mct_aVect_zero(prism_coupler(cplid)%aVect3) |
---|
192 | write(vstring,'(a,i2.2,a)') 'av3loc',prism_coupler(cplid)%trans,'_' |
---|
193 | call oasis_io_read_avfile(rstfile,prism_coupler(cplid)%avect3,& |
---|
194 | prism_part(partid)%gsmap,abort=.false.,nampre=trim(vstring),& |
---|
195 | didread=prism_coupler(cplid)%aVon(3)) |
---|
196 | if (.not. prism_coupler(cplid)%aVon(3)) then |
---|
197 | call mct_aVect_clean(prism_coupler(cplid)%avect3) |
---|
198 | endif |
---|
199 | |
---|
200 | call mct_aVect_init(prism_coupler(cplid)%aVect4,prism_coupler(cplid)%aVect1,lsize) |
---|
201 | call mct_aVect_zero(prism_coupler(cplid)%aVect4) |
---|
202 | write(vstring,'(a,i2.2,a)') 'av4loc',prism_coupler(cplid)%trans,'_' |
---|
203 | call oasis_io_read_avfile(rstfile,prism_coupler(cplid)%avect4,& |
---|
204 | prism_part(partid)%gsmap,abort=.false.,nampre=trim(vstring),& |
---|
205 | didread=prism_coupler(cplid)%aVon(4)) |
---|
206 | if (.not. prism_coupler(cplid)%aVon(4)) then |
---|
207 | call mct_aVect_clean(prism_coupler(cplid)%avect4) |
---|
208 | endif |
---|
209 | |
---|
210 | call mct_aVect_init(prism_coupler(cplid)%aVect5,prism_coupler(cplid)%aVect1,lsize) |
---|
211 | call mct_aVect_zero(prism_coupler(cplid)%aVect5) |
---|
212 | write(vstring,'(a,i2.2,a)') 'av5loc',prism_coupler(cplid)%trans,'_' |
---|
213 | call oasis_io_read_avfile(rstfile,prism_coupler(cplid)%avect5,& |
---|
214 | prism_part(partid)%gsmap,abort=.false.,nampre=trim(vstring),& |
---|
215 | didread=prism_coupler(cplid)%aVon(5)) |
---|
216 | if (.not. prism_coupler(cplid)%aVon(5)) then |
---|
217 | call mct_aVect_clean(prism_coupler(cplid)%avect5) |
---|
218 | endif |
---|
219 | |
---|
220 | if (OASIS_debug >= 20) then |
---|
221 | write(nulprt,*) subname,' DEBUG read loctrans restart',& |
---|
222 | cplid,prism_coupler(cplid)%avcnt |
---|
223 | write(nulprt,*) subname,' DEBUG read loctrans restart',cplid,& |
---|
224 | minval(prism_coupler(cplid)%avect1%rAttr),& |
---|
225 | maxval(prism_coupler(cplid)%avect1%rAttr) |
---|
226 | endif |
---|
227 | endif |
---|
228 | ENDDO ! cplid |
---|
229 | call oasis_timer_stop ('advance_init') |
---|
230 | |
---|
231 | call oasis_debug_exit(subname) |
---|
232 | |
---|
233 | end SUBROUTINE oasis_advance_init |
---|
234 | !--------------------------------------------------------------------- |
---|
235 | SUBROUTINE oasis_advance_run(mop,varid,msec,kinfo,icount,nff,& |
---|
236 | array1din,array1dout,array2dout,readrest,& |
---|
237 | a2on,array2,a3on,array3,a4on,array4,a5on,array5) |
---|
238 | |
---|
239 | IMPLICIT none |
---|
240 | ! ---------------------------------------------------------------- |
---|
241 | integer(kind=ip_i4_p), intent(in) :: mop ! OASIS_Out or OASIS_In |
---|
242 | INTEGER(kind=ip_i4_p), intent(in) :: varid ! prism_var id |
---|
243 | INTEGER(kind=ip_i4_p), intent(in) :: msec ! model time |
---|
244 | INTEGER(kind=ip_i4_p), intent(inout) :: kinfo ! status |
---|
245 | REAL (kind=ip_r8_p), optional :: array1din(:) ! data |
---|
246 | REAL (kind=ip_r8_p), OPTIONAL :: array1dout(:) ! data |
---|
247 | REAL (kind=ip_r8_p), OPTIONAL :: array2dout(:,:) ! data |
---|
248 | logical , optional :: readrest ! special flag to indicate this |
---|
249 | ! is called from the advance_init |
---|
250 | INTEGER(kind=ip_i4_p), OPTIONAL :: nff |
---|
251 | INTEGER(kind=ip_i4_p), optional :: icount |
---|
252 | ! method for restart |
---|
253 | logical , optional :: a2on ! logical for array2 |
---|
254 | REAL (kind=ip_r8_p), optional :: array2(:) ! data |
---|
255 | logical , optional :: a3on ! logical for array3 |
---|
256 | REAL (kind=ip_r8_p), optional :: array3(:) ! data |
---|
257 | logical , optional :: a4on ! logical for array4 |
---|
258 | REAL (kind=ip_r8_p), optional :: array4(:) ! data |
---|
259 | logical , optional :: a5on ! logical for array5 |
---|
260 | REAL (kind=ip_r8_p), optional :: array5(:) ! data |
---|
261 | ! ---------------------------------------------------------------- |
---|
262 | character(len=ic_lvar):: vname |
---|
263 | INTEGER(kind=ip_i4_p) :: cplid,rouid,mapid,partid,namID |
---|
264 | INTEGER(kind=ip_i4_p) :: nfav,nsav,nsa,n,nc,nf |
---|
265 | INTEGER(kind=ip_i4_p) :: lsize,nflds,iicount |
---|
266 | integer(kind=ip_i4_p) :: tag,dt,ltime,lag,getput,maxtime,conserv |
---|
267 | logical :: consbfb |
---|
268 | logical :: sndrcv,output,input,unpack |
---|
269 | logical :: snddiag,rcvdiag |
---|
270 | logical :: arrayon(prism_coupler_avsmax) |
---|
271 | LOGICAL :: a22on,a33on,a44on,a55on |
---|
272 | real(kind=ip_double_p):: sndmult,sndadd,rcvmult,rcvadd |
---|
273 | character(len=ic_xl) :: rstfile ! restart filename |
---|
274 | character(len=ic_xl) :: inpfile ! input filename |
---|
275 | integer(kind=ip_i4_p) :: nx,ny |
---|
276 | integer(kind=ip_i4_p) :: mseclag ! model time + lag |
---|
277 | real(kind=ip_r8_p) :: rcnt ! 1./cnt |
---|
278 | character(len=ic_med) :: tstring ! timer label string |
---|
279 | character(len=ic_med) :: fstring ! output file string |
---|
280 | character(len=ic_med) :: cstring ! temporary string |
---|
281 | character(len=ic_med) :: vstring ! temporary string |
---|
282 | logical :: comm_now ! time to communicate |
---|
283 | logical :: time_now ! coupling time |
---|
284 | logical :: lreadrest ! local readrest |
---|
285 | TYPE(mct_avect) :: avtest ! temporary |
---|
286 | type(mct_avect) :: avtmp ! data read from restart |
---|
287 | type(mct_avect) :: avtmp2 ! data read from restart |
---|
288 | type(mct_avect) :: avtmp3 ! data read from restart |
---|
289 | type(mct_avect) :: avtmp4 ! data read from restart |
---|
290 | type(mct_avect) :: avtmp5 ! data read from restart |
---|
291 | character(len=*),parameter :: subname = 'oasis_advance_run ' |
---|
292 | character(len=*),parameter :: F01 = '(a,i3.3)' |
---|
293 | ! ---------------------------------------------------------------- |
---|
294 | |
---|
295 | call oasis_debug_enter(subname) |
---|
296 | |
---|
297 | if (mpi_comm_local == MPI_COMM_NULL) then |
---|
298 | write(nulprt,*) subname,' ERROR called on non coupling task' |
---|
299 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
300 | CALL oasis_flush(nulprt) |
---|
301 | call oasis_abort_noarg() |
---|
302 | endif |
---|
303 | |
---|
304 | kinfo = OASIS_OK |
---|
305 | vname = prism_var(varid)%name |
---|
306 | |
---|
307 | IF (OASIS_Debug >=15) THEN |
---|
308 | WRITE(nulprt,*) subname,' Variable :',vname,' ncpl :',prism_var(varid)%ncpl |
---|
309 | CALL oasis_flush(nulprt) |
---|
310 | ENDIF |
---|
311 | |
---|
312 | lreadrest = .false. |
---|
313 | if (present(readrest)) then |
---|
314 | lreadrest = readrest |
---|
315 | endif |
---|
316 | if (lreadrest) kinfo = OASIS_fromrest |
---|
317 | |
---|
318 | !------------------------------------------------ |
---|
319 | ! validate mop |
---|
320 | !------------------------------------------------ |
---|
321 | |
---|
322 | if (mop /= OASIS_Out .and. mop /= OASIS_In) then |
---|
323 | write(nulprt,*) subname,' at ',msec,mseclag,' ERROR: ',trim(vname) |
---|
324 | write(nulprt,*) subname,' ERROR mop invalid ',mop |
---|
325 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
326 | CALL oasis_flush(nulprt) |
---|
327 | call oasis_abort_noarg() |
---|
328 | endif |
---|
329 | |
---|
330 | !------------------------------------------------ |
---|
331 | ! for all the couplers associated with this var |
---|
332 | !------------------------------------------------ |
---|
333 | |
---|
334 | call oasis_debug_note(subname//' loop over var ncpl') |
---|
335 | |
---|
336 | IF (PRESENT(icount)) THEN |
---|
337 | icount = 0 |
---|
338 | ENDIF |
---|
339 | DO nc = 1,prism_var(varid)%ncpl |
---|
340 | IF (PRESENT(icount)) THEN |
---|
341 | icount = icount+1 |
---|
342 | ENDIF |
---|
343 | cplid = prism_var(varid)%cpl(nc) |
---|
344 | rouid = prism_coupler(cplid)%routerid |
---|
345 | mapid = prism_coupler(cplid)%mapperid |
---|
346 | tag = prism_coupler(cplid)%tag |
---|
347 | dt = prism_coupler(cplid)%dt |
---|
348 | lag = prism_coupler(cplid)%lag |
---|
349 | ltime = prism_coupler(cplid)%ltime |
---|
350 | getput = prism_coupler(cplid)%getput |
---|
351 | sndrcv = prism_coupler(cplid)%sndrcv |
---|
352 | rstfile =TRIM(prism_coupler(cplid)%rstfile) |
---|
353 | inpfile =TRIM(prism_coupler(cplid)%inpfile) |
---|
354 | maxtime = prism_coupler(cplid)%maxtime |
---|
355 | output = prism_coupler(cplid)%output |
---|
356 | input = prism_coupler(cplid)%input |
---|
357 | partid = prism_coupler(cplid)%partID |
---|
358 | conserv = prism_coupler(cplid)%conserv |
---|
359 | consbfb = .TRUE. |
---|
360 | IF (TRIM(prism_coupler(cplid)%consopt) == "opt") consbfb = .FALSE. |
---|
361 | snddiag = prism_coupler(cplid)%snddiag |
---|
362 | rcvdiag = prism_coupler(cplid)%rcvdiag |
---|
363 | sndadd = prism_coupler(cplid)%sndadd |
---|
364 | sndmult = prism_coupler(cplid)%sndmult |
---|
365 | rcvadd = prism_coupler(cplid)%rcvadd |
---|
366 | rcvmult = prism_coupler(cplid)%rcvmult |
---|
367 | namID = prism_coupler(cplid)%namID |
---|
368 | |
---|
369 | unpack = (sndrcv .OR. input) |
---|
370 | |
---|
371 | CALL oasis_debug_note(subname//' set nx and ny') |
---|
372 | IF (prism_part(partid)%nx >= 1) THEN |
---|
373 | nx = prism_part(partid)%nx |
---|
374 | ny = prism_part(partid)%ny |
---|
375 | ELSE |
---|
376 | nx = prism_part(partid)%gsize |
---|
377 | ny = 1 |
---|
378 | ENDIF |
---|
379 | |
---|
380 | IF (OASIS_debug >= 20) THEN |
---|
381 | WRITE(nulprt,*) subname,' DEBUG nx, ny = ',nx,ny |
---|
382 | CALL oasis_flush(nulprt) |
---|
383 | ENDIF |
---|
384 | |
---|
385 | !------------------------------------------------ |
---|
386 | ! check that lag is reasonable |
---|
387 | !------------------------------------------------ |
---|
388 | |
---|
389 | IF (ABS(lag) > dt) THEN |
---|
390 | WRITE(nulprt,*) subname,' ERROR lag gt dt for cplid',cplid |
---|
391 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
392 | CALL oasis_flush(nulprt) |
---|
393 | CALL oasis_abort_noarg() |
---|
394 | ENDIF |
---|
395 | |
---|
396 | !------------------------------------------------ |
---|
397 | ! read restart and call advance for the current fields |
---|
398 | ! right now, do not know whether any hot map terms are there |
---|
399 | ! assume they are if something is read, otherwise not |
---|
400 | !------------------------------------------------ |
---|
401 | |
---|
402 | call oasis_debug_note(subname//' check for lag restart') |
---|
403 | IF (getput == OASIS3_PUT .AND. lag > 0 .AND. readrest .eqv. .TRUE.) THEN |
---|
404 | ! effective model time of restart : msec |
---|
405 | mseclag = msec + lag |
---|
406 | IF (LEN_TRIM(rstfile) < 1) THEN |
---|
407 | WRITE(nulprt,*) subname,' ERROR restart undefined' |
---|
408 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
409 | CALL oasis_flush(nulprt) |
---|
410 | CALL oasis_abort_noarg() |
---|
411 | ENDIF |
---|
412 | lsize = mct_aVect_lsize(prism_coupler(cplid)%aVect1) |
---|
413 | IF (OASIS_debug >= 2) THEN |
---|
414 | WRITE(nulprt,*) subname,' at ',msec,mseclag,' RRST: ',& |
---|
415 | TRIM(prism_coupler(cplid)%fldlist),' ',TRIM(rstfile) |
---|
416 | ENDIF |
---|
417 | CALL mct_aVect_init(avtmp,rlist=prism_coupler(cplid)%fldlist,lsize=lsize) |
---|
418 | CALL oasis_io_read_avfile(TRIM(rstfile),avtmp,prism_part(partid)%gsmap) |
---|
419 | |
---|
420 | CALL mct_aVect_init(avtmp2,rlist=prism_coupler(cplid)%fldlist,lsize=lsize) |
---|
421 | CALL oasis_io_read_avfile(TRIM(rstfile),avtmp2,prism_part(partid)%gsmap, & |
---|
422 | abort=.FALSE.,nampre='av2_',didread=a22on) |
---|
423 | |
---|
424 | CALL mct_aVect_init(avtmp3,rlist=prism_coupler(cplid)%fldlist,lsize=lsize) |
---|
425 | CALL oasis_io_read_avfile(TRIM(rstfile),avtmp3,prism_part(partid)%gsmap, & |
---|
426 | abort=.FALSE.,nampre='av3_',didread=a33on) |
---|
427 | |
---|
428 | CALL mct_aVect_init(avtmp4,rlist=prism_coupler(cplid)%fldlist,lsize=lsize) |
---|
429 | CALL oasis_io_read_avfile(TRIM(rstfile),avtmp4,prism_part(partid)%gsmap, & |
---|
430 | abort=.FALSE.,nampre='av4_',didread=a44on) |
---|
431 | |
---|
432 | CALL mct_aVect_init(avtmp5,rlist=prism_coupler(cplid)%fldlist,lsize=lsize) |
---|
433 | CALL oasis_io_read_avfile(TRIM(rstfile),avtmp5,prism_part(partid)%gsmap, & |
---|
434 | abort=.false.,nampre='av5_',didread=a55on) |
---|
435 | |
---|
436 | array1din(1:lsize) = avtmp%rAttr(nff,1:lsize) |
---|
437 | IF (a22on) array2(1:lsize) = avtmp2%rAttr(nff,1:lsize) |
---|
438 | IF (a33on) array3(1:lsize) = avtmp3%rAttr(nff,1:lsize) |
---|
439 | IF (a44on) array4(1:lsize) = avtmp4%rAttr(nff,1:lsize) |
---|
440 | IF (a55on) array5(1:lsize) = avtmp5%rAttr(nff,1:lsize) |
---|
441 | |
---|
442 | CALL mct_avect_clean(avtmp) |
---|
443 | IF (a22on) THEN |
---|
444 | CALL mct_avect_clean(avtmp2) |
---|
445 | ENDIF |
---|
446 | IF (a33on) THEN |
---|
447 | CALL mct_avect_clean(avtmp3) |
---|
448 | ENDIF |
---|
449 | IF (a44on) THEN |
---|
450 | CALL mct_avect_clean(avtmp4) |
---|
451 | ENDIF |
---|
452 | IF (a55on) THEN |
---|
453 | CALL mct_avect_clean(avtmp5) |
---|
454 | ENDIF |
---|
455 | ENDIF |
---|
456 | |
---|
457 | |
---|
458 | !------------------------------------------------ |
---|
459 | ! check that model op matches coupler op |
---|
460 | !------------------------------------------------ |
---|
461 | |
---|
462 | if ((mop == OASIS_Out .and. getput == OASIS3_PUT) .or. & |
---|
463 | (mop == OASIS_In .and. getput == OASIS3_GET)) then |
---|
464 | !-continue |
---|
465 | else |
---|
466 | write(nulprt,*) subname,' ERROR model op does not match coupler op',mop,getput |
---|
467 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
468 | CALL oasis_flush(nulprt) |
---|
469 | call oasis_abort_noarg() |
---|
470 | endif |
---|
471 | |
---|
472 | !------------------------------------------------ |
---|
473 | ! compute lag time, only on put side |
---|
474 | ! set time now, is it a coupling period? |
---|
475 | !------------------------------------------------ |
---|
476 | |
---|
477 | call oasis_debug_note(subname//' set mseclag') |
---|
478 | if (getput == OASIS3_PUT) then |
---|
479 | mseclag = msec + lag |
---|
480 | elseif (getput == OASIS3_GET) then |
---|
481 | mseclag = msec |
---|
482 | endif |
---|
483 | |
---|
484 | if (OASIS_debug >= 20) then |
---|
485 | write(nulprt,*) subname,' DEBUG msec,mseclag = ',msec,mseclag |
---|
486 | CALL oasis_flush(nulprt) |
---|
487 | endif |
---|
488 | |
---|
489 | time_now = .false. |
---|
490 | if (mod(mseclag,dt) == 0) time_now = .true. |
---|
491 | |
---|
492 | !------------------------------------------------ |
---|
493 | ! check that model hasn't gone past maxtime |
---|
494 | !------------------------------------------------ |
---|
495 | |
---|
496 | if (msec >= maxtime) then |
---|
497 | write(nulprt,*) subname,' at ',msec,mseclag,' ERROR: ',trim(vname) |
---|
498 | write(nulprt,*) subname,' ERROR model time beyond namcouple maxtime',& |
---|
499 | msec,maxtime |
---|
500 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
501 | CALL oasis_flush(nulprt) |
---|
502 | call oasis_abort_noarg() |
---|
503 | endif |
---|
504 | |
---|
505 | !------------------------------------------------ |
---|
506 | ! check that model isn't going backwards |
---|
507 | ! msec >= 0 does the check only in run mode, not in initialization |
---|
508 | !------------------------------------------------ |
---|
509 | |
---|
510 | if (lcouplertime /= ispval .and. msec >= 0 .and. msec < lcouplertime) then |
---|
511 | write(nulprt,*) subname,' at ',msec,mseclag,' ERROR: ',trim(vname) |
---|
512 | write(nulprt,*) subname,' ERROR model seems to be running backwards',& |
---|
513 | msec,lcouplertime |
---|
514 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
515 | CALL oasis_flush(nulprt) |
---|
516 | call oasis_abort_noarg() |
---|
517 | endif |
---|
518 | |
---|
519 | !------------------------------------------------ |
---|
520 | ! check that varible didn't miss a coupling period |
---|
521 | ! also check that prior sequences weren't missed at this |
---|
522 | ! step if get operation. only done for sndrcv operations. |
---|
523 | ! attempts to trap deadlocks before they happen |
---|
524 | !------------------------------------------------ |
---|
525 | |
---|
526 | do n = 1,prism_ncoupler |
---|
527 | if ((prism_coupler(n)%ltime /= ispval) .and. & |
---|
528 | (sndrcv .and. prism_coupler(n)%sndrcv) .and. & |
---|
529 | (msec > prism_coupler(n)%ltime + prism_coupler(n)%dt)) then |
---|
530 | write(nulprt,*) subname,' ERROR coupling skipped at earlier time, & |
---|
531 | & potential deadlock ' |
---|
532 | write(nulprt,*) subname,' my coupler = ',cplid,' variable = ',& |
---|
533 | trim(prism_coupler(cplid)%fldlist) |
---|
534 | write(nulprt,*) subname,' current time = ',msec,' mseclag = ',mseclag |
---|
535 | write(nulprt,*) subname,' skipped coupler = ',n,' variable = ',& |
---|
536 | trim(prism_coupler(n)%fldlist) |
---|
537 | write(nulprt,*) subname,' skipped coupler last time and dt = ',& |
---|
538 | prism_coupler(n)%ltime,prism_coupler(n)%dt |
---|
539 | WRITE(nulprt,*) subname,' ERROR model timestep does not match coupling timestep' |
---|
540 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
541 | call oasis_flush(nulprt) |
---|
542 | call oasis_abort_noarg() |
---|
543 | endif |
---|
544 | if ((prism_coupler(n)%ltime /= ispval) .and. & |
---|
545 | (sndrcv .and. prism_coupler(n)%sndrcv .and. getput == OASIS3_GET) .and. & |
---|
546 | (prism_coupler(cplid)%seq > prism_coupler(n)%seq) .and. & |
---|
547 | (msec >= prism_coupler(n)%ltime + prism_coupler(n)%dt)) then |
---|
548 | write(nulprt,*) subname,' ERROR coupling sequence out of order, & |
---|
549 | & potential deadlock ' |
---|
550 | write(nulprt,*) subname,' my coupler = ',cplid,' variable = ',& |
---|
551 | trim(prism_coupler(cplid)%fldlist) |
---|
552 | write(nulprt,*) subname,' sequence number = ',prism_coupler(cplid)%seq |
---|
553 | write(nulprt,*) subname,' current time = ',msec,' mseclag = ',mseclag |
---|
554 | write(nulprt,*) subname,' skipped coupler = ',n,' variable = ',& |
---|
555 | trim(prism_coupler(n)%fldlist) |
---|
556 | write(nulprt,*) subname,' skipped coupler last time and dt = ',& |
---|
557 | prism_coupler(n)%ltime,prism_coupler(n)%dt |
---|
558 | write(nulprt,*) subname,' skipped sequence number = ',prism_coupler(n)%seq |
---|
559 | WRITE(nulprt,*) subname,' ERROR model sequence does not match coupling sequence' |
---|
560 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
561 | CALL oasis_flush(nulprt) |
---|
562 | call oasis_abort_noarg() |
---|
563 | endif |
---|
564 | enddo |
---|
565 | |
---|
566 | !------------------------------------------------ |
---|
567 | ! compute field index and check sizes |
---|
568 | !------------------------------------------------ |
---|
569 | |
---|
570 | call oasis_debug_note(subname//' compute field index and sizes') |
---|
571 | nfav = mct_avect_indexra(prism_coupler(cplid)%avect1,trim(vname)) |
---|
572 | nsav = mct_avect_lsize(prism_coupler(cplid)%avect1) |
---|
573 | if (lag > 0 .and. readrest .eqv. .true. ) nsa=size(array1din) |
---|
574 | if (present(array1din )) nsa = size(array1din ) |
---|
575 | if (present(array1dout)) nsa = size(array1dout) |
---|
576 | if (present(array2dout)) nsa = size(array2dout) |
---|
577 | |
---|
578 | if (OASIS_debug >= 20) then |
---|
579 | write(nulprt,*) subname,' DEBUG nfav,nsav,nsa = ',nfav,nsav,nsa |
---|
580 | endif |
---|
581 | |
---|
582 | if (nsav /= nsa) then |
---|
583 | write(nulprt,*) subname,' at ',msec,mseclag,' ERROR: ',trim(vname) |
---|
584 | write(nulprt,*) subname,' ERROR sizes ',nsav,nsa |
---|
585 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
586 | CALL oasis_flush(nulprt) |
---|
587 | call oasis_abort_noarg() |
---|
588 | endif |
---|
589 | |
---|
590 | !------------------------------------------------ |
---|
591 | ! check for higher order coupling fields |
---|
592 | ! and get everything ready |
---|
593 | ! arrayon is what's passed this time |
---|
594 | ! optional args only on put side |
---|
595 | !------------------------------------------------ |
---|
596 | |
---|
597 | IF (readrest .neqv. .TRUE.) THEN |
---|
598 | arrayon = .false. |
---|
599 | arrayon(1) = .true. |
---|
600 | if (present(a2on)) arrayon(2) = a2on |
---|
601 | if (present(a3on)) arrayon(3) = a3on |
---|
602 | if (present(a4on)) arrayon(4) = a4on |
---|
603 | if (present(a5on)) arrayon(5) = a5on |
---|
604 | |
---|
605 | if ((getput == OASIS3_GET) .or. & |
---|
606 | (getput == OASIS3_PUT .and. trim(prism_coupler(cplid)%maploc) == "dst" )) then |
---|
607 | if (arrayon(2) .or. arrayon(3) .or. & |
---|
608 | arrayon(4) .or. arrayon(5)) then |
---|
609 | write(nulprt,*) subname,' at ',msec,mseclag,' ERROR: ',trim(vname) |
---|
610 | write(nulprt,*) subname,' higher order mapping not allowed on get side' |
---|
611 | write(nulprt,*) subname,' consider changing map location from dst to src' |
---|
612 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
613 | CALL oasis_flush(nulprt) |
---|
614 | call oasis_abort_noarg() |
---|
615 | endif |
---|
616 | endif |
---|
617 | |
---|
618 | if ((arrayon(2) .and. .not.present(array2)) .or. & |
---|
619 | (arrayon(3) .and. .not.present(array3)) .or. & |
---|
620 | (arrayon(4) .and. .not.present(array4)) .or. & |
---|
621 | (arrayon(5) .and. .not.present(array5))) then |
---|
622 | write(nulprt,*) subname,' at ',msec,mseclag,' ERROR: ',trim(vname) |
---|
623 | write(nulprt,*) subname,' arrayon true but array not sent' |
---|
624 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
625 | CALL oasis_flush(nulprt) |
---|
626 | call oasis_abort_noarg() |
---|
627 | ! With the current way of using oasis_advance_run, the above test is useless but we keep the test |
---|
628 | ! as someone might be later adding an interface call that would violate the consistency |
---|
629 | endif |
---|
630 | |
---|
631 | ! IF ( (TRIM(namscrmet(namID)) == 'BICUBIC' ) .AND. (.NOT. (arrayon(2))) & |
---|
632 | ! .AND. (.NOT. (arrayon(3))) .AND. (.NOT. (arrayon(4))) ) THEN |
---|
633 | ! WRITE(nulprt,*) subname,' ERROR : BICUBIC GRADIENT is asked ' & |
---|
634 | ! WRITE(nulprt,*) subname,' 2nd, 3rd and 4th arguments must be provided' |
---|
635 | ! CALL oasis_flush(nulprt,*) |
---|
636 | ! CALL oasis_abort_noarg() |
---|
637 | ! ENDIF |
---|
638 | |
---|
639 | ! initialize aVect2-5 here if not already allocated |
---|
640 | |
---|
641 | if (arrayon(2) .and. .not. prism_coupler(cplid)%aVon(2)) then |
---|
642 | call mct_aVect_init(prism_coupler(cplid)%aVect2,prism_coupler(cplid)%aVect1,nsav) |
---|
643 | call mct_aVect_zero(prism_coupler(cplid)%aVect2) |
---|
644 | prism_coupler(cplid)%aVon(2) = .true. |
---|
645 | if (OASIS_debug >= 2) then |
---|
646 | write(nulprt,*) subname,' at ',msec,mseclag,' ALLO: ',& |
---|
647 | trim(vname),' ','aVect2' |
---|
648 | endif |
---|
649 | endif |
---|
650 | |
---|
651 | if (arrayon(3) .and. .not. prism_coupler(cplid)%aVon(3)) then |
---|
652 | call mct_aVect_init(prism_coupler(cplid)%aVect3,prism_coupler(cplid)%aVect1,nsav) |
---|
653 | call mct_aVect_zero(prism_coupler(cplid)%aVect3) |
---|
654 | prism_coupler(cplid)%aVon(3) = .true. |
---|
655 | if (OASIS_debug >= 2) then |
---|
656 | write(nulprt,*) subname,' at ',msec,mseclag,' ALLO: ',& |
---|
657 | trim(vname),' ','aVect3' |
---|
658 | endif |
---|
659 | endif |
---|
660 | |
---|
661 | if (arrayon(4) .and. .not. prism_coupler(cplid)%aVon(4)) then |
---|
662 | call mct_aVect_init(prism_coupler(cplid)%aVect4,prism_coupler(cplid)%aVect1,nsav) |
---|
663 | call mct_aVect_zero(prism_coupler(cplid)%aVect4) |
---|
664 | prism_coupler(cplid)%aVon(4) = .true. |
---|
665 | if (OASIS_debug >= 2) then |
---|
666 | write(nulprt,*) subname,' at ',msec,mseclag,' ALLO: ',& |
---|
667 | trim(vname),' ','aVect4' |
---|
668 | endif |
---|
669 | endif |
---|
670 | |
---|
671 | if (arrayon(5) .and. .not. prism_coupler(cplid)%aVon(5)) then |
---|
672 | call mct_aVect_init(prism_coupler(cplid)%aVect5,prism_coupler(cplid)%aVect1,nsav) |
---|
673 | call mct_aVect_zero(prism_coupler(cplid)%aVect5) |
---|
674 | prism_coupler(cplid)%aVon(5) = .true. |
---|
675 | if (OASIS_debug >= 2) then |
---|
676 | write(nulprt,*) subname,' at ',msec,mseclag,' ALLO: ',& |
---|
677 | trim(vname),' ','aVect5' |
---|
678 | endif |
---|
679 | endif |
---|
680 | endif |
---|
681 | |
---|
682 | !------------------------------------------------ |
---|
683 | ! update avect1-5 on put side, apply appropriate transform |
---|
684 | ! if its coupling time, set status of this var to ready |
---|
685 | ! on restart, treat as instant value |
---|
686 | !------------------------------------------------ |
---|
687 | |
---|
688 | if (getput == OASIS3_PUT) then |
---|
689 | |
---|
690 | call oasis_debug_note(subname//' loctrans operation') |
---|
691 | write(tstring,F01) 'pcpy_',cplid |
---|
692 | call oasis_timer_start(tstring) |
---|
693 | |
---|
694 | cstring = 'none' |
---|
695 | if (lreadrest .or. prism_coupler(cplid)%trans == ip_instant) then |
---|
696 | if (time_now) then |
---|
697 | cstring = 'instant' |
---|
698 | do n = 1,nsav |
---|
699 | prism_coupler(cplid)%avect1%rAttr(nfav,n) = array1din(n) |
---|
700 | if (prism_coupler(cplid)%aVon(2)) then |
---|
701 | if (present(array2)) then |
---|
702 | prism_coupler(cplid)%avect2%rAttr(nfav,n) = array2(n) |
---|
703 | else |
---|
704 | prism_coupler(cplid)%avect2%rAttr(nfav,n) = 0.0 |
---|
705 | endif |
---|
706 | endif |
---|
707 | if (prism_coupler(cplid)%aVon(3)) then |
---|
708 | if (present(array3)) then |
---|
709 | prism_coupler(cplid)%avect3%rAttr(nfav,n) = array3(n) |
---|
710 | else |
---|
711 | prism_coupler(cplid)%avect3%rAttr(nfav,n) = 0.0 |
---|
712 | endif |
---|
713 | endif |
---|
714 | if (prism_coupler(cplid)%aVon(4)) then |
---|
715 | if (present(array4)) then |
---|
716 | prism_coupler(cplid)%avect4%rAttr(nfav,n) = array4(n) |
---|
717 | else |
---|
718 | prism_coupler(cplid)%avect4%rAttr(nfav,n) = 0.0 |
---|
719 | endif |
---|
720 | endif |
---|
721 | if (prism_coupler(cplid)%aVon(5)) then |
---|
722 | if (present(array5)) then |
---|
723 | prism_coupler(cplid)%avect5%rAttr(nfav,n) = array5(n) |
---|
724 | else |
---|
725 | prism_coupler(cplid)%avect5%rAttr(nfav,n) = 0.0 |
---|
726 | endif |
---|
727 | endif |
---|
728 | enddo |
---|
729 | prism_coupler(cplid)%avcnt(nfav) = 1 |
---|
730 | endif |
---|
731 | |
---|
732 | elseif (prism_coupler(cplid)%trans == ip_average) then |
---|
733 | cstring = 'average' |
---|
734 | if (kinfo == OASIS_OK) kinfo = OASIS_LocTrans |
---|
735 | do n = 1,nsav |
---|
736 | prism_coupler(cplid)%avect1%rAttr(nfav,n) = & |
---|
737 | prism_coupler(cplid)%avect1%rAttr(nfav,n) + array1din(n) |
---|
738 | if (prism_coupler(cplid)%aVon(2)) then |
---|
739 | if (present(array2)) then |
---|
740 | prism_coupler(cplid)%avect2%rAttr(nfav,n) = & |
---|
741 | prism_coupler(cplid)%avect2%rAttr(nfav,n) + array2(n) |
---|
742 | endif |
---|
743 | endif |
---|
744 | if (prism_coupler(cplid)%aVon(3)) then |
---|
745 | if (present(array3)) then |
---|
746 | prism_coupler(cplid)%avect3%rAttr(nfav,n) = & |
---|
747 | prism_coupler(cplid)%avect3%rAttr(nfav,n) + array3(n) |
---|
748 | endif |
---|
749 | endif |
---|
750 | if (prism_coupler(cplid)%aVon(4)) then |
---|
751 | if (present(array4)) then |
---|
752 | prism_coupler(cplid)%avect4%rAttr(nfav,n) = & |
---|
753 | prism_coupler(cplid)%avect4%rAttr(nfav,n) + array4(n) |
---|
754 | endif |
---|
755 | endif |
---|
756 | if (prism_coupler(cplid)%aVon(5)) then |
---|
757 | if (present(array5)) then |
---|
758 | prism_coupler(cplid)%avect5%rAttr(nfav,n) = & |
---|
759 | prism_coupler(cplid)%avect5%rAttr(nfav,n) + array5(n) |
---|
760 | endif |
---|
761 | endif |
---|
762 | enddo |
---|
763 | prism_coupler(cplid)%avcnt(nfav) = prism_coupler(cplid)%avcnt(nfav) + 1 |
---|
764 | |
---|
765 | elseif (prism_coupler(cplid)%trans == ip_accumul) then |
---|
766 | cstring = 'accumul' |
---|
767 | if (kinfo == OASIS_OK) kinfo = OASIS_LocTrans |
---|
768 | do n = 1,nsav |
---|
769 | prism_coupler(cplid)%avect1%rAttr(nfav,n) = & |
---|
770 | prism_coupler(cplid)%avect1%rAttr(nfav,n) + array1din(n) |
---|
771 | if (prism_coupler(cplid)%aVon(2)) then |
---|
772 | if (present(array2)) then |
---|
773 | prism_coupler(cplid)%avect2%rAttr(nfav,n) = & |
---|
774 | prism_coupler(cplid)%avect2%rAttr(nfav,n) + array2(n) |
---|
775 | endif |
---|
776 | endif |
---|
777 | if (prism_coupler(cplid)%aVon(3)) then |
---|
778 | if (present(array3)) then |
---|
779 | prism_coupler(cplid)%avect3%rAttr(nfav,n) = & |
---|
780 | prism_coupler(cplid)%avect3%rAttr(nfav,n) + array3(n) |
---|
781 | endif |
---|
782 | endif |
---|
783 | if (prism_coupler(cplid)%aVon(4)) then |
---|
784 | if (present(array4)) then |
---|
785 | prism_coupler(cplid)%avect4%rAttr(nfav,n) = & |
---|
786 | prism_coupler(cplid)%avect4%rAttr(nfav,n) + array4(n) |
---|
787 | endif |
---|
788 | endif |
---|
789 | if (prism_coupler(cplid)%aVon(5)) then |
---|
790 | if (present(array5)) then |
---|
791 | prism_coupler(cplid)%avect5%rAttr(nfav,n) = & |
---|
792 | prism_coupler(cplid)%avect5%rAttr(nfav,n) + array5(n) |
---|
793 | endif |
---|
794 | endif |
---|
795 | enddo |
---|
796 | prism_coupler(cplid)%avcnt(nfav) = 1 |
---|
797 | |
---|
798 | elseif (prism_coupler(cplid)%trans == ip_max) then |
---|
799 | cstring = 'max' |
---|
800 | if (kinfo == OASIS_OK) kinfo = OASIS_LocTrans |
---|
801 | if (prism_coupler(cplid)%aVon(2) .or. prism_coupler(cplid)%aVon(3) .or. & |
---|
802 | prism_coupler(cplid)%aVon(4) .or. prism_coupler(cplid)%aVon(5)) then |
---|
803 | write(nulprt,*) subname,' at ',msec,mseclag,' ERROR: ',trim(vname) |
---|
804 | write(nulprt,*) subname,' higher order mapping with MAX trans not supported' |
---|
805 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
806 | CALL oasis_flush(nulprt) |
---|
807 | call oasis_abort_noarg() |
---|
808 | endif |
---|
809 | do n = 1,nsav |
---|
810 | if (prism_coupler(cplid)%avcnt(nfav) == 0) then |
---|
811 | prism_coupler(cplid)%avect1%rAttr(nfav,n) = array1din(n) |
---|
812 | else |
---|
813 | prism_coupler(cplid)%avect1%rAttr(nfav,n) = & |
---|
814 | max(prism_coupler(cplid)%avect1%rAttr(nfav,n),array1din(n)) |
---|
815 | endif |
---|
816 | enddo |
---|
817 | prism_coupler(cplid)%avcnt(nfav) = 1 |
---|
818 | |
---|
819 | elseif (prism_coupler(cplid)%trans == ip_min) then |
---|
820 | cstring = 'min' |
---|
821 | if (kinfo == OASIS_OK) kinfo = OASIS_LocTrans |
---|
822 | if (prism_coupler(cplid)%aVon(2) .or. prism_coupler(cplid)%aVon(3) .or. & |
---|
823 | prism_coupler(cplid)%aVon(4) .or. prism_coupler(cplid)%aVon(5)) then |
---|
824 | write(nulprt,*) subname,' at ',msec,mseclag,' ERROR: ',trim(vname) |
---|
825 | write(nulprt,*) subname,' higher order mapping with MIN trans not supported' |
---|
826 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
827 | CALL oasis_flush(nulprt) |
---|
828 | call oasis_abort_noarg() |
---|
829 | endif |
---|
830 | do n = 1,nsav |
---|
831 | if (prism_coupler(cplid)%avcnt(nfav) == 0) then |
---|
832 | prism_coupler(cplid)%avect1%rAttr(nfav,n) = array1din(n) |
---|
833 | else |
---|
834 | prism_coupler(cplid)%avect1%rAttr(nfav,n) = & |
---|
835 | min(prism_coupler(cplid)%avect1%rAttr(nfav,n),array1din(n)) |
---|
836 | endif |
---|
837 | enddo |
---|
838 | prism_coupler(cplid)%avcnt(nfav) = 1 |
---|
839 | |
---|
840 | else |
---|
841 | write(nulprt,*) subname,' ERROR: trans not known ',prism_coupler(cplid)%trans |
---|
842 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
843 | CALL oasis_flush(nulprt) |
---|
844 | call oasis_abort_noarg() |
---|
845 | endif |
---|
846 | call oasis_timer_stop(tstring) |
---|
847 | |
---|
848 | if (OASIS_debug >= 2 .and. trim(cstring) /= 'none') then |
---|
849 | write(nulprt,*) subname,' at ',msec,mseclag,' PACK: ',& |
---|
850 | trim(vname),' ',trim(cstring) |
---|
851 | endif |
---|
852 | |
---|
853 | if (OASIS_debug >= 20) then |
---|
854 | write(nulprt,*) subname,' DEBUG loctrans update ',cplid,' ',& |
---|
855 | trim(cstring),prism_coupler(cplid)%avcnt(nfav) |
---|
856 | endif |
---|
857 | |
---|
858 | if (time_now) then |
---|
859 | prism_coupler(cplid)%status(nfav) = OASIS_COMM_READY |
---|
860 | endif |
---|
861 | endif |
---|
862 | |
---|
863 | !------------------------------------------------ |
---|
864 | ! decide if it's time to communicate based on |
---|
865 | ! time. also, on the put side, status of all vars |
---|
866 | ! must be ready which means all vars have called put. |
---|
867 | ! on get side, all ready means all vars have unpacked |
---|
868 | ! from last get. |
---|
869 | !------------------------------------------------ |
---|
870 | |
---|
871 | call oasis_debug_note(subname//' comm_now compute') |
---|
872 | comm_now = .false. |
---|
873 | if (time_now) then |
---|
874 | comm_now = .true. |
---|
875 | do nf = 1,prism_coupler(cplid)%nflds |
---|
876 | if (prism_coupler(cplid)%status(nf) /= OASIS_COMM_READY) then |
---|
877 | comm_now = .false. |
---|
878 | if (OASIS_debug >= 15) then |
---|
879 | write(nulprt,*) subname,' at ',msec,mseclag,' STAT: ',nf,' NOT READY' |
---|
880 | WRITE(nulprt,*) subname,' Status : ',prism_coupler(cplid)%status(nf) |
---|
881 | endif |
---|
882 | kinfo=OASIS_Waitforallingroup |
---|
883 | else |
---|
884 | if (OASIS_debug >= 15) then |
---|
885 | write(nulprt,*) subname,' at ',msec,mseclag,' STAT: ',nf,' READY' |
---|
886 | endif |
---|
887 | endif |
---|
888 | enddo |
---|
889 | endif |
---|
890 | |
---|
891 | if (comm_now) then |
---|
892 | |
---|
893 | call oasis_debug_note(subname//' comm_now') |
---|
894 | |
---|
895 | !------------------------------------------------ |
---|
896 | ! this is the time critical bit, we need to make sure the |
---|
897 | ! model is truly advancing in time when comms are called. |
---|
898 | ! must ignore the initial call, ltime = 0 |
---|
899 | !------------------------------------------------ |
---|
900 | |
---|
901 | if (prism_coupler(cplid)%ltime /= ispval .and. msec <= prism_coupler(cplid)%ltime) then |
---|
902 | write(nulprt,*) subname,' ERROR: model did not advance in time correctly',& |
---|
903 | msec,prism_coupler(cplid)%ltime |
---|
904 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
905 | CALL oasis_flush(nulprt) |
---|
906 | call oasis_abort_noarg() |
---|
907 | endif |
---|
908 | |
---|
909 | !------------------------------------------------ |
---|
910 | ! average as needed (not cache friendly yet) |
---|
911 | !------------------------------------------------ |
---|
912 | |
---|
913 | IF (getput == OASIS3_PUT) THEN |
---|
914 | call oasis_debug_note(subname//' loctrans calc') |
---|
915 | write(tstring,F01) 'pavg_',cplid |
---|
916 | call oasis_timer_start(tstring) |
---|
917 | do nf = 1,prism_coupler(cplid)%nflds |
---|
918 | if (prism_coupler(cplid)%avcnt(nf) > 1) then |
---|
919 | rcnt = 1.0/prism_coupler(cplid)%avcnt(nf) |
---|
920 | do n = 1,nsav |
---|
921 | prism_coupler(cplid)%avect1%rAttr(nf,n) = & |
---|
922 | prism_coupler(cplid)%avect1%rAttr(nf,n) * rcnt |
---|
923 | if (prism_coupler(cplid)%aVon(2)) then |
---|
924 | prism_coupler(cplid)%avect2%rAttr(nf,n) = & |
---|
925 | prism_coupler(cplid)%avect2%rAttr(nf,n) * rcnt |
---|
926 | endif |
---|
927 | if (prism_coupler(cplid)%aVon(3)) then |
---|
928 | prism_coupler(cplid)%avect3%rAttr(nf,n) = & |
---|
929 | prism_coupler(cplid)%avect3%rAttr(nf,n) * rcnt |
---|
930 | endif |
---|
931 | if (prism_coupler(cplid)%aVon(4)) then |
---|
932 | prism_coupler(cplid)%avect4%rAttr(nf,n) = & |
---|
933 | prism_coupler(cplid)%avect4%rAttr(nf,n) * rcnt |
---|
934 | endif |
---|
935 | if (prism_coupler(cplid)%aVon(5)) then |
---|
936 | prism_coupler(cplid)%avect5%rAttr(nf,n) = & |
---|
937 | prism_coupler(cplid)%avect5%rAttr(nf,n) * rcnt |
---|
938 | endif |
---|
939 | enddo |
---|
940 | endif |
---|
941 | if (OASIS_debug >= 20) then |
---|
942 | write(nulprt,*) subname,' DEBUG loctrans calc0 = ',cplid,nf,& |
---|
943 | prism_coupler(cplid)%avcnt(nf) |
---|
944 | write(nulprt,*) subname,' DEBUG loctrans calc1 = ',cplid,nf,& |
---|
945 | minval(prism_coupler(cplid)%avect1%rAttr(nf,:)),& |
---|
946 | maxval(prism_coupler(cplid)%avect1%rAttr(nf,:)) |
---|
947 | call oasis_flush(nulprt) |
---|
948 | if (prism_coupler(cplid)%aVon(2)) & |
---|
949 | write(nulprt,*) subname,' DEBUG loctrans calc2 = ',cplid,nf,& |
---|
950 | minval(prism_coupler(cplid)%avect2%rAttr(nf,:)),& |
---|
951 | maxval(prism_coupler(cplid)%avect2%rAttr(nf,:)) |
---|
952 | if (prism_coupler(cplid)%aVon(3)) & |
---|
953 | write(nulprt,*) subname,' DEBUG loctrans calc3 = ',cplid,nf,& |
---|
954 | minval(prism_coupler(cplid)%avect3%rAttr(nf,:)),& |
---|
955 | maxval(prism_coupler(cplid)%avect3%rAttr(nf,:)) |
---|
956 | if (prism_coupler(cplid)%aVon(4)) & |
---|
957 | write(nulprt,*) subname,' DEBUG loctrans calc4 = ',cplid,nf,& |
---|
958 | minval(prism_coupler(cplid)%avect4%rAttr(nf,:)),& |
---|
959 | maxval(prism_coupler(cplid)%avect4%rAttr(nf,:)) |
---|
960 | if (prism_coupler(cplid)%aVon(5)) & |
---|
961 | write(nulprt,*) subname,' DEBUG loctrans calc5 = ',cplid,nf,& |
---|
962 | minval(prism_coupler(cplid)%avect5%rAttr(nf,:)),& |
---|
963 | maxval(prism_coupler(cplid)%avect5%rAttr(nf,:)) |
---|
964 | endif |
---|
965 | enddo |
---|
966 | call oasis_timer_stop(tstring) |
---|
967 | ENDIF |
---|
968 | |
---|
969 | !------------------------------------------------ |
---|
970 | ! past namcouple runtime (maxtime) no communication |
---|
971 | ! do restart if time+lag = maxtime, this assumes coupling |
---|
972 | ! period and lag and maxtime are all nicely consistent |
---|
973 | !------------------------------------------------ |
---|
974 | |
---|
975 | if (mseclag >= maxtime) then |
---|
976 | sndrcv = .false. ! turn off communication |
---|
977 | unpack = .false. ! nothing to unpack |
---|
978 | if (getput == OASIS3_PUT .and. lag > 0 .and. mseclag == maxtime) then |
---|
979 | kinfo = OASIS_ToRest |
---|
980 | call oasis_debug_note(subname//' lag restart write') |
---|
981 | write(tstring,F01) 'wrst_',cplid |
---|
982 | call oasis_timer_start(tstring) |
---|
983 | call oasis_io_write_avfile(rstfile,prism_coupler(cplid)%avect1, & |
---|
984 | prism_part(partid)%gsmap,nx,ny) |
---|
985 | if (prism_coupler(cplid)%aVon(2)) & |
---|
986 | call oasis_io_write_avfile(rstfile,prism_coupler(cplid)%avect2, & |
---|
987 | prism_part(partid)%gsmap,nx,ny,nampre='av2_') |
---|
988 | if (prism_coupler(cplid)%aVon(3)) & |
---|
989 | call oasis_io_write_avfile(rstfile,prism_coupler(cplid)%avect3, & |
---|
990 | prism_part(partid)%gsmap,nx,ny,nampre='av3_') |
---|
991 | if (prism_coupler(cplid)%aVon(4)) & |
---|
992 | call oasis_io_write_avfile(rstfile,prism_coupler(cplid)%avect4, & |
---|
993 | prism_part(partid)%gsmap,nx,ny,nampre='av4_') |
---|
994 | if (prism_coupler(cplid)%aVon(5)) & |
---|
995 | call oasis_io_write_avfile(rstfile,prism_coupler(cplid)%avect5, & |
---|
996 | prism_part(partid)%gsmap,nx,ny,nampre='av5_') |
---|
997 | call oasis_timer_stop(tstring) |
---|
998 | if (OASIS_debug >= 2) then |
---|
999 | write(nulprt,*) subname,' at ',msec,mseclag,' WRST: ', & |
---|
1000 | trim(mct_avect_exportRList2c(prism_coupler(cplid)%avect1)),' ',& |
---|
1001 | trim(rstfile) |
---|
1002 | call oasis_flush(nulprt) |
---|
1003 | endif |
---|
1004 | endif |
---|
1005 | endif |
---|
1006 | |
---|
1007 | !------------------------------------------------ |
---|
1008 | ! map and communicate operations |
---|
1009 | !------------------------------------------------ |
---|
1010 | |
---|
1011 | if (sndrcv) then |
---|
1012 | if (getput == OASIS3_PUT) then |
---|
1013 | kinfo = OASIS_sent |
---|
1014 | call oasis_debug_note(subname//' put section') |
---|
1015 | if (OASIS_debug >= 2) then |
---|
1016 | write(nulprt,*) subname,' at ',msec,mseclag,' SEND: ', & |
---|
1017 | trim(mct_avect_exportRList2c(prism_coupler(cplid)%avect1)) |
---|
1018 | call oasis_flush(nulprt) |
---|
1019 | endif |
---|
1020 | if (sndadd /= 0.0_ip_double_p .or. sndmult /= 1.0_ip_double_p) then |
---|
1021 | call oasis_debug_note(subname//' apply sndmult sndadd') |
---|
1022 | if (OASIS_debug >= 20) then |
---|
1023 | write(nulprt,*) subname,' DEBUG sndmult,add = ',sndmult,sndadd |
---|
1024 | write(nulprt,*) subname,' DEBUG put b4 sndmult,add = ',cplid,& |
---|
1025 | minval(prism_coupler(cplid)%avect1%rAttr),& |
---|
1026 | maxval(prism_coupler(cplid)%avect1%rAttr) |
---|
1027 | endif |
---|
1028 | prism_coupler(cplid)%avect1%rAttr(:,:) = prism_coupler(cplid)%avect1%rAttr(:,:)*sndmult & |
---|
1029 | + sndadd |
---|
1030 | endif |
---|
1031 | if (snddiag) call oasis_advance_avdiag(prism_coupler(cplid)%avect1,mpi_comm_local) |
---|
1032 | if (mapid > 0) then |
---|
1033 | write(tstring,F01) 'pmap_',cplid |
---|
1034 | call oasis_debug_note(subname//' put map') |
---|
1035 | if (OASIS_debug >= 20) then |
---|
1036 | write(nulprt,*) subname,' DEBUG put av11 b4 map = ',cplid,& |
---|
1037 | minval(prism_coupler(cplid)%avect1%rAttr),& |
---|
1038 | maxval(prism_coupler(cplid)%avect1%rAttr) |
---|
1039 | if (prism_coupler(cplid)%aVon(2)) & |
---|
1040 | write(nulprt,*) subname,' DEBUG put av2 b4 map = ',cplid,& |
---|
1041 | minval(prism_coupler(cplid)%avect2%rAttr),& |
---|
1042 | maxval(prism_coupler(cplid)%avect2%rAttr) |
---|
1043 | if (prism_coupler(cplid)%aVon(3)) & |
---|
1044 | write(nulprt,*) subname,' DEBUG put av3 b4 map = ',cplid,& |
---|
1045 | minval(prism_coupler(cplid)%avect3%rAttr),& |
---|
1046 | maxval(prism_coupler(cplid)%avect3%rAttr) |
---|
1047 | if (prism_coupler(cplid)%aVon(4)) & |
---|
1048 | write(nulprt,*) subname,' DEBUG put av4 b4 map = ',cplid,& |
---|
1049 | minval(prism_coupler(cplid)%avect4%rAttr),& |
---|
1050 | maxval(prism_coupler(cplid)%avect4%rAttr) |
---|
1051 | if (prism_coupler(cplid)%aVon(5)) & |
---|
1052 | write(nulprt,*) subname,' DEBUG put av5 b4 map = ',cplid,& |
---|
1053 | minval(prism_coupler(cplid)%avect5%rAttr),& |
---|
1054 | maxval(prism_coupler(cplid)%avect5%rAttr) |
---|
1055 | endif |
---|
1056 | call oasis_timer_start(tstring) |
---|
1057 | call mct_avect_zero(prism_coupler(cplid)%avect1m) |
---|
1058 | call oasis_advance_map(prism_coupler(cplid)%avect1, & |
---|
1059 | prism_coupler(cplid)%avect1m,prism_mapper(mapid),conserv,consbfb, & |
---|
1060 | prism_coupler(cplid)%aVon ,prism_coupler(cplid)%avect2, & |
---|
1061 | prism_coupler(cplid)%avect3,prism_coupler(cplid)%avect4, & |
---|
1062 | prism_coupler(cplid)%avect5) |
---|
1063 | call oasis_timer_stop(tstring) |
---|
1064 | write(tstring,F01) 'psnd_',cplid |
---|
1065 | call oasis_debug_note(subname//' put send') |
---|
1066 | if (OASIS_debug >= 20) then |
---|
1067 | write(nulprt,*) subname,' DEBUG put av1m b4 send = ',cplid,& |
---|
1068 | minval(prism_coupler(cplid)%avect1m%rAttr),& |
---|
1069 | maxval(prism_coupler(cplid)%avect1m%rAttr) |
---|
1070 | endif |
---|
1071 | call oasis_timer_start(tstring) |
---|
1072 | call mct_waitsend(prism_router(rouid)%router) |
---|
1073 | #if defined balance |
---|
1074 | WRITE(nulprt, FMT='(A,I2.2,A,F16.5)') & |
---|
1075 | 'Balance: ',cplid,' Before MPI put ', MPI_Wtime() |
---|
1076 | #endif |
---|
1077 | call mct_isend(prism_coupler(cplid)%avect1m,prism_router(rouid)%router,tag) |
---|
1078 | #if defined balance |
---|
1079 | WRITE(nulprt, FMT='(A,I2.2,A,F16.5)') & |
---|
1080 | 'Balance: ',cplid,' After MPI put ', MPI_Wtime() |
---|
1081 | #endif |
---|
1082 | call oasis_timer_stop(tstring) |
---|
1083 | ELSE |
---|
1084 | write(tstring,F01) 'psnd_',cplid |
---|
1085 | call oasis_debug_note(subname//' put send') |
---|
1086 | if (OASIS_debug >= 20) then |
---|
1087 | write(nulprt,*) subname,' DEBUG put av1 b4 send = ',cplid,& |
---|
1088 | minval(prism_coupler(cplid)%avect1%rAttr),& |
---|
1089 | maxval(prism_coupler(cplid)%avect1%rAttr) |
---|
1090 | endif |
---|
1091 | call oasis_timer_start(tstring) |
---|
1092 | call mct_waitsend(prism_router(rouid)%router) |
---|
1093 | #if defined balance |
---|
1094 | WRITE(nulprt, FMT='(A,I2.2,A,F16.5)') & |
---|
1095 | 'Balance: ',cplid,' Before MPI put ', MPI_Wtime() |
---|
1096 | #endif |
---|
1097 | call mct_isend(prism_coupler(cplid)%avect1,prism_router(rouid)%router,tag) |
---|
1098 | #if defined balance |
---|
1099 | WRITE(nulprt, FMT='(A,I2.2,A,F16.5)') & |
---|
1100 | 'Balance: ',cplid,' After MPI put ', MPI_Wtime() |
---|
1101 | #endif |
---|
1102 | call oasis_timer_stop(tstring) |
---|
1103 | endif |
---|
1104 | elseif (getput == OASIS3_GET) then |
---|
1105 | call oasis_debug_note(subname//' get section') |
---|
1106 | if (OASIS_debug >= 2 ) then |
---|
1107 | write(nulprt,*) subname,' at ',msec,mseclag,' RECV: ', & |
---|
1108 | trim(mct_avect_exportRList2c(prism_coupler(cplid)%avect1)) |
---|
1109 | call oasis_flush(nulprt) |
---|
1110 | endif |
---|
1111 | if (mapid > 0) then |
---|
1112 | call oasis_debug_note(subname//' get recv') |
---|
1113 | write(tstring,F01) 'grcv_',cplid |
---|
1114 | call oasis_timer_start(tstring) |
---|
1115 | call mct_avect_zero(prism_coupler(cplid)%avect1m) |
---|
1116 | #if defined balance |
---|
1117 | WRITE(nulprt, FMT='(A,I2.2,A,F16.5)') & |
---|
1118 | 'Balance: ',cplid,' Before MPI get ', MPI_Wtime() |
---|
1119 | #endif |
---|
1120 | call mct_recv(prism_coupler(cplid)%avect1m,prism_router(rouid)%router,tag) |
---|
1121 | #if defined balance |
---|
1122 | WRITE(nulprt, FMT='(A,I2.2,A,F16.5)') & |
---|
1123 | 'Balance: ',cplid,' After MPI get ', MPI_Wtime() |
---|
1124 | #endif |
---|
1125 | call oasis_timer_stop(tstring) |
---|
1126 | if (OASIS_debug >= 20) then |
---|
1127 | write(nulprt,*) subname,' DEBUG get af recv = ',cplid,& |
---|
1128 | minval(prism_coupler(cplid)%avect1m%rAttr),& |
---|
1129 | maxval(prism_coupler(cplid)%avect1m%rAttr) |
---|
1130 | endif |
---|
1131 | call oasis_debug_note(subname//' get map') |
---|
1132 | write(tstring,F01) 'gmap_',cplid |
---|
1133 | call oasis_timer_start(tstring) |
---|
1134 | call mct_avect_zero(prism_coupler(cplid)%avect1) |
---|
1135 | call oasis_advance_map(prism_coupler(cplid)%avect1m, & |
---|
1136 | prism_coupler(cplid)%avect1,prism_mapper(mapid),conserv,consbfb) |
---|
1137 | call oasis_timer_stop(tstring) |
---|
1138 | if (OASIS_debug >= 20) then |
---|
1139 | write(nulprt,*) subname,' DEBUG get af map = ',cplid,& |
---|
1140 | minval(prism_coupler(cplid)%avect1%rAttr),& |
---|
1141 | maxval(prism_coupler(cplid)%avect1%rAttr) |
---|
1142 | endif |
---|
1143 | else |
---|
1144 | write(tstring,F01) 'grcv_',cplid |
---|
1145 | call oasis_debug_note(subname//' get recv') |
---|
1146 | call oasis_timer_start(tstring) |
---|
1147 | #if defined balance |
---|
1148 | WRITE(nulprt, FMT='(A,I2.2,A,F16.5)') & |
---|
1149 | 'Balance: ',cplid,' Before MPI get ', MPI_Wtime() |
---|
1150 | #endif |
---|
1151 | call mct_recv(prism_coupler(cplid)%avect1,prism_router(rouid)%router,tag) |
---|
1152 | #if defined balance |
---|
1153 | WRITE(nulprt, FMT='(A,I2.2,A,F16.5)') & |
---|
1154 | 'Balance: ',cplid,' After MPI get ', MPI_Wtime() |
---|
1155 | #endif |
---|
1156 | call oasis_timer_stop(tstring) |
---|
1157 | if (OASIS_debug >= 20) then |
---|
1158 | write(nulprt,*) subname,' DEBUG get af recv = ',cplid,& |
---|
1159 | minval(prism_coupler(cplid)%avect1%rAttr),& |
---|
1160 | maxval(prism_coupler(cplid)%avect1%rAttr) |
---|
1161 | endif |
---|
1162 | endif |
---|
1163 | call oasis_debug_note(subname//' apply rcvmult rcvadd') |
---|
1164 | if (rcvadd /= 0.0_ip_double_p .or. rcvmult /= 1.0_ip_double_p) then |
---|
1165 | prism_coupler(cplid)%avect1%rAttr(:,:) = prism_coupler(cplid)%avect1%rAttr(:,:)*rcvmult & |
---|
1166 | + rcvadd |
---|
1167 | if (OASIS_debug >= 20) then |
---|
1168 | write(nulprt,*) subname,' DEBUG rcvmult,add = ',rcvmult,rcvadd |
---|
1169 | write(nulprt,*) subname,' DEBUG get af rcvmult,add = ',cplid,& |
---|
1170 | minval(prism_coupler(cplid)%avect1%rAttr),& |
---|
1171 | maxval(prism_coupler(cplid)%avect1%rAttr) |
---|
1172 | endif |
---|
1173 | endif |
---|
1174 | if (rcvdiag) call oasis_advance_avdiag(prism_coupler(cplid)%avect1,mpi_comm_local) |
---|
1175 | endif ! getput |
---|
1176 | endif ! sndrcv |
---|
1177 | |
---|
1178 | if (output) then |
---|
1179 | if (kinfo == OASIS_sent) then |
---|
1180 | kinfo = OASIS_sentout |
---|
1181 | elseif (kinfo == OASIS_torest) then |
---|
1182 | kinfo = OASIS_torestout |
---|
1183 | else |
---|
1184 | kinfo = OASIS_output |
---|
1185 | endif |
---|
1186 | call oasis_debug_note(subname//' output') |
---|
1187 | write(tstring,F01) 'wout_',cplid |
---|
1188 | call oasis_timer_start(tstring) |
---|
1189 | if (OASIS_debug >= 2) then |
---|
1190 | write(nulprt,*) subname,' at ',msec,mseclag,' WRIT: ', & |
---|
1191 | trim(mct_avect_exportRList2c(prism_coupler(cplid)%avect1)) |
---|
1192 | call oasis_flush(nulprt) |
---|
1193 | endif |
---|
1194 | write(fstring,'(A,I2.2)') '_'//trim(compnm)//'_',cplid |
---|
1195 | call oasis_io_write_avfbf(prism_coupler(cplid)%avect1,prism_part(partid)%gsmap, & |
---|
1196 | nx,ny,msec,fstring) |
---|
1197 | call oasis_timer_stop(tstring) |
---|
1198 | |
---|
1199 | if (OASIS_debug >= 30) then |
---|
1200 | call mct_avect_init(avtest,prism_coupler(cplid)%avect1,& |
---|
1201 | mct_aVect_lsize(prism_coupler(cplid)%avect1)) |
---|
1202 | write(tstring,F01) 'rinp_',cplid |
---|
1203 | call oasis_timer_start(tstring) |
---|
1204 | call oasis_io_read_avfbf(avtest,prism_part(partid)%gsmap,msec,fstring) |
---|
1205 | write(nulprt,*) subname,' DEBUG write/read test avfbf should be zero ',& |
---|
1206 | sum(prism_coupler(cplid)%avect1%rAttr-avtest%rAttr) |
---|
1207 | call mct_avect_clean(avtest) |
---|
1208 | call oasis_timer_stop(tstring) |
---|
1209 | endif |
---|
1210 | |
---|
1211 | endif |
---|
1212 | |
---|
1213 | !------------------------------------------------ |
---|
1214 | ! set avcnt, avect1, ltime, and status |
---|
1215 | !------------------------------------------------ |
---|
1216 | |
---|
1217 | call oasis_debug_note(subname//' reset status') |
---|
1218 | if (getput == OASIS3_PUT) then |
---|
1219 | prism_coupler(cplid)%ltime = msec |
---|
1220 | prism_coupler(cplid)%status(:) = OASIS_COMM_WAIT |
---|
1221 | prism_coupler(cplid)%avcnt(:) = 0 |
---|
1222 | call mct_avect_zero(prism_coupler(cplid)%avect1) |
---|
1223 | if (prism_coupler(cplid)%aVon(2)) & |
---|
1224 | call mct_avect_zero(prism_coupler(cplid)%avect2) |
---|
1225 | if (prism_coupler(cplid)%aVon(3)) & |
---|
1226 | call mct_avect_zero(prism_coupler(cplid)%avect3) |
---|
1227 | if (prism_coupler(cplid)%aVon(4)) & |
---|
1228 | call mct_avect_zero(prism_coupler(cplid)%avect4) |
---|
1229 | if (prism_coupler(cplid)%aVon(5)) & |
---|
1230 | call mct_avect_zero(prism_coupler(cplid)%avect5) |
---|
1231 | if (OASIS_debug >= 20) then |
---|
1232 | write(nulprt,*) subname,' DEBUG put reset status = ' |
---|
1233 | endif |
---|
1234 | elseif (getput == OASIS3_GET) then |
---|
1235 | prism_coupler(cplid)%ltime = msec |
---|
1236 | prism_coupler(cplid)%status(:) = OASIS_COMM_WAIT |
---|
1237 | if (OASIS_debug >= 20) then |
---|
1238 | write(nulprt,*) subname,' DEBUG get reset status = ' |
---|
1239 | endif |
---|
1240 | endif |
---|
1241 | |
---|
1242 | ELSE |
---|
1243 | |
---|
1244 | !------------------------------------------------ |
---|
1245 | ! no action, document |
---|
1246 | !------------------------------------------------ |
---|
1247 | |
---|
1248 | if (OASIS_debug >= 15) then |
---|
1249 | if (getput == OASIS3_PUT) then |
---|
1250 | write(nulprt,*) subname,' at ',msec,mseclag,' SKIP: ', & |
---|
1251 | trim(mct_avect_exportRList2c(prism_coupler(cplid)%avect1)) |
---|
1252 | elseif (getput == OASIS3_GET) then |
---|
1253 | write(nulprt,*) subname,' at ',msec,mseclag,' SKIP: ', & |
---|
1254 | trim(mct_avect_exportRList2c(prism_coupler(cplid)%avect1)) |
---|
1255 | endif |
---|
1256 | call oasis_flush(nulprt) |
---|
1257 | endif |
---|
1258 | |
---|
1259 | endif ! comm_now |
---|
1260 | |
---|
1261 | !------------------------------------------------ |
---|
1262 | ! sav non-instant loctrans operations for future restart |
---|
1263 | ! at the end of the run only |
---|
1264 | !------------------------------------------------ |
---|
1265 | |
---|
1266 | IF (mseclag + dt >= maxtime .AND. & |
---|
1267 | getput == OASIS3_PUT .and. prism_coupler(cplid)%trans /= ip_instant) then |
---|
1268 | call oasis_debug_note(subname//' loctrans restart write') |
---|
1269 | write(tstring,F01) 'wtrn_',cplid |
---|
1270 | call oasis_timer_start(tstring) |
---|
1271 | WRITE(vstring,'(a,i2.2,a)') 'loc',prism_coupler(cplid)%trans,'_cnt' |
---|
1272 | CALL oasis_io_write_array(rstfile,iarray=prism_coupler(cplid)%avcnt,& |
---|
1273 | ivarname=TRIM(vstring)) |
---|
1274 | write(vstring,'(a,i2.2,a)') 'loc',prism_coupler(cplid)%trans,'_' |
---|
1275 | CALL oasis_io_write_avfile(rstfile,prism_coupler(cplid)%avect1, & |
---|
1276 | prism_part(partid)%gsmap,nx,ny,nampre=TRIM(vstring)) |
---|
1277 | if (prism_coupler(cplid)%aVon(2)) then |
---|
1278 | write(vstring,'(a,i2.2,a)') 'av2loc',prism_coupler(cplid)%trans,'_' |
---|
1279 | CALL oasis_io_write_avfile(rstfile,prism_coupler(cplid)%avect2, & |
---|
1280 | prism_part(partid)%gsmap,nx,ny,nampre=TRIM(vstring)) |
---|
1281 | endif |
---|
1282 | if (prism_coupler(cplid)%aVon(3)) then |
---|
1283 | write(vstring,'(a,i2.2,a)') 'av3loc',prism_coupler(cplid)%trans,'_' |
---|
1284 | CALL oasis_io_write_avfile(rstfile,prism_coupler(cplid)%avect3, & |
---|
1285 | prism_part(partid)%gsmap,nx,ny,nampre=TRIM(vstring)) |
---|
1286 | endif |
---|
1287 | if (prism_coupler(cplid)%aVon(4)) then |
---|
1288 | write(vstring,'(a,i2.2,a)') 'av4loc',prism_coupler(cplid)%trans,'_' |
---|
1289 | CALL oasis_io_write_avfile(rstfile,prism_coupler(cplid)%avect4, & |
---|
1290 | prism_part(partid)%gsmap,nx,ny,nampre=TRIM(vstring)) |
---|
1291 | endif |
---|
1292 | if (prism_coupler(cplid)%aVon(5)) then |
---|
1293 | write(vstring,'(a,i2.2,a)') 'av5loc',prism_coupler(cplid)%trans,'_' |
---|
1294 | CALL oasis_io_write_avfile(rstfile,prism_coupler(cplid)%avect5, & |
---|
1295 | prism_part(partid)%gsmap,nx,ny,nampre=TRIM(vstring)) |
---|
1296 | endif |
---|
1297 | call oasis_timer_stop(tstring) |
---|
1298 | if (OASIS_debug >= 2) then |
---|
1299 | write(nulprt,*) subname,' at ',msec,mseclag,' WTRN: ', & |
---|
1300 | trim(mct_avect_exportRList2c(prism_coupler(cplid)%avect1)),' ',trim(rstfile) |
---|
1301 | call oasis_flush(nulprt) |
---|
1302 | endif |
---|
1303 | if (OASIS_debug >= 20) then |
---|
1304 | write(nulprt,*) subname,' DEBUG write loctrans restart',cplid,& |
---|
1305 | prism_coupler(cplid)%avcnt |
---|
1306 | write(nulprt,*) subname,' DEBUG write loctrans restart',cplid,& |
---|
1307 | minval(prism_coupler(cplid)%avect1%rAttr),& |
---|
1308 | maxval(prism_coupler(cplid)%avect1%rAttr) |
---|
1309 | endif |
---|
1310 | ENDIF |
---|
1311 | |
---|
1312 | !------------------------------------------------ |
---|
1313 | ! GET only, unpack avect1 if its coupling time |
---|
1314 | !------------------------------------------------ |
---|
1315 | |
---|
1316 | IF (getput == OASIS3_GET) THEN |
---|
1317 | IF (time_now .AND. unpack) THEN |
---|
1318 | if (kinfo == OASIS_output) then |
---|
1319 | kinfo = OASIS_recvout |
---|
1320 | elseif (kinfo == OASIS_fromrest) then |
---|
1321 | kinfo = OASIS_fromrestout |
---|
1322 | else |
---|
1323 | kinfo = OASIS_recvd |
---|
1324 | endif |
---|
1325 | if (input) then |
---|
1326 | kinfo = OASIS_input |
---|
1327 | call oasis_debug_note(subname//' input') |
---|
1328 | if (OASIS_debug >= 2) then |
---|
1329 | write(nulprt,*) subname,' at ',msec,mseclag,' READ: ', & |
---|
1330 | trim(mct_avect_exportRList2c(prism_coupler(cplid)%avect1)) |
---|
1331 | call oasis_flush(nulprt) |
---|
1332 | endif |
---|
1333 | write(tstring,F01) 'grin_',cplid |
---|
1334 | call oasis_timer_start(tstring) |
---|
1335 | if (trim(inpfile) /= trim(cspval)) then |
---|
1336 | call oasis_io_read_avfbf(prism_coupler(cplid)%avect1,prism_part(partid)%gsmap,& |
---|
1337 | msec,filename=trim(inpfile)) |
---|
1338 | else |
---|
1339 | fstring = '_'//trim(compnm) |
---|
1340 | call oasis_io_read_avfbf(prism_coupler(cplid)%avect1,prism_part(partid)%gsmap,& |
---|
1341 | msec,f_string=fstring) |
---|
1342 | endif |
---|
1343 | call oasis_timer_stop(tstring) |
---|
1344 | endif |
---|
1345 | if (OASIS_debug >= 2) then |
---|
1346 | write(nulprt,*) subname,' at ',msec,mseclag,' UPCK: ',trim(vname) |
---|
1347 | endif |
---|
1348 | write(tstring,F01) 'gcpy_',cplid |
---|
1349 | call oasis_debug_note(subname//' get copy to array') |
---|
1350 | call oasis_timer_start(tstring) |
---|
1351 | if (present(array1dout)) array1dout(:) = & |
---|
1352 | prism_coupler(cplid)%avect1%rAttr(nfav,:) |
---|
1353 | if (present(array2dout)) array2dout(:,:) = & |
---|
1354 | RESHAPE(prism_coupler(cplid)%avect1%rAttr(nfav,:),SHAPE(array2dout)) |
---|
1355 | call oasis_timer_stop(tstring) |
---|
1356 | if (OASIS_debug >= 20) then |
---|
1357 | if (present(array1dout)) write(nulprt,*) subname,' DEBUG array copy = ',& |
---|
1358 | cplid,minval(array1dout),maxval(array1dout) |
---|
1359 | if (present(array2dout)) write(nulprt,*) subname,' DEBUG array copy = ',& |
---|
1360 | cplid,minval(array2dout),maxval(array2dout) |
---|
1361 | endif |
---|
1362 | ENDIF |
---|
1363 | if (time_now) prism_coupler(cplid)%status(nfav) = OASIS_COMM_READY |
---|
1364 | ENDIF |
---|
1365 | |
---|
1366 | !------------------------------------------------ |
---|
1367 | ! always remember last id and last coupler time |
---|
1368 | !------------------------------------------------ |
---|
1369 | |
---|
1370 | lcouplerid = cplid |
---|
1371 | lcouplertime = msec |
---|
1372 | |
---|
1373 | if (OASIS_debug >= 2) then |
---|
1374 | write(nulprt,*) subname,' at ',msec,mseclag,' KINF: ',trim(vname),kinfo |
---|
1375 | endif |
---|
1376 | |
---|
1377 | enddo ! nc = 1,var%ncpl |
---|
1378 | |
---|
1379 | call oasis_debug_exit(subname) |
---|
1380 | |
---|
1381 | END SUBROUTINE oasis_advance_run |
---|
1382 | |
---|
1383 | |
---|
1384 | !------------------------------------------------------------------- |
---|
1385 | |
---|
1386 | SUBROUTINE oasis_advance_map(av1,avd,mapper,conserv,consbfb,& |
---|
1387 | avon,av2,av3,av4,av5) |
---|
1388 | |
---|
1389 | ! NOTE: mask = 0 is active point according to oasis3 conserv.f |
---|
1390 | |
---|
1391 | implicit none |
---|
1392 | type(mct_aVect) ,intent(in) :: av1 ! source av |
---|
1393 | type(mct_aVect) ,intent(inout) :: avd ! dst av |
---|
1394 | type(prism_mapper_type),intent(inout) :: mapper ! prism_mapper |
---|
1395 | integer(kind=ip_i4_p) ,intent(in),optional :: conserv ! conserv flag |
---|
1396 | logical ,intent(in),optional :: consbfb ! conserv bfb option |
---|
1397 | logical ,intent(in),optional :: avon(:) ! which source avs are on |
---|
1398 | type(mct_aVect) ,intent(in),optional :: av2 ! source av2 |
---|
1399 | type(mct_aVect) ,intent(in),optional :: av3 ! source av2 |
---|
1400 | type(mct_aVect) ,intent(in),optional :: av4 ! source av2 |
---|
1401 | type(mct_aVect) ,intent(in),optional :: av5 ! source av2 |
---|
1402 | |
---|
1403 | integer(kind=ip_i4_p) :: fsize,lsizes,lsized,nf,ni,n,m |
---|
1404 | real(kind=ip_r8_p) :: sumtmp, wts_sums, wts_sumd, zradi, zlagr |
---|
1405 | integer(kind=ip_i4_p),allocatable :: imasks(:),imaskd(:) |
---|
1406 | real(kind=ip_r8_p),allocatable :: areas(:),aread(:) |
---|
1407 | real(kind=ip_r8_p),allocatable :: av_sums(:),av_sumd(:) ! local sums |
---|
1408 | character(len=ic_med) :: tstring ! timer label string |
---|
1409 | type(mct_aVect) :: avdtmp ! for summing multiple mapping weights |
---|
1410 | type(mct_aVect) :: av2g ! for bfb sums |
---|
1411 | logical :: lconsbfb |
---|
1412 | integer(kind=ip_i4_p),parameter :: avsmax = prism_coupler_avsmax |
---|
1413 | logical :: locavon(avsmax) ! local avon |
---|
1414 | integer(kind=ip_i4_p) :: avonsize, nterm |
---|
1415 | character(len=*),parameter :: subname = 'oasis_advance_map' |
---|
1416 | |
---|
1417 | call oasis_debug_enter(subname) |
---|
1418 | |
---|
1419 | lconsbfb = .true. |
---|
1420 | if (present(consbfb)) then |
---|
1421 | lconsbfb = consbfb |
---|
1422 | endif |
---|
1423 | |
---|
1424 | !--- assume avon and av2-5 are not passed but av1 always is --- |
---|
1425 | avonsize = 1 |
---|
1426 | nterm = 1 |
---|
1427 | locavon = .false. |
---|
1428 | locavon(1) = .true. |
---|
1429 | |
---|
1430 | !--- but if avon is passed, use avon flags --- |
---|
1431 | if (present(avon)) then |
---|
1432 | avonsize = size(avon) |
---|
1433 | nterm = min(avsmax,avonsize) |
---|
1434 | locavon(1:nterm) = avon(1:nterm) |
---|
1435 | else |
---|
1436 | !--- if avon is not passed, av2-5 should not be --- |
---|
1437 | if (present(av2) .or. present(av3) .or. present(av4) .or. present(av5)) then |
---|
1438 | WRITE(nulprt,*) subname,' ERROR av2-5 passed but avon not passed' |
---|
1439 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
1440 | CALL oasis_flush(nulprt) |
---|
1441 | CALL oasis_abort_noarg() |
---|
1442 | endif |
---|
1443 | endif |
---|
1444 | |
---|
1445 | ! check consistency between weights and coupling terms |
---|
1446 | do n = 1,nterm |
---|
1447 | if (locavon(n) .and. n > mapper%nwgts) then |
---|
1448 | WRITE(nulprt,*) subname,' ERROR in nwgts and coupling terms',mapper%nwgts,n |
---|
1449 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
1450 | CALL oasis_flush(nulprt) |
---|
1451 | CALL oasis_abort_noarg() |
---|
1452 | endif |
---|
1453 | enddo |
---|
1454 | |
---|
1455 | |
---|
1456 | if (locavon(1)) then |
---|
1457 | if (mct_avect_nRattr(av1) /= mct_avect_nRattr(avd)) then |
---|
1458 | WRITE(nulprt,*) subname,' ERROR in av1 num of flds' |
---|
1459 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
1460 | CALL oasis_flush(nulprt) |
---|
1461 | CALL oasis_abort_noarg() |
---|
1462 | endif |
---|
1463 | call mct_sMat_avMult(av1, mapper%sMatP(1), avd) |
---|
1464 | endif |
---|
1465 | |
---|
1466 | if (locavon(2).or.locavon(3).or.locavon(4).or.locavon(5)) then |
---|
1467 | lsized = mct_avect_lsize(avd) |
---|
1468 | call mct_aVect_init(avdtmp,avd,lsized) |
---|
1469 | |
---|
1470 | if (locavon(2)) then |
---|
1471 | if (mct_avect_nRattr(av2) /= mct_avect_nRattr(avd)) then |
---|
1472 | WRITE(nulprt,*) subname,' ERROR in av2 num of flds' |
---|
1473 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
1474 | CALL oasis_flush(nulprt) |
---|
1475 | CALL oasis_abort_noarg() |
---|
1476 | endif |
---|
1477 | call mct_sMat_avMult(av2, mapper%sMatP(2), avdtmp) |
---|
1478 | avd%rAttr = avd%rAttr + avdtmp%rAttr |
---|
1479 | endif |
---|
1480 | |
---|
1481 | if (locavon(3)) then |
---|
1482 | if (mct_avect_nRattr(av3) /= mct_avect_nRattr(avd)) then |
---|
1483 | WRITE(nulprt,*) subname,' ERROR in av3 num of flds' |
---|
1484 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
1485 | CALL oasis_flush(nulprt) |
---|
1486 | CALL oasis_abort_noarg() |
---|
1487 | endif |
---|
1488 | call mct_sMat_avMult(av3, mapper%sMatP(3), avdtmp) |
---|
1489 | avd%rAttr = avd%rAttr + avdtmp%rAttr |
---|
1490 | endif |
---|
1491 | |
---|
1492 | if (locavon(4)) then |
---|
1493 | if (mct_avect_nRattr(av4) /= mct_avect_nRattr(avd)) then |
---|
1494 | WRITE(nulprt,*) subname,' ERROR in av4 num of flds' |
---|
1495 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
1496 | CALL oasis_flush(nulprt) |
---|
1497 | CALL oasis_abort_noarg() |
---|
1498 | endif |
---|
1499 | call mct_sMat_avMult(av4, mapper%sMatP(4), avdtmp) |
---|
1500 | avd%rAttr = avd%rAttr + avdtmp%rAttr |
---|
1501 | endif |
---|
1502 | |
---|
1503 | if (locavon(5)) then |
---|
1504 | if (mct_avect_nRattr(av5) /= mct_avect_nRattr(avd)) then |
---|
1505 | WRITE(nulprt,*) subname,' ERROR in av5 num of flds' |
---|
1506 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
1507 | CALL oasis_flush(nulprt) |
---|
1508 | CALL oasis_abort_noarg() |
---|
1509 | endif |
---|
1510 | call mct_sMat_avMult(av5, mapper%sMatP(5), avdtmp) |
---|
1511 | avd%rAttr = avd%rAttr + avdtmp%rAttr |
---|
1512 | endif |
---|
1513 | |
---|
1514 | call mct_aVect_clean(avdtmp) |
---|
1515 | endif |
---|
1516 | |
---|
1517 | |
---|
1518 | call oasis_debug_note(subname//' map') |
---|
1519 | |
---|
1520 | IF (PRESENT(conserv)) THEN |
---|
1521 | call oasis_debug_note(subname//' conserv') |
---|
1522 | IF (conserv /= ip_cnone) THEN |
---|
1523 | fsize = mct_avect_nRattr(av1) |
---|
1524 | allocate(av_sums(fsize),av_sumd(fsize)) |
---|
1525 | |
---|
1526 | zradi = 1./(eradius*eradius) |
---|
1527 | |
---|
1528 | !------------------- |
---|
1529 | ! extract mask and area and compute sum of masked area for source |
---|
1530 | !------------------- |
---|
1531 | lsizes = mct_avect_lsize(mapper%av_ms) |
---|
1532 | allocate(imasks(lsizes),areas(lsizes)) |
---|
1533 | nf = mct_aVect_indexIA(mapper%av_ms,'mask') |
---|
1534 | imasks(:) = mapper%av_ms%iAttr(nf,:) |
---|
1535 | nf = mct_aVect_indexRA(mapper%av_ms,'area') |
---|
1536 | areas(:) = mapper%av_ms%rAttr(nf,:)*zradi |
---|
1537 | |
---|
1538 | if (lconsbfb) then |
---|
1539 | call mct_avect_gather(mapper%av_ms,av2g,prism_part(mapper%spart)%gsmap,& |
---|
1540 | 0,mpi_comm_local) |
---|
1541 | wts_sums = 0.0_ip_r8_p |
---|
1542 | if (mpi_rank_local == 0) then |
---|
1543 | ni = mct_aVect_indexIA(av2g,'mask') |
---|
1544 | nf = mct_aVect_indexRA(av2g,'area') |
---|
1545 | do n = 1,mct_avect_lsize(av2g) |
---|
1546 | if (av2g%iAttr(ni,n) == 0) wts_sums = wts_sums + av2g%rAttr(nf,n)*zradi |
---|
1547 | enddo |
---|
1548 | endif |
---|
1549 | call oasis_mpi_bcast(wts_sums,mpi_comm_local,subname//" bcast wts_sums") |
---|
1550 | if (mpi_rank_local == 0) then |
---|
1551 | call mct_avect_clean(av2g) |
---|
1552 | endif |
---|
1553 | else |
---|
1554 | sumtmp = 0.0_ip_r8_p |
---|
1555 | do n = 1,lsizes |
---|
1556 | if (imasks(n) == 0) sumtmp = sumtmp + areas(n) |
---|
1557 | enddo |
---|
1558 | call oasis_mpi_sum(sumtmp,wts_sums,mpi_comm_local,string=subname//':wts_sums',& |
---|
1559 | all=.true.) |
---|
1560 | endif |
---|
1561 | |
---|
1562 | !------------------- |
---|
1563 | ! extract mask and area and compute sum of masked area for destination |
---|
1564 | !------------------- |
---|
1565 | lsized = mct_avect_lsize(mapper%av_md) |
---|
1566 | allocate(imaskd(lsized),aread(lsized)) |
---|
1567 | nf = mct_aVect_indexIA(mapper%av_md,'mask') |
---|
1568 | imaskd(:) = mapper%av_md%iAttr(nf,:) |
---|
1569 | nf = mct_aVect_indexRA(mapper%av_md,'area') |
---|
1570 | aread(:) = mapper%av_md%rAttr(nf,:)*zradi |
---|
1571 | |
---|
1572 | if (lconsbfb) then |
---|
1573 | call mct_avect_gather(mapper%av_md,av2g,prism_part(mapper%dpart)%gsmap,0,mpi_comm_local) |
---|
1574 | wts_sumd = 0.0_ip_r8_p |
---|
1575 | if (mpi_rank_local == 0) then |
---|
1576 | ni = mct_aVect_indexIA(av2g,'mask') |
---|
1577 | nf = mct_aVect_indexRA(av2g,'area') |
---|
1578 | do n = 1,mct_avect_lsize(av2g) |
---|
1579 | if (av2g%iAttr(ni,n) == 0) wts_sumd = wts_sumd + av2g%rAttr(nf,n)*zradi |
---|
1580 | enddo |
---|
1581 | endif |
---|
1582 | call oasis_mpi_bcast(wts_sumd,mpi_comm_local,subname//" bcast wts_sumd") |
---|
1583 | if (mpi_rank_local == 0) then |
---|
1584 | call mct_avect_clean(av2g) |
---|
1585 | endif |
---|
1586 | else |
---|
1587 | sumtmp = 0.0_ip_r8_p |
---|
1588 | do n = 1,lsized |
---|
1589 | if (imaskd(n) == 0) sumtmp = sumtmp + aread(n) |
---|
1590 | enddo |
---|
1591 | call oasis_mpi_sum(sumtmp,wts_sumd,mpi_comm_local,string=subname//':wts_sumd',all=.true.) |
---|
1592 | endif |
---|
1593 | |
---|
1594 | if (OASIS_debug >= 30) then |
---|
1595 | write(nulprt,*) subname,' DEBUG conserve src mask ',minval(imasks),& |
---|
1596 | maxval(imasks),sum(imasks) |
---|
1597 | write(nulprt,*) subname,' DEBUG conserve dst mask ',minval(imaskd),& |
---|
1598 | maxval(imaskd),sum(imaskd) |
---|
1599 | write(nulprt,*) subname,' DEBUG conserve src area ',minval(areas),& |
---|
1600 | maxval(areas),sum(areas) |
---|
1601 | write(nulprt,*) subname,' DEBUG conserve dst area ',minval(aread),& |
---|
1602 | maxval(aread),sum(aread) |
---|
1603 | write(nulprt,*) subname,' DEBUG conserve wts_sum ',wts_sums,wts_sumd |
---|
1604 | endif |
---|
1605 | |
---|
1606 | !------------------- |
---|
1607 | ! compute global sums of av1 |
---|
1608 | ! assume av1 is the thing to be conserved |
---|
1609 | !------------------- |
---|
1610 | call oasis_advance_avsum(av1,av_sums,prism_part(mapper%spart)%gsmap,mpi_comm_local, & |
---|
1611 | mask=imasks,wts=areas,consbfb=lconsbfb) |
---|
1612 | call oasis_advance_avsum(avd,av_sumd,prism_part(mapper%dpart)%gsmap,mpi_comm_local, & |
---|
1613 | mask=imaskd,wts=aread,consbfb=lconsbfb) |
---|
1614 | |
---|
1615 | if (OASIS_debug >= 20) then |
---|
1616 | write(nulprt,*) subname,' DEBUG src sum b4 conserve ',av_sums |
---|
1617 | write(nulprt,*) subname,' DEBUG dst sum b4 conserve ',av_sumd |
---|
1618 | endif |
---|
1619 | |
---|
1620 | if (conserv == ip_cglobal) then |
---|
1621 | if (wts_sumd == 0.0_ip_r8_p) then |
---|
1622 | WRITE(nulprt,*) subname,' ERROR: conserve global wts_sumd/sums zero' |
---|
1623 | WRITE(nulprt,*) subname,' abort by model :',compid,& |
---|
1624 | ' proc :',mpi_rank_local |
---|
1625 | CALL oasis_flush(nulprt) |
---|
1626 | CALL oasis_abort_noarg() |
---|
1627 | endif |
---|
1628 | do m = 1,fsize |
---|
1629 | zlagr = (av_sumd(m) - av_sums(m)) / wts_sumd |
---|
1630 | do n = 1,lsized |
---|
1631 | if (imaskd(n) == 0) avd%rAttr(m,n) = avd%rAttr(m,n) - zlagr |
---|
1632 | enddo |
---|
1633 | enddo |
---|
1634 | elseif (conserv == ip_cglbpos) then |
---|
1635 | do m = 1,fsize |
---|
1636 | if (av_sumd(m) == 0.0_ip_r8_p .and. av_sums(m) /= 0.0_ip_r8_p) then |
---|
1637 | WRITE(nulprt,*) subname,' ERROR: conserve cglbpos av_sumd/sums' |
---|
1638 | WRITE(nulprt,*) subname,' abort by model :',compid,& |
---|
1639 | ' proc :',mpi_rank_local |
---|
1640 | CALL oasis_flush(nulprt) |
---|
1641 | CALL oasis_abort_noarg() |
---|
1642 | elseif (av_sumd(m) /= 0.0_ip_r8_p) then |
---|
1643 | zlagr = av_sums(m) / av_sumd(m) |
---|
1644 | do n = 1,lsized |
---|
1645 | if (imaskd(n) == 0) avd%rAttr(m,n) = avd%rAttr(m,n) * zlagr |
---|
1646 | enddo |
---|
1647 | endif |
---|
1648 | enddo |
---|
1649 | elseif (conserv == ip_cbasbal) then |
---|
1650 | if (wts_sumd == 0.0_ip_r8_p .or. wts_sums == 0.0_ip_r8_p) then |
---|
1651 | WRITE(nulprt,*) subname,' ERROR: conserve wts_sumd/sums zero' |
---|
1652 | WRITE(nulprt,*) subname,' abort by model :',compid,& |
---|
1653 | ' proc :',mpi_rank_local |
---|
1654 | CALL oasis_flush(nulprt) |
---|
1655 | CALL oasis_abort_noarg() |
---|
1656 | endif |
---|
1657 | do m = 1,fsize |
---|
1658 | zlagr = (av_sumd(m) - (av_sums(m)*(wts_sumd/wts_sums))) / wts_sumd |
---|
1659 | do n = 1,lsized |
---|
1660 | if (imaskd(n) == 0) avd%rAttr(m,n) = avd%rAttr(m,n) - zlagr |
---|
1661 | enddo |
---|
1662 | enddo |
---|
1663 | elseif (conserv == ip_cbaspos) then |
---|
1664 | do m = 1,fsize |
---|
1665 | if (av_sumd(m) == 0.0_ip_r8_p .and. av_sums(m) /= 0.0_ip_r8_p) then |
---|
1666 | WRITE(nulprt,*) subname,' ERROR: conserve cglbpos av_sumd/sums' |
---|
1667 | WRITE(nulprt,*) subname,' abort by model :',compid,& |
---|
1668 | ' proc :',mpi_rank_local |
---|
1669 | CALL oasis_flush(nulprt) |
---|
1670 | CALL oasis_abort_noarg() |
---|
1671 | elseif (av_sumd(m) /= 0.0_ip_r8_p) then |
---|
1672 | zlagr = (av_sums(m)/av_sumd(m)) * (wts_sumd/wts_sums) |
---|
1673 | do n = 1,lsized |
---|
1674 | if (imaskd(n) == 0) avd%rAttr(m,n) = avd%rAttr(m,n) * zlagr |
---|
1675 | enddo |
---|
1676 | endif |
---|
1677 | enddo |
---|
1678 | else |
---|
1679 | WRITE(nulprt,*) subname,' ERROR: conserv option' |
---|
1680 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
1681 | CALL oasis_flush(nulprt) |
---|
1682 | CALL oasis_abort_noarg() |
---|
1683 | endif |
---|
1684 | |
---|
1685 | if (OASIS_debug >= 20) then |
---|
1686 | call oasis_advance_avsum(av1,av_sums,prism_part(mapper%spart)%gsmap,mpi_comm_local, & |
---|
1687 | mask=imasks,wts=areas,consbfb=lconsbfb) |
---|
1688 | call oasis_advance_avsum(avd,av_sumd,prism_part(mapper%dpart)%gsmap,mpi_comm_local, & |
---|
1689 | mask=imaskd,wts=aread,consbfb=lconsbfb) |
---|
1690 | write(nulprt,*) subname,' DEBUG src sum af conserve ',av_sums |
---|
1691 | write(nulprt,*) subname,' DEBUG dst sum af conserve ',av_sumd |
---|
1692 | CALL oasis_flush(nulprt) |
---|
1693 | endif |
---|
1694 | |
---|
1695 | deallocate(imasks,imaskd,areas,aread) |
---|
1696 | deallocate(av_sums,av_sumd) |
---|
1697 | ENDIF |
---|
1698 | ENDIF ! present conserve |
---|
1699 | |
---|
1700 | call oasis_debug_exit(subname) |
---|
1701 | |
---|
1702 | END SUBROUTINE oasis_advance_map |
---|
1703 | |
---|
1704 | !------------------------------------------------------------------- |
---|
1705 | |
---|
1706 | SUBROUTINE oasis_advance_avsum(av,sum,gsmap,mpicom,mask,wts,consbfb) |
---|
1707 | |
---|
1708 | implicit none |
---|
1709 | type(mct_aVect) ,intent(in) :: av ! av |
---|
1710 | real(kind=ip_r8_p) ,intent(inout) :: sum(:) ! sum of av fields |
---|
1711 | type(mct_gsMap) ,intent(in) :: gsmap ! gsmap associate with av |
---|
1712 | integer(kind=ip_i4_p),intent(in) :: mpicom ! mpicom |
---|
1713 | integer(kind=ip_i4_p),intent(in),optional :: mask(:) ! mask to apply to av |
---|
1714 | real(kind=ip_r8_p) ,intent(in),optional :: wts(:) ! wts to apply to av |
---|
1715 | logical ,intent(in),optional :: consbfb ! bfb conserve |
---|
1716 | |
---|
1717 | integer(kind=ip_i4_p) :: n,m,ierr,mytask |
---|
1718 | integer(kind=ip_i4_p) :: lsize,fsize ! local size of av, number of flds in av |
---|
1719 | real(kind=ip_r8_p),allocatable :: lsum(:) ! local sums |
---|
1720 | real(kind=ip_r8_p),allocatable :: lwts(:) ! local wts taking into account mask and wts |
---|
1721 | type(mct_aVect) :: av1, av1g ! use av1,av1g for gather and bfb sum |
---|
1722 | logical :: lconsbfb ! local conserve bfb |
---|
1723 | character(len=*),parameter :: subname = 'oasis_advance_avsum' |
---|
1724 | |
---|
1725 | call oasis_debug_enter(subname) |
---|
1726 | |
---|
1727 | lconsbfb = .true. |
---|
1728 | if (present(consbfb)) then |
---|
1729 | lconsbfb = consbfb |
---|
1730 | endif |
---|
1731 | |
---|
1732 | fsize = mct_avect_nRattr(av) |
---|
1733 | lsize = mct_avect_lsize(av) |
---|
1734 | |
---|
1735 | allocate(lsum(fsize)) |
---|
1736 | lsum = 0.0_ip_r8_p |
---|
1737 | allocate(lwts(lsize)) |
---|
1738 | lwts = 1.0_ip_r8_p |
---|
1739 | |
---|
1740 | if (size(sum) /= fsize) then |
---|
1741 | WRITE(nulprt,*) subname,' ERROR: size sum ne size av' |
---|
1742 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
1743 | CALL oasis_flush(nulprt) |
---|
1744 | CALL oasis_abort_noarg() |
---|
1745 | endif |
---|
1746 | |
---|
1747 | if (present(mask)) then |
---|
1748 | if (size(mask) /= lsize) then |
---|
1749 | WRITE(nulprt,*) subname,' ERROR: size mask ne size av' |
---|
1750 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
1751 | CALL oasis_flush(nulprt) |
---|
1752 | CALL oasis_abort_noarg() |
---|
1753 | endif |
---|
1754 | do n = 1,lsize |
---|
1755 | if (mask(n) /= 0) lwts(n) = 0.0_ip_r8_p |
---|
1756 | enddo |
---|
1757 | endif |
---|
1758 | |
---|
1759 | if (present(wts)) then |
---|
1760 | if (size(wts) /= lsize) then |
---|
1761 | WRITE(nulprt,*) subname,' ERROR: size wts ne size av' |
---|
1762 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
1763 | CALL oasis_flush(nulprt) |
---|
1764 | CALL oasis_abort_noarg() |
---|
1765 | endif |
---|
1766 | do n = 1,lsize |
---|
1767 | lwts(n) = lwts(n) * wts(n) |
---|
1768 | enddo |
---|
1769 | endif |
---|
1770 | |
---|
1771 | if (lconsbfb) then |
---|
1772 | call mct_avect_init(av1,av,lsize) |
---|
1773 | do n = 1,lsize |
---|
1774 | do m = 1,fsize |
---|
1775 | av1%rAttr(m,n) = av%rAttr(m,n)*lwts(n) |
---|
1776 | enddo |
---|
1777 | enddo |
---|
1778 | call mct_avect_gather(av1,av1g,gsmap,0,mpicom) |
---|
1779 | call MPI_COMM_RANK(mpicom,mytask,ierr) |
---|
1780 | sum = 0.0_ip_r8_p |
---|
1781 | if (mytask == 0) then |
---|
1782 | do n = 1,mct_avect_lsize(av1g) |
---|
1783 | do m = 1,fsize |
---|
1784 | sum(m) = sum(m) + av1g%rAttr(m,n) |
---|
1785 | enddo |
---|
1786 | enddo |
---|
1787 | endif |
---|
1788 | call oasis_mpi_bcast(sum,mpicom,subname//" bcast sum") |
---|
1789 | call mct_avect_clean(av1) |
---|
1790 | if (mytask == 0) then |
---|
1791 | call mct_avect_clean(av1g) |
---|
1792 | endif |
---|
1793 | else |
---|
1794 | lsum = 0.0_ip_r8_p |
---|
1795 | do n = 1,lsize |
---|
1796 | do m = 1,fsize |
---|
1797 | lsum(m) = lsum(m) + av%rAttr(m,n)*lwts(n) |
---|
1798 | enddo |
---|
1799 | enddo |
---|
1800 | call oasis_mpi_sum(lsum,sum,mpicom,string=trim(subname)//':sum',all=.true.) |
---|
1801 | endif |
---|
1802 | |
---|
1803 | deallocate(lsum) |
---|
1804 | deallocate(lwts) |
---|
1805 | |
---|
1806 | call oasis_debug_exit(subname) |
---|
1807 | |
---|
1808 | END SUBROUTINE oasis_advance_avsum |
---|
1809 | |
---|
1810 | !------------------------------------------------------------------- |
---|
1811 | |
---|
1812 | SUBROUTINE oasis_advance_avdiag(av,mpicom,mask,wts) |
---|
1813 | |
---|
1814 | implicit none |
---|
1815 | type(mct_aVect) ,intent(in) :: av ! av |
---|
1816 | integer(kind=ip_i4_p),intent(in) :: mpicom ! mpicom |
---|
1817 | integer(kind=ip_i4_p),intent(in),optional :: mask(:) ! mask to apply to av |
---|
1818 | real(kind=ip_r8_p) ,intent(in),optional :: wts(:) ! wts to apply to av |
---|
1819 | |
---|
1820 | integer(kind=ip_i4_p) :: n,m,ierr,mype |
---|
1821 | integer(kind=ip_i4_p) :: lsize,fsize ! local size of av, number of flds in av |
---|
1822 | logical :: first_call |
---|
1823 | real(kind=ip_r8_p) :: lval ! local temporary |
---|
1824 | real(kind=ip_r8_p),allocatable :: lsum(:) ! local sum |
---|
1825 | real(kind=ip_r8_p),allocatable :: lmin(:) ! local min |
---|
1826 | real(kind=ip_r8_p),allocatable :: lmax(:) ! local max |
---|
1827 | real(kind=ip_r8_p),allocatable :: gsum(:) ! global sum |
---|
1828 | real(kind=ip_r8_p),allocatable :: gmin(:) ! global min |
---|
1829 | real(kind=ip_r8_p),allocatable :: gmax(:) ! global max |
---|
1830 | real(kind=ip_r8_p),allocatable :: lwts(:) ! local wts taking into account mask and wts |
---|
1831 | type(mct_string) :: mstring ! mct char type |
---|
1832 | character(len=64):: itemc ! string converted to char |
---|
1833 | character(len=*),parameter :: subname = 'oasis_advance_avdiag' |
---|
1834 | |
---|
1835 | call oasis_debug_enter(subname) |
---|
1836 | |
---|
1837 | fsize = mct_avect_nRattr(av) |
---|
1838 | lsize = mct_avect_lsize(av) |
---|
1839 | |
---|
1840 | allocate(lsum(fsize)) |
---|
1841 | allocate(lmin(fsize)) |
---|
1842 | allocate(lmax(fsize)) |
---|
1843 | allocate(gsum(fsize)) |
---|
1844 | allocate(gmin(fsize)) |
---|
1845 | allocate(gmax(fsize)) |
---|
1846 | |
---|
1847 | allocate(lwts(lsize)) |
---|
1848 | lwts = 1.0_ip_r8_p |
---|
1849 | !!$ lmin=HUGE(lwts) |
---|
1850 | !!$ lmax=-lmin |
---|
1851 | if (present(mask)) then |
---|
1852 | if (size(mask) /= lsize) then |
---|
1853 | WRITE(nulprt,*) subname,' ERROR: size mask ne size av' |
---|
1854 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
1855 | CALL oasis_flush(nulprt) |
---|
1856 | CALL oasis_abort_noarg() |
---|
1857 | endif |
---|
1858 | do n = 1,lsize |
---|
1859 | if (mask(n) /= 0) lwts(n) = 0.0_ip_r8_p |
---|
1860 | enddo |
---|
1861 | endif |
---|
1862 | |
---|
1863 | if (present(wts)) then |
---|
1864 | if (size(wts) /= lsize) then |
---|
1865 | WRITE(nulprt,*) subname,' ERROR: size wts ne size av' |
---|
1866 | WRITE(nulprt,*) subname,' abort by model :',compid,' proc :',mpi_rank_local |
---|
1867 | CALL oasis_flush(nulprt) |
---|
1868 | CALL oasis_abort_noarg() |
---|
1869 | endif |
---|
1870 | do n = 1,lsize |
---|
1871 | lwts(n) = lwts(n) * wts(n) |
---|
1872 | enddo |
---|
1873 | endif |
---|
1874 | |
---|
1875 | lsum = 0.0_ip_r8_p |
---|
1876 | do m = 1,fsize |
---|
1877 | first_call = .true. |
---|
1878 | do n = 1,lsize |
---|
1879 | lval = av%rAttr(m,n)*lwts(n) |
---|
1880 | lsum(m) = lsum(m) + lval |
---|
1881 | if (lwts(n) /= 0.0_ip_r8_p) then |
---|
1882 | if (first_call) then |
---|
1883 | lmin(m) = lval |
---|
1884 | lmax(m) = lval |
---|
1885 | first_call = .false. |
---|
1886 | else |
---|
1887 | lmin(m) = min(lmin(m),lval) |
---|
1888 | lmax(m) = max(lmax(m),lval) |
---|
1889 | endif |
---|
1890 | endif |
---|
1891 | enddo |
---|
1892 | enddo |
---|
1893 | |
---|
1894 | mype = -1 |
---|
1895 | if (mpicom /= MPI_COMM_NULL) then |
---|
1896 | call MPI_COMM_RANK(mpicom,mype,ierr) |
---|
1897 | call oasis_mpi_sum(lsum,gsum,mpicom,string=trim(subname)//':sum',all=.false.) |
---|
1898 | call oasis_mpi_min(lmin,gmin,mpicom,string=trim(subname)//':min',all=.false.) |
---|
1899 | call oasis_mpi_max(lmax,gmax,mpicom,string=trim(subname)//':max',all=.false.) |
---|
1900 | endif |
---|
1901 | if (mype == 0) then |
---|
1902 | do m = 1,fsize |
---|
1903 | call mct_aVect_getRList(mstring,m,av) |
---|
1904 | itemc = mct_string_toChar(mstring) |
---|
1905 | call mct_string_clean(mstring) |
---|
1906 | write(nulprt,'(a,a16,3g21.12)') ' diags: ',trim(itemc),gmin(m),gmax(m),gsum(m) |
---|
1907 | enddo |
---|
1908 | endif |
---|
1909 | |
---|
1910 | deallocate(lsum,lmin,lmax) |
---|
1911 | deallocate(gsum,gmin,gmax) |
---|
1912 | deallocate(lwts) |
---|
1913 | |
---|
1914 | call oasis_debug_exit(subname) |
---|
1915 | |
---|
1916 | END SUBROUTINE oasis_advance_avdiag |
---|
1917 | |
---|
1918 | !------------------------------------------------------------------- |
---|
1919 | END MODULE mod_oasis_advance |
---|
1920 | |
---|