MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_viscosity.t
Go to the documentation of this file.
1 !> The module add viscous source terms and check time step
2 !>
3 !> Viscous forces in the momentum equations:
4 !> d m_i/dt += - div (vc_mu * PI)
5 !> !! Viscous work in the energy equation:
6 !> !! de/dt += - div (v . vc_mu * PI)
7 !> where the PI stress tensor is
8 !> PI_i,j = - (dv_j/dx_i + dv_i/dx_j) + (2/3)*Sum_k dv_k/dx_k
9 !> where vc_mu is the dynamic viscosity coefficient (g cm^-1 s^-1).
11  use mod_comm_lib, only: mpistop
12  implicit none
13 
14  !> Viscosity coefficient
15  double precision, public :: vc_mu = 1.d0
16 
17  !> fourth order
18  logical :: vc_4th_order = .false.
19 
20  !> source split or not
21  logical :: vc_split= .false.
22 
23  !> whether to compute the viscous terms as
24  !> fluxes (ie in the div on the LHS), or not (by default)
25  logical :: viscindiv= .false.
26 
27  !> Index of the density (in the w array)
28  integer, private, parameter :: rho_ = 1
29 
30  !> Indices of the momentum density
31  integer, allocatable, private, protected :: mom(:)
32 
33  !> Index of the energy density (-1 if not present)
34  integer, private, protected :: e_
35 
36  ! Public methods
37  public :: visc_get_flux_prim
38 
39 contains
40  !> Read this module"s parameters from a file
41  subroutine vc_params_read(files)
43  character(len=*), intent(in) :: files(:)
44  integer :: n
45 
46  namelist /vc_list/ vc_mu, vc_4th_order, vc_split, viscindiv
47 
48  do n = 1, size(files)
49  open(unitpar, file=trim(files(n)), status="old")
50  read(unitpar, vc_list, end=111)
51 111 close(unitpar)
52  end do
53 
54  end subroutine vc_params_read
55 
56  !> Initialize the module
57  subroutine viscosity_init(phys_wider_stencil,phys_req_diagonal)
59  integer, intent(inout) :: phys_wider_stencil
60  logical, intent(inout) :: phys_req_diagonal
61  integer :: nwx,idir
62 
64 
65  if(vc_split) any_source_split=.true.
66 
67  ! Determine flux variables
68  nwx = 1 ! rho (density)
69 
70  allocate(mom(ndir))
71  do idir = 1, ndir
72  nwx = nwx + 1
73  mom(idir) = nwx ! momentum density
74  end do
75 
76  nwx = nwx + 1
77  e_ = nwx ! energy density
78 
79  if (viscindiv) then
80  ! to compute the derivatives from left and right upwinded values
81  phys_wider_stencil = 1
82  phys_req_diagonal = .true. ! viscInDiv
83  end if
84 
85  end subroutine viscosity_init
86 
87  subroutine viscosity_add_source(qdt,ixI^L,ixO^L,wCT,w,x,&
88  energy,qsourcesplit,active)
89  ! Add viscosity source in isotropic Newtonian fluids to w within ixO
90  ! neglecting bulk viscosity
91  ! dm/dt= +div(mu*[d_j v_i+d_i v_j]-(2*mu/3)* div v * kr)
93  use mod_geometry
94 
95  integer, intent(in) :: ixI^L, ixO^L
96  double precision, intent(in) :: qdt, x(ixI^S,1:ndim), wCT(ixI^S,1:nw)
97  double precision, intent(inout) :: w(ixI^S,1:nw)
98  logical, intent(in) :: energy,qsourcesplit
99  logical, intent(inout) :: active
100 
101  integer:: ix^L,idim,idir,jdir,iw
102  double precision:: lambda(ixI^S,ndir,ndir),tmp(ixI^S),tmp2(ixI^S),v(ixI^S,ndir),vlambda(ixI^S,ndir)
103 
104  if (viscindiv) return
105 
106  if(qsourcesplit .eqv. vc_split) then
107  active = .true.
108  ! standard case, textbook viscosity
109  ! Calculating viscosity sources
110  if(.not.vc_4th_order) then
111  ! involves second derivatives, two extra layers
112  ix^l=ixo^l^ladd2;
113  if({ iximin^d>ixmin^d .or. iximax^d<ixmax^d|.or.})&
114  call mpistop("error for viscous source addition, 2 layers needed")
115  ix^l=ixo^l^ladd1;
116  else
117  ! involves second derivatives, four extra layers
118  ix^l=ixo^l^ladd4;
119  if({ iximin^d>ixmin^d .or. iximax^d<ixmax^d|.or.})&
120  call mpistop("error for viscous source addition"//&
121  "requested fourth order gradients: 4 layers needed")
122  ix^l=ixo^l^ladd2;
123  end if
124 
125  ! get velocity
126  do idir=1,ndir
127  v(ixi^s,idir)=wct(ixi^s,mom(idir))/wct(ixi^s,rho_)
128  end do
129 
130  ! construct lambda tensor: lambda_ij = gradv_ij + gradv_ji
131  ! initialize
132  lambda=zero
133 
134  !next construct
135  do idim=1,ndim; do idir=1,ndir
136  ! Calculate velocity gradient tensor within ixL: gradv= grad v,
137  ! thus gradv_ij=d_j v_i
138  tmp(ixi^s)=v(ixi^s,idir)
139  ! Correction for Christoffel terms in non-cartesian
140  if (coordinate==cylindrical .and. idim==r_ .and. idir==phi_ ) tmp(ixi^s) = tmp(ixi^s)/x(ixi^s,1)
141  if (coordinate==spherical) then
142  if (idim==r_ .and. (idir==2 .or. idir==phi_)) then
143  tmp(ixi^s) = tmp(ixi^s)/x(ixi^s,1)
144 {^nooned
145  elseif (idim==2 .and. idir==phi_) then
146  tmp(ixi^s)=tmp(ixi^s)/dsin(x(ixi^s,2))
147 }
148  endif
149  endif
150  call gradient(tmp,ixi^l,ix^l,idim,tmp2)
151  ! Correction for Christoffel terms in non-cartesian
152  if (coordinate==cylindrical .and. idim==r_ .and. idir==phi_ ) tmp2(ix^s)=tmp2(ix^s)*x(ix^s,1)
153  if (coordinate==cylindrical .and. idim==phi_ .and. idir==phi_ ) tmp2(ix^s)=tmp2(ix^s)+v(ix^s,r_)/x(ix^s,1)
154  if (coordinate==spherical) then
155  if (idim==r_ .and. (idir==2 .or. idir==phi_)) then
156  tmp2(ix^s) = tmp2(ix^s)*x(ix^s,1)
157 {^nooned
158  elseif (idim==2 .and. idir==phi_ ) then
159  tmp2(ix^s)=tmp2(ix^s)*dsin(x(ix^s,2))
160  elseif (idim==2 .and. idir==2 ) then
161  tmp2(ix^s)=tmp2(ix^s)+v(ix^s,r_)/x(ix^s,1)
162  elseif (idim==phi_.and. idir==phi_) then
163  tmp2(ix^s)=tmp2(ix^s)+v(ix^s,r_)/x(ix^s,1)+v(ix^s,2)/(x(ix^s,1)*dtan(x(ix^s,2)))
164 }
165  endif
166  endif
167  lambda(ix^s,idim,idir)= lambda(ix^s,idim,idir)+ tmp2(ix^s)
168  lambda(ix^s,idir,idim)= lambda(ix^s,idir,idim)+ tmp2(ix^s)
169  enddo; enddo;
170 
171  ! Multiply lambda with viscosity coefficient and dt
172  lambda(ix^s,1:ndir,1:ndir)=lambda(ix^s,1:ndir,1:ndir)*vc_mu*qdt
173 
174  !calculate div v term through trace action separately
175  ! rq : it is safe to use the trace rather than compute the divergence
176  ! since we always retrieve the divergence (even with the
177  ! Christoffel terms)
178  tmp=0.d0
179  do idir=1,ndir
180  tmp(ix^s)=tmp(ix^s)+lambda(ix^s,idir,idir)
181  end do
182  tmp(ix^s)=tmp(ix^s)/3.d0
183 
184  !substract trace from diagonal elements
185  do idir=1,ndir
186  lambda(ix^s,idir,idir)=lambda(ix^s,idir,idir)-tmp(ix^s)
187  enddo
188 
189  ! dm/dt= +div(mu*[d_j v_i+d_i v_j]-(2*mu/3)* div v * kr)
190  ! hence m_j=m_j+d_i tensor_ji
191  do idir=1,ndir
192  do idim=1,ndim
193  tmp(ix^s)=lambda(ix^s,idir,idim)
194  ! Correction for divergence of a tensor
195  if (coordinate==cylindrical .and. idim==r_ .and. (idir==r_ .or. idir==z_)) tmp(ix^s) = tmp(ix^s)*x(ix^s,1)
196  if (coordinate==cylindrical .and. idim==r_ .and. idir==phi_ ) tmp(ix^s) = tmp(ix^s)*x(ix^s,1)**two
197  if (coordinate==spherical) then
198  if (idim==r_ .and. idir==r_ ) tmp(ix^s) = tmp(ix^s)*x(ix^s,1)**two
199  if (idim==r_ .and. (idir==2 .or. idir==phi_)) tmp(ix^s) = tmp(ix^s)*x(ix^s,1)**3.d0
200 {^nooned
201  if (idim==2 .and. (idir==r_ .or. idir==2)) tmp(ix^s) = tmp(ix^s)*dsin(x(ix^s,2))
202  if (idim==2 .and. idir==phi_ ) tmp(ix^s) = tmp(ix^s)*dsin(x(ix^s,2))**two
203 }
204  endif
205  call gradient(tmp,ixi^l,ixo^l,idim,tmp2)
206  ! Correction for divergence of a tensor
207  if (coordinate==cylindrical .and. idim==r_ .and. (idir==r_ .or. idir==z_)) tmp2(ixo^s) = tmp2(ixo^s)/x(ixo^s,1)
208  if (coordinate==cylindrical .and. idim==r_ .and. idir==phi_ ) tmp2(ixo^s) = tmp2(ixo^s)/(x(ixo^s,1)**two)
209  if (coordinate==spherical) then
210  if (idim==r_ .and. idir==r_ ) tmp2(ixo^s) = tmp2(ixo^s)/(x(ixo^s,1)**two)
211  if (idim==r_ .and. (idir==2 .or. idir==phi_)) tmp2(ixo^s) = tmp2(ixo^s)/(x(ixo^s,1)**3.d0)
212 {^nooned
213  if (idim==2 .and. (idir==r_ .or. idir==2)) tmp2(ixo^s) = tmp2(ixo^s)/(dsin(x(ixo^s,2)))
214  if (idim==2 .and. idir==phi_ ) tmp2(ixo^s) = tmp2(ixo^s)/(dsin(x(ixo^s,2))**two)
215 }
216  endif
217  w(ixo^s,mom(idir))=w(ixo^s,mom(idir))+tmp2(ixo^s)
218  enddo
219  ! Correction for geometrical terms in the div of a tensor
220  if (coordinate==cylindrical .and. idir==r_ ) w(ixo^s,mom(idir))=w(ixo^s,mom(idir))-lambda(ixo^s,phi_,phi_)/x(ixo^s,1)
221  if (coordinate==spherical .and. idir==r_ ) w(ixo^s,mom(idir))=w(ixo^s,mom(idir))-(lambda(ixo^s,2,2)+lambda(ixo^s,phi_,phi_))/x(ixo^s,1)
222 {^nooned
223  if (coordinate==spherical .and. idir==2 ) w(ixo^s,mom(idir))=w(ixo^s,mom(idir))-lambda(ixo^s,phi_,phi_)/(x(ixo^s,1)/dtan(x(ixo^s,2)))
224 }
225  end do
226 
227  if(energy) then
228  ! de/dt= +div(v.dot.[mu*[d_j v_i+d_i v_j]-(2*mu/3)* div v *kr])
229  ! thus e=e+d_i v_j tensor_ji
230  vlambda=0.d0
231  do idim=1,ndim
232  do idir=1,ndir
233  vlambda(ixi^s,idim)=vlambda(ixi^s,idim)+v(ixi^s,idir)*lambda(ixi^s,idir,idim)
234  end do
235  end do
236  call divvector(vlambda,ixi^l,ixo^l,tmp2)
237  w(ixo^s,e_)=w(ixo^s,e_)+tmp2(ixo^s)
238  end if
239  end if
240 
241  end subroutine viscosity_add_source
242 
243  subroutine viscosity_get_dt(w,ixI^L,ixO^L,dtnew,dx^D,x)
244  ! Check diffusion time limit for dt < dtdiffpar * dx**2 / (mu/rho)
246 
247  integer, intent(in) :: ixI^L, ixO^L
248  double precision, intent(in) :: dx^D, x(ixI^S,1:ndim)
249  double precision, intent(in) :: w(ixI^S,1:nw)
250  double precision, intent(inout) :: dtnew
251 
252  double precision :: tmp(ixI^S)
253  double precision:: dtdiff_visc, dxinv2(1:ndim)
254  integer:: idim
255 
256  ! Calculate the kinematic viscosity tmp=mu/rho
257 
258  tmp(ixo^s)=vc_mu/w(ixo^s,rho_)
259 
260  ^d&dxinv2(^d)=one/dx^d**2;
261  do idim=1,ndim
262  dtdiff_visc=dtdiffpar/maxval(tmp(ixo^s)*dxinv2(idim))
263  ! limit the time step
264  dtnew=min(dtnew,dtdiff_visc)
265  enddo
266 
267  end subroutine viscosity_get_dt
268 
269  ! viscInDiv
270  ! Get the viscous stress tensor terms in the idim direction
271  ! Beware : a priori, won't work for ndir /= ndim
272  ! Rq : we work with primitive w variables here
273  ! Rq : ixO^L is already extended by 1 unit in the direction we work on
274 
275  subroutine visc_get_flux_prim(w, x, ixI^L, ixO^L, idim, f, energy)
277  use mod_geometry
278  integer, intent(in) :: ixi^l, ixo^l, idim
279  double precision, intent(in) :: w(ixi^s, 1:nw), x(ixi^s, 1:^nd)
280  double precision, intent(inout) :: f(ixi^s, nwflux)
281  logical, intent(in) :: energy
282  integer :: idir, i
283  double precision :: v(ixi^s,1:ndir)
284 
285  double precision :: divergence(ixi^s)
286 
287  double precision:: lambda(ixi^s,ndir) !, tmp(ixI^S) !gradV(ixI^S,ndir,ndir)
288 
289  if (.not. viscindiv) return
290 
291  do i=1,ndir
292  v(ixi^s,i)=w(ixi^s,i+1)
293  enddo
294  call divvector(v,ixi^l,ixo^l,divergence)
295 
296  call get_crossgrad(ixi^l,ixo^l,x,w,idim,lambda)
297  lambda(ixo^s,idim) = lambda(ixo^s,idim) - (2.d0/3.d0) * divergence(ixo^s)
298 
299  ! Compute the idim-th row of the viscous stress tensor
300  do idir = 1, ndir
301  f(ixo^s, mom(idir)) = f(ixo^s, mom(idir)) - vc_mu * lambda(ixo^s,idir)
302  if (energy) f(ixo^s, e_) = f(ixo^s, e_) - vc_mu * lambda(ixo^s,idir) * v(ixi^s,idir)
303  enddo
304 
305  end subroutine visc_get_flux_prim
306 
307  ! Compute the cross term ( d_i v_j + d_j v_i in Cartesian BUT NOT IN
308  ! CYLINDRICAL AND SPHERICAL )
309  subroutine get_crossgrad(ixI^L,ixO^L,x,w,idim,cross)
311  use mod_geometry
312  integer, intent(in) :: ixI^L, ixO^L, idim
313  double precision, intent(in) :: w(ixI^S, 1:nw), x(ixI^S, 1:ndim)
314  double precision, intent(out) :: cross(ixI^S,ndir)
315  integer :: idir
316  double precision :: tmp(ixI^S), v(ixI^S)
317 
318  if (ndir/=ndim) call mpistop("This formula are probably wrong for ndim/=ndir")
319  ! Beware also, we work w/ the angle as the 3rd component in cylindrical
320  ! and the colatitude as the 2nd one in spherical
321  cross(ixi^s,:)=zero
322  tmp(ixi^s)=zero
323  select case(coordinate)
325  call cart_cross_grad(ixi^l,ixo^l,x,w,idim,cross)
326  case (cylindrical)
327  if (idim==1) then
328  ! for rr and rz
329  call cart_cross_grad(ixi^l,ixo^l,x,w,idim,cross)
330  ! then we overwrite rth w/ the correct expression
331  {^nooned
332  v(ixi^s)=w(ixi^s,mom(1)) ! v_r
333  call gradient(v,ixi^l,ixo^l,2,tmp) ! d_th (rq : already contains 1/r)
334  cross(ixi^s,2)=tmp(ixi^s)
335  v(ixi^s)=w(ixi^s,mom(2))/x(ixi^s,1) ! v_th / r
336  call gradient(v,ixi^l,ixo^l,1,tmp) ! d_r
337  cross(ixi^s,2)=cross(ixi^s,2)+tmp(ixi^s)*x(ixi^s,1)
338  }
339  elseif (idim==2) then
340  ! thr (idem as above)
341  v(ixi^s)=w(ixi^s,mom(1)) ! v_r
342  {^nooned
343  call gradient(v,ixi^l,ixo^l,2,tmp) ! d_th
344  cross(ixi^s,1)=tmp(ixi^s)
345  v(ixi^s)=w(ixi^s,mom(2))/x(ixi^s,1) ! v_th / r
346  call gradient(v,ixi^l,ixo^l,1,tmp) ! d_r
347  cross(ixi^s,1)=cross(ixi^s,1)+tmp(ixi^s)*x(ixi^s,1)
348  ! thth
349  v(ixi^s)=w(ixi^s,mom(2)) ! v_th
350  call gradient(v,ixi^l,ixo^l,2,tmp) ! d_th
351  cross(ixi^s,2)=two*tmp(ixi^s)
352  v(ixi^s)=w(ixi^s,mom(1)) ! v_r
353  cross(ixi^s,2)=cross(ixi^s,2)+two*v(ixi^s)/x(ixi^s,1) ! + 2 vr/r
354  !thz
355  v(ixi^s)=w(ixi^s,mom(3)) ! v_z
356  call gradient(v,ixi^l,ixo^l,2,tmp) ! d_th
357  }
358  cross(ixi^s,3)=tmp(ixi^s)
359  v(ixi^s)=w(ixi^s,mom(2)) ! v_th
360  {^ifthreed
361  call gradient(v,ixi^l,ixo^l,3,tmp) ! d_z
362  }
363  cross(ixi^s,3)=cross(ixi^s,3)+tmp(ixi^s)
364  {^ifthreed
365  elseif (idim==3) then
366  ! for zz and rz
367  call cart_cross_grad(ixi^l,ixo^l,x,w,idim,cross)
368  ! then we overwrite zth w/ the correct expression
369  !thz
370  v(ixi^s)=w(ixi^s,mom(3)) ! v_z
371  call gradient(v,ixi^l,ixo^l,2,tmp) ! d_th
372  cross(ixi^s,2)=tmp(ixi^s)
373  v(ixi^s)=w(ixi^s,mom(2)) ! v_th
374  call gradient(v,ixi^l,ixo^l,3,tmp) ! d_z
375  cross(ixi^s,2)=cross(ixi^s,2)+tmp(ixi^s)
376  }
377  endif
378  case (spherical)
379  if (idim==1) then
380  ! rr (normal, simply 2 * dr vr)
381  v(ixi^s)=w(ixi^s,mom(1)) ! v_r
382  call gradient(v,ixi^l,ixo^l,1,tmp) ! d_r
383  cross(ixi^s,1)=two*tmp(ixi^s)
384  !rth
385  v(ixi^s)=w(ixi^s,mom(1)) ! v_r
386  {^nooned
387  call gradient(v,ixi^l,ixo^l,2,tmp) ! d_th (rq : already contains 1/r)
388  cross(ixi^s,2)=tmp(ixi^s)
389  v(ixi^s)=w(ixi^s,mom(2))/x(ixi^s,1) ! v_th / r
390  call gradient(v,ixi^l,ixo^l,1,tmp) ! d_r
391  cross(ixi^s,2)=cross(ixi^s,2)+tmp(ixi^s)*x(ixi^s,1)
392  }
393  {^ifthreed
394  !rph
395  v(ixi^s)=w(ixi^s,mom(1)) ! v_r
396  call gradient(v,ixi^l,ixo^l,3,tmp) ! d_phi (rq : contains 1/rsin(th))
397  cross(ixi^s,3)=tmp(ixi^s)
398  v(ixi^s)=w(ixi^s,mom(3))/x(ixi^s,1) ! v_phi / r
399  call gradient(v,ixi^l,ixo^l,1,tmp) ! d_r
400  cross(ixi^s,3)=cross(ixi^s,3)+tmp(ixi^s)*x(ixi^s,1)
401  }
402  elseif (idim==2) then
403  ! thr
404  v(ixi^s)=w(ixi^s,mom(1)) ! v_r
405  {^nooned
406  call gradient(v,ixi^l,ixo^l,2,tmp) ! d_th (rq : already contains 1/r)
407  cross(ixi^s,1)=tmp(ixi^s)
408  v(ixi^s)=w(ixi^s,mom(2))/x(ixi^s,1) ! v_th / r
409  call gradient(v,ixi^l,ixo^l,1,tmp) ! d_r
410  cross(ixi^s,1)=cross(ixi^s,1)+tmp(ixi^s)*x(ixi^s,1)
411  ! thth
412  v(ixi^s)=w(ixi^s,mom(2)) ! v_th
413  call gradient(v,ixi^l,ixo^l,2,tmp) ! d_th
414  cross(ixi^s,2)=two*tmp(ixi^s)
415  v(ixi^s)=w(ixi^s,mom(1)) ! v_r
416  cross(ixi^s,2)=cross(ixi^s,2)+two*v(ixi^s)/x(ixi^s,1) ! + 2 vr/r
417  }
418  {^ifthreed
419  !thph
420  v(ixi^s)=w(ixi^s,mom(2)) ! v_th
421  call gradient(v,ixi^l,ixo^l,3,tmp) ! d_phi (rq : contains 1/rsin(th))
422  cross(ixi^s,3)=tmp(ixi^s)
423  v(ixi^s)=w(ixi^s,mom(3))/dsin(x(ixi^s,2)) ! v_ph / sin(th)
424  call gradient(v,ixi^l,ixo^l,2,tmp) ! d_th
425  cross(ixi^s,3)=cross(ixi^s,3)+tmp(ixi^s)*dsin(x(ixi^s,2))
426  }
427  {^ifthreed
428  elseif (idim==3) then
429  !phr
430  v(ixi^s)=w(ixi^s,mom(1)) ! v_r
431  call gradient(v,ixi^l,ixo^l,3,tmp) ! d_phi
432  cross(ixi^s,1)=tmp(ixi^s)
433  v(ixi^s)=w(ixi^s,mom(3))/x(ixi^s,1) ! v_phi / r
434  call gradient(v,ixi^l,ixo^l,1,tmp) ! d_r
435  cross(ixi^s,1)=cross(ixi^s,1)+tmp(ixi^s)*x(ixi^s,1)
436  !phth
437  v(ixi^s)=w(ixi^s,mom(2)) ! v_th
438  call gradient(v,ixi^l,ixo^l,3,tmp) ! d_phi
439  cross(ixi^s,2)=tmp(ixi^s)
440  v(ixi^s)=w(ixi^s,mom(3))/dsin(x(ixi^s,2)) ! v_ph / sin(th)
441  call gradient(v,ixi^l,ixo^l,2,tmp) ! d_th
442  cross(ixi^s,2)=cross(ixi^s,2)+tmp(ixi^s)*dsin(x(ixi^s,2))
443  !phph
444  v(ixi^s)=w(ixi^s,mom(3)) ! v_ph
445  call gradient(v,ixi^l,ixo^l,3,tmp) ! d_phi
446  cross(ixi^s,3)=two*tmp(ixi^s)
447  v(ixi^s)=w(ixi^s,mom(1)) ! v_r
448  cross(ixi^s,3)=cross(ixi^s,3)+two*v(ixi^s)/x(ixi^s,1) ! + 2 vr/r
449  v(ixi^s)=w(ixi^s,mom(2)) ! v_th
450  cross(ixi^s,3)=cross(ixi^s,3)+two*v(ixi^s)/(x(ixi^s,1)*dtan(x(ixi^s,2))) ! + 2 vth/(rtan(th))
451  }
452  endif
453  case default
454  call mpistop("Unknown geometry specified")
455  end select
456 
457  end subroutine get_crossgrad
458 
459  !> yields d_i v_j + d_j v_i for a given i, OK in Cartesian and for some
460  !> tensor terms in cylindrical (rr & rz) and in spherical (rr)
461  subroutine cart_cross_grad(ixI^L,ixO^L,x,w,idim,cross)
463  use mod_geometry
464  integer, intent(in) :: ixI^L, ixO^L, idim
465  double precision, intent(in) :: w(ixI^S, 1:nw), x(ixI^S, 1:^ND)
466  double precision, intent(out) :: cross(ixI^S,ndir)
467  integer :: idir
468  double precision :: tmp(ixI^S), v(ixI^S)
469 
470  v(ixi^s)=w(ixi^s,mom(idim))
471  do idir=1,ndir
472  call gradient(v,ixi^l,ixo^l,idir,tmp)
473  cross(ixo^s,idir)=tmp(ixo^s)
474  enddo
475  do idir=1,ndir
476  v(ixi^s)=w(ixi^s,mom(idir))
477  call gradient(v,ixi^l,ixo^l,idim,tmp)
478  cross(ixo^s,idir)=cross(ixo^s,idir)+tmp(ixo^s)
479  enddo
480 
481  end subroutine cart_cross_grad
482 
483  subroutine visc_add_source_geom(qdt, ixI^L, ixO^L, wCT, w, x)
485  use mod_geometry
486  ! w and wCT conservative variables here
487  integer, intent(in) :: ixI^L, ixO^L
488  double precision, intent(in) :: qdt, x(ixI^S, 1:ndim)
489  double precision, intent(inout) :: wCT(ixI^S, 1:nw), w(ixI^S, 1:nw)
490  ! to change and to set as a parameter in the parfile once the possibility to
491  ! solve the equations in an angular momentum conserving form has been
492  ! implemented (change tvdlf.t eg)
493  double precision :: v(ixI^S,1:ndir), vv(ixI^S), divergence(ixI^S)
494  double precision :: tmp(ixI^S),tmp1(ixI^S)
495  integer :: i
496 
497  if (.not. viscindiv) return
498 
499  select case (coordinate)
500  case (cylindrical)
501  ! get the velocity components
502  do i=1,ndir
503  v(ixi^s,i)=wct(ixi^s,mom(i))/wct(ixi^s,rho_)
504  enddo
505  ! thth tensor term - - -
506  ! 1st the cross grad term
507 {^nooned
508  vv(ixi^s)=v(ixi^s,2) ! v_th
509  call gradient(vv,ixi^l,ixo^l,2,tmp1) ! d_th
510  tmp(ixo^s)=two*(tmp1(ixo^s)+v(ixi^s,1)/x(ixo^s,1)) ! 2 ( d_th v_th / r + vr/r )
511  ! 2nd the divergence
512  call divvector(v,ixi^l,ixo^l,divergence)
513  tmp(ixo^s) = tmp(ixo^s) - (2.d0/3.d0) * divergence(ixo^s)
514  ! s[mr]=-thth/radius
515  w(ixo^s,mom(1))=w(ixo^s,mom(1))-qdt*vc_mu*tmp(ixo^s)/x(ixo^s,1)
516  ! rth tensor term - - -
517  vv(ixi^s)=v(ixi^s,1) ! v_r
518  call gradient(vv,ixi^l,ixo^l,2,tmp1) ! d_th
519  tmp(ixo^s)=tmp1(ixo^s)
520  vv(ixi^s)=v(ixi^s,2)/x(ixi^s,1) ! v_th / r
521  call gradient(vv,ixi^l,ixo^l,1,tmp1) ! d_r
522  tmp(ixo^s)=tmp(ixo^s)+tmp1(ixo^s)*x(ixo^s,1)
523  ! s[mphi]=+rth/radius
524  w(ixo^s,mom(2))=w(ixo^s,mom(2))+qdt*vc_mu*tmp(ixo^s)/x(ixo^s,1)
525 }
526  case (spherical)
527  ! get the velocity components
528  do i=1,ndir
529  v(ixi^s,i)=wct(ixi^s,mom(i))/wct(ixi^s,rho_)
530  enddo
531  ! thth tensor term - - -
532  ! 1st the cross grad term
533  vv(ixi^s)=v(ixi^s,2) ! v_th
534 {^nooned
535  call gradient(vv,ixi^l,ixo^l,2,tmp1) ! d_th
536  tmp(ixo^s)=two*(tmp1(ixo^s)+v(ixo^s,1)/x(ixo^s,1)) ! 2 ( 1/r * d_th v_th + vr/r )
537  ! 2nd the divergence
538  call divvector(v,ixi^l,ixo^l,divergence)
539  tmp(ixo^s) = tmp(ixo^s) - (2.d0/3.d0) * divergence(ixo^s)
540  ! s[mr]=-thth/radius
541  w(ixo^s,mom(1))=w(ixo^s,mom(1))-qdt*vc_mu*tmp(ixo^s)/x(ixo^s,1)
542 }
543  ! phiphi tensor term - - -
544  ! 1st the cross grad term
545  vv(ixi^s)=v(ixi^s,3) ! v_ph
546 {^ifthreed
547  call gradient(vv,ixi^l,ixo^l,3,tmp1) ! d_phi
548  tmp(ixo^s)=two*(tmp1(ixo^s)+v(ixo^s,1)/x(ixo^s,1)+v(ixo^s,2)/(x(ixo^s,1)*dtan(x(ixo^s,2)))) ! 2 ( 1/rsinth * d_ph v_ph + vr/r + vth/rtanth )
549 }
550  ! 2nd the divergence
551  tmp(ixo^s) = tmp(ixo^s) - (2.d0/3.d0) * divergence(ixo^s)
552  ! s[mr]=-phiphi/radius
553  w(ixo^s,mom(1))=w(ixo^s,mom(1))-qdt*vc_mu*tmp(ixo^s)/x(ixo^s,1)
554  ! s[mth]=-cotanth*phiphi/radius
555 {^nooned
556  w(ixo^s,mom(2))=w(ixo^s,mom(2))-qdt*vc_mu*tmp(ixo^s)/(x(ixo^s,1)*dtan(x(ixo^s,2)))
557 }
558  ! rth tensor term - - -
559  vv(ixi^s)=v(ixi^s,1) ! v_r
560  call gradient(vv,ixi^l,ixo^l,2,tmp) ! d_th (rq : already contains 1/r)
561  vv(ixi^s)=v(ixi^s,2)/x(ixi^s,1) ! v_th / r
562  call gradient(vv,ixi^l,ixo^l,1,tmp1) ! d_r
563  tmp(ixo^s)=tmp(ixo^s)+tmp1(ixo^s)*x(ixo^s,1)
564  ! s[mth]=+rth/radius
565  w(ixo^s,mom(2))=w(ixo^s,mom(2))+qdt*vc_mu*tmp(ixo^s)/x(ixo^s,1)
566  ! rphi tensor term - - -
567  vv(ixi^s)=v(ixi^s,1) ! v_r
568 {^ifthreed
569  call gradient(vv,ixi^l,ixo^l,3,tmp) ! d_phi (rq : contains 1/rsin(th))
570 }
571  vv(ixi^s)=v(ixi^s,3)/x(ixi^s,1) ! v_phi / r
572  call gradient(vv,ixi^l,ixo^l,1,tmp1) ! d_r
573  tmp(ixo^s)=tmp(ixo^s)+tmp1(ixo^s)*x(ixo^s,1)
574  ! s[mphi]=+rphi/radius
575  w(ixo^s,mom(3))=w(ixo^s,mom(3))+qdt*vc_mu*tmp(ixo^s)/x(ixo^s,1)
576  ! phith tensor term - - -
577  vv(ixi^s)=v(ixi^s,2) ! v_th
578 {^ifthreed
579  call gradient(vv,ixi^l,ixo^l,3,tmp) ! d_phi
580 }
581 {^nooned
582  vv(ixi^s)=v(ixi^s,3)/dsin(x(ixi^s,2)) ! v_ph / sin(th)
583  call gradient(vv,ixi^l,ixo^l,2,tmp1) ! d_th
584  tmp(ixo^s)=tmp(ixo^s)+tmp1(ixo^s)*dsin(x(ixo^s,2))
585  ! s[mphi]=+cotanth*phith/radius
586  w(ixo^s,mom(3))=w(ixo^s,mom(3))+qdt*vc_mu*tmp(ixo^s)/(x(ixo^s,1)*dtan(x(ixo^s,2)))
587 }
588  end select
589 
590  end subroutine visc_add_source_geom
591 
592 end module mod_viscosity
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: mod_comm_lib.t:208
Module with geometry-related routines (e.g., divergence, curl)
Definition: mod_geometry.t:2
integer coordinate
Definition: mod_geometry.t:7
integer, parameter spherical
Definition: mod_geometry.t:11
integer, parameter cartesian
Definition: mod_geometry.t:8
integer, parameter cylindrical
Definition: mod_geometry.t:10
integer, parameter cartesian_stretched
Definition: mod_geometry.t:9
subroutine gradient(q, ixIL, ixOL, idir, gradq)
Calculate gradient of a scalar q within ixL in direction idir.
Definition: mod_geometry.t:321
subroutine divvector(qvec, ixIL, ixOL, divq, fourthorder, sixthorder)
Calculate divergence of a vector qvec within ixL.
Definition: mod_geometry.t:516
This module contains definitions of global parameters and variables and some generic functions/subrou...
double precision dtdiffpar
For resistive MHD, the time step is also limited by the diffusion time: .
integer, parameter unitpar
file handle for IO
logical any_source_split
if any normal source term is added in split fasion
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
integer r_
Indices for cylindrical coordinates FOR TESTS, negative value when not used:
The module add viscous source terms and check time step.
Definition: mod_viscosity.t:10
subroutine get_crossgrad(ixIL, ixOL, x, w, idim, cross)
subroutine viscosity_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
subroutine viscosity_add_source(qdt, ixIL, ixOL, wCT, w, x, energy, qsourcesplit, active)
Definition: mod_viscosity.t:89
subroutine viscosity_init(phys_wider_stencil, phys_req_diagonal)
Initialize the module.
Definition: mod_viscosity.t:58
subroutine cart_cross_grad(ixIL, ixOL, x, w, idim, cross)
yields d_i v_j + d_j v_i for a given i, OK in Cartesian and for some tensor terms in cylindrical (rr ...
subroutine visc_add_source_geom(qdt, ixIL, ixOL, wCT, w, x)
double precision, public vc_mu
Viscosity coefficient.
Definition: mod_viscosity.t:15
logical vc_4th_order
fourth order
Definition: mod_viscosity.t:18
logical vc_split
source split or not
Definition: mod_viscosity.t:21
subroutine, public visc_get_flux_prim(w, x, ixIL, ixOL, idim, f, energy)
subroutine vc_params_read(files)
Read this module"s parameters from a file.
Definition: mod_viscosity.t:42
logical viscindiv
whether to compute the viscous terms as fluxes (ie in the div on the LHS), or not (by default)
Definition: mod_viscosity.t:25