Skip to content

Commit 658f4c6

Browse files
committed
fix: work around gfortran bugs
1 parent fa07f94 commit 658f4c6

3 files changed

Lines changed: 108 additions & 25 deletions

File tree

src/matcha/distribution_s.F90

Lines changed: 3 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -39,7 +39,7 @@ pure function monotonically_increasing(f) result(monotonic)
3939
module procedure velocities
4040

4141
double precision, allocatable :: sampled_speeds(:,:), dir(:,:,:)
42-
integer cell, step
42+
integer cell, step, k
4343

4444
call_assert(allocated(self%cumulative_distribution_))
4545
call_assert(allocated(self%vel_))
@@ -48,9 +48,8 @@ pure function monotonically_increasing(f) result(monotonic)
4848
associate(ncells => size(speeds,1), nsteps => size(speeds,2))
4949
allocate(sampled_speeds(ncells,nsteps))
5050
do concurrent(cell = 1:ncells, step = 1:nsteps)
51-
associate(k => findloc(speeds(cell,step) >= self%cumulative_distribution(), value=.false., dim=1)-1)
52-
sampled_speeds(cell,step) = self%vel_(k)
53-
end associate
51+
k = findloc(speeds(cell,step) >= self%cumulative_distribution(), value=.false., dim=1)-1
52+
sampled_speeds(cell,step) = self%vel_(k)
5453
end do
5554

5655
! Create unit vectors

src/matcha/output_s.f90

Lines changed: 9 additions & 13 deletions
Original file line numberDiff line numberDiff line change
@@ -51,8 +51,8 @@ pure function sim_speeds(history) result(speeds)
5151
double precision, allocatable :: speeds(:)
5252

5353
integer, parameter :: nspacedims=3
54-
integer i, j, k
55-
double precision, allocatable :: x(:,:,:)
54+
integer i, j, k, ij
55+
double precision, allocatable :: x(:,:,:), t(:), u(:)
5656

5757
associate( &
5858
npositions => size(history,1), &
@@ -62,17 +62,13 @@ pure function sim_speeds(history) result(speeds)
6262
do concurrent(i=1:npositions)
6363
x(i,:,:) = history(i)%positions()
6464
end do
65-
associate(t => history%time())
66-
allocate(speeds(ncells*(npositions-1)))
67-
do concurrent(i = 1:npositions-1, j = 1:ncells)
68-
associate( &
69-
u => (x(i+1,j,:) - x(i,j,:))/(t(i+1) - t(i)), &
70-
ij => i + (j-1)*(npositions-1) &
71-
)
72-
speeds(ij) = sqrt(sum([(u(k)**2, k=1,nspacedims)]))
73-
end associate
74-
end do
75-
end associate
65+
t = history%time()
66+
allocate(speeds(ncells*(npositions-1)))
67+
do concurrent(i = 1:npositions-1, j = 1:ncells)
68+
u = (x(i+1,j,:) - x(i,j,:))/(t(i+1) - t(i))
69+
ij = i + (j-1)*(npositions-1)
70+
speeds(ij) = sqrt(sum([(u(k)**2, k=1,nspacedims)]))
71+
end do
7672
end associate
7773
end function
7874

test/matcha_test_m.F90

Lines changed: 96 additions & 8 deletions
Original file line numberDiff line numberDiff line change
@@ -5,14 +5,15 @@
55

66
module matcha_test_m
77
use julienne_m, only : &
8-
string_t, &
9-
test_t, &
10-
test_diagnosis_t &
11-
,vector_test_description_t &
12-
#if ! HAVE_PROCEDURE_ACTUAL_FOR_POINTER_DUMMY
13-
,vector_diagnosis_function_i &
8+
string_t &
9+
#if __GFORTRAN__
10+
,diagnosis_function_i &
11+
,test_description_t &
1412
#endif
15-
,test_result_t
13+
,test_diagnosis_t &
14+
,vector_test_description_t &
15+
,test_result_t &
16+
,test_t
1617
use input_m, only : input_t
1718
use output_m, only : output_t
1819
use matcha_m, only : matcha
@@ -35,6 +36,7 @@ pure function subject() result(specimen)
3536
specimen = "A matcha_t object"
3637
end function
3738

39+
#ifndef __GFORTRAN__
3840
function results() result(test_results)
3941
type(test_result_t), allocatable :: test_results(:)
4042
type(vector_test_description_t), allocatable :: vector_test_descriptions(:)
@@ -54,7 +56,26 @@ function results() result(test_results)
5456
test_results = vector_test_descriptions(i)%run()
5557
end do
5658
end function
59+
#else
60+
function results() result(test_results)
61+
type(test_result_t), allocatable :: test_results(:)
62+
type(test_description_t), allocatable :: test_descriptions(:)
63+
procedure(diagnosis_function_i), pointer :: &
64+
compare_image_distributions_ptr &
65+
,compare_global_distributions_ptr
66+
67+
compare_image_distributions_ptr => compare_image_distributions
68+
compare_global_distributions_ptr => compare_global_distributions
69+
70+
test_descriptions = [ &
71+
test_description_t("matching simulated distributions to empirical distribution", compare_image_distributions_ptr) &
72+
,test_description_t("matching simulated global distributions to empirical distribution", compare_global_distributions_ptr) &
73+
]
74+
test_results = test_descriptions%run()
75+
end function
76+
#endif
5777

78+
#ifndef __GFORTRAN__
5879
function compare_image_distributions() result(test_diagnoses)
5980
logical test_passes
6081
type(test_diagnosis_t), allocatable :: test_diagnoses(:)
@@ -87,7 +108,37 @@ function compare_image_distributions() result(test_diagnoses)
87108
end associate
88109
end associate
89110
end function
90-
111+
#else
112+
function compare_image_distributions() result(test_diagnosis)
113+
logical test_passes
114+
type(test_diagnosis_t) test_diagnosis
115+
type(output_t) output
116+
117+
integer, parameter :: speed=1, freq=2 ! subscripts for speeds and frequencies
118+
double precision, parameter :: tolerance = 1.D-02
119+
120+
associate(input => input_t())
121+
output = output_t(input, matcha(input))
122+
associate( &
123+
empirical_distribution => input%sample_distribution(), &
124+
simulated_distribution => output%simulated_distribution() &
125+
)
126+
associate( &
127+
diffmax_speeds=> maxval(abs(empirical_distribution(:,speed)-simulated_distribution(:,speed))), &
128+
diffmax_freqs => maxval(abs(empirical_distribution(:,freq)-simulated_distribution(:,freq))) &
129+
)
130+
test_diagnosis = test_diagnosis_t( &
131+
test_passed = diffmax_freqs < tolerance .and. diffmax_speeds < tolerance &
132+
,diagnostics_string = "diffmax_speeds " // string_t(diffmax_speeds) // ", diffmax_freqs " // string_t(diffmax_freqs) &
133+
// ", expected < " // string_t(tolerance) &
134+
)
135+
end associate
136+
end associate
137+
end associate
138+
end function
139+
#endif
140+
141+
#ifndef __GFORTRAN__
91142
function compare_global_distributions() result(test_diagnoses)
92143
type(test_diagnosis_t), allocatable :: test_diagnoses(:)
93144
type(output_t) output
@@ -128,5 +179,42 @@ function compare_global_distributions() result(test_diagnoses)
128179
end associate
129180

130181
end function
182+
#else
183+
function compare_global_distributions() result(test_diagnosis)
184+
type(test_diagnosis_t) test_diagnosis
185+
type(output_t) output
186+
double precision, allocatable :: simulated_distribution(:,:)
187+
integer num_cells
188+
integer, parameter :: speed=1, freq=2 ! subscripts for speeds and frequencies
189+
double precision, parameter :: tolerance = 1.D-02
190+
191+
associate(input => input_t())
192+
output = output_t(input, matcha(input))
193+
associate(empirical_distribution => input%sample_distribution())
194+
simulated_distribution = output%simulated_distribution()
195+
num_cells = output%my_num_cells()
196+
simulated_distribution(:,freq) = num_cells*simulated_distribution(:,freq)
197+
call co_sum(simulated_distribution(:,freq), result_image=1)
198+
call co_sum(num_cells, result_image=1)
199+
if (this_image()/=1) then
200+
test_diagnosis = test_diagnosis_t(test_passed=.true., diagnostics_string="")
201+
else
202+
simulated_distribution(:,freq) = simulated_distribution(:,freq)/dble(num_cells)
203+
associate( &
204+
diffmax_speeds=> maxval(abs(empirical_distribution(:,speed)-simulated_distribution(:,speed))), &
205+
diffmax_freqs => maxval(abs(empirical_distribution(:,freq)-simulated_distribution(:,freq))) &
206+
)
207+
test_diagnosis = test_diagnosis_t( &
208+
test_passed = diffmax_freqs < tolerance .and. diffmax_speeds < tolerance &
209+
,diagnostics_string = "diffmax_speeds " // string_t(diffmax_speeds) // ", diffmax_freqs " // string_t(diffmax_freqs) &
210+
// ", expected < " // string_t(tolerance) &
211+
)
212+
end associate
213+
end if
214+
end associate
215+
end associate
216+
217+
end function
218+
#endif
131219

132220
end module matcha_test_m

0 commit comments

Comments
 (0)