The TELEMAC-MASCARET system  trunk
lecdon_sisyphe.f
Go to the documentation of this file.
1 ! *************************
2  SUBROUTINE lecdon_sisyphe
3 ! *************************
4 !
5  &(motcar,file_desc,path,ncar,code,cas_file,dico_file)
6 !
7 !***********************************************************************
8 ! SISYPHE V7P2
9 !***********************************************************************
10 !
11 !brief READS THE STEERING FILE BY CALL TO DAMOCLES.
12 !
13 !history C.VILLARET + JMH (EDF-LNHE)
14 !+ 02/05/2012
15 !+ V6P2
16 !+ File for liquid boundaries added
17 !+
18 !
19 !history JWI
20 !+ 31/05/2012
21 !+ V6P2
22 !+ added one increment to include wave orbital velocities
23 !+ (SORLEO(I+28+(NOMBLAY+4)*NSICLA+NOMBLAY).OR.
24 !
25 !history PAT (LNHE)
26 !+ 18/06/2012
27 !+ V6P2
28 !+ updated version with HRW's development
29 !
30 !history Pablo Tassi PAT (EDF-LNHE)
31 !+ 12/02/2013
32 !+ V6P3
33 !+ Preparing for the use of a higher NSICLM value
34 !+ (by Rebekka Kopmann)
35 !
36 !history Pablo Tassi PAT (EDF-LNHE)
37 !+ 12/02/2013
38 !+ V6P3
39 !+ Settling lag: determines choice between Rouse and Miles concentration profile
40 !+ (by Michiel Knaapen HRW)
41 !
42 !history J-M HERVOUET (EDF LAB, LNHE)
43 !+ 18/05/2015
44 !+ V7P1
45 !+ Adding CHECK_MESH for the keyword 'CHECKING THE MESH'
46 !
47 !history R. KOPMANN (BAW)
48 !+ 13/07/2016
49 !+ V7P2
50 !+ Integrating liquid boundary file for QS
51 !
52 !history S.E. BOURBAN (HRW)
53 !+ 20/06/2016
54 !+ V7P2
55 !+ Now taking into account names of differentiators given by user.
56 !
57 !history R.KOPMANN (BAW)
58 !+ 07/12/2017
59 !+ V7P2
60 !+ VSPRES results file for CVSM handled by CAS-file.
61 !
62 !history B.GLANDER (BAW)
63 !+ 06/12/2018
64 !+ V7P2
65 !+ NEW VARIABLE: ZRL REFERENCE LEVEL FOR NESTOR, CHANGE OF SORLEO, SORIMP
66 !
67 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
68 !| FILE_DESC |<--| STORES STRINGS 'SUBMIT' OF DICTIONARY
69 !| MOTCAR |<--| VALUES OF KEY-WORDS OF TYPE CHARACTER
70 !| NCAR |-->| NUMBER OF LETTERS IN STRING PATH
71 !| PATH |-->| FULL PATH TO CODE DICTIONARY
72 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
73 !
77 !
79  IMPLICIT NONE
80 !
81 !+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+
82 !
83  INTEGER, INTENT(IN) :: NCAR
84  CHARACTER(LEN=24), INTENT(IN) :: CODE
85  CHARACTER(LEN=PATH_LEN), INTENT(IN) :: PATH
86  CHARACTER(LEN=PATH_LEN), INTENT(INOUT) :: MOTCAR(maxkeyword)
87  CHARACTER(LEN=PATH_LEN), INTENT(INOUT) :: FILE_DESC(4,maxkeyword)
88 ! API
89  CHARACTER(LEN=PATH_LEN), INTENT(IN) :: CAS_FILE
90  CHARACTER(LEN=PATH_LEN), INTENT(IN) :: DICO_FILE
91 !
92 !+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+
93 !
94  INTEGER :: I,K,ERR
95  INTEGER :: MOTINT(maxkeyword)
96  INTEGER :: TROUVE(4,maxkeyword)
97  INTEGER :: ADRESS(4,maxkeyword)
98  INTEGER :: DIMENS(4,maxkeyword)
99  DOUBLE PRECISION :: SUMAVAI
100  DOUBLE PRECISION :: MOTREA(maxkeyword)
101  LOGICAL :: DOC,EFFPEN
102  LOGICAL :: MOTLOG(maxkeyword)
103  CHARACTER(LEN=PATH_LEN) :: NOM_CAS
104  CHARACTER(LEN=PATH_LEN) :: NOM_DIC
105  CHARACTER(LEN=72) :: MOTCLE(4,maxkeyword,2)
106 
107  CHARACTER(LEN=PATH_LEN) TEMPVAR
108  INTEGER :: ID_DICO, ID_CAS
109 !
110 !-----------------------------------------------------------------------
111 !
112  CHARACTER(LEN=2) CHAR2
113 !
114 !-----------------------------------------------------------------------
115 !
116 ! INITIALISES THE VARIABLES FOR DAMOCLES CALL :
117 !
118  DO k = 1, maxkeyword
119 ! A FILENAME NOT GIVEN BY DAMOCLES WILL BE RECOGNIZED AS A WHITE SPACE
120 ! (IT MAY BE THAT NOT ALL COMPILERS WILL INITIALISE LIKE THAT)
121  motcar(k)(1:1)=' '
122 !
123  dimens(1,k) = 0
124  dimens(2,k) = 0
125  dimens(3,k) = 0
126  dimens(4,k) = 0
127  ENDDO
128 !
129 ! WRITES OUT INFO
130  doc = .false.
131 !
132 !-----------------------------------------------------------------------
133 ! OPENS DICTIONNARY AND STEERING FILES
134 !-----------------------------------------------------------------------
135 !
136  IF(ncar.GT.0) THEN
137 !
138  nom_dic=path(1:ncar)//'SISDICO'
139  nom_cas=path(1:ncar)//'SISCAS'
140 !
141  ELSE
142 !
143  nom_dic='SISDICO'
144  nom_cas='SISCAS'
145 !
146  ENDIF
147  IF((cas_file(1:1).NE.' ').AND.(dico_file(1:1).NE.' ')) THEN
148  WRITE(lu,*) 'FIXED DICO AND STEERING FILE PRESENT'
149  nom_dic=dico_file
150  nom_cas=cas_file
151  WRITE(lu,*) 'NOM_DIC',nom_dic
152  WRITE(lu,*) 'NOM_CAS',nom_cas
153  ENDIF
154 !
155  CALL get_free_id(id_dico)
156  OPEN(id_dico,file=nom_dic,form='FORMATTED',action='READ')
157  CALL get_free_id(id_cas)
158  OPEN(id_cas,file=nom_cas,form='FORMATTED',action='READ')
159 !
160 !-----------------------------------------------------------------------
161 ! CALLS DAMOCLES
162 !-----------------------------------------------------------------------
163 !
164  CALL damocle( adress, dimens ,maxkeyword, doc , lng , lu ,
165  & motint, motrea ,motlog , motcar ,
166  & motcle, trouve ,id_dico, id_cas,.false. ,file_desc)
167 !
168 !-----------------------------------------------------------------------
169 ! CLOSES DICTIONNARY AND STEERING FILES
170 !-----------------------------------------------------------------------
171 !
172  CLOSE(id_dico)
173  CLOSE(id_cas)
174 !
175 ! DECODES 'SUBMIT' CHAINS
176 !
177  CALL read_submit(sis_files,maxlu_sis,file_desc,maxkeyword)
178 !
179 !-----------------------------------------------------------------------
180 !
181 ! RETRIEVES FILES NUMBERS IN TELEMAC-2D FORTRAN PARAMETERS
182 ! AT THIS LEVEL LOGICAL UNITS ARE EQUAL TO THE FILE NUMBER
183 !
184  DO i=1,maxlu_sis
185  IF(sis_files(i)%TELNAME.EQ.'SISHYD') THEN
186  sishyd=i
187  ELSEIF(sis_files(i)%TELNAME.EQ.'SISGEO') THEN
188  sisgeo=i
189  ELSEIF(sis_files(i)%TELNAME.EQ.'SISCLI') THEN
190  siscli=i
191  ELSEIF(sis_files(i)%TELNAME.EQ.'SISPRE') THEN
192  sispre=i
193  ELSEIF(sis_files(i)%TELNAME.EQ.'SISRES') THEN
194  sisres=i
195  ELSEIF(sis_files(i)%TELNAME.EQ.'SISREF') THEN
196  sisref=i
197  ELSEIF(sis_files(i)%TELNAME.EQ.'SISCOU') THEN
198  siscou=i
199  ELSEIF(sis_files(i)%TELNAME.EQ.'SISFON') THEN
200  sisfon=i
201  ELSEIF(sis_files(i)%TELNAME.EQ.'SISSEC') THEN
202  sissec=i
203  ELSEIF(sis_files(i)%TELNAME.EQ.'SISSEO') THEN
204  sisseo=i
205  ELSEIF(sis_files(i)%TELNAME.EQ.'SISLIQ') THEN
206  sisliq=i
207  ELSEIF(sis_files(i)%TELNAME.EQ.'SISFLX') THEN
208  sisflx=i
209 ! === NESTOR FILES
210  ELSEIF(sis_files(i)%TELNAME.EQ.'SINACT') THEN
211  sinact=i
212  ELSEIF(sis_files(i)%TELNAME.EQ.'SINPOL') THEN
213  sinpol=i
214  ELSEIF(sis_files(i)%TELNAME.EQ.'SINREF') THEN
215  sinref=i
216  ELSEIF(sis_files(i)%TELNAME.EQ.'SINRST') THEN
217  sinrst=i
218 ! ===
219  ELSEIF(sis_files(i)%TELNAME.EQ.'VSPRES') THEN
220  vspres=i
221 !
222  ENDIF
223  ENDDO
224 !
225 !-----------------------------------------------------------------------
226 !
227 ! ALLOCATING ARRAYS DEPENDING ON MAXFRO
228 !
229 !-----------------------------------------------------------------------
230 !
231 ! MAXIMUM NUMBER OF BOUNDARIES
232  maxfro = motint( adress(1,58) )
233 !
234  ALLOCATE(soldis(maxfro))
235  ALLOCATE(okcgl(maxfro))
236  ALLOCATE(okqgl(maxfro))
237  ALLOCATE(cbor_classe(nsiclm*maxfro))
238 !
239  DO k=1,maxfro
240  okcgl(k)=.true.
241  okqgl(k)=.true.
242  ENDDO
243 !
244 !-----------------------------------------------------------------------
245 !
246 ! ASSIGNS THE STEERING FILE VALUES TO THE PARAMETER FORTRAN NAME
247 !
248 !-----------------------------------------------------------------------
249 !
250 ! OPTION OF MATRIX ASSEMBLY IS HARD-CODED
251 !
252  optass = 1
253  produc = 1
254 !
255 ! DISCRETISATIONS
256 !
257  ielmt = 11 ! SEDIMENTOLOGICAL VARIABLES
258  ielmh_sis = 11 ! VARIABLES ASSOCIATED WITH WATER DEPTH
259  ielmu_sis = 11 ! VARIABLES ASSOCIATED WITH VELOCITIES
260 !
261 ! FOR NOW PRINTOUTS START AT ZERO
262 !
263  ptinig = 0
264  ptinil = 0
265 !
266 ! NON-EQUILIBIRUM BEDLOAD
267 !
268  loadmeth = 0
269 !
270 ! ICM = MOTINT( ADRESS(1, 1) )
271  icf = motint( adress(1, 2) )
272  npas = motint( adress(1, 3) )
273  nmaree = motint( adress(1, 4) )
274 ! N1 = MOTINT( ADRESS(1, 5) )
275 
276  leopr = motint( adress(1, 6) )
277  lispr = motint( adress(1, 7) )
278  optban = motint( adress(1, 11) )
279  lvmac = motint( adress(1, 12) )
280  nsous = motint( adress(1, 14) )
281 !
282  mardat(1) = motint( adress(1, 15) )
283  mardat(2) = motint( adress(1, 15) + 1 )
284  mardat(3) = motint( adress(1, 15) + 2 )
285  martim(1) = motint( adress(1, 16) )
286  martim(2) = motint( adress(1, 16) + 1 )
287  martim(3) = motint( adress(1, 16) + 2 )
288 !
289  slvsed%SLV = motint( adress(1, 17) )
290  slvsed%KRYLOV = motint( adress(1, 18) )
291  slvsed%PRECON = motint( adress(1, 19) )
292  slvsed%NITMAX = motint( adress(1, 20) )
293  choix = motint( adress(1, 21) )
294  dirflu = motint( adress(1, 22) )
295  npriv = motint( adress(1, 23) )
296  nadvar = motint( adress(1,30) )
297 ! NUMBER OF DIRECTIONS FOR DIFFERENTIATED VARIABLES
298  ad_numofdir = motint( adress(1,59) )
299 !
300 ! NCSIZE = MOTINT( ADRESS(1, 24) )
301 ! NUMBER OF PROCESSORS (ALREADY GIVEN IN INIT_FILES2;
302 ! MUST BE THE SAME, BUT WHEN USING COUPLED MODELS IT CAN
303 ! WRONGLY BE DIFFERENT)
304  IF(ncsize.NE.motint(adress(1,24))) THEN
305  WRITE(lu,*) 'DIFFERENT NUMBER OF PARALLEL PROCESSORS:'
306  WRITE(lu,*) 'DECLARED BEFORE (CASE OF COUPLING ?):',ncsize
307  WRITE(lu,*) 'SISYPHE :',motint(adress(1,24))
308  WRITE(lu,*) 'VALUE ',ncsize,' IS KEPT'
309  ENDIF
310  resol = motint( adress(1, 25) )
311  slvtra%SLV = motint( adress(1, 26) )
312  slvtra%KRYLOV = motint( adress(1, 27) )
313  slvtra%PRECON = motint( adress(1, 28) )
314  slvtra%NITMAX = motint( adress(1, 29) )
315  optdif = motint( adress(1, 31) )
316  optsup = motint( adress(1, 32) )
317  produc = motint( adress(1, 33) )
318  optass = motint( adress(1, 34) )
319  opdtra = motint( adress(1, 35) )
320  kfrot = motint( adress(1, 37) )
321  ncondis = motint( adress(1, 38) )
322  slopeff = motint( adress(1, 39) )
323  devia = motint( adress(1, 40) )
324  nomblay = motint( adress(1,251) )
325  nsicla = motint( adress(1,252) )
326  hidfac = motint( adress(1,253) )
327  icq = motint( adress(1, 41) )
328 ! CONTROL SECTIONS
329  ncp=dimens(1,42)
330  ALLOCATE(ctrlsc(ncp),stat=err)
331  CALL check_allocate(err, 'CTRLSC')
332  IF(ncp.GE.1) THEN
333  DO k=1,ncp
334  ctrlsc(k) = motint( adress(1,42) + k-1 )
335  ENDDO
336  ENDIF
337 ! COORDINATES OF THE ORIGIN
338  i_orig = motint( adress(1,43) )
339  j_orig = motint( adress(1,43)+1 )
340  debug = motint( adress(1,44) )
341 !
342  icr = motint(adress(1,46) )
343 !
344  iks = motint(adress(1,47) )
345 ! CVSM MODEL
346  pro_max_max = max(nsicla*4+4,motint(adress(1,49) ))
347  cvsmpperiod = motint(adress(1,50) )
348  alt_model = motint(adress(1,52) )
349  vsmtype = motint(adress(1,53) )
350 !
351 ! MAXIMUM NUMBER OF ITERATIONS FOR ADVECTION SCHEMES
352 !
353  maxadv = motint(adress(1,54))
354 !
355 ! SCHEME OPTION FOR ADVECTION
356 !
357  optadv=1
358  IF(resol.EQ.1) THEN
359 ! CHARACTERISTICS
360 ! OPTADV=OPTCHA (IN TELEMAC-2D, NOT IN SISYPHE YET)
361  ELSEIF(resol.EQ.2) THEN
362 ! SUPG
363  optadv=optsup
364  ELSEIF(resol.EQ.5) THEN
365 ! PSI SCHEME
366 ! OPTADV_VI=OPTPSI (IN TELEMAC-2D, NOT IN SISYPHE YET)
367  ENDIF
368 ! SCHEME OPTION FOR ADVECTION HAS PRIORITY WHEN PRESENT
369  IF(trouve(1,55).EQ.2) THEN
370  optadv = motint(adress(1,55))
371  ENDIF
372  nco_dist = motint( adress(1,56) )
373  nsp_dist = motint( adress(1,57) )
374 !
375 ! ############### !
376 ! REAL KEYWORDS !
377 ! ############### !
378 !
379  ! NON-EQUILIBRIRUM BEDLOAD
380  ! ------------------------
381  ls0 = 1.d0
382 !
383  rc = motrea( adress(2, 1) )
384  xmve = motrea( adress(2, 2) )
385  xmvs = motrea( adress(2, 3) )
386  DO k=1,nsicla
387  fdm(k) = motrea( adress(2, 4) + k-1 )
388  ENDDO
389 ! IF OLD NAME OF KEYWORD 28 HAS BEEN FOUND
390  IF(trouve(2,28).EQ.2) THEN
391  DO k=1,nsicla
392  fdm(k) = motrea( adress(2,28) + k-1 )
393  ENDDO
394  ENDIF
395  xkv = motrea( adress(2, 5) )
396 ! SHIELDS NUMBERS
397  DO k=1,dimens(2,6)
398  ac(k) = motrea( adress(2, 6) + k-1 )
399  ENDDO
400 ! IF ALL SHIELDS NUMBERS ARE NOT GIVEN, THE LATEST
401 ! ONE IS TAKEN FOR THE FOLLOWING
402 ! FOR EXAMPLE IF ONLY ONE IS GIVEN, ALL WILL HAVE
403 ! THE SAME VALUE
404  IF(dimens(2,6).LT.nsicla) THEN
405  DO k=dimens(2,6)+1,nsicla
406  ac(k) = motrea( adress(2, 6)+max(dimens(2,6),1)-1 )
407  ENDDO
408  ENDIF
409  sfon = motrea( adress(2, 7) )
410  grav = motrea( adress(2, 8) )
411  zero = motrea( adress(2, 9) )
412  slvsed%ZERO = zero
413  vce = motrea( adress(2, 10) )
414  hmin = motrea( adress(2, 11) )
415  delt = motrea( adress(2, 12) )
416  tprec = motrea( adress(2, 13) )
417  pmaree = motrea( adress(2, 14) )
418  teta = motrea( adress(2, 15) )
419  beta = motrea( adress(2, 16) )
420  slvsed%EPS = motrea( adress(2, 17) )
421  teta_susp = motrea( adress(2, 18) )
422  xkx = motrea( adress(2, 19) )
423  xky = motrea( adress(2, 20) )
424  slvtra%EPS = motrea( adress(2, 21) )
425 ! SETTLING VELOCITIES (SAME TREATMENT THAN SHIELDS NUMBERS)
426  DO k=1,dimens(2,22)
427  xwc(k) = motrea( adress(2, 22) + k-1 )
428  ENDDO
429  IF(dimens(2,22).LT.nsicla) THEN
430  DO k=dimens(2,22)+1,nsicla
431  xwc(k) = motrea( adress(2, 22)+dimens(2,22)-1 )
432  ENDDO
433  ENDIF
434  crit_cfd = motrea( adress(2, 23) )
435  kspratio = motrea( adress(2, 24) )
436  phised = motrea( adress(2, 25) )
437  beta2 = motrea( adress(2, 26) )
438  bijk = motrea( adress(2, 27) )
439 !
440 ! INITIAL CONCENTRATIONS
441 ! ++++++++++++++++++++++
442 !
443  DO k=1,nsiclm
444 ! DEFAULT VALUE
445  cs0(k) = 0.d0
446  ENDDO
447  DO k=1,nsicla
448  cs0(k)=motrea( adress(2,30) + k-1 )
449  ENDDO
450  DO k=1,10*maxfro
451  cbor_classe(k)=0.d0
452  ENDDO
453  IF(dimens(2,31).GT.0) THEN
454  DO k=1,dimens(2,31)
455  cbor_classe(k)=motrea( adress(2,31) + k-1 )
456  ENDDO
457  ENDIF
458 !
459 ! COHESIVE SEDIMENT
460 ! +++++++++++++++++
461 !
462  ncouch_tass = motint( adress(1,45) )
463 !
464 ! DEFAULT VALUES, PREVIOUSLY IN THE DICTIONARY
465 ! EVEN IF THERE ARE NOT 10 LAYERS
466  conc_vase(1) = 50.d0
467  conc_vase(2) = 100.d0
468  conc_vase(3) = 150.d0
469  conc_vase(4) = 200.d0
470  conc_vase(5) = 250.d0
471  conc_vase(6) = 300.d0
472  conc_vase(7) = 350.d0
473  conc_vase(8) = 400.d0
474  conc_vase(9) = 450.d0
475  conc_vase(10) = 500.d0
476 !
477  IF(dimens(2,32).GT.0) THEN
478  DO k=1,dimens(2,32)
479  conc_vase(k)=motrea( adress(2,32) + k-1 )
480  ENDDO
481  ENDIF
482 !
483 ! OBSOLETE KEY WORD : CSF_VASE = MOTREA( ADRESS(2, 29) )
484 !
485 ! CSF_VASE = CONC_VASE(1)/XMVS
486 !
487  IF(dimens(2,34).GT.0) THEN
488  DO k=1,dimens(2,34)
489  toce_vase(k)=motrea( adress(2,34) + k-1 )
490  ENDDO
491  ENDIF
492 !
493 ! KRONE AND PARTHENIADES EROSION AND DEPOSITION LAW
494 !
495 ! OBSOLETE KEY WORD : VITCE= MOTREA( ADRESS(2,35))
496 !
497  vitce = sqrt(toce_vase(1)/xmve)
498  vitcd= motrea( adress(2,36))
499 ! PARTHENIADES WITH CONVERSION TO M/S
500  partheniades = motrea( adress(2,37))/xmvs
501 !
502 ! CONSOLIDATION MODEL
503 !
504  tass = motlog(adress(3,23))
505  itass = motint(adress(1,48) )
506 !
507 ! MULTILAYER MODEL (WALTHER, 2008)
508 ! ITASS = 1
509 !
510 ! DEFAULT VALUES, PREVIOUSLY IN THE DICTIONARY
511 ! EVEN IF THERE ARE NOT 10LAYERS
512  trans_mass(1) = 5.d-5
513  trans_mass(2) = 4.5d-5
514  trans_mass(3) = 4.d-5
515  trans_mass(4) = 3.5d-5
516  trans_mass(5) = 3.d-5
517  trans_mass(6) = 2.5d-5
518  trans_mass(7) = 2.d-5
519  trans_mass(8) = 1.5d-5
520  trans_mass(9) = 1.d-5
521  trans_mass(10) = 0.d0
522 !
523  IF(dimens(2,33).GT.0) THEN
524  DO k=1,dimens(2,33)
525  trans_mass(k)=motrea( adress(2,33) + k-1 )
526  ENDDO
527  ENDIF
528 !
529 ! V6P1
530 ! THIEBOT MULTI LAYER MODEL
531 ! ITASS=2
532 !
533  conc_gel=motrea( adress(2,38))
534  coef_n= motrea( adress(2,39))
535  conc_max=motrea( adress(2,50))
536 !
537 ! PRESCRIBED SOLID DISCHARGES
538 !
539  nsoldis=dimens(2,51)
540  IF(nsoldis.GT.0) THEN
541  DO k=1,nsoldis
542  soldis(k)=motrea(adress(2,51)+k-1)
543  ENDDO
544  ENDIF
545 !
546 ! MINIMUM DEPTH FOR BEDLOAD
547 !
548  hmin_bedload=motrea(adress(2,52))
549 !
550 ! HIDING EXPOSURE MULTI GRAIN MODEL
551 !
552  DO k=1,nsicla
553  hidi(k) = motrea( adress(2,253) + k-1 )
554  IF (trouve(2,255).EQ.1) THEN
555  fd90(k)= fdm(k)
556  ELSE
557  fd90(k)= motrea( adress(2,255) + k-1 )
558  ENDIF
559  ava0(k) = motrea( adress(2,258) + k-1 )
560  ENDDO
561  elay0 = motrea( adress(2,259) )
562 !
563 ! UM: MPM-Factor
564  mpm = motrea( adress(2,260) )
565 ! UM: ALPHA-Factor
566  alpha = motrea( adress(2,261) )
567 ! UM: MOFAC-Factor
568  mofac = motrea( adress(2,262) )
569 
570  ! ################## !
571  ! LOGICAL KEYWORDS !
572  ! ################## !
573 ! INDEX 99 IS ALREADY USED FOR KEYWORD 'LIST OF FILES'
574 ! INDEX 54 IS ALREADY USED FOR KEYWORD 'DESCRIPTION OF LIBRARIES'
575 ! INDEX 57 IS ALREADY USED FOR KEYWORD 'DEFAULT EXECUTABLE'
576  ! SPHERICAL EQUATIONS HARD-CODED
577  ! ----------------------------------
578  spheri = .false.
579 
580 
581  ! COMPUTATION OF FALL VELOCITIES
582  ! ------------------------------------------
583  calwc = .false.
584  ! IF TROUVE: VELOCITIES ARE USER-DEFINED
585  ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
586  IF (trouve(2, 22).EQ.2) calwc = .true.
587 ! CV
588  ! SHIELDS PARAMETER
589  ! ------------------------------------------
590  calac = .false.
591  ! IF TROUVE
592  ! ~~~~~~~~~~~~~
593  IF (trouve(2, 6).EQ.2) calac = .true.
594 
595 
596  bilma = motlog( adress(3, 1) )
597  perma = motlog( adress(3, 2) )
598  bandec = motlog( adress(3, 3) )
599  valid = motlog( adress(3, 4) )
600 ! DTVAR = MOTLOG( ADRESS(3, 5) )
601  lumpi = motlog( adress(3, 6) )
602  susp = motlog( adress(3, 7) )
603  charr = motlog( adress(3, 8) )
604  houle = motlog( adress(3, 10) )
605  const_alayer = motlog( adress(3, 11) )
606  lcondis = motlog( adress(3, 12) )
607  lgrafed = motlog( adress(3, 13) )
608 ! USED TO CHECK SIS_FILES(SISPRE)%NAME
609  debu = motlog( adress(3, 14) )
610  imp_inflow_c = motlog( adress(3, 15) )
611  seccurrent = motlog( adress(3, 16) )
612  havesecfile = motlog( adress(3, 59) )
613  IF(code(1:9).EQ.'TELEMAC3D') seccurrent = .false.
614  unit = motlog( adress(3, 17) )
615  vf = motlog( adress(3,253) )
616  corr_conv = motlog( adress(3, 18) )
617  DO k=1,nsicla
618  sedco(k) = .false.
619  ENDDO
620  IF(dimens(3,19).GT.0) THEN
621  DO k=1,dimens(3,19)
622  sedco(k) = motlog( adress(3,19) + k-1 )
623  ENDDO
624  ENDIF
625  slide = motlog( adress(3, 20) )
626  dift = motlog( adress(3, 21) )
627  effpen = motlog( adress(3, 22) )
628  IF(.NOT.effpen) THEN
629  slopeff=0
630  devia=0
631  ENDIF
632 !
633  mixte=motlog(adress(3,24))
634 ! COUPLING WITH NESTOR
635  nestor=motlog(adress(3,25))
636 ! V6P1
637  kspred =motlog(adress(3,26))
638 !
639 ! Settling lag: determines choice between Rouse and Miles concentration profile
640 ! SET_LAG = TRUE : Miles
641 ! = FALSE: Rouse
642 !
643  set_lag = motlog(adress(3,27) )
644 ! STATIONARY MODE: calculate sediment transport without updating the bed.
645  stat_mode = motlog(adress(3,28) )
646 ! Checking the mesh
647  check_mesh = motlog(adress(3,29) )
648 ! NEW IMPLEMENTATION FOR CROSS-SECTION
649  doflux = motlog(adress(3,61) )
650  IF ( code(1:7) .NE. 'TELEMAC' ) THEN
651 ! SYMBOLIC LINEAR SOLVER FOR AD
652  ad_symblinsolv = motlog( adress(3,30) )
653 ! RESET TANGENTS ON ENTRY IN LINEAR SOLVER CG FOR AD
654  ad_linsolv_resetderiv = motlog( adress(3,31) )
655 ! ALLOW ADDITONAL INTERATIONS IN ITERATIVE LINEAR SOLVERS
656 ! SHOULD NOT BE USED IN PARALLEL MODE
657  ad_linsolv_derivative_convergence = motlog( adress(3,32) )
658  ENDIF
659  partel_concat = motlog(adress(3,62))
660  IF(ncsize.LE.1) partel_concat=.false.
661 !
662 ! ################################### !
663 ! CHARACTER STRING KEYWORDS !
664 ! ################################### !
665 !
666  titca = motcar( adress(4, 1) )(1:72)
667  sortis = motcar( adress(4, 2) )(1:72)
668  varim = motcar( adress(4, 3) )(1:72)
669  sis_files(sisgeo)%NAME=motcar( adress(4,6) )
670  IF(sis_files(sisgeo)%NAME(1:1).EQ.' ') THEN
671  WRITE(lu,*) 'THE FOLLOWING KEYWORD IS MANDATORY:'
672  WRITE(lu,*) 'GEOMETRY FILE (FICHIER DE GEOMETRIE)'
673  CALL plante(1)
674  stop
675  ENDIF
676  sis_files(siscli)%NAME=motcar( adress(4,9) )
677  IF(sis_files(siscli)%NAME(1:1).EQ.' ') THEN
678  WRITE(lu,*) 'THE FOLLOWING KEYWORD IS MANDATORY:'
679  WRITE(lu,*) 'BOUNDARY CONDITIONS FILE '//
680  & '(FICHIER DES CONDITIONS AUX LIMITES)'
681  CALL plante(1)
682  stop
683  ENDIF
684  sis_files(sishyd)%NAME=motcar( adress(4,29) )
685  sis_files(sispre)%NAME=motcar( adress(4,11) )
686  sis_files(sisres)%NAME=motcar( adress(4,12) )
687  sis_files(sisfon)%NAME=motcar( adress(4,16) )
688  sis_files(sisres)%FMT = motcar( adress(4,31) )(1:8)
689  IF(sis_files(sisres)%NAME(1:1).EQ.' ') THEN
690  WRITE(lu,*) 'THE FOLLOWING KEYWORD IS MANDATORY:'
691  WRITE(lu,*) 'RESULTS FILE (FICHIER DES RESULTATS)'
692  CALL plante(1)
693  stop
694  ENDIF
695  CALL majus(sis_files(sisres)%FMT)
696 ! RESULT FILE FORMAT FOR PREVIOUS SEDIMENTOLOGICAL
697 ! COMPUTATION...
698  sis_files(sispre)%FMT = motcar( adress(4,34) )(1:8)
699  CALL majus(sis_files(sispre)%FMT)
700 ! REFERENCE FILE FORMAT
701  sis_files(sisref)%FMT = motcar( adress(4,33) )(1:8)
702  CALL majus(sis_files(22)%FMT)
703 ! HYDRODYNAMIC FILE FORMAT
704  sis_files(sishyd)%FMT = motcar( adress(4,32) )(1:8)
705  CALL majus(sis_files(sishyd)%FMT)
706 ! WAVE FILE FORMAT (COUPLING WITH TOMAWAC)
707  sis_files(siscou)%FMT = motcar( adress(4,35) )(1:8)
708  CALL majus(sis_files(siscou)%FMT)
709  bingeosis = motcar( adress(4,18) )(1:3)
710  binhydsis = motcar( adress(4,19) )(1:3)
711  binpresis = motcar( adress(4,20) )(1:3)
712  binressis = motcar( adress(4,21) )(1:3)
713  sis_files(sisref)%NAME=motcar( adress(4,22) )
714  binrefsis = motcar( adress(4,23) )(1:3)
715 ! NESTOR STEERING FILE
716  sis_files(sinact)%NAME = motcar( adress(4,27) ) ! Nestor
717  sis_files(sinpol)%NAME = motcar( adress(4,40) ) ! Nestor
718  sis_files(sinref)%NAME = motcar( adress(4,41) ) ! Nestor
719  sis_files(sinrst)%NAME = motcar( adress(4,64) ) ! Nestor
720 ! ****** = MOTCAR( ADRESS(4,28) )
721 ! WAVE FILE
722  sis_files(siscou)%NAME=motcar( adress(4,30) )
723 ! SECTIONS
724  sis_files(sissec)%NAME=motcar( adress(4,36) )
725  sis_files(sisseo)%NAME=motcar( adress(4,37) )
726 ! FILE FOR LIQUID BOUNDARIES
727  sis_files(sisliq)%NAME=motcar( adress(4,38) )
728 ! GEOMETRY FILE FORMAT
729  sis_files(sisgeo)%FMT = motcar( adress(4,39) )(1:8)
730  CALL majus(sis_files(sisgeo)%FMT)
731 ! NAMES OF PRIVATE VARIABLES
732  n_names_priv = min(4,dimens(4,42))
733  IF(n_names_priv.GT.0) THEN
734  DO k=1,n_names_priv
735  names_prive(k) = motcar(adress(4,42)+k-1)(1:32)
736  ENDDO
737  ENDIF
738 !
739 ! ADDITIONAL DIFFERENTIATED VARIABLES
740 ! NAMES OF THE DIFFERENTIATED VARIABLES
741  n_names_advar = dimens(4,70)
742  nadvar = max( nadvar,n_names_advar ) ! WARNING HERE ?
743  IF(nadvar.GT.0) THEN
744  DO i=1,nadvar
745  WRITE(char2,'(I2)') i
746  names_advar(i) = 'DERIVATIVE '//adjustl(char2)//' '
747  & // '?? '
748  ENDDO
749  ENDIF
750  IF(n_names_advar.GT.0) THEN
751  DO k=1,n_names_advar
752  names_advar(k) = motcar(adress(4,70)+k-1)(1:32)
753  ENDDO
754  ENDIF
755 !
756 ! CVSM, But it's not Beautiful
757  tempvar = motcar(adress(4,51) )
759 ! FLUXLINEFILE
760  sis_files(sisflx)%NAME=motcar( adress(4,69) )
761 ! C-VSM RESULTS FILE
762  sis_files(vspres)%NAME=motcar( adress(4,53) )
763  sis_files(vspres)%FMT =motcar( adress(4,55) )(1:8)
764  CALL majus(sis_files(vspres)%FMT)
765 !
766  WRITE(lu,102)
767 102 FORMAT(1x,/,19x, '********************************************',/,
768  & 19x, '* LECDON: *',/,
769  & 19x, '* AFTER CALLING DAMOCLES *',/,
770  & 19x, '* CHECKING OF DATA READ *',/,
771  & 19x, '* IN THE STEERING FILE *',/,
772  & 19x, '********************************************',/)
773 !
774 !-----------------------------------------------------------------------
775 !
776 ! LOGICALS FOR OUTPUT VARIABLES
777 !-----------------------------------------------------------------------
778 !
781 !
782 ! ADDITIONAL DIFFERENTIATED VARIABLES
785 !
786  CALL sortie(sortis , mnemo , maxvar , sorleo )
787  CALL sortie(varim , mnemo , maxvar , sorimp )
788 !
789  DO i = 1, 4
790  IF ((npriv.LT.i).AND.
791  ! JWI 31/05/2012 - added 1 to include wave orbital velocities
792  & (sorleo(i+29+(nomblay+4)*nsicla+nomblay).OR.
793  & sorimp(i+29+(nomblay+4)*nsicla+nomblay))) THEN
794  npriv=max(npriv,i)
795  ENDIF
796  ENDDO
797 !
798 !-----------------------------------------------------------------------
799 !
800 ! CANCELS OUTPUT OF VARIABLES WHICH ARE NOT BUILT IN THIS CASE
801 !
802  IF(.NOT.susp) THEN
803 !V 7/09/2006 MIGHT WANT TO OUTPUT THE SUSPENDED COMPONENT IN BIJKER
804 ! SORIMP(24+4*NSICLA)=.FALSE.
805 ! SORIMP(25+4*NSICLA)=.FALSE.
806 ! SORIMP(26+4*NSICLA)=.FALSE.
807  ENDIF
808 !
809 ! JWI 31/05/2012 - added 1 to include wave orbital velocities
810  IF(.NOT.charr) THEN
811  sorleo(24+(nomblay+2)*nsicla)=.false.
812  sorleo(25+(nomblay+2)*nsicla)=.false.
813  sorleo(26+(nomblay+2)*nsicla)=.false.
814  sorimp(24+(nomblay+2)*nsicla)=.false.
815  sorimp(25+(nomblay+2)*nsicla)=.false.
816  sorimp(26+(nomblay+2)*nsicla)=.false.
817  ENDIF
818 ! JWI END
819 !
820 !-----------------------------------------------------------------------
821 !
822 ! CHECKS TETA VALUE
823 !
824  IF( teta.LT.0.d0.OR.teta.GT.1.d0) THEN
825  WRITE(lu,51)
826 51 FORMAT(/,1x,'BAD VALUE FOR TETA ! ',/
827  & ,1x,'TETA MUST BE WITHIN 0 AND 1 ')
828  CALL plante(1)
829  stop
830  ENDIF
831 !
832 ! INITIALISES MSK (MASKING VARIABLE)
833 ! FOR NOW MASKING IS ONLY DONE FOR ONE 'OPTION FOR THE TREATMENT
834 ! OF TIDAL FLATS'. SHOULD BE OFFERED AS AN OPTION FOR THE USER TO
835 ! CREATE ISLANDS IN THE FUTURE
836  msk = .false.
837  IF (.NOT.bandec) optban=0
838  IF (optban.EQ.2) msk = .true.
839 !
840 !-----------------------------------------------------------------------
841 !
842 ! CHECKS WHETHER THERE IS A REFERENCE FILE
843 !
844  IF (valid.AND.sis_files(sisref)%NAME.EQ.' ') THEN
845  valid=.false.
846  WRITE(lu,71)
847  WRITE(lu,*)
848 71 FORMAT(/,1x,'VALIDATION IS NOT POSSIBLE : ',/
849  & ,1x,'NO REFERENCE FILE ! ')
850  CALL plante(1)
851  stop
852  ENDIF
853 !
854 ! CHECKS THE NUMBER OF NSICLA
855 ! IF(NSICLA.GT.10) THEN
856  IF(nsicla.GT.nsiclm) THEN
857  WRITE(lu,81) nsiclm
858  WRITE(lu,*)
859 81 FORMAT(/,1x,'THE MAXIMUM NUMBER OF SEDIMENT CLASSES IS', i2)
860  CALL plante(1)
861  stop
862  ENDIF
863 ! CHECKS THE SUM OF INITIAL AVAI
864  sumavai = 0.d0
865  DO i=1,nsicla
866  sumavai = sumavai + ava0(i)
867  ENDDO
868  IF(abs(sumavai-1).GE.1.d-8) THEN
869  WRITE(lu,91)
870  WRITE(lu,*)
871 91 FORMAT(/,1x,'SUM OF SEDIMENT FRACTIONS IS NOT 1 ')
872  CALL plante(1)
873  stop
874  ENDIF
875 !
876 ! WARNING FOR THE CHOICE OF RIGID BED METHOD
877 !
878  IF(choix.GT.0.AND.choix.LT.4.AND.vf) THEN
879  WRITE(lu,201)
880  WRITE(lu,*)
881 201 FORMAT(/,1x,'FINITE VOLUMES CHOSEN: ',/
882  & ,1x,'METHOD 4 FOR RIGID BED WILL BE USED ')
883  choix=4
884  ENDIF
885  IF (choix.EQ.4.AND..NOT.vf) THEN
886  WRITE(lu,301)
887  WRITE(lu,*)
888 301 FORMAT(/,1x,'FINITE ELEMENTS CHOSEN: ',/
889  & ,1x,'METHOD 4 FOR RIGID BED CAN NOT BE USED, METHOD 3 US
890  &ED INSTEAD')
891  choix=3
892  ENDIF
893 !
894 ! CHECKS THAT THE BEDLOAD TRANSPORT FORMULATION AND THE HIDING
895 ! FACTOR FORMULATION CAN GO TOGETHER
896 !
897  IF ((hidfac.EQ.3.AND.icf.NE.6).OR.
898  & (hidfac.EQ.1.AND.icf.NE.1).OR.
899  & (hidfac.EQ.2.AND.icf.NE.1)) THEN
900  WRITE(lu,111)
901  WRITE(lu,*)
902 111 FORMAT(/,1x,'CHOICE OF TRANSPORT FORMULA AND HIDING FACTOR FORMU
903  &LATION NOT ALLOWED ')
904  CALL plante(1)
905  stop
906  ENDIF
907 !
908 ! WITHOUT AND WITH COUPLING, SOME CORRECTIONS
909 !
910  IF(code(1:7).EQ.'TELEMAC'.AND.
911  & sis_files(sishyd)%NAME(1:1).NE.' ') THEN
912  sis_files(sishyd)%NAME(1:1)=' '
913  WRITE(lu,113)
914 113 FORMAT(/,1x,'COUPLING: HYDRODYNAMIC FILE IGNORED')
915  ENDIF
916 !
917 ! COMPUTATION CONTINUED
918 !
919  IF(debu) THEN
920  IF(sis_files(sispre)%NAME(1:1).EQ.' ') THEN
921  WRITE(lu,313)
922 313 FORMAT(/,1x,'COMPUTATION CONTINUED:',/,
923  & 1x,'PREVIOUS SEDIMENTOLOGICAL FILE MISSING')
924  CALL plante(1)
925  stop
926  ENDIF
927  ELSE
928  IF(sis_files(sispre)%NAME(1:1).NE.' ') THEN
929  sis_files(sispre)%NAME(1:1)=' '
930  WRITE(lu,213)
931 213 FORMAT(/,1x,'NO COMPUTATION CONTINUED:',/,
932  & 1x,'PREVIOUS SEDIMENTOLOGICAL FILE IGNORED')
933  ENDIF
934  ENDIF
935 !
936 ! METHODS NOT CODED UP FOR SUSPENSION
937 ! -------------------------------------------
938 !
939  IF(susp) THEN
940  IF(resol.NE.adv_car .AND.resol.NE.adv_sup .AND.
941  & resol.NE.adv_lpo .AND.resol.NE.adv_nsc .AND.
942  & resol.NE.adv_psi .AND.resol.NE.adv_lpo_tf.AND.
943  & resol.NE.adv_nsc_tf.AND.resol.NE.adv_psi_tf ) THEN
944  WRITE(lu,303) resol
945 303 FORMAT(1x,'RESOLVING METHOD NOT IMPLEMENTED : ',1i6)
946  CALL plante(1)
947  stop
948  ENDIF
949  ENDIF
950 !
951  IF(.NOT.houle) sis_files(siscou)%NAME(1:1)=' '
952  IF(houle) THEN
953  IF(icf.NE.4.AND.icf.NE.5.AND.icf.NE.8.AND.icf.NE.9) THEN
954  WRITE(lu,1304) icf
955 1304 FORMAT(' TRANSPORT FORMULA',1i3,1x,
956  & 'DOES NOT TAKE WAVES INTO ACCOUNT,',/,1x,
957  & 'TRY 4, 5, 8 OR 9')
958  CALL plante(1)
959  stop
960  ENDIF
961  ENDIF
962 !
963 ! BEDLOAD AND SUSPENDED TRANSPORT COUPLING
964 ! ---------------------------------
965 !
966  IF((icf==30.OR.icf==3.OR.icf==9).AND.susp.AND.charr) THEN
967  WRITE(lu,1302) icf
968  CALL plante(1)
969  stop
970  ENDIF
971 1302 FORMAT('FOR THE FORMULA',1i3,/,1x,
972  & 'THE SUSPENSION TERM IS CALCULATED TWICE,'
973  & ,' WITH TOTAL LOAD FORMULA AND SUSPENSION ')
974 !
975 ! REFERENCE CONCENTRATION
976 !
977 ! MODIFICATION CV 31/12 IF(ICQ.EQ.2.AND.(PERCOU.NE.1.OR..NOT.CHARR)) THEN
978 !
979  IF(icq.EQ.2.AND.(percou.GT.1.OR..NOT.charr)) THEN
980  WRITE(lu,1402) icq
981 1402 FORMAT('FOR THE BIJKER REFERENCE CONCENTRATION',1i3,/,1x,
982  & 'BEDLOAD MUST BE COMPUTED, CHOOSE:',/,1x,
983  & 'BEDLOAD = YES')
984  CALL plante(1)
985  stop
986  ENDIF
987 !
988 ! CHECKS CONSISTENCY OF BEDLOAD LAWS
989 !
990 ! SOULSBY SLOPE EFFECT : REQUIRES A THRESHOLD FORMULA
991 !
992  IF(slopeff.EQ.2) THEN
993  IF(icf.NE.1) THEN
994  WRITE(lu,1404) icf
995 1404 FORMAT('BED-LOAD TRANSPORT FORMULA, HERE ICF=',1i3,/,1x,
996  & 'MUST HAVE A THRESHOLD',/,1x,
997  & 'IF FORMULA FOR SLOPE EFFECT=2 (SOULSBY)')
998  ENDIF
999  ENDIF
1000 !
1001 ! V6P0 : COHERENCE IF CONSOLIDATION MODEL IS USED
1002 ! VITCE AND CSF_VASE STEM FROM THE FIRST LAYER OF THE MULTI-LAYER MODEL
1003 ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1004 ! Mots cles supprimes vitesse critique d'erosion et concentration du lit
1005 ! CV : si premiere couche est vide cela n'est pas correct
1006 
1007 !
1008  IF(mixte) THEN
1009 !
1010 ! FILLS VOIDS WITH MUD:
1011 ! CV: verifier que la concentration en cohesif est non nulle
1012 !
1013  csf_sable= 1.d0
1014 !V: verrouiller les options
1015  nsicla=2
1016  sedco(1)=.false.
1017  sedco(2)=.true.
1018  charr=.false.
1019  susp=.true.
1020 !V
1021  ELSE
1022  csf_sable= (1.d0-xkv)
1023  ENDIF
1024 !
1025  IF((.NOT.mixte).AND.sedco(1)) THEN
1026  charr=.false.
1027  !SUSP=.TRUE. (In general, but not necessary)
1028  ENDIF
1029 !
1030  IF(nomblay.GT.nlaymax) THEN
1031  WRITE (lu,*) 'NUMBER OF BED LOAD MODEL LAYERS LARGER THAN '
1032  WRITE (lu,*) 'THE MAXIMUM PROGRAMMED VALUE OF ', nlaymax
1033  CALL plante(1)
1034  stop
1035  ENDIF
1036 ! IF(NOMBLAY.LT.2) THEN
1037 ! WRITE (LU,*) 'BEWARE: NUMBER OF BED LOAD MODEL LAYERS'
1038 ! WRITE (LU,*) '======= LOWER THAN THE DEFAULT VALUE OF 2'
1039 ! ENDIF
1040 !
1041 !----------------------------------------------------------------
1042 !
1043 ! V6P1: FOR THE BED FRICTION PREDICTOR USE LAW OF FRICTION 5 (NIKURADSE)
1044 !
1045  IF(kspred) kfrot=5
1046 !
1047 !----------------------------------------------------------------
1048 !
1049  RETURN
1050  END
double precision, dimension(:), allocatable, target cbor_classe
integer, parameter adv_nsc_tf
subroutine nomvar_sisyphe(TEXTE, TEXTPR, MNEMO, NSICLA, UNITE, MAXVAR, NPRIV, NOMBLAY, N_NAMES_PRIV, NAMES_PRIVE, NADVAR, NAMES_ADVAR)
Definition: nomvar_sisyphe.f:8
double precision, target phised
integer, parameter adv_psi
integer, parameter maxvar
integer, parameter adv_lpo
double precision, dimension(nsiclm) hidi
logical, dimension(maxvar) sorleo
integer, parameter adv_sup
double precision hmin_bedload
double precision, dimension(nsiclm) ava0
double precision, target xkv
subroutine lecdon_sisyphe(MOTCAR, FILE_DESC, PATH, NCAR, CODE, CAS_FILE, DICO_FILE)
Definition: lecdon_sisyphe.f:7
subroutine lecdon_split_outputpoints(INT_LIST, POINT_ARRAY, FULLOUTPUT)
logical, dimension(:), allocatable okcgl
subroutine read_submit(FILES, NFILES, SUBMIT, NMOT)
Definition: read_submit.f:7
double precision, dimension(nsiclm), target xwc
integer, dimension(3) mardat
double precision, target tprec
integer, dimension(100) cvsmoutput
integer, parameter maxkeyword
character(len=32), dimension(4) names_prive
integer, parameter nsiclm
integer, parameter adv_nsc
double precision, target partheniades
integer, dimension(:), allocatable ctrlsc
double precision, dimension(nlaymax) trans_mass
double precision, target delt
double precision, dimension(nsiclm) fd90
integer, parameter maxlu_sis
logical, dimension(maxvar) sorimp
double precision, dimension(:), pointer x
integer, parameter adv_psi_tf
double precision, dimension(nsiclm), target fdm
logical, dimension(nsiclm) sedco
subroutine sortie(CHAINE, MNEMO, NBRE, SORLEO)
Definition: sortie.f:7
double precision, dimension(:), allocatable soldis
integer, parameter adv_lpo_tf
logical, dimension(:), allocatable okqgl
subroutine damocle(ADRESS, DIMENS, NMAX, DOC, LLNG, LLU, MOTINT, MOTREA, MOTLOG, MOTCAR, MOTCLE, TROUVE, NFICMO, NFICDA, GESTD, MOTATT)
Definition: damocle.f:9
character(len=32), dimension(maxvar) names_advar
double precision, dimension(nlaymax) toce_vase
double precision, dimension(nlaymax) conc_vase
integer, parameter adv_car
character(len=8), dimension(maxvar) mnemo
double precision, target kspratio
double precision, target alpha
character(len=32), dimension(maxvar) textpr
integer, parameter nlaymax
double precision, dimension(nsiclm) cs0
integer, dimension(3) martim
subroutine majus(CHAINE)
Definition: majus.f:7
character(len=32), dimension(maxvar) texte
double precision, target beta2
double precision, dimension(nsiclm), target ac
double precision, target mpm
double precision, target csf_sable
type(bief_file), dimension(maxlu_sis), target sis_files