FMS  2025.03
Flexible Modeling System
tracer_manager.inc
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 !> @defgroup tracer_manager_mod tracer_manager_mod
20 !> @ingroup tracer_manager
21 
22 
23 !> @brief Subroutine to set the tracer field to the wanted profile.
24 !!
25 !> If the profile type is 'fixed' then the tracer field values are set
26 !! equal to the surface value.
27 !! If the profile type is 'profile' then the top/bottom of model and
28 !! surface values are read and an exponential profile is calculated,
29 !! with the profile being dependent on the number of levels in the
30 !! component model. This should be called from the part of the dynamical
31 !! core where tracer restarts are called in the event that a tracer
32 !! restart file does not exist.
33 !!
34 !! This can be activated by adding a method to the field_table
35 !! e.g.
36 !! @verbose "profile_type","fixed","surface_value = 1e-12" @endverbose
37 !! would return values of surf_value = 1e-12 and a multiplier of 1.0
38 !! One can use these to initialize the entire field with a value of 1e-12.
39 !!
40 !! "profile_type","profile","surface_value = 1e-12, top_value = 1e-15"
41 !! In a 15 layer model this would return values of surf_value = 1e-12 and
42 !! multiplier = 0.6309573 i.e 1e-15 = 1e-12*(0.6309573^15)
43 !! In this case the model should be MODEL_ATMOS as you have a "top" value.
44 !!
45 !! If you wish to initialize the ocean model, one can use bottom_value instead
46 !! of top_value.
47 subroutine set_tracer_profile_(model, n, tracer, err_msg)
48 
49 integer, intent(in) :: model !< Parameter representing component model in use
50 integer, intent(in) :: n !< Tracer number
51 real(FMS_TM_KIND_), intent(inout), dimension(:,:,:) :: tracer !< Initialized tracer array
52 character(len=*), intent(out), optional :: err_msg
53 
54 real(FMS_TM_KIND_) :: surf_value, multiplier
55 integer :: numlevels, k, n1, flag
56 real(FMS_TM_KIND_) :: top_value, bottom_value
57 character(len=80) :: scheme, control,profile_type
58 character(len=128) :: err_msg_local
59 character(len=11) :: chn
60 
61 integer, parameter :: lkind=fms_tm_kind_
62 
63 if(.not.module_is_initialized) call tracer_manager_init
64 
65 if (n < 1 .or. n > total_tracers(model)) then
66  write(chn, '(i11)') n
67  err_msg_local = ' Invalid tracer index. Model name = '//trim(model_names(model))//', Index='//trim(chn)
68  if(error_handler('set_tracer_profile', err_msg_local, err_msg)) return
69 endif
70 n1 = tracer_array(model,n)
71 
72 !default values
73 profile_type = 'Fixed'
74 surf_value = 0.0e+00_lkind
75 top_value = surf_value
76 bottom_value = surf_value
77 multiplier = 1.0_lkind
78 
79 tracer = surf_value
80 
81 if ( query_method( 'profile_type',model,n,scheme,control)) then
82 !Change the tracer_number to the tracer_manager version
83 
84  if(lowercase(trim(scheme(1:5))).eq.'fixed') then
85  profile_type = 'Fixed'
86  flag =parse(control,'surface_value',surf_value)
87  multiplier = 1.0_lkind
88  tracer = surf_value
89  endif
90 
91  if(lowercase(trim(scheme(1:7))).eq.'profile') then
92  profile_type = 'Profile'
93  flag=parse(control,'surface_value',surf_value)
94  if (surf_value .eq. 0.0_lkind) &
95  call mpp_error(fatal,'set_tracer_profile : Cannot have a zero surface value for an exponential profile. Tracer '&
96  //tracers(n1)%tracer_name//" "//control//" "//scheme)
97  select case (tracers(n1)%model)
98  case (model_atmos)
99  flag=parse(control,'top_value',top_value)
100  if(mpp_pe()==mpp_root_pe() .and. flag == 0) &
101  call mpp_error(note,'set_tracer_profile : Parameter top_value needs to be defined for the tracer profile.')
102  case (model_ocean)
103  flag =parse(control,'bottom_value',bottom_value)
104  if(mpp_pe() == mpp_root_pe() .and. flag == 0) &
105  call mpp_error(note, &
106  & 'set_tracer_profile : Parameter bottom_value needs to be defined for the tracer profile.')
107  case default
108 ! Should there be a NOTE or WARNING message here?
109  end select
110 
111 ! If profile type is profile then set the surface value to the input
112 ! value and calculate the vertical multiplier.
113 !
114 ! Assume an exponential decay/increase from the surface to the top level
115 ! C = C0 exp ( -multiplier* level_number)
116 ! => multiplier = exp [ ln(Ctop/Csurf)/number_of_levels]
117 !
118 numlevels = size(tracer,3) -1
119  select case (tracers(n1)%model)
120  case (model_atmos)
121  multiplier = exp( log(top_value/surf_value) /real(numlevels,lkind))
122  tracer(:,:,1) = surf_value
123  do k = 2, size(tracer,3)
124  tracer(:,:,k) = tracer(:,:,k-1) * multiplier
125  enddo
126  case (model_ocean)
127  multiplier = exp( log(bottom_value/surf_value) / real(numlevels,lkind))
128  tracer(:,:,size(tracer,3)) = surf_value
129  do k = size(tracer,3) - 1, 1, -1
130  tracer(:,:,k) = tracer(:,:,k+1) * multiplier
131  enddo
132  case default
133  end select
134  endif !scheme.eq.profile
135 
136  if (mpp_pe() == mpp_root_pe() ) write(*,700) 'Tracer ',trim(tracers(n1)%tracer_name), &
137  ' initialized with surface value of ',surf_value, &
138  ' and vertical multiplier of ',multiplier
139  700 FORMAT (3a,e13.6,a,f13.6)
140 
141 endif ! end of query scheme
142 
143 end subroutine set_tracer_profile_
144 
145 !#######################################################################
146 !> @}
147 ! close documentation grouping
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:407