88module m_hypoelastic
99
1010 use m_derived_types !< Definitions of the derived types
11-
1211 use m_global_parameters !< Definitions of the global parameters
13-
14- use m_mpi_proxy !< Message passing interface (MPI) module proxy
12+ use m_finite_differences
13+ use m_helper
1514
1615 implicit none
1716
1817 private; public :: s_initialize_hypoelastic_module, &
18+ s_finalize_hypoelastic_module, &
1919 s_compute_hypoelastic_rhs
2020
2121 real (wp), allocatable, dimension (:) :: Gs
@@ -29,11 +29,16 @@ module m_hypoelastic
2929 real (wp), allocatable, dimension (:, :, :) :: rho_K_field, G_K_field
3030 !$acc declare create(rho_K_field, G_K_field)
3131
32+ real (wp), allocatable, dimension (:, :) :: fd_coeff_x_h
33+ real (wp), allocatable, dimension (:, :) :: fd_coeff_y_h
34+ real (wp), allocatable, dimension (:, :) :: fd_coeff_z_h
35+ !$acc declare create(fd_coeff_x_h,fd_coeff_y_h,fd_coeff_z_h)
36+
3237contains
3338
3439 subroutine s_initialize_hypoelastic_module
3540
36- integer :: i
41+ integer :: i, k, r
3742
3843 @:ALLOCATE(Gs(1 :num_fluids))
3944 @:ALLOCATE(rho_K_field(0 :m,0 :n,0 :p), G_K_field(0 :m,0 :n,0 :p))
@@ -51,6 +56,29 @@ contains
5156 end do
5257 !$acc update device(Gs)
5358
59+ @:ALLOCATE(fd_coeff_x_h(- fd_number:fd_number, 0 :m))
60+ if (n > 0 ) then
61+ @:ALLOCATE(fd_coeff_y_h(- fd_number:fd_number, 0 :n))
62+ end if
63+ if (p > 0 ) then
64+ @:ALLOCATE(fd_coeff_z_h(- fd_number:fd_number, 0 :p))
65+ end if
66+
67+ ! Computing centered finite difference coefficients
68+ call s_compute_finite_difference_coefficients(m, x_cc, fd_coeff_x_h, buff_size, &
69+ fd_number, fd_order)
70+ !$acc update device(fd_coeff_x_h)
71+ if (n > 0 ) then
72+ call s_compute_finite_difference_coefficients(n, y_cc, fd_coeff_y_h, buff_size, &
73+ fd_number, fd_order)
74+ !$acc update device(fd_coeff_y_h)
75+ end if
76+ if (p > 0 ) then
77+ call s_compute_finite_difference_coefficients(p, z_cc, fd_coeff_z_h, buff_size, &
78+ fd_number, fd_order)
79+ !$acc update device(fd_coeff_z_h)
80+ end if
81+
5482 end subroutine s_initialize_hypoelastic_module
5583
5684 !> The purpose of this procedure is to compute the source terms
@@ -66,7 +94,7 @@ contains
6694
6795 real (wp) :: rho_K, G_K
6896
69- integer :: i, k, l, q !< Loop variables
97+ integer :: i, k, l, q, r !< Loop variables
7098 integer :: ndirs !< Number of coordinate directions
7199
72100 ndirs = 1 ; if (n > 0 ) ndirs = 2 ; if (p > 0 ) ndirs = 3
@@ -79,82 +107,91 @@ contains
79107 do q = 0 , p
80108 do l = 0 , n
81109 do k = 0 , m
82- du_dx(k, l, q) = &
83- (q_prim_vf(momxb)%sf(k - 2 , l, q) &
84- - 8._wp * q_prim_vf(momxb)%sf(k - 1 , l, q) &
85- + 8._wp * q_prim_vf(momxb)%sf(k + 1 , l, q) &
86- - q_prim_vf(momxb)%sf(k + 2 , l, q)) &
87- / (12._wp * dx(k))
110+ du_dx(k, l, q) = 0._wp
88111 end do
89112 end do
90113 end do
114+ !$acc end parallel loop
115+
116+ !$acc parallel loop collapse(3 ) gang vector default(present)
117+ do q = 0 , p
118+ do l = 0 , n
119+ do k = 0 , m
120+ !$acc loop seq
121+ do r = - fd_number, fd_number
122+ du_dx(k, l, q) = du_dx(k, l, q) &
123+ + q_prim_vf(momxb)%sf(k + r, l, q)* fd_coeff_x_h(r, k)
124+ end do
125+
126+ end do
127+ end do
128+ end do
129+ !$acc end parallel loop
91130
92131 if (ndirs > 1 ) then
93132 !$acc parallel loop collapse(3 ) gang vector default(present)
94133 do q = 0 , p
95134 do l = 0 , n
96135 do k = 0 , m
97- du_dy(k, l, q) = &
98- (q_prim_vf(momxb)%sf(k, l - 2 , q) &
99- - 8._wp * q_prim_vf(momxb)%sf(k, l - 1 , q) &
100- + 8._wp * q_prim_vf(momxb)%sf(k, l + 1 , q) &
101- - q_prim_vf(momxb)%sf(k, l + 2 , q)) &
102- / (12._wp * dy(l))
103- dv_dx(k, l, q) = &
104- (q_prim_vf(momxb + 1 )%sf(k - 2 , l, q) &
105- - 8._wp * q_prim_vf(momxb + 1 )%sf(k - 1 , l, q) &
106- + 8._wp * q_prim_vf(momxb + 1 )%sf(k + 1 , l, q) &
107- - q_prim_vf(momxb + 1 )%sf(k + 2 , l, q)) &
108- / (12._wp * dx(k))
109- dv_dy(k, l, q) = &
110- (q_prim_vf(momxb + 1 )%sf(k, l - 2 , q) &
111- - 8._wp * q_prim_vf(momxb + 1 )%sf(k, l - 1 , q) &
112- + 8._wp * q_prim_vf(momxb + 1 )%sf(k, l + 1 , q) &
113- - q_prim_vf(momxb + 1 )%sf(k, l + 2 , q)) &
114- / (12._wp * dy(l))
136+ du_dy(k, l, q) = 0._wp ; dv_dx(k, l, q) = 0._wp ; dv_dy(k, l, q) = 0._wp
137+ end do
138+ end do
139+ end do
140+ !$acc end parallel loop
141+
142+ !$acc parallel loop collapse(3 ) gang vector default(present)
143+ do q = 0 , p
144+ do l = 0 , n
145+ do k = 0 , m
146+ !$acc loop seq
147+ do r = - fd_number, fd_number
148+ du_dy(k, l, q) = du_dy(k, l, q) &
149+ + q_prim_vf(momxb)%sf(k, l + r, q)* fd_coeff_y_h(r, l)
150+ dv_dx(k, l, q) = dv_dx(k, l, q) &
151+ + q_prim_vf(momxb + 1 )%sf(k + r, l, q)* fd_coeff_x_h(r, k)
152+ dv_dy(k, l, q) = dv_dy(k, l, q) &
153+ + q_prim_vf(momxb + 1 )%sf(k, l + r, q)* fd_coeff_y_h(r, l)
154+ end do
115155 end do
116156 end do
117157 end do
158+ !$acc end parallel loop
118159
119160 ! 3D
120161 if (ndirs == 3 ) then
162+
121163 !$acc parallel loop collapse(3 ) gang vector default(present)
122164 do q = 0 , p
123165 do l = 0 , n
124166 do k = 0 , m
125- du_dz(k, l, q) = &
126- (q_prim_vf(momxb)%sf(k, l, q - 2 ) &
127- - 8._wp * q_prim_vf(momxb)%sf(k, l, q - 1 ) &
128- + 8._wp * q_prim_vf(momxb)%sf(k, l, q + 1 ) &
129- - q_prim_vf(momxb)%sf(k, l, q + 2 )) &
130- / (12._wp * dz(q))
131- dv_dz(k, l, q) = &
132- (q_prim_vf(momxb + 1 )%sf(k, l, q - 2 ) &
133- - 8._wp * q_prim_vf(momxb + 1 )%sf(k, l, q - 1 ) &
134- + 8._wp * q_prim_vf(momxb + 1 )%sf(k, l, q + 1 ) &
135- - q_prim_vf(momxb + 1 )%sf(k, l, q + 2 )) &
136- / (12._wp * dz(q))
137- dw_dx(k, l, q) = &
138- (q_prim_vf(momxe)%sf(k - 2 , l, q) &
139- - 8._wp * q_prim_vf(momxe)%sf(k - 1 , l, q) &
140- + 8._wp * q_prim_vf(momxe)%sf(k + 1 , l, q) &
141- - q_prim_vf(momxe)%sf(k + 2 , l, q)) &
142- / (12._wp * dx(k))
143- dw_dy(k, l, q) = &
144- (q_prim_vf(momxe)%sf(k, l - 2 , q) &
145- - 8._wp * q_prim_vf(momxe)%sf(k, l - 1 , q) &
146- + 8._wp * q_prim_vf(momxe)%sf(k, l + 1 , q) &
147- - q_prim_vf(momxe)%sf(k, l + 2 , q)) &
148- / (12._wp * dy(l))
149- dw_dz(k, l, q) = &
150- (q_prim_vf(momxe)%sf(k, l, q - 2 ) &
151- - 8._wp * q_prim_vf(momxe)%sf(k, l, q - 1 ) &
152- + 8._wp * q_prim_vf(momxe)%sf(k, l, q + 1 ) &
153- - q_prim_vf(momxe)%sf(k, l, q + 2 )) &
154- / (12._wp * dz(q))
167+ du_dz(k, l, q) = 0_wp ; dv_dz(k, l, q) = 0_wp ; dw_dx(k, l, q) = 0_wp ;
168+ dw_dy(k, l, q) = 0_wp ; dw_dz(k, l, q) = 0_wp ;
155169 end do
156170 end do
157171 end do
172+ !$acc end parallel loop
173+
174+ !$acc parallel loop collapse(3 ) gang vector default(present)
175+ do q = 0 , p
176+ do l = 0 , n
177+ do k = 0 , m
178+ !$acc loop seq
179+ do r = - fd_number, fd_number
180+ du_dz(k, l, q) = du_dz(k, l, q) &
181+ + q_prim_vf(momxb)%sf(k, l, q + r)* fd_coeff_z_h(r, q)
182+ dv_dz(k, l, q) = dv_dz(k, l, q) &
183+ + q_prim_vf(momxb + 1 )%sf(k, l, q + r)* fd_coeff_z_h(r, q)
184+ dw_dx(k, l, q) = dw_dx(k, l, q) &
185+ + q_prim_vf(momxe)%sf(k + r, l, q)* fd_coeff_x_h(r, k)
186+ dw_dy(k, l, q) = dw_dy(k, l, q) &
187+ + q_prim_vf(momxe)%sf(k, l + r, q)* fd_coeff_y_h(r, l)
188+ dw_dz(k, l, q) = dw_dz(k, l, q) &
189+ + q_prim_vf(momxe)%sf(k, l, q + r)* fd_coeff_z_h(r, q)
190+ end do
191+ end do
192+ end do
193+ end do
194+ !$acc end parallel loop
158195 end if
159196 end if
160197
@@ -171,7 +208,7 @@ contains
171208 G_K_field(k, l, q) = G_K
172209
173210 !TODO: take this out if not needed
174- if (G_K < 1000 ) then
211+ if (G_K < verysmall ) then
175212 G_K_field(k, l, q) = 0
176213 end if
177214 end do
@@ -296,4 +333,21 @@ contains
296333
297334 end subroutine s_compute_hypoelastic_rhs
298335
336+ subroutine s_finalize_hypoelastic_module ()
337+
338+ @:DEALLOCATE(Gs)
339+ @:DEALLOCATE(rho_K_field, G_K_field)
340+ @:DEALLOCATE(du_dx)
341+ @:DEALLOCATE(fd_coeff_x_h)
342+ if (n > 0 ) then
343+ @:DEALLOCATE(du_dy,dv_dx,dv_dy)
344+ @:DEALLOCATE(fd_coeff_y_h)
345+ if (p > 0 ) then
346+ @:DEALLOCATE(du_dz, dv_dz, dw_dx, dw_dy, dw_dz)
347+ @:DEALLOCATE(fd_coeff_z_h)
348+ end if
349+ end if
350+
351+ end subroutine s_finalize_hypoelastic_module
352+
299353end module m_hypoelastic
0 commit comments