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