The TELEMAC-MASCARET system  trunk
lecdon_artemis.f
Go to the documentation of this file.
1 ! *************************
2  SUBROUTINE lecdon_artemis
3 ! *************************
4 !
5  &(file_desc,path,ncar,
6  & cas_file,dico_file)
7 !
8 !***********************************************************************
9 ! ARTEMIS V8P2
10 !***********************************************************************
11 !
12 !brief READS THE STEERING FILE THROUGH A DAMOCLES CALL.
13 !
14 !history J-M HERVOUET (LNH)
15 !+ 17/08/1994
16 !+
17 !+ LINKED TO BIEF 5.0
18 !
19 !history D. AELBRECHT (LNH)
20 !+ 20/04/1999
21 !+ V6P0
22 !+
23 !
24 !history C/PEYRARD (EDF)
25 !+ 18/03/2014
26 !+ V6P1 - V7P0
27 !+ NEW KEY WORDS
28 !
29 !history J-M HERVOUET (EDF LAB, LNHE)
30 !+ 18/05/2015
31 !+ V7P1
32 !+ Adding CHECK_MESH for the keyword 'CHECKING THE MESH'
33 !
34 !history Y AUDOUIN (LNHE)
35 !+ 25/05/2015
36 !+ V7P0
37 !+ Modification to comply with the hermes module
38 !
39 !history N.DURAND (HRW)
40 !+ November 2016
41 !+ V7P2
42 !+ Addition of new keywords in the steering file re: animation of the
43 !+ free surface. ANIMFS and ART_FILES(ARTAMP)%NAME read from file
44 !
45 !history N.DURAND (HRW)
46 !+ August 2017
47 !+ V7P3
48 !+ 1. Consolidation of nesting methods, meaning that CHAINTWC is now a
49 !+ choice of integers, and no longer a logical
50 !+ 2. NDTWC (now NDIR) and NFTWC (now NF) are no longer read in the
51 !+ steering file but directly extracted from the spectrum file
52 !+ 3. Addition of various checks before simulation can start
53 !+ 4. Call to ARTEMIS_CONSTANTS added.
54 !+ 5. Addition of new keywords in the steering file for nesting option 2:
55 !+ X_SFREF, Y_SFREF, ART_FILES(WACRES)%NAME
56 !
57 !history N.DURAND (HRW)
58 !+ December 2018
59 !+ V8P0
60 !+ Addition of new keywords in the steering file for nesting option 2:
61 !+ ART_FILES(WACLQD)%NAME
62 !
63 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
64 !| FILE_DESC |<--| STORES STRINGS 'SUBMIT' OF DICTIONARY
65 !| NCAR |-->| NUMBER OF LETTERS IN STRING PATH
66 !| PATH |-->| FULL PATH TO CODE DICTIONARY
67 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
68 !
72  USE interface_artemis, ex_lecdon_artemis => lecdon_artemis
73 !
75  IMPLICIT NONE
76 !
77 !
78 !+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+
79 !
80  CHARACTER(LEN=PATH_LEN), INTENT(INOUT) :: FILE_DESC(4,maxkeyword)
81  INTEGER, INTENT(IN) :: NCAR
82  CHARACTER(LEN=PATH_LEN), INTENT(IN) :: PATH
83 ! API
84  CHARACTER(LEN=PATH_LEN), INTENT(IN) :: DICO_FILE
85  CHARACTER(LEN=PATH_LEN), INTENT(IN) :: CAS_FILE
86 !
87 !+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+
88 !
89  LOGICAL LSTOP
90  INTEGER I,KK
91 !
92  CHARACTER(LEN=8) MNEMO(maxvar)
93 !
94 !-----------------------------------------------------------------------
95 !
96 ! ARRAYS USED IN THE DAMOCLES CALL
97 !
98  INTEGER ADRESS(4,maxkeyword),DIMENS(4,maxkeyword)
99  DOUBLE PRECISION MOTREA(maxkeyword)
100  INTEGER MOTINT(maxkeyword)
101  LOGICAL MOTLOG(maxkeyword)
102  CHARACTER(LEN=PATH_LEN) MOTCAR(maxkeyword)
103  CHARACTER(LEN=72) MOTCLE(4,maxkeyword,2)
104  INTEGER TROUVE(4,maxkeyword)
105  LOGICAL DOC
106  CHARACTER(LEN=PATH_LEN) :: NOM_CAS
107  CHARACTER(LEN=PATH_LEN) :: NOM_DIC
108  INTEGER :: ID_DICO, ID_CAS
109 ! EXTRA DECLARATIONS FOR OUTPUT OF AMPLITUDE AND PHASE FILE
110  CHARACTER(LEN=2) :: CLT ! PERIOD COMPONENT NUMBER
111 !
112 ! END OF DECLARATIONS FOR DAMOCLES CALL :
113 !
114 !-----------------------------------------------------------------------
115 !
116 ! INITIALISES THE VARIABLES FOR DAMOCLES CALL :
117 !
118  DO kk=1,maxkeyword
119 !
120 ! A FILENAME NOT GIVEN BY DAMOCLES WILL BE RECOGNIZED AS A WHITE SPACE
121 ! (IT MAY BE THAT NOT ALL COMPILERS WILL INITIALISE LIKE THAT)
122 !
123  motcar(kk)(1:1)=' '
124 !
125  dimens(1,kk) = 0
126  dimens(2,kk) = 0
127  dimens(3,kk) = 0
128  dimens(4,kk) = 0
129 !
130  ENDDO
131 !
132 ! WRITES OUT INFO
133 !
134  doc = .false.
135 !
136 !-----------------------------------------------------------------------
137 ! OPENS DICTIONNARY AND STEERING FILES
138 !-----------------------------------------------------------------------
139 !
140  IF(ncar.GT.0) THEN
141 !
142  nom_dic=path(1:ncar)//'ARTDICO'
143  nom_cas=path(1:ncar)//'ARTCAS'
144 !
145  ELSE
146 !
147  nom_dic='ARTDICO'
148  nom_cas='ARTCAS'
149 !
150  ENDIF
151  IF((cas_file(1:1).NE.' ').AND.(dico_file(1:1).NE.' ')) THEN
152  nom_dic=dico_file
153  nom_cas=cas_file
154  ENDIF
155 !
156  CALL get_free_id(id_dico)
157  OPEN(id_dico,file=nom_dic,form='FORMATTED',action='READ')
158  CALL get_free_id(id_cas)
159  OPEN(id_cas,file=nom_cas,form='FORMATTED',action='READ')
160 !
161  CALL damocle( adress , dimens , maxkeyword, doc, lng , lu ,
162  & motint , motrea , motlog , motcar , motcle ,
163  & trouve , id_dico, id_cas , .false. , file_desc)
164 !
165 !-----------------------------------------------------------------------
166 ! CLOSES DICTIONNARY AND STEERING FILES
167 !-----------------------------------------------------------------------
168 !
169  CLOSE(id_dico)
170  CLOSE(id_cas)
171 !
172 ! DECODES 'SUBMIT' CHAINS
173 !
174  CALL read_submit(art_files,maxlu_art,file_desc,maxkeyword)
175 !
176 !-----------------------------------------------------------------------
177 !
178  DO i=1,maxlu_art
179  IF(art_files(i)%TELNAME.EQ.'ARTGEO') THEN
180  artgeo=i
181  ELSEIF(art_files(i)%TELNAME.EQ.'ARTCAS') THEN
182  artcas=i
183  ELSEIF(art_files(i)%TELNAME.EQ.'ARTCLI') THEN
184  artcli=i
185  ELSEIF(art_files(i)%TELNAME.EQ.'ARTFON') THEN
186  artfon=i
187  ELSEIF(art_files(i)%TELNAME.EQ.'ARTRES') THEN
188  artres=i
189  ELSEIF(art_files(i)%TELNAME.EQ.'ARTREF') THEN
190  artref=i
191  ELSEIF(art_files(i)%TELNAME.EQ.'ARTBI1') THEN
192  artbi1=i
193  ELSEIF(art_files(i)%TELNAME.EQ.'ARTBI2') THEN
194  artbi2=i
195  ELSEIF(art_files(i)%TELNAME.EQ.'ARTFO1') THEN
196  artfo1=i
197  ELSEIF(art_files(i)%TELNAME.EQ.'ARTFO2') THEN
198  artfo2=i
199  ELSEIF(art_files(i)%TELNAME.EQ.'ARTRBI') THEN
200  artrbi=i
201  ELSEIF(art_files(i)%TELNAME.EQ.'ARTRFO') THEN
202  artrfo=i
203  ELSEIF(art_files(i)%TELNAME.EQ.'ARTAMP') THEN
204  artamp=i
205  ELSEIF(art_files(i)%TELNAME.EQ.'WACSPE') THEN
206  wacspe=i
207  ELSEIF(art_files(i)%TELNAME.EQ.'WACRES') THEN
208  wacres=i
209  ELSEIF(art_files(i)%TELNAME.EQ.'WACLQD') THEN
210  waclqd=i
211  ENDIF
212  ENDDO
213 !
214 !-----------------------------------------------------------------------
215 !
216 ! SETTING CONSTANTS (PI AND RELATED)
217 !
218  CALL artemis_constants
219 !
220 !-----------------------------------------------------------------------
221 !
222 ! ASSIGNS THE STEERING FILE VALUES TO THE PARAMETER FORTRAN NAME :
223 !
224 !-----------------------------------------------------------------------
225 !
226 ! INTEGER KEYWORDS :
227 !
228  leoprd = motint( adress(1, 1) )
229  lisprd = motint( adress(1, 2) )
230  slvart%NITMAX = motint( adress(1,3) )
231  slvart%PRECON = motint( adress(1,4) )
232 ! DISESP = MOTINT( ADRESS(1, 5) )
233 ! STDGEO = MOTINT( ADRESS(1, 6) )
234 ! STDRES = MOTINT( ADRESS(1, 7) )
235  slvart%SLV = motint( adress(1,8) )
236  lisfon = motint( adress(1, 9) )
237  npale = motint( adress(1,10) )
238  ndale = motint( adress(1,11) )
239  ibreak = motint( adress(1,12) )
240  nitdis = motint( adress(1,13) )
241  regido = motint( adress(1,14) )
242  formfr = motint( adress(1,15) )
243  slvart%KRYLOV = motint( adress(1,16) )
244  lvmac = motint( adress(1,17) )
245 ! KFROT = MOTINT( ADRESS(1,18) )
246  optass = motint( adress(1,19) )
247  produc = motint( adress(1,20) )
248  mardat(1) = motint( adress(1,21) )
249  mardat(2) = motint( adress(1,21) + 1 )
250  mardat(3) = motint( adress(1,21) + 2 )
251  martim(1) = motint( adress(1,22) )
252  martim(2) = motint( adress(1,22) + 1 )
253  martim(3) = motint( adress(1,22) + 2 )
254  npriv = motint( adress(1,23) )
255 ! PLACEHOLDER FOR NCSIZE:
256 ! NCSIZE = MOTINT( ADRESS(1,24) )
257 ! ORIGIN COORDINATES
258  i_orig = motint( adress(1,25) )
259  j_orig = motint( adress(1,25)+1 )
260 ! RAPIDLY VARYING TOPOGRAPHY
261  ipentco = motint( adress(1,26) )
262 ! MAX ITERATION ON TETAP
263  nittp = motint( adress(1,27) )
264 ! DEBUGGER KEYWORD
265  debug = motint( adress(1,28) )
266 ! 2 PLACEHOLDERS FOR NITFS AND NADVAR:
267 ! NITFS = MOTINT( ADRESS(1,29) )
268 ! NADVAR = MOTINT( ADRESS(1,30) )
269  chaintwc= motint( adress(1,31) )
270 
271 ! FOR THE MOMENT WTITES TO FILE FROM THE START OF SIMULATION
272  ptinig = 0
273  ptinil = 0
274 !
275 ! REAL KEYWORDS :
276 !
277  per = motrea( adress(2, 1) )
278  tetah = motrea( adress(2, 2) )
279  grav = motrea( adress(2, 3) )
280  slvart%ZERO = motrea( adress(2, 4) )
281  slvart%EPS = motrea( adress(2, 5) )
282 ! HMIN = MOTREA( ADRESS(2, 6) )
283  cotini = motrea( adress(2, 7) )
284  hautin = motrea( adress(2, 8) )
285  perdeb = motrea( adress(2, 9) )
286  perfin = motrea( adress(2,10) )
287  perpas = motrea( adress(2,11) )
288  perpic = motrea( adress(2,12) )
289  gamma = motrea( adress(2,13) )
290  tetmin = motrea( adress(2,14) )
291  tetmax = motrea( adress(2,15) )
292  expos = motrea( adress(2,16) )
293 ! RELAX = MOTREA( ADRESS(2,17) )
294  epsdis = motrea( adress(2,18) )
295  reldis = motrea( adress(2,19) )
296  alfabj = motrea( adress(2,20) )
297  gammas = motrea( adress(2,21) )
298  kdally = motrea( adress(2,22) )
299  gdally = motrea( adress(2,23) )
300  visco = motrea( adress(2,24) )
301  diam90 = motrea( adress(2,25) )
302  diam50 = motrea( adress(2,26) )
303  mvsed = motrea( adress(2,27) )
304  mveau = motrea( adress(2,28) )
305  fwcoef = motrea( adress(2,29) )
306  ricoef = motrea( adress(2,30) )
307  ffon = motrea( adress(2,31) )
308  pmin = motrea( adress(2,32) )
309  pmax = motrea( adress(2,33) )
310  depref = motrea( adress(2,34) )
311  x_phref = motrea( adress(2,35) )
312  y_phref = motrea( adress(2,35)+1 )
313  epsdir = motrea( adress(2,36) )
314  epstp = motrea( adress(2,37) )
315  reltp = motrea( adress(2,38) )
316  tpstwc = motrea( adress(2,39) )
317 ! PLACEHOLDERS FOR TINIFS AND DTFS:
318 ! TINIFS = MOTREA( ADRESS(2,40) )
319 ! DTFS = MOTREA( ADRESS(2,41) )
320 ! FOR THE F SPECTRUM REFERENCE POINT
321  x_sfref = motrea( adress(2,42) )
322  y_sfref = motrea( adress(2,42)+1 )
323 !
324 ! LOGICAL KEYWORDS :
325 !
326  listin = motlog( adress(3, 1) )
327  infogr = motlog( adress(3, 2) )
328  balaye = motlog( adress(3, 3) )
329  alemon = motlog( adress(3, 4) )
330  alemul = motlog( adress(3, 5) )
331  deferl = motlog( adress(3, 6) )
332  frotte = motlog( adress(3, 7) )
333  entfw = motlog( adress(3, 8) )
334  entreg = motlog( adress(3, 9) )
335  entrug = motlog( adress(3, 10) )
336  lishou = motlog( adress(3, 11) )
337  valid = motlog( adress(3, 12) )
338  courant = motlog( adress(3, 13) )
339  langauto = motlog( adress(3, 14) )
340  lphaseauto= motlog( adress(3, 15) )
341 ! SPHERICAL EQUATIONS, HARD-CODED
342  spheri = .false.
343  check_mesh= motlog( adress(3, 16) )
344  animfs = motlog( adress(3, 17) )
345  partel_concat = motlog(adress(3,18))
346  IF(ncsize.LE.1) partel_concat=.false.!
347 ! STRING KEYWORDS : SOME ARE USED BY THE LAUNCHING
348 ! PROCEDURE
349 !
350  titcas = motcar( adress(4, 1) )(1:72)
351  vardes = motcar( adress(4, 2) )(1:72)
352  CALL majus(vardes)
353  varimp = motcar( adress(4, 3) )(1:72)
354  CALL majus(varimp)
355 ! FROM 4 TO 5: READ AND USED BY PRECOS
356  art_files(artgeo)%NAME=motcar( adress(4,6) )
357  IF(art_files(artgeo)%NAME(1:1).EQ.' ') THEN
358  WRITE(lu,*) 'THE FOLLOWING KEYWORD IS MANDATORY:'
359  WRITE(lu,*) 'GEOMETRY FILE (FICHIER DE GEOMETRIE)'
360  CALL plante(1)
361  stop
362  ENDIF
363 ! PLACEHOLDERS:
364 ! NOMFOR = MOTCAR( ADRESS(4, 7) )
365 ! NOMCAS = MOTCAR( ADRESS(4, 8) )
366 ! ART_FILES(ARTCAS)%NAME=MOTCAR( ADRESS(4,8) )
367  art_files(artcli)%NAME=motcar( adress(4,9) )
368  IF(art_files(artcli)%NAME(1:1).EQ.' ') THEN
369  WRITE(lu,*) 'THE FOLLOWING KEYWORD IS MANDATORY:'
370  WRITE(lu,*) 'BOUNDARY CONDITIONS FILE '//
371  & '(FICHIER DES CONDITIONS AUX LIMITES)'
372  CALL plante(1)
373  stop
374  ENDIF
375 ! WRITE(*,*) 'IN LECDON ',ART_FILES(ARTGEO)%NAME
376  art_files(artres)%NAME=motcar( adress(4,10) )
377  IF(art_files(artres)%NAME(1:1).EQ.' ') THEN
378  WRITE(lu,*) 'THE FOLLOWING KEYWORD IS MANDATORY:'
379  WRITE(lu,*) 'RESULTS FILE (FICHIER DES RESULTATS)'
380  CALL plante(1)
381  stop
382  ENDIF
383 ! FROM 11 TO 14 : READ AND USED BY PRECOS
384  art_files(artfon)%NAME=motcar( adress(4,15) )
385  art_files(artbi1)%NAME=motcar( adress(4,16) )
386  art_files(artbi2)%NAME=motcar( adress(4,17) )
387  art_files(artfo1)%NAME=motcar( adress(4,18) )
388  art_files(artfo2)%NAME=motcar( adress(4,19) )
389  art_files(artrbi)%NAME=motcar( adress(4,20) )
390  art_files(artrfo)%NAME=motcar( adress(4,21) )
391 ! FROM 22 TO 23 : READ AND USED BY PRECOS OR NOT USED
392  cdtini = motcar( adress(4,24) )(1:72)
393  CALL majus(cdtini)
394 ! BINGEO = MOTCAR( ADRESS(4,25) )(1:3)
395 ! CALL MAJUS(BINGEO)
396 ! BINRES = MOTCAR( ADRESS(4,26) )(1:3)
397 ! CALL MAJUS(BINRES)
398 ! ACCOUNTNUMBER = MOTCAR( ADRESS(4, 27) )
399  art_files(artref)%NAME=motcar( adress(4,28) )
400 ! GEOMETRY FILE FORMAT
401  art_files(artgeo)%FMT = motcar( adress(4,29) )(1:8)
402  CALL majus(art_files(artgeo)%FMT)
403 ! RESULTS FILE FORMAT
404  art_files(artres)%FMT = motcar( adress(4,30) )(1:8)
405  CALL majus(art_files(artres)%FMT)
406 ! REFERENCE FILE FORMAT
407  art_files(artref)%FMT = motcar( adress(4,31) )(1:8)
408  CALL majus(art_files(artref)%FMT)
409 ! BINARY FILE 1 FORMAT
410  art_files(artbi1)%FMT = motcar( adress(4,32) )(1:8)
411  CALL majus(art_files(artbi1)%FMT)
412 ! BINARY FILE 2 FORMAT
413  art_files(artbi2)%FMT = motcar( adress(4,33) )(1:8)
414  CALL majus(art_files(artbi2)%FMT)
415  IF(animfs) THEN
416 ! PHASE AND AMPLITUDE FILENAME
417  art_files(artamp)%NAME=motcar( adress(4,34) )
418 ! PHASE AND AMPLITUDE FILE FORMAT
419  art_files(artamp)%FMT = motcar( adress(4,35) )(1:8)
420  CALL majus(art_files(artamp)%FMT)
421  ENDIF
422  IF(chaintwc.GE.1) THEN
423 ! TOMAWAC OUTER SPECTRAL FILENAME
424  art_files(wacspe)%NAME=motcar( adress(4,36) )
425 ! TOMAWAC OUTER SPECTRAL FILE FORMAT
426  art_files(wacspe)%FMT = motcar( adress(4,37) )(1:8)
427  CALL majus(art_files(wacspe)%FMT)
428  ENDIF
429 ! PLACEHOLDERS:
430 ! NOMFST = MOTCAR( ADRESS(4, 38) )
431 ! FREE SURFACE FILE FORMAT = MOTCAR( ADRESS(4, 39) )
432 ! NAME\_ADVAR = MOTCAR( ADRESS(4, 40) )
433  IF(chaintwc.EQ.2) THEN
434 ! TOMAWAC OUTER RESULT FILENAME
435  art_files(wacres)%NAME=motcar( adress(4,41) )
436 ! TOMAWAC OUTER RESULT FILE FORMAT
437  art_files(wacres)%FMT = motcar( adress(4,42) )(1:8)
438  CALL majus(art_files(wacres)%FMT)
439 ! LIQUID BOUNDARY FILENAME, DERIVED FROM TOMAWAC OUTER RESULT FILE
440  art_files(waclqd)%NAME=motcar( adress(4,43) )
441 ! LIQUID BOUNDARY FILE FORMAT
442  art_files(waclqd)%FMT = motcar( adress(4,44) )(1:8)
443  CALL majus(art_files(waclqd)%FMT)
444  ENDIF
445 ! PLACEHOLDERS:
446 ! DICTIONARY = MOTCAR( ADRESS(4, 4) )
447 ! PARTITIONING TOOL = MOTCAR( ADRESS(4, 5) )
448 !
449  IF(listin) THEN
450  WRITE(lu,102)
451  ENDIF
452 102 FORMAT(1x,/,19x, '********************************************',/,
453  & 19x, '* LECDON: *',/,
454  & 19x, '* AFTER CALLING DAMOCLES *',/,
455  & 19x, '* CHECKING OF DATA READ *',/,
456  & 19x, '* IN THE STEERING FILE *',/,
457  & 19x, '********************************************',/)
458 !
459 !-----------------------------------------------------------------------
460 ! NAME OF THE VARIABLES FOR THE RESULTS AND GEOMETRY FILES :
461 !-----------------------------------------------------------------------
462 !
463  mnemo = ' '
464  CALL nomvar_artemis(texte,textpr,mnemo)
465 !
466 ! LOGICAL ARRAY FOR OUTPUT
467 !
468  sorleo = .false.
469  sorimp = .false.
470  CALL sortie(vardes , mnemo , 100 , sorleo )
471  CALL sortie(varimp , mnemo , 100 , sorimp )
472 !
473 !-----------------------------------------------------------------------
474 !
475 ! IN THE CASE OF A PERIOD SWEEPING, THE FIRST PERIOD TO BE COMPUTED
476 ! IS PERDEB
477 !
478  IF (balaye) per = perdeb
479 !
480 !-----------------------------------------------------------------------
481 !
482 ! IF NOTHING IS TO BE WRITTEN TO THE LISTING FILE, THERE SHOULD BE
483 ! NO INFORMATION WRITTEN ABOUT THE SOLVEUR
484 !
485  IF (.NOT.listin) infogr = .false.
486 !
487 !-----------------------------------------------------------------------
488 !
489 ! FOR RANDOM SEA COMPUTATIONS, INHIBITS THE GRAPHICAL AND LISTING
490 ! OUTPUT OF THE PHASES, SPEEDS, FREE SURFACE ELEVATION, AND
491 ! POTENTIAL BECAUSE THEY DO NOT MEAN ANYTHING.
492 !
493 ! BUT WRITES OUT AN AVERAGE WAVE NUMBER, AN AVERAGE PHASE CELERITY
494 ! AND AN AVERAGE GROUP VELOCITY, COMPUTED FROM THE MEAN WAVE
495 ! PERIOD T01. ALSO WRITES OUT AN AVERAGE INCIDENCE
496 !
497 ! MOREOVER, CHECKS THAT THE NUMBER OF DISCRETISED PERIODS AND
498 ! DIRECTIONS IS AT LEAST 1
499 !
500  IF(alemon .OR. alemul) THEN
501 !
502 ! 2 : PHASE
503 !
504  sorleo( 2) = .false.
505  sorimp( 2) = .false.
506 !
507 ! 3 AND 4 : U0 AND V0
508 !
509  sorleo( 3) = .false.
510  sorimp( 3) = .false.
511  sorleo( 4) = .false.
512  sorimp( 4) = .false.
513 !
514 ! 5 : FREE SURFACE
515 !
516  sorleo( 5) = .false.
517  sorimp( 5) = .false.
518 !
519 ! 11 AND 12 : REAL AND IMAGINARY PARTS OF THE POTENTIAL
520 !
521  sorleo(11) = .false.
522  sorimp(11) = .false.
523  sorleo(12) = .false.
524  sorimp(12) = .false.
525 !
526  IF (npale.LE.0) npale = 1
527  IF (ndale.LE.0) ndale = 1
528 !
529  ELSE
530 !
531 ! 17, 18 AND 19 : T01, T02 AND TM
532 !
533  sorleo(17) = .false.
534  sorimp(17) = .false.
535  sorleo(18) = .false.
536  sorimp(18) = .false.
537  sorleo(19) = .false.
538  sorimp(19) = .false.
539  ENDIF
540 !
541 ! NO PERIOD SWEEPING FOR RANDOM SEAS
542 !
543  IF (alemon .OR. alemul) balaye = .false.
544 !
545 ! NDALE=1 FOR MONO-DIRECTIONAL SEAS
546 !
547  IF (alemon .AND. .NOT.alemul) ndale = 1
548 !
549 ! NPALE=1 FOR REGULAR SEAS
550 !
551  IF (.NOT.alemon .AND. .NOT.alemul .AND. .NOT.balaye) npale = 1
552 !
553 ! IF AMPLITUDE AND PHASE FILE IS REQUIRED, NAME OF VARIABLES IN TEXTANIM
554 ! AND LOGICAL ARRAY FOR OUTPUT IN SORNIM
555 !
556  IF (animfs) THEN
557  DO i=1,ndale
558  WRITE(clt,205) i
559 205 FORMAT(i2.2)
560  textanim(2*i-1) = 'WAVE HEIGHT D'//clt//' (M)'
561  textanim(2*i) = 'WAVE PHASE D'//clt//' (RAD)'
562 !
563  sornim(2*i-1) = .true.
564  sornim(2*i) = .true.
565  ENDDO
566  ENDIF
567 !
568 !-----------------------------------------------------------------------
569 !
570 ! DOES NOT WRITE OUT 'PRIVE' VARIABLES IF THEY'VE NOT BEEN ALLOCATED
571 !
572  lstop = .false.
573  DO i=1,4
574  IF ((sorleo(12+i).OR.sorimp(12+i)).AND.(npriv.LT.i)) THEN
575  WRITE(lu,17) i,npriv
576  lstop = .true.
577  ENDIF
578  ENDDO
579  IF (lstop) THEN
580  CALL plante(1)
581  stop
582  ENDIF
583  17 FORMAT(1x,'PRIVATE ARRAY ',1i1,' CANNOT BE USED '
584  & ,1x,'BECAUSE IT WAS NOT ALLOCATED.',/
585  & ,1x,'CHECK ''NPRIV'' (AT THE TIME BEING ',1i1,' ).',/)
586 !
587 !-----------------------------------------------------------------------
588 !
589 ! WRITES OUT THE TITLE
590 !
591  IF(listin) THEN
592  WRITE(lu,3001) titcas
593 3001 FORMAT(/1x,'NAME OF THE STUDY :',1x,a72,/)
594  ENDIF
595 !
596 !-----------------------------------------------------------------------
597 !
598 ! NO TIDAL FLATS ==> MSK = .FALSE.
599 !
600  msk = .false.
601 !
602 !-----------------------------------------------------------------------
603 !
604 ! COMPULSORY CHOICE OF KEYWORD FOR THE DIRECT SOLVEUR
605 !
606  IF( (slvart%SLV.EQ.8).AND.(optass.NE.3) ) THEN
607  WRITE(lu,3003)
608 !
609 3003 FORMAT(1x,'WITH DIRECT SYSTEM SOLVER, EDGE-BASED STORAGE',/,1x,
610  & 'IS MANDATORY',///)
611  CALL plante(1)
612  stop
613  ENDIF
614 !
615 ! USE OF DIRECT SOLVEUR IS NOT POSSIBLE WITH PARALLELISM
616 !
617  IF(ncsize.GT.1) THEN
618  IF(slvart%SLV.EQ.8) THEN
619  WRITE(lu,3005)
620 3005 FORMAT(1x,'WITH PARALLELISM,',/,1x,
621  & 'NO DIRECT SYSTEM SOLVER 8 USE 9',///)
622  CALL plante(1)
623  stop
624  ENDIF
625  ELSE
626  IF(slvart%SLV.EQ.9) THEN
627  WRITE(lu,3006)
628 3006 FORMAT(1x,'WITH NO PARALLELISM,',/,1x,
629  & 'NO DIRECT SYSTEM SOLVER 9 USE 8',///)
630  CALL plante(1)
631  stop
632  ENDIF
633  ENDIF
634 !
635 !-----------------------------------------------------------------------
636 !
637 ! TESTS PMIN%MAX and TETMIN%MAX
638 !
639 ! CHECKS IF TETMIN<TETMAX
640  IF(tetmax.LE.tetmin) THEN
641  WRITE(lu,*) 'LECDON_ARTEMIS:'
642  WRITE(lu,*) 'MAXIMUM ANGLE OF PROPAGATION < MINIMUM ANGLE'
643  WRITE(lu,*) '=> CHANGE MAXIMUM OR MINIMUM ANGLE OF PROPAGATION'
644  WRITE(lu,*) ' IN THE USER PARAMETER FILE. '
645  CALL plante(1)
646  stop
647  ENDIF
648 ! CHECKS IF PMIN<PMAX
649  IF(pmax.LE.pmin) THEN
650  WRITE(lu,*) 'LECDON_ARTEMIS:'
651  WRITE(lu,*) 'MAXIMUM SPECTRAL PERIOD < MINIMUM SPECTRAL PERIOD'
652  WRITE(lu,*) '=> CHANGE MAXIMUM OR MINIMUM SPECTRAL PERIOD '
653  WRITE(lu,*) ' IN THE USER PARAMETER FILE. '
654  CALL plante(1)
655  stop
656  ENDIF
657 !
658 !-----------------------------------------------------------------------
659 !
660 ! TESTS SUPPLY OF REQUIRED INPUT FOR NESTING
661 !
662  IF(chaintwc.GE.1) THEN
663  IF (art_files(wacspe)%NAME .EQ. '') THEN
664  WRITE(lu,3007)
665 3007 FORMAT(/1x,'PLEASE SUPPLY TOMAWAC OUTER SPECTRAL FILE',/)
666  stop
667  ENDIF
668  ENDIF
669  IF(chaintwc.EQ.2) THEN
670  IF (art_files(wacres)%NAME .EQ. '') THEN
671  WRITE(lu,3009)
672 3009 FORMAT(/1x,'PLEASE SUPPLY TOMAWAC OUTER RESULT FILE',/)
673  stop
674  ENDIF
675  IF (art_files(waclqd)%NAME .EQ. '') THEN
676  WRITE(lu,3011)
677 3011 FORMAT(/1x,'PLEASE SUPPLY TOMAWAC LIQUID BOUNDARY FILE',/)
678  stop
679  ENDIF
680  IF (x_sfref .LT. -99999.) THEN
681  WRITE(lu,3013)
682 3013 FORMAT(/1x,'PLEASE SUPPLY X FOR THE REFERENCE F SPECTRUM',/)
683  stop
684  ENDIF
685  IF (y_sfref .LT. -99999.) THEN
686  WRITE(lu,3015)
687 3015 FORMAT(/1x,'PLEASE SUPPLY Y FOR THE REFERENCE F SPECTRUM',/)
688  stop
689  ENDIF
690  ENDIF
691 !
692 !-----------------------------------------------------------------------
693 !
694  RETURN
695  END
logical, dimension(maxvar) sorimp
character(len=32), dimension(maxvar) texte
logical, dimension(maxvar) sornim
integer, dimension(3) martim
integer, parameter maxlu_art
subroutine read_submit(FILES, NFILES, SUBMIT, NMOT)
Definition: read_submit.f:7
integer, parameter maxkeyword
subroutine nomvar_artemis(TEXTE, TEXTPR, MNEMO)
Definition: nomvar_artemis.f:7
integer, parameter maxvar
integer, dimension(3) mardat
logical, dimension(maxvar) sorleo
character(len=32), dimension(maxvar) textanim
character(len=32), dimension(maxvar) textpr
type(bief_file), dimension(maxlu_art), target art_files
subroutine sortie(CHAINE, MNEMO, NBRE, SORLEO)
Definition: sortie.f:7
double precision, dimension(:), pointer x
character(len=72) titcas
subroutine damocle(ADRESS, DIMENS, NMAX, DOC, LLNG, LLU, MOTINT, MOTREA, MOTLOG, MOTCAR, MOTCLE, TROUVE, NFICMO, NFICDA, GESTD, MOTATT)
Definition: damocle.f:9
subroutine lecdon_artemis(FILE_DESC, PATH, NCAR, CAS_FILE, DICO_FILE)
Definition: lecdon_artemis.f:8
subroutine artemis_constants
character(len=72) vardes
character(len=72) cdtini
subroutine majus(CHAINE)
Definition: majus.f:7
character(len=72) varimp