My Project
Functions/Subroutines
swanser.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine readxy (NAMX, NAMY, XX, YY, KONT, XSTA, YSTA)
 
subroutine txpbla (TEXT, IF, IL)
 
character *20 function numstr (IVAL, RVAL, FORM)
 
subroutine swcopi (IARR1, IARR2, LENGTH)
 
subroutine kscip1 (MMT, SIG, D, K, CG, N, ND)
 
subroutine kscip2 (S, D, CG)
 
character *20 function intstr (IVAL)
 
real function gamma (XX)
 
function gammln (XX)
 
subroutine sshape (ACLOC, SPCSIG, SPCDIR, FSHAPL, DSHAPL)
 
real function degcnv (DEGREE)
 
subroutine swanout (DEP2)
 

Function/Subroutine Documentation

◆ degcnv()

real function degcnv ( real  DEGREE)

Definition at line 1225 of file swanser.f90.

1225 ! *
1226 !***********************************************************************
1227 !
1228  USE swcomm3
1229 
1230  IMPLICIT NONE
1231 !
1232 !
1233 ! --|-----------------------------------------------------------|--
1234 ! | Delft University of Technology |
1235 ! | Faculty of Civil Engineering |
1236 ! | Environmental Fluid Mechanics Section |
1237 ! | P.O. Box 5048, 2600 GA Delft, The Netherlands |
1238 ! | |
1239 ! | Programmers: R.C. Ris, N. Booij, |
1240 ! | IJ.G. Haagsma, A.T.M.M. Kieftenburg, |
1241 ! | M. Zijlema, E.E. Kriezi, |
1242 ! | R. Padilla-Hernandez, L.H. Holthuijsen |
1243 ! --|-----------------------------------------------------------|--
1244 !
1245 !
1246 ! SWAN (Simulating WAves Nearshore); a third generation wave model
1247 ! Copyright (C) 2004-2005 Delft University of Technology
1248 !
1249 ! This program is free software; you can redistribute it and/or
1250 ! modify it under the terms of the GNU General Public License as
1251 ! published by the Free Software Foundation; either version 2 of
1252 ! the License, or (at your option) any later version.
1253 !
1254 ! This program is distributed in the hope that it will be useful,
1255 ! but WITHOUT ANY WARRANTY; without even the implied warranty of
1256 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
1257 ! GNU General Public License for more details.
1258 !
1259 ! A copy of the GNU General Public License is available at
1260 ! http://www.gnu.org/copyleft/gpl.html#SEC3
1261 ! or by writing to the Free Software Foundation, Inc.,
1262 ! 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
1263 !
1264 !
1265 ! 0. Authors
1266 !
1267 ! 1. UPDATE
1268 !
1269 ! SEP 1997: New for SWAN 32.01
1270 ! Cor van der Schelde - Delft Hydraulics
1271 ! 30.70, Feb. 98: test output suppressed (causes problem if subr
1272 ! is used during output
1273 !
1274 ! 2. PURPOSE
1275 !
1276 ! Transform degrees from nautical to cartesian or vice versa.
1277 !
1278 ! 3. METHOD
1279 !
1280 ! DEGCNV = 180 + dnorth - degree
1281 !
1282 ! 4. PARAMETERLIST
1283 !
1284 ! DEGCNV direction in cartesian or nautical degrees.
1285 ! DEGREE direction in nautical or cartesian degrees.
1286 !
1287 ! 5. SUBROUTINES CALLING
1288 !
1289 ! ---
1290 !
1291 ! 6. SUBROUTINES USED
1292 !
1293 ! NONE
1294 !
1295 ! 7. ERROR MESSAGES
1296 !
1297 ! NONE
1298 !
1299 ! 8. REMARKS
1300 !
1301 ! Nautical convention Cartesian convention
1302 !
1303 ! 0 90
1304 ! | |
1305 ! | |
1306 ! | |
1307 ! | |
1308 ! 270 --------+-------- 90 180 --------+-------- 0
1309 ! | |
1310 ! | |
1311 ! | |
1312 ! | |
1313 ! 180 270
1314 !
1315 ! 9. STRUCTURE
1316 !
1317 ! ---------------------------------
1318 ! IF (NAUTICAL DEGREES) THEN
1319 ! CONVERT DEGREES
1320 ! IF (DEGREES > 360 OR < 0) THEN
1321 ! CORRECT DEGREES WITHIN 0 - 360
1322 ! ---------------------------------
1323 !
1324 ! 10. SOURCE TEXT
1325 !
1326 !***********************************************************************
1327 !
1328  INTEGER :: IENT
1329  REAL :: DEGREE
1330 
1331  SAVE ient
1332  DATA ient /0/
1333  CALL strace(ient,'DEGCNV')
1334 !
1335  IF ( bnaut ) THEN
1336  degcnv = 180. + dnorth - degree
1337  ELSE
1338  degcnv = degree
1339  ENDIF
1340 !
1341  IF (degcnv >= 360.) THEN
1342  degcnv = mod(degcnv, 360.)
1343  ELSE IF (degcnv < 0.) THEN
1344  degcnv = mod(degcnv, 360.) + 360.
1345  ELSE
1346 ! DEGCNV between 0 and 360; do nothing
1347  ENDIF
1348 !
1349  RETURN
subroutine strace(IENT, SUBNAM)
Definition: ocpmix.f90:468
logical bnaut
Definition: swmod1.f90:2144
real dnorth
Definition: swmod1.f90:1721
real function degcnv(DEGREE)
Definition: swanser.f90:1225
Here is the call graph for this function:
Here is the caller graph for this function:

◆ gamma()

real function gamma ( real  XX)

Definition at line 788 of file swanser.f90.

788 ! *
789 !********************************************************************
790 !
791 ! Updates
792 ! ver 30.70, Oct 1997 by N.Booij: new subroutine
793 !
794 ! Purpose
795 ! Compute the transcendental function Gamma
796 !
797 ! Subroutines used
798 ! GAMMLN (Numerical Recipes)
799 !
800  REAL XX, YY, ABIG
801  SAVE ient, abig
802  DATA ient /0/, abig /30./
803  CALL strace (ient, 'GAMMA')
804  yy = gammln(xx)
805  IF (yy.GT.abig) yy = abig
806  IF (yy.LT.-abig) yy = -abig
807  gamma = exp(yy)
808  RETURN
subroutine strace(IENT, SUBNAM)
Definition: ocpmix.f90:468
real function gamma(XX)
Definition: swanser.f90:788
function gammln(XX)
Definition: swanser.f90:814
Here is the call graph for this function:
Here is the caller graph for this function:

◆ gammln()

function gammln (   XX)

Definition at line 814 of file swanser.f90.

814 ! *
815 !********************************************************************
816 !
817 ! Method:
818 ! function is copied from: Press et al., "Numerical Recipes"
819 !
820  DOUBLE PRECISION COF(6),STP,HALF,ONE,FPF,X,TMP,SER
821  DATA cof,stp/76.18009173d0,-86.50532033d0,24.01409822d0, &
822  -1.231739516d0,.120858003d-2,-.536382d-5,2.50662827465d0/
823  DATA half,one,fpf/0.5d0,1.0d0,5.5d0/
824  x=xx-one
825  tmp=x+fpf
826  tmp=(x+half)*log(tmp)-tmp
827  ser=one
828  DO 11 j=1,6
829  x=x+one
830  ser=ser+cof(j)/x
831 11 CONTINUE
832  gammln=tmp+log(stp*ser)
833  RETURN
function gammln(XX)
Definition: swanser.f90:814
Here is the caller graph for this function:

◆ intstr()

character*20 function intstr ( integer  IVAL)

Definition at line 696 of file swanser.f90.

696 !
697 !****************************************************************
698 !
699  IMPLICIT NONE
700 !
701 !
702 ! --|-----------------------------------------------------------|--
703 ! | Delft University of Technology |
704 ! | Faculty of Civil Engineering and Geosciences |
705 ! | Environmental Fluid Mechanics Section |
706 ! | P.O. Box 5048, 2600 GA Delft, The Netherlands |
707 ! | |
708 ! | Programmer: M. Zijlema |
709 ! --|-----------------------------------------------------------|--
710 !
711 !
712 ! SWAN (Simulating WAves Nearshore); a third generation wave model
713 ! Copyright (C) 2004-2005 Delft University of Technology
714 !
715 ! This program is free software; you can redistribute it and/or
716 ! modify it under the terms of the GNU General Public License as
717 ! published by the Free Software Foundation; either version 2 of
718 ! the License, or (at your option) any later version.
719 !
720 ! This program is distributed in the hope that it will be useful,
721 ! but WITHOUT ANY WARRANTY; without even the implied warranty of
722 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
723 ! GNU General Public License for more details.
724 !
725 ! A copy of the GNU General Public License is available at
726 ! http://www.gnu.org/copyleft/gpl.html#SEC3
727 ! or by writing to the Free Software Foundation, Inc.,
728 ! 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
729 !
730 !
731 ! 0. Authors
732 !
733 ! 40.23: Marcel Zijlema
734 !
735 ! 1. Updates
736 !
737 ! 40.23, Feb. 03: New subroutine
738 !
739 ! 2. Purpose
740 !
741 ! Convert integer to string
742 !
743 ! 4. Argument variables
744 !
745 ! IVAL integer to be converted
746 !
747  INTEGER IVAL
748 !
749 ! 6. Local variables
750 !
751 ! CVAL : character represented an integer of mantisse
752 ! I : counter
753 ! IPOS : position in mantisse
754 ! IQUO : whole quotient
755 !
756  INTEGER I, IPOS, IQUO
757  CHARACTER*1, ALLOCATABLE :: CVAL(:)
758 !
759 ! 12. Structure
760 !
761 ! Trivial.
762 !
763 ! 13. Source text
764 !
765  ipos = 1
766  100 CONTINUE
767  IF (ival/10**ipos.GE.1.) THEN
768  ipos = ipos + 1
769  GO TO 100
770  END IF
771  ALLOCATE(cval(ipos))
772 
773  DO i=ipos,1,-1
774  iquo=ival/10**(i-1)
775  cval(ipos-i+1)=char(int(iquo)+48)
776  ival=ival-iquo*10**(i-1)
777  END DO
778 
779  WRITE (intstr,*) (cval(i), i=1,ipos)
780 
781  RETURN
character *20 function intstr(IVAL)
Definition: swanser.f90:696

◆ kscip1()

subroutine kscip1 ( integer  MMT,
real, dimension(mmt)  SIG,
real  D,
real, dimension(mmt)  K,
real, dimension(mmt)  CG,
real, dimension(mmt)  N,
real, dimension(mmt)  ND 
)

Definition at line 504 of file swanser.f90.

504 ! *
505 !***********************************************************************
506 !
507  USE ocpcomm4
508  USE swcomm3
509 !
510  IMPLICIT NONE
511 !
512 ! --|-----------------------------------------------------------|--
513 ! | Delft University of Technology |
514 ! | Faculty of Civil Engineering |
515 ! | Environmental Fluid Mechanics Section |
516 ! | P.O. Box 5048, 2600 GA Delft, The Netherlands |
517 ! | |
518 ! | Programmers: R.C. Ris, N. Booij, |
519 ! | IJ.G. Haagsma, A.T.M.M. Kieftenburg, |
520 ! | M. Zijlema, E.E. Kriezi, |
521 ! | R. Padilla-Hernandez, L.H. Holthuijsen |
522 ! --|-----------------------------------------------------------|--
523 !
524 ! 0. Authors
525 !
526 ! 1. Updates
527 !
528 ! 2. Purpose
529 !
530 ! Calculation of the wave number, group velocity, group number N
531 ! and the derivative of N w.r.t. depth (=ND)
532 !
533 ! 3. Method
534 !
535 ! --
536 !
537 ! 4. Argument variables
538 !
539 ! MMT input number of frequency points
540 !
541  INTEGER MMT
542 !
543 ! CG output group velocity
544 ! D input local depth
545 ! K output wave number
546 ! N output ratio of group and phase velocity
547 ! ND output derivative of N with respect to D
548 ! SIG input rel. frequency for which wave parameters
549 ! must be determined
550 !
551  REAL CG(MMT), D, &
552  K(MMT), N(MMT), ND(MMT), SIG(MMT)
553 !
554 ! 6. Local variables
555 !
556 ! C phase velocity
557 ! FAC1 auxiliary factor
558 ! FAC2 auxiliary factor
559 ! FAC3 auxiliary factor
560 ! IENT number of entries
561 ! IS counter in frequency (sigma-space)
562 ! KND dimensionless wave number
563 ! ROOTDG square root of D/GRAV
564 ! WGD square root of GRAV*D
565 ! SND dimensionless frequency
566 ! SND2 = SND*SND
567 !
568  INTEGER IENT, IS
569  REAL KND, ROOTDG, SND, WGD, SND2, C, FAC1, FAC2, FAC3
570 !
571 ! 8. Subroutines used
572 !
573 ! --
574 !
575 ! 9. Subroutines calling
576 !
577 ! SWOEXA, SWOEXF (Swan/Output)
578 !
579 ! 10. Error messages
580 !
581 ! --
582 !
583 ! 11. Remarks
584 !
585 ! --
586 !
587 ! 12. Structure
588 !
589 ! -----------------------------------------------------------------
590 ! Compute non-dimensional frequency SND
591 ! IF SND >= 2.5, then
592 ! Compute wave number K, group velocity CGO, ratio of group
593 ! and phase velocity N and its derivative ND according to
594 ! deep water theory
595 ! ELSE IF SND =< 1.e-6
596 ! Compute wave number K, group velocity CGO, ratio of group
597 ! and phase velocity N and its derivative ND
598 ! according to extremely shallow water
599 ! ELSE
600 ! Compute wave number K, group velocity CGO and the ratio of
601 ! group and phase velocity N by Pade and other simple formulas.
602 ! Compute the derivative of N w.r.t. D = ND.
603 ! -----------------------------------------------------------------
604 !
605 ! 13. Source text
606 !
607 ! SAVE IENT
608 ! DATA IENT /0/
609 ! IF (LTRACE) CALL STRACE (IENT, 'KSCIP1')
610 !
611  rootdg = sqrt(d/grav_w)
612  wgd = rootdg*grav_w
613  DO is = 1, mmt
614 ! SND is dimensionless frequency
615  snd = sig(is) * rootdg
616  IF(snd >= 2.5)THEN
617 ! ******* deep water *******
618  k(is) = sig(is) * sig(is) / grav_w
619  cg(is) = 0.5 * grav_w / sig(is)
620  n(is) = 0.5
621  nd(is) = 0.
622  ELSE IF(snd < 1.e-6)THEN
623 ! *** very shallow water ***
624 ! print*,'IN VERY SHALLOW WATER'
625 ! stop
626  k(is) = snd/d
627  cg(is) = wgd
628  n(is) = 1.
629  nd(is) = 0.
630  ELSE
631  snd2 = snd*snd
632  c = sqrt(grav_w*d/(snd2+1./(1.+0.666*snd2+0.445*snd2**2 &
633  -0.105*snd2**3+0.272*snd2**4)))
634  k(is) = sig(is)/c
635  knd = k(is)*d
636  fac1 = 2.*knd/sinh(2.*knd)
637  n(is) = 0.5*(1.+fac1)
638  cg(is)= n(is)*c
639  fac2 = snd2/knd
640  fac3 = 2.*fac2/(1.+fac2*fac2)
641  nd(is)= fac1*(0.5/d - k(is)/fac3)
642  ENDIF
643  END DO
644 
645  RETURN
real grav_w
Definition: swmod1.f90:1721
integer n
Definition: mod_main.f90:55
Here is the caller graph for this function:

◆ kscip2()

subroutine kscip2 ( real, dimension(msc)  S,
real  D,
real, dimension(msc)  CG 
)

Definition at line 650 of file swanser.f90.

650 ! *
651 !***********************************************************************
652 !
653  USE ocpcomm4
654  USE swcomm3
655 !
656  IMPLICIT NONE
657  REAL CG(MSC), D, K,S(MSC)
658  INTEGER IENT, IS
659  REAL KND, ROOTDG, SND, WGD, SND2, C, FAC1,N
660  rootdg = sqrt(d/grav_w)
661  wgd = rootdg*grav_w
662 ! SND is dimensionless frequency
663  DO is=1,msc
664  snd = s(is) * rootdg
665  IF(snd >= 2.5)THEN
666 ! ******* deep water *******
667  k = s(is) * s(is) / grav_w
668  cg(is) = 0.5 * grav_w / s(is)
669  ELSE IF(snd < 1.e-6)THEN
670 ! *** very shallow water ***
671 ! print*,'IN VERY SHALLOW WATER'
672 ! stop
673  k = snd/d
674  cg(is) = wgd
675  ELSE
676  snd2 = snd*snd
677  c = sqrt(grav_w*d/(snd2+1./(1.+0.666*snd2+0.445*snd2**2 &
678  -0.105*snd2**3+0.272*snd2**4)))
679  k = s(is)/c
680  knd = k*d
681  fac1 = 2.*knd/sinh(2.*knd)
682  n = 0.5*(1.+fac1)
683  cg(is)= n*c
684  ENDIF
685  END DO
686  RETURN
real grav_w
Definition: swmod1.f90:1721
integer n
Definition: mod_main.f90:55
integer msc
Definition: swmod1.f90:1673
Here is the caller graph for this function:

◆ numstr()

character*20 function numstr ( integer  IVAL,
real  RVAL,
character  FORM 
)

Definition at line 305 of file swanser.f90.

305 !
306 !****************************************************************
307 !
308  USE ocpcomm4
309 !
310  IMPLICIT NONE
311 !
312 !
313 ! --|-----------------------------------------------------------|--
314 ! | Delft University of Technology |
315 ! | Faculty of Civil Engineering and Geosciences |
316 ! | Environmental Fluid Mechanics Section |
317 ! | P.O. Box 5048, 2600 GA Delft, The Netherlands |
318 ! | |
319 ! | Programmer: M. Zijlema |
320 ! --|-----------------------------------------------------------|--
321 !
322 !
323 ! SWAN (Simulating WAves Nearshore); a third generation wave model
324 ! Copyright (C) 2004-2005 Delft University of Technology
325 !
326 ! This program is free software; you can redistribute it and/or
327 ! modify it under the terms of the GNU General Public License as
328 ! published by the Free Software Foundation; either version 2 of
329 ! the License, or (at your option) any later version.
330 !
331 ! This program is distributed in the hope that it will be useful,
332 ! but WITHOUT ANY WARRANTY; without even the implied warranty of
333 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
334 ! GNU General Public License for more details.
335 !
336 ! A copy of the GNU General Public License is available at
337 ! http://www.gnu.org/copyleft/gpl.html#SEC3
338 ! or by writing to the Free Software Foundation, Inc.,
339 ! 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
340 !
341 !
342 ! 0. Authors
343 !
344 ! 40.23: Marcel Zijlema
345 ! 40.41: Marcel Zijlema
346 !
347 ! 1. Updates
348 !
349 ! 40.23, Feb. 03: New subroutine
350 ! 40.41, Oct. 04: common blocks replaced by modules, include files removed
351 !
352 ! 2. Purpose
353 !
354 ! Convert integer or real to string with given format
355 !
356 ! 4. Argument variables
357 !
358 ! IVAL integer to be converted
359 ! FORM given format
360 ! RVAL real to be converted
361 !
362  INTEGER IVAL
363  REAL RVAL
364  CHARACTER FORM*20
365 !
366 ! 6. Local variables
367 !
368 ! 12. Structure
369 !
370 ! Trivial.
371 !
372 ! 13. Source text
373 !
374  IF(ival /= inan)THEN
375  WRITE(numstr,form) ival
376  ELSE IF( rval /= rnan)THEN
377  WRITE (numstr,form) rval
378  ELSE
379  numstr = ''
380  END IF
381 
382  RETURN
real rnan
Definition: swmod1.f90:485
character *20 function numstr(IVAL, RVAL, FORM)
Definition: swanser.f90:305
integer inan
Definition: swmod1.f90:484
Here is the caller graph for this function:

◆ readxy()

subroutine readxy ( character, dimension(*)  NAMX,
character, dimension(*)  NAMY,
  XX,
  YY,
character, dimension(*)  KONT,
  XSTA,
  YSTA 
)

Definition at line 63 of file swanser.f90.

63 ! *
64 !***********************************************************************
65 !
66  USE ocpcomm1
67  USE swcomm2
68 !
69 !
70 ! --|-----------------------------------------------------------|--
71 ! | Delft University of Technology |
72 ! | Faculty of Civil Engineering |
73 ! | Environmental Fluid Mechanics Section |
74 ! | P.O. Box 5048, 2600 GA Delft, The Netherlands |
75 ! | |
76 ! | Programmers: R.C. Ris, N. Booij, |
77 ! | IJ.G. Haagsma, A.T.M.M. Kieftenburg, |
78 ! | M. Zijlema, E.E. Kriezi, |
79 ! | R. Padilla-Hernandez, L.H. Holthuijsen |
80 ! --|-----------------------------------------------------------|--
81 !
82 !
83 ! SWAN (Simulating WAves Nearshore); a third generation wave model
84 ! Copyright (C) 2004-2005 Delft University of Technology
85 !
86 ! This program is free software; you can redistribute it and/or
87 ! modify it under the terms of the GNU General Public License as
88 ! published by the Free Software Foundation; either version 2 of
89 ! the License, or (at your option) any later version.
90 !
91 ! This program is distributed in the hope that it will be useful,
92 ! but WITHOUT ANY WARRANTY; without even the implied warranty of
93 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
94 ! GNU General Public License for more details.
95 !
96 ! A copy of the GNU General Public License is available at
97 ! http://www.gnu.org/copyleft/gpl.html#SEC3
98 ! or by writing to the Free Software Foundation, Inc.,
99 ! 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
100 !
101 !
102 ! 0. Authors
103 ! 40.22: John Cazes and Tim Campbell
104 ! 40.13: Nico Booij
105 ! 40.51: Marcel Zijlema
106 !
107 ! 1. UPDATE
108 !
109 ! Nov. 1996 offset values are added to standard values
110 ! because they will be subtracted later
111 ! 40.13, Nov. 01: a valid value for YY is required if a valid value
112 ! for XX has been given; ocpcomm1.inc reactivated
113 ! 40.51, Feb. 05: correction to location points equal to offset values
114 !
115 ! 2. PURPOSE
116 !
117 ! Read x and y, initialize offset values XOFFS and YOFFS
118 !
119 ! 3. METHOD
120 !
121 ! ---
122 !
123 ! 4. PARAMETERLIST
124 !
125 ! NAMX, NAMY inp char names of the two coordinates as given in
126 ! the user manual
127 ! XX, YY out real values of x and y taking into account offset
128 ! KONT inp char what to be done if values are missing
129 ! see doc. of INDBLE (Ocean Pack doc.)
130 ! XSTA, YSTA inp real standard values of x and y
131 !
132 ! 5. SUBROUTINES CALLING
133 !
134 !
135 !
136 ! 6. SUBROUTINES USED
137 !
138 ! INDBLE (Ocean Pack)
139  LOGICAL EQREAL
140 !
141 ! 7. ERROR MESSAGES
142 !
143 ! ---
144 !
145 ! 8. REMARKS
146 !
147 ! ---
148 !
149 ! 9. STRUCTURE
150 !
151 ! ----------------------------------------------------------------
152 ! Read x and y in double prec.
153 ! If this is first couple of values
154 ! Then assign values to XOFFS and YOFFS
155 ! make LXOFFS True
156 ! ---------------------------------------------------------------
157 ! make XX and YY equal to x and y taking into account offset
158 ! ----------------------------------------------------------------
159 !
160 ! 10. SOURCE TEXT
161 !
162  DOUBLE PRECISION XTMP, YTMP
163  CHARACTER NAMX *(*), NAMY *(*), KONT *(*)
164  SAVE ient
165  DATA ient /0/
166  CALL strace (ient,'READXY')
167 !
168 !JQI CALL INDBLE (NAMX, XTMP, KONT, DBLE(XSTA)+DBLE(XOFFS))
169  IF(chgval)THEN
170 ! a valid value was given for XX
171 !JQI CALL INDBLE (NAMY, YTMP, 'REQ', DBLE(YSTA)+DBLE(YOFFS))
172  ELSE
173 !JQI CALL INDBLE (NAMY, YTMP, KONT, DBLE(YSTA)+DBLE(YOFFS))
174  ENDIF
175  IF(.NOT.lxoffs)THEN
176  xoffs = real(xtmp)
177  yoffs = real(ytmp)
178  lxoffs = .true.
179  ENDIF
180 !JQI IF(.NOT.EQREAL(XOFFS,REAL(XTMP)))THEN
181 !JQI XX = REAL(XTMP-DBLE(XOFFS))
182 !JQI ELSE IF(OPTG == 3)THEN
183 !JQI XX = 1.E-5
184 !JQI ELSE
185 !JQI XX = 0.
186 !JQI END IF
187 !JQI IF(.NOT.EQREAL(YOFFS,REAL(YTMP)))THEN
188 !JQI YY = REAL(YTMP-DBLE(YOFFS))
189 !JQI ELSE IF(OPTG == 3)THEN
190 !JQI YY = 1.E-5
191 !JQI ELSE
192 !JQI YY = 0.
193 !JQI END IF
194 !
195  RETURN
subroutine strace(IENT, SUBNAM)
Definition: ocpmix.f90:468
real yoffs
Definition: swmod1.f90:1364
logical lxoffs
Definition: swmod1.f90:1365
logical chgval
Definition: swmod1.f90:154
real xoffs
Definition: swmod1.f90:1363
Here is the call graph for this function:

◆ sshape()

subroutine sshape ( real, dimension(mdc,msc)  ACLOC,
real, dimension(msc)  SPCSIG,
real, dimension(mdc,6)  SPCDIR,
integer  FSHAPL,
integer  DSHAPL 
)

Definition at line 840 of file swanser.f90.

840 ! *
841 !***********************************************************************
842 !
843  USE ocpcomm4
844  USE swcomm3
845 !
846 ! 0. Authors
847 !
848 ! 1. Updates
849 !
850 ! 2. Purpose
851 !
852 ! Calculating of energy density at boundary point (x,y,sigma,theta)
853 !
854 ! 3. Method (updated...)
855 !
856 ! see: M. Yamaguchi: Approximate expressions for integral properties
857 ! of the JONSWAP spectrum; Proc. JSCE, No. 345/II-1, pp. 149-152,
858 ! 1984.
859 !
860 ! computation of mean period: see Swan system documentation
861 !
862 ! 4. Argument variables
863 !
864 ! o ACLOC : Energy density at a point in space
865 ! i SPCDIR: (*,1); spectral directions (radians) 30.82
866 ! (*,2); cosine of spectral directions 30.82
867 ! (*,3); sine of spectral directions 30.82
868 ! (*,4); cosine^2 of spectral directions 30.82
869 ! (*,5); cosine*sine of spectral directions 30.82
870 ! (*,6); sine^2 of spectral directions 30.82
871 ! i SPCSIG: Relative frequencies in computational domain in sigma-space 30.82
872 !
873  REAL ACLOC(MDC,MSC)
874  REAL SPCDIR(MDC,6)
875  REAL SPCSIG(MSC)
876 !
877 ! i DSHAPL: Directional distribution
878 ! i FSHAPL: Shape of spectrum:
879 ! =1; Pierson-Moskowitz spectrum
880 ! =2; Jonswap spectrum
881 ! =3; bin
882 ! =4; Gauss curve
883 ! (if >0: period is interpreted as peak per.
884 ! if <0: period is interpreted as mean per.)
885 !
886  INTEGER FSHAPL, DSHAPL
887 !
888 ! 5. Parameter variables
889 !
890 ! 6. Local variables
891 !
892 ! ID counter of directions
893 ! IS counter of frequencies
894 ! LSHAPE absolute value of FSHAPL
895 !
896  INTEGER ID, IS, LSHAPE
897 !
898 ! PKPER peak period 30.80
899 ! APSHAP aux. var. used in computation of spectrum
900 ! AUX1 auxiliary variable
901 ! AUX2 auxiliary variable
902 ! AUX3 auxiliary variable
903 ! COEFF coefficient for behaviour around the peak (Jonswap)
904 ! CPSHAP aux. var. used in computation of spectrum
905 ! CTOT total energy
906 ! CTOTT total energy (used for comparison)
907 ! DIFPER auxiliary variable used to select bin closest
908 ! to given frequency
909 ! MPER
910 ! MS power in directional distribution
911 ! RA action density
912 ! SALPHA
913 ! SF frequency (Hz)
914 ! SF4 SF**4
915 ! SF5 SF**5
916 ! FPK frequency corresponding to peak period (1/PKPER) 30.80
917 ! FPK4 FPK**4
918 ! SYF peakedness parameter
919 !
920  REAL APSHAP, AUX1, AUX2, AUX3
921  REAL COEFF ,SYF ,MPER ,CTOT ,CTOTT,PKPER ,DIFPER
922  REAL MS
923  REAL RA ,SALPHA,SF ,SF4 ,SF5 ,FPK ,FPK4, FAC
924 !
925 ! LOGPM indicates whether peak or mean frequency is used
926 ! DVERIF logical used in verification of incident direction
927 !
928  LOGICAL LOGPM, DVERIF
929 !
930 ! PSHAPE coefficients of spectral distribution (see remarks)
931 ! SPPARM array containing integral wave parameters (see remarks)
932 !
933 ! 8. Subroutines used
934 !
935 ! ---
936 !
937 ! 9. Subroutines calling
938 !
939 ! 10. Error messages
940 !
941 ! 11. Remarks
942 !
943 ! PSHAPE(1): SY0, peak enhancement factor (gamma) in Jonswap spectrum
944 ! PSHAPE(2): spectral width in case of Gauss spectrum in rad/s
945 !
946 ! SPPARM real input incident wave parameters (Hs, Period,
947 ! direction, Ms (dir. spread))
948 ! SPPARM(1): Hs, sign. wave height
949 ! SPPARM(2): Wave period given by the user (either peak or mean) 30.80
950 ! SPPARM(3): average direction
951 ! SPPARM(4): directional spread
952 !
953 ! ---------------------------------------------------------------------
954 !
955 ! In the case of a JONSWAP spectrum the initial conditions are given by
956 ! _ _ _ _ _
957 ! | _ _ -4 | | | S - S |
958 ! 2 | | S | | | | p |
959 ! a g | | _ | | exp|-1/2 * |________ |* 2/pi COS(T-T )
960 ! E(S,D )= ___ exp|-5/4 *| S | | G | | e * S | wi
961 ! wa 5 | | p | | |_ |_ p _|
962 ! S | |_ _| |
963 ! |_ _|
964 !
965 ! where
966 ! S : rel. frequency
967 !
968 ! D : Dir. of wave component
969 ! wa
970 !
971 ! a : equili. range const. (Phillips' constant)
972 ! g : gravity acceleration
973 !
974 ! S : Peak frequency
975 ! p
976 !
977 ! G : Peak enhancement factor
978 ! e : Peak width
979 !
980 ! T : local wind direction
981 ! wi
982 !
983 ! 12. Structure
984 !
985 ! ----------------------------------------------------------------
986 ! case shape
987 ! =1: calculate value of Pierson-Moskowitz spectrum
988 ! =2: calculate value of Jonswap spectrum
989 ! =3: calculate value of bin spectrum
990 ! =4: calculate value of Gauss spectrum
991 ! else: Give error message because of wrong shape
992 ! ----------------------------------------------------------------
993 ! if LOGPM is True
994 ! then calculate average period
995 ! if it differs from given average period
996 ! then recalculate peak period
997 ! restart procedure to compute spectral shape
998 ! ----------------------------------------------------------------
999 ! for all spectral bins do
1000 ! multiply all action densities by directional distribution
1001 ! ----------------------------------------------------------------
1002 !
1003 ! 13. Source text
1004 !
1005  SAVE ient
1006  DATA ient/0/
1007  CALL strace(ient,'SSHAPE')
1008 !
1009  IF (itest >= 80) WRITE (prtest, 8) fshapl, dshapl, &
1010  (spparm(jj), jj = 1,4)
1011  8 FORMAT (' entry SSHAPE ', 2i3, 4e12.4)
1012  IF(fshapl < 0)THEN
1013  lshape = - fshapl
1014  logpm = .false.
1015  ELSE
1016  lshape = fshapl
1017  logpm = .true.
1018  ENDIF
1019 
1020 !
1021  IF(spparm(1) <= 0.) &
1022  CALL msgerr(1,' sign. wave height at boundary is not positive')
1023 !
1024  pkper = spparm(2)
1025  itper = 0
1026  IF(lshape == 3)THEN
1027 ! select bin closest to given period
1028  difper = 1.e10
1029  DO is = 1, msc
1030  IF(abs(pkper - pi2_w/spcsig(is)) < difper)THEN
1031  isp = is
1032  difper = abs(pkper - pi2_w/spcsig(is))
1033  ENDIF
1034  ENDDO
1035  ENDIF
1036 !
1037 ! compute spectral shape using peak period PKPER
1038 !
1039  fac = 1.
1040  100 fpk = (1./pkper)
1041  fpk4 = fpk**4
1042  IF(lshape == 1)THEN
1043  salpha = ((spparm(1) ** 2) * (fpk4)) * 5. / 16.
1044  ELSE IF(lshape == 2)THEN
1045 ! *** SALPHA = alpha*(grav**2)/(2.*pi)**4)
1046  salpha = (spparm(1)**2 * fpk4) / &
1047  ((0.06533*(pshape(1)**0.8015)+0.13467)*16.)
1048  ELSE IF(lshape == 4)THEN
1049  aux1 = spparm(1)**2 / ( 16.* sqrt(pi2_w) * pshape(2))
1050  aux3 = 2. * pshape(2)**2
1051  ENDIF
1052 !
1053  ctott = 0.
1054  DO is = 1, msc
1055 !
1056  IF(lshape == 1)THEN
1057 ! *** LSHAPE = 1 : Pierson and Moskowitz ***
1058  sf = spcsig(is) / pi2_w
1059  sf4 = sf**4
1060  sf5 = sf**5
1061  ra = (salpha/sf5)*exp(-(5.*fpk4)/(4.*sf4))/(pi2_w*spcsig(is))
1062  acloc(mdc,is) = ra
1063  ELSE IF(lshape == 2)THEN
1064 ! *** LSHAPE = 2 : JONSWAP ***
1065  sf = spcsig(is)/(pi2_w)
1066  sf4 = sf**4
1067  sf5 = sf**5
1068  cpshap = 1.25 * fpk4 / sf4
1069  IF(cpshap > 10.)THEN
1070  ra = 0.
1071  ELSE
1072  ra = (salpha/sf5) * exp(-cpshap)
1073  ENDIF
1074  IF(sf < fpk)THEN
1075  coeff = 0.07
1076  ELSE
1077  coeff = 0.09
1078  ENDIF
1079  apshap = 0.5 * ((sf-fpk) / (coeff*fpk)) **2
1080  IF(apshap > 10.)THEN
1081  syf = 1.
1082  ELSE
1083  ppshap = exp(-apshap)
1084  syf = pshape(1)**ppshap
1085  ENDIF
1086  ra = syf*ra/(spcsig(is)*pi2_w)
1087  acloc(mdc,is) = ra
1088  IF(itest >= 120) WRITE (prtest, 112) &
1089  sf, salpha, cpshap, apshap, syf, ra
1090  112 FORMAT (' SSHAPE freq. ', 8e12.4)
1091  ELSE IF(lshape == 3)THEN
1092 !
1093 ! *** all energy concentrated in one BIN ***
1094 !
1095  IF(is == isp)THEN
1096  acloc(mdc,is) = ( spparm(1)**2 ) / &
1097  ( 16. * spcsig(is)**2 * frintf )
1098  ELSE
1099  acloc(mdc,is) = 0.
1100  ENDIF
1101  ELSE IF(lshape == 4)THEN
1102 !
1103 ! *** energy Gaussian distributed (wave-current tests) ***
1104 !
1105  aux2 = ( spcsig(is) - ( pi2_w / pkper ) )**2
1106  ra = aux1 * exp( -1. * aux2 / aux3 ) / spcsig(is)
1107  acloc(mdc,is) = ra
1108  ELSE
1109  IF(is == 1)THEN
1110  CALL msgerr (2,'Wrong type for frequency shape')
1111  WRITE (printf, *) ' -> ', fshapl, lshape
1112  ENDIF
1113  ENDIF
1114  IF (itest >= 10) &
1115  ctott = ctott + frintf * acloc(mdc,is) * spcsig(is)**2
1116  END DO !END DO IS
1117  IF(itest >= 10)THEN
1118  IF(spparm(1) > 0.01)THEN
1119  hstmp = 4. * sqrt(ctott)
1120  IF(abs(hstmp-spparm(1)) > 0.1*spparm(1)) &
1121  WRITE (printf, 303) spparm(1), hstmp
1122  303 FORMAT (' SSHAPE, deviation in Hs, should be ', f8.3, &
1123  ', calculated ', f8.3)
1124  ENDIF
1125  ENDIF
1126 !
1127 ! if mean frequency was given recalculate PKPER and restart
1128 !
1129  IF (.NOT.logpm .AND. itper < 10)THEN
1130  itper = itper + 1
1131 ! calculate average frequency
1132  am0 = 0.
1133  am1 = 0.
1134  DO is = 1, msc
1135  as2 = acloc(mdc,is) * (spcsig(is))**2
1136  as3 = as2 * spcsig(is)
1137  am0 = am0 + as2
1138  am1 = am1 + as3
1139  ENDDO
1140 ! contribution of tail to total energy density
1141  pptail = pwtail(1) - 1.
1142  aptail = 1. / (pptail * (1. + pptail * (frinth-1.)))
1143  am0 = am0 * frintf + aptail * as2
1144  pptail = pwtail(1) - 2.
1145  eptail = 1. / (pptail * (1. + pptail * (frinth-1.)))
1146  am1 = am1 * frintf + eptail * as3
1147 ! Mean period:
1148  IF( am1 /= 0.)THEN
1149  mper = pi2_w * am0 / am1
1150  ELSE
1151  CALL msgerr(3, ' first moment is zero in calculating the')
1152  CALL msgerr(3, ' spectrum at boundary using param. bc.')
1153  END IF
1154  IF(itest >= 80) WRITE (prtest, 72) itper, spparm(2), mper, &
1155  pkper
1156  72 FORMAT (' SSHAPE iter=', i2, ' period values:', 3f7.2)
1157  IF(abs(mper-spparm(2)) > 0.01*spparm(2))THEN
1158 ! modification suggested by Mauro Sclavo
1159  pkper = (spparm(2) / mper) * pkper
1160  GOTO 100
1161  ENDIF
1162  ENDIF
1163 !
1164  IF (itper >= 10)THEN
1165  CALL msgerr(3, 'No convergence calculating the spectrum')
1166  CALL msgerr(3, 'at the boundary using parametric bound. cond.')
1167  ENDIF
1168 !
1169 ! now introduce distribution over directions
1170 !
1171  adir = pi_w * degcnv(spparm(3)) / 180.
1172  IF(dshapl == 1)THEN
1173  dspr = pi_w * spparm(4) / 180.
1174  ms = max(dspr**(-2) - 2., 1.)
1175  ELSE
1176  ms = spparm(4)
1177  ENDIF
1178  IF(ms < 12.)THEN
1179  ctot = (2.**ms) * (gamma(0.5*ms+1.))**2 / (pi_w * gamma(ms+1.))
1180  ELSE
1181  ctot = sqrt(0.5*ms/pi_w) / (1. - 0.25/ms)
1182  ENDIF
1183  IF(itest >= 100)THEN
1184  esom = 0.
1185  DO is = 1, msc
1186  esom = esom + frintf * spcsig(is)**2 * acloc(mdc,is)
1187  ENDDO
1188  WRITE (prtest, *) ' SSHAPE dir ', 4.*sqrt(abs(esom)), &
1189  spparm(1), ctot, ms, gamma(0.5*ms+1.), gamma(ms+1.), &
1190  ctot
1191  ENDIF
1192  dverif = .false.
1193  ctott = 0.
1194  DO id = 1, mdc
1195  acos = cos(spcdir(id,1) - adir)
1196  IF(acos > 0.)THEN
1197  cdir = ctot * max(acos**ms, 1.e-10)
1198  IF(.NOT.fulcir)THEN
1199  IF(acos >= cos(ddir)) dverif = .true.
1200  ENDIF
1201  ELSE
1202  cdir = 0.
1203  ENDIF
1204  IF(itest >= 10) ctott = ctott + cdir * ddir
1205  IF(itest >= 100) WRITE (prtest, 360) id,spcdir(id,1),cdir
1206  360 FORMAT (' ID Spcdir Cdir: ',i3,3(1x,e11.4))
1207  DO is = 1, msc
1208  acloc(id,is) = cdir * acloc(mdc,is)
1209  ENDDO
1210  ENDDO
1211  IF(itest >= 10)THEN
1212  IF(abs(ctott-1.) > 0.1) WRITE (printf, 363) ctott
1213  363 FORMAT (' SSHAPE, integral of Cdir is not 1, but:', f6.3)
1214  ENDIF
1215  IF (.NOT.fulcir .AND. .NOT.dverif) &
1216  CALL msgerr (1, 'incident direction is outside sector')
1217 !
1218  RETURN
1219 
real, dimension(:), allocatable, save spcsig
Definition: swmod2.f90:767
subroutine strace(IENT, SUBNAM)
Definition: ocpmix.f90:468
real, dimension(:,:), allocatable, save spcdir
Definition: swmod2.f90:767
integer prtest
Definition: swmod1.f90:517
real, dimension(mshape) pshape
Definition: swmod1.f90:2134
real, dimension(msppar) spparm
Definition: swmod1.f90:2138
integer mdc
Definition: swmod1.f90:1672
real pi2_w
Definition: swmod1.f90:1722
logical fulcir
Definition: swmod1.f90:1683
real pi_w
Definition: swmod1.f90:1721
integer printf
Definition: swmod1.f90:517
real frintf
Definition: swmod1.f90:1678
integer itest
Definition: swmod1.f90:536
real function gamma(XX)
Definition: swanser.f90:788
real function degcnv(DEGREE)
Definition: swanser.f90:1225
real frinth
Definition: swmod1.f90:1678
subroutine msgerr(LEV, STRING)
Definition: ocpmix.f90:582
real, dimension(10) pwtail
Definition: swmod1.f90:1722
integer msc
Definition: swmod1.f90:1673
real ddir
Definition: swmod1.f90:1676
Here is the call graph for this function:

◆ swanout()

subroutine swanout (   DEP2)

Definition at line 1356 of file swanser.f90.

1356  RETURN

◆ swcopi()

subroutine swcopi ( integer, dimension(length)  IARR1,
integer, dimension(length)  IARR2,
integer  LENGTH 
)

Definition at line 387 of file swanser.f90.

387 !
388 !****************************************************************
389 !
390  USE ocpcomm4
391 !
392  IMPLICIT NONE
393 !
394 !
395 ! --|-----------------------------------------------------------|--
396 ! | Delft University of Technology |
397 ! | Faculty of Civil Engineering and Geosciences |
398 ! | Environmental Fluid Mechanics Section |
399 ! | P.O. Box 5048, 2600 GA Delft, The Netherlands |
400 ! | |
401 ! | Programmer: M. Zijlema |
402 ! --|-----------------------------------------------------------|--
403 !
404 !
405 ! SWAN (Simulating WAves Nearshore); a third generation wave model
406 ! Copyright (C) 2004-2005 Delft University of Technology
407 !
408 ! This program is free software; you can redistribute it and/or
409 ! modify it under the terms of the GNU General Public License as
410 ! published by the Free Software Foundation; either version 2 of
411 ! the License, or (at your option) any later version.
412 !
413 ! This program is distributed in the hope that it will be useful,
414 ! but WITHOUT ANY WARRANTY; without even the implied warranty of
415 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
416 ! GNU General Public License for more details.
417 !
418 ! A copy of the GNU General Public License is available at
419 ! http://www.gnu.org/copyleft/gpl.html#SEC3
420 ! or by writing to the Free Software Foundation, Inc.,
421 ! 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
422 !
423 !
424 ! 0. Authors
425 !
426 ! 40.23: Marcel Zijlema
427 ! 40.41: Marcel Zijlema
428 !
429 ! 1. Updates
430 !
431 ! 40.23, Feb. 03: New subroutine
432 ! 40.41, Oct. 04: common blocks replaced by modules, include files removed
433 !
434 ! 2. Purpose
435 !
436 ! Copies integer array IARR1 to IARR2
437 !
438 ! 3. Method
439 !
440 ! ---
441 !
442 ! 4. Argument variables
443 !
444 ! IARR1 source array
445 ! IARR2 target array
446 ! LENGTH array length
447 !
448  INTEGER LENGTH
449  INTEGER IARR1(LENGTH), IARR2(LENGTH)
450 !
451 ! 6. Local variables
452 !
453 ! I : loop counter
454 ! IENT : number of entries
455 !
456  INTEGER I, IENT
457 !
458 ! 8. Subroutines used
459 !
460 ! MSGERR Writes error message
461 ! STRACE Tracing routine for debugging
462 !
463 ! 9. Subroutines calling
464 !
465 ! ---
466 !
467 ! 10. Error messages
468 !
469 ! ---
470 !
471 ! 11. Remarks
472 !
473 ! ---
474 !
475 ! 12. Structure
476 !
477 ! Trivial.
478 !
479 ! 13. Source text
480 !
481  SAVE ient
482  DATA ient/0/
483  IF(ltrace) CALL strace (ient,'SWCOPI')
484 
485 ! --- check array length
486 
487  IF(length <= 0)THEN
488  CALL msgerr( 3, 'Array length should be positive' )
489  END IF
490 
491 ! --- copy elements of array IARR1 to IARR2
492 
493  DO i = 1, length
494  iarr2(i) = iarr1(i)
495  END DO
496 
497  RETURN
subroutine strace(IENT, SUBNAM)
Definition: ocpmix.f90:468
logical ltrace
Definition: swmod1.f90:537
subroutine msgerr(LEV, STRING)
Definition: ocpmix.f90:582
Here is the call graph for this function:

◆ txpbla()

subroutine txpbla ( character*(*)  TEXT,
integer  IF,
integer  IL 
)

Definition at line 202 of file swanser.f90.

202 !
203 !****************************************************************
204 !
205  IMPLICIT NONE
206 !
207 !
208 ! --|-----------------------------------------------------------|--
209 ! | Delft University of Technology |
210 ! | Faculty of Civil Engineering and Geosciences |
211 ! | Environmental Fluid Mechanics Section |
212 ! | P.O. Box 5048, 2600 GA Delft, The Netherlands |
213 ! | |
214 ! | Programmer: M. Zijlema |
215 ! --|-----------------------------------------------------------|--
216 !
217 !
218 ! SWAN (Simulating WAves Nearshore); a third generation wave model
219 ! Copyright (C) 2004-2005 Delft University of Technology
220 !
221 ! This program is free software; you can redistribute it and/or
222 ! modify it under the terms of the GNU General Public License as
223 ! published by the Free Software Foundation; either version 2 of
224 ! the License, or (at your option) any later version.
225 !
226 ! This program is distributed in the hope that it will be useful,
227 ! but WITHOUT ANY WARRANTY; without even the implied warranty of
228 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
229 ! GNU General Public License for more details.
230 !
231 ! A copy of the GNU General Public License is available at
232 ! http://www.gnu.org/copyleft/gpl.html#SEC3
233 ! or by writing to the Free Software Foundation, Inc.,
234 ! 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
235 !
236 !
237 ! 0. Authors
238 !
239 ! 40.23: Marcel Zijlema
240 !
241 ! 1. Updates
242 !
243 ! 40.23, Feb. 03: New subroutine
244 !
245 ! 2. Purpose
246 !
247 ! determines the position of the first and the last non-blank
248 ! (or non-tabulator) character in the text-string
249 !
250 ! 4. Argument variables
251 !
252 ! IF position of the first non-blank character in TEXT
253 ! IL position of the last non-blank character in TEXT
254 ! TEXT text string
255 !
256  INTEGER IF, IL
257  CHARACTER*(*) TEXT
258 !
259 ! 6. Local variables
260 !
261 ! FOUND : TEXT is found or not
262 ! ITABVL: integer value of tabulator character
263 ! LENTXT: length of TEXT
264 !
265  INTEGER LENTXT, ITABVL
266  LOGICAL FOUND
267 !
268 ! 12. Structure
269 !
270 ! Trivial.
271 !
272 ! 13. Source text
273 !
274 !DOS ITABVL = ICHAR(' ')
275 !UNIX ITABVL = ICHAR(' ')
276  lentxt = len(text)
277  IF = 1
278  found = .false.
279 100 IF(IF <= lentxt .AND. .NOT. found)THEN
280  IF(.NOT. (text(if:if) == ' ' .OR. &
281  ichar(text(if:if)) == itabvl))THEN
282  found = .true.
283  ELSE
284  IF = IF + 1
285  ENDIF
286  GOTO 100
287  ENDIF
288  il = lentxt + 1
289  found = .false.
290 200 IF(il > 1 .AND. .NOT. found)THEN
291  il = il - 1
292  IF(.NOT. (text(il:il) == ' ' .OR. &
293  ichar(text(il:il)) == itabvl))THEN
294  found = .true.
295  ENDIF
296  GOTO 200
297  ENDIF
298 
299  RETURN
Here is the caller graph for this function: