FV3 Bundle
tracer_manager.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 
21 ! <CONTACT EMAIL="William.Cooke@noaa.gov">
22 ! William Cooke
23 ! </CONTACT>
24 
25 ! <REVIEWER EMAIL="Matthew.Harrison@noaa.gov">
26 ! Matt Harrison
27 ! </REVIEWER>
28 
29 ! <REVIEWER EMAIL="Bruce.Wyman@noaa.gov">
30 ! Bruce Wyman
31 ! </REVIEWER>
32 
33 ! <REVIEWER EMAIL="Peter.Phillipps@noaa.gov">
34 ! Peter Phillipps
35 ! </REVIEWER>
36 
37 ! <HISTORY SRC="http://www.gfdl.noaa.gov/fms-cgi-bin/cvsweb.cgi/FMS/"/>
38 
39 ! <OVERVIEW>
40 ! Code to manage the simple addition of tracers to the FMS code.
41 ! This code keeps track of the numbers and names of tracers included
42 ! in a tracer table.
43 ! </OVERVIEW>
44 
45 ! <DESCRIPTION>
46 ! This code is a grouping of calls which will allow the simple
47 ! introduction of tracers into the FMS framework. It is designed to
48 ! allow users of a variety of component models interact easily with
49 ! the dynamical core of the model.
50 !
51 ! In calling the tracer manager routines the user must provide a
52 ! parameter identifying the model that the user is working with. This
53 ! parameter is defined within field_manager as MODEL_X
54 ! where X is one of [ATMOS, OCEAN, LAND, ICE].
55 !
56 ! In many of these calls the argument list includes model and tracer_index. These
57 ! are the parameter corresponding to the component model and the tracer_index N is
58 ! the Nth tracer within the component model. Therefore a call with MODEL_ATMOS and 5
59 ! is different from a call with MODEL_OCEAN and 5.
60 !
61 ! </DESCRIPTION>
62 
63 
64 !----------------------------------------------------------------------
65 
66 use mpp_mod, only : mpp_error, &
67  mpp_pe, &
68  mpp_root_pe, &
69  fatal, &
70  warning, &
71  note, &
72  stdlog
73 use mpp_io_mod, only : mpp_open, &
74  mpp_close, &
75  mpp_ascii, &
76  mpp_append, &
77  mpp_rdonly
78 use fms_mod, only : lowercase, &
79  write_version_number
80 
84  model_atmos, &
85  model_land, &
86  model_ocean, &
87  model_ice, &
88  model_coupler, &
89  num_models, &
90  method_type, &
92  parse, &
93  fm_copy_list, &
97  fm_new_value, &
98  fm_exists, &
100 
101 implicit none
102 private
103 
104 !-----------------------------------------------------------------------
105 
106 public tracer_manager_init, &
112  get_tracer_name, &
113  query_method, &
114  set_tracer_atts, &
118  adjust_mass, &
120  no_tracer, &
122 
123 !-----------------------------------------------------------------------
126 end interface
127 !-----------------------------------------------------------------------
128 
129 integer :: num_tracer_fields = 0
130 integer, parameter :: max_tracer_fields = 150
131 integer, parameter :: max_tracer_method = 20
132 integer, parameter :: no_tracer = 1-huge(1)
133 integer, parameter :: notracer = -huge(1)
134 
136 logical :: model_registered(num_models) = .false.
137 
138 type, private :: tracer_type
139  character(len=32) :: tracer_name, tracer_units
140  character(len=128) :: tracer_longname
141  integer :: num_methods, model, instances
142  logical :: is_prognostic, instances_set
143  logical :: needs_init
144 ! Does tracer need mass or positive definite adjustment?
145 ! (true by default for both)
146  logical :: needs_mass_adjust
147  logical :: needs_positive_adjust
148 end type tracer_type
149 
150 type, private :: tracer_name_type
151  character(len=32) :: model_name, tracer_name, tracer_units
152  character(len=128) :: tracer_longname
153 end type tracer_name_type
154 
155 
156 type, private :: inst_type
157  character(len=128) :: name
158  integer :: instances
159 end type inst_type
160 
163 
164 ! Include variable "version" to be written to log file.
165 #include<file_version.h>
166 
167 logical :: module_is_initialized = .false.
168 
169 logical :: verbose_local
171 
172 contains
173 
174 !
175 !#######################################################################
176 !
177 ! <SUBROUTINE NAME="tracer_manager_init">
178 ! <OVERVIEW>
179 ! It is not necessary to call this routine.
180 ! It is included only for backward compatability.
181 ! </OVERVIEW>
182 ! <DESCRIPTION>
183 ! This routine writes the version to the logfile and
184 ! sets the module initialization flag.
185 ! </DESCRIPTION>
186 ! <TEMPLATE>
187 ! call tracer_manager_init
188 ! </TEMPLATE>
189 subroutine tracer_manager_init
190 integer :: model, num_tracers, num_prog, num_diag
191 
192  if(module_is_initialized) return
193  module_is_initialized = .true.
194 
195  call write_version_number ("TRACER_MANAGER_MOD", version)
196  call field_manager_init()
198  do model=1,num_models
199  call get_tracer_meta_data(model, num_tracers, num_prog, num_diag)
200  enddo
201 
202 end subroutine tracer_manager_init
203 ! </SUBROUTINE>
204 
205 !#######################################################################
206 ! <SUBROUTINE NAME="get_tracer_meta_data">
207 ! <OVERVIEW>
208 ! read tracer table and store tracer information associated with "model"
209 ! in "tracers" array.
210 ! </OVERVIEW>
211 subroutine get_tracer_meta_data(model, num_tracers,num_prog,num_diag)
213 integer, intent(in) :: model ! model being used
214 integer, intent(out) :: num_tracers, num_prog, num_diag
215 character(len=256) :: warnmesg
216 
217 character(len=32) :: name_type, type, name
218 integer :: n, m, mod, num_tracer_methods, nfields, swop
219 integer :: j, log_unit, num_methods
220 logical :: flag_type
221 type(method_type), dimension(MAX_TRACER_METHOD) :: methods
222 integer :: instances, siz_inst,i
223 character(len = 32) :: digit,suffnam
224 
225 character(len=128) :: list_name , control
226 integer :: index_list_name
227 logical :: fm_success
228 
229 ! <ERROR MSG="invalid model type" STATUS="FATAL">
230 ! The index for the model type is invalid.
231 ! </ERROR>
232 if (model .ne. model_atmos .and. model .ne. model_land .and. &
233  model .ne. model_ocean .and. model .ne. model_ice .and. &
234  model .ne. model_coupler) call mpp_error(fatal,'tracer_manager_init : invalid model type')
235 
236 ! One should only call get_tracer_meta_data once for each model type
237 ! Therefore need to set up an array to stop the subroutine being
238 ! unnecssarily called multiple times.
239 
240 if ( model_registered(model) ) then
241 ! This routine has already been called for the component model.
242 ! Fill in the values from the previous registration and return.
243  num_tracers = total_tracers(model)
244  num_prog = prog_tracers(model)
245  num_diag = diag_tracers(model)
246  return
247 endif
248 
249 ! Initialize the number of tracers to zero.
250 num_tracers = 0; num_prog = 0; num_diag = 0
251 
252 call field_manager_init(nfields=nfields)
253 
254 ! <ERROR MSG="No tracers are available to be registered." STATUS="NOTE">
255 ! No tracers are available to be registered. This means that the field
256 ! table does not exist or is empty.
257 ! </ERROR>
258 if (nfields == 0 ) then
259 if (mpp_pe() == mpp_root_pe()) &
260  call mpp_error(note,'tracer_manager_init : No tracers are available to be registered.')
261  return
262 endif
263 
264 ! search through field entries for model tracers
265 total_tracers(model) = 0
266 
267 do n=1,nfields
268  call get_field_info(n,type,name,mod,num_methods)
269 
270  if (mod == model .and. type == 'tracer') then
272  total_tracers(model) = total_tracers(model) + 1
273 ! <ERROR MSG="MAX_TRACER_FIELDS exceeded" STATUS="FATAL">
274 ! The maximum number of tracer fields has been exceeded.
275 ! </ERROR>
276  if(num_tracer_fields > max_tracer_fields) call mpp_error(fatal,'tracer_manager_init: MAX_TRACER_FIELDS exceeded')
278  tracers(num_tracer_fields)%model = model
279  tracers(num_tracer_fields)%tracer_name = name
280  tracers(num_tracer_fields)%tracer_units = 'none'
281  tracers(num_tracer_fields)%tracer_longname = tracers(num_tracer_fields)%tracer_name
282  tracers(num_tracer_fields)%instances_set = .false.
283 ! By default, tracers need mass and positive definite adjustments.
284 ! We hardwire exceptions for compatibility with existing field_tables
285 ! This should ideally be cleaned up.
286  tracers(num_tracer_fields)%needs_mass_adjust = .true.
287  tracers(num_tracer_fields)%needs_positive_adjust = .true.
288  if (name == 'cld_amt') then
289  tracers(num_tracer_fields)%needs_mass_adjust = .false.
290  endif
291  if (name == 'cld_amt' .or. name == 'liq_wat' .or. name == 'ice_wat') then
292  tracers(num_tracer_fields)%needs_positive_adjust = .false.
293  endif
294 
295  num_tracer_methods = 0
296  methods = default_method ! initialize methods array
297  call get_field_methods(n,methods)
298  do j=1,num_methods
299  select case (methods(j)%method_type)
300  case ('units')
301  tracers(num_tracer_fields)%tracer_units = methods(j)%method_name
302  case ('longname')
303  tracers(num_tracer_fields)%tracer_longname = methods(j)%method_name
304  case ('instances')
305 ! tracers(num_tracer_fields)%instances = methods(j)%method_name
306  siz_inst = parse(methods(j)%method_name,"",instances)
307  tracers(num_tracer_fields)%instances = instances
308  tracers(num_tracer_fields)%instances_set = .true.
309  case ('adjust_mass')
310  if (methods(j)%method_name == "false") then
311  tracers(num_tracer_fields)%needs_mass_adjust = .false.
312  endif
313  case ('adjust_positive_def')
314  if (methods(j)%method_name == "false") then
315  tracers(num_tracer_fields)%needs_positive_adjust = .false.
316  endif
317  case default
318  num_tracer_methods = num_tracer_methods+1
319 ! tracers(num_tracer_fields)%methods(num_tracer_methods) = methods(j)
320  end select
321  enddo
322  tracers(num_tracer_fields)%num_methods = num_tracer_methods
323  tracers(num_tracer_fields)%needs_init = .false.
324  flag_type = query_method('tracer_type',model,total_tracers(model),name_type)
325  if (flag_type .and. name_type == 'diagnostic') then
326  tracers(num_tracer_fields)%is_prognostic = .false.
327  else
328  tracers(num_tracer_fields)%is_prognostic = .true.
329  endif
330  if (tracers(num_tracer_fields)%is_prognostic) then
331  num_prog = num_prog+1
332  else
333  num_diag = num_diag+1
334  endif
335  endif
336 enddo
337 
338 ! Now cycle through the tracers and add additional instances of the tracers.
339 
340 do n = 1, num_tracer_fields !{
341 ! call get_field_info(n,type,name,mod,num_methods)
342 
343  if ( model == tracers(n)%model .and. tracers(n)%instances_set ) then !{ We have multiple instances of this tracer
344 
345  if ( num_tracer_fields + tracers(n)%instances > max_tracer_fields ) then
346  write(warnmesg, '("tracer_manager_init: Number of tracers will exceed MAX_TRACER_FIELDS with &
347  &multiple (",I3," instances) setup of tracer ",A)') tracers(n)%instances,tracers(n)%tracer_name
348  call mpp_error(fatal, warnmesg)
349  endif
350 
351  do i = 2, tracers(n)%instances !{
353  total_tracers(model) = total_tracers(model) + 1
355  ! Copy the original tracer type to the multiple instances.
357  if ( query_method('instances', model,model_tracer_number(model,n),name, control)) then !{
358 
359  if (i .lt. 10) then !{
360  write (suffnam,'(''suffix'',i1)') i
361  siz_inst = parse(control, suffnam,digit)
362  if (siz_inst == 0 ) then
363  write (digit,'(''_'',i1)') i
364  else
365  digit = "_"//trim(digit)
366  endif
367  elseif (i .lt. 100) then !}{
368  write (suffnam,'(''suffix'',i2)') i
369  siz_inst = parse(control, suffnam,digit)
370  if (siz_inst == 0 ) then
371  write (digit,'(''_'',i2)') i
372  else
373  digit = "_"//trim(digit)
374  endif
375  else !}{
376  call mpp_error(fatal, 'tracer_manager_init: MULTIPLE_TRACER_SET_UP exceeds 100 for '//tracers(n)%tracer_name )
377  endif !}
378 
379  select case(model)
380  case (model_coupler)
381  list_name = "/coupler_mod/tracer/"//trim(tracers(num_tracer_fields)%tracer_name)
382  case (model_atmos)
383  list_name = "/atmos_mod/tracer/"//trim(tracers(num_tracer_fields)%tracer_name)
384  case (model_ocean)
385  list_name = "/ocean_mod/tracer/"//trim(tracers(num_tracer_fields)%tracer_name)
386  case (model_ice )
387  list_name = "/ice_mod/tracer/"//trim(tracers(num_tracer_fields)%tracer_name)
388  case (model_land )
389  list_name = "/land_mod/tracer/"//trim(tracers(num_tracer_fields)%tracer_name)
390  case default
391  list_name = "/default/tracer/"//trim(tracers(num_tracer_fields)%tracer_name)
392  end select
393 
394  if (mpp_pe() == mpp_root_pe() ) write (*,*) "Creating list name = ",trim(list_name)//trim(digit)
395 
396  index_list_name = fm_copy_list(trim(list_name),digit, create = .true.)
397  tracers(num_tracer_fields)%tracer_name = trim(tracers(num_tracer_fields)%tracer_name)//trim(digit)
398  endif !}
399 
400  if (tracers(num_tracer_fields)%is_prognostic) then !{
401  num_prog = num_prog+1
402  else !}{
403  num_diag = num_diag+1
404  endif !}
405  enddo !}
406  ! Multiple instances of tracers were found so need to rename the original tracer.
407  digit = "_1"
408  siz_inst = parse(control, "suffix1",digit)
409  if (siz_inst > 0 ) then !{
410  digit = "_"//trim(digit)
411  endif !}
412  fm_success = fm_modify_name(trim(list_name), trim(tracers(n)%tracer_name)//trim(digit))
413  tracers(n)%tracer_name = trim(tracers(n)%tracer_name)//trim(digit)
414  endif !}
415 enddo !}
416 
417 ! Find any field entries with the instances keyword.
418 do n=1,nfields
419  call get_field_info(n,type,name,mod,num_methods)
420 
421  if ( mod == model .and. type == 'instances' ) then
422  call get_field_methods(n,methods)
423  do j=1,num_methods
424 
425  if (.not.get_tracer_index(mod,methods(j)%method_type,m)) then
426  call mpp_error(fatal,'tracer_manager_init: The instances keyword was found for undefined tracer '&
427  //trim(methods(j)%method_type))
428  else
429  if ( tracers(m)%instances_set ) &
430  call mpp_error(fatal,'tracer_manager_init: The instances keyword was found for '&
431  //trim(methods(j)%method_type)//' but has previously been defined in the tracer entry')
432  siz_inst = parse(methods(j)%method_name,"",instances)
433  tracers(m)%instances = instances
434  call mpp_error(note,'tracer_manager_init: '//trim(instantiations(j)%name)// &
435  ' will have '//trim(methods(j)%method_name)//' instances')
436  endif
437  if ( num_tracer_fields + instances > max_tracer_fields ) then
438  write(warnmesg, '("tracer_manager_init: Number of tracers will exceed MAX_TRACER_FIELDS with &
439  &multiple (",I3," instances) setup of tracer ",A)') tracers(m)%instances,tracers(m)%tracer_name
440  call mpp_error(fatal, warnmesg)
441  endif
442 ! We have found a valid tracer that has more than one instantiation.
443 ! We need to modify that tracer name to tracer_1 and add extra tracers for the extra instantiations.
444  if (instances .eq. 1) then
445  siz_inst = parse(methods(j)%method_control, 'suffix1',digit)
446  if (siz_inst == 0 ) then
447  digit = '_1'
448  else
449  digit = "_"//trim(digit)
450  endif
451  endif
452  do i = 2, instances
454  total_tracers(model) = total_tracers(model) + 1
457 
458  if (i .lt. 10) then !{
459  write (suffnam,'(''suffix'',i1)') i
460  siz_inst = parse(methods(j)%method_control, suffnam,digit)
461  if (siz_inst == 0 ) then
462  write (digit,'(''_'',i1)') i
463  else
464  digit = "_"//trim(digit)
465  endif
466  elseif (i .lt. 100) then !}{
467  write (suffnam,'(''suffix'',i2)') i
468  siz_inst = parse(methods(j)%method_control, suffnam,digit)
469  if (siz_inst == 0 ) then
470  write (digit,'(''_'',i2)') i
471  else
472  digit = "_"//trim(digit)
473  endif
474  else !}{
475  call mpp_error(fatal, 'tracer_manager_init: MULTIPLE_TRACER_SET_UP exceeds 100 for '&
476  //tracers(num_tracer_fields)%tracer_name )
477  endif !}
478 
479  select case(model)
480  case (model_coupler)
481  list_name = "/coupler_mod/tracer/"//trim(tracers(num_tracer_fields)%tracer_name)
482  case (model_atmos)
483  list_name = "/atmos_mod/tracer/"//trim(tracers(num_tracer_fields)%tracer_name)
484  case (model_ocean)
485  list_name = "/ocean_mod/tracer/"//trim(tracers(num_tracer_fields)%tracer_name)
486  case (model_ice )
487  list_name = "/ice_mod/tracer/"//trim(tracers(num_tracer_fields)%tracer_name)
488  case (model_land )
489  list_name = "/land_mod/tracer/"//trim(tracers(num_tracer_fields)%tracer_name)
490  case default
491  list_name = "/default/tracer/"//trim(tracers(num_tracer_fields)%tracer_name)
492  end select
493 
494  if (mpp_pe() == mpp_root_pe() ) write (*,*) "Creating list name = ",trim(list_name)
495 
496  index_list_name = fm_copy_list(trim(list_name),digit, create = .true.)
497 
498  tracers(num_tracer_fields)%tracer_name = trim(tracers(num_tracer_fields)%tracer_name)//digit
499  if (tracers(num_tracer_fields)%is_prognostic) then
500  num_prog = num_prog+1
501  else
502  num_diag = num_diag+1
503  endif
504  enddo
505 !Now rename the original tracer to tracer_1 (or if suffix1 present to tracer_'value_of_suffix1')
506  siz_inst = parse(methods(j)%method_control, 'suffix1',digit)
507  if (siz_inst == 0 ) then
508  digit = '_1'
509  else
510  digit = "_"//trim(digit)
511  endif
512  fm_success = fm_modify_name(trim(list_name), trim(tracers(m)%tracer_name)//trim(digit))
513  tracers(m)%tracer_name = trim(tracers(m)%tracer_name)//trim(digit)
514  enddo
515  endif
516 enddo
517 
518 num_tracers = num_prog + num_diag
519 ! Make the number of tracers available publicly.
520 total_tracers(model) = num_tracers
521 prog_tracers(model) = num_prog
522 diag_tracers(model) = num_diag
523 model_registered(model) = .true.
524 
525 ! Now sort through the tracer fields and sort them so that the
526 ! prognostic tracers are first.
527 
528 do n=1, num_tracers
529  if (.not.check_if_prognostic(model,n) .and. n.le.num_prog) then
530  ! This is a diagnostic tracer so find a prognostic tracer to swop with
531  do m = n, num_tracers
532  if (check_if_prognostic(model,m) .and. .not.check_if_prognostic(model,n)) then
533  swop = tracer_array(model,n)
534  tracer_array(model,n) = tracer_array(model,m)
535  tracer_array(model,m) = swop
536  cycle
537  endif
538  enddo
539  endif
540 enddo
541 
542 do n=1, num_tracer_fields
543  call print_tracer_info(model,n)
544 enddo
545 
546 log_unit = stdlog()
547 if ( mpp_pe() == mpp_root_pe() ) then
548  write (log_unit,15) trim(model_names(model)),total_tracers(model)
549 endif
550 
551 15 format ('Number of tracers in field table for ',a,' model = ',i4)
552 
553 end subroutine get_tracer_meta_data
554 !</SUBROUTINE>
555 
556 
557 function model_tracer_number(model,n)
558 integer, intent(in) :: model, n
559 integer model_tracer_number
560 
561 integer :: i
562 
564 
565 do i = 1, max_tracer_fields
566  if ( tracer_array(model,i) == n ) then
568  return
569  endif
570 enddo
571 
572 end function model_tracer_number
573 
574 !#######################################################################
575 !
576 ! <SUBROUTINE NAME="register_tracers">
577 
578 ! <OVERVIEW>
579 ! It is not necessary to call this routine.
580 ! It is included only for backward compatability.
581 ! </OVERVIEW>
582 ! <DESCRIPTION>
583 ! This routine returns the total number of valid tracers,
584 ! the number of prognostic and diagnostic tracers.
585 ! </DESCRIPTION>
586 ! <TEMPLATE>
587 ! call register_tracers(model, num_tracers,num_prog,num_diag)
588 ! </TEMPLATE>
589 
590 ! <IN NAME="model" TYPE="integer">
591 ! A parameter to identify which model is being used.
592 ! </IN>
593 ! <OUT NAME="num_tracers" TYPE="integer">
594 ! The total number of valid tracers within the component model.
595 ! </OUT>
596 ! <OUT NAME="num_prog" TYPE="integer">
597 ! The number of prognostic tracers within the component model.
598 ! </OUT>
599 ! <OUT NAME="num_diag" TYPE="integer">
600 ! The number of diagnostic tracers within the component model.
601 ! </OUT>
602 subroutine register_tracers(model, num_tracers, num_prog, num_diag, num_family)
603 integer, intent(in) :: model
604 integer, intent(out) :: num_tracers, num_prog, num_diag
605 integer, intent(out), optional :: num_family
606 
608 
609 call get_number_tracers(model, num_tracers, num_prog, num_diag, num_family)
610 
611 end subroutine register_tracers
612 !</SUBROUTINE>
613 
614 !#######################################################################
615 
616 ! <SUBROUTINE NAME="get_number_tracers">
617 ! <OVERVIEW>
618 ! A routine to return the number of tracers included in a component model.
619 ! </OVERVIEW>
620 ! <DESCRIPTION>
621 ! This routine returns the total number of valid tracers,
622 ! the number of prognostic and diagnostic tracers
623 ! </DESCRIPTION>
624 ! <TEMPLATE>
625 ! call get_number_tracers(model, num_tracers,num_prog,num_diag)
626 ! </TEMPLATE>
627 
628 ! <IN NAME="model" TYPE="integer">
629 ! A parameter to identify which model is being used.
630 ! </IN>
631 ! <OUT NAME="num_tracers" TYPE="integer, optional">
632 ! The total number of valid tracers within the component model.
633 ! </OUT>
634 ! <OUT NAME="num_prog" TYPE="integer, optional">
635 ! The number of prognostic tracers within the component model.
636 ! </OUT>
637 ! <OUT NAME="num_diag" TYPE="integer, optional">
638 ! The number of diagnostic tracers within the component model.
639 ! </OUT>
640 subroutine get_number_tracers(model, num_tracers, num_prog, num_diag, num_family)
642 integer, intent(in) :: model
643 integer, intent(out), optional :: num_tracers, num_prog, num_diag, num_family
644 
646 
647 ! <ERROR MSG="Model number is invalid." STATUS="FATAL">
648 ! The index of the component model is invalid.
649 ! </ERROR>
650 if (model .ne. model_atmos .and. model .ne. model_land .and. &
651  model .ne. model_ocean .and. model .ne. model_ice .and. &
652  model .ne. model_coupler) &
653  call mpp_error(fatal,"get_number_tracers : Model number is invalid.")
654 
655 if (present(num_tracers)) num_tracers = total_tracers(model)
656 if (present(num_prog)) num_prog = prog_tracers(model)
657 if (present(num_diag)) num_diag = diag_tracers(model)
658 if (present(num_family)) num_family = 0 ! Needed only for backward compatability with lima
659 
660 end subroutine get_number_tracers
661 !</SUBROUTINE>
662 
663 
664 ! <SUBROUTINE NAME="get_tracer_indices">
665 
666 ! <OVERVIEW>
667 ! Routine to return the component model tracer indices as defined within
668 ! the tracer manager.
669 ! </OVERVIEW>
670 ! <DESCRIPTION>
671 ! If several models are being used or redundant tracers have been written to
672 ! the tracer_table, then the indices in the component model and the tracer
673 ! manager may not have a one to one correspondence. Therefore the component
674 ! model needs to know what index to pass to calls to tracer_manager routines in
675 ! order that the correct tracer information be accessed.
676 ! </DESCRIPTION>
677 ! <TEMPLATE>
678 ! call get_tracer_indices(model, ind, prog_ind, diag_ind)
679 ! </TEMPLATE>
680 
681 ! <IN NAME="model" TYPE="integer">
682 ! A parameter to identify which model is being used.
683 ! </IN>
684 ! <OUT NAME="ind" TYPE="integer, optional" DIM="(:)" >
685 ! An array containing the tracer manager defined indices for
686 ! all the tracers within the component model.
687 ! </OUT>
688 ! <OUT NAME="prog_ind" TYPE="integer, optional" DIM="(:)" >
689 ! An array containing the tracer manager defined indices for
690 ! the prognostic tracers within the component model.
691 ! </OUT>
692 ! <OUT NAME="diag_ind" TYPE="integer, optional" DIM="(:)" >
693 ! An array containing the tracer manager defined indices for
694 ! the diagnostic tracers within the component model.
695 ! </OUT>
696 subroutine get_tracer_indices(model, ind, prog_ind, diag_ind, fam_ind)
698 integer, intent(in) :: model
699 integer, intent(out), dimension(:), optional :: ind, prog_ind, diag_ind, fam_ind
700 
701 integer :: i, j, np, nd, n
702 
704 
705 nd=0;np=0;n=0
706 
707 ! Initialize arrays with dummy values
708 if (PRESENT(ind)) ind = no_tracer
709 if (PRESENT(prog_ind)) prog_ind = no_tracer
710 if (PRESENT(diag_ind)) diag_ind = no_tracer
711 if (PRESENT(fam_ind)) fam_ind = no_tracer
712 
713 do i = 1, max_tracer_fields
714 j = tracer_array(model,i)
715  if ( j /= notracer) then
716  if ( model == tracers(j)%model) then
717  if (PRESENT(ind)) then
718  n=n+1
719 ! <ERROR MSG="index array size too small in get_tracer_indices" STATUS="FATAL">
720 ! The global index array is too small and cannot contain all the tracer numbers.
721 ! </ERROR>
722  if (n > size(ind(:))) call mpp_error(fatal,'get_tracer_indices : index array size too small in get_tracer_indices')
723  ind(n) = i
724  endif
725 
726  if (tracers(j)%is_prognostic.and.PRESENT(prog_ind)) then
727  np=np+1
728 ! <ERROR MSG="prognostic array size too small in get_tracer_indices" STATUS="FATAL">
729 ! The prognostic index array is too small and cannot contain all the tracer numbers.
730 ! </ERROR>
731  if ( np > size( prog_ind(:)))call mpp_error(fatal,&
732  'get_tracer_indices : prognostic array size too small in get_tracer_indices')
733  prog_ind(np) = i
734  else if (.not.tracers(j)%is_prognostic .and. PRESENT(diag_ind)) then
735  nd = nd+1
736 ! <ERROR MSG="diagnostic array size too small in get_tracer_indices" STATUS="FATAL">
737 ! The diagnostic index array is too small and cannot contain all the tracer numbers.
738 ! </ERROR>
739  if (nd > size(diag_ind(:))) call mpp_error(fatal,&
740  'get_tracer_indices : diagnostic array size too small in get_tracer_indices')
741  diag_ind(nd) = i
742  endif
743  endif
744  endif
745 enddo
746 
747 return
748 end subroutine get_tracer_indices
749 !</SUBROUTINE>
750 
751 !<FUNCTION NAME= "get_tracer_index">
752 ! <OVERVIEW>
753 ! Function which returns the number assigned to the tracer name.
754 ! </OVERVIEW>
755 ! <DESCRIPTION>
756 ! This is a function which returns the index, as implied within the component model.
757 ! There are two overloaded interfaces: one of type integer, one logical.
758 ! </DESCRIPTION>
759 ! <TEMPLATE>
760 ! integer: index = get_tracer_index(model, name, indices, verbose)
761 ! logical: if ( get_tracer_index(model, name, index, indices, verbose) ) then
762 ! </TEMPLATE>
763 ! <IN NAME="model" TYPE="integer">
764 ! A parameter to identify which model is being used.
765 ! </IN>
766 ! <IN NAME="name" TYPE="character">
767 ! The name of the tracer (as assigned in the field table).
768 ! </IN>
769 ! <IN NAME="indices" TYPE="integer, optional" DIM="(:)">
770 ! An array indices.
771 ! When present, the returned index will limit the search for the tracer
772 ! to those tracers whos indices are amoung those in array "indices".
773 ! This would be useful when it is desired to limit the search to a subset
774 ! of the tracers. Such a subset might be the diagnostic or prognostic tracers.
775 ! (Note that subroutine get_tracer_indices returns these subsets)
776 ! </IN>
777 ! <IN NAME="verbose" TYPE="logical, optional">
778 ! A flag to allow the message saying that a tracer with this name has not
779 ! been found. This should only be used for debugging purposes.
780 ! </IN>
781 ! <OUT NAME="get_tracer_index" TYPE="integer">
782 ! integer function:
783 ! The index of the tracer named "name".
784 ! If no tracer by that name exists then the returned value is NO_TRACER.
785 ! logical function:
786 ! If no tracer by that name exists then the returned value is .false.,
787 ! otherwise the returned value is .true.
788 ! </OUT>
789 function get_tracer_index_integer(model, name, indices, verbose)
791 integer, intent(in) :: model
792 character(len=*), intent(in) :: name
793 integer, intent(in), dimension(:), optional :: indices
794 logical, intent(in), optional :: verbose
795 integer :: get_tracer_index_integer
796 
797 integer :: i
798 
800 
802 
803 if (PRESENT(indices)) then
804  do i = 1, size(indices(:))
805  if (model == tracers(indices(i))%model .and. lowercase(trim(name)) == trim(tracers(indices(i))%tracer_name)) then
807  exit
808  endif
809  enddo
810 else
811  do i=1, num_tracer_fields
812  if(tracer_array(model,i) == notracer) cycle
813  if (lowercase(trim(name)) == trim(tracers(tracer_array(model,i))%tracer_name)) then
814  get_tracer_index_integer = i!TRACER_ARRAY(model,i)
815  exit
816  endif
817  enddo
818 end if
819 
820 verbose_local=.false.
821 if (present(verbose)) verbose_local=verbose
822 
823 if (verbose_local) then
824 ! <ERROR MSG="tracer with this name not found: X" STATUS="NOTE">
825  if (get_tracer_index_integer == no_tracer ) then
826  call mpp_error(note,'get_tracer_index : tracer with this name not found: '//trim(name))
827  endif
828 ! </ERROR>
829 endif
830 
831 return
832 
833 end function get_tracer_index_integer
834 
835 !#######################################################################
836 function get_tracer_index_logical(model, name, index, indices, verbose)
838 integer, intent(in) :: model
839 character(len=*), intent(in) :: name
840 integer, intent(out) :: index
841 integer, intent(in), dimension(:), optional :: indices
842 logical, intent(in), optional :: verbose
843 logical :: get_tracer_index_logical
844 
845 index = get_tracer_index_integer(model, name, indices, verbose)
846 if(index == no_tracer) then
847  get_tracer_index_logical = .false.
848 else
849  get_tracer_index_logical = .true.
850 endif
851 
852 end function get_tracer_index_logical
853 !</FUNCTION>
854 
855 !#######################################################################
856 ! <SUBROUTINE NAME="tracer_manager_end" >
857 ! <OVERVIEW>
858 ! Routine to write to the log file that the tracer manager is ending.
859 ! </OVERVIEW>
860 ! <DESCRIPTION>
861 ! Routine to write to the log file that the tracer manager is ending.
862 ! </DESCRIPTION>
863 ! <TEMPLATE>
864 ! call tracer_manager_end
865 ! </TEMPLATE>
866 subroutine tracer_manager_end
868 integer :: log_unit
869 
870 log_unit = stdlog()
871 if ( mpp_pe() == mpp_root_pe() ) then
872  write (log_unit,'(/,(a))') 'Exiting tracer_manager, have a nice day ...'
873 endif
874 
875 module_is_initialized = .false.
876 
877 end subroutine tracer_manager_end
878 !</SUBROUTINE>
879 
880 !#######################################################################
881 !
882 subroutine print_tracer_info(model,n)
883 !
884 ! Routine to print out the components of the tracer.
885 ! This is useful for informational purposes.
886 ! Used in get_tracer_meta_data.
887 !
888 ! Arguments:
889 ! INTENT IN
890 ! i : index of the tracer that is being printed.
891 !
892 integer, intent(in) :: model,n
893 integer :: i,log_unit
894 
896 
897 if(mpp_pe()==mpp_root_pe() .and. tracer_array(model,n)> 0 ) then
898  i = tracer_array(model,n)
899  log_unit = stdlog()
900  write(log_unit, *)'----------------------------------------------------'
901  write(log_unit, *) 'Contents of tracer entry ', i
902  write(log_unit, *) 'Model type and field name'
903  write(log_unit, *) 'Model : ', tracers(i)%model
904  write(log_unit, *) 'Field name : ', trim(tracers(i)%tracer_name)
905  write(log_unit, *) 'Tracer units : ', trim(tracers(i)%tracer_units)
906  write(log_unit, *) 'Tracer longname : ', trim(tracers(i)%tracer_longname)
907  write(log_unit, *) 'Tracer is_prognostic : ', tracers(i)%is_prognostic
908  write(log_unit, *)'----------------------------------------------------'
909 endif
910 
911 900 FORMAT(a,2(1x,e12.6))
912 901 FORMAT(e12.6,1x,e12.6)
913 
914 
915 end subroutine print_tracer_info
916 
917 !#######################################################################
918 !
919 ! <SUBROUTINE NAME="get_tracer_names" >
920 ! <OVERVIEW>
921 ! Routine to find the names associated with a tracer number.
922 ! </OVERVIEW>
923 ! <DESCRIPTION>
924 ! This routine can return the name, long name and units associated
925 ! with a tracer.
926 ! </DESCRIPTION>
927 ! <TEMPLATE>
928 ! call get_tracer_names(model,n,name,longname, units)
929 ! </TEMPLATE>
930 
931 ! <IN NAME="model" TYPE="integer">
932 ! A parameter representing the component model in use.
933 ! </IN>
934 ! <IN NAME="n" TYPE="integer">
935 ! Tracer number.
936 ! </IN>
937 ! <OUT NAME="name" TYPE="character" >
938 ! Field name associated with tracer number.
939 ! </OUT>
940 ! <OUT NAME="longname" TYPE="character, optional" >
941 ! The long name associated with tracer number.
942 ! </OUT>
943 ! <OUT NAME="units" TYPE="character, optional" >
944 ! The units associated with tracer number.
945 ! </OUT>
946 
947 subroutine get_tracer_names(model,n,name,longname, units, err_msg)
949 integer, intent(in) :: model, n
950 character (len=*),intent(out) :: name
951 character (len=*), intent(out), optional :: longname, units, err_msg
952 character (len=128) :: err_msg_local
953 integer :: n1
954 character(len=11) :: chn
955 
957 
958  if (n < 1 .or. n > total_tracers(model)) then
959  write(chn, '(i11)') n
960  err_msg_local = ' Invalid tracer index. Model name = '//trim(model_names(model))//', Index='//trim(chn)
961  if(error_handler('get_tracer_names', err_msg_local, err_msg)) return
962  endif
963  n1 = tracer_array(model,n)
964 
965 name = trim(tracers(n1)%tracer_name)
966 if (PRESENT(longname)) longname = trim(tracers(n1)%tracer_longname)
967 if (PRESENT(units)) units = trim(tracers(n1)%tracer_units)
968 
969 end subroutine get_tracer_names
970 !</SUBROUTINE>
971 !
972 !#######################################################################
973 !
974 ! <FUNCTION NAME="get_tracer_name" >
975 ! <OVERVIEW>
976 ! Routine to find the names associated with a tracer number.
977 ! </OVERVIEW>
978 ! <DESCRIPTION>
979 ! This routine can return the name, long name and units associated with a tracer.
980 ! The return value of get_tracer_name is .false. when a FATAL error condition is
981 ! detected, otherwise the return value is .true.
982 ! </DESCRIPTION>
983 ! <TEMPLATE>
984 ! if(.not.get_tracer_name(model,n,name,longname, units, err_msg)) call mpp_error(.....
985 ! </TEMPLATE>
986 
987 ! <IN NAME="model" TYPE="integer">
988 ! A parameter representing the component model in use.
989 ! </IN>
990 ! <IN NAME="n" TYPE="integer">
991 ! Tracer number.
992 ! </IN>
993 ! <OUT NAME="name" TYPE="character" >
994 ! Field name associated with tracer number.
995 ! </OUT>
996 ! <OUT NAME="longname" TYPE="character, optional" >
997 ! The long name associated with tracer number.
998 ! </OUT>
999 ! <OUT NAME="units" TYPE="character, optional" >
1000 ! The units associated with tracer number.
1001 ! </OUT>
1002 ! <OUT NAME="err_msg" TYPE="character, optional" >
1003 ! When present:
1004 ! If a FATAL error condition is detected then err_msg will contain an error message
1005 ! and the return value of get_tracer_name will be .false.
1006 ! If no FATAL error is detected err_msg will be filled with space characters and
1007 ! and the return value of get_tracer_name will be .true.
1008 ! When not present:
1009 ! A FATAL error will result in termination inside get_tracer_name without returning.
1010 ! If no FATAL error is detected the return value of get_tracer_name will be .true.
1011 ! </OUT>
1012 
1013 function get_tracer_name(model,n,name,longname, units, err_msg)
1015 logical :: get_tracer_name
1016 integer, intent(in) :: model, n
1017 character (len=*),intent(out) :: name
1018 character (len=*), intent(out), optional :: longname, units, err_msg
1019 character (len=128) :: err_msg_local
1020 integer :: n1
1021 character(len=11) :: chn
1022 
1024 
1025  if (n < 1 .or. n > total_tracers(model)) then
1026  write(chn, '(i11)') n
1027  err_msg_local = ' Invalid tracer index. Model name = '//trim(model_names(model))//', Index='//trim(chn)
1028  if(error_handler('get_tracer_name', err_msg_local, err_msg)) then
1029  get_tracer_name = .false.
1030  return
1031  endif
1032  else
1033  get_tracer_name = .true.
1034  endif
1035  n1 = tracer_array(model,n)
1036 
1037 name = trim(tracers(n1)%tracer_name)
1038 if (PRESENT(longname)) longname = trim(tracers(n1)%tracer_longname)
1039 if (PRESENT(units)) units = trim(tracers(n1)%tracer_units)
1040 
1041 end function get_tracer_name
1042 !</FUNCTION>
1043 !
1044 !#######################################################################
1045 !
1046 !<FUNCTION NAME= "check_if_prognostic">
1047 ! <OVERVIEW>
1048 ! Function to see if a tracer is prognostic or diagnostic.
1049 ! </OVERVIEW>
1050 ! <DESCRIPTION>
1051 ! All tracers are assumed to be prognostic when read in from the field_table
1052 ! However a tracer can be changed to a diagnostic tracer by adding the line
1053 ! "tracer_type","diagnostic"
1054 ! to the tracer description in field_table.
1055 ! </DESCRIPTION>
1056 ! <TEMPLATE>
1057 ! logical =check_if_prognostic(model, n)
1058 ! </TEMPLATE>
1059 
1060 ! <IN NAME="model" TYPE="integer">
1061 ! A parameter representing the component model in use.
1062 ! </IN>
1063 ! <IN NAME="n" TYPE="integer">
1064 ! Tracer number
1065 ! </IN>
1066 ! <OUT NAME="check_if_prognostic" TYPE="logical">
1067 ! A logical flag set TRUE if the tracer is
1068 ! prognostic.
1069 ! </OUT>
1070 function check_if_prognostic(model, n, err_msg)
1072 integer, intent(in) :: model, n
1073 logical :: check_if_prognostic
1074 character(len=*), intent(out), optional :: err_msg
1075 character(len=128) :: err_msg_local
1076 character(len=11) :: chn
1077 
1079 
1080 if (n < 1 .or. n > total_tracers(model)) then
1081  write(chn, '(i11)') n
1082  err_msg_local = ' Invalid tracer index. Model name = '//trim(model_names(model))//', Index='//trim(chn)
1083  check_if_prognostic = .true.
1084  if(error_handler('check_if_prognostic', err_msg_local, err_msg)) return
1085 endif
1086 
1087 !Convert local model index to tracer_manager index
1088 
1089 check_if_prognostic = tracers(tracer_array(model,n))%is_prognostic
1090 
1091 end function check_if_prognostic
1092 !</FUNCTION>
1093 
1094 ! Does tracer need mass or positive definite adjustments?
1095 !#######################################################################
1096 ! Function to check whether tracer should have its mass adjusted
1097 function adjust_mass(model, n, err_msg)
1099 integer, intent(in) :: model, n
1100 logical :: adjust_mass
1101 character(len=*), intent(out), optional :: err_msg
1102 character(len=128) :: err_msg_local
1103 character(len=11) :: chn
1104 
1106 
1107 if (n < 1 .or. n > total_tracers(model)) then
1108  write(chn, '(i11)') n
1109  err_msg_local = ' Invalid tracer index. Model name = '//trim(model_names(model))//', Index='//trim(chn)
1110  adjust_mass = .true.
1111  if(error_handler('adjust_mass', err_msg_local, err_msg)) return
1112 endif
1113 
1114 !Convert local model index to tracer_manager index
1115 
1116 adjust_mass = tracers(tracer_array(model,n))%needs_mass_adjust
1117 
1118 end function adjust_mass
1119 
1120 ! Function to check whether tracer should be adjusted to remain positive definite
1121 function adjust_positive_def(model, n, err_msg)
1123 integer, intent(in) :: model, n
1124 logical :: adjust_positive_def
1125 character(len=*), intent(out), optional :: err_msg
1126 character(len=128) :: err_msg_local
1127 character(len=11) :: chn
1128 
1130 
1131 if (n < 1 .or. n > total_tracers(model)) then
1132  write(chn, '(i11)') n
1133  err_msg_local = ' Invalid tracer index. Model name = '//trim(model_names(model))//', Index='//trim(chn)
1134  adjust_positive_def = .true.
1135  if(error_handler('adjust_positive_def', err_msg_local, err_msg)) return
1136 endif
1137 
1138 !Convert local model index to tracer_manager index
1139 
1140 adjust_positive_def = tracers(tracer_array(model,n))%needs_positive_adjust
1141 
1142 end function adjust_positive_def
1143 
1144 !
1145 !#######################################################################
1146 !
1147 ! <SUBROUTINE NAME="set_tracer_profile" >
1148 ! <OVERVIEW>
1149 ! Subroutine to set the tracer field to the wanted profile.
1150 ! </OVERVIEW>
1151 ! <DESCRIPTION>
1152 ! If the profile type is 'fixed' then the tracer field values are set
1153 ! equal to the surface value.
1154 ! If the profile type is 'profile' then the top/bottom of model and
1155 ! surface values are read and an exponential profile is calculated,
1156 ! with the profile being dependent on the number of levels in the
1157 ! component model. This should be called from the part of the dynamical
1158 ! core where tracer restarts are called in the event that a tracer
1159 ! restart file does not exist.
1160 !
1161 ! This can be activated by adding a method to the field_table
1162 ! e.g.
1163 ! "profile_type","fixed","surface_value = 1e-12"
1164 ! would return values of surf_value = 1e-12 and a multiplier of 1.0
1165 ! One can use these to initialize the entire field with a value of 1e-12.
1166 !
1167 ! "profile_type","profile","surface_value = 1e-12, top_value = 1e-15"
1168 ! In a 15 layer model this would return values of surf_value = 1e-12 and
1169 ! multiplier = 0.6309573 i.e 1e-15 = 1e-12*(0.6309573^15)
1170 ! In this case the model should be MODEL_ATMOS as you have a "top" value.
1171 !
1172 ! If you wish to initialize the ocean model, one can use bottom_value instead
1173 ! of top_value.
1174 
1175 ! </DESCRIPTION>
1176 ! <TEMPLATE>
1177 ! call set_tracer_profile(model, n, tracer)
1178 ! </TEMPLATE>
1179 
1180 ! <IN NAME="model" TYPE="integer">
1181 ! A parameter representing the component model in use.
1182 ! </IN>
1183 ! <IN NAME="n" TYPE="integer">
1184 ! Tracer number.
1185 ! </IN>
1186 ! <INOUT NAME="tracer_array" TYPE="real">
1187 ! The initialized tracer array.
1188 ! </INOUT>
1189 
1190 subroutine set_tracer_profile(model, n, tracer, err_msg)
1192 integer, intent(in) :: model, n
1193  real, intent(inout), dimension(:,:,:) :: tracer
1194 character(len=*), intent(out), optional :: err_msg
1195 
1196 real :: surf_value, multiplier
1197 integer :: numlevels, k, n1, flag
1198 real :: top_value, bottom_value
1199 character(len=80) :: scheme, control,profile_type
1200 character(len=128) :: err_msg_local
1201 character(len=11) :: chn
1202 
1204 
1205 if (n < 1 .or. n > total_tracers(model)) then
1206  write(chn, '(i11)') n
1207  err_msg_local = ' Invalid tracer index. Model name = '//trim(model_names(model))//', Index='//trim(chn)
1208  if(error_handler('set_tracer_profile', err_msg_local, err_msg)) return
1209 endif
1210 n1 = tracer_array(model,n)
1211 
1212 !default values
1213 profile_type = 'Fixed'
1214 surf_value = 0.0e+00
1215 top_value = surf_value
1216 bottom_value = surf_value
1217 multiplier = 1.0
1218 
1219 tracer = surf_value
1220 
1221 if ( query_method( 'profile_type',model,n,scheme,control)) then
1222 !Change the tracer_number to the tracer_manager version
1223 
1224  if(lowercase(trim(scheme(1:5))).eq.'fixed') then
1225  profile_type = 'Fixed'
1226  flag =parse(control,'surface_value',surf_value)
1227  multiplier = 1.0
1228  tracer = surf_value
1229  endif
1230 
1231  if(lowercase(trim(scheme(1:7))).eq.'profile') then
1232  profile_type = 'Profile'
1233  flag=parse(control,'surface_value',surf_value)
1234  if (surf_value .eq. 0.0) &
1235  call mpp_error(fatal,'set_tracer_profile : Cannot have a zero surface value for an exponential profile. Tracer '&
1236  //tracers(n1)%tracer_name//" "//control//" "//scheme)
1237  select case (tracers(n1)%model)
1238  case (model_atmos)
1239  flag=parse(control,'top_value',top_value)
1240  if(mpp_pe()==mpp_root_pe() .and. flag == 0) &
1241  call mpp_error(note,'set_tracer_profile : Parameter top_value needs to be defined for the tracer profile.')
1242  case (model_ocean)
1243  flag =parse(control,'bottom_value',bottom_value)
1244  if(mpp_pe() == mpp_root_pe() .and. flag == 0) &
1245  call mpp_error(note,'set_tracer_profile : Parameter bottom_value needs to be defined for the tracer profile.')
1246  case default
1247 ! Should there be a NOTE or WARNING message here?
1248  end select
1249 
1250 ! If profile type is profile then set the surface value to the input
1251 ! value and calculate the vertical multiplier.
1252 !
1253 ! Assume an exponential decay/increase from the surface to the top level
1254 ! C = C0 exp ( -multiplier* level_number)
1255 ! => multiplier = exp [ ln(Ctop/Csurf)/number_of_levels]
1256 !
1257 numlevels = size(tracer,3) -1
1258  select case (tracers(n1)%model)
1259  case (model_atmos)
1260  multiplier = exp( log(top_value/surf_value) /numlevels)
1261  tracer(:,:,1) = surf_value
1262  do k = 2, size(tracer,3)
1263  tracer(:,:,k) = tracer(:,:,k-1) * multiplier
1264  enddo
1265  case (model_ocean)
1266  multiplier = exp( log(bottom_value/surf_value) /numlevels)
1267  tracer(:,:,size(tracer,3)) = surf_value
1268  do k = size(tracer,3) - 1, 1, -1
1269  tracer(:,:,k) = tracer(:,:,k+1) * multiplier
1270  enddo
1271  case default
1272  end select
1273  endif !scheme.eq.profile
1274 
1275  if (mpp_pe() == mpp_root_pe() ) write(*,700) 'Tracer ',trim(tracers(n1)%tracer_name), &
1276  ' initialized with surface value of ',surf_value, &
1277  ' and vertical multiplier of ',multiplier
1278  700 FORMAT (3a,e12.6,a,f10.6)
1279 
1280 endif ! end of query scheme
1281 
1282 end subroutine set_tracer_profile
1283 !</SUBROUTINE>
1284 
1285 !
1286 !#######################################################################
1287 !
1288 ! <FUNCTION NAME="query_method" >
1289 ! <OVERVIEW>
1290 ! A function to query the "methods" associated with each tracer.
1291 ! </OVERVIEW>
1292 ! <DESCRIPTION>
1293 ! A function to query the "methods" associated with each tracer. The
1294 ! "methods" are the parameters of the component model that can be
1295 ! adjusted by user by placing formatted strings, associated with a
1296 ! particular tracer, within the field table.
1297 ! These methods can control the advection, wet deposition, dry
1298 ! deposition or initial profile of the tracer in question. Any
1299 ! parametrization can use this function as long as a routine for parsing
1300 ! the name and control strings are provided by that routine.
1301 ! </DESCRIPTION>
1302 ! <TEMPLATE>
1303 ! logical =query_method (method_type, model, n, name, control)
1304 ! </TEMPLATE>
1305 
1306 ! <IN NAME="method_type" TYPE="character">
1307 ! The method that is being requested.
1308 ! </IN>
1309 ! <IN NAME="model" TYPE="integer">
1310 ! A parameter representing the component model in use.
1311 ! </IN>
1312 ! <IN NAME="n" TYPE="integer">
1313 ! Tracer number
1314 ! </IN>
1315 ! <OUT NAME="name" TYPE="character">
1316 ! A string containing the modified name to be used with
1317 ! method_type. i.e. "2nd_order" might be the default for
1318 ! advection. One could use "4th_order" here to modify
1319 ! that behaviour.
1320 ! </OUT>
1321 ! <OUT NAME="control" TYPE="character, optional">
1322 ! A string containing the modified parameters that are
1323 ! associated with the method_type and name.
1324 ! </OUT>
1325 ! <OUT NAME="query_method" TYPE="logical">
1326 ! A flag to show whether method_type exists with regard to
1327 ! tracer n. If method_type is not present then one must
1328 ! have default values.
1329 ! </OUT>
1330 
1331 !<NOTE>
1332 ! At present the tracer manager module allows the initialization of a tracer
1333 ! profile if a restart does not exist for that tracer.
1334 ! Options for this routine are as follows
1335 !
1336 ! Tracer profile setup
1337 ! ==================================================================
1338 ! |method_type |method_name |method_control |
1339 ! ==================================================================
1340 ! |profile_type |fixed |surface_value = X |
1341 ! |profile_type |profile |surface_value = X, top_value = Y |(atmosphere)
1342 ! |profile_type |profile |surface_value = X, bottom_value = Y |(ocean)
1343 ! ==================================================================
1344 !
1345 !</NOTE>
1346  function query_method (method_type, model, n, name, control, err_msg)
1348 ! A function to query the schemes associated with each tracer.
1349 !
1350 ! INTENT IN
1351 ! method_type : The method that is being requested.
1352 ! model : The model that you are calling this function from.
1353 ! n : The tracer number.
1354 ! INTENT OUT
1355 ! name : A string containing the modified name to be used with
1356 ! method_type. i.e. "2nd_order" might be the default for
1357 ! advection. One could use "4th_order" here to modify
1358 ! that behaviour.
1359 ! control : A string containing the modified parameters that are
1360 ! associated with the method_type and name.
1361 ! query_method : A flag to show whether method_type exists with regard
1362 ! to tracer n. If method_type is not present then one
1363 ! must have default values.
1364 
1365  character(len=*), intent(in) :: method_type
1366  integer , intent(in) :: model, n
1367  character(len=*), intent(out) :: name
1368  character(len=*), intent(out), optional :: control, err_msg
1369  logical :: query_method
1370 
1371  integer :: n1
1372  character(len=256) :: list_name
1373  character(len=1024):: control_tr
1374  character(len=16) :: chn,chn1
1375  character(len=128) :: err_msg_local
1376 
1378 
1379 !Convert the local model tracer number to the tracer_manager version.
1380 
1381  if (n < 1 .or. n > total_tracers(model)) then
1382  write(chn, '(i11)') n
1383  err_msg_local = ' Invalid tracer index. Model name = '//trim(model_names(model))//', Index='//trim(chn)
1384  if(error_handler('query_method', err_msg_local, err_msg)) return
1385  endif
1386 
1387  n1 = tracer_array(model,n)
1388 
1389  select case(model)
1390  case (model_coupler)
1391  list_name = "/coupler_mod/tracer/"//trim(tracers(n1)%tracer_name)//"/"//trim(method_type)
1392  case (model_atmos)
1393  list_name = "/atmos_mod/tracer/"//trim(tracers(n1)%tracer_name)//"/"//trim(method_type)
1394  case (model_ocean)
1395  list_name = "/ocean_mod/tracer/"//trim(tracers(n1)%tracer_name)//"/"//trim(method_type)
1396  case (model_ice )
1397  list_name = "/ice_mod/tracer/"//trim(tracers(n1)%tracer_name)//"/"//trim(method_type)
1398  case (model_land )
1399  list_name = "/land_mod/tracer/"//trim(tracers(n1)%tracer_name)//"/"//trim(method_type)
1400  case default
1401  list_name = "/default/tracer/"//trim(tracers(n1)%tracer_name)//"/"//trim(method_type)
1402  end select
1403 
1404  name = ''
1405  control_tr = ''
1406  query_method = fm_query_method(list_name, name, control_tr)
1407 
1408  if ( present(control) ) then
1409  if ( len_trim(control_tr)>len(control) ) then
1410  write(chn,*)len(control)
1411  write(chn1,*)len_trim(control_tr)
1412  if(error_handler('query_method', &
1413  ' Output string length ('//trim(adjustl(chn)) &
1414  // ') is not enough to return all "control" parameters ("'//trim(control_tr) &
1415  // '", length='//trim(adjustl(chn1))//')', &
1416  err_msg)) return
1417  endif
1418  control = trim(control_tr)
1419  endif
1420 
1421  end function query_method
1422 !</FUNCTION>
1423 
1424 !<SUBROUTINE NAME="set_tracer_atts">
1425 ! <OVERVIEW>
1426 ! A subroutine to allow the user set the tracer longname and units from the
1427 ! tracer initialization routine.
1428 ! </OVERVIEW>
1429 ! <DESCRIPTION>
1430 ! A function to allow the user set the tracer longname and units from the
1431 ! tracer initialization routine. It seems sensible that the user who is
1432 ! coding the tracer code will know what units they are working in and it
1433 ! is probably safer to set the value in the tracer code rather than in
1434 ! the field table.
1435 ! </DESCRIPTION>
1436 ! <TEMPLATE>
1437 ! call set_tracer_atts(model, name, longname, units)
1438 ! </TEMPLATE>
1439 
1440 ! <IN NAME="model" TYPE="integer">
1441 ! A parameter representing the component model in use.
1442 ! </IN>
1443 ! <IN NAME="name" TYPE="character">
1444 ! Tracer name.
1445 ! </IN>
1446 ! <OUT NAME="longname" TYPE="character, optional">
1447 ! A string describing the longname of the tracer for output to NetCDF files
1448 ! </OUT>
1449 ! <OUT NAME="units" TYPE="character, optional">
1450 ! A string describing the units of the tracer for output to NetCDF files
1451 ! </OUT>
1452 subroutine set_tracer_atts(model, name, longname, units)
1454 integer, intent(in) :: model
1455 character(len=*), intent(in) :: name
1456 character(len=*), intent(in), optional :: longname, units
1457 
1458 integer :: n, index
1459 logical :: success
1460 character(len=128) :: list_name
1461 
1462 if ( get_tracer_index(model,name,n) ) then
1463  tracers(tracer_array(model,n))%tracer_units = units
1464  tracers(tracer_array(model,n))%tracer_longname = longname
1465  select case(model)
1466  case(model_coupler)
1467  list_name = "/coupler_mod/tracer/"//trim(name)
1468  case(model_atmos)
1469  list_name = "/atmos_mod/tracer/"//trim(name)
1470  case(model_ocean)
1471  list_name = "/ocean_mod/tracer/"//trim(name)
1472  case(model_land)
1473  list_name = "/land_mod/tracer/"//trim(name)
1474  case(model_ice)
1475  list_name = "/ice_mod/tracer/"//trim(name)
1476  case DEFAULT
1477  list_name = "/"//trim(name)
1478  end select
1479 
1480 ! Method_type is a list, method_name is a name of a parameter and method_control has the value.
1481 ! list_name = trim(list_name)//"/longname"
1482  if ( fm_exists(list_name)) then
1483  success = fm_change_list(list_name)
1484  if ( present(longname) ) then
1485  if ( longname .ne. "" ) index = fm_new_value('longname',longname)
1486  endif
1487  if ( present(units) ) then
1488  if (units .ne. "" ) index = fm_new_value('units',units)
1489  endif
1490  endif
1491 
1492 else
1493  call mpp_error(note,'set_tracer_atts : Trying to set longname and/or units for non-existent tracer : '//trim(name))
1494 endif
1495 
1496 end subroutine set_tracer_atts
1497 !</SUBROUTINE>
1498 
1499 !<SUBROUTINE NAME="set_tracer_method">
1500 ! <OVERVIEW>
1501 ! A subroutine to allow the user to set some tracer specific methods.
1502 ! </OVERVIEW>
1503 ! <DESCRIPTION>
1504 ! A subroutine to allow the user to set methods for a specific tracer.
1505 ! </DESCRIPTION>
1506 ! <TEMPLATE>
1507 ! call set_tracer_method(model, name, method_type, method_name, method_control)
1508 ! </TEMPLATE>
1509 
1510 ! <IN NAME="model" TYPE="integer">
1511 ! A parameter representing the component model in use.
1512 ! </IN>
1513 ! <IN NAME="name" TYPE="character">
1514 ! Tracer name.
1515 ! </IN>
1516 ! <IN NAME="method_type" TYPE="character">
1517 ! The type of the method to be set.
1518 ! </IN>
1519 ! <IN NAME="method_name" TYPE="character">
1520 ! The name of the method to be set.
1521 ! </IN>
1522 ! <IN NAME="method_control" TYPE="character">
1523 ! The control parameters of the method to be set.
1524 ! </IN>
1525 
1526 subroutine set_tracer_method(model, name, method_type, method_name, method_control)
1528 integer, intent(in) :: model
1529 character(len=*), intent(in) :: name
1530 character(len=*), intent(in) :: method_type
1531 character(len=*), intent(in) :: method_name
1532 character(len=*), intent(in) :: method_control
1533 
1534 integer :: n, num_method, index
1535 logical :: success
1536 character(len=128) :: list_name
1537 
1538 if ( get_tracer_index(model,name,n) ) then
1539  tracers(n)%num_methods = tracers(n)%num_methods + 1
1540  num_method = tracers(n)%num_methods
1541 
1542  select case(model)
1543  case(model_coupler)
1544  list_name = "/coupler_mod/tracer/"//trim(name)
1545  case(model_atmos)
1546  list_name = "/atmos_mod/tracer/"//trim(name)
1547  case(model_ocean)
1548  list_name = "/ocean_mod/tracer/"//trim(name)
1549  case(model_land)
1550  list_name = "/land_mod/tracer/"//trim(name)
1551  case(model_ice)
1552  list_name = "/ice_mod/tracer/"//trim(name)
1553  case DEFAULT
1554  list_name = "/"//trim(name)
1555  end select
1556 
1557  if ( method_control .ne. "" ) then
1558 ! Method_type is a list, method_name is a name of a parameter and method_control has the value.
1559  list_name = trim(list_name)//"/"//trim(method_type)
1560  if ( fm_exists(list_name)) then
1561  success = fm_change_list(list_name)
1562  index = fm_new_value(method_type,method_control)
1563  endif
1564  else
1565  call mpp_error(note,'set_tracer_method : Trying to set a method for non-existent tracer : '//trim(name))
1566  endif
1567 endif
1568 
1569 end subroutine set_tracer_method
1570 !</SUBROUTINE>
1571 
1572 function error_handler(routine, err_msg_local, err_msg)
1573 logical :: error_handler
1574 character(len=*), intent(in) :: routine, err_msg_local
1575 character(len=*), intent(out), optional :: err_msg
1576 
1577 if(present(err_msg)) then
1578  err_msg = err_msg_local
1579  error_handler = .true.
1580 else
1581  call mpp_error(fatal,trim(routine)//': '//trim(err_msg_local))
1582 endif
1583 
1584 end function error_handler
1585 
1586 end module tracer_manager_mod
Definition: fms.F90:20
integer, parameter, public model_ice
subroutine print_tracer_info(model, n)
integer, parameter, public model_atmos
logical function, public get_tracer_name(model, n, name, longname, units, err_msg)
subroutine get_tracer_meta_data(model, num_tracers, num_prog, num_diag)
void error_handler(const char *msg)
Definition: mosaic_util.c:55
subroutine set_tracer_method(model, name, method_type, method_name, method_control)
integer, dimension(num_models) total_tracers
logical function, public adjust_positive_def(model, n, err_msg)
integer, parameter, public model_ocean
integer, parameter, public max_tracer_fields
subroutine, public get_tracer_indices(model, ind, prog_ind, diag_ind, fam_ind)
integer, dimension(num_models) prog_tracers
subroutine, public get_field_info(n, fld_type, fld_name, model, num_methods)
subroutine, public create(self, geom, vars)
Linked list implementation.
Definition: qg_fields.F90:76
integer, parameter, public model_land
logical function, public fm_change_list(name)
logical, dimension(num_models) model_registered
subroutine, public get_number_tracers(model, num_tracers, num_prog, num_diag, num_family)
type(method_type), public default_method
Definition: mpp.F90:39
logical function, public query_method(method_type, model, n, name, control, err_msg)
subroutine, public tracer_manager_end
subroutine, public set_tracer_atts(model, name, longname, units)
logical function, public adjust_mass(model, n, err_msg)
subroutine, public set_tracer_profile(model, n, tracer, err_msg)
type(inst_type), dimension(max_tracer_fields), save instantiations
integer function get_tracer_index_integer(model, name, indices, verbose)
integer, dimension(num_models) diag_tracers
logical function, public check_if_prognostic(model, n, err_msg)
integer, parameter notracer
logical function, public fm_exists(name)
logical function, public fm_query_method(name, method_name, method_control)
integer, parameter max_tracer_method
logical function get_tracer_index_logical(model, name, index, indices, verbose)
character(len=11), dimension(num_models), parameter, public model_names
logical function, public fm_modify_name(oldname, newname)
integer, parameter, public num_models
type(tracer_type), dimension(max_tracer_fields), save tracers
integer, dimension(num_models, max_tracer_fields) tracer_array
subroutine, public get_tracer_names(model, n, name, longname, units, err_msg)
subroutine, public tracer_manager_init
integer function, public fm_copy_list(list_name, suffix, create)
integer, parameter, public model_coupler
integer, parameter, public no_tracer
integer function model_tracer_number(model, n)
subroutine, public get_field_methods(n, methods)
subroutine, public field_manager_init(nfields, table_name)
subroutine, public register_tracers(model, num_tracers, num_prog, num_diag, num_family)