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