SFEMaNS  version 5.3
Reference documentation for SFEMaNS
condlim_test_22.f90
Go to the documentation of this file.
2  USE my_util
3  USE def_type_mesh
4  USE input_data
5  USE test_22
6  USE bessel
7 !!$ATTENTION
8 !!$Some subroutines have been commented to avoid warning messages when compiling executable.
9 !!$It can not be done in the module boundary_generic that expects all subroutines to be present.
10 !!$END ATTENTION
11 !!$ PUBLIC :: init_velocity_pressure
12 !!$ PUBLIC :: init_temperature
13 !!$ PUBLIC :: init_level_set
14 !!$ PUBLIC :: source_in_NS_momentum
15 !!$ PUBLIC :: source_in_temperature
16 !!$ PUBLIC :: source_in_level_set
17 !!$ PUBLIC :: vv_exact
18 !!$ PUBLIC :: imposed_velocity_by_penalty
19 !!$ PUBLIC :: pp_exact
20 !!$ PUBLIC :: temperature_exact
21 !!$ PUBLIC :: level_set_exact
22 !!$ PUBLIC :: penal_in_real_space
23 !!$ PUBLIC :: extension_velocity
24  PUBLIC :: vexact
25 !!$ PUBLIC :: H_B_quasi_static
26  PUBLIC :: hexact
27  PUBLIC :: phiexact
28  PUBLIC :: jexact_gauss
29  PUBLIC :: eexact_gauss
30  PUBLIC :: init_maxwell
31  PUBLIC :: mu_bar_in_fourier_space
33  PUBLIC :: mu_in_real_space
34  PRIVATE
35 
36 CONTAINS
37  !===============================================================================
38  ! Boundary conditions for Navier-Stokes
39  !===============================================================================
40 
41 !!$ !===Initialize velocity, pressure
42 !!$ SUBROUTINE init_velocity_pressure(mesh_f, mesh_c, time, dt, list_mode, &
43 !!$ un_m1, un, pn_m1, pn, phin_m1, phin)
44 !!$ IMPLICIT NONE
45 !!$ TYPE(mesh_type) :: mesh_f, mesh_c
46 !!$ REAL(KIND=8), INTENT(OUT):: time
47 !!$ REAL(KIND=8), INTENT(IN) :: dt
48 !!$ INTEGER, DIMENSION(:), INTENT(IN) :: list_mode
49 !!$ REAL(KIND=8), DIMENSION(:,:,:), INTENT(OUT):: un_m1, un
50 !!$ REAL(KIND=8), DIMENSION(:,:,:), INTENT(OUT):: pn_m1, pn, phin_m1, phin
51 !!$ INTEGER :: mode, i, j
52 !!$ REAL(KIND=8), DIMENSION(mesh_c%np) :: pn_m2
53 !!$
54 !!$ time = 0.d0
55 !!$ DO i= 1, SIZE(list_mode)
56 !!$ mode = list_mode(i)
57 !!$ DO j = 1, 6
58 !!$ !===velocity
59 !!$ un_m1(:,j,i) = vv_exact(j,mesh_f%rr,mode,time-dt)
60 !!$ un (:,j,i) = vv_exact(j,mesh_f%rr,mode,time)
61 !!$ END DO
62 !!$ DO j = 1, 2
63 !!$ !===pressure
64 !!$ pn_m2(:) = pp_exact(j,mesh_c%rr,mode,time-2*dt)
65 !!$ pn_m1 (:,j,i) = pp_exact(j,mesh_c%rr,mode,time-dt)
66 !!$ pn (:,j,i) = pp_exact(j,mesh_c%rr,mode,time)
67 !!$ phin_m1(:,j,i) = pn_m1(:,j,i) - pn_m2(:)
68 !!$ phin (:,j,i) = Pn (:,j,i) - pn_m1(:,j,i)
69 !!$ ENDDO
70 !!$ ENDDO
71 !!$ END SUBROUTINE init_velocity_pressure
72 
73 !!$ !===Initialize temperature
74 !!$ SUBROUTINE init_temperature(mesh, time, dt, list_mode, tempn_m1, tempn)
75 !!$ IMPLICIT NONE
76 !!$ TYPE(mesh_type) :: mesh
77 !!$ REAL(KIND=8), INTENT(OUT):: time
78 !!$ REAL(KIND=8), INTENT(IN) :: dt
79 !!$ INTEGER, DIMENSION(:), INTENT(IN) :: list_mode
80 !!$ REAL(KIND=8), DIMENSION(:,:,:), INTENT(OUT):: tempn_m1, tempn
81 !!$ INTEGER :: mode, i, j
82 !!$
83 !!$ time = 0.d0
84 !!$ DO i= 1, SIZE(list_mode)
85 !!$ mode = list_mode(i)
86 !!$ DO j = 1, 2
87 !!$ tempn_m1(:,j,i) = temperature_exact(j, mesh%rr, mode, time-dt)
88 !!$ tempn (:,j,i) = temperature_exact(j, mesh%rr, mode, time)
89 !!$ ENDDO
90 !!$ ENDDO
91 !!$ END SUBROUTINE init_temperature
92 
93 !!$ !===Initialize level_set
94 !!$ SUBROUTINE init_level_set(vv_mesh, time, &
95 !!$ dt, list_mode, level_set_m1, level_set)
96 !!$ IMPLICIT NONE
97 !!$ TYPE(mesh_type) :: vv_mesh
98 !!$ REAL(KIND=8), INTENT(OUT):: time
99 !!$ REAL(KIND=8), INTENT(IN) :: dt
100 !!$ INTEGER, DIMENSION(:), INTENT(IN) :: list_mode
101 !!$ REAL(KIND=8), DIMENSION(:,:,:,:), INTENT(OUT):: level_set, level_set_m1
102 !!$ INTEGER :: mode, i, j, n
103 !!$
104 !!$ time = 0.d0
105 !!$ DO i= 1, SIZE(list_mode)
106 !!$ mode = list_mode(i)
107 !!$ DO j = 1, 2
108 !!$ !===level_set
109 !!$ DO n = 1, inputs%nb_fluid -1
110 !!$ level_set_m1(n,:,j,i) = level_set_exact(n,j,vv_mesh%rr,mode,time-dt)
111 !!$ level_set (n,:,j,i) = level_set_exact(n,j,vv_mesh%rr,mode,time)
112 !!$ END DO
113 !!$ END DO
114 !!$ END DO
115 !!$
116 !!$ END SUBROUTINE init_level_set
117 
118 !!$ !===Source in momemtum equation. Always called.
119 !!$ FUNCTION source_in_NS_momentum(TYPE, rr, mode, i, time, Re, ty, opt_density, opt_tempn) RESULT(vv)
120 !!$ IMPLICIT NONE
121 !!$ INTEGER , INTENT(IN) :: TYPE
122 !!$ REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
123 !!$ INTEGER , INTENT(IN) :: mode, i
124 !!$ REAL(KIND=8), INTENT(IN) :: time
125 !!$ REAL(KIND=8), INTENT(IN) :: Re
126 !!$ CHARACTER(LEN=2), INTENT(IN) :: ty
127 !!$ REAL(KIND=8), DIMENSION(:,:,:), OPTIONAL, INTENT(IN) :: opt_density
128 !!$ REAL(KIND=8), DIMENSION(:,:,:), OPTIONAL, INTENT(IN) :: opt_tempn
129 !!$ REAL(KIND=8), DIMENSION(SIZE(rr,2)) :: vv
130 !!$
131 !!$ vv = 0.d0
132 !!$ CALL error_petsc('source_in_NS_momentum: should not be called for this test')
133 !!$ RETURN
134 !!$ END FUNCTION source_in_NS_momentum
135 
136 !!$ !===Extra source in temperature equation. Always called.
137 !!$ FUNCTION source_in_temperature(TYPE, rr, m, t)RESULT(vv)
138 !!$ IMPLICIT NONE
139 !!$ INTEGER , INTENT(IN) :: TYPE
140 !!$ REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
141 !!$ INTEGER , INTENT(IN) :: m
142 !!$ REAL(KIND=8), INTENT(IN) :: t
143 !!$ REAL(KIND=8), DIMENSION(SIZE(rr,2)) :: vv
144 !!$
145 !!$ vv = 0.d0
146 !!$ CALL error_petsc('source_in_temperature: should not be called for this test')
147 !!$ RETURN
148 !!$ END FUNCTION source_in_temperature
149 
150 !!$ !===Extra source in level set equation. Always called.
151 !!$ FUNCTION source_in_level_set(interface_nb,TYPE, rr, m, t)RESULT(vv)
152 !!$ IMPLICIT NONE
153 !!$ INTEGER , INTENT(IN) :: TYPE
154 !!$ REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
155 !!$ INTEGER , INTENT(IN) :: m, interface_nb
156 !!$ REAL(KIND=8), INTENT(IN) :: t
157 !!$ REAL(KIND=8), DIMENSION(SIZE(rr,2)) :: vv
158 !!$
159 !!$ vv=0.d0
160 !!$ CALL error_petsc('sourece_in_temperature: should not be called for this test')
161 !!$ END FUNCTION source_in_level_set
162 
163 !!$ !===Velocity for boundary conditions in Navier-Stokes.
164 !!$ !===Can be used also to initialize velocity in: init_velocity_pressure_temperature
165 !!$ FUNCTION vv_exact(TYPE,rr,m,t) RESULT(vv)
166 !!$
167 !!$ IMPLICIT NONE
168 !!$ INTEGER , INTENT(IN) :: TYPE
169 !!$ REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
170 !!$ INTEGER, INTENT(IN) :: m
171 !!$ REAL(KIND=8), INTENT(IN) :: t
172 !!$ REAL(KIND=8), DIMENSION(SIZE(rr,2)) :: vv
173 !!$
174 !!$ vv(:) = 0.d0
175 !!$ CALL error_petsc('vv_exact: should not be called for this test')
176 !!$ RETURN
177 !!$ END FUNCTION vv_exact
178 
179 !!$ !===Solid velocity imposed when using penalty technique
180 !!$ !===Defined in Fourier space on mode 0 only.
181 !!$ FUNCTION imposed_velocity_by_penalty(rr,t) RESULT(vv)
182 !!$ IMPLICIT NONE
183 !!$ REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
184 !!$ REAL(KIND=8), INTENT(IN) :: t
185 !!$ REAL(KIND=8), DIMENSION(SIZE(rr,2),6) :: vv
186 !!$
187 !!$ vv=0.d0
188 !!$ RETURN
189 !!$ END FUNCTION imposed_velocity_by_penalty
190 
191 !!$ !===Pressure for boundary conditions in Navier-Stokes.
192 !!$ !===Can be used also to initialize pressure in the subroutine init_velocity_pressure.
193 !!$ !===Use this routine for outflow BCs only.
194 !!$ !===CAUTION: Do not enfore BCs on pressure where normal component
195 !!$ ! of velocity is prescribed.
196 !!$ FUNCTION pp_exact(TYPE,rr,m,t) RESULT (vv)
197 !!$ IMPLICIT NONE
198 !!$ INTEGER , INTENT(IN) :: TYPE
199 !!$ REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
200 !!$ INTEGER , INTENT(IN) :: m
201 !!$ REAL(KIND=8), INTENT(IN) :: t
202 !!$ REAL(KIND=8), DIMENSION(SIZE(rr,2)) :: vv
203 !!$
204 !!$ vv=0.d0
205 !!$ CALL error_petsc('pp_exact: should not be called for this test')
206 !!$ RETURN
207 !!$ END FUNCTION pp_exact
208 
209 !!$ !===Temperature for boundary conditions in temperature equation.
210 !!$ FUNCTION temperature_exact(TYPE,rr,m,t) RESULT (vv)
211 !!$ IMPLICIT NONE
212 !!$ INTEGER , INTENT(IN) :: TYPE
213 !!$ REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
214 !!$ INTEGER , INTENT(IN) :: m
215 !!$ REAL(KIND=8), INTENT(IN) :: t
216 !!$ REAL(KIND=8), DIMENSION(SIZE(rr,2)) :: vv
217 !!$
218 !!$ vv = 0.d0
219 !!$ CALL error_petsc('temperature_exact: should not be called for this test')
220 !!$ RETURN
221 !!$ END FUNCTION temperature_exact
222 
223 !!$ !===Can be used to initialize level set in the subroutine init_level_set.
224 !!$ FUNCTION level_set_exact(interface_nb,TYPE,rr,m,t) RESULT (vv)
225 !!$ IMPLICIT NONE
226 !!$ INTEGER , INTENT(IN) :: TYPE
227 !!$ REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
228 !!$ INTEGER , INTENT(IN) :: m, interface_nb
229 !!$ REAL(KIND=8), INTENT(IN) :: t
230 !!$ REAL(KIND=8), DIMENSION(SIZE(rr,2)) :: vv
231 !!$
232 !!$ vv = 0.d0
233 !!$ CALL error_petsc('level_set_exact: should not be called for this test')
234 !!$ RETURN
235 !!$
236 !!$ END FUNCTION level_set_exact
237 
238 !!$ !===Penalty coefficient (if needed)
239 !!$ !===This coefficient is equal to zero in subdomain
240 !!$ !===where penalty is applied (penalty is zero in solid)
241 !!$ FUNCTION penal_in_real_space(mesh,rr_gauss,angles,nb_angles,nb,ne,time) RESULT(vv)
242 !!$ IMPLICIT NONE
243 !!$ TYPE(mesh_type) :: mesh
244 !!$ REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr_gauss
245 !!$ REAL(KIND=8), DIMENSION(:), INTENT(IN) :: angles
246 !!$ INTEGER, INTENT(IN) :: nb_angles
247 !!$ INTEGER, INTENT(IN) :: nb, ne
248 !!$ REAL(KIND=8), INTENT(IN) :: time
249 !!$ REAL(KIND=8), DIMENSION(nb_angles,ne-nb+1) :: vv
250 !!$
251 !!$ vv = 1.d0
252 !!$ CALL error_petsc('penal_in_real_space: should not be called for this test')
253 !!$ RETURN
254 !!$ END FUNCTION penal_in_real_space
255 
256 !!$ !===Extension of the velocity field in the solid.
257 !!$ !===Used when temperature or Maxwell equations are solved.
258 !!$ !===It extends the velocity field on the Navier-Stokes domain to a
259 !!$ !===velocity field on the temperature and the Maxwell domain.
260 !!$ !===It is also used if problem type=mxw and restart velocity
261 !!$ !===is set to true in data (type problem denoted mxx in the code).
262 !!$ FUNCTION extension_velocity(TYPE, H_mesh, mode, t, n_start) RESULT(vv)
263 !!$ IMPLICIT NONE
264 !!$ TYPE(mesh_type), INTENT(IN) :: H_mesh
265 !!$ INTEGER , INTENT(IN) :: TYPE, n_start
266 !!$ INTEGER, INTENT(IN) :: mode
267 !!$ REAL(KIND=8), INTENT(IN) :: t
268 !!$ REAL(KIND=8), DIMENSION(H_Mesh%np) :: vv
269 !!$
270 !!$ vv = 0.d0
271 !!$ RETURN
272 !!$
273 !!$ END FUNCTION extension_velocity
274 
275  !===============================================================================
276  ! Boundary conditions for Maxwell
277  !===============================================================================
278  !===Velocity used in the induction equation.
279  !===Used only if problem type is mxw and restart velocity is false
280  FUNCTION vexact(m, H_mesh) RESULT(vv) !Set uniquement a l'induction
281  IMPLICIT NONE
282  TYPE(mesh_type), INTENT(IN) :: H_mesh
283  INTEGER, INTENT(IN) :: m
284  REAL(KIND=8), DIMENSION(H_mesh%np,6) :: vv
285  INTEGER :: n
286 
287  vv = 0.d0
288  RETURN
289 
290  !===Dummies variables to avoid warning
291  n=h_mesh%np; n=m
292  !===Dummies variables to avoid warning
293  END FUNCTION vexact
294 
295 !!$ !===Magnetic field and magnetic induction for quasi-static approximation
296 !!$ !===if needed
297 !!$ FUNCTION H_B_quasi_static(char_h_b, rr, m) RESULT(vv)
298 !!$ IMPLICIT NONE
299 !!$ CHARACTER(LEN=1), INTENT(IN) :: char_h_b
300 !!$ REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
301 !!$ INTEGER, INTENT(IN) :: m
302 !!$ REAL(KIND=8), DIMENSION(SIZE(rr,2),6) :: vv
303 !!$
304 !!$ vv = 0.d0
305 !!$ RETURN
306 !!$ END FUNCTION H_B_quasi_static
307 
308  !===Magnetic field for boundary conditions in the Maxwell equations.
309  FUNCTION hexact(H_mesh,TYPE, rr, m, mu_H_field, t) RESULT(vv)
310  IMPLICIT NONE
311  TYPE(mesh_type), INTENT(IN) :: H_mesh
312  INTEGER , INTENT(IN) :: TYPE
313  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
314  INTEGER , INTENT(IN) :: m
315  REAL(KIND=8), INTENT(IN) :: t
316  REAL(KIND=8), DIMENSION(:), INTENT(IN) :: mu_H_field
317  REAL(KIND=8), DIMENSION(SIZE(rr,2)) :: vv
318  INTEGER :: n
319  REAL(KIND=8) :: r
320 
321  IF (m==0) THEN
322  IF (type==1) THEN !Bessel functions defined for reals
323  vv = cosh(rr(2,:))
324  DO n = 1, SIZE(rr,2)
325  vv(n) = -bessj1(rr(1,n))*vv(n)
326  END DO
327  ELSE IF (type==5) THEN
328  vv = sinh(rr(2,:))
329  DO n = 1, SIZE(rr,2)
330  vv(n) = bessj0(rr(1,n))*vv(n)
331  END DO
332  ELSE
333  vv = 0.d0
334  END IF
335  ELSE IF (m==mode_mu_t22) THEN
336  IF (type==1) THEN !Bessel functions defined for reals
337  vv = cosh(rr(2,:))*f_test_t22(rr(1,:),rr(2,:))
338  DO n = 1, SIZE(rr,2)
339  vv(n) = -bessj1(rr(1,n))*vv(n)
340  END DO
341  ELSE IF (type==5) THEN
342  vv = sinh(rr(2,:))*f_test_t22(rr(1,:),rr(2,:))
343  DO n = 1, SIZE(rr,2)
344  vv(n) = bessj0(rr(1,n))*vv(n)
345  END DO
346  ELSE
347  vv = 0.d0
348  END IF
349  ELSE
350  vv(:) = 0.d0
351  END IF
352  RETURN
353 
354  !===Dummies variables to avoid warning
355  n=h_mesh%np; r=mu_h_field(1); r=t
356  !===Dummies variables to avoid warning
357  END FUNCTION hexact
358 
359  !===Scalar potential for boundary conditions in the Maxwell equations.
360  FUNCTION phiexact(TYPE, rr, m, mu_phi,t) RESULT(vv)
361  IMPLICIT NONE
362  INTEGER , INTENT(IN) :: TYPE
363  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
364  INTEGER , INTENT(IN) :: m
365  REAL(KIND=8), INTENT(IN) :: mu_phi, t
366  REAL(KIND=8), DIMENSION(SIZE(rr,2)) :: vv
367  REAL(KIND=8), DIMENSION(SIZE(rr,2)) :: r, z
368  INTEGER :: n
369 
370  IF (m/=0) THEN
371  vv = 0.d0
372  RETURN
373  END IF
374  r = rr(1,:)
375  z = rr(2,:)
376  IF (type==1) THEN !Bessel functions defined for reals
377  DO n = 1, SIZE(rr,2)
378  vv(n) = bessj0(r(n))*cosh(z(n))
379  END DO
380  ELSE
381  vv = 0.d0
382  END IF
383  RETURN
384 
385  !===Dummies variables to avoid warning
386  r=mu_phi; r=t
387  !===Dummies variables to avoid warning
388  END FUNCTION phiexact
389 
390  !===Current in Ohm's law. Curl(H) = sigma(E + uxB) + current
391  FUNCTION jexact_gauss(TYPE, rr, m, mu_phi, sigma, mu_H, t, mesh_id, opt_B_ext) RESULT(vv)
392  IMPLICIT NONE
393  INTEGER , INTENT(IN) :: TYPE
394  REAL(KIND=8), DIMENSION(:), INTENT(IN) :: rr
395  INTEGER , INTENT(IN) :: m
396  REAL(KIND=8), INTENT(IN) :: mu_phi, sigma, mu_H, t
397  INTEGER , INTENT(IN) :: mesh_id
398  REAL(KIND=8), DIMENSION(6), OPTIONAL,INTENT(IN) :: opt_B_ext
399  REAL(KIND=8) :: vv
400  REAL(KIND=8) :: r, z
401  REAL(KIND=8), DIMENSION(1) :: dummy
402  INTEGER :: n
403 
404  IF (m/=mode_mu_t22) THEN
405  vv = 0.d0
406  RETURN
407  END IF
408 
409  r = rr(1)
410  z = rr(2)
411  dummy = f_test_t22(rr(1:1),rr(2:2))
412  IF (type==2) THEN
413  vv = -(mode_mu_t22/r)*dummy(1)*bessj0(r)*sinh(z)
414  ELSE IF (type==3) THEN
415  vv = -dfdr_test_t22(r,z)*bessj0(r)*sinh(z)-dfdz_test_t22(r,z)*bessj1(r)*cosh(z)
416  ELSE IF (type==6) THEN
417  vv = -(mode_mu_t22/r)*dummy(1)*bessj1(r)*cosh(z)
418  ELSE
419  vv = 0.d0
420  END IF
421  RETURN
422 
423  !===Dummies variables to avoid warning
424  r=mu_phi; r=sigma; r=mu_h; r=t; n=mesh_id
425  IF (PRESENT(opt_b_ext)) r=opt_b_ext(1)
426  !===Dummies variables to avoid warning
427  END FUNCTION jexact_gauss
428 
429  !===Electric field for Neumann BC (cf. doc)
430  FUNCTION eexact_gauss(TYPE, rr, m, mu_phi, sigma, mu_H, t) RESULT(vv)
431  IMPLICIT NONE
432  INTEGER, INTENT(IN) :: TYPE
433  REAL(KIND=8), DIMENSION(:), INTENT(IN) :: rr
434  INTEGER, INTENT(IN) :: m
435  REAL(KIND=8), INTENT(IN) :: mu_phi, sigma, mu_H, t
436  REAL(KIND=8) :: vv
437  REAL(KIND=8) :: r
438  INTEGER :: n
439 
440  vv = 0.d0
441  RETURN
442 
443  !===Dummies variables to avoid warning
444  r=rr(1); r=mu_phi; r=sigma; r=mu_h; r=t; n=type; n=m
445  !===Dummies variables to avoid warning
446  END FUNCTION eexact_gauss
447 
448  !===Initialization of magnetic field and scalar potential (if present)
449  SUBROUTINE init_maxwell(H_mesh, phi_mesh, time, dt, mu_H_field, mu_phi, &
450  list_mode, hn1, hn, phin1, phin)
451  IMPLICIT NONE
452  TYPE(mesh_type) :: H_mesh, phi_mesh
453  REAL(KIND=8), INTENT(OUT):: time
454  REAL(KIND=8), INTENT(IN) :: dt
455  REAL(KIND=8), DIMENSION(:), INTENT(IN) :: mu_H_field
456  REAL(KIND=8), INTENT(IN) :: mu_phi
457  INTEGER, DIMENSION(:), INTENT(IN) :: list_mode
458  REAL(KIND=8), DIMENSION(:,:,:), INTENT(OUT):: Hn, Hn1
459  REAL(KIND=8), DIMENSION(:,:,:), INTENT(OUT):: phin, phin1
460  REAL(KIND=8) :: r
461  INTEGER :: n
462 
463  hn1 = 0.d0
464  hn = 0.d0
465  phin1 = 0.d0
466  phin = 0.d0
467  time=0.d0
468  RETURN
469 
470  !===Dummies variables to avoid warning
471  r=h_mesh%rr(1,1); r=phi_mesh%rr(1,1); r=dt; r=mu_h_field(1); r=mu_phi; n=SIZE(list_mode)
472  !===Dummies variables to avoid warning
473  END SUBROUTINE init_maxwell
474 
475  !===Analytical permeability (if needed)
476  !===This function is not needed unless the flag
477  !=== ===Use FEM Interpolation for magnetic permeability (true/false)
478  !===is activated and set to .FALSE. in the data data file. Default is .TRUE.
479  FUNCTION mu_bar_in_fourier_space(H_mesh,nb,ne,pts,pts_ids) RESULT(vv)
480  IMPLICIT NONE
481  TYPE(mesh_type), INTENT(IN) :: H_mesh
482  REAL(KIND=8), DIMENSION(ne-nb+1) :: vv
483  INTEGER, INTENT(IN) :: nb, ne
484  REAL(KIND=8),DIMENSION(2,ne-nb+1),OPTIONAL :: pts
485  INTEGER, DIMENSION(ne-nb+1), OPTIONAL :: pts_ids
486 
487  IF( PRESENT(pts) .AND. PRESENT(pts_ids) ) THEN
488  vv=mu_bar_in_fourier_space_anal_t22(h_mesh,nb,ne,pts,pts_ids)
489  ELSE
490  vv=mu_bar_in_fourier_space_anal_t22(h_mesh,nb,ne,pts)
491  END IF
492  RETURN
493  END FUNCTION mu_bar_in_fourier_space
494 
495  !===Analytical mu_in_fourier_space (if needed)
496  !===This function is not needed unless the flag
497  !=== ===Use FEM Interpolation for magnetic permeability (true/false)
498  !===is activated and set to .FALSE. in the data data file. Default is .TRUE.
499  FUNCTION grad_mu_bar_in_fourier_space(pt,pt_id) RESULT(vv)
500  IMPLICIT NONE
501  REAL(KIND=8),DIMENSION(2), INTENT(in):: pt
502  INTEGER,DIMENSION(1), INTENT(in) :: pt_id
503  REAL(KIND=8),DIMENSION(2) :: vv
504 
506  RETURN
507  END FUNCTION grad_mu_bar_in_fourier_space
508 
509  !===Analytical permeability, mu in real space (if needed)
510  FUNCTION mu_in_real_space(H_mesh,angles,nb_angles,nb,ne,time) RESULT(vv)
511  IMPLICIT NONE
512  TYPE(mesh_type), INTENT(IN) :: H_mesh
513  REAL(KIND=8), DIMENSION(:), INTENT(IN) :: angles
514  INTEGER, INTENT(IN) :: nb_angles
515  INTEGER, INTENT(IN) :: nb, ne
516  REAL(KIND=8), INTENT(IN) :: time
517  REAL(KIND=8), DIMENSION(nb_angles,ne-nb+1) :: vv
518  REAL(KIND=8) :: r
519 
520  vv = mu_in_real_space_anal_t22(h_mesh,angles,nb_angles,nb,ne)
521  RETURN
522 
523  !===Dummies variables to avoid warning
524  r=time
525  !===Dummies variables to avoid warning
526  END FUNCTION mu_in_real_space
527 
528 END MODULE boundary_test_22
real(kind=8) function, dimension(nb_angles, ne-nb+1), public mu_in_real_space(H_mesh, angles, nb_angles, nb, ne, time)
real(kind=8) function, dimension(nb_angles, ne-nb+1) mu_in_real_space_anal_t22(H_mesh, angles, nb_angles, nb, ne)
Definition: test_22.f90:89
Definition: bessel.f90:4
real(kind=8) function, dimension(ne-nb+1), public mu_bar_in_fourier_space(H_mesh, nb, ne, pts, pts_ids)
integer, public mode_mu_t22
Definition: test_22.f90:6
real(kind=8) function, public bessj0(x)
Definition: bessel.f90:206
real(kind=8) function, dimension(size(r)) f_test_t22(r, z)
Definition: test_22.f90:11
real(kind=8) function, dimension(size(rr, 2)), public phiexact(TYPE, rr, m, mu_phi, t)
real(kind=8) function dfdz_test_t22(r, z)
Definition: test_22.f90:27
real(kind=8) function, dimension(h_mesh%np, 6), public vexact(m, H_mesh)
real(kind=8) function, dimension(ne-nb+1) mu_bar_in_fourier_space_anal_t22(H_mesh, nb, ne, pts, pts_ids)
Definition: test_22.f90:36
real(kind=8) function, dimension(2), public grad_mu_bar_in_fourier_space(pt, pt_id)
real(kind=8) function, dimension(size(rr, 2)), public hexact(H_mesh, TYPE, rr, m, mu_H_field, t)
real(kind=8) function dfdr_test_t22(r, z)
Definition: test_22.f90:19
real(kind=8) function, public eexact_gauss(TYPE, rr, m, mu_phi, sigma, mu_H, t)
real(kind=8) function, public jexact_gauss(TYPE, rr, m, mu_phi, sigma, mu_H, t, mesh_id, opt_B_ext)
real(kind=8) function, dimension(2) grad_mu_bar_in_fourier_space_anal_t22(pt, pt_id)
Definition: test_22.f90:61
real(kind=8) function, public bessj1(x)
Definition: bessel.f90:236
subroutine, public init_maxwell(H_mesh, phi_mesh, time, dt, mu_H_field, mu_phi, list_mode, Hn1, Hn, phin1, phin)