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