SFEMaNS  version 4.1 (work in progress)
Reference documentation for SFEMaNS
 All Classes Files Functions Variables Groups Pages
computing_rhs_gauss.f90
Go to the documentation of this file.
3 PRIVATE
4 CONTAINS
5 
6  SUBROUTINE rhs_ns_gauss_3x3(vv_mesh, pp_mesh, communicator, list_mode, time, V1m, pn, pn_inc, rotv_v, &
7  rhs_gauss, opt_tempn)
8  !=================================
9  !RHS for Navier-Stokes
10  USE def_type_mesh
11  USE my_util
12  USE input_data
13  USE sft_parallele
14  USE boundary
15  IMPLICIT NONE
16  TYPE(mesh_type) :: vv_mesh, pp_mesh
17  REAL(KIND=8), INTENT(IN) :: time
18  REAL(KIND=8), DIMENSION(:,:,:), INTENT(IN) :: rotv_v
19  REAL(KIND=8), DIMENSION(:,:,:), INTENT(IN) :: v1m
20  REAL(KIND=8), DIMENSION(:,:,:), INTENT(IN) :: pn, pn_inc
21  INTEGER, DIMENSION(:), INTENT(IN) :: list_mode
22  REAL(KIND=8), DIMENSION(:,:,:), INTENT(IN), OPTIONAL :: opt_tempn
23  REAL(KIND=8), DIMENSION(vv_mesh%gauss%l_G*vv_mesh%dom_me,6,SIZE(list_mode)), INTENT(OUT) :: rhs_gauss
24  REAL(KIND=8), DIMENSION(6) :: fs, ft, fp_inc
25  INTEGER, DIMENSION(vv_mesh%gauss%n_w) :: j_loc
26  REAL(KIND=8), DIMENSION(vv_mesh%np,6) :: ff, imposed_vel
27  REAL(KIND=8), DIMENSION(vv_mesh%np,2) :: p, p_inc
28  REAL(KIND=8), DIMENSION(vv_mesh%gauss%k_d,vv_mesh%gauss%n_w) :: dw_loc
29  REAL(KIND=8), DIMENSION(vv_mesh%dom_me*vv_mesh%gauss%l_G,6) :: imposed_vel_gauss
30  REAL(KIND=8), DIMENSION(vv_mesh%dom_me*vv_mesh%gauss%l_G,6,SIZE(list_mode)) :: fp, rhs_gauss_penal
31  REAL(KIND=8), DIMENSION(2,vv_mesh%gauss%l_G*vv_mesh%dom_me) :: rr_gauss
32  REAL(KIND=8) :: ray
33  INTEGER :: m, l , i, k, index, type
34  INTEGER :: nb_procs, m_max_pad, bloc_size
35 #include "petsc/finclude/petsc.h"
36  petscerrorcode :: ierr
37  mpi_comm :: communicator
38 
39  DO i = 1, SIZE(list_mode)
40  DO k= 1, 6
41  IF (present(opt_tempn)) THEN
42  ff(:,k) = source_in_ns_momentum(k, vv_mesh%rr, list_mode(i), i, time, inputs%Re, 'ns', &
43  opt_tempn=opt_tempn)
44  ELSE
45  ff(:,k) = source_in_ns_momentum(k, vv_mesh%rr, list_mode(i), i, time, inputs%Re, 'ns')
46  END IF
47  END DO
48  DO k = 1, 2
49  CALL inject(pp_mesh%jj, vv_mesh%jj, pn(:,k,i), p(:,k))
50  CALL inject(pp_mesh%jj, vv_mesh%jj, pn_inc(:,k,i), p_inc(:,k))
51  ENDDO
52 
53  index = 0
54  DO m = 1, vv_mesh%dom_me
55  j_loc = vv_mesh%jj(:,m)
56 
57  DO l = 1, vv_mesh%gauss%l_G
58  index = index +1
59  dw_loc = vv_mesh%gauss%dw(:,:,l,m)
60 
61  !===Compute radius of Gauss point
62  ray = sum(vv_mesh%rr(1,j_loc)*vv_mesh%gauss%ww(:,l))
63  rr_gauss(1,index) = ray
64  rr_gauss(2,index) = sum(vv_mesh%rr(2,j_loc)*vv_mesh%gauss%ww(:,l))
65 
66  !--------calcul de la premiere composante : u0(1,:) <--> f(r,m,c)
67  fs(1) = sum(ff(j_loc,1) * vv_mesh%gauss%ww(:,l))
68  ft(1) = sum(v1m(j_loc,1,i) * vv_mesh%gauss%ww(:,l))
69  fp(index,1,i) = -sum(p(j_loc,1)*dw_loc(1,:))
70  fp_inc(1) = -sum(p_inc(j_loc,1)*dw_loc(1,:))
71  !--------calcul de la seconde composante : u0(2,:) <--> f(r,m,s)
72  fs(2) = sum(ff(j_loc,2) * vv_mesh%gauss%ww(:,l))
73  ft(2) = sum(v1m(j_loc,2,i) * vv_mesh%gauss%ww(:,l))
74  fp(index,2,i) = -sum(p(j_loc,2)*dw_loc(1,:))
75  fp_inc(2) = -sum(p_inc(j_loc,2)*dw_loc(1,:))
76  !--------calcul de la troisieme composante : u0(3,:) <--> f(th,m,c)
77  fs(3) = sum(ff(j_loc,3) * vv_mesh%gauss%ww(:,l))
78  ft(3) = sum(v1m(j_loc,3,i) * vv_mesh%gauss%ww(:,l))
79  fp(index,3,i) = -sum(p(j_loc,2)*vv_mesh%gauss%ww(:,l))/ray*list_mode(i)
80  fp_inc(3) = -sum(p_inc(j_loc,2)*vv_mesh%gauss%ww(:,l))/ray*list_mode(i)
81  !--------calcul de la quatrieme composante : u0(4,:) <--> f(th,m,s)
82  fs(4) = sum(ff(j_loc,4) * vv_mesh%gauss%ww(:,l))
83  ft(4) = sum(v1m(j_loc,4,i) * vv_mesh%gauss%ww(:,l))
84  fp(index,4,i) = sum(p(j_loc,1)*vv_mesh%gauss%ww(:,l))/ray*list_mode(i)
85  fp_inc(4) = sum(p_inc(j_loc,1)*vv_mesh%gauss%ww(:,l))/ray*list_mode(i)
86  !--------calcul de la cinquieme composante : u0(5,:) <--> f(z,m,c)
87  fs(5) = sum(ff(j_loc,5) * vv_mesh%gauss%ww(:,l))
88  ft(5) = sum(v1m(j_loc,5,i) * vv_mesh%gauss%ww(:,l))
89  fp(index,5,i) = -sum(p(j_loc,1)*dw_loc(2,:))
90  fp_inc(5) = -sum(p_inc(j_loc,1)*dw_loc(2,:))
91  !--------calcul de la sixieme composante : u0(6,:) <--> f(z,m,s)
92  fs(6) = sum(ff(j_loc,6) * vv_mesh%gauss%ww(:,l))
93  ft(6) = sum(v1m(j_loc,6,i) * vv_mesh%gauss%ww(:,l))
94  fp(index,6,i) = -sum(p(j_loc,2)*dw_loc(2,:))
95  fp_inc(6) = -sum(p_inc(j_loc,2)*dw_loc(2,:))
96 
97  rhs_gauss(index,:,i) = (ft+fp_inc+fs-rotv_v(index,:,i))
98 
99  ENDDO
100  ENDDO
101  IF (inputs%if_ns_penalty) THEN
102  IF(inputs%if_impose_vel_in_solids) THEN
103  IF (list_mode(i)==0) THEN
104  !Velocity imposed by penalty in solids (using BDF2).
105  imposed_vel(:,:) = 3.d0*imposed_velocity_by_penalty(vv_mesh%rr,time)/(2*inputs%dt)
106  index = 0
107  DO m = 1, vv_mesh%dom_me
108  j_loc = vv_mesh%jj(:,m)
109  DO l = 1, vv_mesh%gauss%l_G
110  index = index +1
111  DO TYPE = 1, 6
112  imposed_vel_gauss(index,type) = sum(imposed_vel(j_loc,type) &
113  * vv_mesh%gauss%ww(:,l))
114  END DO
115  END DO
116  END DO
117  rhs_gauss(:,:,i) = rhs_gauss(:,:,i) - imposed_vel_gauss(:,:)
118  ENDIF
119  END IF
120  END IF
121  END DO
122 
123  IF (inputs%if_ns_penalty) THEN
124  CALL mpi_comm_size(communicator, nb_procs, ierr)
125  m_max_pad = 3*SIZE(list_mode)*nb_procs/2
126  bloc_size = SIZE(rhs_gauss,1)/nb_procs+1
127  CALL fft_par_var_eta_prod_gauss_dcl(communicator, penal_in_real_space, vv_mesh, &
128  rhs_gauss, rhs_gauss_penal, nb_procs, bloc_size, m_max_pad, rr_gauss, time)
129 
130  rhs_gauss = rhs_gauss_penal
131 
132  IF(inputs%if_impose_vel_in_solids) THEN
133  DO i = 1, SIZE(list_mode)
134  IF (list_mode(i)==0) THEN
135  rhs_gauss(:,:,i) = rhs_gauss_penal(:,:,i) + imposed_vel_gauss(:,:)
136  END IF
137  END DO
138  END IF
139  END IF
140 
141  rhs_gauss = rhs_gauss + fp
142 
143  END SUBROUTINE rhs_ns_gauss_3x3
144 
145  SUBROUTINE rhs_residual_ns_gauss_3x3(vv_mesh, pp_mesh, communicator, list_mode, time, du_dt,&
146  pn, rotv_v, rhs_gauss, opt_tempn)
147  !=================================
148  !RHS for Residual of Navier-Stokes
149  USE def_type_mesh
150  USE my_util
151  USE input_data
152  USE sft_parallele
153  USE boundary
154  IMPLICIT NONE
155  TYPE(mesh_type) :: vv_mesh, pp_mesh
156  INTEGER, DIMENSION(:), INTENT(IN) :: list_mode
157  REAL(KIND=8), INTENT(IN) :: time
158  REAL(KIND=8), DIMENSION(:,:,:), INTENT(IN) :: du_dt
159  REAL(KIND=8), DIMENSION(:,:,:), INTENT(IN) :: pn
160  REAL(KIND=8), DIMENSION(:,:,:), INTENT(IN) :: rotv_v
161  REAL(KIND=8), DIMENSION(:,:,:), INTENT(IN), OPTIONAL :: opt_tempn
162  REAL(KIND=8), DIMENSION(vv_mesh%gauss%l_G*vv_mesh%dom_me,6,SIZE(list_mode)), INTENT(OUT) :: rhs_gauss
163  REAL(KIND=8), DIMENSION(6) :: fs
164  INTEGER, DIMENSION(vv_mesh%gauss%n_w) :: j_loc
165  REAL(KIND=8), DIMENSION(vv_mesh%np,6) :: ff
166  REAL(KIND=8), DIMENSION(vv_mesh%np,2) :: p
167  REAL(KIND=8), DIMENSION(vv_mesh%gauss%k_d,vv_mesh%gauss%n_w) :: dw_loc
168  REAL(KIND=8), DIMENSION(vv_mesh%dom_me*vv_mesh%gauss%l_G,6,SIZE(list_mode)) :: rhs_gauss_penal, fp, ft
169  REAL(KIND=8), DIMENSION(2,vv_mesh%gauss%l_G*vv_mesh%dom_me) :: rr_gauss
170  REAL(KIND=8) :: ray
171  INTEGER :: m, l , i, k, index
172  INTEGER :: nb_procs, m_max_pad, bloc_size
173 #include "petsc/finclude/petsc.h"
174  petscerrorcode :: ierr
175  mpi_comm :: communicator
176 
177  DO i = 1, SIZE(list_mode)
178  DO k = 1, 6
179  IF (present(opt_tempn)) THEN
180  ff(:,k) = source_in_ns_momentum(k, vv_mesh%rr, list_mode(i), i, time, inputs%Re, 'ns', &
181  opt_tempn=opt_tempn)
182  ELSE
183  ff(:,k) = source_in_ns_momentum(k, vv_mesh%rr, list_mode(i), i, time, inputs%Re, 'ns')
184  END IF
185  END DO
186  DO k = 1, 2
187  CALL inject(pp_mesh%jj, vv_mesh%jj, pn(:,k,i), p(:,k))
188  ENDDO
189 
190  index = 0
191  DO m = 1, vv_mesh%dom_me
192  j_loc = vv_mesh%jj(:,m)
193 
194  DO l = 1, vv_mesh%gauss%l_G
195  index = index +1
196  dw_loc = vv_mesh%gauss%dw(:,:,l,m)
197 
198  !===Compute radius of Gauss point
199  ray = sum(vv_mesh%rr(1,j_loc)*vv_mesh%gauss%ww(:,l))
200  rr_gauss(1,index) = ray
201  rr_gauss(2,index) = sum(vv_mesh%rr(2,j_loc)*vv_mesh%gauss%ww(:,l))
202 
203  !--------calcul de la premiere composante : u0(1,:) <--> f(r,m,c)
204  fs(1) = sum(ff(j_loc,1) * vv_mesh%gauss%ww(:,l))
205  ft(index,1,i) = sum(du_dt(j_loc,1,i) * vv_mesh%gauss%ww(:,l))
206  fp(index,1,i) = sum(p(j_loc,1)*dw_loc(1,:))
207  !--------calcul de la seconde composante : u0(2,:) <--> f(r,m,s)
208  fs(2) = sum(ff(j_loc,2) * vv_mesh%gauss%ww(:,l))
209  ft(index,2,i) = sum(du_dt(j_loc,2,i) * vv_mesh%gauss%ww(:,l))
210  fp(index,2,i) = sum(p(j_loc,2)*dw_loc(1,:))
211  !--------calcul de la troisieme composante : u0(3,:) <--> f(th,m,c)
212  fs(3) = sum(ff(j_loc,3) * vv_mesh%gauss%ww(:,l))
213  ft(index,3,i) = sum(du_dt(j_loc,3,i) * vv_mesh%gauss%ww(:,l))
214  fp(index,3,i) = sum(p(j_loc,2)*vv_mesh%gauss%ww(:,l))/ray*list_mode(i)
215  !--------calcul de la quatrieme composante : u0(4,:) <--> f(th,m,s)
216  fs(4) = sum(ff(j_loc,4) * vv_mesh%gauss%ww(:,l))
217  ft(index,4,i) = sum(du_dt(j_loc,4,i) * vv_mesh%gauss%ww(:,l))
218  fp(index,4,i) = -sum(p(j_loc,1)*vv_mesh%gauss%ww(:,l))/ray*list_mode(i)
219  !--------calcul de la cinquieme composante : u0(5,:) <--> f(z,m,c)
220  fs(5) = sum(ff(j_loc,5) * vv_mesh%gauss%ww(:,l))
221  ft(index,5,i) = sum(du_dt(j_loc,5,i) * vv_mesh%gauss%ww(:,l))
222  fp(index,5,i) = sum(p(j_loc,1)*dw_loc(2,:))
223  !--------calcul de la sixieme composante : u0(6,:) <--> f(z,m,s)
224  fs(6) = sum(ff(j_loc,6) * vv_mesh%gauss%ww(:,l))
225  ft(index,6,i) = sum(du_dt(j_loc,6,i) * vv_mesh%gauss%ww(:,l))
226  fp(index,6,i) = sum(p(j_loc,2)*dw_loc(2,:))
227 
228  rhs_gauss(index,:,i) = -fs+rotv_v(index,:,i)
229  ENDDO
230  ENDDO
231  END DO
232 
233  IF (inputs%if_ns_penalty) THEN
234  CALL mpi_comm_size(communicator, nb_procs, ierr)
235  m_max_pad = 3*SIZE(list_mode)*nb_procs/2
236  bloc_size = SIZE(rhs_gauss,1)/nb_procs+1
237  CALL fft_par_var_eta_prod_gauss_dcl(communicator, penal_in_real_space, vv_mesh, &
238  rhs_gauss, rhs_gauss_penal, nb_procs, bloc_size, m_max_pad, rr_gauss, time)
239  rhs_gauss = rhs_gauss_penal
240  END IF
241 
242  !===Add term time derivative and grad pressure
243  rhs_gauss = rhs_gauss + ft + fp
244 
245  END SUBROUTINE rhs_residual_ns_gauss_3x3
246 
247  SUBROUTINE rhs_residual_ns_gauss_3x3_mom(vv_mesh, pp_mesh, list_mode, time, du_dt, pn, &
248  density, rotb_b, rhs_gauss, opt_tempn)
249  !=================================
250  !RHS for Navier-Stokes
251  USE def_type_mesh
252  USE my_util
253  USE input_data
254  USE sft_parallele
255  USE boundary
256  IMPLICIT NONE
257  TYPE(mesh_type) :: vv_mesh, pp_mesh
258  INTEGER, DIMENSION(:), INTENT(IN) :: list_mode
259  REAL(KIND=8), INTENT(IN) :: time
260  REAL(KIND=8), DIMENSION(:,:,:), INTENT(IN) :: du_dt
261  REAL(KIND=8), DIMENSION(:,:,:), INTENT(IN) :: pn
262  REAL(KIND=8), DIMENSION(:,:,:), INTENT(IN) :: density
263  REAL(KIND=8), DIMENSION(:,:,:), INTENT(IN) :: rotb_b
264  REAL(KIND=8), DIMENSION(:,:,:), INTENT(IN), OPTIONAL :: opt_tempn
265  REAL(KIND=8), DIMENSION(vv_mesh%gauss%l_G*vv_mesh%dom_me,6,SIZE(list_mode)), INTENT(OUT) :: rhs_gauss
266  REAL(KIND=8), DIMENSION(6) :: fs, ft, fp
267  INTEGER, DIMENSION(vv_mesh%gauss%n_w) :: j_loc
268  REAL(KIND=8), DIMENSION(vv_mesh%np,6) :: ff
269  REAL(KIND=8), DIMENSION(vv_mesh%np,2) :: p
270  REAL(KIND=8), DIMENSION(vv_mesh%gauss%k_d,vv_mesh%gauss%n_w) :: dw_loc
271  REAL(KIND=8), DIMENSION(2,vv_mesh%gauss%l_G*vv_mesh%dom_me) :: rr_gauss
272  REAL(KIND=8) :: ray
273  INTEGER :: m, l , i, k, index
274 
275  DO i = 1, SIZE(list_mode)
276  DO k = 1, 6
277  IF (inputs%if_temperature) THEN
278  ff(:,k) = source_in_ns_momentum(k, vv_mesh%rr, list_mode(i), i, time, inputs%Re, 'ns', &
279  opt_density=density, opt_tempn=opt_tempn)
280  ELSE
281  ff(:,k) = source_in_ns_momentum(k, vv_mesh%rr, list_mode(i), i, time, inputs%Re, 'ns', &
282  opt_density=density)
283  END IF
284  END DO
285  DO k = 1, 2
286  CALL inject(pp_mesh%jj, vv_mesh%jj, pn(:,k,i), p(:,k))
287  ENDDO
288 
289  index = 0
290  DO m = 1, vv_mesh%dom_me
291  j_loc = vv_mesh%jj(:,m)
292 
293  DO l = 1, vv_mesh%gauss%l_G
294  index = index +1
295  dw_loc = vv_mesh%gauss%dw(:,:,l,m)
296 
297  !===Compute radius of Gauss point
298  ray = sum(vv_mesh%rr(1,j_loc)*vv_mesh%gauss%ww(:,l))
299  rr_gauss(1,index) = ray
300  rr_gauss(2,index) = sum(vv_mesh%rr(2,j_loc)*vv_mesh%gauss%ww(:,l))
301 
302  !--------calcul de la premiere composante : u0(1,:) <--> f(r,m,c)
303  fs(1) = sum(ff(j_loc,1) * vv_mesh%gauss%ww(:,l))
304  ft(1) = sum(du_dt(j_loc,1,i) * vv_mesh%gauss%ww(:,l))
305  fp(1) = sum(p(j_loc,1)*dw_loc(1,:))
306  !--------calcul de la seconde composante : u0(2,:) <--> f(r,m,s)
307  fs(2) = sum(ff(j_loc,2) * vv_mesh%gauss%ww(:,l))
308  ft(2) = sum(du_dt(j_loc,2,i) * vv_mesh%gauss%ww(:,l))
309  fp(2) = sum(p(j_loc,2)*dw_loc(1,:))
310  !--------calcul de la troisieme composante : u0(3,:) <--> f(th,m,c)
311  fs(3) = sum(ff(j_loc,3) * vv_mesh%gauss%ww(:,l))
312  ft(3) = sum(du_dt(j_loc,3,i) * vv_mesh%gauss%ww(:,l))
313  fp(3) = sum(p(j_loc,2)*vv_mesh%gauss%ww(:,l))/ray*list_mode(i)
314  !--------calcul de la quatrieme composante : u0(4,:) <--> f(th,m,s)
315  fs(4) = sum(ff(j_loc,4) * vv_mesh%gauss%ww(:,l))
316  ft(4) = sum(du_dt(j_loc,4,i) * vv_mesh%gauss%ww(:,l))
317  fp(4) = -sum(p(j_loc,1)*vv_mesh%gauss%ww(:,l))/ray*list_mode(i)
318  !--------calcul de la cinquieme composante : u0(5,:) <--> f(z,m,c)
319  fs(5) = sum(ff(j_loc,5) * vv_mesh%gauss%ww(:,l))
320  ft(5) = sum(du_dt(j_loc,5,i) * vv_mesh%gauss%ww(:,l))
321  fp(5) = sum(p(j_loc,1)*dw_loc(2,:))
322  !--------calcul de la sixieme composante : u0(6,:) <--> f(z,m,s)
323  fs(6) = sum(ff(j_loc,6) * vv_mesh%gauss%ww(:,l))
324  ft(6) = sum(du_dt(j_loc,6,i) * vv_mesh%gauss%ww(:,l))
325  fp(6) = sum(p(j_loc,2)*dw_loc(2,:))
326 
327  rhs_gauss(index,:,i) = ft+fp-fs-rotb_b(index,:,i)
328  ENDDO
329  ENDDO
330  END DO
331 
332  END SUBROUTINE rhs_residual_ns_gauss_3x3_mom
333 
334 !===Local subroutine
335  SUBROUTINE inject(jj_c, jj_f, pp_c, pp_f)
336  IMPLICIT NONE
337  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj_c, jj_f
338  REAL(KIND=8), DIMENSION(:), INTENT(IN) :: pp_c
339  REAL(KIND=8), DIMENSION(:), INTENT(OUT) :: pp_f
340  REAL(KIND=8) :: half = 0.5
341  INTEGER:: m
342  IF (SIZE(jj_c,1)==3) THEN
343  DO m = 1, SIZE(jj_f,2)
344  pp_f(jj_f(1:3,m)) = pp_c(jj_c(:,m))
345  pp_f(jj_f(4,m)) = (pp_c(jj_c(2,m)) + pp_c(jj_c(3,m)))*half
346  pp_f(jj_f(5,m)) = (pp_c(jj_c(3,m)) + pp_c(jj_c(1,m)))*half
347  pp_f(jj_f(6,m)) = (pp_c(jj_c(1,m)) + pp_c(jj_c(2,m)))*half
348  END DO
349 
350  ELSE
351  DO m = 1, SIZE(jj_f,2)
352  pp_f(jj_f(1:4,m)) = pp_c(jj_c(:,m))
353  END DO
354  pp_f(jj_f(5,:)) = (pp_c(jj_c(3,:)) + pp_c(jj_c(4,:)))*half
355  pp_f(jj_f(6,:)) = (pp_c(jj_c(4,:)) + pp_c(jj_c(2,:)))*half
356  pp_f(jj_f(7,:)) = (pp_c(jj_c(2,:)) + pp_c(jj_c(3,:)))*half
357  pp_f(jj_f(8,:)) = (pp_c(jj_c(1,:)) + pp_c(jj_c(4,:)))*half
358  pp_f(jj_f(9,:)) = (pp_c(jj_c(3,:)) + pp_c(jj_c(1,:)))*half
359  pp_f(jj_f(10,:)) = (pp_c(jj_c(1,:)) + pp_c(jj_c(2,:)))*half
360 
361  END IF
362 
363  END SUBROUTINE inject
364 END MODULE rhs_gauss_computing
real(kind=8) function, dimension(size(rr, 2)), public source_in_ns_momentum(TYPE, rr, mode, i, time, Re, ty, opt_density, opt_tempn)
subroutine, public rhs_residual_ns_gauss_3x3_mom(vv_mesh, pp_mesh, list_mode, time, du_dt, pn, density, rotb_b, rhs_gauss, opt_tempn)
subroutine, public rhs_ns_gauss_3x3(vv_mesh, pp_mesh, communicator, list_mode, time, V1m, pn, pn_inc, rotv_v, rhs_gauss, opt_tempn)
subroutine inject(jj_c, jj_f, pp_c, pp_f)
real(kind=8) function, dimension(nb_angles, ne-nb+1), public penal_in_real_space(mesh, rr_gauss, angles, nb_angles, nb, ne, time)
real(kind=8) function, dimension(size(rr, 2), 6), public imposed_velocity_by_penalty(rr, t)
subroutine, public fft_par_var_eta_prod_gauss_dcl(communicator, eta, H_mesh, c_in, c_out, nb_procs, bloc_size, m_max_pad, rr_gauss, time, temps)
subroutine, public rhs_residual_ns_gauss_3x3(vv_mesh, pp_mesh, communicator, list_mode, time, du_dt, pn, rotv_v, rhs_gauss, opt_tempn)