MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
mod_twofl_hllc.t
Go to the documentation of this file.
2 #include "amrvac.h"
5  use mod_functions_bfield, only: mag
6 
7  implicit none
8  private
9 
10 
11  public :: twofl_hllc_init
12 
13  integer :: i_rho, i_e
14  integer, allocatable :: i_mom(:)
15 
16 
17 contains
18 
19  subroutine twofl_hllc_init()
22 
23  allocate(i_mom(1:ndir))
24 
26 
27  end subroutine twofl_hllc_init
28 
29  subroutine twofl_hllc_init_species(ii, rho_, mom, e_)
31  integer, intent(in) :: ii
32  integer, intent(out) :: rho_, mom(1:ndir),e_
33 
34  if(ii==1) then
38  i_rho = rho_c_
39  i_mom(1:ndir) = mom_c(1:ndir)
40  i_e = e_c_
41  else
45  i_rho = rho_n_
46  i_mom(1:ndir) = mom_n(1:ndir)
47  i_e = e_n_
48  end if
49 
50  rho_ = i_rho
51  mom(:) = i_mom(:)
52  e_ = i_e
53 
54  end subroutine twofl_hllc_init_species
55 
56  subroutine twofl_diffuse_hllcd_n(ixI^L,ixO^L,idim,wLC,wRC,fLC,fRC,patchf)
57 
58  ! when method is hllcd or hllcd1 then:
59 
60  ! this subroutine is to enforce regions where we AVOID HLLC
61  ! and use TVDLF instead: this is achieved by setting patchf to 4 in
62  ! certain regions. An additional input parameter is nxdiffusehllc
63  ! which sets the size of the fallback region.
64 
66 
67  integer, intent(in) :: ixI^L,ixO^L,idim
68  double precision, dimension(ixI^S,1:nw), intent(in) :: wRC,wLC
69  double precision, dimension(ixI^S,1:nwflux),intent(in) :: fLC, fRC
70  integer , dimension(ixI^S), intent(inout) :: patchf
71 
72  integer :: ixOO^D,TxOO^L
73 
74  ! In a user-controlled region around any point with flux sign change between
75  ! left and right, ensure fallback to TVDLF
76  {do ixoo^d= ixo^lim^d\}
77  {
78  txoomin^d= max(ixoo^d - nxdiffusehllc*kr(idim,^d), ixomin^d);
79  txoomax^d= min(ixoo^d + nxdiffusehllc*kr(idim,^d), ixomax^d);
80  \}
81  if(abs(patchf(ixoo^d)) == 1 .or. abs(patchf(ixoo^d)) == 4)Then
82  if(any(frc(ixoo^d,1:nwflux)*flc(ixoo^d,1:nwflux)<-smalldouble))Then
83  where(abs(patchf(txoo^s))==1)
84  patchf(txoo^s) = 4
85  endwhere
86  endif
87  endif
88  {enddo^d&\}
89 
90  end subroutine twofl_diffuse_hllcd_n
91 
92  subroutine twofl_get_lcd_n(wLC,wRC,fLC,fRC,cmin,cmax,idim,ixI^L,ixO^L, &
93  whll,Fhll,lambdaCD,patchf)
94  ! Calculate lambda at CD and set the patchf to know the orientation
95  ! of the riemann fan and decide on the flux choice
96  ! We also compute here the HLL flux and w value, for fallback strategy
97 
98  ! In this nul version, we simply compute nothing and ensure TVDLF fallback
100 
101  integer, intent(in) :: ixI^L,ixO^L,idim
102  double precision, dimension(ixI^S,1:nw), intent(in) :: wLC,wRC
103  double precision, dimension(ixI^S,1:nwflux), intent(in) :: fLC,fRC
104  double precision, dimension(ixI^S), intent(in) :: cmax,cmin
105  integer , dimension(ixI^S), intent(inout) :: patchf
106  double precision, dimension(ixI^S,1:nwflux), intent(out) :: Fhll,whll
107  double precision, dimension(ixI^S), intent(out) :: lambdaCD
108 
109  logical , dimension(ixI^S) :: Cond_patchf
110  double precision :: Epsilon
111  integer :: iw,ix^D
112 
113  !--------------------------------------------
114  ! on entry, patch is preset to contain values from -2,1,2,4
115  ! -2: take left flux, no computation here
116  ! +2: take right flux, no computation here
117  ! +4: take TVDLF flux, no computation here
118  ! 1: compute the characteristic speed for the CD
119 
120  cond_patchf(ixo^s)=(abs(patchf(ixo^s))==1)
121  lambdacd=0.d0
122 
123  do iw=1,nwflux
124  where(cond_patchf(ixo^s))
125  !============= compute HLL flux ==============!
126  fhll(ixo^s,iw)= (cmax(ixo^s)*flc(ixo^s,iw)-cmin(ixo^s)*frc(ixo^s,iw) &
127  + cmin(ixo^s)*cmax(ixo^s)*(wrc(ixo^s,iw)-wlc(ixo^s,iw)))&
128  /(cmax(ixo^s)-cmin(ixo^s))
129  !======== compute intermediate HLL state =======!
130  whll(ixo^s,iw) = (cmax(ixo^s)*wrc(ixo^s,iw)-cmin(ixo^s)*wlc(ixo^s,iw)&
131  +flc(ixo^s,iw)-frc(ixo^s,iw))/(cmax(ixo^s)-cmin(ixo^s))
132  end where
133  end do
134 
135  ! deduce the characteristic speed at the CD
136  where(cond_patchf(ixo^s))
137  lambdacd(ixo^s)=whll(ixo^s,i_mom(idim))/whll(ixo^s,i_rho)
138  end where
139 
140  {do ix^db=ixomin^db,ixomax^db\}
141  if(cond_patchf(ix^d)) then
142  ! double check whether obtained speed is in between min and max speeds given
143  ! and identify in which part of the Riemann fan the time-axis is
144  if(cmin(ix^d) < zero .and. lambdacd(ix^d)>zero&
145  .and.lambdacd(ix^d)<cmax(ix^d)) then
146  patchf(ix^d) = -1
147  else if(cmax(ix^d) > zero .and. lambdacd(ix^d) < zero&
148  .and.lambdacd(ix^d)>cmin(ix^d)) then
149  patchf(ix^d) = 1
150  else if(lambdacd(ix^d) >= cmax(ix^d) .or. &
151  lambdacd(ix^d) <= cmin(ix^d)) then
152  lambdacd(ix^d) = zero
153  ! we will fall back to HLL flux case in this degeneracy
154  patchf(ix^d) = 3
155  end if
156  end if
157  {end do\}
158 
159  where(patchf(ixo^s)== 3)
160  cond_patchf(ixo^s)=.false.
161  end where
162 
163  ! handle the specific case where the time axis is exactly on the CD
164  if(any(cond_patchf(ixo^s).and.lambdacd(ixo^s)==zero))then
165  ! determine which sector (forward or backward) of the Riemann fan is smallest
166  ! and select left or right flux accordingly
167  {do ix^db=ixomin^db,ixomax^db\}
168  if(lambdacd(ix^d)==zero.and.cond_patchf(ix^d)) then
169  if(-cmin(ix^d)>=cmax(ix^d)) then
170  patchf(ix^d) = 1
171  else
172  patchf(ix^d) = -1
173  end if
174  end if
175  {end do\}
176  end if
177 
178  ! eigenvalue lambda for contact is near zero: decrease noise by this trick
179  if(flathllc)then
180  epsilon=1.0d-6
181  where(cond_patchf(ixo^s).and. &
182  dabs(lambdacd(ixo^s))/max(cmax(ixo^s),epsilon)< epsilon .and. &
183  dabs(lambdacd(ixo^s))/max(dabs(cmin(ixo^s)),epsilon)< epsilon)
184  lambdacd(ixo^s) = zero
185  end where
186  end if
187 
188  end subroutine twofl_get_lcd_n
189 
190  subroutine twofl_get_wcd_n(wLC,wRC,whll,fRC,fLC,Fhll,patchf,lambdaCD,cmin,cmax,&
191  ixI^L,ixO^L,idim,f)
192  ! compute the intermediate state U*
193  ! only needed where patchf=-1/1
194 
195  ! reference Li S., JCP, 203, 2005, 344-357
196  ! reference T. Miyoski, Kusano JCP, 2008, 2005.
197 
199  use mod_physics
200 
201  integer, intent(in) :: ixI^L,ixO^L,idim
202  double precision, dimension(ixI^S,1:nw), intent(in) :: wRC,wLC
203  double precision, dimension(ixI^S,1:nwflux), intent(in):: whll, Fhll
204  double precision, dimension(ixI^S), intent(in) :: lambdaCD
205  double precision, dimension(ixI^S), intent(in) :: cmax,cmin
206  double precision, dimension(ixI^S,1:nwflux), intent(in):: fRC,fLC
207  double precision, dimension(ixI^S,1:nwflux),intent(out):: f
208  double precision, dimension(ixI^S,1:nw) :: wCD,wSub
209  double precision, dimension(ixI^S,1:nwflux) :: fSub
210  double precision, dimension(ixI^S) :: vSub,cspeed,pCD
211  integer , dimension(ixI^S), intent(in) :: patchf
212 
213  double precision :: csmls
214  integer :: n, iw, ix^D
215 
216  !-------------- auxiliary Speed and array-------------!
217  {do ix^db=ixomin^db,ixomax^db\}
218  if(patchf(ix^d)==1) then
219  cspeed(ix^d)=cmax(ix^d)
220  vsub(ix^d)=wrc(ix^d,i_mom(idim))/wrc(ix^d,i_rho)
221  wsub(ix^d,:)=wrc(ix^d,:)
222  fsub(ix^d,:)=frc(ix^d,:)
223  else if(patchf(ix^d)==-1) then
224  cspeed(ix^d)=cmin(ix^d)
225  vsub(ix^d)=wlc(ix^d,i_mom(idim))/wlc(ix^d,i_rho)
226  wsub(ix^d,:)=wlc(ix^d,:)
227  fsub(ix^d,:)=flc(ix^d,:)
228  end if
229  {end do\}
230 
231  {do ix^db=ixomin^db,ixomax^db\}
232  if(abs(patchf(ix^d))==1) then
233  csmls=one/(cspeed(ix^d)-lambdacd(ix^d))
234  wcd(ix^d,i_rho) = wsub(ix^d,i_rho)&
235  *(cspeed(ix^d)-vsub(ix^d))*csmls
236 
237  !------- Momentum ------!
238  do iw=1, ndir
239  if(iw /= idim)then
240  ! eq. 21 22
241  wcd(ix^d,i_mom(iw))=(cspeed(ix^d)*wsub(ix^d,i_mom(iw))-fsub(ix^d,i_mom(iw)))*csmls
242  else
243  ! eq. 20
244  wcd(ix^d,i_mom(iw)) = wcd(ix^d,i_rho) * lambdacd(ix^d)
245  endif
246  enddo
247  if(phys_energy) then
248  pcd(ix^d) = wsub(ix^d,i_rho)*(cspeed(ix^d)-vsub(ix^d))&
249  *(lambdacd(ix^d)-vsub(ix^d))&
250  +fsub(ix^d,i_mom(idim))-wsub(ix^d,i_mom(idim))*vsub(ix^d)
251  ! Eq 31
252  wcd(ix^d,i_e) = (cspeed(ix^d)*wsub(ix^d,i_e) &
253  -fsub(ix^d,i_e)+lambdacd(ix^d)*pcd(ix^d))*csmls
254  end if
255  do iw=1,nwflux
256  ! f_i=fsub+lambda (wCD-wSub)
257  f(ix^d,iw)=fsub(ix^d,iw)+cspeed(ix^d)*(wcd(ix^d,iw)-wsub(ix^d,iw))
258  end do
259  end if
260  {end do\}
261 
262  end subroutine twofl_get_wcd_n
263 
264  subroutine twofl_diffuse_hllcd_c(ixI^L,ixO^L,idim,wLC,wRC,fLC,fRC,patchf)
265  ! when method is hllcd or hllcd1 then:
266  ! this subroutine is to enforce regions where we AVOID HLLC
267  ! and use TVDLF instead: this is achieved by setting patchf to 4 in
268  ! certain regions. An additional input parameter is nxdiffusehllc
269  ! which sets the size of the fallback region.
271 
272  integer, intent(in) :: ixI^L,ixO^L,idim
273  double precision, dimension(ixI^S,1:nw), intent(in) :: wRC,wLC
274  double precision, dimension(ixI^S,1:nwflux),intent(in) :: fLC, fRC
275  integer, dimension(ixI^S), intent(inout) :: patchf
276 
277  integer :: ixOO^D,TxOO^L
278 
279 
280  ! In a user-controlled region around any point with flux sign change between
281  ! left and right, ensure fallback to TVDLF
282  {do ixoo^d= ixo^lim^d\}
283  {
284  txoomin^d= max(ixoo^d - nxdiffusehllc*kr(idim,^d), ixomin^d);
285  txoomax^d= min(ixoo^d + nxdiffusehllc*kr(idim,^d), ixomax^d);
286  \}
287  if(abs(patchf(ixoo^d)) == 1 .or. abs(patchf(ixoo^d)) == 4)Then
288  if(any(frc(ixoo^d,1:nwflux)*flc(ixoo^d,1:nwflux)<-smalldouble))Then
289  where(abs(patchf(txoo^s))==1)
290  patchf(txoo^s) = 4
291  endwhere
292  endif
293  endif
294  {enddo^d&\}
295 
296  end subroutine twofl_diffuse_hllcd_c
297 
298  subroutine twofl_get_lcd_c(wLC,wRC,fLC,fRC,cmin,cmax,idim,ixI^L,ixO^L, &
299  whll,Fhll,lambdaCD,patchf)
300 
301  ! Calculate lambda at CD and set the patchf to know the orientation
302  ! of the riemann fan and decide on the flux choice
303  ! We also compute here the HLL flux and w value, for fallback strategy
304 
306 
307  integer, intent(in) :: ixI^L,ixO^L,idim
308  double precision, dimension(ixI^S,1:nw), intent(in) :: wLC,wRC
309  double precision, dimension(ixI^S,1:nwflux), intent(in):: fLC,fRC
310  double precision, dimension(ixI^S), intent(in) :: cmax,cmin
311  integer , dimension(ixI^S), intent(inout) :: patchf
312  double precision, dimension(ixI^S,1:nwflux), intent(out) :: Fhll,whll
313  double precision, dimension(ixI^S), intent(out) :: lambdaCD
314 
315  logical , dimension(ixI^S) :: Cond_patchf
316  double precision :: Epsilon
317  integer :: iw
318 
319  ! on entry, patch is preset to contain values from -2,1,2,4
320  ! -2: take left flux, no computation here
321  ! +2: take right flux, no computation here
322  ! +4: take TVDLF flux, no computation here
323  ! 1: compute the characteristic speed for the CD
324 
325  cond_patchf(ixo^s)=(abs(patchf(ixo^s))==1)
326 
327  do iw=1,nwflux
328  where(cond_patchf(ixo^s))
329  !============= compute HLL flux ==============!
330  fhll(ixo^s,iw)= (cmax(ixo^s)*flc(ixo^s,iw)-cmin(ixo^s)*frc(ixo^s,iw) &
331  + cmin(ixo^s)*cmax(ixo^s)*(wrc(ixo^s,iw)-wlc(ixo^s,iw)))&
332  /(cmax(ixo^s)-cmin(ixo^s))
333  !======== compute intermediate HLL state =======!
334  whll(ixo^s,iw) = (cmax(ixo^s)*wrc(ixo^s,iw)-cmin(ixo^s)*wlc(ixo^s,iw)&
335  +flc(ixo^s,iw)-frc(ixo^s,iw))/(cmax(ixo^s)-cmin(ixo^s))
336  endwhere
337  enddo
338 
339  ! deduce the characteristic speed at the CD
340  where(cond_patchf(ixo^s))
341  lambdacd(ixo^s)=whll(ixo^s,i_mom(idim))/whll(ixo^s,i_rho)
342  end where
343 
344 
345  where(cond_patchf(ixo^s))
346  ! double check whether obtained speed is in between min and max speeds given
347  ! and identify in which part of the Riemann fan the time-axis is
348  where(cmin(ixo^s)<zero.and.lambdacd(ixo^s)>zero&
349  .and.lambdacd(ixo^s)<cmax(ixo^s))
350  patchf(ixo^s) = -1
351  elsewhere(cmax(ixo^s)>zero.and.lambdacd(ixo^s)<zero&
352  .and.lambdacd(ixo^s)>cmin(ixo^s))
353  patchf(ixo^s) = 1
354  elsewhere(lambdacd(ixo^s)>=cmax(ixo^s).or.lambdacd(ixo^s) <= cmin(ixo^s))
355  lambdacd(ixo^s) = zero
356  ! we will fall back to HLL flux case in this degeneracy
357  patchf(ixo^s) = 3
358  endwhere
359  endwhere
360 
361 
362  where(patchf(ixo^s)== 3)
363  cond_patchf(ixo^s)=.false.
364  end where
365 
366 
367  ! handle the specific case where the time axis is exactly on the CD
368  if(any(lambdacd(ixo^s)==zero.and.cond_patchf(ixo^s)))then
369  ! determine which sector (forward or backward) of the Riemann fan is smallest
370  ! and select left or right flux accordingly
371  where(lambdacd(ixo^s)==zero.and.cond_patchf(ixo^s))
372  where(-cmin(ixo^s)>=cmax(ixo^s))
373  patchf(ixo^s) = 1
374  elsewhere
375  patchf(ixo^s) = -1
376  endwhere
377  endwhere
378  endif
379 
380  ! eigenvalue lambda for contact is near zero: decrease noise by this trick
381  if(flathllc)then
382  epsilon=1.0d-6
383  where(cond_patchf(ixo^s).and. &
384  dabs(lambdacd(ixo^s))/max(cmax(ixo^s),epsilon)< epsilon .and. &
385  dabs(lambdacd(ixo^s))/max(dabs(cmin(ixo^s)),epsilon)< epsilon)
386  lambdacd(ixo^s) = zero
387  end where
388  end if
389 
390  end subroutine twofl_get_lcd_c
391 
392  subroutine twofl_get_wcd_c(wLC,wRC,whll,fRC,fLC,Fhll,patchf,lambdaCD,cmin,cmax,&
393  ixI^L,ixO^L,idim,f)
394  ! compute the intermediate state U*
395  ! only needed where patchf=-1/1
396 
397  ! reference Li S., JCP, 203, 2005, 344-357
398  ! reference T. Miyoski, Kusano JCP, 2008, 2005.
400  use mod_physics
401  integer, intent(in) :: ixI^L,ixO^L,idim
402  double precision, dimension(ixI^S,1:nw), intent(in) :: wRC,wLC
403  double precision, dimension(ixI^S,1:nwflux), intent(in):: whll, Fhll
404  double precision, dimension(ixI^S), intent(in) :: lambdaCD
405  double precision, dimension(ixI^S), intent(in) :: cmax,cmin
406  double precision, dimension(ixI^S,1:nwflux), intent(in):: fRC,fLC
407  double precision, dimension(ixI^S,1:nwflux),intent(out):: f
408  double precision, dimension(ixI^S,1:nw) :: wCD,wSub
409  double precision, dimension(ixI^S,1:nwflux) :: fSub
410  double precision, dimension(ixI^S) :: vSub,cspeed,pCD,VdotBCD
411  integer , dimension(ixI^S), intent(in) :: patchf
412 
413  integer :: n, iw, idir,ix^D
414  double precision, dimension(ixI^S) :: rho
415 
416 
417  !-------------- auxiliary Speed and array-------------!
418  {do ix^db=ixomin^db,ixomax^db\}
419  if(patchf(ix^d)==1) then
420  cspeed(ix^d)=cmax(ix^d)
421  vsub(ix^d)=wrc(ix^d,i_mom(idim))/wrc(ix^d,i_rho)
422  wsub(ix^d,:)=wrc(ix^d,:)
423  fsub(ix^d,:)=frc(ix^d,:)
424  else if(patchf(ix^d)==-1) then
425  cspeed(ix^d)=cmin(ix^d)
426  vsub(ix^d)=wlc(ix^d,i_mom(idim))/wlc(ix^d,i_rho)
427  wsub(ix^d,:)=wlc(ix^d,:)
428  fsub(ix^d,:)=flc(ix^d,:)
429  end if
430  {end do\}
431 
432  {do ix^db=ixomin^db,ixomax^db\}
433  if(abs(patchf(ix^d))==1) then
434  wcd(ix^d,i_rho) = wsub(ix^d,i_rho)&
435  *(cspeed(ix^d)-vsub(ix^d))/(cspeed(ix^d)-lambdacd(ix^d))
436  !TODO tracer
437 ! do n=1,twofl_n_tracer
438 ! iw = tracer(n)
439 ! wCD(ix^D,iw) = wSub(ix^D,iw)*(cspeed(ix^D)-vSub(ix^D))&
440 ! /(cspeed(ix^D)-lambdaCD(ix^D))
441 ! end do
442  !==== Magnetic field ====!
443  do idir=1,ndir
444  ! case from eq 31
445  wcd(ix^d,mag(idir)) = whll(ix^d,mag(idir))
446  end do
447  !------- Momentum ------!
448  do iw=1, ndir
449  if(iw /= idim)then
450  ! eq. 21 22
451  wcd(ix^d,i_mom(iw))=(cspeed(ix^d)*wsub(ix^d,i_mom(iw))-fsub(ix^d,i_mom(iw))&
452  -wcd(ix^d,mag(idim))*wcd(ix^d,mag(iw)))/&
453  (cspeed(ix^d)-lambdacd(ix^d))
454  else
455  ! eq. 20
456  wcd(ix^d,i_mom(iw)) = wcd(ix^d,i_rho) * lambdacd(ix^d)
457  endif
458  enddo
459  if(phys_energy) then
460  vdotbcd(ix^d) = sum(whll(ix^d,i_mom(:))*whll(ix^d,mag(:)))/whll(ix^d,i_rho)
461  ! Eq 17
462  pcd(ix^d) = wsub(ix^d,i_rho)*(cspeed(ix^d)-vsub(ix^d))&
463  *(lambdacd(ix^d)-vsub(ix^d))&
464  +fsub(ix^d,i_mom(idim))-wsub(ix^d,i_mom(idim))*vsub(ix^d)&
465  + wcd(ix^d,mag(idim))**2
466  ! Eq 31
467  wcd(ix^d,i_e) = (cspeed(ix^d)*wsub(ix^d,i_e) &
468  -fsub(ix^d,i_e)+lambdacd(ix^d)*pcd(ix^d)&
469  -vdotbcd(ix^d)*wcd(ix^d,mag(idim)))&
470  /(cspeed(ix^d)-lambdacd(ix^d))
471  end if
472  end if
473  {end do\}
474 
475  do iw=1,nwflux
476  if(iw == mag(idim)) then
477  f(ixo^s,iw)=zero
478  else if(twofl_glm .and. iw == psi_) then
479  f(ixo^s,iw)=zero
480  else
481  where(abs(patchf(ixo^s))==1)
482  ! f_i=fsub+lambda (wCD-wSub)
483  f(ixo^s,iw)=fsub(ixo^s,iw)+cspeed(ixo^s)*(wcd(ixo^s,iw)-wsub(ixo^s,iw))
484  endwhere
485  end if
486  end do
487 
488  end subroutine twofl_get_wcd_c
489 
490 end module mod_twofl_hllc
integer, dimension(:), allocatable, public mag
Indices of the magnetic field.
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, dimension(3, 3) kr
Kronecker delta tensor.
double precision, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
procedure(sub_get_wcd), pointer phys_get_wcd
procedure(sub_hllc_init_species), pointer phys_hllc_init_species
procedure(sub_get_lcd), pointer phys_get_lcd
procedure(sub_diffuse_hllcd), pointer phys_diffuse_hllcd
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition: mod_physics.t:4
subroutine twofl_get_wcd_c(wLC, wRC, whll, fRC, fLC, Fhll, patchf, lambdaCD, cmin, cmax, ixIL, ixOL, idim, f)
subroutine twofl_get_wcd_n(wLC, wRC, whll, fRC, fLC, Fhll, patchf, lambdaCD, cmin, cmax, ixIL, ixOL, idim, f)
subroutine twofl_hllc_init_species(ii, rho_, mom, e_)
subroutine twofl_get_lcd_n(wLC, wRC, fLC, fRC, cmin, cmax, idim, ixIL, ixOL, whll, Fhll, lambdaCD, patchf)
subroutine twofl_get_lcd_c(wLC, wRC, fLC, fRC, cmin, cmax, idim, ixIL, ixOL, whll, Fhll, lambdaCD, patchf)
subroutine, public twofl_hllc_init()
subroutine twofl_diffuse_hllcd_c(ixIL, ixOL, idim, wLC, wRC, fLC, fRC, patchf)
subroutine twofl_diffuse_hllcd_n(ixIL, ixOL, idim, wLC, wRC, fLC, fRC, patchf)
Magneto-hydrodynamics module.
Definition: mod_twofl_phys.t:2
integer, public e_n_
integer, public e_c_
Index of the energy density (-1 if not present)
integer, dimension(:), allocatable, public mom_c
Indices of the momentum density.
integer, dimension(:), allocatable, public mom_n
integer, public rho_c_
Index of the density (in the w array)
integer, public rho_n_