SICOPOLIS V3.2
 All Classes Files Functions Variables Macros
sico_init.F90
Go to the documentation of this file.
1 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 !
3 ! Subroutine : s i c o _ i n i t
4 !
5 !> @file
6 !!
7 !! Initialisations for SICOPOLIS.
8 !!
9 !! @section Copyright
10 !!
11 !! Copyright 2009-2016 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 !> Initialisations for SICOPOLIS.
34 !<------------------------------------------------------------------------------
35 subroutine sico_init(delta_ts, glac_index, &
36  mean_accum, &
37  dtime, dtime_temp, dtime_wss, dtime_out, dtime_ser, &
38  time, time_init, time_end, time_output, &
39  dxi, deta, dzeta_c, dzeta_t, dzeta_r, &
40  z_sl, dzsl_dtau, z_mar, &
41  ii, jj, nn, &
42  ndat2d, ndat3d, n_output, &
43  runname)
44 
45 use sico_types
47 use sico_vars
51 
52 implicit none
53 
54 integer(i4b), intent(out) :: ii((imax+1)*(jmax+1)), &
55  jj((imax+1)*(jmax+1)), &
56  nn(0:jmax,0:imax)
57 integer(i4b), intent(out) :: ndat2d, ndat3d
58 integer(i4b), intent(out) :: n_output
59 real(dp), intent(out) :: delta_ts, glac_index
60 real(dp), intent(out) :: mean_accum
61 real(dp), intent(out) :: dtime, dtime_temp, dtime_wss, &
62  dtime_out, dtime_ser
63 real(dp), intent(out) :: time, time_init, time_end, time_output(100)
64 real(dp), intent(out) :: dxi, deta, dzeta_c, dzeta_t, dzeta_r
65 real(dp), intent(out) :: z_sl, dzsl_dtau, z_mar
66 character(len=100), intent(out) :: runname
67 
68 integer(i4b) :: i, j, kc, kt, kr, m, n, ir, jr
69 integer(i4b) :: ios, ios1, ios2, ios3, ios4
70 integer(i4b) :: ierr
71 real(dp) :: dtime0, dtime_temp0, dtime_wss0, dtime_out0, dtime_ser0
72 real(dp) :: time_init0, time_end0, time_output0(100)
73 real(dp) :: d_dummy
74 character (len=100) :: anfdatname
75 character (len=256) :: shell_command
76 character :: ch_dummy
77 logical :: flag_3d_output
78 
79 character(len=64), parameter :: fmt1 = '(a)', &
80  fmt2 = '(a,i4)', &
81  fmt3 = '(a,es12.4)'
82 
83 character(len= 8) :: ch_imax
84 character(len=128) :: fmt4
85 
86 write(ch_imax, fmt='(i8)') imax
87 write(fmt4, fmt='(a)') '('//trim(adjustl(ch_imax))//'(i1),i1)'
88 
89 !-------- Name of the computational domain --------
90 
91 #if defined(ANT)
92 ch_domain_long = 'Antarctica'
93 ch_domain_short = 'ant'
94 
95 #elif defined(ASF)
96 ch_domain_long = 'Austfonna'
97 ch_domain_short = 'asf'
98 
99 #elif defined(EMTP2SGE)
100 ch_domain_long = 'EISMINT Phase 2 Simplified Geometry Experiment'
101 ch_domain_short = 'emtp2sge'
102 
103 #elif defined(GRL)
104 ch_domain_long = 'Greenland'
105 ch_domain_short = 'grl'
106 
107 #elif defined(NHEM)
108 ch_domain_long = 'Northern hemisphere'
109 ch_domain_short = 'nhem'
110 
111 #elif defined(SCAND)
112 ch_domain_long = 'Scandinavia and Eurasia'
113 ch_domain_short = 'scand'
114 
115 #elif defined(TIBET)
116 ch_domain_long = 'Tibet'
117 ch_domain_short = 'tibet'
118 
119 #elif defined(NMARS)
120 ch_domain_long = 'North polar cap of Mars'
121 ch_domain_short = 'nmars'
122 
123 #elif defined(SMARS)
124 ch_domain_long = 'South polar cap of Mars'
125 ch_domain_short = 'smars'
126 
127 #elif (defined(XYZ))
128 ch_domain_long = 'XYZ'
129 ch_domain_short = 'xyz'
130 #if (defined(HEINO))
131 ch_domain_long = trim(ch_domain_long)//'/ISMIP HEINO'
132 #endif
133 
134 #else
135 stop ' sico_init: No valid domain specified!'
136 #endif
137 
138 !-------- Some initial values --------
139 
140 n_output = 0
141 
142 dtime = 0.0_dp
143 dtime_temp = 0.0_dp
144 dtime_wss = 0.0_dp
145 dtime_out = 0.0_dp
146 dtime_ser = 0.0_dp
147 
148 time = 0.0_dp
149 time_init = 0.0_dp
150 time_end = 0.0_dp
151 time_output = 0.0_dp
152 
153 !-------- Initialisation of the Library of Iterative Solvers Lis,
154 ! if required --------
155 
156 #if (CALCZS==4 || MARGIN==3 || DYNAMICS==2)
157 
158 #include "lisf.h"
159 call lis_initialize(ierr)
160 
161 #endif
162 
163 !-------- Read physical parameters --------
164 
165 call phys_para()
166 
167 call ice_mat_eqs_pars(rf, r_t, kappa, c, -190, 10)
168 
169 temp_min = temp_min
170 s_t = s_t *1.0e-09_dp ! K/km3 -> K/m3
171 x_hat = x_hat *1.0e+03_dp ! km -> m
172 y_hat = y_hat *1.0e+03_dp ! km -> m
173 rad = rad *1.0e+03_dp ! km -> m
174 b_min = b_min /year_sec ! m/a -> m/s
175 b_max = b_max /year_sec ! m/a -> m/s
176 
177 ! ------ Some auxiliary quantities required for the enthalpy method
178 
179 call calc_c_int_table(c, -190, 10, l)
181 
182 !-------- Compatibility check of the SICOPOLIS version with the header file
183 
184 if ( trim(version) /= trim(sico_version) ) &
185  stop ' sico_init: SICOPOLIS version not compatible with header file!'
186 
187 !-------- Check whether the dynamics and thermodynamics modes are defined
188 
189 #if (!defined(DYNAMICS))
190 stop ' sico_init: DYNAMICS not defined in the header file!'
191 #endif
192 
193 #if (!defined(CALCMOD))
194 stop ' sico_init: CALCMOD not defined in the header file!'
195 #endif
196 
197 #if (defined(ENTHMOD))
198 write(6, fmt='(a)') ' sico_init: ENTHMOD must not be defined any more.'
199 write(6, fmt='(a)') ' Please update your header file!'
200 stop
201 #endif
202 
203 !-------- Compatibility check of the horizontal resolution with the
204 ! number of grid points --------
205 
206 #if (GRID==0)
207 
208 if (abs(dx-50.0_dp) < eps) then
209 
210  if ((imax /= 80).or.(jmax /= 80)) then
211  stop ' sico_init: IMAX and/or JMAX wrong!'
212  end if
213 
214 else if (abs(dx-25.0_dp) < eps) then
215 
216  if ((imax /= 160).or.(jmax /= 160)) then
217  stop ' sico_init: IMAX and/or JMAX wrong!'
218  end if
219 
220 else
221  stop ' sico_init: DX wrong!'
222 end if
223 
224 #elif (GRID==1)
225 
226 stop ' sico_init: GRID==1 not allowed for this application!'
227 
228 #elif (GRID==2)
229 
230 stop ' sico_init: GRID==2 not allowed for this application!'
231 
232 #endif
233 
234 !-------- Compatibility check of the thermodynamics mode
235 ! (cold vs. polythermal vs. enthalpy method)
236 ! and the number of grid points in the lower (kt) ice domain --------
237 
238 #if (CALCMOD==0 || CALCMOD==2 || CALCMOD==3 || CALCMOD==-1)
239 
240 if (ktmax /= 1) then
241  write(6, fmt='(a)') ' sico_init: For options CALCMOD==0, 2, 3 or -1,'
242  write(6, fmt='(a)') ' the separate kt domain is redundant.'
243  write(6, fmt='(a)') ' Therefore, consider setting KTMAX to 1.'
244  write(6, fmt='(a)') ' '
245 end if
246 
247 #endif
248 
249 !-------- Compatibility check of discretization schemes for the horizontal and
250 ! vertical advection terms in the temperature and age equations --------
251 
252 #if (ADV_HOR==1)
253 stop ' sico_init: Option ADV_HOR==1 (central differences) not defined!'
254 #endif
255 
256 !-------- Check whether for the shallow shelf
257 ! or shelfy stream approximation
258 ! the chosen grid is Cartesian coordinates
259 ! without distortion correction (GRID==0) --------
260 
261 #if ((MARGIN==3 && DYNAMICS==1) || DYNAMICS==2) /* requires SSA or SStA */
262 #if (GRID != 0)
263 write(6, fmt='(a)') ' sico_init: Distortion correction not implemented'
264 write(6, fmt='(a)') ' for the shallow shelf approximation (SSA)'
265 write(6, fmt='(a)') ' or the shelfy stream approximation (SStA)'
266 write(6, fmt='(a)') ' -> GRID==0 required!'
267 stop
268 #endif
269 #endif
270 
271 !-------- Setting of forcing flag --------
272 
273 forcing_flag = 1 ! forcing by delta_ts
274 
275 !-------- Initialization of numerical time steps --------
276 
277 dtime0 = dtime0
278 dtime_temp0 = dtime_temp0
279 
280 !-------- Further initializations --------
281 
282 dzeta_c = 1.0_dp/real(kcmax,dp)
283 dzeta_t = 1.0_dp/real(ktmax,dp)
284 dzeta_r = 1.0_dp/real(krmax,dp)
285 
286 ndat2d = 1
287 ndat3d = 1
288 
289 !-------- General abbreviations --------
290 
291 ! ------ kc domain
292 
293 if (deform >= eps) then
294 
295  flag_aa_nonzero = .true. ! non-equidistant grid
296 
297  aa = deform
298  ea = exp(aa)
299 
300  kc=0
301  zeta_c(kc) = 0.0_dp
302  eaz_c(kc) = 1.0_dp
303  eaz_c_quotient(kc) = 0.0_dp
304 
305  do kc=1, kcmax-1
306  zeta_c(kc) = kc*dzeta_c
307  eaz_c(kc) = exp(aa*zeta_c(kc))
308  eaz_c_quotient(kc) = (eaz_c(kc)-1.0_dp)/(ea-1.0_dp)
309  end do
310 
311  kc=kcmax
312  zeta_c(kc) = 1.0_dp
313  eaz_c(kc) = exp(aa)
314  eaz_c_quotient(kc) = 1.0_dp
315 
316 else
317 
318  flag_aa_nonzero = .false. ! equidistant grid
319 
320  aa = 0.0_dp
321  ea = 1.0_dp
322 
323  kc=0
324  zeta_c(kc) = 0.0_dp
325  eaz_c(kc) = 1.0_dp
326  eaz_c_quotient(kc) = 0.0_dp
327 
328  do kc=1, kcmax-1
329  zeta_c(kc) = kc*dzeta_c
330  eaz_c(kc) = 1.0_dp
331  eaz_c_quotient(kc) = zeta_c(kc)
332  end do
333 
334  kc=kcmax
335  zeta_c(kc) = 1.0_dp
336  eaz_c(kc) = 1.0_dp
337  eaz_c_quotient(kc) = 1.0_dp
338 
339 end if
340 
341 ! ------ kt domain
342 
343 kt=0
344 zeta_t(kt) = 0.0_dp
345 
346 do kt=1, ktmax-1
347  zeta_t(kt) = kt*dzeta_t
348 end do
349 
350 kt=ktmax
351 zeta_t(kt) = 1.0_dp
352 
353 ! ------ kr domain
354 
355 kr=0
356 zeta_r(kr) = 0.0_dp
357 
358 do kr=1, krmax-1
359  zeta_r(kr) = kr*dzeta_r
360 end do
361 
362 kr=krmax
363 zeta_r(kr) = 1.0_dp
364 
365 !-------- Construction of a vector (with index n) from a 2-d array
366 ! (with indices i, j) by diagonal numbering --------
367 
368 n=1
369 
370 do m=0, imax+jmax
371  do i=m, 0, -1
372  j = m-i
373  if ((i <= imax).and.(j <= jmax)) then
374  ii(n) = i
375  jj(n) = j
376  nn(j,i) = n
377  n=n+1
378  end if
379  end do
380 end do
381 
382 !-------- Specification of current simulation --------
383 
384 runname = runname
385 anfdatname = anfdatname
386 
387 #if (defined(YEAR_ZERO))
388 year_zero = year_zero
389 #else
390 year_zero = 2000.0_dp ! default value 2000 CE
391 #endif
392 
393 time_init0 = time_init0
394 time_end0 = time_end0
395 dtime_ser0 = dtime_ser0
396 
397 #if (OUTPUT==1 || OUTPUT==3)
398 dtime_out0 = dtime_out0
399 #endif
400 
401 #if (OUTPUT==2 || OUTPUT==3)
402 n_output = n_output
403 time_output0( 1) = time_out0_01
404 time_output0( 2) = time_out0_02
405 time_output0( 3) = time_out0_03
406 time_output0( 4) = time_out0_04
407 time_output0( 5) = time_out0_05
408 time_output0( 6) = time_out0_06
409 time_output0( 7) = time_out0_07
410 time_output0( 8) = time_out0_08
411 time_output0( 9) = time_out0_09
412 time_output0(10) = time_out0_10
413 time_output0(11) = time_out0_11
414 time_output0(12) = time_out0_12
415 time_output0(13) = time_out0_13
416 time_output0(14) = time_out0_14
417 time_output0(15) = time_out0_15
418 time_output0(16) = time_out0_16
419 time_output0(17) = time_out0_17
420 time_output0(18) = time_out0_18
421 time_output0(19) = time_out0_19
422 time_output0(20) = time_out0_20
423 #endif
424 
425 !-------- Write log file --------
426 
427 shell_command = 'if [ ! -d'
428 shell_command = trim(shell_command)//' '//outpath
429 shell_command = trim(shell_command)//' '//'] ; then mkdir'
430 shell_command = trim(shell_command)//' '//outpath
431 shell_command = trim(shell_command)//' '//'; fi'
432 call system(trim(shell_command))
433  ! Check whether directory OUTPATH exists. If not, it is created.
434 
435 open(10, iostat=ios, &
436  file=outpath//'/'//trim(runname)//'.log', &
437  status='new')
438 if (ios /= 0) stop ' sico_init: Error when opening the log file!'
439 
440 write(10, fmt=trim(fmt1)) 'Computational domain:'
441 write(10, fmt=trim(fmt1)) trim(ch_domain_long)
442 write(10, fmt=trim(fmt1)) ' '
443 
444 write(10, fmt=trim(fmt2)) 'imax =', imax
445 write(10, fmt=trim(fmt2)) 'jmax =', jmax
446 write(10, fmt=trim(fmt2)) 'kcmax =', kcmax
447 write(10, fmt=trim(fmt2)) 'ktmax =', ktmax
448 write(10, fmt=trim(fmt2)) 'krmax =', krmax
449 write(10, fmt=trim(fmt1)) ' '
450 
451 write(10, fmt=trim(fmt3)) 'a =', aa
452 write(10, fmt=trim(fmt1)) ' '
453 
454 #if (GRID==0 || GRID==1)
455 write(10, fmt=trim(fmt3)) 'x0 =', x0
456 write(10, fmt=trim(fmt3)) 'y0 =', y0
457 write(10, fmt=trim(fmt3)) 'dx =', dx
458 #elif (GRID==2)
459 write(10, fmt=trim(fmt3)) 'dlambda =', dlambda
460 write(10, fmt=trim(fmt3)) 'dphi =', dphi
461 #endif
462 write(10, fmt=trim(fmt1)) ' '
463 
464 write(10, fmt=trim(fmt3)) 'year_zero =', year_zero
465 write(10, fmt=trim(fmt3)) 'time_init =', time_init0
466 write(10, fmt=trim(fmt3)) 'time_end =', time_end0
467 write(10, fmt=trim(fmt3)) 'dtime =', dtime0
468 write(10, fmt=trim(fmt3)) 'dtime_temp =', dtime_temp0
469 #if (OUTPUT==1 || OUTPUT==3)
470 write(10, fmt=trim(fmt3)) 'dtime_out =', dtime_out0
471 #endif
472 write(10, fmt=trim(fmt3)) 'dtime_ser =', dtime_ser0
473 write(10, fmt=trim(fmt1)) ' '
474 
475 write(10, fmt=trim(fmt2)) 'ANF_DAT =', anf_dat
476 write(10, fmt=trim(fmt1)) 'mask_present file = '//mask_present_file
477 #if (ANF_DAT==1 && defined(TEMP_INIT))
478 write(10, fmt=trim(fmt2)) 'TEMP_INIT =', temp_init
479 #endif
480 #if (ANF_DAT==3)
481 write(10, fmt=trim(fmt1)) 'Initial-value file = '//anfdatname
482 write(10, fmt=trim(fmt1)) 'Path to initial-value file = '//anfdatpath
483 #endif
484 write(10, fmt=trim(fmt1)) ' '
485 
486 write(10, fmt=trim(fmt1)) 'Physical-parameter file = '//phys_para_file
487 write(10, fmt=trim(fmt1)) ' '
488 
489 #if (CALCZS==3 || CALCZS==4)
490 write(10, fmt=trim(fmt3)) 'ovi_weight =', ovi_weight
491 #if (CALCZS==3)
492 write(10, fmt=trim(fmt3)) 'omega_sor =', omega_sor
493 #endif
494 write(10, fmt=trim(fmt1)) ' '
495 #endif
496 
497 write(10, fmt=trim(fmt3)) 'temp_min =', temp_min
498 write(10, fmt=trim(fmt3)) 's_t =', s_t
499 write(10, fmt=trim(fmt3)) 'x_hat =', x_hat
500 write(10, fmt=trim(fmt3)) 'y_hat =', y_hat
501 write(10, fmt=trim(fmt3)) 'rad =', rad
502 write(10, fmt=trim(fmt3)) 'b_min =', b_min
503 write(10, fmt=trim(fmt3)) 'b_max =', b_max
504 write(10, fmt=trim(fmt1)) ' '
505 
506 #if (TSURFACE==1)
507 write(10, fmt=trim(fmt3)) 'delta_ts0 =', delta_ts0
508 #elif (TSURFACE==3)
509 write(10, fmt=trim(fmt3)) 'sine_amplit =', sine_amplit
510 write(10, fmt=trim(fmt3)) 'sine_period =', sine_period
511 #elif (TSURFACE==4)
512 write(10, fmt=trim(fmt1)) 'GRIP file = '//grip_temp_file
513 write(10, fmt=trim(fmt3)) 'grip_temp_fact =', grip_temp_fact
514 #endif
515 
516 #if (SEA_LEVEL==1)
517 write(10, fmt=trim(fmt3)) 'z_sl0 =', z_sl0
518 #elif (SEA_LEVEL==3)
519 write(10, fmt=trim(fmt1)) 'sea-level file = '//sea_level_file
520 #endif
521 write(10, fmt=trim(fmt1)) ' '
522 
523 #if (MARGIN==2)
524 #if (MARINE_ICE_CALVING==2 || MARINE_ICE_CALVING==3)
525 write(10, fmt=trim(fmt3)) 'z_mar =', z_mar
526 write(10, fmt=trim(fmt1)) ' '
527 #elif (MARINE_ICE_CALVING==4 || MARINE_ICE_CALVING==5 \
528  || marine_ice_calving==6 || marine_ice_calving==7)
529 write(10, fmt=trim(fmt3)) 'fact_z_mar =', fact_z_mar
530 write(10, fmt=trim(fmt1)) ' '
531 #elif (MARINE_ICE_FORMATION==2 && MARINE_ICE_CALVING==9)
532 write(10, fmt=trim(fmt3)) 'calv_uw_coeff =', calv_uw_coeff
533 write(10, fmt=trim(fmt3)) 'r1_calv_uw =', r1_calv_uw
534 write(10, fmt=trim(fmt3)) 'r2_calv_uw =', r2_calv_uw
535 write(10, fmt=trim(fmt1)) ' '
536 #endif
537 #elif (MARGIN==3)
538 #if (ICE_SHELF_CALVING==2)
539 write(10, fmt=trim(fmt3)) 'H_calv =', h_calv
540 write(10, fmt=trim(fmt1)) ' '
541 #endif
542 #endif
543 
544 write(10, fmt=trim(fmt1)) 'Sediment-mask file = '//mask_sedi_file
545 write(10, fmt=trim(fmt1)) ' '
546 
547 write(10, fmt=trim(fmt3)) 'c_slide =', c_slide
548 write(10, fmt=trim(fmt3)) 'gamma_slide =', gamma_slide
549 write(10, fmt=trim(fmt2)) 'p_weert =', p_weert
550 write(10, fmt=trim(fmt2)) 'q_weert =', q_weert
551 #if (SLIDE_LAW==2)
552 write(10, fmt=trim(fmt3)) 'red_pres_limit_fact =', red_pres_limit_fact
553 #endif
554 write(10, fmt=trim(fmt1)) ' '
555 
556 write(10, fmt=trim(fmt3)) 'c_slide_sedi =', c_slide_sedi
557 write(10, fmt=trim(fmt3)) 'gamma_slide_sedi =', gamma_slide_sedi
558 write(10, fmt=trim(fmt2)) 'p_weert_sedi =', p_weert_sedi
559 write(10, fmt=trim(fmt2)) 'q_weert_sedi =', q_weert_sedi
560 write(10, fmt=trim(fmt1)) ' '
561 
562 #if (Q_GEO_MOD==1)
563 write(10, fmt=trim(fmt3)) 'q_geo =', q_geo
564 write(10, fmt=trim(fmt1)) ' '
565 #endif
566 
567 write(10, fmt=trim(fmt2)) 'REBOUND =', rebound
568 #if (REBOUND==1)
569 write(10, fmt=trim(fmt3)) 'frac_llra =', frac_llra
570 #endif
571 #if (REBOUND==1 || REBOUND==2)
572 write(10, fmt=trim(fmt2)) 'TIME_LAG_MOD =', time_lag_mod
573 #if (TIME_LAG_MOD==1)
574 write(10, fmt=trim(fmt3)) 'time_lag =', time_lag
575 #elif (TIME_LAG_MOD==2)
576 write(10, fmt=trim(fmt1)) 'time_lag_file = '//time_lag_file
577 #else
578 stop ' sico_init: TIME_LAG_MOD must be either 1 or 2!'
579 #endif
580 #endif
581 #if (REBOUND==2)
582 write(10, fmt=trim(fmt2)) 'FLEX_RIG_MOD =', flex_rig_mod
583 #if (FLEX_RIG_MOD==1)
584 write(10, fmt=trim(fmt3)) 'flex_rig =', flex_rig
585 #elif (FLEX_RIG_MOD==2)
586 write(10, fmt=trim(fmt1)) 'flex_rig_file = '//flex_rig_file
587 #else
588 stop ' sico_init: FLEX_RIG_MOD must be either 1 or 2!'
589 #endif
590 #endif
591 write(10, fmt=trim(fmt2)) 'Q_LITHO =', q_litho
592 write(10, fmt=trim(fmt1)) ' '
593 
594 #if (FLOW_LAW==2)
595 write(10, fmt=trim(fmt3)) 'gr_size =', gr_size
596 write(10, fmt=trim(fmt1)) ' '
597 #endif
598 #if (FIN_VISC==2)
599 write(10, fmt=trim(fmt3)) 'sigma_res =', sigma_res
600 write(10, fmt=trim(fmt1)) ' '
601 #endif
602 
603 write(10, fmt=trim(fmt2)) 'ENHMOD =', enhmod
604 #if (ENHMOD==1 || ENHMOD==2 || ENHMOD==3)
605 write(10, fmt=trim(fmt3)) 'enh_fact =', enh_fact
606 #endif
607 #if (ENHMOD==2 || ENHMOD==3)
608 write(10, fmt=trim(fmt3)) 'enh_intg =', enh_intg
609 #endif
610 #if (ENHMOD==2)
611 write(10, fmt=trim(fmt3)) 'age_trans =', age_trans_0
612 #endif
613 #if (ENHMOD==3)
614 write(10, fmt=trim(fmt3)) 'date_trans1 =', date_trans1_0
615 write(10, fmt=trim(fmt3)) 'date_trans2 =', date_trans2_0
616 write(10, fmt=trim(fmt3)) 'date_trans3 =', date_trans3_0
617 #endif
618 #if (ENHMOD==4 || ENHMOD==5)
619 write(10, fmt=trim(fmt3)) 'enh_compr =', enh_compr
620 write(10, fmt=trim(fmt3)) 'enh_shear =', enh_shear
621 #endif
622 #if ((ENHMOD==1 || ENHMOD==2 || ENHMOD==3 || ENHMOD==4) && MARGIN==3)
623 write(10, fmt=trim(fmt3)) 'enh_shelf =', enh_shelf
624 #endif
625 write(10, fmt=trim(fmt1)) ' '
626 
627 write(10, fmt=trim(fmt3)) 'numdiff_H_t =', numdiff_h_t
628 write(10, fmt=trim(fmt3)) 'tau_cts =', tau_cts
629 write(10, fmt=trim(fmt3)) 'H_min =', h_min
630 write(10, fmt=trim(fmt3)) 'vh_max =', vh_max
631 write(10, fmt=trim(fmt3)) 'hd_min =', hd_min
632 write(10, fmt=trim(fmt3)) 'hd_max =', hd_max
633 #if (defined(QBM_MIN))
634 write(10, fmt=trim(fmt3)) 'qbm_min =', qbm_min
635 #elif (defined(QB_MIN)) /* obsolete */
636 write(10, fmt=trim(fmt3)) 'qb_min =', qb_min
637 #endif
638 #if (defined(QBM_MAX))
639 write(10, fmt=trim(fmt3)) 'qbm_max =', qbm_max
640 #elif (defined(QB_MAX)) /* obsolete */
641 write(10, fmt=trim(fmt3)) 'qb_max =', qb_max
642 #endif
643 write(10, fmt=trim(fmt3)) 'age_min =', age_min
644 write(10, fmt=trim(fmt3)) 'age_max =', age_max
645 write(10, fmt=trim(fmt3)) 'mean_accum =', mean_accum
646 #if (ADV_VERT==1)
647 write(10, fmt=trim(fmt3)) 'age_diff =', agediff
648 #endif
649 write(10, fmt=trim(fmt1)) ' '
650 
651 write(10, fmt=trim(fmt2)) 'GRID =', grid
652 write(10, fmt=trim(fmt2)) 'DYNAMICS =', dynamics
653 #if (DYNAMICS==2 && defined(RATIO_SL_THRESH))
654 write(10, fmt=trim(fmt3)) 'ratio_sl_thresh =', ratio_sl_thresh
655 #endif
656 write(10, fmt=trim(fmt2)) 'CALCMOD =', calcmod
657 #if (CALCMOD==-1 && defined(TEMP_CONST))
658 write(10, fmt=trim(fmt3)) 'TEMP_CONST =', temp_const
659 #endif
660 #if (CALCMOD==-1 && defined(AGE_CONST))
661 write(10, fmt=trim(fmt3)) 'AGE_CONST =', age_const
662 #endif
663 #if (CALCMOD==1 && defined(CTS_MELTING_FREEZING))
664 write(10, fmt=trim(fmt2)) 'CTS_MELTING_FREEZING =', cts_melting_freezing
665 #endif
666 write(10, fmt=trim(fmt2)) 'FLOW_LAW =', flow_law
667 write(10, fmt=trim(fmt2)) 'FIN_VISC =', fin_visc
668 write(10, fmt=trim(fmt2)) 'MARGIN =', margin
669 #if (MARGIN==2)
670 write(10, fmt=trim(fmt2)) 'MARINE_ICE_FORMATION =', marine_ice_formation
671 write(10, fmt=trim(fmt2)) 'MARINE_ICE_CALVING =', marine_ice_calving
672 #elif (MARGIN==3)
673 write(10, fmt=trim(fmt2)) 'ICE_SHELF_CALVING =', ice_shelf_calving
674 #endif
675 write(10, fmt=trim(fmt2)) 'ZS_EVOL =', zs_evol
676 write(10, fmt=trim(fmt2)) 'CALCZS =', calczs
677 write(10, fmt=trim(fmt2)) 'ADV_HOR =', adv_hor
678 write(10, fmt=trim(fmt2)) 'ADV_VERT =', adv_vert
679 write(10, fmt=trim(fmt2)) 'TOPOGRAD =', topograd
680 write(10, fmt=trim(fmt2)) 'TSURFACE =', tsurface
681 write(10, fmt=trim(fmt2)) 'SEA_LEVEL =', sea_level
682 write(10, fmt=trim(fmt2)) 'SLIDE_LAW =', slide_law
683 write(10, fmt=trim(fmt2)) 'Q_GEO_MOD =', q_geo_mod
684 write(10, fmt=trim(fmt1)) ' '
685 
686 #if (defined(OUT_TIMES))
687 write(10, fmt=trim(fmt2)) 'OUT_TIMES =', out_times
688 #endif
689 write(10, fmt=trim(fmt2)) 'OUTPUT =', output
690 write(10, fmt=trim(fmt2)) 'OUTSER =', outser
691 #if (defined(OUTSER_V_A))
692 write(10, fmt=trim(fmt2)) 'OUTSER_V_A =', outser_v_a
693 #endif
694 #if (OUTPUT==1 || OUTPUT==2)
695 write(10, fmt=trim(fmt2)) 'ERGDAT =', ergdat
696 #endif
697 write(10, fmt=trim(fmt2)) 'NETCDF =', netcdf
698 write(10, fmt=trim(fmt1)) ' '
699 #if (OUTPUT==2 || OUTPUT==3)
700 write(10, fmt=trim(fmt2)) 'n_output =', n_output
701 do n=1, n_output
702  write(10, fmt=trim(fmt3)) 'time_output =', time_output0(n)
703 end do
704 write(10, fmt=trim(fmt1)) ' '
705 #endif
706 
707 write(10, fmt=trim(fmt1)) 'Program version and date: '//version//' / '//date
708 
709 close(10, status='keep')
710 
711 !-------- Conversion of time quantities --------
712 
713 year_zero = year_zero*year_sec ! a --> s
714 time_init = time_init0*year_sec ! a --> s
715 time_end = time_end0*year_sec ! a --> s
716 dtime = dtime0*year_sec ! a --> s
717 dtime_temp = dtime_temp0*year_sec ! a --> s
718 dtime_ser = dtime_ser0*year_sec ! a --> s
719 #if (OUTPUT==1 || OUTPUT==3)
720 dtime_out = dtime_out0*year_sec ! a --> s
721 #endif
722 #if (OUTPUT==2 || OUTPUT==3)
723 do n=1, n_output
724  time_output(n) = time_output0(n)*year_sec ! a --> s
725 end do
726 #endif
727 
728 if (abs(dtime_temp/dtime-nint(dtime_temp/dtime)) > eps) &
729  stop ' sico_init: dtime_temp must be a multiple of dtime!'
730 
731 time = time_init
732 
733 !-------- Mean accumulation --------
734 
735 mean_accum = mean_accum*(1.0e-03_dp/year_sec)*(rho_w/rho)
736 ! ! mm/a water equiv. --> m/s ice equiv.
737 
738 !-------- Read present topography mask --------
739 
740 open(24, iostat=ios, &
741  file=inpath//'/'//trim(ch_domain_short)//'/'//mask_present_file, &
742  recl=2048, status='old')
743 
744 if (ios /= 0) stop ' sico_init: Error when opening the mask file!'
745 
746 do n=1, 6; read(24, fmt='(a)') ch_dummy; end do
747 
748 do j=jmax, 0, -1
749  read(24, fmt=trim(fmt4)) (maske_help(j,i), i=0,imax)
750 end do
751 
752 close(24, status='keep')
753 
754 !-------- Read rock/sediment mask --------
755 
756 open(24, iostat=ios, &
757  file=inpath//'/'//trim(ch_domain_short)//'/'//mask_sedi_file, &
758  recl=2048, status='old')
759 
760 if (ios /= 0) stop ' sico_init: Error when opening the rock/sediment mask file!'
761 
762 do n=1, 6; read(24, fmt='(a)') ch_dummy; end do
763 
764 do j=jmax, 0, -1
765  read(24, fmt=trim(fmt4)) (maske_sedi(j,i), i=0,imax)
766 end do
767 
768 close(24, status='keep')
769 
770 !-------- Read data for delta_ts --------
771 
772 #if (TSURFACE==4)
773 
774 open(21, iostat=ios, &
775  file=inpath//'/general/'//grip_temp_file, &
776  status='old')
777 if (ios /= 0) stop ' sico_init: Error when opening the data file for delta_ts!'
778 
779 read(21, fmt=*) ch_dummy, grip_time_min, grip_time_stp, grip_time_max
780 
781 if (ch_dummy /= '#') then
782  write(6, fmt=*) ' sico_init: grip_time_min, grip_time_stp, grip_time_max'
783  write(6, fmt=*) ' not defined in data file for delta_ts!'
784  stop
785 end if
786 
787 ndata_grip = (grip_time_max-grip_time_min)/grip_time_stp
788 
789 allocate(griptemp(0:ndata_grip))
790 
791 do n=0, ndata_grip
792  read(21, fmt=*) d_dummy, griptemp(n)
793 end do
794 
795 close(21, status='keep')
796 
797 #endif
798 
799 !-------- Read data for z_sl --------
800 
801 #if (SEA_LEVEL==3)
802 
803 open(21, iostat=ios, &
804  file=inpath//'/general/'//sea_level_file, &
805  status='old')
806 if (ios /= 0) stop ' sico_init: Error when opening the data file for z_sl!'
807 
808 read(21, fmt=*) ch_dummy, specmap_time_min, specmap_time_stp, specmap_time_max
809 
810 if (ch_dummy /= '#') then
811  write(6, fmt=*) ' sico_init: specmap_time_min, specmap_time_stp, specmap_time_max'
812  write(6, fmt=*) ' not defined in data file for z_sl!'
813  stop
814 end if
815 
816 ndata_specmap = (specmap_time_max-specmap_time_min)/specmap_time_stp
817 
818 allocate(specmap_zsl(0:ndata_specmap))
819 
820 do n=0, ndata_specmap
821  read(21, fmt=*) d_dummy, specmap_zsl(n)
822 end do
823 
824 close(21, status='keep')
825 
826 #endif
827 
828 !-------- Determination of the geothermal heat flux --------
829 
830 #if (Q_GEO_MOD==1)
831 
832 ! ------ Constant value
833 
834 do i=0, imax
835 do j=0, jmax
836  q_geo(j,i) = q_geo *1.0e-03_dp ! mW/m2 -> W/m2
837 end do
838 end do
839 
840 #elif (Q_GEO_MOD==2)
841 
842 stop ' sico_init: Option Q_GEO_MOD==2 not available for this application!'
843 
844 #endif
845 
846 !-------- Determination of the time lag
847 ! of the relaxing asthenosphere --------
848 
849 #if (REBOUND==1 || REBOUND==2)
850 
851 #if (TIME_LAG_MOD==1)
852 
853 time_lag_asth = time_lag*year_sec ! a -> s
854 
855 #elif (TIME_LAG_MOD==2)
856 
857 open(29, iostat=ios, &
858  file=inpath//'/'//trim(ch_domain_short)//'/'//time_lag_file, &
859  recl=8192, status='old')
860 
861 if (ios /= 0) stop ' sico_init: Error when opening the time-lag file!'
862 
863 do n=1, 6; read(29, fmt='(a)') ch_dummy; end do
864 
865 do j=jmax, 0, -1
866  read(29, fmt=*) (time_lag_asth(j,i), i=0,imax)
867 end do
868 
869 close(29, status='keep')
870 
871 time_lag_asth = time_lag_asth*year_sec ! a -> s
872 
873 #endif
874 
875 #elif (REBOUND==0)
876 
877 time_lag_asth = 0.0_dp ! dummy values
878 
879 #endif
880 
881 !-------- Determination of the flexural rigidity of the lithosphere --------
882 
883 #if (REBOUND==2)
884 
885 #if (FLEX_RIG_MOD==1)
886 
887 flex_rig_lith = flex_rig
888 
889 #elif (FLEX_RIG_MOD==2)
890 
891 open(29, iostat=ios, &
892  file=inpath//'/'//trim(ch_domain_short)//'/'//flex_rig_file, &
893  recl=8192, status='old')
894 
895 if (ios /= 0) stop ' sico_init: Error when opening the flex-rig file!'
896 
897 do n=1, 6; read(29, fmt='(a)') ch_dummy; end do
898 
899 do j=jmax, 0, -1
900  read(29, fmt=*) (flex_rig_lith(j,i), i=0,imax)
901 end do
902 
903 close(29, status='keep')
904 
905 #endif
906 
907 #elif (REBOUND==0 || REBOUND==1)
908 
909 flex_rig_lith = 0.0_dp ! dummy values
910 
911 #endif
912 
913 !-------- Definition of initial values --------
914 
915 ! ------ Initial topography with a thin ice layer everywhere
916 
917 #if (ANF_DAT==1)
918 
919 call topography1(dxi, deta)
920 
921 z_sl = -1.11e+11_dp ! dummy value for initial call of subroutine boundary
922 
923 call boundary(time_init, dtime, dxi, deta, &
924  delta_ts, glac_index, z_sl, dzsl_dtau, z_mar)
925 
926 q_bm = 0.0_dp
927 q_tld = 0.0_dp
928 q_b_tot = 0.0_dp
929 
930 p_b_w = 0.0_dp
931 h_w = 0.0_dp
932 
933 call init_temp_water_age()
934 
935 call calc_enhance(time_init)
936 
937 ! ------ Ice-free, relaxed bedrock
938 
939 #elif (ANF_DAT==2)
940 
941 call topography2(dxi, deta)
942 
943 z_sl = -1.11e+11_dp ! dummy value for initial call of subroutine boundary
944 
945 call boundary(time_init, dtime, dxi, deta, &
946  delta_ts, glac_index, z_sl, dzsl_dtau, z_mar)
947 
948 q_bm = 0.0_dp
949 q_tld = 0.0_dp
950 q_b_tot = 0.0_dp
951 
952 p_b_w = 0.0_dp
953 h_w = 0.0_dp
954 
955 call init_temp_water_age()
956 
957 call calc_enhance(time_init)
958 
959 ! ------ Read initial state from output of previous simulation
960 
961 #elif (ANF_DAT==3)
962 
963 #if (NETCDF==1)
964 call topography3(dxi, deta, z_sl, anfdatname)
965 #elif (NETCDF==2)
966 call topography3_nc(dxi, deta, z_sl, anfdatname)
967 #else
968 stop ' sico_init: Parameter NETCDF must be either 1 or 2!'
969 #endif
970 
971 call boundary(time_init, dtime, dxi, deta, &
972  delta_ts, glac_index, z_sl, dzsl_dtau, z_mar)
973 
974 q_b_tot = q_bm + q_tld
975 
976 #endif
977 
978 !-------- Grounding line and calving front flags --------
979 
980 flag_grounding_line_1 = .false.
981 flag_grounding_line_2 = .false.
982 
983 flag_calving_front_1 = .false.
984 flag_calving_front_2 = .false.
985 
986 #if (MARGIN==1 || MARGIN==2)
987 
988 continue
989 
990 #elif (MARGIN==3)
991 
992 do i=1, imax-1
993 do j=1, jmax-1
994 
995  if ( (maske(j,i)==0_i2b) & ! grounded ice
996  .and. &
997  ( (maske(j,i+1)==3_i2b) & ! with
998  .or.(maske(j,i-1)==3_i2b) & ! one
999  .or.(maske(j+1,i)==3_i2b) & ! neighbouring
1000  .or.(maske(j-1,i)==3_i2b) ) & ! floating ice point
1001  ) &
1002  flag_grounding_line_1(j,i) = .true.
1003 
1004  if ( (maske(j,i)==3_i2b) & ! floating ice
1005  .and. &
1006  ( (maske(j,i+1)==0_i2b) & ! with
1007  .or.(maske(j,i-1)==0_i2b) & ! one
1008  .or.(maske(j+1,i)==0_i2b) & ! neighbouring
1009  .or.(maske(j-1,i)==0_i2b) ) & ! grounded ice point
1010  ) &
1011  flag_grounding_line_2(j,i) = .true.
1012 
1013  if ( (maske(j,i)==3_i2b) & ! floating ice
1014  .and. &
1015  ( (maske(j,i+1)==2_i2b) & ! with
1016  .or.(maske(j,i-1)==2_i2b) & ! one
1017  .or.(maske(j+1,i)==2_i2b) & ! neighbouring
1018  .or.(maske(j-1,i)==2_i2b) ) & ! ocean point
1019  ) &
1020  flag_calving_front_1(j,i) = .true.
1021 
1022  if ( (maske(j,i)==2_i2b) & ! ocean
1023  .and. &
1024  ( (maske(j,i+1)==3_i2b) & ! with
1025  .or.(maske(j,i-1)==3_i2b) & ! one
1026  .or.(maske(j+1,i)==3_i2b) & ! neighbouring
1027  .or.(maske(j-1,i)==3_i2b) ) & ! floating ice point
1028  ) &
1029  flag_calving_front_2(j,i) = .true.
1030 
1031 end do
1032 end do
1033 
1034 do i=1, imax-1
1035 
1036  j=0
1037  if ( (maske(j,i)==2_i2b) & ! ocean
1038  .and. (maske(j+1,i)==3_i2b) & ! with one neighbouring
1039  ) & ! floating ice point
1040  flag_calving_front_2(j,i) = .true.
1041 
1042  j=jmax
1043  if ( (maske(j,i)==2_i2b) & ! ocean
1044  .and. (maske(j-1,i)==3_i2b) & ! with one neighbouring
1045  ) & ! floating ice point
1046  flag_calving_front_2(j,i) = .true.
1047 
1048 end do
1049 
1050 do j=1, jmax-1
1051 
1052  i=0
1053  if ( (maske(j,i)==2_i2b) & ! ocean
1054  .and. (maske(j,i+1)==3_i2b) & ! with one neighbouring
1055  ) & ! floating ice point
1056  flag_calving_front_2(j,i) = .true.
1057 
1058  i=imax
1059  if ( (maske(j,i)==2_i2b) & ! ocean
1060  .and. (maske(j,i-1)==3_i2b) & ! with one neighbouring
1061  ) & ! floating ice point
1062  flag_calving_front_2(j,i) = .true.
1063 
1064 end do
1065 
1066 #else
1067 stop ' sico_init: MARGIN must be either 1, 2 or 3!'
1068 #endif
1069 
1070 !-------- Distance between grid points with delta_i=ir, delta_j=jr --------
1071 
1072 #if (GRID==0 || GRID==1)
1073 
1074 do ir=-imax, imax
1075 do jr=-jmax, jmax
1076  dist_dxdy(jr,ir) = sqrt( (real(ir,dp)*dxi)**2 + (real(jr,dp)*deta)**2 )
1077  ! distortion due to stereographic projection not accounted for
1078 end do
1079 end do
1080 
1081 #endif
1082 
1083 !-------- Initial velocities --------
1084 
1085 call calc_temp_melt()
1086 
1087 #if (DYNAMICS==1 || DYNAMICS==2)
1088 
1089 call calc_vxy_b_sia(time, z_sl)
1090 call calc_vxy_sia(dzeta_c, dzeta_t)
1091 
1092 #if (MARGIN==3 || DYNAMICS==2)
1093 call calc_vxy_ssa(z_sl, dxi, deta, dzeta_c, dzeta_t, ii, jj, nn)
1094 #endif
1095 
1096 call calc_vz_sia(dxi, deta, dzeta_c, dzeta_t)
1097 
1098 #if (MARGIN==3)
1099 call calc_vz_ssa(z_sl, dxi, deta, dzeta_c)
1100 #endif
1101 
1102 #elif (DYNAMICS==0)
1103 
1104 call calc_vxy_static()
1105 call calc_vz_static()
1106 
1107 #else
1108  stop ' sico_init: DYNAMICS must be either 0, 1 or 2!'
1109 #endif
1110 
1111 call calc_dxyz(dxi, deta, dzeta_c, dzeta_t)
1112 
1113 !-------- Initial enthalpies --------
1114 
1115 #if (CALCMOD==0 || CALCMOD==-1)
1116 
1117 do i=0, imax
1118 do j=0, jmax
1119 
1120  do kc=0, kcmax
1121  enth_c(kc,j,i) = enth_fct_temp_omega(temp_c(kc,j,i), 0.0_dp)
1122  end do
1123 
1124  do kt=0, ktmax
1125  enth_t(kt,j,i) = enth_c(0,j,i)
1126  end do
1127 
1128 end do
1129 end do
1130 
1131 #elif (CALCMOD==1)
1132 
1133 do i=0, imax
1134 do j=0, jmax
1135 
1136  do kc=0, kcmax
1137  enth_c(kc,j,i) = enth_fct_temp_omega(temp_c(kc,j,i), 0.0_dp)
1138  end do
1139 
1140  if ( (maske(j,i) == 0_i2b).and.(n_cts(j,i) == 1_i2b) ) then
1141  do kt=0, ktmax
1142  enth_t(kt,j,i) = enth_fct_temp_omega(temp_t_m(kt,j,i), omega_t(kt,j,i))
1143  end do
1144  else
1145  do kt=0, ktmax
1146  enth_t(kt,j,i) = enth_c(0,j,i)
1147  end do
1148  end if
1149 
1150 end do
1151 end do
1152 
1153 #elif (CALCMOD==2 || CALCMOD==3)
1154 
1155 do i=0, imax
1156 do j=0, jmax
1157 
1158  do kc=0, kcmax
1159  enth_c(kc,j,i) = enth_fct_temp_omega(temp_c(kc,j,i), omega_c(kc,j,i))
1160  end do
1161 
1162  do kt=0, ktmax
1163  enth_t(kt,j,i) = enth_c(0,j,i)
1164  end do
1165 
1166 end do
1167 end do
1168 
1169 #else
1170 
1171 stop ' sico_init: Parameter CALCMOD must be either -1, 0, 1, 2 or 3!'
1172 
1173 #endif
1174 
1175 !-------- Initialize time-series files --------
1176 
1177 ! ------ Time-series file for the ice sheet on the whole
1178 
1179 open(12, iostat=ios, &
1180  file=outpath//'/'//trim(runname)//'.ser', &
1181  status='new')
1182 if (ios /= 0) stop ' sico_init: Error when opening the ser file!'
1183 
1184 write(12,1102)
1185 write(12,1103)
1186 
1187 #if (!defined(OUTSER_V_A) || OUTSER_V_A==1)
1188 
1189  1102 format(' t(a) D_Ts(deg C) z_sl(m)',/, &
1190  ' H_max(m) zs_max(m) V_g(m^3)', &
1191  ' V_t(m^3) V_fw(m^3/a) V_sle(m)',/, &
1192  ' A_g(m^2) A_t(m^2) V_bm(m^3/a)', &
1193  ' V_tld(m^3/a) H_t_max(m) vs_max(m/a)')
1194  1103 format('----------------------------------------------------', &
1195  '---------------------------------------')
1196 
1197 #elif (OUTSER_V_A==2)
1198 
1199  1102 format(' t(a) D_Ts(deg C) z_sl(m)',/, &
1200  ' V(m^3) V_g(m^3) V_f(m^3)', &
1201  ' A(m^2) A_g(m^2) A_f(m^2)',/, &
1202  ' H_max(m) zs_max(m)', &
1203  ' V_t(m^3) V_fw(m^3/a) V_sle(m)',/, &
1204  ' A_t(m^2) V_bm(m^3/a)', &
1205  ' V_tld(m^3/a) H_t_max(m) vs_max(m/a)')
1206  1103 format('----------------------------------------------------', &
1207  '---------------------------------------')
1208 
1209 #else
1210  stop ' sico_init: OUTSER_V_A must be either 1 or 2!'
1211 #endif
1212 
1213 ! ------ Time-series file for the sediment area ("Hudson Bay, Hudson Strait")
1214 
1215 open(15, iostat=ios, &
1216  file=outpath//'/'//trim(runname)//'.sed', &
1217  status='new')
1218 if (ios /= 0) stop ' sico_init: Error when opening the sed file!'
1219 
1220 write(15,1108)
1221 write(15,1109)
1222 
1223 1108 format(' t(a) D_Ts(deg C) z_sl(m)',/, &
1224  ' H_ave(m) Tbh_ave(C) Atb(m^2)')
1225 1109 format('----------------------------------------------------')
1226 
1227 ! ------ Time-series file for the ice-thickness field
1228 
1229 #if (OUTSER==2)
1230 
1231 open(13, iostat=ios, &
1232  file=outpath//'/'//trim(runname)//'.thk', &
1233  status='new')
1234 if (ios /= 0) stop ' sico_init: Error when opening the thk file!'
1235 
1236 write(13,1104)
1237 write(13,1105)
1238 
1239 1104 format(' t(a) D_Ts(deg C) z_sl(m)',/, &
1240  ' H(m): JMAX+1 records [j = JMAX (-1) 0] with',/, &
1241  ' IMAX+1 values [i = 0 (1) IMAX] in each record')
1242 1105 format('----------------------------------------------------')
1243 
1244 #endif
1245 
1246 ! ------ Time-series file for selected positions ("deep boreholes")
1247 
1248 #if (OUTSER==3)
1249 
1250 n_core = 7 ! Points P1 - P7
1251 
1252 allocate(lambda_core(n_core), phi_core(n_core), &
1253  x_core(n_core), y_core(n_core))
1254 
1255 lambda_core(1) = 0.0_dp ! dummy
1256 phi_core(1) = 0.0_dp ! dummy
1257 x_core(1) = 3900.0_dp *1.0e+03_dp ! Point P1,
1258 y_core(1) = 2000.0_dp *1.0e+03_dp ! conversion km -> m
1259 
1260 lambda_core(2) = 0.0_dp ! dummy
1261 phi_core(2) = 0.0_dp ! dummy
1262 x_core(2) = 3800.0_dp *1.0e+03_dp ! Point P2,
1263 y_core(2) = 2000.0_dp *1.0e+03_dp ! conversion km -> m
1264 
1265 lambda_core(3) = 0.0_dp ! dummy
1266 phi_core(3) = 0.0_dp ! dummy
1267 x_core(3) = 3700.0_dp *1.0e+03_dp ! Point P3,
1268 y_core(3) = 2000.0_dp *1.0e+03_dp ! conversion km -> m
1269 
1270 lambda_core(4) = 0.0_dp ! dummy
1271 phi_core(4) = 0.0_dp ! dummy
1272 x_core(4) = 3500.0_dp *1.0e+03_dp ! Point P4,
1273 y_core(4) = 2000.0_dp *1.0e+03_dp ! conversion km -> m
1274 
1275 lambda_core(5) = 0.0_dp ! dummy
1276 phi_core(5) = 0.0_dp ! dummy
1277 x_core(5) = 3200.0_dp *1.0e+03_dp ! Point P5,
1278 y_core(5) = 2000.0_dp *1.0e+03_dp ! conversion km -> m
1279 
1280 lambda_core(6) = 0.0_dp ! dummy
1281 phi_core(6) = 0.0_dp ! dummy
1282 x_core(6) = 2900.0_dp *1.0e+03_dp ! Point P6,
1283 y_core(6) = 2000.0_dp *1.0e+03_dp ! conversion km -> m
1284 
1285 lambda_core(7) = 0.0_dp ! dummy
1286 phi_core(7) = 0.0_dp ! dummy
1287 x_core(7) = 2600.0_dp *1.0e+03_dp ! Point P7,
1288 y_core(7) = 2000.0_dp *1.0e+03_dp ! conversion km -> m
1289 
1290 open(14, iostat=ios, &
1291  file=outpath//'/'//trim(runname)//'.core', &
1292  status='new')
1293 if (ios /= 0) stop ' sico_init: Error when opening the core file!'
1294 
1295 if (forcing_flag == 1) then
1296 
1297  write(14,1106)
1298  write(14,1107)
1299 
1300  1106 format(' t(a) D_Ts(C) z_sl(m)',/, &
1301  ' H_P1(m) H_P2(m) H_P3(m) H_P4(m)', &
1302  ' H_P5(m) H_P6(m) H_P7(m)',/, &
1303  ' v_P1(m/a) v_P2(m/a) v_P3(m/a) v_P4(m/a)', &
1304  ' v_P5(m/a) v_P6(m/a) v_P7(m/a)',/, &
1305  ' T_P1(C) T_P2(C) T_P3(C) T_P4(C)', &
1306  ' T_P5(C) T_P6(C) T_P7(C)',/, &
1307  ' Rb_P1(W/m2) Rb_P2(W/m2) Rb_P3(W/m2) Rb_P4(W/m2)', &
1308  ' Rb_P5(W/m2) Rb_P6(W/m2) Rb_P7(W/m2)')
1309  1107 format('-----------------------------------------------------------------', &
1310  '---------------------------------------')
1311 
1312 end if
1313 
1314 #endif
1315 
1316 !-------- Output of the initial state --------
1317 
1318 #if (OUTPUT==1)
1319 
1320 #if (ERGDAT==0)
1321  flag_3d_output = .false.
1322 #elif (ERGDAT==1)
1323  flag_3d_output = .true.
1324 #else
1325  stop ' sico_init: ERGDAT must be either 0 or 1!'
1326 #endif
1327 
1328 #if (NETCDF==1)
1329  call output1(runname, time_init, delta_ts, glac_index, z_sl, &
1330  flag_3d_output, ndat2d, ndat3d)
1331 #elif (NETCDF==2)
1332  call output_nc(runname, time_init, delta_ts, glac_index, z_sl, &
1333  flag_3d_output, ndat2d, ndat3d)
1334 #else
1335  stop ' sico_init: Parameter NETCDF must be either 1 or 2!'
1336 #endif
1337 
1338 #elif (OUTPUT==2)
1339 
1340 if (time_output(1) <= time_init+eps) then
1341 
1342 #if (ERGDAT==0)
1343  flag_3d_output = .false.
1344 #elif (ERGDAT==1)
1345  flag_3d_output = .true.
1346 #else
1347  stop ' sico_init: ERGDAT must be either 0 or 1!'
1348 #endif
1349 
1350 #if (NETCDF==1)
1351  call output1(runname, time_init, delta_ts, glac_index, z_sl, &
1352  flag_3d_output, ndat2d, ndat3d)
1353 #elif (NETCDF==2)
1354  call output_nc(runname, time_init, delta_ts, glac_index, z_sl, &
1355  flag_3d_output, ndat2d, ndat3d)
1356 #else
1357  stop ' sico_init: Parameter NETCDF must be either 1 or 2!'
1358 #endif
1359 
1360 end if
1361 
1362 #elif (OUTPUT==3)
1363 
1364  flag_3d_output = .false.
1365 
1366 #if (NETCDF==1)
1367  call output1(runname, time_init, delta_ts, glac_index, z_sl, &
1368  flag_3d_output, ndat2d, ndat3d)
1369 #elif (NETCDF==2)
1370  call output_nc(runname, time_init, delta_ts, glac_index, z_sl, &
1371  flag_3d_output, ndat2d, ndat3d)
1372 #else
1373  stop ' sico_init: Parameter NETCDF must be either 1 or 2!'
1374 #endif
1375 
1376 if (time_output(1) <= time_init+eps) then
1377 
1378  flag_3d_output = .true.
1379 
1380 #if (NETCDF==1)
1381  call output1(runname, time_init, delta_ts, glac_index, z_sl, &
1382  flag_3d_output, ndat2d, ndat3d)
1383 #elif (NETCDF==2)
1384  call output_nc(runname, time_init, delta_ts, glac_index, z_sl, &
1385  flag_3d_output, ndat2d, ndat3d)
1386 #else
1387  stop ' sico_init: Parameter NETCDF must be either 1 or 2!'
1388 #endif
1389 
1390 end if
1391 
1392 #else
1393  stop ' sico_init: OUTPUT must be either 1, 2 or 3!'
1394 #endif
1395 
1396 call output2(time_init, dxi, deta, delta_ts, glac_index, z_sl)
1397 #if (OUTSER==2)
1398 call output3(time_init, delta_ts, glac_index, z_sl)
1399 #endif
1400 #if (OUTSER==3)
1401 call output4(time_init, dxi, deta, delta_ts, glac_index, z_sl)
1402 #endif
1403 
1404 end subroutine sico_init
1405 !
subroutine calc_vz_ssa(z_sl, dxi, deta, dzeta_c)
Computation of the vertical velocity vz in the shallow shelf approximation.
Definition: calc_vz_ssa.F90:35
subroutine calc_vxy_static()
Computation of the horizontal velocity vx, vy, the horizontal volume flux qx, qy etc. for static ice.
subroutine calc_vxy_sia(dzeta_c, dzeta_t)
Computation of the shear stresses txz, tyz, the effective shear stress sigma, the depth-averaged flui...
subroutine topography1(dxi, deta)
Definition of the initial surface and bedrock topography (including gradients) and of the horizontal ...
Definition: topography1.F90:39
subroutine, public ice_mat_eqs_pars(RF_table, R_T_val, KAPPA_table, C_table, n_temp_min, n_temp_max, RHO_I_val, RHO_C_val, KAPPA_C_val, C_C_val)
Setting of required physical parameters.
subroutine output2(time, dxi, deta, delta_ts, glac_index, z_sl)
Writing of time-series data on file in ASCII format.
Definition: output2.F90:35
subroutine phys_para()
Reading of physical parameters.
Definition: phys_para.F90:35
Declarations of kind types for SICOPOLIS.
Definition: sico_types.F90:35
subroutine topography3(dxi, deta, z_sl, anfdatname)
Definition of the initial surface and bedrock topography (including gradients) and of the horizontal ...
Definition: topography3.F90:39
subroutine output3(time, delta_ts, glac_index, z_sl)
Writing of time-series data of the ice thickness field on file in ASCII format.
Definition: output3.F90:38
subroutine topography3_nc(dxi, deta, z_sl, anfdatname)
Definition of the initial surface and bedrock topography (including gradients) and of the horizontal ...
subroutine calc_temp_melt()
Computation of the melting temperatures.
subroutine output_nc(runname, time, delta_ts, glac_index, z_sl, flag_3d_output, ndat2d, ndat3d)
Writing of time-slice files in NetCDF format.
Definition: output_nc.F90:35
subroutine calc_vxy_ssa(z_sl, dxi, deta, dzeta_c, dzeta_t, ii, jj, nn)
Computation of the horizontal velocity vx, vy, the horizontal volume flux qx, qy and the flux across ...
Declarations of global variables for SICOPOLIS (for the ANT domain).
Definition: sico_vars.F90:35
subroutine topography2(dxi, deta)
Definition of the initial surface and bedrock topography (including gradients) and of the horizontal ...
Definition: topography2.F90:39
subroutine calc_vz_static()
Computation of the vertical velocity vz for static ice.
subroutine calc_dxyz(dxi, deta, dzeta_c, dzeta_t)
Computation of all components of the strain-rate tensor, the full effective strain rate and the shear...
Definition: calc_dxyz.F90:37
subroutine calc_enhance(time)
Computation of the flow enhancement factor.
subroutine, public calc_c_int_table(c_table, n_tmp_min, n_tmp_max, L_val)
Computation of the temperature integral of the specific heat of ice as a table (c_int_table). Further, definition of the latent heat of ice.
subroutine calc_vxy_b_sia(time, z_sl)
Computation of the basal horizontal velocity vx_b, vy_b in the shallow ice approximation.
subroutine sico_init(delta_ts, glac_index, mean_accum, dtime, dtime_temp, dtime_wss, dtime_out, dtime_ser, time, time_init, time_end, time_output, dxi, deta, dzeta_c, dzeta_t, dzeta_r, z_sl, dzsl_dtau, z_mar, ii, jj, nn, ndat2d, ndat3d, n_output, runname)
Initialisations for SICOPOLIS.
Definition: sico_init.F90:35
subroutine boundary(time, dtime, dxi, deta, delta_ts, glac_index, z_sl, dzsl_dtau, z_mar)
Computation of the surface temperature (must be less than 0 deg C!) and of the accumulation-ablation ...
Definition: boundary.F90:37
subroutine init_temp_water_age()
Initial temperature, water content and age.
Material quantities of ice: Rate factor, heat conductivity, specific heat (heat capacity).
subroutine output4(time, dxi, deta, delta_ts, glac_index, z_sl)
Writing of time-series data of the deep ice cores on file in ASCII format.
Definition: output4.F90:35
Declarations of global variables for SICOPOLIS.
subroutine, public calc_c_int_inv_table()
Computation of the inverse of the temperature integral of the specific heat of ice as a table (c_int_...
Conversion from temperature (temp) and water content (omega) to enthalpy (enth) and vice versa...
subroutine calc_vz_sia(dxi, deta, dzeta_c, dzeta_t)
Computation of the vertical velocity vz in the shallow ice approximation.
Definition: calc_vz_sia.F90:35
subroutine output1(runname, time, delta_ts, glac_index, z_sl, flag_3d_output, ndat2d, ndat3d)
Writing of time-slice files in binary format.
Definition: output1.F90:35
real(dp) function, public enth_fct_temp_omega(temp_val, omega_val)
Enthalpy as a function of temperature and water content.