ROMS
Loading...
Searching...
No Matches
esmf_regcm_mod Module Reference

Functions/Subroutines

subroutine, public atm_setservices (model, rc)
 
subroutine, private regcm_setinitializep1 (model, importstate, exportstate, clock, rc)
 
subroutine, private regcm_setinitializep2 (model, importstate, exportstate, clock, rc)
 
subroutine, private regcm_datainit (model, rc)
 
subroutine, private regcm_setclock (model, rc)
 
subroutine, private regcm_setrunclock (model, rc)
 
subroutine, private regcm_checkimport (model, rc)
 
subroutine, private regcm_setgridarrays (ng, model, rc)
 
subroutine, private regcm_setstates (ng, model, rc)
 
subroutine, private regcm_modeladvance (model, rc)
 
subroutine, private regcm_setfinalize (model, importstate, exportstate, clock, rc)
 
subroutine, private regcm_import (ng, model, rc)
 
subroutine, private regcm_export (ng, model, rc)
 
subroutine, private regcm_uvrot (u, v)
 

Function/Subroutine Documentation

◆ atm_setservices()

subroutine, public esmf_regcm_mod::atm_setservices ( type (esmf_gridcomp) model,
integer, intent(out) rc )

Definition at line 115 of file esmf_atm_regcm.h.

116!
117!=======================================================================
118! !
119! Sets RegCM component shared-object entry points for "initialize", !
120! "run", and "finalize" by using NUOPC generic methods. !
121! !
122!=======================================================================
123!
124! Imported variable declarations.
125!
126 integer, intent(out) :: rc
127!
128 TYPE (ESMF_GridComp) :: model
129!
130! Local variable declarations.
131!
132 character (len=*), parameter :: MyFile = &
133 & __FILE__//", ATM_SetServices"
134!
135!-----------------------------------------------------------------------
136! Initialize return code flag to success state (no error).
137!-----------------------------------------------------------------------
138!
139 IF (esm_track) THEN
140 WRITE (trac,'(a,a,i0)') '==> Entering ATM_SetServices', &
141 & ', PET', petrank
142 FLUSH (trac)
143 END IF
144 rc=esmf_success
145!
146!-----------------------------------------------------------------------
147! Register NUOPC generic routines.
148!-----------------------------------------------------------------------
149!
150 CALL nuopc_compderive (model, &
151 & nuopc_setservices, &
152 & rc=rc)
153 IF (esmf_logfounderror(rctocheck=rc, &
154 & msg=esmf_logerr_passthru, &
155 & line=__line__, &
156 & file=myfile)) THEN
157 RETURN
158 END IF
159!
160!-----------------------------------------------------------------------
161! Register initialize routines.
162!-----------------------------------------------------------------------
163!
164! Set routine for Phase 1 initialization (import and export fields).
165!
166 CALL nuopc_compsetentrypoint (model, &
167 & methodflag=esmf_method_initialize, &
168 & phaselabellist=(/"IPDv00p1"/), &
169 & userroutine=regcm_initializep1, &
170 & rc=rc)
171 IF (esmf_logfounderror(rctocheck=rc, &
172 & msg=esmf_logerr_passthru, &
173 & line=__line__, &
174 & file=myfile)) THEN
175 RETURN
176 END IF
177!
178! Set routine for Phase 2 initialization.
179!
180 CALL nuopc_compsetentrypoint (model, &
181 & methodflag=esmf_method_initialize, &
182 & phaselabellist=(/"IPDv00p2"/), &
183 & userroutine=regcm_initializep2, &
184 & rc=rc)
185 IF (esmf_logfounderror(rctocheck=rc, &
186 & msg=esmf_logerr_passthru, &
187 & line=__line__, &
188 & file=myfile)) THEN
189 RETURN
190 END IF
191!
192!-----------------------------------------------------------------------
193! Attach RegCM component phase independent specializing methods.
194!-----------------------------------------------------------------------
195!
196! Set routine for export initial/restart fields.
197!
198 CALL nuopc_compspecialize (model, &
199 & speclabel=nuopc_label_datainitialize, &
200 & specroutine=regcm_datainit, &
201 & rc=rc)
202 IF (esmf_logfounderror(rctocheck=rc, &
203 & msg=esmf_logerr_passthru, &
204 & line=__line__, &
205 & file=myfile)) THEN
206 RETURN
207 END IF
208!
209! Set routine for setting RegCM clock.
210!
211 CALL nuopc_compspecialize (model, &
212 & speclabel=nuopc_label_setclock, &
213 & specroutine=regcm_setclock, &
214 & rc=rc)
215 IF (esmf_logfounderror(rctocheck=rc, &
216 & msg=esmf_logerr_passthru, &
217 & line=__line__, &
218 & file=myfile)) THEN
219 RETURN
220 END IF
221
222# ifdef ESM_SETRUNCLOCK
223!
224! Set routine for setting COAMPS run clock manually. First, remove the
225! default.
226!
227 CALL esmf_methodremove (model, &
228 & nuopc_label_setrunclock, &
229 & rc=rc)
230 IF (esmf_logfounderror(rctocheck=rc, &
231 & msg=esmf_logerr_passthru, &
232 & line=__line__, &
233 & file=myfile)) THEN
234 RETURN
235 END IF
236!
237 CALL nuopc_compspecialize (model, &
238 & speclabel=nuopc_label_setrunclock, &
239 & specroutine=regcm_setrunclock, &
240 & rc=rc)
241 IF (esmf_logfounderror(rctocheck=rc, &
242 & msg=esmf_logerr_passthru, &
243 & line=__line__, &
244 & file=myfile)) THEN
245 RETURN
246 END IF
247# endif
248!
249! Set routine for checking import state.
250!
251 CALL nuopc_compspecialize (model, &
252 & speclabel=nuopc_label_checkimport, &
253 & specphaselabel="RunPhase1", &
254 & specroutine=regcm_checkimport, &
255 & rc=rc)
256 IF (esmf_logfounderror(rctocheck=rc, &
257 & msg=esmf_logerr_passthru, &
258 & line=__line__, &
259 & file=myfile)) THEN
260 RETURN
261 END IF
262!
263! Set routine for time-stepping RegCM component.
264!
265 CALL nuopc_compspecialize (model, &
266 & speclabel=nuopc_label_advance, &
267 & specroutine=regcm_modeladvance, &
268 & rc=rc)
269 IF (esmf_logfounderror(rctocheck=rc, &
270 & msg=esmf_logerr_passthru, &
271 & line=__line__, &
272 & file=myfile)) THEN
273 RETURN
274 END IF
275!
276!-----------------------------------------------------------------------
277! Register RegCM finalize routine.
278!-----------------------------------------------------------------------
279!
280 CALL esmf_gridcompsetentrypoint (model, &
281 & methodflag=esmf_method_finalize, &
282 & userroutine=regcm_setfinalize, &
283 & rc=rc)
284 IF (esmf_logfounderror(rctocheck=rc, &
285 & msg=esmf_logerr_passthru, &
286 & line=__line__, &
287 & file=myfile)) THEN
288 RETURN
289 END IF
290!
291 IF (esm_track) THEN
292 WRITE (trac,'(a,a,i0)') '<== Exiting RegCM_SetServices', &
293 & ', PET', petrank
294 FLUSH (trac)
295 END IF
296!
297 RETURN

References mod_esmf_esm::esm_track, mod_esmf_esm::petrank, regcm_checkimport(), regcm_datainit(), regcm_modeladvance(), regcm_setclock(), regcm_setfinalize(), regcm_setrunclock(), and mod_esmf_esm::trac.

Here is the call graph for this function:

◆ regcm_checkimport()

subroutine, private esmf_regcm_mod::regcm_checkimport ( type (esmf_gridcomp) model,
integer, intent(out) rc )
private

Definition at line 1108 of file esmf_atm_regcm.h.

1109!
1110!=======================================================================
1111! !
1112! Checks if RegCM component import field is at the correct time. !
1113! !
1114!=======================================================================
1115!
1116! Imported variable declarations.
1117!
1118 integer, intent(out) :: rc
1119!
1120 TYPE (ESMF_GridComp) :: model
1121!
1122! Local variable declarations.
1123!
1124 logical :: IsValid, atCorrectTime
1125!
1126 integer :: ImportCount, i, is, localPET, ng
1127!
1128 real (dp) :: TcurrentInSeconds
1129!
1130 character (len=22) :: DriverTimeString, FieldTimeString
1131
1132 character (len=*), parameter :: MyFile = &
1133 & __FILE__//", RegCM_CheckImport"
1134!
1135 character (ESMF_MAXSTR) :: string, FieldName
1136 character (ESMF_MAXSTR), allocatable :: ImportNameList(:)
1137!
1138 TYPE (ESMF_Clock) :: DriverClock
1139 TYPE (ESMF_Field) :: field
1140 TYPE (ESMF_Time) :: StartTime, CurrentTime
1141 TYPE (ESMF_Time) :: DriverTime, FieldTime
1142 TYPE (ESMF_TimeInterval) :: TimeStep
1143 TYPE (ESMF_VM) :: vm
1144!
1145!-----------------------------------------------------------------------
1146! Initialize return code flag to success state (no error).
1147!-----------------------------------------------------------------------
1148!
1149 IF (esm_track) THEN
1150 WRITE (trac,'(a,a,i0)') '==> Entering RegCM_CheckImport', &
1151 & ', PET', petrank
1152 FLUSH (trac)
1153 END IF
1154 rc=esmf_success
1155!
1156!-----------------------------------------------------------------------
1157! Query component.
1158!-----------------------------------------------------------------------
1159!
1160 CALL nuopc_modelget (model, &
1161 & driverclock=driverclock, &
1162 & rc=rc)
1163 IF (esmf_logfounderror(rctocheck=rc, &
1164 & msg=esmf_logerr_passthru, &
1165 & line=__line__, &
1166 & file=myfile)) THEN
1167 RETURN
1168 END IF
1169!
1170 CALL esmf_gridcompget (model, &
1171 & localpet=localpet, &
1172 & vm=vm, &
1173 & rc=rc)
1174 IF (esmf_logfounderror(rctocheck=rc, &
1175 & msg=esmf_logerr_passthru, &
1176 & line=__line__, &
1177 & file=myfile)) THEN
1178 RETURN
1179 END IF
1180!
1181!-----------------------------------------------------------------------
1182! Get the start time and current time from driver clock.
1183!-----------------------------------------------------------------------
1184!
1185 CALL esmf_clockget (driverclock, &
1186 & timestep=timestep, &
1187 & starttime=starttime, &
1188 & currtime=drivertime, &
1189 & rc=rc)
1190 IF (esmf_logfounderror(rctocheck=rc, &
1191 & msg=esmf_logerr_passthru, &
1192 & line=__line__, &
1193 & file=myfile)) THEN
1194 RETURN
1195 END IF
1196!
1197 CALL esmf_timeget (drivertime, &
1198 & s_r8=tcurrentinseconds, &
1199 & timestringisofrac=drivertimestring, &
1200 & rc=rc)
1201 IF (esmf_logfounderror(rctocheck=rc, &
1202 & msg=esmf_logerr_passthru, &
1203 & line=__line__, &
1204 & file=myfile)) THEN
1205 RETURN
1206 END IF
1207 is=index(drivertimestring, 'T') ! remove 'T' in
1208 IF (is.gt.0) drivertimestring(is:is)=' ' ! ISO 8601 format
1209!
1210!-----------------------------------------------------------------------
1211! Get list of import fields.
1212!-----------------------------------------------------------------------
1213!
1214 IF (nimport(iatmos).gt.0) THEN
1215 nested_loop : DO ng=1,models(iatmos)%Ngrids
1216 IF (any(coupled(iatmos)%LinkedGrid(ng,:))) THEN
1217 CALL esmf_stateget (models(iatmos)%ImportState(ng), &
1218 & itemcount=importcount, &
1219 & rc=rc)
1220 IF (esmf_logfounderror(rctocheck=rc, &
1221 & msg=esmf_logerr_passthru, &
1222 & line=__line__, &
1223 & file=myfile)) THEN
1224 RETURN
1225 END IF
1226!
1227 IF (.not.allocated(importnamelist)) THEN
1228 allocate ( importnamelist(importcount) )
1229 END IF
1230!
1231 CALL esmf_stateget (models(iatmos)%ImportState(ng), &
1232 & itemnamelist=importnamelist, &
1233 & rc=rc)
1234 IF (esmf_logfounderror(rctocheck=rc, &
1235 & msg=esmf_logerr_passthru, &
1236 & line=__line__, &
1237 & file=myfile)) THEN
1238 RETURN
1239 END IF
1240!
1241!-----------------------------------------------------------------------
1242! Only check fields in the ImportState object.
1243!-----------------------------------------------------------------------
1244!
1245 field_loop : DO i=1,importcount
1246 fieldname=trim(importnamelist(i))
1247 CALL esmf_stateget (models(iatmos)%ImportState(ng), &
1248 & itemname=trim(fieldname), &
1249 & field=field, &
1250 & rc=rc)
1251 IF (esmf_logfounderror(rctocheck=rc, &
1252 & msg=esmf_logerr_passthru, &
1253 & line=__line__, &
1254 & file=myfile)) THEN
1255 RETURN
1256 END IF
1257!
1258! If debugging, report field timestamp.
1259!
1260 IF (debuglevel.gt.1) THEN
1261 CALL nuopc_gettimestamp (field, &
1262 & isvalid = isvalid, &
1263 & time = fieldtime, &
1264 & rc = rc)
1265 IF (esmf_logfounderror(rctocheck=rc, &
1266 & msg=esmf_logerr_passthru, &
1267 & line=__line__, &
1268 & file=myfile)) THEN
1269 RETURN
1270 END IF
1271!
1272 IF (isvalid) THEN
1273 CALL esmf_timeget (fieldtime, &
1274 & timestringisofrac = fieldtimestring, &
1275 & rc=rc)
1276 IF (esmf_logfounderror(rctocheck=rc, &
1277 & msg=esmf_logerr_passthru, &
1278 & line=__line__, &
1279 & file=myfile)) THEN
1280 RETURN
1281 END IF
1282 is=index(fieldtimestring, 'T') ! remove 'T'
1283 IF (is.gt.0) fieldtimestring(is:is)=' '
1284!
1285 IF (localpet.eq.0) THEN
1286 WRITE (cplout,10) trim(fieldname), &
1287 & trim(fieldtimestring), &
1288 & trim(drivertimestring)
1289 END IF
1290 END IF
1291 END IF
1292!
1293! Check if import field is at the correct time.
1294!
1295 string='COAMPS_CheckImport - '//trim(fieldname)//' field'
1296 currenttime=drivertime
1297!
1298 atcorrecttime=nuopc_isattime(field, &
1299 & currenttime, &
1300 & rc=rc)
1301 IF (esmf_logfounderror(rctocheck=rc, &
1302 & msg=esmf_logerr_passthru, &
1303 & line=__line__, &
1304 & file=myfile)) THEN
1305 RETURN
1306 END IF
1307!
1308 IF (.not.atcorrecttime) THEN
1309 CALL report_timestamp (field, currenttime, &
1310 & localpet, trim(string), rc)
1311!
1312 string='NUOPC INCOMPATIBILITY DETECTED: Import '// &
1313 & 'Fields not at correct time'
1314 CALL esmf_logseterror(esmf_rc_not_valid, &
1315 & msg=trim(string), &
1316 & line=__line__, &
1317 & file=myfile, &
1318 & rctoreturn=rc)
1319 RETURN
1320 END IF
1321 END DO field_loop
1322 IF (allocated(importnamelist)) deallocate (importnamelist)
1323 END IF
1324 END DO nested_loop
1325 END IF
1326!
1327 IF (esm_track) THEN
1328 WRITE (trac,'(a,a,i0)') '<== Exiting RegCM_CheckImport', &
1329 & ', PET', petrank
1330 FLUSH (trac)
1331 END IF
1332!
1333 10 FORMAT (1x,'RegCM_CheckImport - ',a,':',t32,'TimeStamp = ',a, &
1334 & ', DriverTime = ',a)
1335!
1336 RETURN

References mod_esmf_esm::coupled, mod_esmf_esm::cplout, mod_esmf_esm::debuglevel, mod_esmf_esm::esm_track, mod_esmf_esm::iatmos, mod_esmf_esm::models, mod_esmf_esm::nimport, mod_esmf_esm::petrank, mod_esmf_esm::report_timestamp(), mod_esmf_esm::timestep, and mod_esmf_esm::trac.

Referenced by atm_setservices().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ regcm_datainit()

subroutine, private esmf_regcm_mod::regcm_datainit ( type (esmf_gridcomp) model,
integer, intent(out) rc )
private

Definition at line 555 of file esmf_atm_regcm.h.

556!
557!=======================================================================
558! !
559! Exports RegCM component fields during initialization or restart. !
560! !
561!=======================================================================
562!
563 USE mod_runparams, ONLY : dtsrf
564!
565! Imported variable declarations.
566!
567 integer, intent(out) :: rc
568!
569 TYPE (ESMF_GridComp) :: model
570!
571! Local variable declarations.
572!
573 integer :: is, ng
574 integer :: localPET, PETcount, phase
575!
576 real (dp) :: tstr
577!
578 character (len=*), parameter :: MyFile = &
579 & __FILE__//", RegCM_DataInit"
580
581 character (ESMF_MAXSTR) :: str1, str2
582!
583 TYPE (ESMF_Clock) :: clock
584 TYPE (ESMF_Time) :: CurrentTime
585 TYPE (ESMF_TimeInterval) :: TimeStep
586!
587!-----------------------------------------------------------------------
588! Initialize return code flag to success state (no error).
589!-----------------------------------------------------------------------
590!
591 IF (esm_track) THEN
592 WRITE (trac,'(a,a,i0)') '==> Entering RegCM_DataInit', &
593 & ', PET', petrank
594 FLUSH (trac)
595 END IF
596 rc=esmf_success
597!
598!-----------------------------------------------------------------------
599! Get gridded component clock.
600!-----------------------------------------------------------------------
601!
602 CALL esmf_gridcompget (model, &
603 & clock=clock, &
604 & localpet=localpet, &
605 & petcount=petcount, &
606 & currentphase=phase, &
607 & rc=rc)
608 IF (esmf_logfounderror(rctocheck=rc, &
609 & msg=esmf_logerr_passthru, &
610 & line=__line__, &
611 & file=myfile)) THEN
612 RETURN
613 END IF
614!
615 CALL esmf_clockget (clock, &
616 & currtime=currenttime, &
617 & timestep=timestep, &
618 & rc=rc)
619 IF (esmf_logfounderror(rctocheck=rc, &
620 & msg=esmf_logerr_passthru, &
621 & line=__line__, &
622 & file=myfile)) THEN
623 RETURN
624 END IF
625!
626!-----------------------------------------------------------------------
627! If initialization or restart, put export fields.
628!-----------------------------------------------------------------------
629!
630 IF (restarted.and.(currenttime.eq.esm_restarttime)) THEN
631!
632! Debugging: report time information.
633!
634 IF ((debuglevel.ge.0).and.(localpet.eq.0)) THEN
635 CALL esmf_timeget (currenttime, &
636 & timestringisofrac=str1, &
637 & rc=rc)
638 IF (esmf_logfounderror(rctocheck=rc, &
639 & msg=esmf_logerr_passthru, &
640 & line=__line__, &
641 & file=myfile)) THEN
642 RETURN
643 END IF
644 is=index(str1, 'T') ! remove 'T' in
645 IF (is.gt.0) str1(is:is)=' ' ! ISO 8601 format
646!
647 CALL esmf_timeget (currenttime+timestep, &
648 & timestringisofrac=str2, &
649 & rc=rc)
650 IF (esmf_logfounderror(rctocheck=rc, &
651 & msg=esmf_logerr_passthru, &
652 & line=__line__, &
653 & file=myfile)) THEN
654 RETURN
655 END IF
656 is=index(str2, 'T') ! remove 'T' in
657 IF (is.gt.0) str2(is:is)=' ' ! ISO 8601 format
658!
659 IF (debuglevel.eq.0) THEN
660 WRITE (cplout,10) trim(str1), trim(str2), phase
661 ELSE
662 WRITE (cplout,20) trim(str1), trim(str2), phase, 0.0_r8, &
663 & dtsrf
664 END IF
665 END IF
666!
667! Run RegCM component only for one time-step to fill variables.
668!
669 CALL regcm_run (0.0_r8, dtsrf)
670!
671! Put export fields.
672!
673 IF (nexport(iatmos).gt.0) THEN
674 DO ng=1,models(iatmos)%Ngrids
675 CALL regcm_export (ng, model, rc=rc)
676 IF (esmf_logfounderror(rctocheck=rc, &
677 & msg=esmf_logerr_passthru, &
678 & line=__line__, &
679 & file=myfile)) THEN
680 RETURN
681 END IF
682 END DO
683 END IF
684 END IF
685!
686 IF (esm_track) THEN
687 WRITE (trac,'(a,a,i0)') '<== Exiting RegCM_DataInit', &
688 & ', PET', petrank
689 FLUSH (trac)
690 END IF
691 IF (debuglevel.gt.0) FLUSH (cplout)
692!
693 10 FORMAT (/,' ESMF, Running RegCM: ',a,' --> ',a,' Phase: ',i1)
694 20 FORMAT (/,' ESMF, Running RegCM: ',a,' --> ',a,' Phase: ',i1, &
695 & ' [',f15.2, '-', f15.2, ' s]')
696
697 RETURN

References mod_esmf_esm::cplout, mod_esmf_esm::debuglevel, mod_esmf_esm::esm_track, mod_esmf_esm::iatmos, mod_esmf_esm::models, mod_esmf_esm::nexport, mod_esmf_esm::petrank, regcm_export(), mod_esmf_esm::timestep, and mod_esmf_esm::trac.

Referenced by atm_setservices().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ regcm_export()

subroutine, private esmf_regcm_mod::regcm_export ( integer, intent(in) ng,
type (esmf_gridcomp) model,
integer, intent(out) rc )
private

Definition at line 2769 of file esmf_atm_regcm.h.

2770!
2771!=======================================================================
2772! !
2773! Exports RegCM fields to other coupled gridded components. !
2774! !
2775!=======================================================================
2776!
2777 USE mod_update, ONLY : exportfields
2778 USE mod_dynparam, ONLY : ici1, ici2, jci1, jci2
2779!
2780! Imported variable declarations.
2781!
2782 integer, intent(in) :: ng
2783 integer, intent(out) :: rc
2784
2785 TYPE (ESMF_GridComp) :: model
2786!
2787! Local variable declarations.
2788!
2789 integer :: ifld, i, is, j
2790 integer :: Istr, Iend, Jstr, Jend
2791 integer :: year, month, day, hour, minutes, seconds, sN, SD
2792 integer :: ExportCount
2793 integer :: localDE, localDEcount, localPET, PETcount
2794!
2795 real (dp), parameter :: pi = 3.14159265358979323846_dp
2796 real (dp) :: Fseconds, TimeInDays, Time_Current
2797
2798 real (dp) :: MyFmax(1), MyFmin(1), Fmin(1), Fmax(1), Fval
2799!
2800 real (dp), pointer :: ptr2d(:,:) => null()
2801!
2802 character (len=22) :: Time_CurrentString
2803
2804 character (len=*), parameter :: MyFile = &
2805 & __FILE__//", RegCM_Export"
2806!
2807 character (ESMF_MAXSTR) :: cname, ofile
2808 character (ESMF_MAXSTR), allocatable :: ExportNameList(:)
2809!
2810 TYPE (ESMF_Clock) :: clock
2811 TYPE (ESMF_Field) :: field
2812 TYPE (ESMF_Time) :: CurrentTime
2813 TYPE (ESMF_VM) :: vm
2814!
2815!-----------------------------------------------------------------------
2816! Initialize return code flag to success state (no error).
2817!-----------------------------------------------------------------------
2818!
2819 IF (esm_track) THEN
2820 WRITE (trac,'(a,a,i0)') '==> Entering RegCM_Export', &
2821 & ', PET', petrank
2822 FLUSH (trac)
2823 END IF
2824 rc=esmf_success
2825!
2826!-----------------------------------------------------------------------
2827! Get information about the gridded component.
2828!-----------------------------------------------------------------------
2829!
2830 CALL esmf_gridcompget (model, &
2831 & clock=clock, &
2832 & localpet=localpet, &
2833 & petcount=petcount, &
2834 & vm=vm, &
2835 & name=cname, &
2836 & rc=rc)
2837 IF (esmf_logfounderror(rctocheck=rc, &
2838 & msg=esmf_logerr_passthru, &
2839 & line=__line__, &
2840 & file=myfile)) THEN
2841 RETURN
2842 END IF
2843!
2844! Get number of local decomposition elements (DEs). Usually, a single
2845! DE is associated with each Persistent Execution Thread (PETs). Thus,
2846! localDEcount=1.
2847!
2848 CALL esmf_gridget (models(iatmos)%grid(ng), &
2849 & localdecount=localdecount, &
2850 & rc=rc)
2851 IF (esmf_logfounderror(rctocheck=rc, &
2852 & msg=esmf_logerr_passthru, &
2853 & line=__line__, &
2854 & file=myfile)) THEN
2855 RETURN
2856 END IF
2857!
2858!-----------------------------------------------------------------------
2859! Get current time.
2860!-----------------------------------------------------------------------
2861!
2862 CALL esmf_clockget (clock, &
2863 & currtime=currenttime, &
2864 & rc=rc)
2865 IF (esmf_logfounderror(rctocheck=rc, &
2866 & msg=esmf_logerr_passthru, &
2867 & line=__line__, &
2868 & file=myfile)) THEN
2869 RETURN
2870 END IF
2871!
2872 CALL esmf_timeget (currenttime, &
2873 & yy=year, &
2874 & mm=month, &
2875 & dd=day, &
2876 & h =hour, &
2877 & m =minutes, &
2878 & s =seconds, &
2879 & sn=sn, &
2880 & sd=sd, &
2881 & rc=rc)
2882 IF (esmf_logfounderror(rctocheck=rc, &
2883 & msg=esmf_logerr_passthru, &
2884 & line=__line__, &
2885 & file=myfile)) THEN
2886 RETURN
2887 END IF
2888!
2889 CALL esmf_timeget (currenttime, &
2890 & s_r8=time_current, &
2891 & timestring=time_currentstring, &
2892 & rc=rc)
2893 IF (esmf_logfounderror(rctocheck=rc, &
2894 & msg=esmf_logerr_passthru, &
2895 & line=__line__, &
2896 & file=myfile)) THEN
2897 RETURN
2898 END IF
2899 fseconds=real(seconds,dp)+real(sn,dp)/real(sd,dp)
2900 timeindays=time_current/86400.0_dp
2901 is=index(time_currentstring, 'T') ! remove 'T' in
2902 IF (is.gt.0) time_currentstring(is:is)=' ' ! ISO 8601 format
2903!
2904!-----------------------------------------------------------------------
2905! Get list of export fields.
2906!-----------------------------------------------------------------------
2907!
2908 CALL esmf_stateget (models(iatmos)%ExportState(ng), &
2909 & itemcount=exportcount, &
2910 & rc=rc)
2911 IF (esmf_logfounderror(rctocheck=rc, &
2912 & msg=esmf_logerr_passthru, &
2913 & line=__line__, &
2914 & file=myfile)) THEN
2915 RETURN
2916 END IF
2917!
2918 IF (.not. allocated(exportnamelist)) THEN
2919 allocate ( exportnamelist(exportcount) )
2920 END IF
2921 CALL esmf_stateget (models(iatmos)%ExportState(ng), &
2922 & itemnamelist=exportnamelist, &
2923 & rc=rc)
2924 IF (esmf_logfounderror(rctocheck=rc, &
2925 & msg=esmf_logerr_passthru, &
2926 & line=__line__, &
2927 & file=myfile)) THEN
2928 RETURN
2929 END IF
2930!
2931!-----------------------------------------------------------------------
2932! Rotate vector components (wind, wind stress) to a rectangular
2933! north-south and east-west oriented grid.
2934!-----------------------------------------------------------------------
2935!
2936 CALL regcm_uvrot (exportfields%wndu, exportfields%wndv)
2937 CALL regcm_uvrot (exportfields%taux, exportfields%tauy)
2938!
2939!-----------------------------------------------------------------------
2940! Load export fields.
2941!-----------------------------------------------------------------------
2942!
2943 fld_loop : DO ifld=1,exportcount
2944!
2945! Get field from export state.
2946!
2947 CALL esmf_stateget (models(iatmos)%ExportState(ng), &
2948 & trim(exportnamelist(ifld)), &
2949 & field, &
2950 & rc=rc)
2951 IF (esmf_logfounderror(rctocheck=rc, &
2952 & msg=esmf_logerr_passthru, &
2953 & line=__line__, &
2954 & file=myfile)) THEN
2955 RETURN
2956 END IF
2957!
2958! Get field pointer. Usually, the DO-loop is executed once since
2959! localDEcount=1.
2960!
2961 de_loop : DO localde=0,localdecount-1
2962 CALL esmf_fieldget (field, &
2963 & localde=localde, &
2964 & farrayptr=ptr2d, &
2965 & rc=rc)
2966 IF (esmf_logfounderror(rctocheck=rc, &
2967 & msg=esmf_logerr_passthru, &
2968 & line=__line__, &
2969 & file=myfile)) THEN
2970 RETURN
2971 END IF
2972 istr=ici1
2973 iend=ici2
2974 jstr=jci1
2975 jend=jci2
2976!
2977! Initialize pointer to missing value.
2978!
2979 ptr2d=missing_dp
2980!
2981! Load field data into export state. Notice that all export fields
2982! are kept as computed by RegCM. The imported component does the
2983! proper scaling, physical units conversion, and other manipulations.
2984! It is done to avoid applying such transformations twice.
2985!
2986 SELECT CASE (trim(adjustl(exportnamelist(ifld))))
2987!
2988! Surface air pressure (hPa or mb).
2989!
2990 CASE ('psfc', 'Pair')
2991 fval=exportfields%psfc(jstr,istr)
2992 myfmin(1)=fval
2993 myfmax(1)=fval
2994 DO i=istr,iend
2995 DO j=jstr,jend
2996 fval=exportfields%psfc(j,i)
2997 myfmin(1)=min(myfmin(1),fval)
2998 myfmax(1)=max(myfmax(1),fval)
2999 ptr2d(i,j)=fval
3000 END DO
3001 END DO
3002!
3003! 2 meter surface air temperature (K).
3004!
3005 CASE ('tsfc', 'Tair')
3006 fval=exportfields%tsfc(jstr,istr)
3007 myfmin(1)=fval
3008 myfmax(1)=fval
3009 DO i=istr,iend
3010 DO j=jstr,jend
3011 fval=exportfields%tsfc(j,i)
3012 myfmin(1)=min(myfmin(1),fval)
3013 myfmax(1)=max(myfmax(1),fval)
3014 ptr2d(i,j)=fval
3015 END DO
3016 END DO
3017!
3018! 2 meter surface specific humidity (kg/kg)
3019!
3020 CASE ('qsfc', 'Hair')
3021 fval=exportfields%qsfc(jstr,istr)
3022 myfmin(1)=fval
3023 myfmax(1)=fval
3024 DO i=istr,iend
3025 DO j=jstr,jend
3026 fval=exportfields%qsfc(j,i)
3027 myfmin(1)=min(myfmin(1),fval)
3028 myfmax(1)=max(myfmax(1),fval)
3029 ptr2d(i,j)=fval
3030 END DO
3031 END DO
3032!
3033! Surface net longwave radiation (W m-2)
3034!
3035 CASE ('lwrd', 'LWrad')
3036 fval=exportfields%lwrd(jstr,istr)
3037 myfmin(1)=fval
3038 myfmax(1)=fval
3039 DO i=istr,iend
3040 DO j=jstr,jend
3041 fval=exportfields%lwrd(j,i)
3042 myfmin(1)=min(myfmin(1),fval)
3043 myfmax(1)=max(myfmax(1),fval)
3044 ptr2d(i,j)=fval
3045 END DO
3046 END DO
3047!
3048! Surface downward longwave radiation (W m-2).
3049!
3050 CASE ('dlwrd', 'dLWrad', 'lwrad_down')
3051 fval=exportfields%dlwr(jstr,istr)
3052 myfmin(1)=fval
3053 myfmax(1)=fval
3054 DO i=istr,iend
3055 DO j=jstr,jend
3056 fval=exportfields%dlwr(j,i)
3057 myfmin(1)=min(myfmin(1),fval)
3058 myfmax(1)=max(myfmax(1),fval)
3059 ptr2d(i,j)=fval
3060 END DO
3061 END DO
3062!
3063! Surface latent heat flux (W m-2)
3064!
3065 CASE ('lhfx', 'LHfx')
3066 fval=exportfields%lhfx(jstr,istr)
3067 myfmin(1)=fval
3068 myfmax(1)=fval
3069 DO i=istr,iend
3070 DO j=jstr,jend
3071 fval=exportfields%lhfx(j,i)
3072 myfmin(1)=min(myfmin(1),fval)
3073 myfmax(1)=max(myfmax(1),fval)
3074 ptr2d(i,j)=fval
3075 END DO
3076 END DO
3077!
3078! Surface sensible heat flux (W m-2)
3079!
3080 CASE ('shfx')
3081 fval=exportfields%shfx(jstr,istr)
3082 myfmin(1)=fval
3083 myfmax(1)=fval
3084 DO i=istr,iend
3085 DO j=jstr,jend
3086 fval=exportfields%shfx(j,i)
3087 myfmin(1)=min(myfmin(1),fval)
3088 myfmax(1)=max(myfmax(1),fval)
3089 ptr2d(i,j)=fval
3090 END DO
3091 END DO
3092!
3093! Total precipitation rate (m s-1).
3094!
3095 CASE ('prec')
3096 fval=exportfields%prec(jstr,istr)
3097 myfmin(1)=fval
3098 myfmax(1)=fval
3099 DO i=istr,iend
3100 DO j=jstr,jend
3101 fval=exportfields%prec(j,i)
3102 myfmin(1)=min(myfmin(1),fval)
3103 myfmax(1)=max(myfmax(1),fval)
3104 ptr2d(i,j)=fval
3105 END DO
3106 END DO
3107!
3108! Surface eastward wind component (m s-1).
3109!
3110 CASE ('Uwind', 'u10', 'wndu')
3111 fval=exportfields%wndu(jstr,istr)
3112 myfmin(1)=fval
3113 myfmax(1)=fval
3114 DO i=istr,iend
3115 DO j=jstr,jend
3116 fval=exportfields%wndu(j,i)
3117 myfmin(1)=min(myfmin(1),fval)
3118 myfmax(1)=max(myfmax(1),fval)
3119 ptr2d(i,j)=fval
3120 END DO
3121 END DO
3122!
3123! Surface northward wind component (m s-1).
3124!
3125 CASE ('wndv')
3126 fval=exportfields%wndv(jstr,istr)
3127 myfmin(1)=fval
3128 myfmax(1)=fval
3129 DO i=istr,iend
3130 DO j=jstr,jend
3131 fval=exportfields%wndv(j,i)
3132 myfmin(1)=min(myfmin(1),fval)
3133 myfmax(1)=max(myfmax(1),fval)
3134 ptr2d(i,j)=fval
3135 END DO
3136 END DO
3137!
3138! Surface net solar shortwave radiation (W m-2).
3139!
3140 CASE ('swrd')
3141 fval=exportfields%swrd(jstr,istr)
3142 myfmin(1)=fval
3143 myfmax(1)=fval
3144 DO i=istr,iend
3145 DO j=jstr,jend
3146 fval=exportfields%swrd(j,i)
3147 myfmin(1)=min(myfmin(1),fval)
3148 myfmax(1)=max(myfmax(1),fval)
3149 ptr2d(i,j)=fval
3150 END DO
3151 END DO
3152!
3153! Downward shortwave radiation (W m-2).
3154!
3155 CASE ('dswr')
3156 fval=exportfields%dswr(jstr,istr)
3157 myfmin(1)=fval
3158 myfmax(1)=fval
3159 DO i=istr,iend
3160 DO j=jstr,jend
3161 fval=exportfields%dswr(j,i)
3162 myfmin(1)=min(myfmin(1),fval)
3163 myfmax(1)=max(myfmax(1),fval)
3164 ptr2d(i,j)=fval
3165 END DO
3166 END DO
3167!
3168! Surface runoff over land (m s-1).
3169!
3170 CASE ('rnof')
3171 fval=exportfields%rnof(jstr,istr)
3172 myfmin(1)=fval
3173 myfmax(1)=fval
3174 DO i=istr,iend
3175 DO j=jstr,jend
3176 fval=exportfields%rnof(j,i)
3177 myfmin(1)=min(myfmin(1),fval)
3178 myfmax(1)=max(myfmax(1),fval)
3179 ptr2d(i,j)=fval
3180 END DO
3181 END DO
3182!
3183! Sub-surface runoff over land (m s-1).
3184!
3185 CASE ('snof')
3186 fval=exportfields%snof(jstr,istr)
3187 myfmin(1)=fval
3188 myfmax(1)=fval
3189 DO i=istr,iend
3190 DO j=jstr,jend
3191 fval=exportfields%snof(j,i)
3192 myfmin(1)=min(myfmin(1),fval)
3193 myfmax(1)=max(myfmax(1),fval)
3194 ptr2d(i,j)=fval
3195 END DO
3196 END DO
3197!
3198! Surface eastward wind stress component (N m-2 or Pa).
3199!
3200 CASE ('taux')
3201 fval=exportfields%taux(jstr,istr)
3202 myfmin(1)=fval
3203 myfmax(1)=fval
3204 DO i=istr,iend
3205 DO j=jstr,jend
3206 fval=exportfields%taux(j,i)
3207 myfmin(1)=min(myfmin(1),fval)
3208 myfmax(1)=max(myfmax(1),fval)
3209 ptr2d(i,j)=fval
3210 END DO
3211 END DO
3212!
3213! Surface eastward wind stress component (N m-2 or Pa).
3214!
3215 CASE ('tauy')
3216 fval=exportfields%tauy(jstr,istr)
3217 myfmin(1)=fval
3218 myfmax(1)=fval
3219 DO i=istr,iend
3220 DO j=jstr,jend
3221 fval=exportfields%tauy(j,i)
3222 myfmin(1)=min(myfmin(1),fval)
3223 myfmax(1)=max(myfmax(1),fval)
3224 ptr2d(i,j)=fval
3225 END DO
3226 END DO
3227!
3228! Surface wind speed ( m s-1).
3229!
3230 CASE ('wspd')
3231 fval=exportfields%wspd(jstr,istr)
3232 myfmin(1)=fval
3233 myfmax(1)=fval
3234 DO i=istr,iend
3235 DO j=jstr,jend
3236 fval=exportfields%wspd(j,i)
3237 myfmin(1)=min(myfmin(1),fval)
3238 myfmax(1)=max(myfmax(1),fval)
3239 ptr2d(i,j)=fval
3240 END DO
3241 END DO
3242!
3243! Surface wind direction (radians).
3244!
3245 CASE ('wdir')
3246 fval=atan2(exportfields%wndu(jstr,istr), &
3247 & exportfields%wndv(jstr,istr))
3248 myfmin(1)=fval
3249 myfmax(1)=fval
3250 DO i=istr,iend
3251 DO j=jstr,jend
3252 fval=atan2(exportfields%wndu(j,i), &
3253 & exportfields%wndv(j,i))
3254 IF (dd.lt.0.0_r8) fval=fval+2.0_r8*pi
3255 myfmin(1)=min(myfmin(1),fval)
3256 myfmax(1)=max(myfmax(1),fval)
3257 ptr2d(i,j)=fval
3258 END DO
3259 END DO
3260!
3261! Surface net heat flux (W m-2).
3262!
3263 CASE ('nflx')
3264 fval=exportfields%nflx(jstr,istr)
3265 myfmin(1)=fval
3266 myfmax(1)=fval
3267 DO i=istr,iend
3268 DO j=jstr,jend
3269 fval=exportfields%nflx(j,i)
3270 myfmin(1)=min(myfmin(1),fval)
3271 myfmax(1)=max(myfmax(1),fval)
3272 ptr2d(i,j)=fval
3273 END DO
3274 END DO
3275!
3276! Surface net freshwater flux: E-P (m s-1).
3277!
3278 CASE ('sflx')
3279 fval=exportfields%sflx(jstr,istr)
3280 myfmin(1)=fval
3281 myfmax(1)=fval
3282 DO i=istr,iend
3283 DO j=jstr,jend
3284 fval=exportfields%sflx(j,i)
3285 myfmin(1)=min(myfmin(1),fval)
3286 myfmax(1)=max(myfmax(1),fval)
3287 ptr2d(i,j)=fval
3288 END DO
3289 END DO
3290!
3291! Surface snow concentration.
3292!
3293 CASE ('snow')
3294 fval=exportfields%snow(jstr,istr)
3295 myfmin(1)=fval
3296 myfmax(1)=fval
3297 DO i=istr,iend
3298 DO j=jstr,jend
3299 fval=exportfields%snow(j,i)
3300 myfmin(1)=min(myfmin(1),fval)
3301 myfmax(1)=max(myfmax(1),fval)
3302 ptr2d(i,j)=fval
3303 END DO
3304 END DO
3305!
3306! Export field not found.
3307!
3308 CASE DEFAULT
3309 IF (localpet.eq.0) THEN
3310 WRITE (cplout,10) trim(adjustl(exportnamelist(ifld))), &
3311 & trim(cinpname)
3312 END IF
3313 rc=esmf_rc_not_found
3314 IF (esmf_logfounderror(rctocheck=rc, &
3315 & msg=esmf_logerr_passthru, &
3316 & line=__line__, &
3317 & file=myfile)) THEN
3318 RETURN
3319 END IF
3320 END SELECT
3321!
3322! Nullify pointer to make sure that it does not point on a random
3323! part in the memory.
3324!
3325 IF (associated(ptr2d)) nullify (ptr2d)
3326 END DO de_loop
3327!
3328! Get export field minimun and maximum values.
3329!
3330 CALL esmf_vmallreduce (vm, &
3331 & senddata=myfmin, &
3332 & recvdata=fmin, &
3333 & count=1, &
3334 & reduceflag=esmf_reduce_min, &
3335 & rc=rc)
3336 IF (esmf_logfounderror(rctocheck=rc, &
3337 & msg=esmf_logerr_passthru, &
3338 & line=__line__, &
3339 & file=myfile)) THEN
3340 RETURN
3341 END IF
3342!
3343 CALL esmf_vmallreduce (vm, &
3344 & senddata=myfmax, &
3345 & recvdata=fmax, &
3346 & count=1, &
3347 & reduceflag=esmf_reduce_max, &
3348 & rc=rc)
3349 IF (esmf_logfounderror(rctocheck=rc, &
3350 & msg=esmf_logerr_passthru, &
3351 & line=__line__, &
3352 & file=myfile)) THEN
3353 RETURN
3354 END IF
3355!
3356! Report export field information.
3357!
3358 IF ((debuglevel.ge.0).and.(localpet.eq.0)) THEN
3359 WRITE (cplout,20) trim(exportnamelist(ifld)), &
3360 & trim(time_currentstring), ng, &
3361 & fmin(1), fmax(1)
3362 END IF
3363!
3364! Debugging: write out export field into a NetCDF file.
3365!
3366 IF ((debuglevel.ge.3).and. &
3367 & models(iatmos)%ExportField(ifld)%debug_write) THEN
3368 WRITE (ofile,30) ng, trim(exportnamelist(ifld)), &
3369 & year, month, day, hour, minutes, seconds
3370 CALL esmf_fieldwrite (field, &
3371 & trim(ofile), &
3372 & overwrite=.true., &
3373 & rc=rc)
3374 IF (esmf_logfounderror(rctocheck=rc, &
3375 & msg=esmf_logerr_passthru, &
3376 & line=__line__, &
3377 & file=myfile)) THEN
3378 RETURN
3379 END IF
3380 END IF
3381 END DO fld_loop
3382!
3383! Deallocate local arrays.
3384!
3385 IF (allocated(exportnamelist)) deallocate(exportnamelist)
3386!
3387! Update RegCM export calls counter.
3388!
3389 IF (exportcount.gt.0) THEN
3390 models(iatmos)%ExportCalls=models(iatmos)%ExportCalls+1
3391 END IF
3392!
3393 IF (esm_track) THEN
3394 WRITE (trac,'(a,a,i0)') '<== Exiting RegCM_Export', &
3395 & ', PET', petrank
3396 FLUSH (trac)
3397 END IF
3398 IF (debuglevel.gt.0) FLUSH (cplout)
3399!
3400 10 FORMAT (/,2x,'RegCM_Export - unable to find option to export: ', &
3401 & a,/,18x,'check ''Export(atmos)'' in input script: ',a)
3402 20 FORMAT (2x,'RegCM_Export - ESMF: exporting field ''',a,'''', &
3403 & t72,a,2x,'Grid ',i2.2,/, &
3404 & 19x,'(Cmin = ', 1p,e15.8,0p,' Cmax = ',1p,e15.8,0p,')')
3405 30 FORMAT ('regcm',i2.2,'_export_',a,'_',i4.4,2('-',i2.2),'_', &
3406 & i2.2,2('.',i2.2),'.nc')
3407
3408 RETURN

References mod_esmf_esm::cinpname, mod_esmf_esm::cplout, mod_esmf_esm::debuglevel, mod_esmf_esm::esm_track, mod_esmf_esm::iatmos, mod_esmf_esm::missing_dp, mod_esmf_esm::models, mod_esmf_esm::petrank, regcm_uvrot(), and mod_esmf_esm::trac.

Referenced by regcm_datainit(), and regcm_modeladvance().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ regcm_import()

subroutine, private esmf_regcm_mod::regcm_import ( integer, intent(in) ng,
type (esmf_gridcomp) model,
integer, intent(out) rc )
private

Definition at line 2349 of file esmf_atm_regcm.h.

2350!
2351!=======================================================================
2352! !
2353! Imports fields into RegCM array structure from other coupled !
2354! gridded components. !
2355! !
2356!=======================================================================
2357!
2358 USE mod_update, ONLY : importfields
2359 USE mod_dynparam, ONLY : ici1, ici2, jci1, jci2
2360!
2361! Imported variable declarations.
2362!
2363 integer, intent(in) :: ng
2364 integer, intent(out) :: rc
2365!
2366 TYPE (ESMF_GridComp) :: model
2367!
2368! Local variable declarations.
2369!
2370 integer :: id, ifld, i, is, j
2371 integer :: year, month, day, hour, minutes, seconds, sN, SD
2372 integer :: ImportCount
2373 integer :: localDE, localDEcount, localPET, PETcount
2374 integer :: LBi, UBi, LBj, UBj
2375 integer :: IminP, ImaxP, JminP, JmaxP
2376!
2377 real (dp) :: Fseconds, TimeInDays, Time_Current
2378
2379 real (dp) :: MyFmax(2), MyFmin(2), Fmin(2), Fmax(2), Fval
2380 real (dp) :: scale, add_offset
2381!
2382 real (dp), pointer :: ptr2d(:,:) => null()
2383!
2384 character (len=22 ) :: Time_CurrentString
2385 character (len=100) :: FieldName
2386
2387 character (len=*), parameter :: MyFile = &
2388 & __FILE__//", RegCM_Import"
2389
2390 character (ESMF_MAXSTR) :: cname, ofile
2391 character (ESMF_MAXSTR), allocatable :: ImportNameList(:)
2392!
2393 TYPE (ESMF_Clock) :: clock
2394 TYPE (ESMF_Field) :: field
2395 TYPE (ESMF_State) :: ImportState
2396 TYPE (ESMF_Time) :: CurrentTime
2397 TYPE (ESMF_VM) :: vm
2398!
2399!-----------------------------------------------------------------------
2400! Initialize return code flag to success state (no error).
2401!-----------------------------------------------------------------------
2402!
2403 IF (esm_track) THEN
2404 WRITE (trac,'(a,a,i0)') '==> Entering RegCM_Import', &
2405 & ', PET', petrank
2406 FLUSH (trac)
2407 END IF
2408 rc=esmf_success
2409!
2410!-----------------------------------------------------------------------
2411! Compute RegCM lower and upper bounds (non-overlapping) for physical
2412! area per nested grid and tile.
2413!-----------------------------------------------------------------------
2414!
2415 iminp=ice1
2416 imaxp=ice2
2417 jminp=jce1
2418 jmaxp=jce2
2419!
2420!-----------------------------------------------------------------------
2421! Get information about the gridded component.
2422!-----------------------------------------------------------------------
2423!
2424 CALL esmf_gridcompget (model, &
2425 & importstate=importstate, &
2426 & clock=clock, &
2427 & localpet=localpet, &
2428 & petcount=petcount, &
2429 & vm=vm, &
2430 & name=cname, &
2431 & rc=rc)
2432 IF (esmf_logfounderror(rctocheck=rc, &
2433 & msg=esmf_logerr_passthru, &
2434 & line=__line__, &
2435 & file=myfile)) THEN
2436 RETURN
2437 END IF
2438!
2439! Get number of local decomposition elements (DEs). Usually, a single
2440! DE is associated with each Persistent Execution Thread (PETs). Thus,
2441! localDEcount=1.
2442!
2443 CALL esmf_gridget (models(iatmos)%grid(ng), &
2444 & localdecount=localdecount, &
2445 & rc=rc)
2446 IF (esmf_logfounderror(rctocheck=rc, &
2447 & msg=esmf_logerr_passthru, &
2448 & line=__line__, &
2449 & file=myfile)) THEN
2450 RETURN
2451 END IF
2452!
2453!-----------------------------------------------------------------------
2454! Get current time.
2455!-----------------------------------------------------------------------
2456!
2457 CALL esmf_clockget (clock, &
2458 & currtime=currenttime, &
2459 & rc=rc)
2460 IF (esmf_logfounderror(rctocheck=rc, &
2461 & msg=esmf_logerr_passthru, &
2462 & line=__line__, &
2463 & file=myfile)) THEN
2464 RETURN
2465 END IF
2466!
2467 CALL esmf_timeget (currenttime, &
2468 & yy=year, &
2469 & mm=month, &
2470 & dd=day, &
2471 & h =hour, &
2472 & m =minutes, &
2473 & s =seconds, &
2474 & sn=sn, &
2475 & sd=sd, &
2476 & rc=rc)
2477 IF (esmf_logfounderror(rctocheck=rc, &
2478 & msg=esmf_logerr_passthru, &
2479 & line=__line__, &
2480 & file=myfile)) THEN
2481 RETURN
2482 END IF
2483!
2484 CALL esmf_timeget (currenttime, &
2485 & s_r8=time_current, &
2486 & timestring=time_currentstring, &
2487 & rc=rc)
2488 IF (esmf_logfounderror(rctocheck=rc, &
2489 & msg=esmf_logerr_passthru, &
2490 & line=__line__, &
2491 & file=myfile)) THEN
2492 RETURN
2493 END IF
2494 fseconds=real(seconds,dp)+real(sn,dp)/real(sd,dp)
2495 timeindays=time_current/86400.0_dp
2496 is=index(time_currentstring, 'T') ! remove 'T' in
2497 IF (is.gt.0) time_currentstring(is:is)=' ' ! ISO 8601 format
2498!
2499!-----------------------------------------------------------------------
2500! Get list of import fields.
2501!-----------------------------------------------------------------------
2502!
2503 CALL esmf_stateget (models(iatmos)%ImportState(ng), &
2504 & itemcount=importcount, &
2505 & rc=rc)
2506 IF (esmf_logfounderror(rctocheck=rc, &
2507 & msg=esmf_logerr_passthru, &
2508 & line=__line__, &
2509 & file=myfile)) THEN
2510 RETURN
2511 END IF
2512!
2513 IF (.not.allocated(importnamelist)) THEN
2514 allocate ( importnamelist(importcount) )
2515 END IF
2516 CALL esmf_stateget (models(iatmos)%ImportState(ng), &
2517 & itemnamelist=importnamelist, &
2518 & rc=rc)
2519 IF (esmf_logfounderror(rctocheck=rc, &
2520 & msg=esmf_logerr_passthru, &
2521 & line=__line__, &
2522 & file=myfile)) THEN
2523 RETURN
2524 END IF
2525!
2526!-----------------------------------------------------------------------
2527! Get import fields.
2528!-----------------------------------------------------------------------
2529!
2530 fld_loop : DO ifld=1,importcount
2531 id=field_index(models(iatmos)%ImportField, importnamelist(ifld))
2532!
2533! Get field from import state.
2534!
2535 CALL esmf_stateget (models(iatmos)%ImportState(ng), &
2536 & trim(importnamelist(ifld)), &
2537 & field, &
2538 & rc=rc)
2539 IF (esmf_logfounderror(rctocheck=rc, &
2540 & msg=esmf_logerr_passthru, &
2541 & line=__line__, &
2542 & file=myfile)) THEN
2543 RETURN
2544 END IF
2545!
2546! Get field pointer. Usually, the DO-loop is executed once since
2547! localDEcount=1.
2548!
2549 de_loop : DO localde=0,localdecount-1
2550 CALL esmf_fieldget (field, &
2551 & localde=localde, &
2552 & farrayptr=ptr2d, &
2553 & rc=rc)
2554 IF (esmf_logfounderror(rctocheck=rc, &
2555 & msg=esmf_logerr_passthru, &
2556 & line=__line__, &
2557 & file=myfile)) THEN
2558 RETURN
2559 END IF
2560 lbi=lbound(ptr2d,1)
2561 ubi=ubound(ptr2d,1)
2562 lbj=lbound(ptr2d,2)
2563 ubj=ubound(ptr2d,2)
2564!
2565! Initialize import field parameters. Set "scale" and "add_offset"
2566! values need to convert imported fields to RegCM requirements.
2567!
2568 scale =models(iatmos)%ImportField(id)%scale_factor
2569 add_offset=models(iatmos)%ImportField(id)%add_offset
2570!
2571 fval=ptr2d(iminp,jminp)
2572 myfmin(1)=fval
2573 myfmax(1)=fval
2574 myfmin(2)=fval
2575 myfmax(2)=fval
2576!
2577! Load import data into RegCM component variable.
2578!
2579 SELECT CASE (trim(adjustl(importnamelist(ifld))))
2580!
2581! Sea surface temperature from OCN component (C).
2582!
2583 CASE ('sst', 'SST')
2584 fieldname=trim(importnamelist(ifld))
2585 DO i=iminp,imaxp
2586 DO j=jminp,jmaxp
2587 myfmin(1)=min(myfmin(1),ptr2d(i,j))
2588 myfmax(1)=max(myfmax(1),ptr2d(i,j))
2589 fval=scale*ptr2d(i,j)+add_offset
2590 myfmin(2)=min(myfmin(2),fval)
2591 myfmax(2)=max(myfmax(2),fval)
2592 importfields%sst(j,i)=fval
2593 END DO
2594 END DO
2595
2596# if defined SEA_ICE || defined OCNICE
2597!
2598! Sea ice thickness (m), from either OCN or ICE.
2599!
2600 CASE ('sit')
2601 DO i=iminp,imaxp
2602 DO j=jminp,jmaxp
2603 myfmin(1)=min(myfmin(1),ptr2d(i,j))
2604 myfmax(1)=max(myfmax(1),ptr2d(i,j))
2605 fval=scale*ptr2d(i,j)+add_offset
2606 myfmin(2)=min(myfmin(2),fval)
2607 myfmax(2)=max(myfmax(2),fval)
2608 importfields%sit(j,i)=fval
2609 END DO
2610 END DO
2611# endif
2612!
2613! Land-sea or wet-dry mask (nondimensional), from OCN.
2614!
2615 CASE ('msk')
2616 DO i=iminp,imaxp
2617 DO j=jminp,jmaxp
2618 myfmin(1)=min(myfmin(1),ptr2d(i,j))
2619 myfmax(1)=max(myfmax(1),ptr2d(i,j))
2620 fval=scale*ptr2d(i,j)+add_offset
2621 myfmin(2)=min(myfmin(2),fval)
2622 myfmax(2)=max(myfmax(2),fval)
2623 importfields%msk(j,i)=fval
2624 END DO
2625 END DO
2626!
2627! Surface roughness length scale (m), from WAV.
2628!
2629 CASE ('zo', 'Zo')
2630 DO i=iminp,imaxp
2631 DO j=jminp,jmaxp
2632 myfmin(1)=min(myfmin(1),ptr2d(i,j))
2633 myfmax(1)=max(myfmax(1),ptr2d(i,j))
2634 fval=scale*ptr2d(i,j)+add_offset
2635 myfmin(2)=min(myfmin(2),fval)
2636 myfmax(2)=max(myfmax(2),fval)
2637 importfields%zo(j,i)=fval
2638 END DO
2639 END DO
2640!
2641! Friction velocity (m s-1), from WAV.
2642!
2643 CASE ('ustar')
2644 DO i=iminp,imaxp
2645 DO j=jminp,jmaxp
2646 myfmin(1)=min(myfmin(1),ptr2d(i,j))
2647 myfmax(1)=max(myfmax(1),ptr2d(i,j))
2648 fval=scale*ptr2d(i,j)+add_offset
2649 myfmin(2)=min(myfmin(2),fval)
2650 myfmax(2)=max(myfmax(2),fval)
2651 importfields%ustar(j,i)=fval
2652 END DO
2653 END DO
2654!
2655! Import field not found.
2656!
2657 CASE DEFAULT
2658 IF (localpet.eq.0) THEN
2659 WRITE (cplout,10) trim(importnamelist(ifld)), &
2660 & trim(cinpname)
2661 END IF
2662 rc=esmf_rc_not_found
2663 IF (esmf_logfounderror(rctocheck=rc, &
2664 & msg=esmf_logerr_passthru, &
2665 & line=__line__, &
2666 & file=myfile)) THEN
2667 RETURN
2668 END IF
2669 END SELECT
2670!
2671! Nullify pointer to make sure that it does not point on a random
2672! part in the memory.
2673!
2674 IF (associated(ptr2d)) nullify (ptr2d)
2675 END DO de_loop
2676!
2677! Get import field minimun and maximum values.
2678!
2679 CALL esmf_vmallreduce (vm, &
2680 & senddata=myfmin, &
2681 & recvdata=fmin, &
2682 & count=2, &
2683 & reduceflag=esmf_reduce_min, &
2684 & rc=rc)
2685 IF (esmf_logfounderror(rctocheck=rc, &
2686 & msg=esmf_logerr_passthru, &
2687 & line=__line__, &
2688 & file=myfile)) THEN
2689 RETURN
2690 END IF
2691!
2692 CALL esmf_vmallreduce (vm, &
2693 & senddata=myfmax, &
2694 & recvdata=fmax, &
2695 & count=2, &
2696 & reduceflag=esmf_reduce_max, &
2697 & rc=rc)
2698 IF (esmf_logfounderror(rctocheck=rc, &
2699 & msg=esmf_logerr_passthru, &
2700 & line=__line__, &
2701 & file=myfile)) THEN
2702 RETURN
2703 END IF
2704!
2705! Report import field information.
2706!
2707 IF ((debuglevel.ge.0).and.(localpet.eq.0)) THEN
2708 WRITE (cplout,20) trim(importnamelist(ifld)), &
2709 & trim(time_currentstring), ng, &
2710 & fmin(1), fmax(1)
2711 IF (scale.ne.1.0_dp) THEN
2712 WRITE (cplout,30) fmin(2), fmax(2), &
2713 & ' regcmScale = ', scale
2714 ELSE IF (add_offset.ne.0.0_dp) THEN
2715 WRITE (cplout,30) fmin(2), fmax(2), &
2716 & ' AddOffset = ', add_offset
2717 END IF
2718 END IF
2719!
2720! Debugging: write out import field into a NetCDF file.
2721!
2722 IF ((debuglevel.ge.3).and. &
2723 & models(iatmos)%ImportField(id)%debug_write) THEN
2724 WRITE (ofile,40) ng, trim(importnamelist(ifld)), &
2725 & year, month, day, hour, minutes, seconds
2726 CALL esmf_fieldwrite (field, &
2727 & trim(ofile), &
2728 & overwrite=.true., &
2729 & rc=rc)
2730 IF (esmf_logfounderror(rctocheck=rc, &
2731 & msg=esmf_logerr_passthru, &
2732 & line=__line__, &
2733 & file=myfile)) THEN
2734 RETURN
2735 END IF
2736 END IF
2737 END DO fld_loop
2738!
2739! Deallocate local arrays.
2740!
2741 IF (allocated(importnamelist)) deallocate (importnamelist)
2742!
2743! Update RegCM import calls counter.
2744!
2745 IF (importcount.gt.0) THEN
2746 models(iatmos)%ImportCalls=models(iatmos)%ImportCalls+1
2747 END IF
2748!
2749 IF (esm_track) THEN
2750 WRITE (trac,'(a,a,i0)') '<== Exiting RegCM_Import', &
2751 & ', PET', petrank
2752 FLUSH (trac)
2753 END IF
2754 IF (debuglevel.gt.0) FLUSH (cplout)
2755!
2756 10 FORMAT (/,2x,'RegCM_Import - unable to find option to import: ', &
2757 & a,/,18x,'check ''Import(atmos)'' in input script: ', a)
2758 20 FORMAT (2x,'RegCM_Import - ESMF: importing field ''',a,'''', &
2759 & t72,a,2x,'Grid ',i2.2,/, &
2760 & 19x,'(Dmin = ', 1p,e15.8,0p,' Dmax = ',1p,e15.8,0p,')')
2761 30 FORMAT (19x,'(Cmin = ', 1p,e15.8,0p,' Cmax = ',1p,e15.8,0p, &
2762 & a,1p,e15.8,0p,')')
2763 40 FORMAT ('regcm_',i2.2,'_import_',a,'_',i4.4,2('-',i2.2),'_', &
2764 & i2.2,2('.',i2.2),'.nc')
2765
2766 RETURN

References mod_esmf_esm::cinpname, mod_esmf_esm::cplout, mod_esmf_esm::debuglevel, mod_esmf_esm::esm_track, mod_esmf_esm::field_index(), mod_esmf_esm::iatmos, mod_esmf_esm::models, mod_esmf_esm::petrank, and mod_esmf_esm::trac.

Referenced by regcm_modeladvance().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ regcm_modeladvance()

subroutine, private esmf_regcm_mod::regcm_modeladvance ( type (esmf_gridcomp) model,
integer, intent(out) rc )
private

Definition at line 2079 of file esmf_atm_regcm.h.

2080!
2081!=======================================================================
2082! !
2083! Advance RegCM component for a coupling interval (seconds) using !
2084! "RegCM_run". It also calls "RegCM_Import" and "RegCM_Export" to !
2085! import and export coupling fields, respectively. !
2086! !
2087!=======================================================================
2088!
2089 USE mod_runparams, ONLY : ifrest, ktau, dtsrf
2090!
2091! Imported variable declarations.
2092!
2093 integer, intent(out) :: rc
2094
2095 TYPE (ESMF_GridComp) :: model
2096!
2097! Local variable declarations.
2098!
2099 integer :: is, ng
2100 integer :: localPET, PETcount, phase
2101!
2102 real (dp) :: tstr, tend
2103!
2104 character (len=*), parameter :: MyFile = &
2105 & __FILE__//", RegCM_SetModelAdvance"
2106!
2107 character (ESMF_MAXSTR) :: str1, str2
2108!
2109 TYPE (ESMF_Clock) :: clock
2110 TYPE (ESMF_State) :: ExportState, ImportState
2111 TYPE (ESMF_TimeInterval) :: TimeFrom, TimeTo, TimeStep
2112 TYPE (ESMF_Time) :: ReferenceTime
2113 TYPE (ESMF_Time) :: CurrentTime, StartTime, StopTime
2114 TYPE (ESMF_VM) :: vm
2115!
2116!-----------------------------------------------------------------------
2117! Initialize return code flag to success state (no error).
2118!-----------------------------------------------------------------------
2119!
2120 IF (esm_track) THEN
2121 WRITE (trac,'(a,a,i0)') '==> Entering RegCM_ModelAdvance', &
2122 & ', PET', petrank
2123 FLUSH (trac)
2124 END IF
2125 rc=esmf_success
2126!
2127!-----------------------------------------------------------------------
2128! Get information about the gridded component.
2129!-----------------------------------------------------------------------
2130!
2131! Inquire about RegCM component.
2132!
2133 CALL esmf_gridcompget (model, &
2134 & importstate=importstate, &
2135 & exportstate=exportstate, &
2136 & clock=clock, &
2137 & localpet=localpet, &
2138 & petcount=petcount, &
2139 & currentphase=phase, &
2140 & vm=vm, &
2141 & rc=rc)
2142 IF (esmf_logfounderror(rctocheck=rc, &
2143 & msg=esmf_logerr_passthru, &
2144 & line=__line__, &
2145 & file=myfile)) THEN
2146 RETURN
2147 END IF
2148!
2149! Get driver time step interval, start time, stopping time, reference
2150! time, and current time.
2151!
2152 CALL esmf_clockget (clock, &
2153 & timestep=timestep, &
2154 & starttime=starttime, &
2155 & stoptime=stoptime, &
2156 & reftime=referencetime, &
2157 & currtime=clockinfo(iatmos)%CurrentTime, &
2158 & rc=rc)
2159 IF (esmf_logfounderror(rctocheck=rc, &
2160 & msg=esmf_logerr_passthru, &
2161 & line=__line__, &
2162 & file=myfile)) THEN
2163 RETURN
2164 END IF
2165!
2166!-----------------------------------------------------------------------
2167! Calculate run time.
2168!-----------------------------------------------------------------------
2169!
2170 IF (clockinfo(idriver)%Restarted) THEN
2171 timefrom=currenttime-clockinfo(iatmos)%RestartTime
2172 ELSE
2173 timefrom=currenttime-clockinfo(iatmos)%StartTime
2174 END IF
2175!
2176 CALL esmf_timeintervalget (timefrom, &
2177 & s_r8=tstr, &
2178 & rc=rc)
2179 IF (esmf_logfounderror(rctocheck=rc, &
2180 & msg=esmf_logerr_passthru, &
2181 & line=__line__, &
2182 & file=myfile)) THEN
2183 RETURN
2184 END IF
2185!
2186 timeto=timefrom+timestep
2187 CALL esmf_timeintervalget (timeto, &
2188 & s_r8=tend, &
2189 & rc=rc)
2190 IF (esmf_logfounderror(rctocheck=rc, &
2191 & msg=esmf_logerr_passthru, &
2192 & line=__line__, &
2193 & file=myfile)) THEN
2194 RETURN
2195 END IF
2196!
2197 IF (clockinfo(idriver)%Restarted.and. &
2198 & (starttime.eq.currenttime)) THEN
2199 tstr=tstr+dtsrf
2200 END IF
2201!
2202!-----------------------------------------------------------------------
2203! Debugging: report time information string (YYYY-MM-DD hh:mm:ss).
2204!-----------------------------------------------------------------------
2205!
2206 IF ((debuglevel.ge.0).and.(localpet.eq.0)) THEN
2207!
2208! Current driver time.
2209!
2210 CALL esmf_timeget (currenttime, &
2211 & timestringisofrac=str1, &
2212 & rc=rc)
2213 IF (esmf_logfounderror(rctocheck=rc, &
2214 & msg=esmf_logerr_passthru, &
2215 & line=__line__, &
2216 & file=myfile)) THEN
2217 RETURN
2218 END IF
2219!
2220! Next driver coupling time.
2221!
2222 CALL esmf_timeget (currenttime+timestep, &
2223 & timestringisofrac=str2, &
2224 & rc=rc)
2225 IF (esmf_logfounderror(rctocheck=rc, &
2226 & msg=esmf_logerr_passthru, &
2227 & line=__line__, &
2228 & file=myfile)) THEN
2229 RETURN
2230 END IF
2231!
2232 IF (debuglevel.eq.0) THEN
2233 WRITE (cplout,10) trim(str1), trim(str2), phase
2234 ELSE
2235 WRITE (cplout,20) trim(str1), trim(str2), phase, tstr, tend
2236 END IF
2237 END IF
2238!
2239!-----------------------------------------------------------------------
2240! Get import fields.
2241!-----------------------------------------------------------------------
2242!
2243 IF ((currenttime.ne.reftime).or.restarted) THEN
2244 IF (nimport(iatmos).gt.0) THEN
2245 DO ng=1,models(iatmos)%Ngrids
2246 IF (any(coupled(iatmos)%LinkedGrid(ng,:))) THEN
2247 CALL regcm_import (ng, model, rc=rc)
2248 IF (esmf_logfounderror(rctocheck=rc, &
2249 & msg=esmf_logerr_passthru, &
2250 & line=__line__, &
2251 & file=myfile)) THEN
2252 RETURN
2253 END IF
2254 END IF
2255 END DO
2256 END IF
2257 END IF
2258!
2259!-----------------------------------------------------------------------
2260! Run RegCM component. Notice that atmosphere component is advanced
2261! when ng=1. In nested application, its numerical kernel will advance
2262! all the nested grids in their logical order.
2263!-----------------------------------------------------------------------
2264!
2265 CALL regcm_run (tstr, tend)
2266!
2267!-----------------------------------------------------------------------
2268! Put export fields.
2269!-----------------------------------------------------------------------
2270!
2271 IF (nexport(iatmos).gt.0) THEN
2272 DO ng=1,models(iatmos)%Ngrids
2273 IF (any(coupled(iatmos)%LinkedGrid(ng,:))) THEN
2274 CALL regcm_export (ng, model, rc)
2275 IF (esmf_logfounderror(rctocheck=rc, &
2276 & msg=esmf_logerr_passthru, &
2277 & line=__line__, &
2278 & file=myfile)) THEN
2279 RETURN
2280 END IF
2281 END IF
2282 END DO
2283 END IF
2284!
2285 IF (esm_track) THEN
2286 WRITE (trac,'(a,a,i0)') '<== Exiting RegCM_ModelAdvance', &
2287 & ', PET', petrank
2288 FLUSH (trac)
2289 END IF
2290!
2291 10 FORMAT (' Running RegCM Component: ',a,' --> ',a,' Phase: ',i1)
2292 20 FORMAT (' Running RegCM Component: ',a,' --> ',a,' Phase: ',i1, &
2293 & ' [',f15.2, '-', f15.2, ']')
2294
2295 RETURN

References mod_esmf_esm::clockinfo, mod_esmf_esm::coupled, mod_esmf_esm::cplout, mod_esmf_esm::debuglevel, mod_esmf_esm::esm_track, mod_esmf_esm::iatmos, mod_esmf_esm::idriver, mod_esmf_esm::models, mod_esmf_esm::nexport, mod_esmf_esm::nimport, mod_esmf_esm::petrank, regcm_export(), regcm_import(), mod_esmf_esm::timestep, and mod_esmf_esm::trac.

Referenced by atm_setservices().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ regcm_setclock()

subroutine, private esmf_regcm_mod::regcm_setclock ( type (esmf_gridcomp) model,
integer, intent(out) rc )
private

Definition at line 700 of file esmf_atm_regcm.h.

701!
702!=======================================================================
703! !
704! Sets RegCM component date calendar, start and stop time, and !
705! coupling interval. !
706! !
707!=======================================================================
708!
709 USE mod_dynparam , ONLY : calendar
710 USE mod_runparams, ONLY : idate0, idate1, idate2, dtsec
711 USE mod_date, ONLY : split_idate
712!
713! Imported variable declarations.
714!
715 integer, intent(out) :: rc
716!
717 TYPE (ESMF_GridComp) :: model
718!
719! Local variable declarations.
720!
721 integer :: ng
722 integer :: ref_year, start_year, stop_year
723 integer :: ref_month, start_month, stop_month
724 integer :: ref_day, start_day, stop_day
725 integer :: ref_hour, start_hour, stop_hour
726 integer :: ref_minute, start_minute, stop_minute
727 integer :: ref_second, start_second, stop_second
728 integer :: TimeFrac
729!
730 character (len= 20) :: Calendar
731 character (len=160) :: message
732
733 character (len=*), parameter :: MyFile = &
734 & __FILE__//", RegCM_SetClock"
735!
736 TYPE (ESMF_CalKind_Flag) :: CalType
737 TYPE (ESMF_Time) :: StartTime
738!
739!-----------------------------------------------------------------------
740! Initialize return code flag to success state (no error).
741!-----------------------------------------------------------------------
742!
743 IF (esm_track) THEN
744 WRITE (trac,'(a,a,i0)') '==> Entering RegCM_SetClock', &
745 & ', PET', petrank
746 FLUSH (trac)
747 END IF
748 rc=esmf_success
749!
750!-----------------------------------------------------------------------
751! Create RegCM component clock.
752!-----------------------------------------------------------------------
753!
754 calendar=trim(clockinfo(iatmos)%CalendarString)
755 IF (trim(calendar).eq.'gregorian') THEN
756 caltype=esmf_calkind_gregorian
757 ELSE IF ((calendar.eq.'noleap').or.(calendar.eq.'365_day')) THEN
758 caltype=esmf_calkind_noleap
759 ELSE IF (calendar.eq.'360_day') THEN
760 caltype=esmf_calkind_360day
761 ELSE
762 caltype=esmf_calkind_gregorian
763 END IF
764!
765 clockinfo(iatmos)%Calendar=esmf_calendarcreate(caltype, &
766 & name=trim(calendar), &
767 & rc=rc)
768 IF (esmf_logfounderror(rctocheck=rc, &
769 & msg=esmf_logerr_passthru, &
770 & line=__line__, &
771 & file=myfile)) THEN
772 RETURN
773 END IF
774!
775! Set reference time.
776!
777 CALL split_idate (idate0, ref_year, ref_month, ref_day, ref_hour)
778 ref_minute=0
779 ref_second=0
780!
781 CALL esmf_timeset (clockinfo(iatmos)%ReferenceTime, &
782 & yy=ref_year, &
783 & mm=ref_month, &
784 & dd=ref_day, &
785 & h =ref_hour, &
786 & m =ref_minute, &
787 & s =ref_second, &
788 & calendar=clockinfo(iatmos)%Calendar, &
789 & rc=rc)
790 IF (esmf_logfounderror(rctocheck=rc, &
791 & msg=esmf_logerr_passthru, &
792 & line=__line__, &
793 & file=myfile)) THEN
794 RETURN
795 END IF
796!
797! Set start time.
798!
799 CALL split_idate (idate1, start_year, start_month, start_day, &
800 & start_hour)
801 start_minute=0
802 start_second=0
803!
804 CALL esmf_timeset (clockinfo(iatmos)%StartTime, &
805 & yy=start_year, &
806 & mm=start_month, &
807 & dd=start_day, &
808 & h =start_hour, &
809 & m =start_minute, &
810 & s =start_second, &
811 & calendar=clockinfo(iatmos)%Calendar, &
812 & rc=rc)
813 IF (esmf_logfounderror(rctocheck=rc, &
814 & msg=esmf_logerr_passthru, &
815 & line=__line__, &
816 & file=myfile)) THEN
817 RETURN
818 END IF
819!
820! Set stop time.
821!
822 CALL split_idate (idate2, stop_year, stop_month, stop_day, &
823 & stop_hour)
824 stop_minute=0
825 stop_second=0
826!
827 CALL esmf_timeset (clockinfo(iatmos)%StopTime, &
828 & yy=stop_year, &
829 & mm=stop_month, &
830 & dd=stop_day, &
831 & h =stop_hour, &
832 & m =stop_minute, &
833 & s =stop_second, &
834 & calendar=clockinfo(iatmos)%Calendar, &
835 & rc=rc)
836 IF (esmf_logfounderror(rctocheck=rc, &
837 & msg=esmf_logerr_passthru, &
838 & line=__line__, &
839 & file=myfile)) THEN
840 RETURN
841 END IF
842!
843!-----------------------------------------------------------------------
844! Get component clock.
845!-----------------------------------------------------------------------
846!
847 CALL esmf_gridcompget (model, &
848 & clock=clockinfo(iatmos)%Clock, &
849 & rc=rc)
850 IF (esmf_logfounderror(rctocheck=rc, &
851 & msg=esmf_logerr_passthru, &
852 & line=__line__, &
853 & file=myfile)) THEN
854 RETURN
855 END IF
856!
857 CALL esmf_clockget (clockinfo(iatmos)%Clock, &
858 & timestep=clockinfo(iatmos)%TimeStep, &
859 & currtime=clockinfo(iatmos)%CurrentTime, &
860 & rc=rc)
861 IF (esmf_logfounderror(rctocheck=rc, &
862 & msg=esmf_logerr_passthru, &
863 & line=__line__, &
864 & file=myfile)) THEN
865 RETURN
866 END IF
867!
868!-----------------------------------------------------------------------
869! Compare driver time against RegCM component time.
870!-----------------------------------------------------------------------
871!
872 IF (clockinfo(idriver)%Restarted) THEN
873 starttime=clockinfo(idriver)%RestartTime
874 ELSE
875 starttime=clockinfo(idriver)%StartTime
876 END IF
877!
878 IF (clockinfo(iatmos)%StartTime.ne.starttime) THEN
879 CALL esmf_timeprint (clockinfo(iatmos)%StartTime, &
880 & options="string", &
881 & rc=rc)
882 IF (esmf_logfounderror(rctocheck=rc, &
883 & msg=esmf_logerr_passthru, &
884 & line=__line__, &
885 & file=myfile)) THEN
886 RETURN
887 END IF
888!
889 CALL esmf_timeprint (starttime, &
890 & options="string", &
891 & rc=rc)
892 IF (esmf_logfounderror(rctocheck=rc, &
893 & msg=esmf_logerr_passthru, &
894 & line=__line__, &
895 & file=myfile)) THEN
896 RETURN
897 END IF
898!
899 message='Driver and RegCM start times do not match: '// &
900 & 'please check the config files.'
901 CALL esmf_logseterror (esmf_failure, rctoreturn=rc, &
902 & msg=trim(message))
903 RETURN
904 END IF
905!
906 IF (clockinfo(iatmos )%StopTime.ne. &
907 & clockinfo(idriver)%StopTime) THEN
908 CALL esmf_timeprint (clockinfo(iatmos)%StopTime, &
909 & options="string", &
910 & rc=rc)
911 IF (esmf_logfounderror(rctocheck=rc, &
912 & msg=esmf_logerr_passthru, &
913 & line=__line__, &
914 & file=myfile)) THEN
915 RETURN
916 END IF
917!
918 CALL esmf_timeprint (clockinfo(idriver)%StopTime, &
919 & options="string", &
920 & rc=rc)
921 IF (esmf_logfounderror(rctocheck=rc, &
922 & msg=esmf_logerr_passthru, &
923 & line=__line__, &
924 & file=myfile)) THEN
925 RETURN
926 END IF
927!
928 message='Driver and RegCM stop times do not match: '// &
929 & 'please check the config files.'
930 CALL esmf_logseterror (esmf_failure, rctoreturn=rc, &
931 & msg=trim(message))
932 RETURN
933 END IF
934!
935 IF (clockinfo(iatmos )%Calendar.ne. &
936 & clockinfo(idriver)%Calendar) THEN
937 CALL esmf_calendarprint (clockinfo(iatmos)%Calendar, &
938 & options="calkindflag", &
939 & rc=rc)
940 IF (esmf_logfounderror(rctocheck=rc, &
941 & msg=esmf_logerr_passthru, &
942 & line=__line__, &
943 & file=myfile)) THEN
944 RETURN
945 END IF
946!
947 CALL esmf_calendarprint (clockinfo(idriver)%Calendar, &
948 & options="calkindflag", &
949 & rc=rc)
950 IF (esmf_logfounderror(rctocheck=rc, &
951 & msg=esmf_logerr_passthru, &
952 & line=__line__, &
953 & file=myfile)) THEN
954 RETURN
955 END IF
956!
957 message='Driver and RegCM calendars do not match: '// &
958 & 'please check the config files.'
959 CALL esmf_logseterror (esmf_failure, rctoreturn=rc, &
960 & msg=trim(message))
961 RETURN
962 END IF
963!
964!-----------------------------------------------------------------------
965! Modify component clock time step.
966!-----------------------------------------------------------------------
967!
968 timefrac=0
969 DO ng=1,models(iatmos)%Ngrids
970 IF (any(coupled(iatmos)%LinkedGrid(ng,:))) THEN
971 timefrac=max(timefrac, &
972 & maxval(models(iatmos)%TimeFrac(ng,:), &
973 & mask=models(:)%IsActive))
974 END IF
975 END DO
976 IF (timefrac.lt.1) THEN ! needs to be 1 or greater
977 rc=esmf_rc_not_set ! cannot be 0
978 IF (esmf_logfounderror(rctocheck=rc, &
979 & msg=esmf_logerr_passthru, &
980 & line=__line__, &
981 & file=myfile)) THEN
982 RETURN
983 END IF
984 END IF
985 clockinfo(iatmos)%TimeStep=clockinfo(idriver)%TimeStep/timefrac
986!
987 clockinfo(iatmos)%Name='RegCM_clock'
988 CALL esmf_clockset (clockinfo(iatmos)%Clock, &
989 & name=trim(clockinfo(iatmos)%Name), &
990 & reftime =clockinfo(iatmos)%ReferenceTime, &
991 & timestep =clockinfo(iatmos)%TimeStep, &
992 & starttime=clockinfo(iatmos)%StartTime, &
993 & stoptime =clockinfo(iatmos)%StopTime, &
994 rc=rc)
995 IF (esmf_logfounderror(rctocheck=rc, &
996 & msg=esmf_logerr_passthru, &
997 & line=__line__, &
998 & file=myfile)) THEN
999 RETURN
1000 END IF
1001!
1002 IF (esm_track) THEN
1003 WRITE (trac,'(a,a,i0)') '<== Exiting RegCM_SetClock', &
1004 & ', PET', petrank
1005 FLUSH (trac)
1006 END IF
1007!
1008 RETURN

References mod_esmf_esm::clockinfo, mod_esmf_esm::coupled, mod_esmf_esm::esm_track, mod_esmf_esm::iatmos, mod_esmf_esm::idriver, mod_esmf_esm::models, mod_esmf_esm::petrank, mod_esmf_esm::timestep, and mod_esmf_esm::trac.

Referenced by atm_setservices().

Here is the caller graph for this function:

◆ regcm_setfinalize()

subroutine, private esmf_regcm_mod::regcm_setfinalize ( type (esmf_gridcomp) model,
type (esmf_state) importstate,
type (esmf_state) exportstate,
type (esmf_clock) clock,
integer, intent(out) rc )
private

Definition at line 2298 of file esmf_atm_regcm.h.

2301!
2302!=======================================================================
2303! !
2304! Finalize RegCM component execution. It calls RegCM_finalize. !
2305! !
2306!=======================================================================
2307!
2308! Imported variable declarations.
2309!
2310 integer, intent(out) :: rc
2311!
2312 TYPE (ESMF_Clock) :: clock
2313 TYPE (ESMF_GridComp) :: model
2314 TYPE (ESMF_State) :: ExportState
2315 TYPE (ESMF_State) :: ImportState
2316!
2317! Local variable declarations.
2318!
2319 character (len=*), parameter :: MyFile = &
2320 & __FILE__//", RegCM_SetFinalize"
2321!
2322!-----------------------------------------------------------------------
2323! Initialize return code flag to success state (no error).
2324!-----------------------------------------------------------------------
2325!
2326 IF (esm_track) THEN
2327 WRITE (trac,'(a,a,i0)') '==> Entering RegCM_SetFinalize', &
2328 & ', PET', petrank
2329 FLUSH (trac)
2330 END IF
2331 rc=esmf_success
2332!
2333!-----------------------------------------------------------------------
2334! Finalize RegCM component.
2335!-----------------------------------------------------------------------
2336!
2337 CALL regcm_finalize ()
2338 FLUSH (6) ! flush standard output buffer
2339!
2340 IF (esm_track) THEN
2341 WRITE (trac,'(a,a,i0)') '<== Exiting RegCM_SetFinalize', &
2342 & ', PET', petrank
2343 FLUSH (trac)
2344 END IF
2345!
2346 RETURN

References mod_esmf_esm::esm_track, mod_esmf_esm::petrank, and mod_esmf_esm::trac.

Referenced by atm_setservices().

Here is the caller graph for this function:

◆ regcm_setgridarrays()

subroutine, private esmf_regcm_mod::regcm_setgridarrays ( integer, intent(in) ng,
type (esmf_gridcomp), intent(inout) model,
integer, intent(out) rc )
private

Definition at line 1339 of file esmf_atm_regcm.h.

1340!
1341!=======================================================================
1342! !
1343! Sets RegCM component staggered, horizontal grids arrays, grid area, !
1344! and land/sea mask. !
1345! !
1346!=======================================================================
1347!
1348 USE mod_mppparam, ONLY : ma
1349 USE mod_runparams, ONLY : dxsq
1350 USE mod_atm_interface, ONLY : mddom
1351 USE mod_dynparam, ONLY : iy, jx, nproc, &
1352 & ide1, ide2, jde1, jde2, &
1353 & idi1, idi2, jdi1, jdi2, &
1354 & ice1, ice2, jce1, jce2
1355!
1356! Imported variable declarations.
1357!
1358 integer, intent(in) :: ng
1359 integer, intent(out) :: rc
1360!
1361 TYPE (ESMF_GridComp), intent(inout) :: model
1362!
1363! Local variable declarations.
1364!
1365 integer :: gtype, i, ivar, j
1366 integer :: localDE, localDEcount, localPET, PETcount
1367 integer :: cpus_per_dim(2)
1368!
1369 integer (i4b), pointer :: ptrM(:,:) => null()
1370!
1371 real (dp), pointer :: ptrA(:,:) => null()
1372 real (dp), pointer :: ptrX(:,:) => null()
1373 real (dp), pointer :: ptrY(:,:) => null()
1374!
1375 character (len=40) :: name
1376
1377 character (len=*), parameter :: MyFile = &
1378 & __FILE__//", RegCM_SetGridArrays"
1379!
1380 TYPE (ESMF_DistGrid) :: distGrid
1381 TYPE (ESMF_StaggerLoc) :: staggerLoc
1382 TYPE (ESMF_VM) :: vm
1383!
1384!-----------------------------------------------------------------------
1385! Initialize return code flag to success state (no error).
1386!-----------------------------------------------------------------------
1387!
1388 IF (esm_track) THEN
1389 WRITE (trac,'(a,a,i0)') '==> Entering RegCM_SetGridArrays', &
1390 & ', PET', petrank
1391 FLUSH (trac)
1392 END IF
1393 rc=esmf_success
1394!
1395!-----------------------------------------------------------------------
1396! Querry the Virtual Machine (VM) parallel environmemt for the MPI
1397! communicator handle and current node rank.
1398!-----------------------------------------------------------------------
1399!
1400 CALL esmf_gridcompget (model, &
1401 & localpet=localpet, &
1402 & petcount=petcount, &
1403 & vm=vm, &
1404 & rc=rc)
1405 IF (esmf_logfounderror(rctocheck=rc, &
1406 & msg=esmf_logerr_passthru, &
1407 & line=__line__, &
1408 & file=myfile)) THEN
1409 RETURN
1410 END IF
1411!
1412!-----------------------------------------------------------------------
1413! Calculate number of CPUs in each direction.
1414!-----------------------------------------------------------------------
1415!
1416 IF (nproc.lt.4) THEN
1417 cpus_per_dim(2)=1
1418 cpus_per_dim(1)=nproc
1419 ELSE IF (nproc.ge.4) THEN
1420 cpus_per_dim(2)=(nint(sqrt(dble(nproc)))/2)*2
1421 IF (iy.gt.int(1.5*dble(jx))) THEN
1422 cpus_per_dim(2)=cpus_per_dim(2)-1
1423 DO WHILE (mod(nproc,cpus_per_dim(2)).ne.0)
1424 cpus_per_dim(2)=cpus_per_dim(2)-1
1425 end do
1426 ELSE IF (jx.gt.int(1.5*dble(iy))) THEN
1427 cpus_per_dim(2)=cpus_per_dim(2)+1
1428 DO WHILE (mod(nproc,cpus_per_dim(2)).ne.0)
1429 cpus_per_dim(2)=cpus_per_dim(2)+1
1430 END DO
1431 ELSE
1432 DO WHILE (mod(nproc,cpus_per_dim(2)).ne.0)
1433 cpus_per_dim(2)=cpus_per_dim(2)+1
1434 END DO
1435 END IF
1436 cpus_per_dim(1)=nproc/cpus_per_dim(2)
1437 END IF
1438!
1439!-----------------------------------------------------------------------
1440! Create DistGrid based on model domain decomposition.
1441!
1442! ESMF is basically using a right handed coordinate system, and
1443! using the Fortran way of using the smallest stride to the first
1444! dimension but RegCM not. The order of dimension is reversed
1445! because of this limitation.
1446!-----------------------------------------------------------------------
1447!
1448 distgrid=esmf_distgridcreate(minindex=(/ 1, 1 /), &
1449 & maxindex=(/ iy, jx /), &
1450 & regdecomp=cpus_per_dim, &
1451 & rc=rc)
1452 IF (esmf_logfounderror(rctocheck=rc, &
1453 & msg=esmf_logerr_passthru, &
1454 & line=__line__, &
1455 & file=myfile)) THEN
1456 RETURN
1457 END IF
1458!
1459!-----------------------------------------------------------------------
1460! Set component grid coordinates.
1461!-----------------------------------------------------------------------
1462!
1463! Define component grid location type: Arakawa B-grid.
1464!
1465 IF (.not.allocated(models(iatmos)%mesh)) THEN
1466 allocate ( models(iatmos)%mesh(2) )
1467 models(iatmos)%mesh(1)%gtype=icenter
1468 models(iatmos)%mesh(2)%gtype=icorner
1469 END IF
1470!
1471! Create ESMF Grid.
1472!
1473 models(iatmos)%grid(ng)=esmf_gridcreate(distgrid=distgrid, &
1474 & indexflag=esmf_index_global, &
1475 & name=trim(models(iatmos)%name), &
1476 & rc=rc)
1477 IF (esmf_logfounderror(rctocheck=rc, &
1478 & msg=esmf_logerr_passthru, &
1479 & line=__line__, &
1480 & file=myfile)) THEN
1481 RETURN
1482 END IF
1483!
1484! Get number of local decomposition elements (DEs). Usually, a single
1485! DE is associated with each Persistent Execution Thread (PETs). Thus,
1486! localDEcount=1.
1487!
1488 CALL esmf_gridget (models(iatmos)%grid(ng), &
1489 & localdecount=localdecount, &
1490 & rc=rc)
1491 IF (esmf_logfounderror(rctocheck=rc, &
1492 & msg=esmf_logerr_passthru, &
1493 & line=__line__, &
1494 & file=myfile)) THEN
1495 RETURN
1496 END IF
1497!
1498! Mesh coordinates for each variable type.
1499!
1500 mesh_loop : DO ivar=1,ubound(models(iatmos)%mesh, dim=1)
1501!
1502! Set staggering type, Arakawa B-grid.
1503!
1504 SELECT CASE (models(iatmos)%mesh(ivar)%gtype)
1505 CASE (icenter)
1506 staggerloc=esmf_staggerloc_center
1507 CASE (icorner)
1508 staggerloc=esmf_staggerloc_corner
1509 END SELECT
1510!
1511! Allocate coordinate storage associated with staggered grid type.
1512! No coordinate values are set yet.
1513!
1514 CALL esmf_gridaddcoord (models(iatmos)%grid(ng), &
1515 & staggerloc=staggerloc, &
1516 & rc=rc)
1517 IF (esmf_logfounderror(rctocheck=rc, &
1518 & msg=esmf_logerr_passthru, &
1519 & line=__line__, &
1520 & file=myfile)) THEN
1521 RETURN
1522 END IF
1523!
1524! Allocate storage for masking.
1525!
1526 CALL esmf_gridadditem (models(iatmos)%grid(ng), &
1527 & staggerloc=staggerloc, &
1528 & itemflag=esmf_griditem_mask, &
1529 & rc=rc)
1530 IF (esmf_logfounderror(rctocheck=rc, &
1531 & msg=esmf_logerr_passthru, &
1532 & line=__line__, &
1533 & file=myfile)) THEN
1534 RETURN
1535 END IF
1536!
1537! The RegCM masking is as follows, 0: ocean
1538! 2: land
1539!
1540 models(iatmos)%LandValue=2
1541 models(iatmos)%SeaValue=0
1542!
1543! Allocate storage for grid area.
1544!
1545 CALL esmf_gridadditem (models(iatmos)%grid(ng), &
1546 & staggerloc=staggerloc, &
1547 & itemflag=esmf_griditem_area, &
1548 & rc=rc)
1549 IF (esmf_logfounderror(rctocheck=rc, &
1550 & msg=esmf_logerr_passthru, &
1551 & line=__line__, &
1552 & file=myfile)) THEN
1553 RETURN
1554 END IF
1555!
1556! Get pointers and set coordinates for the grid. Usually, the DO-loop
1557! is executed once since localDEcount=1.
1558!
1559 de_loop : DO localde=0,localdecount-1
1560 CALL esmf_gridgetcoord (models(iatmos)%grid(ng), &
1561 & coorddim=1, &
1562 & staggerloc=staggerloc, &
1563 & localde=localde, &
1564 & farrayptr=ptrx, &
1565 & rc=rc)
1566 IF (esmf_logfounderror(rctocheck=rc, &
1567 & msg=esmf_logerr_passthru, &
1568 & line=__line__, &
1569 & file=myfile)) THEN
1570 RETURN
1571 END IF
1572!
1573 CALL esmf_gridgetcoord (models(iatmos)%grid(ng), &
1574 & coorddim=2, &
1575 & staggerloc=staggerloc, &
1576 & localde=localde, &
1577 & farrayptr=ptry, &
1578 & rc=rc)
1579 IF (esmf_logfounderror(rctocheck=rc, &
1580 & msg=esmf_logerr_passthru, &
1581 & line=__line__, &
1582 & file=myfile)) THEN
1583 RETURN
1584 END IF
1585!
1586 CALL esmf_gridgetitem (models(iatmos)%grid(ng), &
1587 & itemflag=esmf_griditem_mask, &
1588 & staggerloc=staggerloc, &
1589 & localde=localde, &
1590 & farrayptr=ptrm, &
1591 & rc=rc)
1592 IF (esmf_logfounderror(rctocheck=rc, &
1593 & msg=esmf_logerr_passthru, &
1594 & line=__line__, &
1595 & file=myfile)) THEN
1596 RETURN
1597 END IF
1598!
1599 CALL esmf_gridgetitem (models(iatmos)%grid(ng), &
1600 & itemflag=esmf_griditem_area, &
1601 & staggerloc=staggerloc, &
1602 & localde=localde, &
1603 & farrayptr=ptra, &
1604 & rc=rc)
1605 IF (esmf_logfounderror(rctocheck=rc, &
1606 & msg=esmf_logerr_passthru, &
1607 & line=__line__, &
1608 & file=myfile)) THEN
1609 RETURN
1610 END IF
1611!
1612! Fill grid pointers.
1613!
1614 SELECT CASE (models(iatmos)%mesh(ivar)%gtype)
1615 CASE (icorner)
1616 DO i=ide1,ide2
1617 DO j=jde1,jde2
1618 ptrx(i,j)=mddom%dlon(j,i)
1619 ptry(i,j)=mddom%dlat(j,i)
1620 END DO
1621 END DO
1622 ptra=dxsq
1623!
1624 IF (ma%has_bdyright) THEN
1625 ptrx(:,jde2+1)=ptrx(:,jde2)+ &
1626 & (ptrx(:,jde2)-ptrx(:,jde2-1))
1627 ptry(:,jde2+1)=ptry(:,jde2)+ &
1628 & (ptry(:,jde2)-ptry(:,jde2-1))
1629 END IF
1630!
1631 IF (ma%has_bdytop) THEN
1632 ptrx(ide2+1,:)=ptrx(ide2,:)+ &
1633 & (ptrx(ide2,:)-ptrx(ide2-1,:))
1634 ptry(ide2+1,:)=ptry(ide2,:)+ &
1635 & (ptry(ide2,:)-ptry(ide2-1,:))
1636 END IF
1637 CASE (icenter)
1638 DO i=ice1,ice2
1639 DO j=jce1,jce2
1640 ptrx(i,j)=mddom%xlon(j,i)
1641 ptry(i,j)=mddom%xlat(j,i)
1642 ptrm(i,j)=int(mddom%mask(j,i))
1643 END DO
1644 END DO
1645 ptra=dxsq
1646!
1647 IF (ma%has_bdyright) THEN
1648 ptrx(:,jce2+1)=ptrx(:,jce2)+ &
1649 & (ptrx(:,jce2)-ptrx(:,jce2-1))
1650 ptry(:,jce2+1)=ptry(:,jce2)+ &
1651 & (ptry(:,jce2)-ptry(:,jce2-1))
1652 END IF
1653!
1654 IF (ma%has_bdytop) THEN
1655 ptrx(ice2+1,:)=ptrx(ice2,:)+ &
1656 & (ptrx(ice2,:)-ptrx(ice2-1,:))
1657 ptry(ice2+1,:)=ptry(ice2,:)+ &
1658 & (ptry(ice2,:)-ptry(ice2-1,:))
1659 END IF
1660 END SELECT
1661!
1662! Nullify pointers.
1663!
1664 IF ( associated(ptrx) ) nullify (ptrx)
1665 IF ( associated(ptry) ) nullify (ptry)
1666 IF ( associated(ptrm) ) nullify (ptrm)
1667 IF ( associated(ptra) ) nullify (ptra)
1668 END DO de_loop
1669!
1670! Debugging: write out component grid in VTK format.
1671!
1672 IF (debuglevel.ge.4) THEN
1673 gtype=models(iatmos)%mesh(ivar)%gtype
1674 CALL esmf_gridwritevtk (models(iatmos)%grid(ng), &
1675 & filename="regcm_"// &
1676 & trim(gridtype(gtype))// &
1677 & "_point", &
1678 & staggerloc=staggerloc, &
1679 & rc=rc)
1680 IF (esmf_logfounderror(rctocheck=rc, &
1681 & msg=esmf_logerr_passthru, &
1682 & line=__line__, &
1683 & file=myfile)) THEN
1684 RETURN
1685 END IF
1686 END IF
1687 END DO mesh_loop
1688!
1689! Assign grid to gridded component.
1690!
1691 CALL esmf_gridcompset (model, &
1692 & grid=models(iatmos)%grid(ng), &
1693 & rc=rc)
1694 IF (esmf_logfounderror(rctocheck=rc, &
1695 & msg=esmf_logerr_passthru, &
1696 & line=__line__, &
1697 & file=myfile)) THEN
1698 RETURN
1699 END IF
1700!
1701 IF (esm_track) THEN
1702 WRITE (trac,'(a,a,i0)') '<== Exiting RegCM_SetGridArrays', &
1703 & ', PET', petrank
1704 FLUSH (trac)
1705 END IF
1706!
1707 RETURN

References mod_esmf_esm::debuglevel, mod_esmf_esm::esm_track, mod_esmf_esm::gridtype, mod_esmf_esm::iatmos, mod_esmf_esm::icenter, mod_esmf_esm::icorner, mod_esmf_esm::models, mod_esmf_esm::petrank, and mod_esmf_esm::trac.

Referenced by regcm_setinitializep2().

Here is the caller graph for this function:

◆ regcm_setinitializep1()

subroutine, private esmf_regcm_mod::regcm_setinitializep1 ( type (esmf_gridcomp) model,
type (esmf_state) importstate,
type (esmf_state) exportstate,
type (esmf_clock) clock,
integer, intent(out) rc )
private

Definition at line 300 of file esmf_atm_regcm.h.

303!
304!=======================================================================
305! !
306! RegCM component Phase 1 initialization: sets import and export !
307! fields long and short names into its respective state. !
308! !
309!=======================================================================
310!
311! Imported variable declarations.
312!
313 integer, intent(out) :: rc
314!
315 TYPE (ESMF_GridComp) :: model
316 TYPE (ESMF_State) :: ImportState
317 TYPE (ESMF_State) :: ExportState
318 TYPE (ESMF_Clock) :: clock
319!
320! Local variable declarations.
321!
322 integer :: i, ng
323!
324 character (len=100) :: CoupledSet, StateLabel
325 character (len=240) :: StandardName, ShortName
326
327 character (len=*), parameter :: MyFile = &
328 & __FILE__//", RegCM_SetInitializeP1"
329!
330!-----------------------------------------------------------------------
331! Initialize return code flag to success state (no error).
332!-----------------------------------------------------------------------
333!
334 IF (esm_track) THEN
335 WRITE (trac,'(a,a,i0)') '==> Entering RegCM_SetInitializeP1', &
336 & ', PET', petrank
337 FLUSH (trac)
338 END IF
339 rc=esmf_success
340!
341!-----------------------------------------------------------------------
342! Set RegCM import state and fields.
343!-----------------------------------------------------------------------
344!
345 importing : IF (nimport(iatmos).gt.0) THEN
346 DO ng=1,models(iatmos)%Ngrids
347 IF (any(coupled(iatmos)%LinkedGrid(ng,:))) THEN
348 coupledset=trim(coupled(iatmos)%SetLabel(ng))
349 statelabel=trim(coupled(iatmos)%ImpLabel(ng))
350 CALL nuopc_addnestedstate (importstate, &
351 & cplset=trim(coupledset), &
352 & nestedstatename=trim(statelabel),&
353 & nestedstate=models(iatmos)% &
354 & importstate(ng), &
355 rc=rc)
356 IF (esmf_logfounderror(rctocheck=rc, &
357 & msg=esmf_logerr_passthru, &
358 & line=__line__, &
359 & file=myfile)) THEN
360 RETURN
361 END IF
362!
363! Add fields import state.
364!
365 DO i=1,nimport(iatmos)
366 standardname=models(iatmos)%ImportField(i)%standard_name
367 shortname =models(iatmos)%ImportField(i)%short_name
368 CALL nuopc_advertise (models(iatmos)%ImportState(ng), &
369 & standardname=trim(standardname), &
370 & name=trim(shortname), &
371 & rc=rc)
372 IF (esmf_logfounderror(rctocheck=rc, &
373 & msg=esmf_logerr_passthru, &
374 & line=__line__, &
375 & file=myfile)) THEN
376 RETURN
377 END IF
378 END DO
379 END IF
380 END DO
381 END IF importing
382!
383!-----------------------------------------------------------------------
384! Set RegCM export state and fields.
385!-----------------------------------------------------------------------
386!
387 exporting : IF (nexport(iatmos).gt.0) THEN
388 DO ng=1,models(iatmos)%Ngrids
389 IF (any(coupled(iatmos)%LinkedGrid(ng,:))) THEN
390 coupledset=trim(coupled(iatmos)%SetLabel(ng))
391 statelabel=trim(coupled(iatmos)%ExpLabel(ng))
392 CALL nuopc_addnestedstate (exportstate, &
393 & cplset=trim(coupledset), &
394 & nestedstatename=trim(statelabel),&
395 & nestedstate=models(iatmos)% &
396 & exportstate(ng), &
397 rc=rc)
398 IF (esmf_logfounderror(rctocheck=rc, &
399 & msg=esmf_logerr_passthru, &
400 & line=__line__, &
401 & file=myfile)) THEN
402 RETURN
403 END IF
404!
405! Add fields to export state.
406!
407 DO i=1,nexport(iatmos)
408 standardname=models(iatmos)%ExportField(i)%standard_name
409 shortname =models(iatmos)%ExportField(i)%short_name
410 CALL nuopc_advertise (models(iatmos)%ExportState(ng), &
411 & standardname=trim(standardname), &
412 & name=trim(shortname), &
413 & rc=rc)
414 IF (esmf_logfounderror(rctocheck=rc, &
415 & msg=esmf_logerr_passthru, &
416 & line=__line__, &
417 & file=myfile)) THEN
418 RETURN
419 END IF
420 END DO
421 END IF
422 END DO
423 END IF exporting
424!
425 IF (esm_track) THEN
426 WRITE (trac,'(a,a,i0)') '<== Exiting RegCM_SetInitializeP1', &
427 & ', PET', petrank
428 FLUSH (trac)
429 END IF
430!
431 RETURN

References mod_esmf_esm::coupled, mod_esmf_esm::esm_track, mod_esmf_esm::iatmos, mod_esmf_esm::models, mod_esmf_esm::nexport, mod_esmf_esm::nimport, mod_esmf_esm::petrank, and mod_esmf_esm::trac.

◆ regcm_setinitializep2()

subroutine, private esmf_regcm_mod::regcm_setinitializep2 ( type (esmf_gridcomp) model,
type (esmf_state) importstate,
type (esmf_state) exportstate,
type (esmf_clock) clock,
integer, intent(out) rc )
private

Definition at line 434 of file esmf_atm_regcm.h.

437!
438!=======================================================================
439! !
440! RegCM component Phase 2 initialization: Initializes RegCM, sets !
441! component grid, and adds import and export fields to respective !
442! states. !
443! !
444!=======================================================================
445!
446 USE mod_runparams, ONLY : dtsrf
447 USE mod_update, ONLY : regcm_allocate => rcm_allocate
448!
449! Imported variable declarations.
450!
451 integer, intent(out) :: rc
452!
453 TYPE (ESMF_GridComp) :: model
454 TYPE (ESMF_State) :: ImportState
455 TYPE (ESMF_State) :: ExportState
456 TYPE (ESMF_Clock) :: clock
457!
458! Local variable declarations.
459!
460 integer :: ng
461 integer :: MyComm, localPET, PETcount
462!
463 character (len=*), parameter :: MyFile = &
464 & __FILE__//", RegCM_SetInitializeP2"
465!
466 TYPE (ESMF_Time) :: CurrentTime, StartTime
467 TYPE (ESMF_VM) :: vm
468!
469!-----------------------------------------------------------------------
470! Initialize return code flag to success state (no error).
471!-----------------------------------------------------------------------
472!
473 IF (esm_track) THEN
474 WRITE (trac,'(a,a,i0)') '==> Entering RegCM_SetInitializeP2', &
475 & ', PET', petrank
476 FLUSH (trac)
477 END IF
478 rc=esmf_success
479!
480!-----------------------------------------------------------------------
481! Querry the Virtual Machine (VM) parallel environmemt for the MPI
482! communicator handle and current node rank.
483!-----------------------------------------------------------------------
484!
485 CALL esmf_gridcompget (model, &
486 & vm=vm, &
487 & rc=rc)
488 IF (esmf_logfounderror(rctocheck=rc, &
489 & msg=esmf_logerr_passthru, &
490 & line=__line__, &
491 & file=myfile)) THEN
492 RETURN
493 END IF
494!
495 CALL esmf_vmget (vm, &
496 & localpet=localpet, &
497 & petcount=petcount, &
498 & mpicommunicator=mycomm, &
499 & rc=rc)
500 IF (esmf_logfounderror(rctocheck=rc, &
501 & msg=esmf_logerr_passthru, &
502 & line=__line__, &
503 & file=myfile)) THEN
504 RETURN
505 END IF
506!
507!-----------------------------------------------------------------------
508! Initialize RegCM component.
509!-----------------------------------------------------------------------
510!
511 CALL regcm_initialize (mpicommunicator=mycomm)
512 CALL regcm_allocate ()
513!
514!-----------------------------------------------------------------------
515! Set-up grid and load coordinate data.
516!-----------------------------------------------------------------------
517!
518 DO ng=1,models(iatmos)%Ngrids
519 IF (any(coupled(iatmos)%LinkedGrid(ng,:))) THEN
520 CALL regcm_setgridarrays (ng, model, rc)
521 IF (esmf_logfounderror(rctocheck=rc, &
522 & msg=esmf_logerr_passthru, &
523 & line=__line__, &
524 & file=myfile)) THEN
525 RETURN
526 END IF
527 END IF
528 END DO
529!
530!-----------------------------------------------------------------------
531! Set-up fields and register to import/export states.
532!-----------------------------------------------------------------------
533!
534 DO ng=1,models(iatmos)%Ngrids
535 IF (any(coupled(iatmos)%LinkedGrid(ng,:))) THEN
536 CALL regcm_setstates (ng, model, rc)
537 IF (esmf_logfounderror(rctocheck=rc, &
538 & msg=esmf_logerr_passthru, &
539 & line=__line__, &
540 & file=myfile)) THEN
541 RETURN
542 END IF
543 END IF
544 END DO
545!
546 IF (esm_track) THEN
547 WRITE (trac,'(a,a,i0)') '<== Exiting RegCM_SetInitializeP2', &
548 & ', PET', petrank
549 FLUSH (trac)
550 END IF
551!
552 RETURN

References mod_esmf_esm::coupled, mod_esmf_esm::esm_track, mod_esmf_esm::iatmos, mod_esmf_esm::models, mod_esmf_esm::petrank, regcm_setgridarrays(), regcm_setstates(), and mod_esmf_esm::trac.

Here is the call graph for this function:

◆ regcm_setrunclock()

subroutine, private esmf_regcm_mod::regcm_setrunclock ( type (esmf_gridcomp) model,
integer, intent(out) rc )
private

Definition at line 1013 of file esmf_atm_regcm.h.

1014!
1015!=======================================================================
1016! !
1017! Sets RegCM run clock manually to avoid getting zero time stamps at !
1018! the first regridding call. !
1019! !
1020!=======================================================================
1021!
1022! Imported variable declarations.
1023!
1024 integer, intent(out) :: rc
1025!
1026 TYPE (ESMF_GridComp) :: model
1027!
1028! Local variable declarations.
1029!
1030 character (len=*), parameter :: MyFile = &
1031 & __FILE__//", COAMPS_SetRunClock"
1032!
1033 TYPE (ESMF_Clock) :: driverClock, modelClock
1034 TYPE (ESMF_Time) :: currTime
1035!
1036!-----------------------------------------------------------------------
1037! Initialize return code flag to success state (no error).
1038!-----------------------------------------------------------------------
1039!
1040 IF (esm_track) THEN
1041 WRITE (trac,'(a,a,i0)') '==> Entering RegCM_SetRunClock', &
1042 & ', PET', petrank
1043 FLUSH (trac)
1044 END IF
1045 rc=esmf_success
1046!
1047!-----------------------------------------------------------------------
1048! Set ROMS run clock manually.
1049!-----------------------------------------------------------------------
1050!
1051! Inquire driver and model clock.
1052!
1053 CALL nuopc_modelget (model, &
1054 & driverclock=driverclock, &
1055 & modelclock=modelclock, &
1056 & rc=rc)
1057 IF (esmf_logfounderror(rctocheck=rc, &
1058 & msg=esmf_logerr_passthru, &
1059 & line=__line__, &
1060 & file=myfile)) THEN
1061 RETURN
1062 END IF
1063!
1064! Set model clock to have the current start time as the driver clock.
1065!
1066 CALL esmf_clockget (driverclock, &
1067 & currtime=currtime, &
1068 & rc=rc)
1069 IF (esmf_logfounderror(rctocheck=rc, &
1070 & msg=esmf_logerr_passthru, &
1071 & line=__line__, &
1072 & file=myfile)) THEN
1073 RETURN
1074 END IF
1075!
1076 CALL esmf_clockset (modelclock, &
1077 & currtime=currtime, &
1078 & rc=rc)
1079 IF (esmf_logfounderror(rctocheck=rc, &
1080 & msg=esmf_logerr_passthru, &
1081 & line=__line__, &
1082 & file=myfile)) THEN
1083 RETURN
1084 END IF
1085!
1086! Check and set the component clock against the driver clock.
1087!
1088 CALL nuopc_compchecksetclock (model, &
1089 & driverclock, &
1090 & rc=rc)
1091 IF (esmf_logfounderror(rctocheck=rc, &
1092 & msg=esmf_logerr_passthru, &
1093 & line=__line__, &
1094 & file=myfile)) THEN
1095 RETURN
1096 END IF
1097!
1098 IF (esm_track) THEN
1099 WRITE (trac,'(a,a,i0)') '<== Exiting RegCM_SetRunClock', &
1100 & ', PET', petrank
1101 FLUSH (trac)
1102 END IF
1103!
1104 RETURN

References mod_esmf_esm::esm_track, mod_esmf_esm::petrank, and mod_esmf_esm::trac.

Referenced by atm_setservices().

Here is the caller graph for this function:

◆ regcm_setstates()

subroutine, private esmf_regcm_mod::regcm_setstates ( integer, intent(in) ng,
type (esmf_gridcomp) model,
integer, intent(out) rc )
private

Definition at line 1710 of file esmf_atm_regcm.h.

1711!
1712!=======================================================================
1713! !
1714! Adds RegCM component export and import fields into its respective !
1715! state. !
1716! !
1717!=======================================================================
1718!
1719! Imported variable declarations.
1720!
1721 integer, intent(in) :: ng
1722 integer, intent(out) :: rc
1723!
1724 TYPE (ESMF_GridComp) :: model
1725!
1726! Local variable declarations.
1727!
1728 integer :: i, id
1729 integer :: localDE, localDEcount
1730 integer :: localPET, PETcount
1731 integer :: ExportCount, ImportCount
1732!
1733 real (dp), dimension(:,:), pointer :: ptr2d => null()
1734!
1735 character (len=*), parameter :: MyFile = &
1736 & __FILE__//", RegCM_SetStates"
1737!
1738 character (ESMF_MAXSTR), allocatable :: ExportNameList(:)
1739 character (ESMF_MAXSTR), allocatable :: ImportNameList(:)
1740!
1741 TYPE (ESMF_ArraySpec) :: arraySpec2d
1742 TYPE (ESMF_Field) :: field
1743 TYPE (ESMF_StaggerLoc) :: staggerLoc
1744 TYPE (ESMF_VM) :: vm
1745!
1746!-----------------------------------------------------------------------
1747! Initialize return code flag to success state (no error).
1748!-----------------------------------------------------------------------
1749!
1750 IF (esm_track) THEN
1751 WRITE (trac,'(a,a,i0)') '==> Entering RegCM_SetStates', &
1752 & ', PET', petrank
1753 FLUSH (trac)
1754 END IF
1755 rc=esmf_success
1756!
1757!-----------------------------------------------------------------------
1758! Get gridded component information.
1759!-----------------------------------------------------------------------
1760!
1761! Get import and export states.
1762!
1763 CALL esmf_gridcompget (model, &
1764 & localpet=localpet, &
1765 & petcount=petcount, &
1766 & vm=vm, &
1767 & rc=rc)
1768 IF (esmf_logfounderror(rctocheck=rc, &
1769 & msg=esmf_logerr_passthru, &
1770 & line=__line__, &
1771 & file=myfile)) THEN
1772 RETURN
1773 END IF
1774!
1775! Get number of local decomposition elements (DEs). Usually, a single
1776! Decomposition Element (DE) is associated with each Persistent
1777! Execution Thread (PETs). Thus, localDEcount=1.
1778!
1779 CALL esmf_gridget (models(iatmos)%grid(ng), &
1780 & localdecount=localdecount, &
1781 & rc=rc)
1782 IF (esmf_logfounderror(rctocheck=rc, &
1783 & msg=esmf_logerr_passthru, &
1784 & line=__line__, &
1785 & file=myfile)) THEN
1786 RETURN
1787 END IF
1788!
1789!-----------------------------------------------------------------------
1790! Set a 2D floating-point array descriptor.
1791!-----------------------------------------------------------------------
1792!
1793 CALL esmf_arrayspecset (arrayspec2d, &
1794 & typekind=esmf_typekind_r8, &
1795 & rank=2, &
1796 & rc=rc)
1797 IF (esmf_logfounderror(rctocheck=rc, &
1798 & msg=esmf_logerr_passthru, &
1799 & line=__line__, &
1800 & file=myfile)) THEN
1801 RETURN
1802 END IF
1803!
1804!-----------------------------------------------------------------------
1805! Add export fields into export state.
1806!-----------------------------------------------------------------------
1807!
1808 exporting : IF (nexport(iatmos).gt.0) THEN
1809!
1810! Get number of fields to export.
1811!
1812 CALL esmf_stateget (models(iatmos)%ExportState(ng), &
1813 & itemcount=exportcount, &
1814 & rc=rc)
1815 IF (esmf_logfounderror(rctocheck=rc, &
1816 & msg=esmf_logerr_passthru, &
1817 & line=__line__, &
1818 & file=myfile)) THEN
1819 RETURN
1820 END IF
1821!
1822! Get a list of export fields names.
1823!
1824 IF (.not.allocated(exportnamelist)) THEN
1825 allocate ( exportnamelist(exportcount) )
1826 END IF
1827 CALL esmf_stateget (models(iatmos)%ExportState(ng), &
1828 & itemnamelist=exportnamelist, &
1829 & rc=rc)
1830 IF (esmf_logfounderror(rctocheck=rc, &
1831 & msg=esmf_logerr_passthru, &
1832 & line=__line__, &
1833 & file=myfile)) THEN
1834 RETURN
1835 END IF
1836!
1837! Set export field(s).
1838!
1839 DO i=1,exportcount
1840 id=field_index(models(iatmos)%ExportField, exportnamelist(i))
1841!
1842 IF (nuopc_isconnected(models(iatmos)%ExportState(ng), &
1843 & fieldname=trim(exportnamelist(i)), &
1844 & rc=rc)) THEN
1845!
1846! Set staggering type.
1847!
1848 SELECT CASE (models(iatmos)%ExportField(id)%gtype)
1849 CASE (icenter)
1850 staggerloc=esmf_staggerloc_center
1851 CASE (icorner)
1852 staggerloc=esmf_staggerloc_corner
1853 END SELECT
1854!
1855! Create 2D field from the Grid and arraySpec.
1856!
1857 field=esmf_fieldcreate(models(iatmos)%grid(ng), &
1858 & arrayspec2d, &
1859 & staggerloc=staggerloc, &
1860 & name=trim(exportnamelist(i)), &
1861 & rc=rc)
1862 IF (esmf_logfounderror(rctocheck=rc, &
1863 & msg=esmf_logerr_passthru, &
1864 & line=__line__, &
1865 & file=myfile)) THEN
1866 RETURN
1867 END IF
1868!
1869! Put data into state. Usually, the DO-loop is executed once since
1870! localDEcount=1.
1871!
1872 DO localde=0,localdecount-1
1873!
1874! Get pointer to DE-local memory allocation within field.
1875!
1876 CALL esmf_fieldget (field, &
1877 & localde=localde, &
1878 & farrayptr=ptr2d, &
1879 & rc=rc)
1880 IF (esmf_logfounderror(rctocheck=rc, &
1881 & msg=esmf_logerr_passthru, &
1882 & line=__line__, &
1883 & file=myfile)) THEN
1884 RETURN
1885 END IF
1886!
1887! Initialize pointer.
1888!
1889 ptr2d=missing_dp
1890!
1891! Nullify pointer to make sure that it does not point on a random part
1892! in the memory.
1893!
1894 IF ( associated(ptr2d) ) nullify (ptr2d)
1895 END DO
1896!
1897! Add field export state.
1898!
1899 CALL nuopc_realize (models(iatmos)%ExportState(ng), &
1900 & field=field, &
1901 & rc=rc)
1902 IF (esmf_logfounderror(rctocheck=rc, &
1903 & msg=esmf_logerr_passthru, &
1904 & line=__line__, &
1905 & file=myfile)) THEN
1906 RETURN
1907 END IF
1908!
1909! Remove field from export state because it is not connected.
1910!
1911 ELSE
1912 IF (localpet.eq.0) THEN
1913 WRITE (cplout,10) trim(exportnamelist(i)), &
1914 & 'Export State: ', &
1915 & trim(coupled(iatmos)%ExpLabel(ng))
1916 END IF
1917 CALL esmf_stateremove (models(iatmos)%ExportState(ng), &
1918 & (/ trim(exportnamelist(i)) /), &
1919 & rc=rc)
1920 IF (esmf_logfounderror(rctocheck=rc, &
1921 & msg=esmf_logerr_passthru, &
1922 & line=__line__, &
1923 & file=myfile)) THEN
1924 RETURN
1925 END IF
1926 END IF
1927 END DO
1928!
1929! Deallocate arrays.
1930!
1931 IF ( allocated(exportnamelist) ) deallocate (exportnamelist)
1932!
1933 END IF exporting
1934!
1935!-----------------------------------------------------------------------
1936! Add import fields into import state.
1937!-----------------------------------------------------------------------
1938!
1939 importing : IF (nimport(iatmos).gt.0) THEN
1940!
1941! Get number of fields to import.
1942!
1943 CALL esmf_stateget (models(iatmos)%ImportState(ng), &
1944 & itemcount=importcount, &
1945 & rc=rc)
1946 IF (esmf_logfounderror(rctocheck=rc, &
1947 & msg=esmf_logerr_passthru, &
1948 & line=__line__, &
1949 & file=myfile)) THEN
1950 RETURN
1951 END IF
1952!
1953! Get a list of import fields names.
1954!
1955 IF (.not.allocated(importnamelist)) THEN
1956 allocate (importnamelist(importcount))
1957 END IF
1958 CALL esmf_stateget (models(iatmos)%ImportState(ng), &
1959 & itemnamelist=importnamelist, &
1960 & rc=rc)
1961 IF (esmf_logfounderror(rctocheck=rc, &
1962 & msg=esmf_logerr_passthru, &
1963 & line=__line__, &
1964 & file=myfile)) THEN
1965 RETURN
1966 END IF
1967!
1968! Set import field(s).
1969!
1970 DO i=1,importcount
1971 id=field_index(models(iatmos)%ImportField, importnamelist(i))
1972!
1973 IF (nuopc_isconnected(models(iatmos)%ImportState(ng), &
1974 & fieldname=trim(importnamelist(i)), &
1975 & rc=rc)) THEN
1976!
1977! Set staggering type.
1978!
1979 SELECT CASE (models(iatmos)%ImportField(id)%gtype)
1980 CASE (icenter)
1981 staggerloc=esmf_staggerloc_center
1982 CASE (icorner)
1983 staggerloc=esmf_staggerloc_corner
1984 END SELECT
1985!
1986! Create 2D field from the Grid, arraySpec.
1987!
1988 field=esmf_fieldcreate(models(iatmos)%grid(ng), &
1989 & arrayspec2d, &
1990 & staggerloc=staggerloc, &
1991 & name=trim(importnamelist(i)), &
1992 & rc=rc)
1993 IF (esmf_logfounderror(rctocheck=rc, &
1994 & msg=esmf_logerr_passthru, &
1995 & line=__line__, &
1996 & file=myfile)) THEN
1997 RETURN
1998 END IF
1999!
2000! Put data into state. Usually, the DO-loop is executed once since
2001! localDEcount=1.
2002!
2003 DO localde=0,localdecount-1
2004!
2005! Get pointer to DE-local memory allocation within field.
2006!
2007 CALL esmf_fieldget (field, &
2008 & localde=localde, &
2009 & farrayptr=ptr2d, &
2010 & rc=rc)
2011 IF (esmf_logfounderror(rctocheck=rc, &
2012 & msg=esmf_logerr_passthru, &
2013 & line=__line__, &
2014 & file=myfile)) THEN
2015 RETURN
2016 END IF
2017!
2018! Initialize pointer.
2019!
2020 ptr2d=missing_dp
2021!
2022! Nullify pointer to make sure that it does not point on a random
2023! part in the memory.
2024!
2025 IF (associated(ptr2d)) nullify (ptr2d)
2026 END DO
2027!
2028! Add field import state.
2029!
2030 CALL nuopc_realize (models(iatmos)%ImportState(ng), &
2031 & field=field, &
2032 & rc=rc)
2033 IF (esmf_logfounderror(rctocheck=rc, &
2034 & msg=esmf_logerr_passthru, &
2035 & line=__line__, &
2036 & file=myfile)) THEN
2037 RETURN
2038 END IF
2039!
2040! Remove field from import state because it is not connected.
2041!
2042 ELSE
2043 IF (localpet.eq.0) THEN
2044 WRITE (cplout,10) trim(importnamelist(i)), &
2045 & 'Import State: ', &
2046 & trim(coupled(iatmos)%ImpLabel(ng))
2047 END IF
2048 CALL esmf_stateremove (models(iatmos)%ImportState(ng), &
2049 & (/ trim(importnamelist(i)) /), &
2050 & rc=rc)
2051 IF (esmf_logfounderror(rctocheck=rc, &
2052 & msg=esmf_logerr_passthru, &
2053 & line=__line__, &
2054 & file=myfile)) THEN
2055 RETURN
2056 END IF
2057 END IF
2058 END DO
2059!
2060! Deallocate arrays.
2061!
2062 IF (allocated(importnamelist)) deallocate (importnamelist)
2063!
2064 END IF importing
2065!
2066 IF (esm_track) THEN
2067 WRITE (trac,'(a,a,i0)') '<== Exiting RegCM_SetStates', &
2068 & ', PET', petrank
2069 FLUSH (trac)
2070 END IF
2071!
2072!
2073 10 FORMAT ('RegCM_SetStates - Removing field ''',a,''' from ',a, &
2074 & '''',a,'''',/,18x,'because it is not connected.')
2075!
2076 RETURN

References mod_esmf_esm::coupled, mod_esmf_esm::cplout, mod_esmf_esm::esm_track, mod_esmf_esm::field_index(), mod_esmf_esm::iatmos, mod_esmf_esm::icenter, mod_esmf_esm::icorner, mod_esmf_esm::missing_dp, mod_esmf_esm::models, mod_esmf_esm::nexport, mod_esmf_esm::nimport, mod_esmf_esm::petrank, and mod_esmf_esm::trac.

Referenced by regcm_setinitializep2().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ regcm_uvrot()

subroutine, private esmf_regcm_mod::regcm_uvrot ( real (r8), dimension(jci1:jci2,ici1:ici2), intent(inout) u,
real (r8), dimension(jci1:jci2,ici1:ici2), intent(inout) v )
private

Definition at line 3411 of file esmf_atm_regcm.h.

3412!
3413!=======================================================================
3414! !
3415! Rotates RegCM vector components to a rectangular north-south and !
3416! east-west oriented grid. !
3417! !
3418!=======================================================================
3419!
3420 USE mod_constants, ONLY : degrad
3421 USE mod_atm_interface, ONLY : mddom
3422 USE mod_dynparam, ONLY : iproj
3423 USE mod_dynparam, ONLY : clon, clat, plon, plat, xcone
3424 USE mod_dynparam, ONLY : ici1, ici2, jci1, jci2
3425!
3426! Imported variable declarations.
3427!
3428 real (r8), intent(inout) :: u(jci1:jci2,ici1:ici2)
3429 real (r8), intent(inout) :: v(jci1:jci2,ici1:ici2)
3430!
3431! Local variable declarations.
3432!
3433 integer :: i, j
3434!
3435 real (r8) :: x, xs, xc, d, us, vs, sindel, cosdel
3436 real (r8) :: pollam, polphi, polcphi, polsphi
3437 real (r8) :: zarg1, zarg2, znorm, zphi, zrla, zrlap
3438!
3439 character (len=*), parameter :: MyFile = &
3440 & __FILE__//", RegCM_uvrot"
3441!
3442!-----------------------------------------------------------------------
3443! Rotate vector components to true EAST and true NORTH.
3444!-----------------------------------------------------------------------
3445!
3446 IF (esm_track) THEN
3447 WRITE (trac,'(a,a,i0)') '==> Entering RegCM_uvrot', &
3448 & ', PET', petrank
3449 FLUSH (trac)
3450 END IF
3451!
3452! Rotated Mercator (ROTMER) or Normal Mercator (NORMER) projections.
3453!
3454 IF ((iproj.eq.'ROTMER').or.(iproj.eq.'NORMER')) THEN
3455 IF (plat.gt.0.0_r8) THEN
3456 pollam=plon+180.0_r8
3457 polphi=90.0_r8-plat
3458 ELSE
3459 polphi=90.0_r8+plat
3460 pollam=plon
3461 END if
3462 IF (pollam.gt.180.0_r8) pollam=pollam-360.0_r8
3463!
3464 polcphi=dcos(degrad*polphi)
3465 polsphi=dsin(degrad*polphi)
3466!
3467 DO j=jci1,jci2
3468 DO i=ici1,ici2
3469 zphi=mddom%dlat(j,i)*degrad
3470 zrla=mddom%dlon(j,i)*degrad
3471 IF (mddom%dlat(j,i).gt.89.999999_r8) zrla=0.0_r8
3472 zrlap=pollam*degrad-zrla
3473 zarg1=polcphi*dsin(zrlap)
3474 zarg2=polsphi*cos(zphi)-polcphi*sin(zphi)*dcos(zrlap)
3475 znorm=1.0_r8/sqrt(zarg1**2+zarg2**2)
3476 sindel=zarg1*znorm
3477 cosdel=zarg2*znorm
3478!
3479 us= u(j,i)*cosdel+v(j,i)*sindel
3480 vs=-u(j,i)*sindel+v(j,i)*cosdel
3481 u(j,i)=us
3482 v(j,i)=vs
3483 END DO
3484 END DO
3485!
3486! Otherwise, rotate from Lambert conformal (LAMCON) projection.
3487!
3488 ELSE
3489 DO i=ici1,ici2
3490 DO j=jci1,jci2
3491 IF (((clon.ge.0.0_r8).and.(mddom%xlon(j,i).ge.0.0_r8)).or. &
3492 & ((clon.lt.0.0_r8).and.(mddom%xlon(j,i).lt.0.0_r8))) THEN
3493 x=(clon-mddom%xlon(j,i))*degrad*xcone
3494 ELSE
3495 IF (clon.ge.0.0_r8) THEN
3496 IF (abs(clon-(mddom%xlon(j,i)+360.0_r8)).lt. &
3497 & abs(clon-mddom%xlon(j,i))) THEN
3498 x=(clon-(mddom%xlon(j,i)+360.0d0))*degrad*xcone
3499 ELSE
3500 x=(clon-mddom%xlon(j,i))*degrad*xcone
3501 END IF
3502 ELSE
3503 IF (abs(clon-(mddom%xlon(j,i)-360.0_r8)).lt. &
3504 & abs(clon-mddom%xlon(j,i))) THEN
3505 x=(clon-(mddom%xlon(j,i)-360.0_r8))*degrad*xcone
3506 ELSE
3507 x=(clon-mddom%xlon(j,i))*degrad*xcone
3508 END IF
3509 END IF
3510 END IF
3511!
3512 xs=sin(x)
3513 xc=cos(x)
3514!
3515 IF (clat.ge.0.0_r8) THEN
3516 d=u(j,i)*xc-v(j,i)*xs
3517 v(j,i)=u(j,i)*xs+v(j,i)*xc
3518 u(j,i)=d
3519 ELSE
3520 d=u(j,i)*xc+v(j,i)*xs
3521 v(j,i)=v(j,i)*xc-u(j,i)*xs
3522 u(j,i)=d
3523 END IF
3524 END DO
3525 END DO
3526 END IF
3527!
3528 IF (esm_track) THEN
3529 WRITE (trac,'(a,a,i0)') '<== Exiting RegCM_uvrot', &
3530 & ', PET', petrank
3531 FLUSH (trac)
3532 END IF
3533!
3534 RETURN

References mod_esmf_esm::esm_track, mod_esmf_esm::petrank, and mod_esmf_esm::trac.

Referenced by regcm_export().

Here is the caller graph for this function: