Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
ElmerCSC
GitHub Repository: ElmerCSC/elmerfem
Path: blob/devel/fem/src/GeneralUtils.F90
3203 views
1
!/*****************************************************************************/
2
! *
3
! * Elmer, A Finite Element Software for Multiphysical Problems
4
! *
5
! * Copyright 1st April 1995 - , CSC - IT Center for Science Ltd., Finland
6
! *
7
! * This library is free software; you can redistribute it and/or
8
! * modify it under the terms of the GNU Lesser General Public
9
! * License as published by the Free Software Foundation; either
10
! * version 2.1 of the License, or (at your option) any later version.
11
! *
12
! * This library is distributed in the hope that it will be useful,
13
! * but WITHOUT ANY WARRANTY; without even the implied warranty of
14
! * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
15
! * Lesser General Public License for more details.
16
! *
17
! * You should have received a copy of the GNU Lesser General Public
18
! * License along with this library (in file ../LGPL-2.1); if not, write
19
! * to the Free Software Foundation, Inc., 51 Franklin Street,
20
! * Fifth Floor, Boston, MA 02110-1301 USA
21
! *
22
! *****************************************************************************/
23
!
24
!/******************************************************************************
25
! *
26
! * Authors: Juha Ruokolainen
27
! * Email: [email protected]
28
! * Web: http://www.csc.fi/elmer
29
! * Address: CSC - IT Center for Science Ltd.
30
! * Keilaranta 14
31
! * 02101 Espoo, Finland
32
! *
33
! * Original Date: 01 Oct 1996
34
! *
35
! *****************************************************************************/
36
37
#include "../config.h"
38
39
!> \ingroup ElmerLib
40
!> \}
41
42
!-----------------------------------------------------------------------------
43
!> Miscellaneous utilities.
44
!-----------------------------------------------------------------------------
45
MODULE GeneralUtils
46
47
USE LoadMod
48
49
#ifdef HAVE_LUA
50
USE, INTRINSIC :: ISO_C_BINDING
51
#endif
52
53
IMPLICIT NONE
54
55
INTERFACE AllocateVector
56
MODULE PROCEDURE AllocateRealVector, AllocateIntegerVector, &
57
AllocateComplexVector, AllocateLogicalVector, &
58
AllocateElementVector
59
END INTERFACE
60
61
INTERFACE AllocateArray
62
MODULE PROCEDURE AllocateRealArray, AllocateIntegerArray, &
63
AllocateComplexArray, AllocateLogicalArray
64
END INTERFACE
65
66
INTERFACE ComponentName
67
MODULE PROCEDURE ComponentNameStr, ComponentNameVar
68
END INTERFACE
69
70
REAL(KIND=dp), PRIVATE :: AdvanceTime1, AdvanceTime2
71
72
CONTAINS
73
74
!------------------------------------------------------------------------------
75
SUBROUTINE StartAdvanceOutput( SolverName, OutputType )
76
!------------------------------------------------------------------------------
77
CHARACTER(LEN=*) :: SolverName, OutputType
78
!------------------------------------------------------------------------------
79
AdvanceTime1 = RealTime()
80
AdvanceTime2 = RealTime()
81
CALL Info( SolverName, OutputType, Level=5 )
82
!------------------------------------------------------------------------------
83
END SUBROUTINE StartAdvanceOutput
84
!------------------------------------------------------------------------------
85
86
87
!------------------------------------------------------------------------------
88
SUBROUTINE AdvanceOutput(t,n,dot_t,percent_t)
89
!------------------------------------------------------------------------------
90
IMPLICIT NONE
91
INTEGER :: t,n
92
REAL(KIND=dp), OPTIONAL :: dot_t,percent_t
93
!------------------------------------------------------------------------------
94
INTEGER :: i
95
REAL(KIND=dp) :: d_t, p_t
96
!------------------------------------------------------------------------------
97
d_t = 1._dp
98
p_t = 20._dp
99
IF ( PRESENT(dot_t) ) d_t = dot_t
100
IF ( PRESENT(percent_t) ) p_t = percent_t
101
102
IF ( RealTime() - AdvanceTime1 > d_t ) THEN
103
CALL Info( '', '.', Level=5, noAdvance=.TRUE. )
104
105
IF ( RealTime() - AdvanceTime2 > p_t ) THEN
106
i = NINT(t*100.0/n)
107
WRITE(Message, '(i3,a)' ) i, '%'
108
CALL Info( '', Message, Level=5 )
109
AdvanceTime2 = RealTime()
110
END IF
111
AdvanceTime1 = RealTime()
112
END IF
113
!------------------------------------------------------------------------------
114
END SUBROUTINE AdvanceOutput
115
!------------------------------------------------------------------------------
116
117
118
!------------------------------------------------------------------------------
119
PURE FUNCTION lentrim(str) RESULT(n)
120
!------------------------------------------------------------------------------
121
CHARACTER(LEN=*), INTENT(IN) :: str
122
INTEGER :: n
123
DO n=LEN(str),1,-1
124
IF ( str(n:n) /= ' ' ) EXIT
125
END DO
126
!------------------------------------------------------------------------------
127
END FUNCTION lentrim
128
!------------------------------------------------------------------------------
129
130
131
!------------------------------------------------------------------------------
132
!> Compare equality of start of s1 to (in most uses string literal) s2.
133
!------------------------------------------------------------------------------
134
!------------------------------------------------------------------------------
135
PURE FUNCTION SEQL(s1,s2) RESULT(L)
136
!------------------------------------------------------------------------------
137
LOGICAL :: L
138
CHARACTER(LEN=*), INTENT(IN) :: s1,s2
139
!------------------------------------------------------------------------------
140
INTEGER :: n
141
!------------------------------------------------------------------------------
142
L = .FALSE.
143
n = LEN(s2)
144
IF(LEN(s1) < n) RETURN
145
IF (s1(1:n)==s2) L=.TRUE.
146
!------------------------------------------------------------------------------
147
END FUNCTION SEQL
148
!------------------------------------------------------------------------------
149
150
151
!------------------------------------------------------------------------------
152
!> Converts integer to string. Handy when writing output with integer data.
153
!------------------------------------------------------------------------------
154
PURE FUNCTION i2s(ival) RESULT(s)
155
!------------------------------------------------------------------------------
156
INTEGER, INTENT(in) :: ival
157
CHARACTER(:), ALLOCATABLE :: s
158
!------------------------------------------------------------------------------
159
INTEGER :: i,j,n,t,v,len
160
INTEGER(8) :: m
161
CHARACTER, PARAMETER :: DIGITS(0:9)=['0','1','2','3','4','5','6','7','8','9']
162
!------------------------------------------------------------------------------
163
164
IF(ival>=0) THEN
165
v = ival
166
IF (v<10) THEN
167
s = DIGITS(v)
168
ELSE IF (ival<100) THEN
169
i = v/10
170
s = DIGITS(i)//DIGITS(v-10*i)
171
ELSE
172
n=3
173
m=100
174
DO WHILE(10*m<=v)
175
n=n+1
176
m=m*10
177
END DO
178
179
ALLOCATE(CHARACTER(n)::s)
180
DO i=1,n
181
t = v / m
182
s(i:i) = DIGITS(t)
183
v = v - t*m
184
m = m / 10
185
END DO
186
END IF
187
ELSE
188
v = -ival
189
IF (v<10) THEN
190
s = '-'//DIGITS(v)
191
ELSE IF (v<100) THEN
192
i = v/10
193
s = '-'//DIGITS(i)//DIGITS(v-10*i)
194
ELSE
195
n=3
196
m=100
197
DO WHILE(10*m<=v)
198
n=n+1
199
m=m*10
200
END DO
201
202
ALLOCATE(CHARACTER(n+1)::s)
203
s(1:1) = '-'
204
DO i=2,n+1
205
t = v / m
206
s(i:i) = DIGITS(t)
207
v = v - t*m
208
m = m / 10
209
END DO
210
END IF
211
END IF
212
!------------------------------------------------------------------------------
213
END FUNCTION i2s
214
!------------------------------------------------------------------------------
215
216
217
218
!------------------------------------------------------------------------------
219
!> Converts string of length n to an integer number.
220
!------------------------------------------------------------------------------
221
PURE FUNCTION s2i(str,n) RESULT(ival)
222
!------------------------------------------------------------------------------
223
INTEGER, INTENT(IN) :: n
224
INTEGER :: ival
225
CHARACTER(LEN=n), INTENT(IN) :: str
226
!------------------------------------------------------------------------------
227
LOGICAL :: neg
228
INTEGER :: j,k
229
INTEGER, PARAMETER :: ic0 = ICHAR('0')
230
231
neg = str(1:1)=='-'
232
k=1
233
IF ( neg ) k=2
234
235
ival = 0
236
DO j=k,n
237
ival = 10*ival + ICHAR(str(j:j)) - ic0
238
END DO
239
IF(neg) ival=-ival
240
END FUNCTION s2i
241
!------------------------------------------------------------------------------
242
243
244
!------------------------------------------------------------------------------
245
!> Converts a string into a number of integer numbers
246
!> It is assumed that the integers may also be separated by
247
!> the given separator.
248
!------------------------------------------------------------------------------
249
FUNCTION str2ints(str,ints,sep) RESULT(n)
250
!------------------------------------------------------------------------------
251
INTEGER, INTENT(out) :: ints(:)
252
CHARACTER(LEN=*), INTENT(in) :: str
253
CHARACTER, OPTIONAL, INTENT(in) :: sep
254
255
INTEGER :: i,k,l,m,n,ic, icsep
256
INTEGER, PARAMETER :: ic0 = ICHAR('0'), ic9 = ICHAR('9'), icm = ICHAR('-'), &
257
ics = ICHAR(' ')
258
259
IF( PRESENT( sep ) ) THEN
260
icsep = ICHAR(sep)
261
ELSE
262
icsep = ics
263
END IF
264
265
266
k = LEN_TRIM(str)
267
l = 1
268
n = 0
269
DO WHILE(l<=k.AND.n<SIZE(ints))
270
DO WHILE(l<=k)
271
ic = ICHAR(str(l:l))
272
IF( ic == ics .OR. ic == icsep ) THEN
273
CONTINUE
274
ELSE
275
EXIT
276
END IF
277
l=l+1
278
END DO
279
IF(l>k) EXIT
280
IF(.NOT.(ic==icm .OR. ic>=ic0 .AND. ic<=ic9)) EXIT
281
282
m = l+1
283
DO WHILE(m<=k)
284
ic = ICHAR(str(m:m))
285
IF(ic<ic0 .OR. ic>ic9) EXIT
286
m=m+1
287
END DO
288
289
n = n + 1
290
ints(n) = s2i(str(l:m-1),m-l)
291
l = m
292
END DO
293
END FUNCTION str2ints
294
!------------------------------------------------------------------------------
295
296
SUBROUTINE WaitSec(t)
297
REAL(KIND=dp) :: t,t0,t1
298
299
t0 = RealTime()
300
DO WHILE(.TRUE.)
301
t1 = RealTime()
302
IF(t1-t0 > t) EXIT
303
END DO
304
305
END SUBROUTINE WaitSec
306
307
!------------------------------------------------------------------------------
308
SUBROUTINE SystemCommand( cmd )
309
!------------------------------------------------------------------------------
310
CHARACTER(LEN=*) :: cmd
311
CALL SystemC( TRIM(cmd) // CHAR(0) )
312
!------------------------------------------------------------------------------
313
END SUBROUTINE SystemCommand
314
!------------------------------------------------------------------------------
315
316
317
!------------------------------------------------------------------------------
318
SUBROUTINE RenameF(old,new)
319
!------------------------------------------------------------------------------
320
CHARACTER(LEN=*) :: old,new
321
!------------------------------------------------------------------------------
322
INTERFACE
323
SUBROUTINE rename_c(old,new) bind(c,name="rename_c")
324
USE, INTRINSIC :: iso_c_binding
325
CHARACTER(KIND=C_CHAR) :: old(*), new(*)
326
END SUBROUTINE rename_c
327
END INTERFACE
328
329
CALL rename_c(TRIM(old)//CHAR(0), TRIM(new)//CHAR(0))
330
!------------------------------------------------------------------------------
331
END SUBROUTINE RenameF
332
!------------------------------------------------------------------------------
333
334
335
!------------------------------------------------------------------------------
336
PURE FUNCTION FileNameQualified(file) RESULT(L)
337
!------------------------------------------------------------------------------
338
LOGICAL :: L
339
CHARACTER(*), INTENT(IN) :: file
340
!------------------------------------------------------------------------------
341
L = INDEX(file,':')>0 .OR. file(1:1)=='/' .OR. file(1:1)==Backslash
342
!------------------------------------------------------------------------------
343
END FUNCTION FileNameQualified
344
!------------------------------------------------------------------------------
345
346
!------------------------------------------------------------------------------
347
PURE FUNCTION LittleEndian() RESULT(L)
348
!------------------------------------------------------------------------------
349
LOGICAL :: L
350
!------------------------------------------------------------------------------
351
INTEGER(1) :: s(2)
352
INTEGER(2), PARAMETER :: t = 256*7+8
353
!------------------------------------------------------------------------------
354
s = TRANSFER(t,s)
355
L=s(1)==8
356
!------------------------------------------------------------------------------
357
END FUNCTION LittleEndian
358
!------------------------------------------------------------------------------
359
360
!------------------------------------------------------------------------------
361
FUNCTION FormatDate() RESULT( date )
362
!------------------------------------------------------------------------------
363
CHARACTER(20) :: date
364
INTEGER :: dates(8)
365
366
CALL DATE_AND_TIME( VALUES=dates )
367
WRITE( date, &
368
'(I4,"/",I2.2,"/",I2.2," ",I2.2,":",I2.2,":",I2.2)' ) &
369
dates(1),dates(2),dates(3),dates(5),dates(6),dates(7)
370
!------------------------------------------------------------------------------
371
END FUNCTION FormatDate
372
!------------------------------------------------------------------------------
373
374
375
!------------------------------------------------------------------------------
376
!> Sort an array of integer values.
377
!------------------------------------------------------------------------------
378
PURE SUBROUTINE Sort( n,a )
379
!------------------------------------------------------------------------------
380
INTEGER, INTENT(in) :: n
381
INTEGER, INTENT(inout) :: a(:)
382
!------------------------------------------------------------------------------
383
384
INTEGER :: i,j,l,ir,ra
385
!------------------------------------------------------------------------------
386
387
IF ( n <= 1 ) RETURN
388
389
l = n / 2 + 1
390
ir = n
391
DO WHILE( .TRUE. )
392
IF ( l > 1 ) THEN
393
l = l - 1
394
ra = a(l)
395
ELSE
396
ra = a(ir)
397
a(ir) = a(1)
398
ir = ir - 1
399
IF ( ir == 1 ) THEN
400
a(1) = ra
401
RETURN
402
END IF
403
END IF
404
i = l
405
j = l + l
406
DO WHILE( j <= ir )
407
IF ( j<ir ) THEN
408
IF ( a(j)<a(j+1) ) j = j+1
409
END IF
410
411
IF ( ra<a(j) ) THEN
412
a(i) = a(j)
413
i = j
414
j = j + i
415
ELSE
416
j = ir + 1
417
END IF
418
a(i) = ra
419
END DO
420
END DO
421
422
!------------------------------------------------------------------------------
423
END SUBROUTINE Sort
424
!------------------------------------------------------------------------------
425
426
!------------------------------------------------------------------------------
427
!> Sort an integer array a, together with an another integer array.
428
!------------------------------------------------------------------------------
429
PURE SUBROUTINE SortI( n,a,b )
430
!------------------------------------------------------------------------------
431
INTEGER, INTENT(in) :: n
432
INTEGER, INTENT(inout) :: a(:),b(:)
433
!------------------------------------------------------------------------------
434
435
INTEGER :: i,j,l,ir,ra,rb
436
!------------------------------------------------------------------------------
437
438
IF ( n <= 1 ) RETURN
439
440
l = n / 2 + 1
441
ir = n
442
DO WHILE( .TRUE. )
443
IF ( l > 1 ) THEN
444
l = l - 1
445
ra = a(l)
446
rb = b(l)
447
ELSE
448
ra = a(ir)
449
rb = b(ir)
450
a(ir) = a(1)
451
b(ir) = b(1)
452
ir = ir - 1
453
IF ( ir == 1 ) THEN
454
a(1) = ra
455
b(1) = rb
456
RETURN
457
END IF
458
END IF
459
i = l
460
j = l + l
461
DO WHILE( j <= ir )
462
IF ( j<ir ) THEN
463
IF ( a(j)<a(j+1) ) j = j+1
464
END IF
465
IF ( ra<a(j) ) THEN
466
a(i) = a(j)
467
b(i) = b(j)
468
i = j
469
j = j + i
470
ELSE
471
j = ir + 1
472
END IF
473
a(i) = ra
474
b(i) = rb
475
END DO
476
END DO
477
478
!------------------------------------------------------------------------------
479
END SUBROUTINE SortI
480
!------------------------------------------------------------------------------
481
482
483
!------------------------------------------------------------------------------
484
!> Sort an index array, and change the order of an real array accordingly.
485
!------------------------------------------------------------------------------
486
PURE SUBROUTINE SortF( n,a,b )
487
!------------------------------------------------------------------------------
488
INTEGER, INTENT(in) :: n
489
INTEGER, INTENT(inout) :: a(:)
490
REAL(KIND=dp), INTENT(inout) :: b(:)
491
!------------------------------------------------------------------------------
492
493
INTEGER :: i,j,l,ir,ra
494
REAL(KIND=dp) :: rb
495
!------------------------------------------------------------------------------
496
497
IF ( n <= 1 ) RETURN
498
499
l = n / 2 + 1
500
ir = n
501
DO WHILE( .TRUE. )
502
503
IF ( l > 1 ) THEN
504
l = l - 1
505
ra = a(l)
506
rb = b(l)
507
ELSE
508
ra = a(ir)
509
rb = b(ir)
510
a(ir) = a(1)
511
b(ir) = b(1)
512
ir = ir - 1
513
IF ( ir == 1 ) THEN
514
a(1) = ra
515
b(1) = rb
516
RETURN
517
END IF
518
END IF
519
i = l
520
j = l + l
521
DO WHILE( j <= ir )
522
IF ( j<ir ) THEN
523
IF ( a(j)<a(j+1) ) j = j+1
524
END IF
525
IF ( ra<a(j) ) THEN
526
a(i) = a(j)
527
b(i) = b(j)
528
i = j
529
j = j + i
530
ELSE
531
j = ir + 1
532
END IF
533
a(i) = ra
534
b(i) = rb
535
END DO
536
END DO
537
538
!------------------------------------------------------------------------------
539
END SUBROUTINE SortF
540
!------------------------------------------------------------------------------
541
542
543
!------------------------------------------------------------------------------
544
!> Sort an real array, and change the order of an index array accordingly.
545
!------------------------------------------------------------------------------
546
PURE SUBROUTINE SortD( n,a,b )
547
!------------------------------------------------------------------------------
548
INTEGER, INTENT(in) :: n
549
INTEGER, INTENT(inout) :: b(:)
550
REAL(KIND=dp), INTENT(inout) :: a(:)
551
!------------------------------------------------------------------------------
552
553
INTEGER :: i,j,l,ir,rb
554
REAL(KIND=dp) :: ra
555
!------------------------------------------------------------------------------
556
557
IF ( n <= 1 ) RETURN
558
559
l = n / 2 + 1
560
ir = n
561
DO WHILE( .TRUE. )
562
563
IF ( l > 1 ) THEN
564
l = l - 1
565
ra = a(l)
566
rb = b(l)
567
ELSE
568
ra = a(ir)
569
rb = b(ir)
570
a(ir) = a(1)
571
b(ir) = b(1)
572
ir = ir - 1
573
IF ( ir == 1 ) THEN
574
a(1) = ra
575
b(1) = rb
576
RETURN
577
END IF
578
END IF
579
i = l
580
j = l + l
581
DO WHILE( j <= ir )
582
IF ( j<ir ) THEN
583
IF ( a(j)<a(j+1) ) j = j+1
584
END IF
585
IF ( ra<a(j) ) THEN
586
a(i) = a(j)
587
b(i) = b(j)
588
i = j
589
j = j + i
590
ELSE
591
j = ir + 1
592
END IF
593
a(i) = ra
594
b(i) = rb
595
END DO
596
END DO
597
598
!------------------------------------------------------------------------------
599
END SUBROUTINE SortD
600
!------------------------------------------------------------------------------
601
602
603
!------------------------------------------------------------------------------
604
!> Sort an complex array, and organize an index table accordingly.
605
!------------------------------------------------------------------------------
606
PURE SUBROUTINE SortC( n,a,b )
607
!------------------------------------------------------------------------------
608
INTEGER, INTENT(in) :: n
609
INTEGER, INTENT(inout) :: b(:)
610
COMPLEX(KIND=dp), INTENT(inout) :: a(:)
611
!------------------------------------------------------------------------------
612
613
INTEGER :: i,j,l,ir,rb
614
COMPLEX(KIND=dp) :: ra
615
!------------------------------------------------------------------------------
616
617
IF ( n <= 1 ) RETURN
618
619
l = n / 2 + 1
620
ir = n
621
DO WHILE( .TRUE. )
622
IF ( l > 1 ) THEN
623
l = l - 1
624
ra = a(l)
625
rb = b(l)
626
ELSE
627
ra = a(ir)
628
rb = b(ir)
629
a(ir) = a(1)
630
b(ir) = b(1)
631
ir = ir - 1
632
IF ( ir == 1 ) THEN
633
a(1) = ra
634
b(1) = rb
635
RETURN
636
END IF
637
END IF
638
i = l
639
j = l + l
640
DO WHILE( j <= ir )
641
IF ( j<ir ) THEN
642
IF ( ABS(a(j))<ABS(a(j+1)) ) j = j+1
643
END IF
644
IF ( ABS(ra)<ABS(a(j)) ) THEN
645
a(i) = a(j)
646
b(i) = b(j)
647
i = j
648
j = j + i
649
ELSE
650
j = ir + 1
651
END IF
652
a(i) = ra
653
b(i) = rb
654
END DO
655
END DO
656
657
!------------------------------------------------------------------------------
658
END SUBROUTINE SortC
659
!------------------------------------------------------------------------------
660
661
662
!------------------------------------------------------------------------------
663
!> Order real components in b in a decreasing order and return the new order
664
!> of indexes in a.
665
!------------------------------------------------------------------------------
666
PURE SUBROUTINE SortR( n,a,b )
667
!------------------------------------------------------------------------------
668
INTEGER, INTENT(in) :: n
669
INTEGER, INTENT(inout) :: a(:)
670
REAL(KIND=dp), INTENT(inout) :: b(:)
671
!------------------------------------------------------------------------------
672
673
INTEGER :: i,j,l,ir,ra
674
REAL(KIND=dp) :: rb
675
!------------------------------------------------------------------------------
676
677
IF ( n <= 1 ) RETURN
678
679
l = n / 2 + 1
680
ir = n
681
DO WHILE( .TRUE. )
682
683
IF ( l > 1 ) THEN
684
l = l - 1
685
ra = a(l)
686
rb = b(l)
687
ELSE
688
ra = a(ir)
689
rb = b(ir)
690
a(ir) = a(1)
691
b(ir) = b(1)
692
ir = ir - 1
693
IF ( ir == 1 ) THEN
694
a(1) = ra
695
b(1) = rb
696
RETURN
697
END IF
698
END IF
699
i = l
700
j = l + l
701
DO WHILE( j <= ir )
702
IF ( j<ir ) THEN
703
IF ( b(j) > b(j+1) ) j = j+1
704
END IF
705
IF ( rb > b(j) ) THEN
706
a(i) = a(j)
707
b(i) = b(j)
708
i = j
709
j = j + i
710
ELSE
711
j = ir + 1
712
END IF
713
a(i) = ra
714
b(i) = rb
715
END DO
716
END DO
717
718
!------------------------------------------------------------------------------
719
END SUBROUTINE SortR
720
!------------------------------------------------------------------------------
721
722
!------------------------------------------------------------------------------
723
!> Search an integer value in an ordered array.
724
!------------------------------------------------------------------------------
725
PURE FUNCTION SearchI( N,Array,Val ) RESULT ( Idx )
726
!------------------------------------------------------------------------------
727
INTEGER, INTENT(in) :: N,Val,Array(:)
728
!------------------------------------------------------------------------------
729
INTEGER :: Lower, Upper,Lou,Idx
730
!------------------------------------------------------------------------------
731
732
Idx = 0
733
Upper = N
734
Lower = 1
735
736
! Handle the special case
737
738
IF ( Upper == 0 ) RETURN
739
740
DO WHILE( .TRUE. )
741
IF ( Array(Lower) == Val) THEN
742
Idx = Lower
743
EXIT
744
ELSE IF ( Array(Upper) == Val ) THEN
745
Idx = Upper
746
EXIT
747
END IF
748
749
IF ( (Upper-Lower)>1 ) THEN
750
Lou = ISHFT((Upper + Lower), -1)
751
IF ( Array(Lou) < Val ) THEN
752
Lower = Lou
753
ELSE
754
Upper = Lou
755
END IF
756
ELSE
757
EXIT
758
END IF
759
END DO
760
761
RETURN
762
763
!------------------------------------------------------------------------------
764
END FUNCTION SearchI
765
!------------------------------------------------------------------------------
766
767
768
769
!------------------------------------------------------------------------------
770
!> Search a real value in an ordered array.
771
!------------------------------------------------------------------------------
772
PURE FUNCTION SearchR( N,Array,Val ) RESULT ( Idx )
773
!------------------------------------------------------------------------------
774
775
INTEGER, INTENT(in) :: N
776
INTEGER :: Idx
777
REAL(KIND=dP), INTENT(in) :: Val,Array(:)
778
!------------------------------------------------------------------------------
779
INTEGER :: Lower, Upper,Lou
780
!------------------------------------------------------------------------------
781
782
Idx = 0
783
Upper = N
784
Lower = 1
785
786
! Handle the special case
787
IF ( Upper == 0 ) RETURN
788
789
DO WHILE( .TRUE. )
790
IF ( ABS( Array(Lower) - Val) < TINY(Val) ) THEN
791
Idx = Lower
792
EXIT
793
ELSE IF ( ABS( Array(Upper) - Val ) < TINY(Val) ) THEN
794
Idx = Upper
795
EXIT
796
END IF
797
798
IF ( (Upper-Lower) > 1 ) THEN
799
Lou = ISHFT((Upper + Lower), -1)
800
IF ( Array(Lou) < Val ) THEN
801
Lower = Lou
802
ELSE
803
Upper = Lou
804
END IF
805
ELSE
806
EXIT
807
END IF
808
END DO
809
810
RETURN
811
812
!------------------------------------------------------------------------------
813
END FUNCTION SearchR
814
!------------------------------------------------------------------------------
815
816
817
!------------------------------------------------------------------------------
818
SUBROUTINE OpenIncludeFile( Unit, FileName, IncludePath )
819
!------------------------------------------------------------------------------
820
INTEGER :: Unit
821
CHARACTER(LEN=*) :: FileName, IncludePath
822
!------------------------------------------------------------------------------
823
INTEGER :: i,j,k,k0,k1,l,iostat
824
CHARACTER(LEN=1024) :: name, TmpName
825
!------------------------------------------------------------------------------
826
827
i = 1
828
name = FileName
829
DO WHILE( name(i:i) == ' ' .OR. name(i:i)=='"')
830
i = i + 1
831
END DO
832
j = LEN_TRIM(name)
833
IF ( name(j:j) == '"' ) j=j-1
834
name = TRIM(name(i:j))
835
836
IF ( INDEX(name,':') == 0 .AND. name(1:1) /= '/' .AND. &
837
name(1:1) /= Backslash ) THEN
838
k0 = 1
839
DO WHILE( IncludePath(k0:k0) == '"' )
840
k0 = k0+1
841
END DO
842
k1 = INDEX( IncludePath, ';' )
843
844
DO WHILE( k1 >= k0 )
845
DO k = k1-1,k0,-1
846
IF ( IncludePath(k:k) /= ' ' .AND. IncludePath(k:k)/='"' ) EXIT
847
END DO
848
IF ( IncludePath(k:k) == '"' ) k=k-1
849
IF ( k >= k0 ) THEN
850
WRITE( tmpName,'(a,a,a)' ) IncludePath(k0:k), '/', TRIM(name)
851
OPEN( Unit, FILE=TRIM(tmpName), STATUS='OLD',ERR=10 )
852
RETURN
853
END IF
854
10 CONTINUE
855
k0 = k1+1
856
k1 = INDEX( IncludePath(k0:), ';' ) + k0 - 1
857
END DO
858
859
IF ( LEN_TRIM(IncludePath(k0:))>0 ) THEN
860
k1 = INDEX( IncludePath(k0:), '"' ) + k0 - 2
861
IF ( k1 < k0 ) k1=LEN_TRIM(IncludePath)
862
tmpName = TRIM(IncludePath(k0:k1)) // '/' // TRIM(name)
863
OPEN( Unit, FILE=TRIM(TmpName), STATUS='OLD',ERR=20 )
864
RETURN
865
END IF
866
867
20 CONTINUE
868
OPEN( Unit, FILE=TRIM(name), STATUS='OLD',IOSTAT=iostat )
869
ELSE
870
OPEN( Unit, FILE=TRIM(name), STATUS='OLD',IOSTAT=iostat )
871
END IF
872
873
IF( iostat /= 0 ) THEN
874
CALL Fatal('OpenIncludeFile','Cannot open include file: '//TRIM(Name))
875
END IF
876
877
878
!------------------------------------------------------------------------------
879
END SUBROUTINE OpenIncludeFile
880
!------------------------------------------------------------------------------
881
882
883
!------------------------------------------------------------------------------
884
!> Read a (logical) line from FORTRAN device Unit and remove leading, trailing,
885
!> and multiple blanks between words. Also convert uppercase characters to
886
!> lowercase.The logical line can continue the several physical lines by adding
887
!> the backslash (\) mark at the end of a physical line.
888
!------------------------------------------------------------------------------
889
RECURSIVE FUNCTION ReadAndTrim( Unit,str,echo,literal,noeval ) RESULT(l)
890
!------------------------------------------------------------------------------
891
INTEGER :: Unit !< Fortran unit number to read from
892
CHARACTER(LEN=:), ALLOCATABLE :: str !< The string read from the file
893
LOGICAL, OPTIONAL :: Echo
894
LOGICAL, OPTIONAL :: literal
895
LOGICAL, OPTIONAL :: noeval
896
LOGICAL :: l !< Success of the read operation
897
!------------------------------------------------------------------------------
898
INTEGER, PARAMETER :: MAXLEN = 163840
899
900
CHARACTER(LEN=:), ALLOCATABLE :: temp
901
CHARACTER(LEN=12) :: tmpstr
902
CHARACTER(LEN=MAXLEN) :: readstr = ' ', copystr = ' ', matcstr=' ' , IncludePath=' '
903
904
LOGICAL :: InsideQuotes, OpenSection=.FALSE., DoEval
905
INTEGER :: i,j,k,m,ValueStarts=0,inlen,ninlen,outlen,IncludeUnit=28,IncludeUnitBase=28
906
907
CHARACTER(LEN=MAX_NAME_LEN) :: Prefix = ' '
908
909
INTEGER, PARAMETER :: A=ICHAR('A'),Z=ICHAR('Z'),U2L=ICHAR('a')-ICHAR('A'),Tab=9
910
CHARACTER(LEN=MAXLEN) :: tmatcstr, tcmdstr
911
INTEGER :: tninlen
912
913
SAVE ReadStr, ValueStarts, Prefix, OpenSection
914
915
IF ( PRESENT(literal) ) literal=.FALSE.
916
l = .TRUE.
917
918
! Optionally do not expand the MATC and LUA expressions.
919
DoEval = .TRUE.
920
IF( PRESENT( NoEval ) ) THEN
921
DoEval = .NOT. NoEval
922
END IF
923
924
IF(.NOT.ALLOCATED(str)) ALLOCATE(CHARACTER(512)::str)
925
outlen = LEN(str)
926
927
IF ( ValueStarts==0 .AND. OpenSection ) THEN
928
str = 'end'
929
ValueStarts = 0
930
OpenSection = .FALSE.
931
RETURN
932
END IF
933
934
IF ( ValueStarts == 0 ) THEN
935
tmpstr = ' '
936
DO WHILE( .TRUE. )
937
IF ( IncludeUnit < IncludeUnitBase ) THEN
938
READ( IncludeUnit,'(A)',END=1,ERR=1 ) readstr
939
GO TO 2
940
1 CLOSE(IncludeUnit)
941
IncludeUnit = IncludeUnit+1
942
READ( IncludeUnit,'(A)',END=10,ERR=10 ) readstr
943
2 CONTINUE
944
ELSE
945
READ( Unit,'(A)',END=10,ERR=10 ) readstr
946
END IF
947
948
readstr = ADJUSTL(readstr)
949
950
DO k=1,12
951
j = ICHAR(readstr(k:k))
952
IF ( j >= A .AND. j<= Z ) THEN
953
Tmpstr(k:k) = CHAR(j+U2L)
954
ELSE
955
tmpstr(k:k) = readstr(k:k)
956
END IF
957
END DO
958
959
IF ( SEQL(Tmpstr, 'include path') ) THEN
960
k = LEN_TRIM(readstr)
961
IncludePath(1:k-13) = readstr(14:k)
962
Tmpstr = ''
963
ELSE
964
EXIT
965
END IF
966
END DO
967
968
IF ( SEQL(tmpstr, 'include ') ) THEN
969
IncludeUnit = IncludeUnit-1
970
971
CALL Info('OpenIncludeFile','Trying to include file: '//TRIM(readstr(9:)),Level=10)
972
973
inlen = LEN_TRIM(readstr)
974
975
i = INDEX( readstr(1:inlen), '$' )
976
IF ( i>0 .AND. i<inlen ) THEN
977
CALL TrimMatcExpression()
978
CALL Info('ReadAndTrim','Include file after MATC trimming: '&
979
//TRIM(readstr(9:)),Level=6)
980
END IF
981
982
i = INDEX( readstr(1:inlen),'#')
983
IF( i>0 .AND. i<inlen ) THEN
984
#ifdef HAVE_LUA
985
CALL TrimLuaExpression()
986
CALL Info('ReadAndTrim','Include file after LUA trimming: '&
987
//TRIM(readstr(9:)),Level=6)
988
#else
989
CALL Fatal('ReadAndTrim','LUA not included, cannot continue')
990
#endif
991
END IF
992
993
CALL OpenIncludeFile( IncludeUnit, TRIM(readstr(9:)), IncludePath )
994
995
READ( IncludeUnit,'(A)',END=3,ERR=3 ) readstr
996
GO TO 4
997
3 CLOSE(IncludeUnit)
998
IncludeUnit = IncludeUnit+1
999
READ( Unit,'(A)',END=10,ERR=10 ) readstr
1000
4 CONTINUE
1001
END IF
1002
ninlen = LEN_TRIM(readstr)
1003
ELSE
1004
inlen = LEN_TRIM(readstr)
1005
ninlen = inlen-ValueStarts+1
1006
IF ( Prefix == ' ' ) THEN
1007
readstr = readstr(ValueStarts:inlen)
1008
ELSE IF ( Prefix == '::' ) THEN
1009
readstr = readstr(ValueStarts:inlen)
1010
OpenSection = .TRUE.
1011
Prefix = ' '
1012
ELSE
1013
DO i=ValueStarts,inlen
1014
IF ( readstr(i:i) == ')' ) THEN
1015
readstr(i:i) = ' '
1016
EXIT
1017
ELSE IF ( readstr(i:i) == ',' ) THEN
1018
readstr(i:i) = ' '
1019
END IF
1020
END DO
1021
ninlen = ninlen + LEN_TRIM(Prefix) + 1
1022
readstr = TRIM(Prefix) // ' ' // readstr(ValueStarts:inlen)
1023
END IF
1024
END IF
1025
1026
ValueStarts = 0
1027
InsideQuotes = .FALSE.
1028
1029
i = INDEX( readstr(1:ninlen), '!' )
1030
IF ( i>0 ) ninlen=i-1
1031
1032
i = 1
1033
inlen = ninlen
1034
DO WHILE( i <= inlen )
1035
IF ( readstr(i:i) == '"' ) InsideQuotes = .NOT.InsideQuotes
1036
IF ( .NOT. InsideQuotes .AND. readstr(i:i) == CHAR(92) .AND. i==inlen ) THEN
1037
readstr(i:i) = ' '
1038
IF ( IncludeUnit < IncludeUnitBase ) THEN
1039
READ( IncludeUnit,'(A)',END=10,ERR=10 ) readstr(i+1:MAXLEN)
1040
ELSE
1041
READ( Unit,'(A)',END=10,ERR=10 ) readstr(i+1:MAXLEN)
1042
END IF
1043
DO j=LEN(readstr),i+1,-1
1044
IF ( readstr(j:j) /= ' ' ) EXIT
1045
END DO
1046
inlen = inlen + j-i
1047
END IF
1048
i = i + 1
1049
END DO
1050
1051
i = INDEX( readstr(1:inlen), '#' )
1052
IF ( i>0 .AND. i<inlen .AND. DoEval ) THEN
1053
#ifdef HAVE_LUA
1054
CALL TrimLuaExpression()
1055
#else
1056
CALL Fatal('ReadAndTrim','LUA not included, cannot continue')
1057
#endif
1058
END IF
1059
1060
i = INDEX( readstr(1:inlen), '$' )
1061
IF ( i>0 .AND. i<inlen .AND. DoEval ) THEN
1062
CALL TrimMatcExpression()
1063
END IF
1064
1065
IF ( PRESENT( Echo ) ) THEN
1066
IF ( Echo .AND. inlen > 0 ) WRITE( 6, '(a)' ) readstr(1:inlen)
1067
END IF
1068
1069
i = 1
1070
DO WHILE(i <= inlen )
1071
IF (readstr(i:i) /= ' ' .AND. ICHAR(readstr(i:i))/=Tab ) EXIT
1072
i = i + 1
1073
END DO
1074
1075
InsideQuotes = .FALSE.
1076
1077
IF ( PRESENT(literal) ) THEN
1078
IF ( readstr(i:i) == '"' ) literal=.TRUE.
1079
END IF
1080
1081
k = 1
1082
DO WHILE( i<=inlen )
1083
IF ( readstr(i:i) == '"' ) THEN
1084
InsideQuotes = .NOT.InsideQuotes
1085
i=i+1
1086
IF ( i>inlen ) EXIT
1087
END IF
1088
1089
IF ( .NOT.InsideQuotes ) THEN
1090
IF ( readstr(i:i) == '!' .OR. readstr(i:i) == '#' .OR. &
1091
readstr(i:i) == '=' .OR. readstr(i:i) == '(' .OR. &
1092
readstr(i:i) == ';' .OR. readstr(i:i+1) == '::' ) EXIT
1093
IF (ICHAR( readstr(i:i))<32.AND.ICHAR(readstr(i:i))/=Tab) EXIT
1094
END IF
1095
1096
DO WHILE( i <= inlen )
1097
IF ( readstr(i:i) == '"' ) THEN
1098
InsideQuotes = .NOT.InsideQuotes
1099
i=i+1
1100
IF ( i>inlen ) EXIT
1101
END IF
1102
1103
IF ( .NOT.InsideQuotes ) THEN
1104
IF ( readstr(i:i) == ' ' .OR. readstr(i:i) == '=' .OR. &
1105
readstr(i:i) == ';' .OR. readstr(i:i) == '(' .OR. &
1106
readstr(i:i+1) == '::' ) EXIT
1107
IF ( ICHAR( readstr(i:i))<32 ) EXIT
1108
END IF
1109
1110
IF ( k>outlen ) THEN
1111
temp = str
1112
DEALLOCATE(str)
1113
outlen=LEN(temp)+512
1114
ALLOCATE(CHARACTER(outlen)::str)
1115
str(1:LEN(temp))=temp; str(LEN(temp)+1:)=''
1116
DEALLOCATE(temp)
1117
END IF
1118
1119
j = ICHAR( readstr(i:i) )
1120
IF ( .NOT.InsideQuotes .AND. j>=A .AND. j<=Z ) THEN
1121
str(k:k) = CHAR(j+U2L)
1122
ELSE IF ( .NOT.InsideQuotes .AND. j==Tab ) THEN
1123
str(k:k) = ' '
1124
ELSE
1125
str(k:k) = readstr(i:i)
1126
ENDIF
1127
1128
i = i + 1
1129
k = k + 1
1130
END DO
1131
1132
IF ( k <= outlen ) str(k:k) = ' '
1133
k = k + 1
1134
1135
DO WHILE( i<=inlen )
1136
IF ( readstr(i:i) /= ' ' .AND. ICHAR(readstr(i:i))/=Tab ) EXIT
1137
i = i + 1
1138
END DO
1139
END DO
1140
str(k:)=' '
1141
1142
IF ( i <= inlen ) THEN
1143
Prefix = ' '
1144
IF ( ReadStr(i:i) == '=' ) THEN
1145
ValueStarts = i + 1
1146
ELSE IF ( ReadStr(i:i) == ';' ) THEN
1147
ValueStarts = i + 1
1148
ELSE IF ( ReadStr(i:i) == '(' ) THEN
1149
ValueStarts = i + 1
1150
Prefix = 'Size'
1151
ELSE IF ( ReadStr(i:i+1) == '::' ) THEN
1152
ValueStarts = i + 2
1153
Prefix = '::'
1154
ELSE IF ( ICHAR(readstr(i:i)) < 32 ) THEN
1155
DO WHILE( i <= inlen )
1156
IF ( ICHAR(readstr(i:i)) >= 32 ) EXIT
1157
i = i + 1
1158
END DO
1159
IF ( i <= inlen ) ValueStarts = i
1160
END IF
1161
END IF
1162
RETURN
1163
10 CONTINUE
1164
l = .FALSE.
1165
!------------------------------------------------------------------------------
1166
1167
CONTAINS
1168
1169
SUBROUTINE TrimMatcExpression()
1170
1171
i = INDEX( readstr(1:inlen), '$' )
1172
IF ( i>0 .AND. i<inlen ) THEN
1173
m = i
1174
copystr(i:inlen) = readstr(i:inlen)
1175
DO WHILE(i<=inlen)
1176
IF ( copystr(i:i) == '$' ) THEN
1177
DO j=i+1,inlen-1
1178
IF ( copystr(j:j) == '$' ) EXIT
1179
END DO
1180
ninlen = j - i
1181
! Initialize variables for each copy of MATC separately
1182
1183
!$OMP PARALLEL DEFAULT(NONE) &
1184
!$OMP SHARED(copystr, i, matcstr, ninlen, inlen) &
1185
!$OMP PRIVATE(tcmdstr, tmatcstr, tninlen)
1186
1187
tninlen = ninlen
1188
tcmdstr = copystr(i+1:inlen)
1189
tninlen = MATC(tcmdstr, tmatcstr, tninlen)
1190
!$OMP BARRIER
1191
1192
!$OMP SINGLE
1193
matcstr(1:tninlen) = tmatcstr(1:tninlen)
1194
ninlen = tninlen
1195
!$OMP END SINGLE
1196
1197
!$OMP END PARALLEL
1198
1199
DO k=1,ninlen
1200
readstr(m:m) = matcstr(k:k)
1201
m = m + 1
1202
END DO
1203
i = j+1
1204
ELSE
1205
readstr(m:m) = copystr(i:i)
1206
i = i + 1
1207
m = m + 1
1208
END IF
1209
END DO
1210
inlen = m-1
1211
readstr(inlen+1:) = ' '
1212
END IF
1213
1214
END SUBROUTINE TrimMatcExpression
1215
1216
#ifdef HAVE_LUA
1217
SUBROUTINE TrimLuaExpression()
1218
1219
INTEGER :: lstat
1220
character(kind=c_char, len=:), pointer :: lua_result
1221
integer :: result_len
1222
logical :: closed_region, first_bang
1223
closed_region = .false.
1224
first_bang = .true.
1225
m = i
1226
copystr(i:inlen) = readstr(i:inlen)
1227
DO WHILE(i<=inlen)
1228
IF ( copystr(i:i) == '#' ) THEN
1229
DO j=i+1,inlen-1
1230
IF ( copystr(j:j) == '#' ) EXIT
1231
END DO
1232
ninlen = j - i
1233
1234
! Initialize variables for each copy of Lua interpreter separately
1235
1236
!$OMP PARALLEL DEFAULT(NONE) &
1237
!$OMP SHARED(copystr, i, matcstr, ninlen, inlen, closed_region, first_bang, j) &
1238
!$OMP PRIVATE(tcmdstr, tninlen, lstat, result_len, lua_result)
1239
1240
tninlen = ninlen
1241
tcmdstr = copystr(i+1:inlen)
1242
1243
IF(tcmdstr(tninlen:tninlen) == '#') then
1244
closed_region = .TRUE.
1245
ELSE
1246
closed_region = .FALSE.
1247
END IF
1248
1249
IF(closed_region) THEN
1250
lstat = lua_dostring( LuaState, &
1251
'return tostring('// tcmdstr(1:tninlen-1) // ')'//c_null_char, 1)
1252
ELSE
1253
IF (i == 1 .and. first_bang .and. j == inlen) THEN ! ' # <luacode>' case, do not do 'return tostring(..)'.
1254
! Instead, just execute the line in the lua interpreter
1255
lstat = lua_dostring( LuaState, tcmdstr(1:tninlen) // c_null_char, 1)
1256
ELSE ! 'abc = # <luacode>' case, oneliners only
1257
lstat = lua_dostring( LuaState, &
1258
'return tostring('// tcmdstr(1:tninlen) // ')'//c_null_char, 1)
1259
END IF
1260
END IF
1261
lua_result => lua_popstring(LuaState, result_len)
1262
1263
!$OMP SINGLE
1264
matcstr(1:result_len) = lua_result(1:result_len)
1265
ninlen = result_len
1266
!$OMP END SINGLE
1267
1268
!$OMP END PARALLEL
1269
1270
DO k=1,ninlen
1271
readstr(m:m) = matcstr(k:k)
1272
m = m + 1
1273
END DO
1274
i = j+1
1275
ELSE
1276
readstr(m:m) = copystr(i:i)
1277
i = i + 1
1278
m = m + 1
1279
END IF
1280
first_bang = .false.
1281
END DO
1282
inlen = m-1
1283
readstr(inlen+1:) = ' '
1284
END SUBROUTINE TrimLuaExpression
1285
#endif
1286
1287
END FUNCTION ReadAndTrim
1288
!------------------------------------------------------------------------------
1289
1290
1291
!------------------------------------------------------------------------------
1292
PURE FUNCTION GetVarName(Var) RESULT(str)
1293
!------------------------------------------------------------------------------
1294
TYPE(Variable_t), INTENT(in) :: Var
1295
CHARACTER(LEN=Var % NameLen) :: str
1296
!------------------------------------------------------------------------------
1297
str = Var % Name(1:Var % NameLen)
1298
!------------------------------------------------------------------------------
1299
END FUNCTION GetVarName
1300
!------------------------------------------------------------------------------
1301
1302
1303
!------------------------------------------------------------------------------
1304
FUNCTION ComponentNameVar( Var, Component ) RESULT(str)
1305
!------------------------------------------------------------------------------
1306
TYPE(Variable_t),INTENT(in) :: Var
1307
INTEGER, OPTIONAL,INTENT(in) :: Component
1308
!------------------------------------------------------------------------------
1309
CHARACTER(:), ALLOCATABLE :: str
1310
!------------------------------------------------------------------------------
1311
IF ( Var % Name(1:Var % NameLen) == 'flow solution' ) THEN
1312
str='flow solution'
1313
IF ( .NOT. PRESENT(Component) ) RETURN
1314
IF ( Component == Var % DOFs ) THEN
1315
str = 'pressure'
1316
RETURN
1317
ELSE
1318
str = 'velocity ' // i2s(Component)
1319
END IF
1320
ELSE
1321
str = ComponentName(Var % Name, Component)
1322
END IF
1323
!------------------------------------------------------------------------------
1324
END FUNCTION ComponentNameVar
1325
!------------------------------------------------------------------------------
1326
1327
1328
!------------------------------------------------------------------------------
1329
FUNCTION ComponentNameStr( BaseName, Component_arg ) RESULT(str)
1330
!------------------------------------------------------------------------------
1331
INTEGER, OPTIONAL, INTENT(in) :: Component_arg
1332
CHARACTER(:), ALLOCATABLE :: str
1333
CHARACTER(LEN=*), INTENT(in) :: BaseName
1334
!------------------------------------------------------------------------------
1335
INTEGER :: ind, ind1, DOFsTot, DOFs, Component
1336
!------------------------------------------------------------------------------
1337
ind = INDEX( BaseName,'[' )
1338
1339
Component = 0
1340
IF ( PRESENT(Component_arg) ) Component=Component_arg
1341
1342
IF ( ind<=0 ) THEN
1343
str = TRIM(BaseName)
1344
IF ( Component > 0 ) THEN
1345
str = str // ' ' // i2s(Component)
1346
END IF
1347
ELSE IF( Component == 0 ) THEN
1348
str = BaseName(1:ind-1)
1349
ELSE
1350
DOFsTot = 0
1351
DO WHILE( .TRUE. )
1352
ind1 = INDEX( BaseName(ind+1:),':' )+ind
1353
IF ( ind1 <= ind ) THEN
1354
CALL Fatal( 'ComponentName', 'Syntax error in variable definition.' )
1355
END IF
1356
READ(BaseName(ind1+1:),'(i1)') DOFs
1357
DOFsTot = DOFsTot+DOFs
1358
IF ( DOFsTot>=Component ) EXIT
1359
ind = ind1+2
1360
END DO
1361
str = BaseName(ind+1:ind1-1)
1362
IF ( DOFs>1 ) THEN
1363
DOFs = Component - DOFsTot + DOFs
1364
str = str // ' ' //i2s(DOFs)
1365
END IF
1366
END IF
1367
!------------------------------------------------------------------------------
1368
END FUNCTION ComponentNameStr
1369
!------------------------------------------------------------------------------
1370
1371
1372
!------------------------------------------------------------------------------
1373
!> Solves a tridiagonal linear system.
1374
!------------------------------------------------------------------------------
1375
PURE SUBROUTINE SolveTriDiag( n, y, h, r )
1376
!------------------------------------------------------------------------------
1377
INTEGER, INTENT(in) :: n
1378
REAL(KIND=dp), INTENT(out) :: r(:)
1379
REAL(KIND=dp), INTENT(in) :: y(:), h(:)
1380
1381
REAL(KIND=dp) :: s,b(n)
1382
INTEGER :: i
1383
1384
DO i=2,n-1
1385
b(i) = 2 * ( h(i-1) + h(i) )
1386
r(i) = 3 * ( h(i) * ( y(i)-y(i-1) ) / h(i-1) + &
1387
h(i-1) * ( y(i+1)-y(i) ) / h(i) )
1388
END DO
1389
1390
r(2) = r(2) - h(2) * r(1)
1391
DO i=2,n-2
1392
s = -h(i+1) / b(i)
1393
r(i+1) = r(i+1) + s * r(i)
1394
b(i+1) = b(i+1) + s * h(i-1)
1395
END DO
1396
1397
DO i=n-1,2,-1
1398
r(i) = (r(i) - h(i-1) * r(i+1)) / b(i)
1399
END DO
1400
!------------------------------------------------------------------------------
1401
END SUBROUTINE SolveTriDiag
1402
!------------------------------------------------------------------------------
1403
1404
1405
FUNCTION CheckMonotone(n,x) RESULT ( Monotone )
1406
REAL(KIND=dp), INTENT(in) :: x(:)
1407
INTEGER, INTENT(in) :: n
1408
LOGICAL :: Monotone
1409
1410
INTEGER :: i
1411
1412
Monotone = .TRUE.
1413
DO i=1,n-1
1414
IF( x(i+1) <= x(i) ) THEN
1415
Monotone = .FALSE.
1416
WRITE (Message,'(E14.7,A,E14.7)') x(i),'>=',x(i+1)
1417
CALL WARN('CheckMonotone', Message)
1418
EXIT
1419
END IF
1420
END DO
1421
1422
END FUNCTION CheckMonotone
1423
1424
!------------------------------------------------------------------------------
1425
!> Solver for the coefficients of a cubic spline.
1426
!------------------------------------------------------------------------------
1427
PURE SUBROUTINE CubicSpline( n,x,y,r, monotone )
1428
!------------------------------------------------------------------------------
1429
REAL(KIND=dp), INTENT(in) :: x(:),y(:)
1430
REAL(KIND=dp), INTENT(out) :: r(:)
1431
INTEGER, INTENT(in) :: n
1432
LOGICAL, OPTIONAL, INTENT(in) :: monotone
1433
1434
REAL(KIND=dp) :: t,h(n),tau, alpha, beta
1435
INTEGER :: i
1436
LOGICAL :: mono
1437
1438
DO i=1,n-1
1439
h(i) = x(i+1) - x(i)
1440
END DO
1441
1442
r(1) = (y(2) - y(1) ) / h(1)
1443
r(n) = (y(n) - y(n-1) ) / h(n-1)
1444
1445
mono = .FALSE.
1446
IF(PRESENT(monotone)) mono = Monotone
1447
1448
IF (mono) THEN
1449
DO i=1,n-1
1450
h(i) = (y(i+1) - y(i) ) / h(i)
1451
END DO
1452
1453
DO i=2,n-1
1454
r(i) = (h(i-1) + h(i))/2
1455
END DO
1456
1457
DO i=1,n-1
1458
IF(ABS(h(i))<10*AEPS) THEN
1459
r(i) = 0._dp; r(i+1) = 0._dp
1460
CYCLE
1461
END IF
1462
1463
alpha = r(i) / h(i); beta = r(i+1) / h(i)
1464
IF ( alpha < 0._dp .OR. beta < 0._dp ) THEN
1465
r(i) = 0._dp;
1466
CYCLE
1467
END IF
1468
1469
tau = SQRT(alpha**2 + beta**2)
1470
IF(tau > 3) THEN
1471
tau = 3._dp / tau
1472
r(i) = alpha*tau*h(i); r(i+1)=beta*tau*h(i)
1473
END IF
1474
END DO
1475
ELSE
1476
CALL SolveTriDiag( n,y,h,r )
1477
END IF
1478
!------------------------------------------------------------------------------
1479
END SUBROUTINE CubicSpline
1480
!------------------------------------------------------------------------------
1481
1482
!------------------------------------------------------------------------------
1483
!> Evalulate a cubic spline.
1484
!------------------------------------------------------------------------------
1485
PURE FUNCTION CubicSplineVal(x,y,r,t) RESULT(s)
1486
!------------------------------------------------------------------------------
1487
REAL(KIND=dp) :: s
1488
REAL(KIND=dp), INTENT(in) :: x(:),y(:),r(:),t
1489
!------------------------------------------------------------------------------
1490
REAL(KIND=dp) :: a,b,c,d,h,lt
1491
1492
h = x(2)-x(1)
1493
a = -2 * ( y(2) - y(1) ) + ( r(1) + r(2) ) * h
1494
b = 3 * ( y(2) - y(1) ) - ( 2*r(1) + r(2) ) * h
1495
c = r(1) * h
1496
d = y(1)
1497
1498
lt = (t - x(1)) / h
1499
s = ((a*lt + b) * lt + c) * lt + d
1500
!------------------------------------------------------------------------------
1501
END FUNCTION CubicSplineVal
1502
!------------------------------------------------------------------------------
1503
1504
1505
!------------------------------------------------------------------------------
1506
!> Evalulate derivative of cubic spline.
1507
!------------------------------------------------------------------------------
1508
PURE FUNCTION CubicSplinedVal(x,y,r,t) RESULT(s)
1509
!------------------------------------------------------------------------------
1510
REAL(KIND=dp) :: s
1511
REAL(KIND=dp), INTENT(in) :: x(:),y(:),r(:),t
1512
!------------------------------------------------------------------------------
1513
REAL(KIND=dp) :: a,b,c,h,lt
1514
1515
h = x(2)-x(1)
1516
a = -2 * ( y(2) - y(1) ) + ( r(1) + r(2) ) * h
1517
b = 3 * ( y(2) - y(1) ) - ( 2*r(1) + r(2) ) * h
1518
c = r(1) * h
1519
1520
lt = (t - x(1)) / h
1521
s = ((3*a*lt + 2*b) * lt + c)/h
1522
!------------------------------------------------------------------------------
1523
END FUNCTION CubicSplinedVal
1524
!------------------------------------------------------------------------------
1525
1526
1527
!------------------------------------------------------------------------------
1528
!> Search array index such that tval(i) <= t < tval(i+1)
1529
!------------------------------------------------------------------------------
1530
PURE FUNCTION SearchInterval( tval, t ) RESULT(i)
1531
!------------------------------------------------------------------------------
1532
INTEGER :: i
1533
REAL(KIND=dp), INTENT(in) :: tval(:), t
1534
!------------------------------------------------------------------------------
1535
INTEGER :: n,n0,n1
1536
!------------------------------------------------------------------------------
1537
1538
n = SIZE(tval)
1539
1540
IF (t < tval(2)) THEN
1541
i = 1
1542
ELSE IF (t>=tval(n-1)) THEN
1543
i = n-1
1544
ELSE
1545
n0 = 1
1546
n1 = n
1547
i = (n0+n1)/2
1548
DO WHILE(.TRUE.)
1549
IF ( tval(i) <= t .AND. tval(i+1)>t ) EXIT
1550
1551
IF ( tval(i) > t ) THEN
1552
n1 = i-1
1553
ELSE
1554
n0 = i+1
1555
END IF
1556
i = (n0+n1)/2
1557
END DO
1558
END IF
1559
IF(i>n-1) i=n-1
1560
1561
!------------------------------------------------------------------------------
1562
END FUNCTION SearchInterval
1563
!------------------------------------------------------------------------------
1564
1565
!------------------------------------------------------------------------------
1566
!> As SearchInterval, but doesn't assume we'll find the value in the interval
1567
!------------------------------------------------------------------------------
1568
PURE FUNCTION SearchIntPosition( tval, t ) RESULT(i)
1569
!------------------------------------------------------------------------------
1570
INTEGER :: i
1571
INTEGER, INTENT(in) :: tval(:), t
1572
!------------------------------------------------------------------------------
1573
INTEGER :: n,n0,n1
1574
!------------------------------------------------------------------------------
1575
1576
n = SIZE(tval)
1577
1578
IF (t < tval(1)) THEN
1579
i = 0
1580
ELSE IF (t>=tval(n)) THEN
1581
i = n
1582
ELSE
1583
n0 = 1
1584
n1 = n
1585
i = (n0+n1)/2
1586
DO WHILE(.TRUE.)
1587
IF ( tval(i) <= t .AND. tval(i+1)>t ) EXIT
1588
1589
IF ( tval(i) > t ) THEN
1590
n1 = i-1
1591
ELSE
1592
n0 = i+1
1593
END IF
1594
i = (n0+n1)/2
1595
END DO
1596
END IF
1597
IF(i>n) i=n
1598
1599
!------------------------------------------------------------------------------
1600
END FUNCTION SearchIntPosition
1601
!------------------------------------------------------------------------------
1602
1603
1604
!------------------------------------------------------------------------------
1605
!> Interpolate values in a curve given by linear table or splines.
1606
!------------------------------------------------------------------------------
1607
PURE FUNCTION InterpolateCurve( TValues,FValues,T, CubicCoeff) RESULT( F )
1608
!------------------------------------------------------------------------------
1609
REAL(KIND=dp) :: F
1610
REAL(KIND=dp), INTENT(iN) :: TValues(:),FValues(:),T
1611
REAL(KIND=dp), OPTIONAL, POINTER, INTENT(in) :: CubicCoeff(:)
1612
!------------------------------------------------------------------------------
1613
INTEGER :: i,n
1614
LOGICAL :: Cubic
1615
!------------------------------------------------------------------------------
1616
1617
n = SIZE(TValues)
1618
1619
! This is a misuse of the interpolation in case of standard dependency
1620
! of type y=a*x.
1621
IF( n == 1 ) THEN
1622
F = FValues(1) * T
1623
RETURN
1624
END IF
1625
1626
i = SearchInterval( Tvalues, t )
1627
1628
Cubic = PRESENT(CubicCoeff)
1629
Cubic = Cubic .AND. T>=Tvalues(1) .AND. T<=Tvalues(n)
1630
IF ( Cubic ) Cubic = Cubic.AND.ASSOCIATED(CubicCoeff)
1631
1632
IF ( Cubic ) THEN
1633
F = CubicSplineVal(Tvalues(i:i+1),FValues(i:i+1),CubicCoeff(i:i+1),T)
1634
ELSE
1635
F = (T-TValues(i)) / (TValues(i+1)-TValues(i))
1636
F = (1-F)*FValues(i) + F*FValues(i+1)
1637
END IF
1638
END FUNCTION InterpolateCurve
1639
!------------------------------------------------------------------------------
1640
1641
1642
!------------------------------------------------------------------------------
1643
!> Derivate a curve given by linear table or splines.
1644
!------------------------------------------------------------------------------
1645
PURE FUNCTION DerivateCurve( TValues,FValues,T,CubicCoeff ) RESULT( F )
1646
!------------------------------------------------------------------------------
1647
REAL(KIND=dp) :: F
1648
REAL(KIND=dp), INTENT(in) :: TValues(:),FValues(:),T
1649
REAL(KIND=dp), OPTIONAL, POINTER, INTENT(in) :: CubicCoeff(:)
1650
!------------------------------------------------------------------------------
1651
INTEGER :: i,n
1652
LOGICAL :: Cubic
1653
!------------------------------------------------------------------------------
1654
n = SIZE(TValues)
1655
1656
i = SearchInterval( Tvalues, t )
1657
1658
Cubic = PRESENT(CubicCoeff)
1659
Cubic = Cubic .AND. T>=Tvalues(1) .AND. T<=Tvalues(n)
1660
IF ( Cubic ) Cubic = Cubic.AND.ASSOCIATED(CubicCoeff)
1661
1662
IF (Cubic) THEN
1663
F = CubicSplinedVal(Tvalues(i:i+1),FValues(i:i+1),CubicCoeff(i:i+1),T)
1664
ELSE
1665
F = (FValues(i+1)-FValues(i)) / (TValues(i+1)-TValues(i))
1666
END IF
1667
!------------------------------------------------------------------------------
1668
END FUNCTION DerivateCurve
1669
!------------------------------------------------------------------------------
1670
1671
1672
!------------------------------------------------------------------------------
1673
!> Integrate a curve given by linear table or splines.
1674
!------------------------------------------------------------------------------
1675
PURE SUBROUTINE CumulativeIntegral(TValues,FValues,CubicCoeff,Cumulative)
1676
!------------------------------------------------------------------------------
1677
REAL(KIND=dp), INTENT(in) :: TValues(:),FValues(:)
1678
REAL(KIND=dp), INTENT(out) :: Cumulative(:)
1679
REAL(KIND=dp), OPTIONAL, POINTER, INTENT(in) :: CubicCoeff(:)
1680
!------------------------------------------------------------------------------
1681
INTEGER :: i,n
1682
LOGICAL :: Cubic
1683
REAL(KIND=dp) :: t(2), y(2), r(2), h, a, b, c, d
1684
!------------------------------------------------------------------------------
1685
n = SIZE(TValues)
1686
1687
Cubic = PRESENT(CubicCoeff)
1688
IF ( Cubic ) Cubic = Cubic.AND.ASSOCIATED(CubicCoeff)
1689
1690
! here only complete intervals:
1691
! -----------------------------
1692
Cumulative(1) = 0._dp
1693
IF ( Cubic ) THEN
1694
DO i=1,n-1
1695
t(1) = Tvalues(i)
1696
t(2) = Tvalues(i+1)
1697
1698
y(1) = FValues(i)
1699
y(2) = FValues(i+1)
1700
1701
r(1) = CubicCoeff(i)
1702
r(2) = CubicCoeff(i+1)
1703
1704
h = t(2) - t(1)
1705
1706
a = (-2 * (y(2) - y(1)) + ( r(1) + r(2)) * h)/4
1707
b = ( 3 * (y(2) - y(1)) - (2*r(1) + r(2)) * h)/3
1708
c = (r(1) * h)/2
1709
d = y(1)
1710
Cumulative(i+1) = Cumulative(i) + h*(a+b+c+d)
1711
END DO
1712
ELSE
1713
DO i=1,n-1
1714
t(1) = Tvalues(i)
1715
t(2) = Tvalues(i+1)
1716
1717
y(1) = FValues(i)
1718
y(2) = FValues(i+1)
1719
1720
h = t(2) - t(1)
1721
c = (y(2)-y(1))/2
1722
d = y(1)
1723
Cumulative(i+1) = Cumulative(i) + h*(c+d)
1724
END DO
1725
END IF
1726
!------------------------------------------------------------------------------
1727
END SUBROUTINE CumulativeIntegral
1728
!------------------------------------------------------------------------------
1729
1730
!------------------------------------------------------------------------------
1731
!> Integrate a curve given by linear table or splines.
1732
!------------------------------------------------------------------------------
1733
FUNCTION IntegrateCurve(TValues,FValues,CubicCoeff,T0,T1,Cumulative,Found) RESULT(sumf)
1734
!------------------------------------------------------------------------------
1735
REAL(KIND=dp) :: sumf
1736
1737
REAL(KIND=dp), INTENT(in) :: TValues(:),FValues(:)
1738
REAL(KIND=dp), OPTIONAL, INTENT(in) :: T0, T1
1739
REAL(KIND=dp), OPTIONAL, INTENT(in) :: Cumulative(:)
1740
REAL(KIND=dp), OPTIONAL, POINTER, INTENT(in) :: CubicCoeff(:)
1741
LOGICAL, OPTIONAL :: Found
1742
!------------------------------------------------------------------------------
1743
INTEGER :: i,n,i0,i1
1744
LOGICAL :: Cubic
1745
REAL(KIND=dp) :: t(2), y(2), r(2), h, a, b, c, d, s0, s1, tt0, tt1
1746
!------------------------------------------------------------------------------
1747
1748
sumf = 0.0_dp
1749
IF(PRESENT(Found)) Found = .FALSE.
1750
1751
n = SIZE(TValues)
1752
IF(n<2) RETURN
1753
1754
IF(SIZE(FValues) /= n ) THEN
1755
CALL Warn('IntegrateCurve','TValues and Fvalues should be of same size!')
1756
RETURN
1757
END IF
1758
1759
tt0 = TValues(1)
1760
IF(PRESENT(t0)) tt0=t0
1761
1762
tt1 = TValues(n)
1763
IF(PRESENT(t1)) tt1=t1
1764
1765
IF(tt0>=tt1) RETURN
1766
1767
IF(PRESENT(Found)) Found = .TRUE.
1768
1769
1770
! t0 < first, t1 <= first
1771
IF(tt1<=Tvalues(1)) THEN
1772
t(1) = Tvalues(1)
1773
t(2) = Tvalues(2)
1774
1775
y(1) = FValues(1)
1776
y(2) = FValues(2)
1777
1778
h = t(2) - t(1)
1779
s0 = (tt0 - t(1)) / h
1780
s1 = (tt1 - t(1)) / h
1781
c = (y(2) - y(1)) / 2
1782
d = y(1)
1783
sumf = sumf + h * ((c*s1 + d)*s1 - (c*s0 + d)*s0)
1784
RETURN
1785
END IF
1786
1787
! t0 >= last, t1 > last
1788
IF(tt0>=Tvalues(n)) THEN
1789
t(1) = Tvalues(n-1)
1790
t(2) = Tvalues(n)
1791
1792
y(1) = FValues(n-1)
1793
y(2) = FValues(n)
1794
1795
h = t(2) - t(1)
1796
s0 = (tt0 - t(1)) / h
1797
s1 = (tt1 - t(1)) / h
1798
c = (y(2) - y(1)) / 2
1799
d = y(1)
1800
sumf = sumf + h * ((c*s1 + d)*s1 - (c*s0 + d)*s0)
1801
RETURN
1802
END IF
1803
1804
! first interval outside
1805
IF(tt0<Tvalues(1)) THEN
1806
t(1) = Tvalues(1)
1807
t(2) = Tvalues(2)
1808
1809
y(1) = FValues(1)
1810
y(2) = FValues(2)
1811
1812
h = t(2) - t(1)
1813
s0 = (tt0 - t(1)) / h
1814
c = (y(2) - y(1)) / 2
1815
d = y(1)
1816
sumf = sumf - h * (c*s0 + d)*s0
1817
tt0 = Tvalues(1)
1818
END IF
1819
1820
! last interval outside
1821
IF(tt1>Tvalues(n)) THEN
1822
t(1) = Tvalues(n-1)
1823
t(2) = Tvalues(n)
1824
1825
y(1) = FValues(n-1)
1826
y(2) = FValues(n)
1827
1828
h = t(2) - t(1)
1829
s1 = (tt1 - t(1)) / h
1830
c = (y(2) - y(1)) / 2
1831
d = y(1)
1832
sumf = sumf + h * ( (c*s1 + d)*s1 - (c+d) )
1833
tt1 = Tvalues(n)
1834
END IF
1835
1836
IF(tt0 >= tt1) RETURN
1837
1838
Cubic = PRESENT(CubicCoeff)
1839
IF ( Cubic ) Cubic = Cubic.AND.ASSOCIATED(CubicCoeff)
1840
1841
i0 = SearchInterval( Tvalues, tt0 )
1842
1843
! first (possibly partial) interval:
1844
! -------------------------------------
1845
t(1) = Tvalues(i0)
1846
t(2) = Tvalues(i0+1)
1847
1848
h = t(2) - t(1)
1849
s0 = (tt0-t(1))/h
1850
s1 = MIN((tt1-t(1))/h,1._dp)
1851
1852
IF(s0>0 .OR. s1<1) THEN
1853
y(1) = FValues(i0)
1854
y(2) = FValues(i0+1)
1855
1856
IF(Cubic) THEN
1857
r(1) = CubicCoeff(i0)
1858
r(2) = CubicCoeff(i0+1)
1859
1860
a = (-2 * (y(2) - y(1)) + ( r(1) + r(2)) * h)/4
1861
b = ( 3 * (y(2) - y(1)) - (2*r(1) + r(2)) * h)/3
1862
c = (r(1) * h)/2
1863
d = y(1)
1864
sumf = sumf + h * ( (((a*s1 + b)*s1 + c)*s1 + d)*s1 - &
1865
(((a*s0 + b)*s0 + c)*s0 + d)*s0 )
1866
1867
ELSE
1868
c = (y(2)-y(1))/2
1869
d = y(1)
1870
sumf = sumf + h * ( (c*s1 + d)*s1 - (c*s0 + d)*s0 )
1871
END IF
1872
i0 = i0 + 1
1873
tt0 = Tvalues(i0)
1874
IF(tt0 >= tt1) RETURN
1875
END IF
1876
1877
i1 = SearchInterval( Tvalues, tt1 )
1878
1879
! last (possibly partial) interval:
1880
! ------------------------------------
1881
t(1) = Tvalues(i1)
1882
t(2) = Tvalues(i1+1)
1883
1884
h = t(2) - t(1)
1885
1886
s0 = MAX((tt0-t(1))/h, 0.0_dp)
1887
s1 = (tt1-t(1))/h
1888
1889
IF(s0>0 .OR. s1<1) THEN
1890
y(1) = FValues(i1)
1891
y(2) = FValues(i1+1)
1892
1893
IF(Cubic) THEN
1894
r(1) = CubicCoeff(i1)
1895
r(2) = CubicCoeff(i1+1)
1896
1897
a = (-2 * (y(2) - y(1)) + ( r(1) + r(2)) * h)/4
1898
b = ( 3 * (y(2) - y(1)) - (2*r(1) + r(2)) * h)/3
1899
c = (r(1) * h)/2
1900
d = y(1)
1901
sumf = sumf + h * ( (((a*s1 + b)*s1 + c)*s1 + d)*s1 - &
1902
(((a*s0 + b)*s0 + c)*s0 + d)*s0 )
1903
ELSE
1904
c = (y(2)-y(1))/2
1905
d = y(1)
1906
sumf = sumf + h * ( (c*s1 + d)*s1 - (c*s0 + d)*s0 )
1907
END IF
1908
i1 = i1 - 1
1909
tt1 = Tvalues(i1+1)
1910
IF(tt0 >= tt1) RETURN
1911
END IF
1912
1913
! here only complete intervals:
1914
! -----------------------------
1915
1916
IF(PRESENT(Cumulative)) THEN
1917
sumf = sumf + Cumulative(i1+1) - Cumulative(i0)
1918
RETURN
1919
END IF
1920
1921
IF ( Cubic ) THEN
1922
DO i=i0,i1
1923
t(1) = Tvalues(i)
1924
t(2) = Tvalues(i+1)
1925
1926
y(1) = FValues(i)
1927
y(2) = FValues(i+1)
1928
1929
r(1) = CubicCoeff(i)
1930
r(2) = CubicCoeff(i+1)
1931
1932
h = t(2) - t(1)
1933
1934
a = (-2 * (y(2) - y(1)) + ( r(1) + r(2)) * h)/4
1935
b = ( 3 * (y(2) - y(1)) - (2*r(1) + r(2)) * h)/3
1936
c = (r(1) * h)/2
1937
d = y(1)
1938
sumf = sumf + h * (a+b+c+d)
1939
END DO
1940
ELSE
1941
DO i=i0,i1
1942
t(1) = Tvalues(i)
1943
t(2) = Tvalues(i+1)
1944
1945
y(1) = FValues(i)
1946
y(2) = FValues(i+1)
1947
1948
h = t(2) - t(1)
1949
c = (y(2)-y(1))/2
1950
d = y(1)
1951
sumf = sumf + h * (c+d)
1952
END DO
1953
END IF
1954
!------------------------------------------------------------------------------
1955
END FUNCTION IntegrateCurve
1956
!------------------------------------------------------------------------------
1957
1958
1959
!------------------------------------------------------------------------------
1960
SUBROUTINE ClearMatrix( Matrix )
1961
#if defined(ELMER_HAVE_MPI_MODULE)
1962
USE mpi
1963
#endif
1964
TYPE(Matrix_t), POINTER, INTENT(in) :: Matrix
1965
#if defined(ELMER_HAVE_MPIF_HEADER)
1966
INCLUDE "mpif.h"
1967
#endif
1968
1969
Matrix % FORMAT = MATRIX_CRS
1970
1971
NULLIFY( Matrix % Child )
1972
NULLIFY( Matrix % Parent )
1973
NULLIFY( Matrix % EMatrix )
1974
NULLIFY( Matrix % ConstraintMatrix )
1975
1976
NULLIFY( Matrix % Perm )
1977
NULLIFY( Matrix % InvPerm )
1978
1979
NULLIFY( Matrix % Cols )
1980
NULLIFY( Matrix % Rows )
1981
NULLIFY( Matrix % Diag )
1982
1983
NULLIFY( Matrix % RHS )
1984
NULLIFY( Matrix % Force )
1985
NULLIFY( Matrix % RHS_im )
1986
1987
NULLIFY( Matrix % Values )
1988
NULLIFY( Matrix % ILUValues )
1989
NULLIFY( Matrix % MassValues )
1990
NULLIFY( Matrix % DampValues )
1991
1992
NULLIFY( Matrix % BulkRHS )
1993
NULLIFY( Matrix % BulkValues )
1994
1995
NULLIFY( Matrix % BulkResidual )
1996
1997
NULLIFY( Matrix % ILUCols )
1998
NULLIFY( Matrix % ILURows )
1999
NULLIFY( Matrix % ILUDiag )
2000
2001
NULLIFY( Matrix % CRHS )
2002
NULLIFY( Matrix % CForce )
2003
2004
NULLIFY( Matrix % ParMatrix )
2005
2006
NULLIFY( Matrix % CValues )
2007
NULLIFY( Matrix % CILUValues )
2008
NULLIFY( Matrix % CMassValues )
2009
NULLIFY( Matrix % CDampValues )
2010
2011
! NULLIFY( Matrix % GRows )
2012
! NULLIFY( Matrix % RowOwner )
2013
NULLIFY( Matrix % GOrder )
2014
NULLIFY( Matrix % EPerm )
2015
2016
NULLIFY( Matrix % ParMatrix )
2017
2018
NULLIFY( Matrix % ParallelInfo )
2019
#ifdef HAVE_UMFPACK
2020
Matrix % UMFPack_Numeric = 0
2021
#endif
2022
2023
Matrix % Cholesky = .FALSE.
2024
Matrix % Lumped = .FALSE.
2025
Matrix % Ordered = .FALSE.
2026
Matrix % COMPLEX = .FALSE.
2027
Matrix % Symmetric = .FALSE.
2028
Matrix % SolveCount = 0
2029
Matrix % NumberOfRows = 0
2030
Matrix % Ndeg = -1
2031
Matrix % ProjectorBC = 0
2032
Matrix % ProjectorType = PROJECTOR_TYPE_DEFAULT
2033
2034
Matrix % Solver => NULL()
2035
2036
Matrix % DGMatrix = .FALSE.
2037
Matrix % Comm = ELMER_COMM_WORLD
2038
2039
END SUBROUTINE ClearMatrix
2040
2041
2042
!------------------------------------------------------------------------------
2043
FUNCTION AllocateMatrix() RESULT(Matrix)
2044
!------------------------------------------------------------------------------
2045
TYPE(Matrix_t), POINTER :: Matrix
2046
!------------------------------------------------------------------------------
2047
ALLOCATE( Matrix )
2048
CALL ClearMatrix( Matrix )
2049
!------------------------------------------------------------------------------
2050
END FUNCTION AllocateMatrix
2051
!------------------------------------------------------------------------------
2052
2053
2054
!------------------------------------------------------------------------------
2055
RECURSIVE SUBROUTINE FreeQuadrantTree( Root )
2056
!------------------------------------------------------------------------------
2057
TYPE(Quadrant_t), POINTER :: Root
2058
2059
INTEGER :: i
2060
2061
IF ( .NOT. ASSOCIATED( Root ) ) RETURN
2062
2063
IF ( ASSOCIATED(Root % Elements) ) DEALLOCATE( Root % Elements )
2064
2065
IF ( ASSOCIATED( Root % ChildQuadrants ) ) THEN
2066
DO i=1,SIZE(Root % ChildQuadrants)
2067
CALL FreeQuadrantTree( Root % ChildQuadrants(i) % Quadrant )
2068
END DO
2069
DEALLOCATE( Root % ChildQuadrants )
2070
NULLIFY( Root % ChildQuadrants )
2071
END IF
2072
2073
DEALLOCATE( Root )
2074
!------------------------------------------------------------------------------
2075
END SUBROUTINE FreeQuadrantTree
2076
!------------------------------------------------------------------------------
2077
2078
2079
!------------------------------------------------------------------------------
2080
SUBROUTINE AllocateRealVector( F, n, From, FailureMessage )
2081
!------------------------------------------------------------------------------
2082
REAL(KIND=dp), POINTER :: F(:)
2083
INTEGER :: n
2084
CHARACTER(LEN=*), OPTIONAL :: From, FailureMessage
2085
!------------------------------------------------------------------------------
2086
INTEGER :: istat
2087
!------------------------------------------------------------------------------
2088
2089
ALLOCATE( F(n), STAT=istat )
2090
IF ( istat /= 0 ) THEN
2091
WRITE( Message, * )'Unable to allocate ', n, ' element real array.'
2092
CALL Error( 'AllocateRealVector', Message )
2093
IF ( PRESENT( From ) ) THEN
2094
WRITE( Message, * )'Requested From: ', TRIM(From)
2095
CALL Error( 'AllocateRealVector', Message )
2096
END IF
2097
IF ( PRESENT( FailureMessage ) ) THEN
2098
CALL Fatal( 'AllocateRealVector', FailureMessage )
2099
END IF
2100
END IF
2101
!------------------------------------------------------------------------------
2102
END SUBROUTINE AllocateRealVector
2103
!------------------------------------------------------------------------------
2104
2105
2106
!------------------------------------------------------------------------------
2107
SUBROUTINE AllocateComplexVector( f, n, From, FailureMessage )
2108
!------------------------------------------------------------------------------
2109
COMPLEX(KIND=dp), POINTER :: f(:)
2110
INTEGER :: n
2111
CHARACTER(LEN=*), OPTIONAL :: From, FailureMessage
2112
!------------------------------------------------------------------------------
2113
INTEGER :: istat
2114
!------------------------------------------------------------------------------
2115
2116
ALLOCATE( f(n), STAT=istat )
2117
IF ( istat /= 0 ) THEN
2118
WRITE( Message, * )'Unable to allocate ', n, ' element complex array.'
2119
CALL Error( 'AllocateComplexVector', Message )
2120
IF ( PRESENT( From ) ) THEN
2121
WRITE( Message, * )'Requested From: ', TRIM(From)
2122
CALL Error( 'AllocateComplexVector', Message )
2123
END IF
2124
IF ( PRESENT( FailureMessage ) ) THEN
2125
CALL Fatal( 'AllocateComplexVector', FailureMessage )
2126
END IF
2127
END IF
2128
!------------------------------------------------------------------------------
2129
END SUBROUTINE AllocateComplexVector
2130
!------------------------------------------------------------------------------
2131
2132
2133
!------------------------------------------------------------------------------
2134
SUBROUTINE AllocateIntegerVector( f, n, From, FailureMessage )
2135
!------------------------------------------------------------------------------
2136
INTEGER, POINTER :: f(:)
2137
INTEGER :: n
2138
CHARACTER(LEN=*), OPTIONAL :: From, FailureMessage
2139
!------------------------------------------------------------------------------
2140
INTEGER :: istat
2141
!------------------------------------------------------------------------------
2142
2143
ALLOCATE( f(n), STAT=istat )
2144
IF ( istat /= 0 ) THEN
2145
WRITE( Message, * )'Unable to allocate ', n, ' element integer array.'
2146
CALL Error( 'AllocateIntegerVector', Message )
2147
IF ( PRESENT( From ) ) THEN
2148
WRITE( Message, * )'Requested From: ', TRIM(From)
2149
CALL Error( 'AllocateIntegerVector', Message )
2150
END IF
2151
IF ( PRESENT( FailureMessage ) ) THEN
2152
CALL Fatal( 'AllocateIntegerVector', FailureMessage )
2153
END IF
2154
END IF
2155
!------------------------------------------------------------------------------
2156
END SUBROUTINE AllocateIntegerVector
2157
!------------------------------------------------------------------------------
2158
2159
2160
!------------------------------------------------------------------------------
2161
SUBROUTINE AllocateLogicalVector( f, n, From, FailureMessage )
2162
!------------------------------------------------------------------------------
2163
LOGICAL, POINTER :: f(:)
2164
INTEGER :: n
2165
CHARACTER(LEN=*), OPTIONAL :: From, FailureMessage
2166
!------------------------------------------------------------------------------
2167
INTEGER :: istat
2168
!------------------------------------------------------------------------------
2169
2170
ALLOCATE( f(n), STAT=istat )
2171
IF ( istat /= 0 ) THEN
2172
WRITE( Message, * )'Unable to allocate ', n, ' element logical array.'
2173
CALL Error( 'AllocateLogicalVector', Message )
2174
IF ( PRESENT( From ) ) THEN
2175
WRITE( Message, * )'Requested From: ', TRIM(From)
2176
CALL Error( 'AllocateLogicalVector', Message )
2177
END IF
2178
IF ( PRESENT( FailureMessage ) ) THEN
2179
CALL Fatal( 'AllocateLogicalVector', FailureMessage )
2180
END IF
2181
END IF
2182
!------------------------------------------------------------------------------
2183
END SUBROUTINE AllocateLogicalVector
2184
!------------------------------------------------------------------------------
2185
2186
2187
!------------------------------------------------------------------------------
2188
SUBROUTINE AllocateElementVector( f, n, From, FailureMessage )
2189
!------------------------------------------------------------------------------
2190
TYPE(Element_t), POINTER :: f(:)
2191
INTEGER :: n
2192
CHARACTER(LEN=*), OPTIONAL :: From, FailureMessage
2193
!------------------------------------------------------------------------------
2194
INTEGER :: istat
2195
!------------------------------------------------------------------------------
2196
2197
ALLOCATE( f(n), STAT=istat )
2198
IF ( istat /= 0 ) THEN
2199
WRITE( Message, * )'Unable to allocate structured ', n, ' element array.'
2200
CALL Error( 'AllocateElementVector', Message )
2201
IF ( PRESENT( From ) ) THEN
2202
WRITE( Message, * )'Requested From: ', TRIM(From)
2203
CALL Error( 'AllocateElementVector', Message )
2204
END IF
2205
IF ( PRESENT( FailureMessage ) ) THEN
2206
CALL Fatal( 'AllocateElementVector', FailureMessage )
2207
END IF
2208
END IF
2209
!------------------------------------------------------------------------------
2210
END SUBROUTINE AllocateElementVector
2211
!------------------------------------------------------------------------------
2212
2213
2214
!------------------------------------------------------------------------------
2215
SUBROUTINE AllocateRealArray( f, n1, n2, From, FailureMessage )
2216
!------------------------------------------------------------------------------
2217
REAL(KIND=dp), POINTER :: f(:,:)
2218
INTEGER :: n1,n2
2219
CHARACTER(LEN=*), OPTIONAL :: From, FailureMessage
2220
!------------------------------------------------------------------------------
2221
INTEGER :: istat
2222
!------------------------------------------------------------------------------
2223
2224
ALLOCATE( f(n1,n2), STAT=istat )
2225
IF ( istat /= 0 ) THEN
2226
WRITE( Message, * )'Unable to allocate ', n1, ' by ', n2, ' element real matrix.'
2227
CALL Error( 'AllocateRealArray', Message )
2228
IF ( PRESENT( From ) ) THEN
2229
WRITE( Message, * )'Requested From: ', TRIM(From)
2230
CALL Error( 'AllocateRealArray', Message )
2231
END IF
2232
IF ( PRESENT( FailureMessage ) ) THEN
2233
CALL Fatal( 'AllocateRealArray', FailureMessage )
2234
END IF
2235
END IF
2236
!------------------------------------------------------------------------------
2237
END SUBROUTINE AllocateRealArray
2238
!------------------------------------------------------------------------------
2239
2240
!------------------------------------------------------------------------------
2241
SUBROUTINE AllocateComplexArray( f, n1, n2, From, FailureMessage )
2242
!------------------------------------------------------------------------------
2243
COMPLEX(KIND=dp), POINTER :: f(:,:)
2244
INTEGER :: n1,n2
2245
CHARACTER(LEN=*), OPTIONAL :: From, FailureMessage
2246
!------------------------------------------------------------------------------
2247
INTEGER :: istat
2248
!------------------------------------------------------------------------------
2249
2250
ALLOCATE( f(n1,n2), STAT=istat )
2251
IF ( istat /= 0 ) THEN
2252
WRITE( Message, * )'Unable to allocate ', n1, ' by ', n2, ' element complex matrix.'
2253
CALL Error( 'AllocateComplexArray', Message )
2254
IF ( PRESENT( From ) ) THEN
2255
WRITE( Message, * )'Requested From: ', TRIM(From)
2256
CALL Error( 'AllocateComplexArray', Message )
2257
END IF
2258
IF ( PRESENT( FailureMessage ) ) THEN
2259
CALL Fatal( 'AllocateComplexArray', FailureMessage )
2260
END IF
2261
END IF
2262
!------------------------------------------------------------------------------
2263
END SUBROUTINE AllocateComplexArray
2264
!------------------------------------------------------------------------------
2265
2266
2267
!------------------------------------------------------------------------------
2268
SUBROUTINE AllocateIntegerArray( f, n1, n2, From, FailureMessage )
2269
!------------------------------------------------------------------------------
2270
INTEGER, POINTER :: f(:,:)
2271
INTEGER :: n1,n2
2272
CHARACTER(LEN=*), OPTIONAL :: From, FailureMessage
2273
!------------------------------------------------------------------------------
2274
INTEGER :: istat
2275
!------------------------------------------------------------------------------
2276
2277
istat = -1
2278
IF ( n1 > 0 .AND. n2 > 0 ) THEN
2279
ALLOCATE( f(n1,n2), STAT=istat )
2280
END IF
2281
IF ( istat /= 0 ) THEN
2282
WRITE( Message, * )'Unable to allocate ', n1, ' by ', n2, ' element integer matrix.'
2283
CALL Error( 'AllocateIntegerArray', Message )
2284
IF ( PRESENT( From ) ) THEN
2285
WRITE( Message, * )'Requested From: ', TRIM(From)
2286
CALL Error( 'AllocateIntegerArray', Message )
2287
END IF
2288
IF ( PRESENT( FailureMessage ) ) THEN
2289
CALL Fatal( 'AllocateIntegerArray', FailureMessage )
2290
END IF
2291
END IF
2292
!------------------------------------------------------------------------------
2293
END SUBROUTINE AllocateIntegerArray
2294
!------------------------------------------------------------------------------
2295
2296
2297
2298
!------------------------------------------------------------------------------
2299
SUBROUTINE AllocateLogicalArray( f, n1, n2, From, FailureMessage )
2300
!------------------------------------------------------------------------------
2301
LOGICAL, POINTER :: f(:,:)
2302
INTEGER :: n1,n2
2303
CHARACTER(LEN=*), OPTIONAL :: From, FailureMessage
2304
!------------------------------------------------------------------------------
2305
INTEGER :: istat
2306
!------------------------------------------------------------------------------
2307
2308
istat = -1
2309
IF ( n1 > 0 .AND. n2 > 0 ) THEN
2310
ALLOCATE( f(n1,n2), STAT=istat )
2311
END IF
2312
IF ( istat /= 0 ) THEN
2313
WRITE( Message, * )'Unable to allocate ', n1, ' by ', n2, ' element logical matrix.'
2314
CALL Error( 'AllocateLogicalArray', Message )
2315
IF ( PRESENT( From ) ) THEN
2316
WRITE( Message, * )'Requested From: ', TRIM(From)
2317
CALL Error( 'AllocateLogicalArray', Message )
2318
END IF
2319
IF ( PRESENT( FailureMessage ) ) THEN
2320
CALL Fatal( 'AllocateLogicalArray', FailureMessage )
2321
END IF
2322
END IF
2323
!------------------------------------------------------------------------------
2324
END SUBROUTINE AllocateLogicalArray
2325
!------------------------------------------------------------------------------
2326
2327
! Pad given integer value to be the next largest multiple of nbyte
2328
FUNCTION IntegerNBytePad(val, nbyte) RESULT(padval)
2329
IMPLICIT NONE
2330
2331
INTEGER, INTENT(IN) :: val, nbyte
2332
INTEGER :: padval
2333
! Parameters and variables
2334
INTEGER, PARAMETER :: bytesinint = KIND(val)
2335
INTEGER :: nbytesinint
2336
2337
! Compute number of nbytes in int
2338
nbytesinint = nbyte/bytesinint
2339
! Compute value padded to multiples of n-byte
2340
padval=((val-1)/nbytesinint)*nbytesinint+nbytesinint
2341
END FUNCTION IntegerNBytePad
2342
2343
! Pad given value to be the next largest multiple of nbyte
2344
FUNCTION NBytePad(val, bytesinelem, nbyte) RESULT(padval)
2345
IMPLICIT NONE
2346
2347
INTEGER, INTENT(IN) :: val, bytesinelem, nbyte
2348
INTEGER :: padval
2349
! Variables
2350
INTEGER :: nbytesinelem
2351
2352
! Compute number of nbytes in a single element
2353
nbytesinelem = nbyte/bytesinelem
2354
! Compute value padded to multiples of n-byte
2355
padval=((val-1)/nbytesinelem)*nbytesinelem+nbytesinelem
2356
END FUNCTION NBytePad
2357
2358
!------------------------------------------------------------------------------
2359
!> Given the filename0 (and suffix0) find the 1st free filename
2360
!> that does not exist in the current working directory
2361
!------------------------------------------------------------------------------
2362
2363
FUNCTION NextFreeFilename(Filename0,Suffix0,LastExisting) RESULT (Filename)
2364
2365
CHARACTER(LEN=*) :: Filename0
2366
CHARACTER(LEN=*), OPTIONAL :: Suffix0
2367
LOGICAL, OPTIONAL :: LastExisting
2368
CHARACTER(:), ALLOCATABLE :: Filename,Prefix,Suffix,PrevFilename
2369
LOGICAL :: FileIs
2370
INTEGER :: No, ind, len
2371
2372
ind = INDEX( FileName0,'.',.TRUE. )
2373
len = LEN_TRIM(Filename0)
2374
IF(ind > 0) THEN
2375
Prefix = Filename0(1:ind-1)
2376
Suffix = Filename0(ind:len)
2377
ELSE
2378
Prefix = Filename0(1:len)
2379
IF(PRESENT(Suffix0)) THEN
2380
Suffix = '.'//TRIM(Suffix0)
2381
ELSE
2382
Suffix = '.dat'
2383
END IF
2384
END IF
2385
2386
DO No = 1,9999
2387
IF( No > 0 ) PrevFilename = Filename
2388
FileName = TRIM(Prefix)//I2S(No)//TRIM(Suffix)
2389
2390
INQUIRE( FILE=Filename, EXIST=FileIs )
2391
IF(.NOT. FileIs) EXIT
2392
END DO
2393
2394
IF( PRESENT(LastExisting)) THEN
2395
IF( LastExisting ) Filename = PrevFilename
2396
END IF
2397
2398
CALL Info('NextFreeFilename','Next Free filename is: '//TRIM(Filename),Level=12)
2399
2400
!------------------------------------------------------------------------------
2401
END FUNCTION NextFreeFilename
2402
!------------------------------------------------------------------------------
2403
2404
!------------------------------------------------------------------------------
2405
!> Given the filename0 add a string related to the partitioning.
2406
!------------------------------------------------------------------------------
2407
2408
FUNCTION AddFilenameParSuffix(Filename0,Suffix0,Parallel,MyPe,NumWidth,PeMax,PeSeparator) RESULT (Filename)
2409
2410
CHARACTER(LEN=*) :: Filename0
2411
CHARACTER(LEN=*), OPTIONAL :: Suffix0
2412
CHARACTER(LEN=*), OPTIONAL :: PeSeparator
2413
LOGICAL :: Parallel
2414
INTEGER :: MyPe
2415
INTEGER, OPTIONAL :: NumWidth
2416
INTEGER, OPTIONAL :: PeMax
2417
CHARACTER(LEN=MAX_NAME_LEN) :: Filename
2418
!------------------------------------------------------------------------------
2419
CHARACTER(LEN=MAX_NAME_LEN) :: OutStyle
2420
CHARACTER(:), ALLOCATABLE :: Prefix, Suffix
2421
INTEGER :: No, ind, len, NumW, NoLim
2422
2423
ind = INDEX( FileName0,'.',.TRUE. )
2424
len = LEN_TRIM(Filename0)
2425
2426
! If the only dot is the first one it only related to the current working directory.
2427
IF(ind > 1) THEN
2428
Prefix = Filename0(1:ind-1)
2429
Suffix = Filename0(ind:len)
2430
ELSE
2431
Prefix = Filename0(1:len)
2432
IF(PRESENT(Suffix0)) THEN
2433
Suffix = '.'//TRIM(Suffix0)
2434
ELSE
2435
Suffix = '.dat'
2436
END IF
2437
END IF
2438
2439
IF( Parallel ) THEN
2440
No = MyPe + 1
2441
2442
IF( PRESENT(NumWidth) ) THEN
2443
NumW = NumWidth
2444
ELSE IF( PRESENT( PeMax ) ) THEN
2445
NumW = CEILING( LOG10( 1.0_dp * PeMax ) )
2446
ELSE
2447
NumW = 4
2448
END IF
2449
NoLim = 10**NumW
2450
2451
IF( PRESENT( PeSeparator ) ) THEN
2452
Prefix = TRIM(Prefix)//TRIM(PeSeparator)
2453
END IF
2454
2455
IF( No >= NoLim ) THEN
2456
FileName = TRIM(Prefix)//I2S(No)//TRIM(Suffix)
2457
ELSE
2458
WRITE( OutStyle,'(A,I1,A,I1,A)') '(A,I',NumW,'.',NumW,',A)'
2459
WRITE( FileName,OutStyle) TRIM(Prefix),No,TRIM(Suffix)
2460
END IF
2461
ELSE
2462
FileName = TRIM(Prefix)//TRIM(Suffix)
2463
END IF
2464
2465
!------------------------------------------------------------------------------
2466
END FUNCTION AddFilenameParSuffix
2467
!------------------------------------------------------------------------------
2468
2469
2470
! This takes union of two integer vectors
2471
! and returns the number of common values.
2472
!---------------------------------------------
2473
FUNCTION CountSameIntegers(v1,v2,vsame) RESULT ( n )
2474
INTEGER, POINTER :: v1(:), v2(:)
2475
INTEGER, POINTER, OPTIONAL :: vsame(:)
2476
INTEGER :: n
2477
2478
INTEGER :: i1,i2
2479
2480
n = 0
2481
IF(.NOT. ASSOCIATED(v1)) RETURN
2482
IF(.NOT. ASSOCIATED(v2)) RETURN
2483
2484
DO i1=1,SIZE(v1)
2485
DO i2=1,SIZE(v2)
2486
IF( v1(i1) == v2(i2) ) n = n+1
2487
END DO
2488
END DO
2489
2490
IF(n==0) RETURN
2491
2492
IF( PRESENT(vsame) ) THEN
2493
IF(.NOT. ASSOCIATED(vsame) ) THEN
2494
ALLOCATE(vsame(n) )
2495
END IF
2496
vsame = 0
2497
n = 0
2498
2499
DO i1=1,SIZE(v1)
2500
DO i2=1,SIZE(v2)
2501
IF( v1(i1) == v2(i2) ) THEN
2502
n = n+1
2503
vsame(n) = v1(i1)
2504
END IF
2505
END DO
2506
END DO
2507
END IF
2508
2509
END FUNCTION CountSameIntegers
2510
2511
2512
2513
!---------------------------------------------------------<
2514
!> Returns values from a normal distribution to be used in
2515
!> thermal velocity distribution, for example.
2516
!---------------------------------------------------------
2517
FUNCTION NormalRandom() RESULT ( normalrand )
2518
2519
REAL(KIND=dp) :: normalrand,mean
2520
INTEGER :: flag = 0
2521
REAL(KIND=dp) :: fac,gsave,rsq,r1,r2
2522
2523
SAVE flag,gsave
2524
2525
IF (flag == 0) THEN
2526
rsq=2.0_dp
2527
2528
DO WHILE(rsq >= 1.0_dp .OR. rsq == 0.0_dp )
2529
CALL RANDOM_NUMBER(r1)
2530
CALL RANDOM_NUMBER(r2)
2531
r1 = 2.0_dp * r1 - 1.0_dp
2532
r2 = 2.0_dp * r2 - 1.0_dp
2533
rsq = r1*r1 + r2*r2
2534
ENDDO
2535
2536
fac = SQRT(-2.0_dp * LOG(rsq) / rsq)
2537
gsave = r1 * fac
2538
normalrand = r2 * fac
2539
flag = 1
2540
ELSE
2541
normalrand = gsave
2542
flag = 0
2543
ENDIF
2544
2545
END FUNCTION NormalRandom
2546
2547
2548
!---------------------------------------------------------
2549
!> Returns values from a even distribution [0,1]
2550
!---------------------------------------------------------
2551
FUNCTION EvenRandom() RESULT ( rand )
2552
REAL(KIND=dp) :: rand
2553
CALL RANDOM_NUMBER(rand)
2554
END FUNCTION EvenRandom
2555
2556
2557
!-----------------------------------------------------
2558
! Convert to effective BH-curve for harmonic analysis
2559
!-----------------------------------------------------
2560
SUBROUTINE ConvertTableToHarmonic(n,bVal,hVal)
2561
INTEGER :: n
2562
REAL(KIND=dp) :: bval(:), hval(:)
2563
2564
INTEGER :: i,j,n_int = 200
2565
REAL(KIND=dp) :: alpha,b,h,nu_eff, hOrig(n)
2566
2567
hOrig = hVal(1:n)
2568
DO i=1,n
2569
nu_eff = 0._dp
2570
DO j=1,n_int
2571
alpha = PI/2._dp*j/(1._dp*N_int)
2572
b = sin(alpha)*bVal(i)
2573
h = linterpolate(n,b,bVal,hOrig)
2574
IF(b>0.0_dp) nu_eff = nu_eff + h/b
2575
END DO
2576
hVal(i) = nu_eff * bVal(i) / N_int
2577
END DO
2578
2579
CONTAINS
2580
2581
! Evaluates y=f(x) for a piecewise linear function f defined by x_points and y_points
2582
FUNCTION linterpolate(n, x, xp, yp) RESULT(y)
2583
INTEGER :: n
2584
REAL(KIND=dp) :: x, y, xp(:), yp(:)
2585
2586
INTEGER :: i
2587
REAL(KIND=dp) :: x0,x1,y0,y1,t
2588
2589
y = 0._dp
2590
DO i=2,n
2591
x0=xp(i-1); x1=xp(i)
2592
IF((x >= x0) .AND. (x <= x1)) THEN
2593
y0=yp(i-1); y1=yp(i)
2594
t=(x-x0)/(x1-x0);
2595
y=(1-t)*y0 + t*y1;
2596
EXIT
2597
END IF
2598
END DO
2599
END FUNCTION linterpolate
2600
END SUBROUTINE ConvertTableToHarmonic
2601
2602
2603
SUBROUTINE ForceLoad
2604
CALL MPI_SEND()
2605
END SUBROUTINE ForceLoad
2606
2607
END MODULE GeneralUtils
2608
2609
2610
!---------------------------------------------------------
2611
!> Module mainly for writing xml based vtk files.
2612
!> The idea is that same routines save both the ascii
2613
!> and binary format.
2614
!---------------------------------------------------------
2615
MODULE AscBinOutputUtils
2616
2617
2618
USE Types
2619
IMPLICIT NONE
2620
2621
LOGICAL, PRIVATE :: AsciiOutput, SinglePrec, CalcSum = .FALSE.
2622
INTEGER, PRIVATE :: VtuUnit = 0, BufferSize = 0
2623
REAL, POINTER, PRIVATE :: FVals(:)
2624
REAL(KIND=dp), POINTER, PRIVATE :: DVals(:)
2625
INTEGER, POINTER, PRIVATE :: IVals(:)
2626
INTEGER, PRIVATE :: INoVals, NoVals
2627
REAL(KIND=dp) :: RSum = 0.0_dp, Isum = 0.0_dp
2628
INTEGER :: Rcount = 0, Icount = 0, Scount = 0, Ssum = 0
2629
2630
SAVE :: AsciiOutput, SinglePrec, VtuUnit, BufferSize, &
2631
FVals, DVals, IVals, CalcSum, Rsum, Isum, Ssum, RCount, Icount, Scount
2632
2633
2634
2635
CONTAINS
2636
2637
2638
! Initialize the buffer for writing, choose mode etc.
2639
!-----------------------------------------------------------------
2640
SUBROUTINE AscBinWriteInit( IsAscii, IsSingle, UnitNo, BufSize )
2641
2642
LOGICAL :: IsAscii, IsSingle
2643
INTEGER :: UnitNo, BufSize
2644
2645
AsciiOutput = IsAscii
2646
SinglePrec = IsSingle
2647
VtuUnit = UnitNo
2648
BufferSize = BufSize
2649
2650
CALL Info('AscBinWriteInit','Initializing buffered ascii/binary writing',Level=8)
2651
IF( AsciiOutput ) THEN
2652
CALL Info('AscBinWriteInit','Writing in ascii',Level=10)
2653
ELSE
2654
CALL Info('AscBinWriteInit','Writing in binary',Level=10)
2655
END IF
2656
2657
IF( SinglePrec ) THEN
2658
CALL Info('AscBinWriteInit','Writing in single precision',Level=10)
2659
ELSE
2660
CALL Info('AscBinWriteInit','Writing in double precision',Level=10)
2661
END IF
2662
2663
WRITE(Message,'(A,I0)') 'Writing to unit number: ',VtuUnit
2664
CALL Info('AscBinWriteInit',Message,Level=10)
2665
2666
IF(.NOT. AsciiOutput ) THEN
2667
WRITE(Message,'(A,I0)') 'Size of buffer is: ',BufferSize
2668
CALL Info('AscBinWriteInit',Message,Level=10)
2669
2670
ALLOCATE( Ivals( BufferSize ) )
2671
IF( SinglePrec ) THEN
2672
ALLOCATE( FVals( BufferSize ) )
2673
ELSE
2674
ALLOCATE( Dvals( BufferSize ) )
2675
END IF
2676
2677
INoVals = 0
2678
NoVals = 0
2679
END IF
2680
2681
END SUBROUTINE AscBinWriteInit
2682
2683
2684
! Free the buffer, next buffer can be different in size and type
2685
!---------------------------------------------------------------
2686
SUBROUTINE AscBinWriteFree()
2687
2688
CALL Info('AscBinWriteFree','Terminating buffered ascii/binary writing',Level=10)
2689
2690
IF( AsciiOutput ) RETURN
2691
2692
IF( SinglePrec ) THEN
2693
DEALLOCATE( FVals )
2694
ELSE
2695
DEALLOCATE( DVals )
2696
END IF
2697
DEALLOCATE( IVals )
2698
2699
BufferSize = 0
2700
VtuUnit = 0
2701
2702
END SUBROUTINE AscBinWriteFree
2703
2704
2705
2706
! The writing of xml strings is done here to allow easier modification
2707
! of output strategies.
2708
!-------------------------------------------------------------------------
2709
SUBROUTINE AscBinStrWrite( Str )
2710
2711
CHARACTER(LEN=1024) :: Str
2712
INTEGER, PARAMETER :: VtuUnit = 58
2713
2714
WRITE( VtuUnit ) TRIM(Str)
2715
IF( CalcSum ) THEN
2716
Scount = Scount + 1
2717
Ssum = Ssum + len_trim( Str )
2718
END IF
2719
2720
2721
END SUBROUTINE AscBinStrWrite
2722
2723
2724
! Write a binary value, either in single or double precision
2725
!------------------------------------------------------------------------
2726
SUBROUTINE AscBinRealWrite( val, EmptyBuffer )
2727
2728
INTEGER, PARAMETER :: VtuUnit = 58
2729
REAL(KIND=dp) :: val
2730
LOGICAL, OPTIONAL :: EmptyBuffer
2731
LOGICAL :: Empty
2732
CHARACTER(LEN=1024) :: Str
2733
2734
IF( VtuUnit == 0 ) THEN
2735
CALL Fatal('AscBinRealWrite','Buffer not initialized for writing')
2736
END IF
2737
2738
IF( PRESENT( EmptyBuffer ) ) THEN
2739
Empty = EmptyBuffer
2740
ELSE
2741
Empty = .FALSE.
2742
END IF
2743
2744
! Ascii output is not buffered
2745
IF( AsciiOutput ) THEN
2746
IF( Empty ) RETURN
2747
IF( ABS( val ) <= TINY ( val ) ) THEN
2748
WRITE(Str,'(A)') " 0.0"
2749
ELSE IF( SinglePrec ) THEN
2750
WRITE( Str,'(ES12.3E3)') val
2751
ELSE
2752
WRITE( Str,'(ES16.7E3)') val
2753
END IF
2754
WRITE( VtuUnit ) TRIM(Str)
2755
IF( CalcSum ) THEN
2756
Rcount = Rcount + 1
2757
Rsum = Rsum + val
2758
END IF
2759
RETURN
2760
END IF
2761
2762
! Buffered binary output
2763
IF( Empty .OR. NoVals == BufferSize ) THEN
2764
IF( NoVals == 0 ) THEN
2765
RETURN
2766
ELSE IF( SinglePrec ) THEN
2767
WRITE( VtuUnit ) Fvals(1:NoVals)
2768
ELSE
2769
WRITE( VtuUnit ) DVals(1:NoVals)
2770
END IF
2771
2772
IF( CalcSum ) THEN
2773
Rcount = Rcount + NoVals
2774
IF( SinglePrec ) THEN
2775
Rsum = 1.0_dp * SUM( Fvals(1:NoVals) )
2776
ELSE
2777
Rsum = SUM( DVals(1:NoVals) )
2778
END IF
2779
END IF
2780
2781
NoVals = 0
2782
IF( Empty ) RETURN
2783
END IF
2784
2785
! Save values in the buffer (either single or double prec.)
2786
NoVals = NoVals + 1
2787
IF( SinglePrec ) THEN
2788
Fvals(NoVals) = val
2789
ELSE
2790
DVals(NoVals) = val
2791
END IF
2792
2793
2794
END SUBROUTINE AscBinRealWrite
2795
2796
2797
! Write an integer value
2798
!-------------------------------------------------
2799
SUBROUTINE AscBinIntegerWrite( ival, EmptyBuffer )
2800
2801
INTEGER, PARAMETER :: VtuUnit = 58
2802
INTEGER :: ival
2803
LOGICAL, OPTIONAL :: EmptyBuffer
2804
LOGICAL :: Empty
2805
CHARACTER(LEN=1024) :: Str
2806
2807
IF( VtuUnit == 0 ) THEN
2808
CALL Fatal('AscBinIntegerWrite','Buffer not initialized for writing')
2809
END IF
2810
2811
IF( PRESENT( EmptyBuffer ) ) THEN
2812
Empty = EmptyBuffer
2813
ELSE
2814
Empty = .FALSE.
2815
END IF
2816
2817
IF( AsciiOutput ) THEN
2818
IF( Empty ) RETURN
2819
WRITE( Str, '(" ",I0)') ival
2820
WRITE( VtuUnit ) TRIM(Str)
2821
IF( CalcSum ) Isum = Isum + ival
2822
RETURN
2823
END IF
2824
2825
IF( Empty .OR. INoVals == BufferSize ) THEN
2826
IF( INoVals == 0 ) THEN
2827
RETURN
2828
ELSE
2829
WRITE( VtuUnit ) Ivals(1:INoVals)
2830
IF( CalcSum ) THEN
2831
Icount = Icount + 1
2832
Isum = Isum + SUM( Ivals(1:INoVals) )
2833
END IF
2834
END IF
2835
INoVals = 0
2836
IF( Empty ) RETURN
2837
END IF
2838
2839
INoVals = INoVals + 1
2840
Ivals(INoVals) = ival
2841
2842
END SUBROUTINE AscBinIntegerWrite
2843
2844
2845
SUBROUTINE AscBinInitNorm(CalcNorm)
2846
LOGICAL :: CalcNorm
2847
2848
CalcSum = CalcNorm
2849
RSum = 0.0_dp
2850
ISum = 0.0_dp
2851
Ssum = 0
2852
Rcount = 0
2853
Icount = 0
2854
Scount = 0
2855
2856
END SUBROUTINE AscBinInitNorm
2857
2858
2859
FUNCTION AscBinCompareNorm(RefResults,ExtResults) RESULT ( RelativeNorm )
2860
REAL(KIND=dp), DIMENSION(*) :: RefResults
2861
REAL(KIND=dp), DIMENSION(*), OPTIONAL :: ExtResults
2862
REAL(KIND=dp) :: RelativeNorm
2863
REAL(KIND=dp), POINTER :: ThisResults(:)
2864
REAL(KIND=dp) :: c
2865
INTEGER :: i, n
2866
2867
n = 6 !SIZE(RefResults)
2868
ALLOCATE(ThisResults(n))
2869
2870
IF( PRESENT( ExtResults ) ) THEN
2871
ThisResults(1:n) = ExtResults(1:n)
2872
ELSE
2873
ThisResults(1) = scount
2874
ThisResults(2) = icount
2875
ThisResults(3) = rcount
2876
ThisResults(4) = ssum
2877
ThisResults(5) = isum ! we use real for isum since it could be huge too!
2878
ThisResults(6) = rsum
2879
END IF
2880
2881
PRINT *,'Checksums for file output:'
2882
PRINT *,'RefResults:',RefResults(1:n)
2883
PRINT *,'ThisResults:',ThisResults(1:n)
2884
2885
! We have a special relative pseunonorm that should be one!
2886
RelativeNorm = 0.0
2887
DO i=1,n
2888
IF( ABS(RefResults(i) ) > EPSILON(c) ) THEN
2889
c = ThisResults(i)/RefResults(i)
2890
IF( ABS(ThisResults(i) ) > EPSILON(c) ) THEN
2891
c = MAX( c, 1.0_dp /c )
2892
END IF
2893
ELSE
2894
c = 1.0_dp + ABS(ThisResults(i))
2895
END IF
2896
RelativeNorm = RelativeNorm + c
2897
END DO
2898
RelativeNorm = RelativeNorm / n
2899
2900
END FUNCTION AscBinCompareNorm
2901
2902
2903
END MODULE AscBinOutputUtils
2904
2905
2906
!> \}
2907
2908
2909