SICOPOLIS V3.0
 All Classes Files Functions Variables Macros
output4.F90
Go to the documentation of this file.
1 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 !
3 ! Subroutine : o u t p u t 4
4 !
5 !> @file
6 !!
7 !! Writing of time-series data of the deep ice cores on file in ASCII format.
8 !!
9 !! @section Copyright
10 !!
11 !! Copyright 2009-2013 Ralf Greve
12 !!
13 !! @section License
14 !!
15 !! This file is part of SICOPOLIS.
16 !!
17 !! SICOPOLIS is free software: you can redistribute it and/or modify
18 !! it under the terms of the GNU General Public License as published by
19 !! the Free Software Foundation, either version 3 of the License, or
20 !! (at your option) any later version.
21 !!
22 !! SICOPOLIS is distributed in the hope that it will be useful,
23 !! but WITHOUT ANY WARRANTY; without even the implied warranty of
24 !! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
25 !! GNU General Public License for more details.
26 !!
27 !! You should have received a copy of the GNU General Public License
28 !! along with SICOPOLIS. If not, see <http://www.gnu.org/licenses/>.
29 !<
30 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
31 
32 !-------------------------------------------------------------------------------
33 !> Writing of time-series data of the deep ice cores on file in ASCII format.
34 !<------------------------------------------------------------------------------
35 subroutine output4(time, dxi, deta, delta_ts, glac_index, z_sl)
36 
37 use sico_types
39 
40 implicit none
41 real(dp), intent(in) :: time, dxi, deta, delta_ts, glac_index, z_sl
42 
43 integer(i4b) :: i, j, n
44 real(dp) :: time_val
45 real(dp), dimension(0:JMAX,0:IMAX) :: field
46 real(dp), dimension(:), allocatable :: h_core, temp_core, &
47  vx_core, vy_core, v_core, &
48  rbx_core, rby_core, rb_core
49 
50 allocate(h_core(n_core), temp_core(n_core), &
51  vx_core(n_core), vy_core(n_core), v_core(n_core), &
52  rbx_core(n_core), rby_core(n_core), rb_core(n_core))
53 
54 !-------- Determination of ice-core data --------
55 
56 do n=1, n_core
57 
58 ! ------ Ice thickness
59 
60  field = h_c + h_t
61  call borehole(field, x_core(n), y_core(n), dxi, deta, 'grid', h_core(n))
62 
63 ! ------ Surface velocity
64 
65  do i=0, imax
66  do j=0, jmax
67  field(j,i) = vx_c(kcmax,j,i)
68  end do
69  end do
70 
71  call borehole(field, x_core(n), y_core(n), dxi, deta, 'sg_x', vx_core(n))
72 
73  do i=0, imax
74  do j=0, jmax
75  field(j,i) = vy_c(kcmax,j,i)
76  end do
77  end do
78 
79  call borehole(field, x_core(n), y_core(n), dxi, deta, 'sg_y', vy_core(n))
80 
81  v_core(n) = sqrt(vx_core(n)**2+vy_core(n)**2)
82 
83 ! ------ Basal temperature
84 
85  do i=0, imax
86  do j=0, jmax
87  field(j,i) = temp_r(krmax,j,i)
88  end do
89  end do
90 
91  call borehole(field, x_core(n), y_core(n), dxi, deta, 'grid', temp_core(n))
92 
93 ! ------ Basal frictional heating
94 
95  do i=0, imax
96  do j=0, jmax
97  field(j,i) = vx_t(0,j,i)*txz_t(0,j,i)
98  end do
99  end do
100 
101  call borehole(field, x_core(n), y_core(n), dxi, deta, 'sg_x', rbx_core(n))
102 
103  do i=0, imax
104  do j=0, jmax
105  field(j,i) = vy_t(0,j,i)*tyz_t(0,j,i)
106  end do
107  end do
108 
109  call borehole(field, x_core(n), y_core(n), dxi, deta, 'sg_y', rby_core(n))
110 
111  rb_core(n) = rbx_core(n) + rby_core(n)
112 
113 end do
114 
115 ! ------ Conversion
116 
117 #if ( !defined(OUT_TIMES) || OUT_TIMES==1 )
118 time_val = time /year_sec ! s -> a
119 #elif OUT_TIMES == 2
120 time_val = (time+year_zero) /year_sec ! s -> a
121 #else
122 stop ' output4: OUT_TIMES must be either 1 or 2!'
123 #endif
124 
125 do n=1, n_core
126  v_core(n) = v_core(n) *year_sec ! m/s -> m/a
127 end do
128 
129 !-------- Writing of data on file --------
130 
131 if (forcing_flag == 1) then
132  write(14,'(1pe13.6,2(1pe13.4))') time_val, delta_ts, z_sl
133 else if (forcing_flag == 2) then
134  write(14,'(1pe13.6,2(1pe13.4))') time_val, glac_index, z_sl
135 else if (forcing_flag == 3) then
136  write(14,'(1pe13.6,2(1pe13.4))') time_val, 1.11e+11_dp, z_sl
137 end if
138 
139 n=1
140 write(14,'(13x,1pe13.4)',advance='no') h_core(n)
141 do n=2, n_core-1
142  write(14,'(1pe13.4)',advance='no') h_core(n)
143 end do
144 n=n_core
145 write(14,'(1pe13.4)') h_core(n)
146 
147 n=1
148 write(14,'(13x,1pe13.4)',advance='no') v_core(n)
149 do n=2, n_core-1
150  write(14,'(1pe13.4)',advance='no') v_core(n)
151 end do
152 n=n_core
153 write(14,'(1pe13.4)') v_core(n)
154 
155 n=1
156 write(14,'(13x,1pe13.4)',advance='no') temp_core(n)
157 do n=2, n_core-1
158  write(14,'(1pe13.4)',advance='no') temp_core(n)
159 end do
160 n=n_core
161 write(14,'(1pe13.4)') temp_core(n)
162 
163 n=1
164 write(14,'(13x,1pe13.4)',advance='no') rb_core(n)
165 do n=2, n_core-1
166  write(14,'(1pe13.4)',advance='no') rb_core(n)
167 end do
168 n=n_core
169 write(14,'(1pe13.4,/)') rb_core(n)
170 
171 deallocate(h_core, vx_core, vy_core, v_core, temp_core, &
172  rbx_core, rby_core, rb_core)
173 
174 end subroutine output4
175 !