ROMS
Loading...
Searching...
No Matches
lbc_mod::lbc_getatt Interface Reference

Public Member Functions

subroutine lbc_getatt_nf90 (ng, model, ncid, ncname, aname, s)
 
subroutine lbc_getatt_pio (ng, model, piofile, ncname, aname, s)
 

Detailed Description

Definition at line 28 of file lbc.F.

Member Function/Subroutine Documentation

◆ lbc_getatt_nf90()

subroutine lbc_mod::lbc_getatt::lbc_getatt_nf90 ( integer, intent(in) ng,
integer, intent(in) model,
integer, intent(in) ncid,
character (*), intent(in) ncname,
character (*), intent(in) aname,
type(t_lbc), dimension(4,nlbcvar,ngrids), intent(in) s )

Definition at line 45 of file lbc.F.

46!***********************************************************************
47! !
48! If restart, this routine reads lateral boundary conditions !
49! keywords strings from restart NetCDF file global attribute. !
50! Then, it checks their values with those provided from input !
51! script for consistency. !
52! !
53! On Input: !
54! !
55! ng Nested grid number (integer). !
56! model Calling model identifier (integer). !
57! ncid NetCDF file ID (integer). !
58! ncname NetCDF file name (character). !
59! aname NetCDF global attribute name (character). !
60! S Derived type structure, TYPE(T_LBC) !
61! !
62! On Output: !
63! !
64! exit_flag Error flag (integer) stored in MOD_SCALARS !
65! ioerror NetCDF return code (integer) stored in MOD_IOUNITS !
66! !
67!***********************************************************************
68!
69 USE mod_param
70 USE mod_parallel
71 USE mod_iounits
72 USE mod_ncparam
73 USE mod_netcdf
74 USE mod_scalars
75!
76#if !defined PARALLEL_IO && defined DISTRIBUTE
78#endif
79 USE strings_mod, ONLY : founderror
80!
81 implicit none
82!
83! Imported variable declarations.
84!
85 integer, intent(in) :: ng, model, ncid
86 character (*), intent(in) :: ncname
87 character (*), intent(in) :: aname
88!
89 TYPE(T_LBC), intent(in) :: S(4,nLBCvar,Ngrids)
90!
91! Local variable declarations
92!
93 integer :: i, ibry, ie, ifield, is, ne, lstr, lvar, status
94#if !defined PARALLEL_IO && defined DISTRIBUTE
95 integer, dimension(2) :: ibuffer
96#endif
97!
98 character (len= 7) :: string(4)
99 character (len= 8) :: B(4)
100 character (len= 40) :: BryVar1, BryVar2
101 character (len= 70) :: Bstring, line
102 character (len=2816) :: lbc_att
103
104 character (len=*), parameter :: MyFile = &
105 & __FILE__//", lbc_getatt_nf90"
106!
107!-----------------------------------------------------------------------
108! If restart, read and check lateral boundary conditions global
109! attribute.
110!-----------------------------------------------------------------------
111!
112! Read in global attribute.
113!
114 IF (inpthread) THEN
115 status=nf90_get_att(ncid, nf90_global, trim(aname), lbc_att)
116 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
117 WRITE (stdout,10) trim(aname), trim(ncname), &
118 & trim(sourcefile)
119 exit_flag=3
120 ioerror=status
121 END IF
122 END IF
123
124#if !defined PARALLEL_IO && defined DISTRIBUTE
125!
126! Broadcast error flags to all processors in the group.
127!
128 ibuffer(1)=exit_flag
129 ibuffer(2)=ioerror
130 CALL mp_bcasti (ng, model, ibuffer)
131 exit_flag=ibuffer(1)
132 ioerror=ibuffer(2)
133 IF (exit_flag.eq.noerror) THEN
134 CALL mp_bcasts (ng, model, lbc_att)
135 END IF
136#endif
137!
138! Check keyword values from global attribute and compare against
139! logical switches in lateral boundary condition structure.
140!
141 b(iwest )='western '
142 b(isouth)='southern'
143 b(ieast )='eastern '
144 b(inorth)='northern'
145 DO i=1,len(bstring)
146 bstring(i:i)=' '
147 END DO
148 DO i=1,len(line)
149 line(i:i)=' '
150 END DO
151!
152 DO ifield=1,nlbcvar
153 bryvar1=trim(vname(1,idbvar(ifield))) // ':'
154 is=index(lbc_att, trim(bryvar1))
155 IF (ifield.lt.nlbcvar) THEN
156 bryvar2=trim(vname(1,idbvar(ifield+1))) // ':'
157 ie=index(lbc_att, trim(bryvar2))-1
158 ELSE
159 ie=len_trim(lbc_att)
160 END IF
161 IF ((is.gt.0).and.(ie.gt.0).and.(ie.gt.is)) THEN
162 line=lbc_att(is:ie)
163 is=index(line, ':')+1 ! find colon location
164 ie=index(line, char(10))-1 ! find Line Feed location
165 IF (ie.le.0) THEN
166 ie=len_trim(line) ! last attribute entry
167 END IF
168 bstring=trim(adjustl(line(is:ie))) ! extract boundary string
169 ne=min(len_trim(bstring), 28) ! north keyword < 7 chars
170 string(1)=bstring( 1: 7) ! west keyword
171 string(2)=bstring( 8:14) ! south keyword
172 string(3)=bstring(15:21) ! east keyword
173 string(4)=bstring(22:ne) ! north keyword
174 DO ibry=1,4
175 SELECT CASE (trim(string(ibry)))
176 CASE ('Cha')
177 IF (.not.s(ibry,ifield,ng)%Chapman_implicit) THEN
178 IF (master) THEN
179 WRITE (stdout,20) b(ibry), &
180 & trim(vname(1,idbvar(ifield))), &
181 & trim(string(ibry)), &
182 & 'S(',ibry,ifield,ng,')%Chapman_implicit', &
183 & s(ibry,ifield,ng)%Chapman_implicit, &
184 & trim(ncname)
185 END IF
186 exit_flag=5
187 END IF
188 CASE ('Che')
189 IF (.not.s(ibry,ifield,ng)%Chapman_explicit) THEN
190 IF (master) THEN
191 WRITE (stdout,20) b(ibry), &
192 & trim(vname(1,idbvar(ifield))), &
193 & trim(string(ibry)), &
194 & 'S(',ibry,ifield,ng,')%Chapman_explicit', &
195 & s(ibry,ifield,ng)%Chapman_explicit, &
196 & trim(ncname)
197 END IF
198 exit_flag=5
199 END IF
200 CASE ('Cla')
201 IF (.not.s(ibry,ifield,ng)%clamped) THEN
202 IF (master) THEN
203 WRITE (stdout,20) b(ibry), &
204 & trim(vname(1,idbvar(ifield))), &
205 & trim(string(ibry)), &
206 & 'S(',ibry,ifield,ng,')%clamped', &
207 & s(ibry,ifield,ng)%clamped, &
208 & trim(ncname)
209 END IF
210 exit_flag=5
211 END IF
212 CASE ('Clo')
213 IF (.not.s(ibry,ifield,ng)%closed) THEN
214 IF (master) THEN
215 WRITE (stdout,20) b(ibry), &
216 & trim(vname(1,idbvar(ifield))), &
217 & trim(string(ibry)), &
218 & 'S(',ibry,ifield,ng,')%closed', &
219 & s(ibry,ifield,ng)%closed, &
220 & trim(ncname)
221 END IF
222 exit_flag=5
223 END IF
224 CASE ('Fla')
225 IF (.not.s(ibry,ifield,ng)%Flather) THEN
226 IF (master) THEN
227 WRITE (stdout,20) b(ibry), &
228 & trim(vname(1,idbvar(ifield))), &
229 & trim(string(ibry)), &
230 & 'S(',ibry,ifield,ng,')%Flather', &
231 & s(ibry,ifield,ng)%Flather, &
232 & trim(ncname)
233 END IF
234 exit_flag=5
235 END IF
236 CASE ('Gra')
237 IF (.not.s(ibry,ifield,ng)%gradient) THEN
238 IF (master) THEN
239 WRITE (stdout,20) b(ibry), &
240 & trim(vname(1,idbvar(ifield))), &
241 & trim(string(ibry)), &
242 & 'S(',ibry,ifield,ng,')%gradient', &
243 & s(ibry,ifield,ng)%gradient, &
244 & trim(ncname)
245 END IF
246 exit_flag=5
247 END IF
248 CASE ('Mix')
249 IF (.not.s(ibry,ifield,ng)%mixed) THEN
250 IF (master) THEN
251 WRITE (stdout,20) b(ibry), &
252 & trim(vname(1,idbvar(ifield))), &
253 & trim(string(ibry)), &
254 & 'S(',ibry,ifield,ng,')%mixed', &
255 & s(ibry,ifield,ng)%mixed, &
256 & trim(ncname)
257 END IF
258 exit_flag=5
259 END IF
260 CASE ('Nes')
261 IF (.not.s(ibry,ifield,ng)%nested) THEN
262 IF (master) THEN
263 WRITE (stdout,20) b(ibry), &
264 & trim(vname(1,idbvar(ifield))), &
265 & trim(string(ibry)), &
266 & 'S(',ibry,ifield,ng,')%nested', &
267 & s(ibry,ifield,ng)%nested, &
268 & trim(ncname)
269 END IF
270 exit_flag=5
271 END IF
272 CASE ('Per')
273 IF (.not.s(ibry,ifield,ng)%periodic) THEN
274 IF (master) THEN
275 WRITE (stdout,20) b(ibry), &
276 & trim(vname(1,idbvar(ifield))), &
277 & trim(string(ibry)), &
278 & 'S(',ibry,ifield,ng,')%periodic', &
279 & s(ibry,ifield,ng)%periodic, &
280 & trim(ncname)
281 END IF
282 exit_flag=5
283 END IF
284 CASE ('Rad')
285 IF (.not.s(ibry,ifield,ng)%radiation) THEN
286 IF (master) THEN
287 WRITE (stdout,20) b(ibry), &
288 & trim(vname(1,idbvar(ifield))), &
289 & trim(string(ibry)), &
290 & 'S(',ibry,ifield,ng,')%radiation', &
291 & s(ibry,ifield,ng)%radiation, &
292 & trim(ncname)
293 END IF
294 exit_flag=5
295 END IF
296 CASE ('RadNud')
297 IF (.not.(s(ibry,ifield,ng)%radiation.and. &
298 & s(ibry,ifield,ng)%nudging)) THEN
299 IF (master) THEN
300 WRITE (stdout,20) b(ibry), &
301 & trim(vname(1,idbvar(ifield))), &
302 & trim(string(ibry)), &
303 & 'S(',ibry,ifield,ng,')%radiation', &
304 & s(ibry,ifield,ng)%radiation, &
305 & trim(ncname)
306 END IF
307 exit_flag=5
308 END IF
309 CASE ('Red')
310 IF (.not.s(ibry,ifield,ng)%reduced) THEN
311 IF (master) THEN
312 WRITE (stdout,20) b(ibry), &
313 & trim(vname(1,idbvar(ifield))), &
314 & trim(string(ibry)), &
315 & 'S(',ibry,ifield,ng,')%reduced', &
316 & s(ibry,ifield,ng)%reduced, &
317 & trim(ncname)
318 END IF
319 exit_flag=5
320 END IF
321 CASE ('Shc')
322 IF (.not.s(ibry,ifield,ng)%Shchepetkin) THEN
323 IF (master) THEN
324 WRITE (stdout,20) b(ibry), &
325 & trim(vname(1,idbvar(ifield))), &
326 & trim(string(ibry)), &
327 & 'S(',ibry,ifield,ng,')%Shchepetkin', &
328 & s(ibry,ifield,ng)%Shchepetkin, &
329 & trim(ncname)
330 END IF
331 exit_flag=5
332 END IF
333 CASE DEFAULT
334 IF (master) THEN
335 WRITE (stdout,30) b(ibry), &
336 & trim(vname(1,idbvar(ifield))), &
337 & trim(string(ibry)), trim(ncname)
338 END IF
339 exit_flag=5
340 END SELECT
341 END DO
342 END IF
343 END DO
344!
345 10 FORMAT (/,' LBC_GETATT_NF90 - error while reading global ', &
346 & 'attribute:',2x,a,/,19x,'in restart file:',2x,a,/, &
347 & 19x,'call from:',2x,a, &
348 & /,19x,'Probably global attribute was not found ...', &
349 & /,19x,'restart file needs to be generated by ROMS ', &
350 & 'version 3.6 or higher', &
351 & /,19x,'Alternatively, you may use NO_LBC_ATT at your ', &
352 & 'own risk!')
353 20 FORMAT (/,' LBC_GETATT_NF90 - inconsistent ',a,' lateral', &
354 & 'boundary condition for variable: ',2x,a, &
355 & /,19x,'restart file LBC keyword = ',1x,a, &
356 & /,19x,'but assigned structure switch: ', &
357 & 1x,a,i1,',',i2,',',i1,a,' = ',l1, &
358 & /,19x,'check input script LBC keyword for consitency ...',&
359 & /,19x,'restart file:',2x,a)
360 30 FORMAT (/,' LBC_GETATT_NF90 - inconsistent ',a,' boundary for ', &
361 & 'variable: ',a,2x,'Keyword = ',a,/,19x,'in input file:', &
362 & 2x,a)
363!
364 RETURN
integer ioerror
integer stdout
character(len=256) sourcefile
integer, dimension(:), allocatable idbvar
character(len=maxlen), dimension(6, 0:nv) vname
logical inpthread
logical master
integer nlbcvar
Definition mod_param.F:355
integer, parameter iwest
integer exit_flag
integer, parameter isouth
integer, parameter ieast
integer, parameter inorth
integer noerror
logical function, public founderror(flag, noerr, line, routine)
Definition strings.F:52

References mod_scalars::exit_flag, strings_mod::founderror(), mod_ncparam::idbvar, mod_scalars::ieast, mod_scalars::inorth, mod_parallel::inpthread, mod_iounits::ioerror, mod_scalars::isouth, mod_scalars::iwest, mod_parallel::master, mod_scalars::noerror, mod_iounits::sourcefile, mod_iounits::stdout, and mod_ncparam::vname.

Here is the call graph for this function:

◆ lbc_getatt_pio()

subroutine lbc_mod::lbc_getatt::lbc_getatt_pio ( integer, intent(in) ng,
integer, intent(in) model,
type (file_desc_t), intent(in) piofile,
character (*), intent(in) ncname,
character (*), intent(in) aname,
type(t_lbc), dimension(4,nlbcvar,ngrids), intent(in) s )

Definition at line 370 of file lbc.F.

371!***********************************************************************
372! !
373! If restart, this routine reads lateral boundary conditions !
374! keywords strings from restart NetCDF file global attribute. !
375! Then, it checks their values with those provided from input !
376! script for consistency. !
377! !
378! On Input: !
379! !
380! ng Nested grid number (integer) !
381! model Calling model identifier (integer) !
382! pioFile PIO file descriptor structure, TYPE(File_desc_t) !
383! pioFile%fh file handler !
384! pioFile%iosystem IO system descriptor (struct) !
385! ncname PIO filename (character) !
386! aname PIO global attribute name (character) !
387! S Derived type structure, TYPE(T_LBC) !
388! !
389! On Output: !
390! !
391! exit_flag Error flag (integer) stored in MOD_SCALARS !
392! ioerror NetCDF return code (integer) stored in MOD_IOUNITS !
393! !
394!***********************************************************************
395!
396 USE mod_param
397 USE mod_parallel
398 USE mod_iounits
399 USE mod_ncparam
400 USE mod_scalars
401 USE pio
402!
403 USE strings_mod, ONLY : founderror
404!
405 implicit none
406!
407! Imported variable declarations.
408!
409 integer, intent(in) :: ng, model
410 character (*), intent(in) :: ncname
411 character (*), intent(in) :: aname
412!
413 TYPE (File_desc_t), intent(in) :: pioFile
414 TYPE(T_LBC), intent(in) :: S(4,nLBCvar,Ngrids)
415!
416! Local variable declarations
417!
418 integer :: i, ibry, ie, ifield, is, ne, lstr, lvar, status
419!
420 character (len= 7) :: string(4)
421 character (len= 8) :: B(4)
422 character (len= 40) :: BryVar1, BryVar2
423 character (len= 70) :: Bstring, line
424 character (len=2816) :: lbc_att
425
426 character (len=*), parameter :: MyFile = &
427 & __FILE__//", lbc_getatt_pio"
428!
429!-----------------------------------------------------------------------
430! If restart, read and check lateral boundary conditions global
431! attribute.
432!-----------------------------------------------------------------------
433!
434! Read in global attribute.
435!
436 status=pio_get_att(piofile, pio_global, trim(aname), lbc_att)
437 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
438 WRITE (stdout,10) trim(aname), trim(ncname), &
439 & trim(sourcefile)
440 exit_flag=3
441 ioerror=status
442 END IF
443!
444! Check keyword values from global attribute and compare against
445! logical switches in lateral boundary condition structure.
446!
447 b(iwest )='western '
448 b(isouth)='southern'
449 b(ieast )='eastern '
450 b(inorth)='northern'
451 DO i=1,len(bstring)
452 bstring(i:i)=' '
453 END DO
454 DO i=1,len(line)
455 line(i:i)=' '
456 END DO
457!
458 DO ifield=1,nlbcvar
459 bryvar1=trim(vname(1,idbvar(ifield))) // ':'
460 is=index(lbc_att, trim(bryvar1))
461 IF (ifield.lt.nlbcvar) THEN
462 bryvar2=trim(vname(1,idbvar(ifield+1))) // ':'
463 ie=index(lbc_att, trim(bryvar2))-1
464 ELSE
465 ie=len_trim(lbc_att)
466 END IF
467 IF ((is.gt.0).and.(ie.gt.0).and.(ie.gt.is)) THEN
468 line=lbc_att(is:ie)
469 is=index(line, ':')+1 ! find colon location
470 ie=index(line, char(10))-1 ! find Line Feed location
471 IF (ie.le.0) THEN
472 ie=len_trim(line) ! last attribute entry
473 END IF
474 bstring=trim(adjustl(line(is:ie))) ! extract boundary string
475 ne=min(len_trim(bstring), 28) ! north keyword < 7 chars
476 string(1)=bstring( 1: 7) ! west keyword
477 string(2)=bstring( 8:14) ! south keyword
478 string(3)=bstring(15:21) ! east keyword
479 string(4)=bstring(22:ne) ! north keyword
480 DO ibry=1,4
481 SELECT CASE (trim(string(ibry)))
482 CASE ('Cha')
483 IF (.not.s(ibry,ifield,ng)%Chapman_implicit) THEN
484 IF (master) THEN
485 WRITE (stdout,20) b(ibry), &
486 & trim(vname(1,idbvar(ifield))), &
487 & trim(string(ibry)), &
488 & 'S(',ibry,ifield,ng,')%Chapman_implicit', &
489 & s(ibry,ifield,ng)%Chapman_implicit, &
490 & trim(ncname)
491 END IF
492 exit_flag=5
493 END IF
494 CASE ('Che')
495 IF (.not.s(ibry,ifield,ng)%Chapman_explicit) THEN
496 IF (master) THEN
497 WRITE (stdout,20) b(ibry), &
498 & trim(vname(1,idbvar(ifield))), &
499 & trim(string(ibry)), &
500 & 'S(',ibry,ifield,ng,')%Chapman_explicit', &
501 & s(ibry,ifield,ng)%Chapman_explicit, &
502 & trim(ncname)
503 END IF
504 exit_flag=5
505 END IF
506
507 CASE ('Cla')
508 IF (.not.s(ibry,ifield,ng)%clamped) THEN
509 IF (master) THEN
510 WRITE (stdout,20) b(ibry), &
511 & trim(vname(1,idbvar(ifield))), &
512 & trim(string(ibry)), &
513 & 'S(',ibry,ifield,ng,')%clamped', &
514 & s(ibry,ifield,ng)%clamped, &
515 & trim(ncname)
516 END IF
517 exit_flag=5
518 END IF
519 CASE ('Clo')
520 IF (.not.s(ibry,ifield,ng)%closed) THEN
521 IF (master) THEN
522 WRITE (stdout,20) b(ibry), &
523 & trim(vname(1,idbvar(ifield))), &
524 & trim(string(ibry)), &
525 & 'S(',ibry,ifield,ng,')%closed', &
526 & s(ibry,ifield,ng)%closed, &
527 & trim(ncname)
528 END IF
529 exit_flag=5
530 END IF
531 CASE ('Fla')
532 IF (.not.s(ibry,ifield,ng)%Flather) THEN
533 IF (master) THEN
534 WRITE (stdout,20) b(ibry), &
535 & trim(vname(1,idbvar(ifield))), &
536 & trim(string(ibry)), &
537 & 'S(',ibry,ifield,ng,')%Flather', &
538 & s(ibry,ifield,ng)%Flather, &
539 & trim(ncname)
540 END IF
541 exit_flag=5
542 END IF
543 CASE ('Gra')
544 IF (.not.s(ibry,ifield,ng)%gradient) THEN
545 IF (master) THEN
546 WRITE (stdout,20) b(ibry), &
547 & trim(vname(1,idbvar(ifield))), &
548 & trim(string(ibry)), &
549 & 'S(',ibry,ifield,ng,')%gradient', &
550 & s(ibry,ifield,ng)%gradient, &
551 & trim(ncname)
552 END IF
553 exit_flag=5
554 END IF
555 CASE ('Mix')
556 IF (.not.s(ibry,ifield,ng)%mixed) THEN
557 IF (master) THEN
558 WRITE (stdout,20) b(ibry), &
559 & trim(vname(1,idbvar(ifield))), &
560 & trim(string(ibry)), &
561 & 'S(',ibry,ifield,ng,')%mixed', &
562 & s(ibry,ifield,ng)%mixed, &
563 & trim(ncname)
564 END IF
565 exit_flag=5
566 END IF
567 CASE ('Nes')
568 IF (.not.s(ibry,ifield,ng)%nested) THEN
569 IF (master) THEN
570 WRITE (stdout,20) b(ibry), &
571 & trim(vname(1,idbvar(ifield))), &
572 & trim(string(ibry)), &
573 & 'S(',ibry,ifield,ng,')%nested', &
574 & s(ibry,ifield,ng)%nested, &
575 & trim(ncname)
576 END IF
577 exit_flag=5
578 END IF
579 CASE ('Per')
580 IF (.not.s(ibry,ifield,ng)%periodic) THEN
581 IF (master) THEN
582 WRITE (stdout,20) b(ibry), &
583 & trim(vname(1,idbvar(ifield))), &
584 & trim(string(ibry)), &
585 & 'S(',ibry,ifield,ng,')%periodic', &
586 & s(ibry,ifield,ng)%periodic, &
587 & trim(ncname)
588 END IF
589 exit_flag=5
590 END IF
591 CASE ('Rad')
592 IF (.not.s(ibry,ifield,ng)%radiation) THEN
593 IF (master) THEN
594 WRITE (stdout,20) b(ibry), &
595 & trim(vname(1,idbvar(ifield))), &
596 & trim(string(ibry)), &
597 & 'S(',ibry,ifield,ng,')%radiation', &
598 & s(ibry,ifield,ng)%radiation, &
599 & trim(ncname)
600 END IF
601 exit_flag=5
602 END IF
603 CASE ('RadNud')
604 IF (.not.(s(ibry,ifield,ng)%radiation.and. &
605 & s(ibry,ifield,ng)%nudging)) THEN
606 IF (master) THEN
607 WRITE (stdout,20) b(ibry), &
608 & trim(vname(1,idbvar(ifield))), &
609 & trim(string(ibry)), &
610 & 'S(',ibry,ifield,ng,')%radiation', &
611 & s(ibry,ifield,ng)%radiation, &
612 & trim(ncname)
613 END IF
614 exit_flag=5
615 END IF
616 CASE ('Red')
617 IF (.not.s(ibry,ifield,ng)%reduced) THEN
618 IF (master) THEN
619 WRITE (stdout,20) b(ibry), &
620 & trim(vname(1,idbvar(ifield))), &
621 & trim(string(ibry)), &
622 & 'S(',ibry,ifield,ng,')%reduced', &
623 & s(ibry,ifield,ng)%reduced, &
624 & trim(ncname)
625 END IF
626 exit_flag=5
627 END IF
628 CASE ('Shc')
629 IF (.not.s(ibry,ifield,ng)%Shchepetkin) THEN
630 IF (master) THEN
631 WRITE (stdout,20) b(ibry), &
632 & trim(vname(1,idbvar(ifield))), &
633 & trim(string(ibry)), &
634 & 'S(',ibry,ifield,ng,')%Shchepetkin', &
635 & s(ibry,ifield,ng)%Shchepetkin, &
636 & trim(ncname)
637 END IF
638 exit_flag=5
639 END IF
640 CASE DEFAULT
641 IF (master) THEN
642 WRITE (stdout,30) b(ibry), &
643 & trim(vname(1,idbvar(ifield))), &
644 & trim(string(ibry)), trim(ncname)
645 END IF
646 exit_flag=5
647 END SELECT
648 END DO
649 END IF
650 END DO
651!
652 10 FORMAT (/,' LBC_GETATT_PIO - error while reading global ', &
653 & 'attribute:',2x,a,/,18x,'in restart file:',2x,a,/, &
654 & 18x,'call from:',2x,a, &
655 & /,18x,'Probably global attribute was not found ...', &
656 & /,18x,'restart file needs to be generated by ROMS ', &
657 & 'version 3.6 or higher', &
658 & /,18x,'Alternatively, you may use NO_LBC_ATT at your ', &
659 & 'own risk!')
660 20 FORMAT (/,' LBC_GETATT_PIO - inconsistent ',a,' lateral', &
661 & 'boundary condition for variable: ',2x,a, &
662 & /,18x,'restart file LBC keyword = ',1x,a, &
663 & /,18x,'but assigned structure switch: ', &
664 & 1x,a,i1,',',i2,',',i1,a,' = ',l1, &
665 & /,18x,'check input script LBC keyword for consitency ...',&
666 & /,18x,'restart file:',2x,a)
667 30 FORMAT (/,' LBC_GETATT_PIO - inconsistent ',a,' boundary for ', &
668 & 'variable: ',a,2x,'Keyword = ',a,/,18x,'in input file:', &
669 & 2x,a)
670!
671 RETURN

References mod_scalars::exit_flag, strings_mod::founderror(), mod_ncparam::idbvar, mod_scalars::ieast, mod_scalars::inorth, mod_iounits::ioerror, mod_scalars::isouth, mod_scalars::iwest, mod_parallel::master, mod_iounits::sourcefile, mod_iounits::stdout, and mod_ncparam::vname.

Here is the call graph for this function:

The documentation for this interface was generated from the following file: