updated top-level README and version_decl for V4.4.2 (#1795)
[WRF.git] / external / atm_ocn / atm_comm.F
blob839b045f05d5d378cf33cea8ac65574ed7ef26dc
1       MODULE ATM_cc
3       USE CMP_COMM, ONLY: &
5      &   MPI_COMM_Atmos => COMM_local, &
7      &   Coupler_id, &
8      &   component_master_rank_local, &
9      &   process_rank_local, &
10      &   component_nprocs, &
11      &   ibuffer, &
13      &   MPI_INTEGER,MPI_STATUS_SIZE, &
14      &   kind_REAL,kind_alt_REAL, &
15      &   MPI_kind_REAL,MPI_kind_alt_REAL
17       implicit none
19       integer,parameter:: ND=3
20       integer Ocean_spec /-1/, WM_id /-10/
21       integer NSF
22       integer NSF_WM
23       real dtc,              &  !<- Coupling period
24      &     dta,              &  !<- AM time step ("physical")
25      &     dta2dtc              !<- AM time step / Coupling period
26       integer i_dtc2dta /100/   !<- Coupling period / AM time step
27       integer & !,dimension(ND)::
28      &ims,ime,jms,jme,its,ite,jts,jte,ids,idf,jds,jdf,  NGP
29       integer kms,kme,kts,kte,kds,kde
30       integer,parameter:: kind_R=kind_alt_REAL
31 !c     integer,parameter:: kind_tiling=kind_R
32       integer,parameter:: kind_sfcflux=kind_R, &
33      &                    kind_SST=kind_R, &
34      &                    kind_SLM=kind_R, &
35      &                    kind_lonlat=kind_REAL
36       INTEGER, PARAMETER  :: kind_cur = kind_r, kind_wstate = kind_r, &
37                              kind_windp = kind_r
38       integer MPI_kind_R, &
39      &MPI_kind_sfcflux,MPI_kind_SST,MPI_kind_SLM,MPI_kind_lonlat
40       INTEGER :: MPI_kind_cur, MPI_kind_wstate, MPI_kind_windp
41       integer n_ts(ND) /0,0,0/, gid
42       integer rc /5/
43       real,parameter:: &
44      &   SLM_OS_value=1., & !<-must be real open sea mask value in AM
45      &   unrealistically_low_SST=0.01, &  ! <- must be unreal low but >=0.,
46                                        ! see interp. --- check!
47      &   unrealistically_low_SV=-1.E30, &
48                      ! <- must be negative unreal low surface flux
49                      ! or other surface value to be sent
50                      ! to Coupler, see Coupler code
51      &   unrealistically_low_SF=unrealistically_low_SV, & !<- same thing
52      &   unrealistically_low_SVp=0.99*unrealistically_low_SV
54       logical initialized /.false./
55       logical PHYS,zeroSF,nrmSF,sendSF,getSST
57       TYPE SST_ARRAY
58         real(kind=kind_SST),dimension(:,:),pointer:: a 
59       END TYPE SST_ARRAY
60       TYPE SF_ARRAY
61         real(kind=kind_sfcflux),dimension(:,:,:),pointer:: a
62       END TYPE SF_ARRAY
64       TYPE cur_array
65         REAL(KIND = kind_cur), dimension(:,:), pointer :: a
66       END TYPE cur_array
67       TYPE wstate_array
68         REAL(KIND = kind_wstate), dimension(:,:), pointer :: a
69       END TYPE wstate_array
70       TYPE windp_array
71         REAL(KIND = kind_windp), dimension(:,:,:), pointer :: a
72       END TYPE windp_array
74       TYPE (SST_ARRAY), dimension(ND):: SST_cc
75       TYPE (SF_ARRAY), dimension(min(ND,2)):: sf
76       TYPE (cur_array), dimension(nd) ::  ucur_cc, vcur_cc
77       TYPE (wstate_array), dimension(nd) :: alpha_cc, gamma_cc
78       TYPE (windp_array), dimension(nd) :: wwinp
80       character*12 sgid
82 !Controls:
83       integer nunit_announce /6/, VerbLev /3/
84 ! To control awo couplings
85       integer ia2o /1/, &
86      &        io2a /1/, &
87      &        ia2w /1/, &
88      &        iw2a /0/
90       SAVE
92       END MODULE ATM_cc
94 !C***********************************************************************
96       SUBROUTINE ATM_SET_COMM(new_atm_comm)
97         USE ATM_cc
98         integer, intent(in) :: new_atm_comm
100 !       This routine is called when the atmospheric model wants to
101 !       remove processors from taking part in coupling so that they
102 !       can perform I/O or diagnostics.  Any processors that will
103 !       continue to be in MPI_COMM_Atmos must call this routine, and
104 !       any processors that will leave MPI_COMM_Atmos must call
105 !       ATM_LEAVE_COUPLING.
107         MPI_COMM_Atmos=new_atm_comm
108       end
110 !C***********************************************************************
112       SUBROUTINE ATM_LEAVE_COUPLING()
113         USE ATM_cc
115 !       This routine is called when the atmospheric model wants to
116 !       remove processors from taking part in coupling so that they
117 !       can perform I/O or diagnostics.  Any processors that will
118 !       continue to be in MPI_COMM_Atmos must call ATM_SET_COMM, and
119 !       any processors that will leave MPI_COMM_Atmos must call
120 !       this routine.
122 !       Currently, there is nothing we have to do here.
123         return
125       end
127 !C***********************************************************************
129       SUBROUTINE ATM_CMP_START(atm_comm)
131       USE ATM_cc
133       implicit none
135       integer atm_comm
137 !     integer Atmos_id /1/, Atmos_master_rank_local /0/, Atmos_spec /1/
138       integer Atmos_id /1/, Atmos_master_rank_local /0/, Atmos_spec /0/
139 !                                                            <- D.S.
140 !        Atmos_spec=1 for the case when the only field AM receives
141 !        from Coupler is SST. Atmos_spec=0 allows receiving additional
142 !        fields from C., originating from both OM, WM. (Atmos_spec does
143 !        not control receiving in AM but is sent to C. thus transferring
144 !        the control to C.)
145       integer ibuf(1),ierr
146       character*20 s
149                       !<-id of OM as a component of the coupled system
150       call CMP_INIT(Atmos_id,1)
151                              !<-"flexibility level"
152       if (Coupler_id.ge.0) VerbLev=min(VerbLev,ibuffer(4))
153       write(s,'(i2)') VerbLev
155       call CMP_INTRO(Atmos_master_rank_local)
156       call ATM_ANNOUNCE('back from CMP_INTRO, VerbLev='//s,2)
158       initialized=.true.
160       call CMP_INTEGER_SEND(Atmos_spec,1)
162       call CMP_gnr_RECV(Ocean_spec,1,MPI_INTEGER)
163       write(s,'(i2)') Ocean_spec
164       call ATM_ANNOUNCE('back from CMP_INTEGER_RECV, OM spec is '//s,2)
165       call MPI_BCAST(Ocean_spec,1,MPI_INTEGER, &
166      &component_master_rank_local,MPI_COMM_Atmos,ierr)
167       call ATM_ANNOUNCE('ATM_CMP_START: Ocean_spec broadcast',2)
169       call CMP_gnr_RECV(WM_id,1,MPI_INTEGER)
170       write(s,'(i4)') WM_id
171       call ATM_ANNOUNCE('back from CMP_INTEGER_RECV, WM id is '//s,2)
172       call MPI_BCAST(WM_id,1,MPI_INTEGER, &
173      &component_master_rank_local,MPI_COMM_Atmos,ierr)
174       call ATM_ANNOUNCE('ATM_CMP_START: WM_id broadcast',2)
175       if (WM_id.gt.0) then
176         NSF_WM=4
177       else
178         NSF_WM=0
179       end if
181       if (Ocean_spec.eq.1) then
182         NSF=4
183       else if (Ocean_spec.eq.2) then
184         NSF=8
185       else if (Ocean_spec.eq.0) then
186         NSF=1
187       else if (Coupler_id.ge.0) then
188         call GLOB_ABORT(Ocean_spec-1, &
189      &  'ATM_CMP_START received wrong Ocean_spec value, aborted',rc)
190       else
191         Ocean_spec=1
192         NSF=4
193         call ATM_ANNOUNCE('AM is standalone: Ocean_spec=1, NSF=4'// &
194      &  ' assigned (as if for POM coupling)',2)
195       end if
197       call CMP_gnr_RECV(ia2o,1,MPI_INTEGER)
198       write(s,'(i4)') ia2o
199       call ATM_ANNOUNCE('back from CMP_INTEGER_RECV, ia2o is '//s,2)
200       call MPI_BCAST(ia2o,1,MPI_INTEGER, &
201      &component_master_rank_local,MPI_COMM_Atmos,ierr)
202       call ATM_ANNOUNCE('ATM_CMP_START: ia2o broadcast',2)
204       call CMP_gnr_RECV(io2a,1,MPI_INTEGER)
205       write(s,'(i4)') io2a
206       call ATM_ANNOUNCE('back from CMP_INTEGER_RECV, io2a is '//s,2)
207       call MPI_BCAST(io2a,1,MPI_INTEGER, &
208      &component_master_rank_local,MPI_COMM_Atmos,ierr)
209       call ATM_ANNOUNCE('ATM_CMP_START: io2a broadcast',2)
211       call CMP_gnr_RECV(ia2w,1,MPI_INTEGER)
212       write(s,'(i4)') ia2w
213       call ATM_ANNOUNCE('back from CMP_INTEGER_RECV, ia2w is '//s,2)
214       call MPI_BCAST(ia2w,1,MPI_INTEGER, &
215      &component_master_rank_local,MPI_COMM_Atmos,ierr)
216       call ATM_ANNOUNCE('ATM_CMP_START: ia2w broadcast',2)
218       call CMP_gnr_RECV(iw2a,1,MPI_INTEGER)
219       write(s,'(i4)') iw2a
220       call ATM_ANNOUNCE('back from CMP_INTEGER_RECV, iw2a is '//s,2)
221       call MPI_BCAST(iw2a,1,MPI_INTEGER, &
222      &component_master_rank_local,MPI_COMM_Atmos,ierr)
223       call ATM_ANNOUNCE('ATM_CMP_START: iw2a broadcast',2)
225       if (kind_R.eq.kind_REAL) then
226         MPI_kind_R=MPI_kind_REAL
227       else 
228         MPI_kind_R=MPI_kind_alt_REAL
229       end if
230       if (kind_sfcflux.eq.kind_REAL) then
231         MPI_kind_sfcflux=MPI_kind_REAL
232       else 
233         MPI_kind_sfcflux=MPI_kind_alt_REAL
234       end if
235       if (kind_SST.eq.kind_REAL) then
236         MPI_kind_SST=MPI_kind_REAL
237       else 
238         MPI_kind_SST=MPI_kind_alt_REAL
239       end if
240       if (kind_SLM.eq.kind_REAL) then
241         MPI_kind_SLM=MPI_kind_REAL
242       else 
243         MPI_kind_SLM=MPI_kind_alt_REAL
244       end if
245       if (kind_lonlat.eq.kind_REAL) then
246         MPI_kind_lonlat=MPI_kind_REAL
247       else 
248         MPI_kind_lonlat=MPI_kind_alt_REAL
249       end if
251       IF (kind_cur == kind_real) THEN
252          MPI_kind_cur = MPI_kind_real
253       ELSE
254          MPI_kind_cur = MPI_kind_alt_real
255       END IF
256       IF (kind_wstate == kind_real) THEN
257          MPI_kind_wstate = MPI_kind_real
258       ELSE
259          MPI_kind_wstate = MPI_kind_alt_real
260       END IF
261       IF (kind_windp == kind_real) THEN
262          MPI_kind_windp = MPI_kind_real
263       ELSE
264          MPI_kind_windp = MPI_kind_alt_real
265       END IF
267       atm_comm=MPI_COMM_Atmos
269       return
270       END
272 !C***********************************************************************
274       SUBROUTINE ATM_INIT_CHECK(s)
276       USE ATM_cc, ONLY: initialized,rc
278       implicit none
280       character*(*) s
282       if (.not. initialized) call GLOB_ABORT(1,s,rc)
284       return
285       END
287 !C***********************************************************************
289       subroutine ATM_TSTEP_INIT(NTSD,NPHS,gid_,dta_, &
290      &ids_,idf_,jds_,jdf_,its_,ite_,jts_,jte_,ims_,ime_,jms_,jme_, &
291        !<-"domain"         !<-"tile"           !<-"memory" (tile+halo)
292      &kds_,kde_,kts_,kte_,kms_,kme_,&
293      &HLON,HLAT,VLON,VLAT, &
294      &SLM, &
295      &i_parent_start,j_parent_start,&
296      &guessdtc,dtc_)
298       USE ATM_cc
300       implicit none
302       real, intent(in) :: guessdtc
303       real, intent(out) :: dtc_
305       integer NTSD,NPHS,gid_
306       real dta_
307       integer ids_,idf_,jds_,jdf_,its_,ite_,jts_,jte_, &
308      &ims_,ime_,jms_,jme_,kds_,kde_,kts_,kte_,kms_,kme_
309       real(kind=kind_lonlat),dimension(ims_:ime_,jms_:jme_):: &
310      &HLON,HLAT,VLON,VLAT
311       real(kind=kind_SLM),dimension(ims_:ime_,jms_:jme_):: SLM
312       integer i_parent_start,j_parent_start
314       integer KDT,buf(2) /0,0/
315       character*24 s
316       character*80 s1
317       character*255 message
318       SAVE
321       gid=gid_
322       call GLOB_ABORT((gid-1)*(gid-2)*(gid-3), &
323      &'Abort: in ATM_TSTEP_INIT gid is neither 1 nor 2 nor 3',rc)
324       KDT=NTSD/NPHS+1
325       PHYS=MOD(NTSD,NPHS).eq.0 ! .and. gid.eq.1 <-removed to bring MG in
326       dta=dta_ 
328       write(s1,'("gid=",i1," NTSD=",i5," NPHS=",i3," KDT=",i5,'// &
329      &'" PHYS=",L1)') gid,NTSD,NPHS,KDT,PHYS
330       call ATM_ANNOUNCE('ATM_TSTEP_INIT entered: '//trim(s1),3)
332 !c     IF (n_ts.eq.-1 .and. PHYS) THEN
333 !c       PHYS=.false.
334 !c       n_ts=0   ! init. value must be -1 . But if PHYS does not need
335 !c                ! this correction, init. value must be 0 (whereas this
336 !c                ! IF statement may stay)
337 !c     END IF
338       if (.not.PHYS) then
339         zeroSF=.false.
340         nrmSF=.false.
341         sendSF=.false.
342         RETURN
343       end if
345       n_ts(gid)=n_ts(gid)+1  ! init. value must be 0   ***0***
346       write(s,'(2i8)') KDT,n_ts(gid)
347       write(sgid,'(" grid id = ",i1)') gid
348       call ATM_ANNOUNCE('ATM_TSTEP_INIT working:'// &
349      &sgid//'; KDT, n_ts: '//s,3)
350       call GLOB_ABORT(KDT-n_ts(gid), &
351      &'Abort: in ATM_TSTEP_INIT KDT, n_ts(gid) differ '//s,rc)
353       call ATM_RECVdtc(guessdtc)
354       dtc_=dtc
356       zeroSF=((n_ts(gid)-1)/i_dtc2dta)*i_dtc2dta .eq. n_ts(gid)-1
357       nrmSF=(n_ts(gid)/i_dtc2dta)*i_dtc2dta .eq. n_ts(gid)
358       sendSF=(n_ts(gid)/i_dtc2dta)*i_dtc2dta .eq. n_ts(gid)
359                                     !<-check, this depends
360                                     ! on where ATM_SENDFLUXES is called.
361                                     ! MOD(n_ts,i_dtc2dta).eq.0 should
362                                     ! be good for calling it after
363                                     ! ATM_DOFLUXES at the same t.s.
365       ids=ids_
366       idf=idf_
367       jds=jds_
368       jdf=jdf_
369       its=its_
370       ite=ite_
371       jts=jts_
372       jte=jte_
373       ims=ims_
374       ime=ime_
375       jms=jms_
376       jme=jme_
378       kds=kds_
379       kde=kde_
380       kts=kts_
381       kms=kms_
382       kme=kme_
383       kte=kte_
385       NGP=(idf-ids+1)*(jdf-jds+1)
387       IF (n_ts(gid).eq.1) THEN
389       call ATM_ANNOUNCE('ATM_TSTEP_INIT to allocate sf, SST_cc',3)
391       IF (gid.le.2) THEN !** innermost grid not active in coupling **
392         allocate(sf(gid)%a(ims:ime,jms:jme,NSF))
393         ALLOCATE(wwinp(gid)%a(ims:ime,jms:jme,NSF_WM))
394       END IF !** innermost grid not active in coupling **
395         allocate(SST_cc(gid)%a(ims:ime,jms:jme))
396       ALLOCATE(ucur_cc(gid)%a(ims:ime,jms:jme), vcur_cc(gid)%a(ims:ime,jms:jme))
397       ALLOCATE(alpha_cc(gid)%a(ims:ime,jms:jme), gamma_cc(gid)%a(ims:ime,jms:jme))
399       END IF
401       if (gid.eq.2) then
402         write(s,'(2i8)') i_parent_start,j_parent_start
403         if (zeroSF) then
404           buf(1)=i_parent_start
405           buf(2)=j_parent_start
406           call CMP_INTEGER_SEND(buf,2)
407           call ATM_ANNOUNCE( &
408      &    'ATM_TSTEP_INIT: i_parent_start, j_parent_start sent '//s,3)
409         else
410           call GLOB_ABORT(abs(i_parent_start-buf(1))+abs(j_parent_start- &
411      &    buf(2)),'NESTED GRID MOVED DURING C TIME STEP: ABORTED '// &
412      &    s,rc)
413         end if
414       end if
416       IF (gid.le.2) THEN !** innermost grid not active in coupling **
418       CALL ATM_SENDGRIDS(HLON,HLAT,VLON,VLAT)
420       CALL ATM_SENDSLM(SLM)
422       END IF !** innermost grid not active in coupling **
424       if (VerbLev.ge.2) then
425          write(message,*) 'AM: ATM_TSTEP_INIT: returning ',gid,         &
426      &n_ts(gid),ids,idf,jds,jdf,its,ite,jts,jte,ims,ime,jms,jme,NGP,NSF
427          call wrf_debug(2,message)
428       endif
430       RETURN
431       end
433 !C***********************************************************************
435       SUBROUTINE ATM_RECVdtc(guessdtc)
437       USE ATM_cc
439       implicit none
441       real,intent(in) :: guessdtc
442       real(kind=kind_R) buf(1),dtc2dta
443       integer ierr,i
444       logical first/.true./
445       character*20 s
446       SAVE
449       write(s,'(1pe20.12)') dta
450       call ATM_ANNOUNCE('ATM_RECVdtc: AM time step dta='//s,3)
452       IF (first) THEN
453         call ATM_ANNOUNCE( &
454      &  'ATM_RECVdtc: to receive C time step; AM time step dta='//s,2)
456         call CMP_gnr_RECV(buf,1,MPI_kind_R)
458         call MPI_BCAST(buf,1,MPI_kind_R, &
459      &  component_master_rank_local,MPI_COMM_Atmos,ierr)
460         call ATM_ANNOUNCE('ATM_RECVdtc: C time step broadcast',2)
461         dtc=buf(1)
463         if (Coupler_id.lt.0) then
464           dtc=guessdtc
465           write(s,'(1pe20.12)') dtc
466           call ATM_ANNOUNCE('ATM_RECVdtc: C time step assigned '// &
467      &    trim(s)//' : standalone mode',2)
468         else
469           write(s,'(1pe20.12)') buf
470           call ATM_ANNOUNCE( &
471      &    'ATM_RECVdtc: C time step dtc='//s//' received',2)
472         end if
473       END IF
475       dtc2dta=dtc/dta
476       i_dtc2dta=nint(dtc2dta)
477       if (abs(i_dtc2dta-dtc2dta).gt.1.E-5) call GLOB_ABORT(1, &
478      &'AM: ABORTED: dtc is not a multiple of dta',1)
480       i=3
481       if (n_ts(gid).eq.1) i=2
482       if (i_dtc2dta.eq.0) then
483         i_dtc2dta=4
484         call ATM_ANNOUNCE('ratio of C/AM time steps =0, assigned 4 .'// &
485      &  ' This should only occur in standalone mode and ONLY IF dtc '// &
486      &  'HAS NOT BEEN ASSIGNED A POSITIVE VALUE: ** ATTENTION **',i)
487       else
488         write(s,'(i2)') i_dtc2dta
489         call ATM_ANNOUNCE('ratio of C/AM time steps: '//trim(s),i)
490       end if
492       dta2dtc=1./i_dtc2dta
494       first=.false.
496       RETURN
497       END
499 !C***********************************************************************
501       SUBROUTINE ATM_SENDGRIDS(HLON,HLAT,VLON,VLAT)
503       USE ATM_cc
505       implicit none
507       real(kind=kind_lonlat),dimension(ims:ime,jms:jme):: &
508      &HLON,HLAT,VLON,VLAT 
510       real(kind=kind_lonlat),dimension(ids:idf,jds:jdf):: &
511      &ALONt,ALATt,ALONv,ALATv
513       integer buf(2)
516 !c     IF (gid.ne.1) RETURN ! only "parent grid" dim. and coor. are sent
518       IF (.not.PHYS .or. n_ts(gid).ne.1) RETURN
520       IF (gid.gt.2) RETURN ! innermost grid's dim. / coor. are not sent
521       
522 !temporarily excluded      if (Coupler_id.lt.0) return    !   <- standalone mode
524       buf(1)=idf-ids+1
525       buf(2)=jdf-jds+1
526       call ATM_ANNOUNCE('to send grid dimensions,'//sgid,1)
527       call CMP_INTEGER_SEND(buf,2)
528       call ATM_ANNOUNCE('grid dimensions sent,'//sgid,1)
530 !c     IF (gid.eq.1) THEN    !  only "parent grid" coordinates are sent
532         if (kind_lonlat.eq.4) then
533           call ASSEMBLE(ALONt,HLON,kind_lonlat)
534           call ASSEMBLE(ALATt,HLAT,kind_lonlat)
535           call ASSEMBLE(ALONv,VLON,kind_lonlat)
536           call ASSEMBLE(ALATv,VLAT,kind_lonlat)
537         else if (kind_lonlat.eq.8) then
538           call ASSEMBLE_R8(ALONt,HLON,kind_lonlat)
539           call ASSEMBLE_R8(ALATt,HLAT,kind_lonlat)
540           call ASSEMBLE_R8(ALONv,VLON,kind_lonlat)
541           call ASSEMBLE_R8(ALATv,VLAT,kind_lonlat)
542         else
543           call GLOB_ABORT(1,'wrong value of kind_lonlat in ATM_SENDGRIDS',1)
544         end if
546         call ATM_ANNOUNCE('(BP) to send grid arrays (4 MPI calls)',2)
548         call CMP_gnr_SEND(ALONt,NGP,MPI_kind_lonlat)
549         call CMP_gnr_SEND(ALATt,NGP,MPI_kind_lonlat)
550         call CMP_gnr_SEND(ALONv,NGP,MPI_kind_lonlat)
551         call CMP_gnr_SEND(ALATv,NGP,MPI_kind_lonlat)
553         call ATM_ANNOUNCE('the 4 grid arrays sent',1)
555 !c     END IF
557       call ATM_ANNOUNCE('(BP) ATM_SENDGRIDS: returning',2)
559       return
560       END
562 !C***********************************************************************
564       SUBROUTINE ATM_SENDSLM(SLM)
566       USE ATM_cc
568       implicit none
570       real(kind=kind_SLM),dimension(ims:ime,jms:jme):: SLM
572       real(kind=kind_SLM),dimension(ids:idf,jds:jdf):: SLM_g
573       integer buf(2)
576 !c     IF (gid.ne.1) RETURN  !  only "parent grid" mask is sent
577       IF (.not.PHYS .or. n_ts(gid).ne.1) RETURN
579       IF (gid.gt.2) RETURN ! innermost grid's mask is not sent
580       
581 !temporarily excluded      if (Coupler_id.lt.0) return    !   <- standalone mode
583       call ASSEMBLE(SLM_g,SLM,kind_SLM)
585       call ATM_ANNOUNCE('(BP) to send SLM',2)
587       call CMP_gnr_SEND(SLM_g,NGP,MPI_kind_SLM)
588       call CMP_gnr_SEND(SLM_g,NGP,MPI_kind_SLM)
589            ! Coupler requires as many copies of mask as there are grids
590            ! [and mask array is the same for H- (=t-) and V- grids]
592       call ATM_ANNOUNCE('(BP) ATM_SENDSLM: returning',2)
594       return
595       END
597 !C***********************************************************************
599       SUBROUTINE ATM_GETSST(SST,SLM)
601       USE ATM_cc
603       implicit none
605       real(kind=kind_SST) SST(ims:ime,jms:jme)
606       real(kind=kind_SLM) SLM(ims:ime,jms:jme)
608       integer i,j
609       real(kind=kind_SST) SST_g(ids:idf,jds:jdf)
611       IF ( io2a .LT. 1 ) RETURN
613       IF (.not.PHYS) RETURN
615       IF (gid.gt.2) RETURN ! nothing is done to get innermost grid's
616                            ! SST ** IN THIS PRELIMINARY VERSION **
618       call ATM_ANNOUNCE('ATM_GETSST entered (PHYS=.true.)',3)
620       getSST=((n_ts(gid)-1)/i_dtc2dta)*i_dtc2dta .eq. n_ts(gid)-1
621       if (getSST.neqv.zeroSF) then
622         call GLOB_ABORT(1,'getSST differs from zeroSF, which screws'// &
623      &  ' up the design for exchanges with C',rc)
624       end if
626       if (getSST) then
627         if (n_ts(gid).eq.1 .and. gid.eq.1) then
628           call ATM_ANNOUNCE('ATM_GETSST: to send ref. SST'//sgid,2)
629           call ASSEMBLE(SST_g,SST,kind_SST)
630           call CMP_gnr_SEND(SST_g,NGP,MPI_kind_SST)
631           call ATM_ANNOUNCE('ATM_GETSST: ref. SST sent'//sgid,2)
632         end if
633         call ATM_ANNOUNCE('ATM_GETSST: to receive SST',3)
634         call CMP_gnr_RECV(SST_g,NGP,MPI_kind_SST)
635         call DISASSEMBLE(SST_g,SST_cc(gid)%a,kind_SST)
636         call ATM_ANNOUNCE('ATM_GETSST: SST received',3)
637       end if
638       
639       if (Coupler_id.lt.0) return    !   <- standalone mode
641       do j=jts,jte
642       do i=its,ite
643         if (abs(SLM(i,j)-SLM_OS_value).lt.0.01) then 
644                                   ! i.e. if it is OS (open sea) AMGP
645                                   !
646           if (SST_cc(gid)%a(i,j).gt.unrealistically_low_SST) &
647                                           ! i.e. if there is a valid
648                                           ! result of interpolation from
649                                           ! OMG for this AMGP
650      &       SST(i,j)=SST_cc(gid)%a(i,j)
651         end if
652       end do
653       end do
655       return
656       END
658 !C***********************************************************************
660       SUBROUTINE ATM_DOFLUXES(TWBS,QWBS,RLWIN,RSWIN,RADOT,RSWOUT, &
661 !c    &USTAR,U10,V10,PINT,PREC)
662      &TX,TY,PINT,PREC)
664       USE ATM_cc
666       implicit none
668       real(kind=kind_sfcflux),dimension(ims:ime,jms:jme,kms:kme):: PINT
670       real(kind=kind_sfcflux),dimension(ims:ime,jms:jme):: &
671      &TWBS,QWBS,RLWIN,RSWIN,RADOT,RSWOUT,TX,TY,PREC
672 !c    &TWBS,QWBS,RLWIN,RSWIN,RADOT,RSWOUT,USTAR,U10,V10,PINT,PREC
673 !       Act. arg. for PINT is a 3d array - so this only is OK if
674 !       Ps=Act.arg.(:,:.1) - actually, Ps=PINT(:,1,:)
676       real(kind=kind_sfcflux),dimension(ims:ime,jms:jme):: SWR,R
677       real dtainv
680       IF ( ia2o .LT. 1 .and. ia2w .LT. 1 ) RETURN
681       IF (.not.PHYS) RETURN
683       IF (gid.gt.2) RETURN
685 ! Debug insertion:->
686 !c     if (PREC(ims+3,jms+3).ne.0 .or. PREC(ims+5,jms+5).ne.0) then
687 !c       print*,'ATM_DOFLUXES,gid,n_ts(gid),PREC(3,3),PREC(5,5): ',
688 !c    &  gid,n_ts(gid),PREC(ims+3,jms+3),PREC(ims+5,jms+5)
689 !c     end if
690 ! <-:Debug insertion
692       call ATM_ANNOUNCE('ATM_DOFLUXES entered',3)
694       dtainv=1./dta
696       if (zeroSF) sf(gid)%a=0.
698       SWR(its:ite,jts:jte)=-RSWIN(its:ite,jts:jte)+RSWOUT(its:ite,jts:jte)          ! Check sign! here SWR is meant to be
699                                  ! positive upward
700 !c     sf(gid)%a(:,:,NSF-1)=sf(gid)%a(:,:,NSF-1)-TX
701 !c     sf(gid)%a(:,:,NSF)=sf(gid)%a(:,:,NSF)-TY
702 !c                    ! <- signs for stress components are     changed
703 !c                    ! so it is -stress
705 !c     R=SWR+RADOT-RLWIN          ! Check sign! here R (net radiation)
706                                  ! is meant to be positive upward
708 !oooooooooooooooooooooooooooooo
709       IF (Ocean_spec.eq.1) THEN
710 !oooooooooooooooooooooooooooooo
711         sf(gid)%a(its:ite,jts:jte,NSF-3)=sf(gid)%a(its:ite,jts:jte,NSF-3)-TWBS(its:ite,jts:jte)-QWBS(its:ite,jts:jte)+RADOT(its:ite,jts:jte)-RLWIN(its:ite,jts:jte)
712                                        ! -TWBS (-QWBS) is supposed to
713                                        ! be sensible (latent) heat flux,
714                                        ! positive upward
715         sf(gid)%a(its:ite,jts:jte,NSF-2)=sf(gid)%a(its:ite,jts:jte,NSF-2)+SWR(its:ite,jts:jte)
716         sf(gid)%a(its:ite,jts:jte,NSF-1)=sf(gid)%a(its:ite,jts:jte,NSF-1)-TX(its:ite,jts:jte)
717         sf(gid)%a(its:ite,jts:jte,NSF)=sf(gid)%a(its:ite,jts:jte,NSF)-TY(its:ite,jts:jte)
718                      ! <- signs for stress components are changed
719 !ooooooooooooooooooooooooooooooooooo
720       ELSE IF (Ocean_spec.eq.2) THEN
721 !ooooooooooooooooooooooooooooooooooo
722         sf(gid)%a(its:ite,jts:jte,NSF-7)=sf(gid)%a(its:ite,jts:jte,NSF-7)+PREC(its:ite,jts:jte)
723         sf(gid)%a(its:ite,jts:jte,NSF-6)=sf(gid)%a(its:ite,jts:jte,NSF-6)-TWBS(its:ite,jts:jte)
724         sf(gid)%a(its:ite,jts:jte,NSF-5)=sf(gid)%a(its:ite,jts:jte,NSF-5)-QWBS(its:ite,jts:jte)
725         sf(gid)%a(its:ite,jts:jte,NSF-4)=sf(gid)%a(its:ite,jts:jte,NSF-4)+PINT(its:ite,jts:jte,1)-101300.
726         sf(gid)%a(its:ite,jts:jte,NSF-3)=sf(gid)%a(its:ite,jts:jte,NSF-3)-SWR(its:ite,jts:jte)-RADOT(its:ite,jts:jte)+RLWIN(its:ite,jts:jte)
727         sf(gid)%a(its:ite,jts:jte,NSF-2)=sf(gid)%a(its:ite,jts:jte,NSF-2)-SWR(its:ite,jts:jte)
729         sf(gid)%a(its:ite,jts:jte,NSF-1)=sf(gid)%a(its:ite,jts:jte,NSF-1)+TX(its:ite,jts:jte)
730         sf(gid)%a(its:ite,jts:jte,NSF)=sf(gid)%a(its:ite,jts:jte,NSF)+TY(its:ite,jts:jte)
731                      ! <- signs for stress components are NOT changed
732         if (nrmSF) then
733           sf(gid)%a(its:ite,jts:jte,1)=sf(gid)%a(its:ite,jts:jte,1)*dtainv
734                         ! so this will be m/s; check what OM wants
735         end if
736 !ooooooooooo
737       END IF
738 !ooooooooooo
740       if (nrmSF) then
741         sf(gid)%a=sf(gid)%a*dta2dtc
742       end if
744       call ATM_ANNOUNCE('ATM_DOFLUXES to return',3)
746       return
747       END
749 !C***********************************************************************
751       SUBROUTINE ATM_SENDFLUXES
753       USE ATM_cc
755       implicit none
757       real(kind=kind_sfcflux) F(ids:idf,jds:jdf)
758       integer n
761       IF ( ia2o .LT. 1 .and. ia2w .LT. 1 ) RETURN
762       if (.not.PHYS) RETURN
764       IF (gid.gt.2) RETURN
766       if (.not.sendSF) then
767         call ATM_ANNOUNCE( &
768      &  'ATM_SENDLUXES entered with PHYS but not sendSF: returning'// &
769      &  sgid,3)
770         RETURN
771       end if
773       call ATM_ANNOUNCE('In ATM_SENDLUXES'//sgid,3)
775       do n=1,NSF
776         call ASSEMBLE(F,sf(gid)%a(:,:,n),kind_sfcflux)
777         call CMP_gnr_SEND(F,NGP,MPI_kind_sfcflux)
778       end do
780       call ATM_ANNOUNCE('ATM_SENDFLUXES to return'//sgid,3)
782       return
783       END
785 !C***********************************************************************
787       SUBROUTINE ATM_ANNOUNCE(s,DbgLev)
789       USE ATM_cc, ONLY: nunit_announce,VerbLev,MPI_COMM_Atmos
791       implicit none
793       character*(*) s
794       integer DbgLev
796       integer ierr
798       if (DbgLev.le.VerbLev) then
799         if (s(1:5).eq.'(BP) ') then
800           call MPI_BARRIER(MPI_COMM_Atmos,ierr)
801         end if
802         CALL CMP_ANNOUNCE(nunit_announce,'AM: '//s)
803       end if
805       return
806       END
808 SUBROUTINE atm_getcur(ucur,vcur)
810 ! Bringing ocean currents .. 
811 ! Biju Thomas,  GSO/URI  on 4/8/2015
813    USE atm_cc
814    IMPLICIT NONE
815    REAL(KIND = kind_cur), DIMENSION(ims:ime,jms:jme) :: ucur, vcur
816    REAL(KIND = kind_cur), DIMENSION(ids:idf,jds:jdf) :: ucur_g, vcur_g
817    REAL, PARAMETER :: cur_ll = 0._kind_cur, cur_ul = 5._kind_cur,  &
818                       cur_k = 0._kind_cur
819                       
820    INTEGER :: i, j
821    LOGICAL :: getcur
823       IF ( io2a .LT. 2 ) RETURN
824    IF (.NOT. phys .OR. gid > 2) RETURN
825    IF (ocean_spec /= 1) CALL atm_announce('Warn: ocean currents needed',3)
826    CALL atm_announce('atm_getcur entered (phys = .true.)',3)
827    
828    getcur = ((n_ts(gid)-1)/i_dtc2dta)*i_dtc2dta == n_ts(gid)-1
829    
830    IF (getcur .NEQV. zerosf) THEN
831       CALL glob_abort(1, 'Warn: getcur does not match zerosf', rc)
832    END IF
834    IF (getcur) THEN
835       CALL atm_announce('atm_getcur: to receive CUR',3)
836       CALL cmp_gnr_recv(ucur_g, ngp, mpi_kind_cur)
837       CALL cmp_gnr_recv(vcur_g, ngp, mpi_kind_cur)
838       CALL disassemble(ucur_g, ucur_cc(gid)%a, kind_cur)
839       CALL disassemble(vcur_g, vcur_cc(gid)%a, kind_cur)
840       CALL atm_announce('atm_getcur: CUR received',3)
841    END IF  
842    
843    IF ( coupler_id .LT. 0 ) RETURN
846    DO j = jts,jte
847    DO i = its,ite
848      IF ( ABS(ucur_cc(gid)%a(i,j)) .GE. cur_ll .AND.                &
849           ABS(ucur_cc(gid)%a(i,j)) .LE. cur_ul ) THEN
850         ucur(i,j) = ucur_cc(gid)%a(i,j)
851      ELSE
852         ucur(i,j) = cur_k
853      ENDIF
854      IF ( ABS(vcur_cc(gid)%a(i,j)) .GE. cur_ll .AND.                &
855           ABS(vcur_cc(gid)%a(i,j)) .LE. cur_ul ) THEN
856         vcur(i,j) = vcur_cc(gid)%a(i,j)
857      ELSE
858         vcur(i,j) = cur_k
859      ENDIF
860    ENDDO
861    ENDDO
863 END SUBROUTINE atm_getcur 
865 SUBROUTINE atm_getwstate(alpha,gamma)
867 ! Bringing Wave state (Charnok coeff & misalignment Angle) 
868 ! Biju Thomas,  GSO/URI  on 4/8/2015
870    USE atm_cc
871    IMPLICIT NONE
872    REAL(KIND = kind_wstate), DIMENSION(ims:ime,jms:jme) :: alpha, gamma
873    REAL(KIND = kind_wstate), DIMENSION(ids:idf,jds:jdf) :: alpha_g, &
874                                                            gamma_g
875    REAL, PARAMETER :: deg2rad=3.1415926_kind_wstate/180_kind_wstate 
876    REAL, PARAMETER :: alpha_ll = 0.0_kind_wstate, &
877                       alpha_ul = 0.2_kind_wstate, &
878                       alpha_k = 0.0185_kind_wstate, &
879                       gamma_ll = -20.0_kind_wstate*deg2rad, &
880                       gamma_ul = 20.0_kind_wstate*deg2rad, &
881                       gamma_k = 0.0_kind_wstate
884    INTEGER :: i, j
885    LOGICAL :: getwstate
887       IF ( iw2a .LT. 1 ) RETURN
888    IF (wm_id <= 0)  RETURN
889    IF (.NOT. phys .OR. gid > 2) RETURN
890    CALL atm_announce('atm_getstate entered (phys = .true.)',3)
892    getwstate = ((n_ts(gid)-1)/i_dtc2dta)*i_dtc2dta == n_ts(gid)-1
894    IF (getwstate .NEQV. zerosf) THEN
895       CALL glob_abort(1, 'Warn: getwstate does not match zerosf', rc)
896    END IF
898    IF (getwstate) THEN
899       CALL atm_announce('atm_getwsate: to receive WSTATE',3)
900       CALL cmp_gnr_recv(alpha_g, ngp, mpi_kind_wstate)
901       CALL cmp_gnr_recv(gamma_g, ngp, mpi_kind_wstate)
902       CALL disassemble(alpha_g, alpha_cc(gid)%a, kind_wstate)
903       CALL disassemble(gamma_g, gamma_cc(gid)%a, kind_wstate)
904       CALL atm_announce('atm_getwstate: WSTATE received',3)
905    END IF
907    IF ( coupler_id .LT. 0 ) RETURN
908    
909    
910    DO j = jts,jte
911    DO i = its,ite
912 !     IF ( alpha_cc(gid)%a(i,j) .GT. alpha_ll .AND.                &
913 !          alpha_cc(gid)%a(i,j) .LT. alpha_ul ) THEN
914 !        alpha(i,j) = alpha_cc(gid)%a(i,j)
915 !     ELSE
916 !        alpha(i,j) = alpha_k
917 !     ENDIF 
918      alpha(i,j) = alpha_cc(gid)%a(i,j)
919      IF ( gamma_cc(gid)%a(i,j) .GT. gamma_ll .AND.                &
920           gamma_cc(gid)%a(i,j) .LT. gamma_ul ) THEN
921         gamma(i,j) =  gamma_cc(gid)%a(i,j)
922      ELSE
923         gamma(i,j) = gamma_k
924      ENDIF
925    ENDDO
926    ENDDO
928 END SUBROUTINE atm_getwstate
930 SUBROUTINE atm_prepwindp(ulowl, vlowl, richn, zlowl)
932 ! Preparing Wind and adjusting variables (Height of lowest model level 
933 !                                          Richarson number)
934 ! Biju Thomas,  GSO/URI  on 4/8/2015
936    USE atm_cc
937    IMPLICIT NONE
939    REAL(KIND = kind_windp), DIMENSION(ims:ime,jms:jme) :: ulowl, vlowl, &
940                                                           richn, zlowl
941       IF ( ia2w .LT. 1 ) RETURN
942    IF (wm_id <= 0)  RETURN
943    IF (.NOT. phys .OR. gid > 2) RETURN
945    CALL atm_announce('atm_atm_prepwindp: entered',3)
947    IF (zerosf) wwinp(gid)%a = 0.0
948    
949    wwinp(gid)%a(its:ite,jts:jte,NSF_WM-1) =    &  ! D.S.
950                               wwinp(gid)%a(its:ite,jts:jte,NSF_WM-1) + &
951                               ulowl(its:ite,jts:jte)
952    wwinp(gid)%a(its:ite,jts:jte,NSF_WM) =      &  ! D.S.
953                               wwinp(gid)%a(its:ite,jts:jte,NSF_WM)   + &
954                               vlowl(its:ite,jts:jte) 
955    wwinp(gid)%a(its:ite,jts:jte,NSF_WM-3) =    &  ! D.S.
956                               wwinp(gid)%a(its:ite,jts:jte,NSF_WM-3) + &
957                               richn(its:ite,jts:jte)
958    wwinp(gid)%a(its:ite,jts:jte,NSF_WM-2) =    &  ! D.S.
959                               wwinp(gid)%a(its:ite,jts:jte,NSF_WM-2) + &
960                               zlowl(its:ite,jts:jte)
961    IF (nrmsf) THEN
962       wwinp(gid)%a = wwinp(gid)%a*dta2dtc
963    END IF  
965    CALL atm_announce('atm_atm_prepwindp: returned',3)
967 END SUBROUTINE atm_prepwindp
970 SUBROUTINE atm_sendwindp
972 ! Sending wind and it adjustment fields (U1, V1, Charnok coeff & misalignment Angle)
973 ! Biju Thomas,  GSO/URI  on 4/8/2015
975    USE atm_cc
976    IMPLICIT NONE
978    INTEGER :: n
979    REAL(KIND = kind_windp), DIMENSION(ids:idf,jds:jdf) :: field
981       IF ( ia2w .LT. 1 ) RETURN
982    IF (wm_id <= 0)  RETURN
983    IF (.NOT. phys .OR. gid > 2) RETURN
985    IF (.NOT. sendsf) THEN
986       CALL atm_announce('atm_sendwindp entered with PHYS but not sendSF: returning'// &
987       sgid,3)
988       RETURN
989    END IF
991    CALL atm_announce('atm_prepwindp: entered'//sgid,3)
993    DO n = 1, NSF_WM
994       CALL assemble(field, wwinp(gid)%a(:,:,n), kind_windp)
995       CALL cmp_gnr_send(field, ngp, mpi_kind_windp) 
996    END DO
998    CALL atm_announce('atm_prepwindp: reterned'//sgid,3)
1000 END SUBROUTINE atm_sendwindp