MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
mod_usr_methods.t
Go to the documentation of this file.
1 !> Module with all the methods that users can customize in AMRVAC
2 !>
3 !> Each procedure pointer can be initialized in a user's mod_usr.t
5 
6  implicit none
7  public
8 
9  !> Initialize the user's settings (after initializing amrvac)
10  procedure(p_no_args), pointer :: usr_set_parameters => null()
11  !> Initialize earch grid block data
12  procedure(init_one_grid), pointer :: usr_init_one_grid => null()
13 
14  ! Boundary condition related
15  procedure(special_bc), pointer :: usr_special_bc => null()
16  procedure(special_mg_bc), pointer :: usr_special_mg_bc => null()
17 
18  procedure(internal_bc), pointer :: usr_internal_bc => null()
19 
20  ! Output related
21  procedure(p_no_args), pointer :: usr_print_log => null()
22  procedure(p_no_args), pointer :: usr_write_analysis => null()
23  procedure(transform_w), pointer :: usr_transform_w => null()
24  procedure(aux_output), pointer :: usr_aux_output => null()
25  procedure(add_aux_names), pointer :: usr_add_aux_names => null()
26  procedure(sub_modify_io), pointer :: usr_modify_output => null()
27  procedure(special_convert), pointer :: usr_special_convert => null()
28 
29  ! Called at the beginning of every time step (after determining dt)
30  procedure(process_grid), pointer :: usr_process_grid => null()
31  procedure(process_global), pointer :: usr_process_global => null()
32 
33  ! Called every time step just after advance (with w^(n+1), it^n, global_time^n)
34  procedure(process_adv_grid), pointer :: usr_process_adv_grid => null()
35  procedure(process_adv_global), pointer :: usr_process_adv_global => null()
36 
37  ! Called after initial condition before the start of the simulation
38  procedure(p_no_args), pointer :: usr_improve_initial_condition => null()
39 
40  ! Called before the start of the simulation
41  procedure(p_no_args), pointer :: usr_before_main_loop => null()
42 
43  ! Source terms
44  procedure(source), pointer :: usr_source => null()
45  procedure(get_dt), pointer :: usr_get_dt => null()
46  procedure(phys_gravity), pointer :: usr_gravity => null()
47 
48  ! Usr defined dust drag force
49  procedure(phys_dust_get_dt), pointer :: usr_dust_get_dt => null()
50  procedure(phys_dust_get_3d_dragforce), pointer :: usr_get_3d_dragforce => null()
51 
52  ! Usr defined space varying viscosity
53  procedure(phys_visco), pointer :: usr_setvisco => null()
54 
55  ! Usr defined thermal pressure for hydro & energy=.False.
56  procedure(hd_pthermal), pointer :: usr_set_pthermal => null()
57 
58  ! Refinement related procedures
59  procedure(refine_grid), pointer :: usr_refine_grid => null()
60  procedure(var_for_errest), pointer :: usr_var_for_errest => null()
61  procedure(a_refine_threshold), pointer :: usr_refine_threshold => null()
62  procedure(flag_grid), pointer :: usr_flag_grid => null()
63 
64  ! Set time-independent magnetic field for B0 splitting
65  procedure(set_b0), pointer :: usr_set_b0 => null()
66  ! Set time-independent variables for equilibrium splitting, except for B0
67  procedure(set_equi_vars), pointer :: usr_set_equi_vars => null()
68  ! Set time-independent current density for B0 splitting
69  procedure(set_j0), pointer :: usr_set_j0 => null()
70  procedure(special_resistivity), pointer :: usr_special_resistivity => null()
71 
72  ! Particle module related
73  procedure(update_payload), pointer :: usr_update_payload => null()
74  procedure(create_particles), pointer :: usr_create_particles => null()
75  procedure(check_particle), pointer :: usr_check_particle => null()
76  procedure(particle_fields), pointer :: usr_particle_fields => null()
77  procedure(particle_analytic), pointer :: usr_particle_analytic => null()
78  procedure(particle_position), pointer :: usr_particle_position => null()
79 
80  ! Radiation quantity related
81  procedure(special_opacity), pointer :: usr_special_opacity => null()
82  procedure(special_aniso_opacity), pointer :: usr_special_aniso_opacity => null()
83  procedure(special_opacity_qdot), pointer :: usr_special_opacity_qdot => null()
84  procedure(special_fluxlimiter), pointer :: usr_special_fluxlimiter => null()
85  procedure(special_diffcoef), pointer :: usr_special_diffcoef => null()
86 
87  ! Called after the mesh has been adjuste
88  procedure(after_refine), pointer :: usr_after_refine => null()
89 
90  ! initialize vector potential on cell edges for magnetic field
91  procedure(init_vector_potential), pointer :: usr_init_vector_potential => null()
92 
93  ! allow user to change inductive electric field, especially for boundary driven applications
94  procedure(set_electric_field), pointer :: usr_set_electric_field => null()
95 
96  ! allow user to specify variables at physical boundaries
97  procedure(set_wlr), pointer :: usr_set_wlr => null()
98 
99  ! allow user to specify the expansion function for the surface of a cross sectional
100  ! area of a 1D prominence, along with the analytical derivative of that function and its
101  ! primitive shape evaluated in the boundaries \int_(x_i-dx_i/2)^(x_i+dx_i/2) A(s) ds
102  procedure(set_surface), pointer :: usr_set_surface => null()
103 
104  ! for tracing field. allow user to specify variables and field
105  procedure(set_field_w), pointer :: usr_set_field_w => null()
106  procedure(set_field), pointer :: usr_set_field => null()
107 
108  ! allow user to specify R factor in ideal gas law with partial ionization
109  procedure(rfactor), pointer :: usr_rfactor => null()
110 
111  abstract interface
112 
113  subroutine p_no_args()
114  end subroutine p_no_args
115 
116  !> Initialize one grid
117  subroutine init_one_grid(ixI^L,ixO^L,w,x)
119  integer, intent(in) :: ixI^L, ixO^L
120  double precision, intent(in) :: x(ixI^S,1:ndim)
121  double precision, intent(inout) :: w(ixI^S,1:nw)
122  end subroutine init_one_grid
123 
124  !> special boundary types, users must assign conservative
125  !> variables in boundaries
126  subroutine special_bc(qt,ixI^L,ixO^L,iB,w,x)
128  !> Shape of input arrays
129  integer, intent(in) :: ixI^L
130  !> Region where boundary values have to be set
131  integer, intent(in) :: ixO^L
132  !> Integer indicating direction of boundary
133  integer, intent(in) :: iB
134  double precision, intent(in) :: qt, x(ixI^S,1:ndim)
135  double precision, intent(inout) :: w(ixI^S,1:nw)
136  end subroutine special_bc
137 
138  !> Special boundary type for radiation hydrodynamics module, only used to
139  !> set the boundary conditions for the radiation energy.
140  subroutine special_mg_bc(iB)
142  integer, intent(in) :: iB
143  end subroutine special_mg_bc
144 
145  !> internal boundary, user defined
146  !> This subroutine can be used to artificially overwrite ALL conservative
147  !> variables in a user-selected region of the mesh, and thereby act as
148  !> an internal boundary region. It is called just before external (ghost cell)
149  !> boundary regions will be set by the BC selection. Here, you could e.g.
150  !> want to introduce an extra variable (nwextra, to be distinguished from nwaux)
151  !> which can be used to identify the internal boundary region location.
152  !> Its effect should always be local as it acts on the mesh.
153  subroutine internal_bc(level,qt,ixI^L,ixO^L,w,x)
155  integer, intent(in) :: ixI^L,ixO^L,level
156  double precision, intent(in) :: qt
157  double precision, intent(inout) :: w(ixI^S,1:nw)
158  double precision, intent(in) :: x(ixI^S,1:ndim)
159  end subroutine internal_bc
160 
161  !> this subroutine is ONLY to be used for computing auxiliary variables
162  !> which happen to be non-local (like div v), and are in no way used for
163  !> flux computations. As auxiliaries, they are also not advanced
164  subroutine process_grid(igrid,level,ixI^L,ixO^L,qt,w,x)
166  integer, intent(in) :: igrid,level,ixI^L,ixO^L
167  double precision, intent(in) :: qt,x(ixI^S,1:ndim)
168  double precision, intent(inout) :: w(ixI^S,1:nw)
169  end subroutine process_grid
170 
171  !> If defined, this routine is called before writing output, and it can
172  !> set/modify the variables in the w array.
173  subroutine sub_modify_io(ixI^L,ixO^L,qt,w,x)
175  integer, intent(in) :: ixI^L,ixO^L
176  double precision, intent(in) :: qt,x(ixI^S,1:ndim)
177  double precision, intent(inout) :: w(ixI^S,1:nw)
178  end subroutine sub_modify_io
179 
180  !> This subroutine is called at the beginning of each time step
181  !> by each processor. No communication is specified, so the user
182  !> has to implement MPI routines if information has to be shared
183  subroutine process_global(iit,qt)
185  integer, intent(in) :: iit
186  double precision, intent(in) :: qt
187  end subroutine process_global
188 
189  !> for processing after the advance (PIC-MHD, e.g.)
190  subroutine process_adv_grid(igrid,level,ixI^L,ixO^L,qt,w,x)
192  integer, intent(in) :: igrid,level,ixI^L,ixO^L
193  double precision, intent(in) :: qt,x(ixI^S,1:ndim)
194  double precision, intent(inout) :: w(ixI^S,1:nw)
195  end subroutine process_adv_grid
196 
197  !> for processing after the advance (PIC-MHD, e.g.)
198  subroutine process_adv_global(iit,qt)
200  integer, intent(in) :: iit
201  double precision, intent(in) :: qt
202  end subroutine process_adv_global
203 
204  !> this subroutine can be used in convert, to add auxiliary variables to the
205  !> converted output file, for further analysis using tecplot, paraview, ....
206  !> these auxiliary values need to be stored in the nw+1:nw+nwauxio slots
207  !
208  !> the array normconv can be filled in the (nw+1:nw+nwauxio) range with
209  !> corresponding normalization values (default value 1)
210  subroutine aux_output(ixI^L,ixO^L,w,x,normconv)
212  integer, intent(in) :: ixI^L,ixO^L
213  double precision, intent(in) :: x(ixI^S,1:ndim)
214  double precision :: w(ixI^S,nw+nwauxio)
215  double precision :: normconv(0:nw+nwauxio)
216  end subroutine aux_output
217 
218  !> Add names for the auxiliary variables
219  subroutine add_aux_names(varnames)
221  character(len=*) :: varnames
222  end subroutine add_aux_names
223 
224  !> Calculate w(iw)=w(iw)+qdt*SOURCE[wCT,qtC,x] within ixO for all indices
225  !> iw=iwmin...iwmax. wCT is at time qCT
226  subroutine source(qdt,ixI^L,ixO^L,iw^LIM,qtC,wCT,qt,w,x)
228  integer, intent(in) :: ixI^L, ixO^L, iw^LIM
229  double precision, intent(in) :: qdt, qtC, qt
230  double precision, intent(in) :: wCT(ixI^S,1:nw), x(ixI^S,1:ndim)
231  double precision, intent(inout) :: w(ixI^S,1:nw)
232  end subroutine source
233 
234  !> Limit "dt" further if necessary, e.g. due to the special source terms.
235  !> The getdt_courant (CFL condition) and the getdt subroutine in the AMRVACPHYS
236  !> module have already been called.
237  subroutine get_dt(w,ixI^L,ixO^L,dtnew,dx^D,x)
239  integer, intent(in) :: ixI^L, ixO^L
240  double precision, intent(in) :: dx^D, x(ixI^S,1:ndim)
241  double precision, intent(in) :: w(ixI^S,1:nw)
242  double precision, intent(inout) :: dtnew
243  end subroutine get_dt
244 
245  !> Calculate gravitational acceleration in each dimension
246  subroutine phys_gravity(ixI^L,ixO^L,wCT,x,gravity_field)
248  integer, intent(in) :: ixI^L, ixO^L
249  double precision, intent(in) :: x(ixI^S,1:ndim)
250  double precision, intent(in) :: wCT(ixI^S,1:nw)
251  double precision, intent(out) :: gravity_field(ixI^S,ndim)
252  end subroutine phys_gravity
253 
254  !> Calculate the 3d drag force of gas onto dust
255  subroutine phys_dust_get_3d_dragforce(ixI^L, ixO^L, w, x, fdrag, ptherm, vgas,dust_n_species)
257  integer, intent(in) :: ixI^L, ixO^L, dust_n_species
258  double precision, intent(in) :: x(ixI^S, 1:ndim)
259  double precision, intent(in) :: w(ixI^S, 1:nw)
260  double precision, intent(out) :: &
261  fdrag(ixI^S, 1:ndir, 1:dust_n_species)
262  double precision, intent(in) :: ptherm(ixI^S), vgas(ixI^S, ndir)
263  end subroutine phys_dust_get_3d_dragforce
264 
265  !> Calculate the time step associated with the usr drag force
266  subroutine phys_dust_get_dt(w, ixI^L, ixO^L, dtdust, dx^D, x, dust_n_species)
268  integer, intent(in) :: ixI^L, ixO^L, dust_n_species
269  double precision, intent(in) :: dx^D, x(ixI^S,1:ndim)
270  double precision, intent(in) :: w(ixI^S,1:nw)
271  double precision, intent(inout) :: dtdust(1:dust_n_species)
272  end subroutine phys_dust_get_dt
273 
274  !>Calculation anormal viscosity depending on space
275  subroutine phys_visco(ixI^L,ixO^L,x,w,mu)
277  integer, intent(in) :: ixI^L, ixO^L
278  double precision, intent(in) :: x(ixI^S,1:ndim)
279  double precision, intent(in) :: w(ixI^S,1:nw)
280  double precision, intent(out) :: mu(ixI^S)
281  end subroutine phys_visco
282 
283  !>Calculation anormal pressure for hd & energy=.False.
284  subroutine hd_pthermal(w,x,ixI^L,ixO^L,pth)
286  integer, intent(in) :: ixI^L, ixO^L
287  double precision, intent(in) :: x(ixI^S,1:ndim)
288  double precision, intent(in) :: w(ixI^S,1:nw)
289  double precision, intent(out) :: pth(ixI^S)
290  end subroutine hd_pthermal
291 
292  !>Calculation R factor for ideal gas law with partial ionization
293  subroutine rfactor(w,x,ixI^L,ixO^L,pth)
295  integer, intent(in) :: ixI^L, ixO^L
296  double precision, intent(in) :: x(ixI^S,1:ndim)
297  double precision, intent(in) :: w(ixI^S,1:nw)
298  double precision, intent(out) :: pth(ixI^S)
299  end subroutine rfactor
300 
301  !> Set the "eta" array for resistive MHD based on w or the
302  !> "current" variable which has components between idirmin and 3.
303  subroutine special_resistivity(w,ixI^L,ixO^L,idirmin,x,current,eta)
305  integer, intent(in) :: ixI^L, ixO^L, idirmin
306  double precision, intent(in) :: w(ixI^S,nw), x(ixI^S,1:ndim)
307  double precision :: current(ixI^S,7-2*ndir:3), eta(ixI^S)
308  end subroutine special_resistivity
309 
310 
311  !> Set user defined opacity for use in diffusion coeff, heating and cooling, and radiation force
312  subroutine special_opacity(ixI^L,ixO^L,w,x,kappa)
314  integer, intent(in) :: ixI^L, ixO^L
315  double precision, intent(in) :: w(ixI^S,1:nw), x(ixI^S,1:ndim)
316  double precision, intent(out):: kappa(ixO^S)
317  end subroutine special_opacity
318 
319  !> Set user defined, anisotropic opacity for use in diffusion coeff, heating and cooling, and radiation force
320  subroutine special_aniso_opacity(ixI^L,ixO^L,w,x,kappa,idir)
322  integer, intent(in) :: ixI^L, ixO^L, idir
323  double precision, intent(in) :: w(ixI^S,1:nw), x(ixI^S,1:ndim)
324  double precision, intent(out):: kappa(ixO^S)
325  end subroutine special_aniso_opacity
326 
327  !> Set user defined opacity for use in diffusion coeff, heating and cooling, and radiation force. Overwrites special_opacity
328  subroutine special_opacity_qdot(ixI^L,ixO^L,w,x,kappa)
330  integer, intent(in) :: ixI^L, ixO^L
331  double precision, intent(in) :: w(ixI^S,1:nw), x(ixI^S,1:ndim)
332  double precision, intent(out):: kappa(ixO^S)
333  end subroutine special_opacity_qdot
334 
335  !> Set user defined FLD flux limiter, lambda
336  subroutine special_fluxlimiter(ixI^L,ixO^L,w,x,fld_lambda,fld_R)
338  integer, intent(in) :: ixI^L, ixO^L
339  double precision, intent(in) :: w(ixI^S,1:nw), x(ixI^S,1:ndim)
340  double precision, intent(out):: fld_lambda(ixI^S),fld_R(ixI^S)
341  end subroutine special_fluxlimiter
342 
343  !> Set user defined FLD diffusion coefficient
344  subroutine special_diffcoef(w, wCT, x, ixI^L, ixO^L)
346  integer, intent(in) :: ixI^L, ixO^L
347  double precision, intent(inout) :: w(ixI^S, 1:nw)
348  double precision, intent(in) :: wCT(ixI^S, 1:nw)
349  double precision, intent(in) :: x(ixI^S, 1:ndim)
350  end subroutine special_diffcoef
351 
352  !> Enforce additional refinement or coarsening
353  !> One can use the coordinate info in x and/or time qt=t_n and w(t_n) values w.
354  !> you must set consistent values for integers refine/coarsen:
355  !> refine = -1 enforce to not refine
356  !> refine = 0 doesn't enforce anything
357  !> refine = 1 enforce refinement
358  !> coarsen = -1 enforce to not coarsen
359  !> coarsen = 0 doesn't enforce anything
360  !> coarsen = 1 enforce coarsen
361  !> e.g. refine for negative first coordinate x < 0 as
362  !> if (any(x(ix^S,1) < zero)) refine=1
363  subroutine refine_grid(igrid,level,ixI^L,ixO^L,qt,w,x,refine,coarsen)
365  integer, intent(in) :: igrid, level, ixI^L, ixO^L
366  double precision, intent(in) :: qt, w(ixI^S,1:nw), x(ixI^S,1:ndim)
367  integer, intent(inout) :: refine, coarsen
368  end subroutine refine_grid
369 
370  !> this is the place to compute a local auxiliary variable to be used
371  !> as refinement criterion for the Lohner error estimator only
372  !> -->it is then requiring and iflag>nw
373  !> note that ixO=ixI=ixG, hence the term local (gradients need special attention!)
374  subroutine var_for_errest(ixI^L,ixO^L,iflag,w,x,var)
376  integer, intent(in) :: ixI^L,ixO^L,iflag
377  double precision, intent(in) :: w(ixI^S,1:nw), x(ixI^S,1:ndim)
378  double precision, intent(out) :: var(ixI^S)
379  end subroutine var_for_errest
380 
381  !> Here one can add a steady (time-independent) potential background field
382  subroutine set_b0(ixI^L,ixO^L,x,wB0)
384  integer, intent(in) :: ixI^L,ixO^L
385  double precision, intent(in) :: x(ixI^S,1:ndim)
386  double precision, intent(inout) :: wB0(ixI^S,1:ndir)
387  end subroutine set_b0
388 
389  !> Here one can add a time-independent background current density
390  subroutine set_j0(ixI^L,ixO^L,x,wJ0)
392  integer, intent(in) :: ixI^L,ixO^L
393  double precision, intent(in) :: x(ixI^S,1:ndim)
394  double precision, intent(inout) :: wJ0(ixI^S,7-2*ndir:ndir)
395  end subroutine set_j0
396 
397  !> Here one can add a steady (time-independent) equi vars
398  subroutine set_equi_vars(ixI^L,ixO^L,x,w0)
400  integer, intent(in) :: ixI^L,ixO^L
401  double precision, intent(in) :: x(ixI^S,1:ndim)
402  double precision, intent(inout) :: w0(ixI^S,1:number_equi_vars)
403  end subroutine set_equi_vars
404 
405  !> adjust w when restart from dat file with different w variables
406  subroutine transform_w(ixI^L,ixO^L,nw_in,w_in,x,w_out)
408  integer, intent(in) :: ixI^L, ixO^L, nw_in
409  double precision, intent(in) :: w_in(ixI^S,1:nw_in)
410  double precision, intent(in) :: x(ixI^S, 1:ndim)
411  double precision, intent(out) :: w_out(ixI^S,1:nw)
412  end subroutine transform_w
413 
414  !> use different threshold in special regions for AMR to
415  !> reduce/increase resolution there where nothing/something interesting happens.
416  subroutine a_refine_threshold(wlocal,xlocal,threshold,qt,level)
418  double precision, intent(in) :: wlocal(1:nw),xlocal(1:ndim),qt
419  double precision, intent(inout) :: threshold
420  integer, intent(in) :: level
421  end subroutine a_refine_threshold
422 
423  !> Allow user to use their own data-postprocessing procedures
424  subroutine special_convert(qunitconvert)
426  integer, intent(in) :: qunitconvert
427  character(len=20) :: userconvert_type
428  end subroutine special_convert
429 
430  !> flag=-1 : Treat all cells active, omit deactivation (onentry, default)
431  !> flag=0 : Treat as normal domain
432  !> flag=1 : Treat as passive, but reduce by safety belt
433  !> flag=2 : Always treat as passive
434  subroutine flag_grid(qt,ixI^L,ixO^L,w,x,flag)
436  integer, intent(in) :: ixI^L, ixO^L
437  integer, intent(inout) :: flag
438  double precision, intent(in) :: qt
439  double precision, intent(inout) :: w(ixI^S,1:nw)
440  double precision, intent(in) :: x(ixI^S,1:ndim)
441  end subroutine flag_grid
442 
443  !> Update payload of particles
444  subroutine update_payload(igrid,x,u,q,m,mypayload,mynpayload,particle_time)
446  integer, intent(in) :: igrid,mynpayload
447  double precision, intent(in) :: x(1:ndir),u(1:ndir),q,m,particle_time
448  double precision, intent(out) :: mypayload(mynpayload)
449  end subroutine update_payload
450 
451  !> Create particles
452  subroutine create_particles(n_particles, x, v, q, m, follow)
453  integer, intent(in) :: n_particles
454  double precision, intent(out) :: x(3, n_particles)
455  double precision, intent(out) :: v(3, n_particles)
456  double precision, intent(out) :: q(n_particles)
457  double precision, intent(out) :: m(n_particles)
458  logical, intent(out) :: follow(n_particles)
459  end subroutine create_particles
460 
461  !> Check arbitrary particle conditions or modifications
462  subroutine check_particle(igrid,x,v,q,m,follow,check)
464  integer, intent(in) :: igrid
465  double precision, intent(inout) :: x(1:ndir)
466  double precision, intent(inout) :: v(1:ndir),q,m
467  logical, intent(inout) :: follow
468  logical, intent(out) :: check
469  end subroutine check_particle
470 
471  !> Associate fields to particle
472  subroutine particle_fields(w, x, E, B)
474  double precision, intent(in) :: w(ixG^T,1:nw)
475  double precision, intent(in) :: x(ixG^T,1:ndim)
476  double precision, intent(out) :: E(ixG^T, ndir)
477  double precision, intent(out) :: B(ixG^T, ndir)
478  end subroutine particle_fields
479 
480  subroutine particle_analytic(ix, x, tloc, vec)
482  integer, intent(in) :: ix(ndir) !< Indices in gridvars
483  double precision, intent(in) :: x(ndir)
484  double precision, intent(in) :: tloc
485  double precision, intent(out) :: vec(ndir)
486  end subroutine particle_analytic
487 
488  !> User-defined particle movement
489  subroutine particle_position(x, n, tloc, tlocnew)
491  integer, intent(in) :: n
492  double precision, intent(inout) :: x(3)
493  double precision, intent(in) :: tloc, tlocnew
494  end subroutine particle_position
495 
496  subroutine after_refine(n_coarsen, n_refine)
497  integer, intent(in) :: n_coarsen
498  integer, intent(in) :: n_refine
499  end subroutine after_refine
500 
501  !> initialize vector potential on cell edges for magnetic field
502  subroutine init_vector_potential(ixI^L, ixC^L, xC, A, idir)
504 
505  integer, intent(in) :: ixI^L, ixC^L, idir
506  double precision, intent(in) :: xC(ixI^S,1:ndim)
507  double precision, intent(out) :: A(ixI^S)
508 
509  end subroutine init_vector_potential
510 
511  ! allow user to change inductive electric field, especially for boundary driven applications
512  subroutine set_electric_field(ixI^L,ixO^L,qt,qdt,fE,s)
514  integer, intent(in) :: ixI^L, ixO^L
515  double precision, intent(in) :: qt, qdt
516  type(state) :: s
517  double precision, intent(inout) :: fE(ixI^S,sdim:3)
518 
519  !integer :: ixC^L,ixA^L
520  ! For example, to set inductive electric field at bottom boundary in a 3D box for induction equation
521  ! v and b are from observational data for data-driven application
522 
523  !associate(w=>s%w,ws=>s%ws)
524 
525  !if(s%is_physical_boundary(5)) then
526  ! ixCmin^D=ixOmin^D-1;
527  ! ixCmax^D=ixOmax^D;
528  ! ixAmin^D=ixCmin^D;
529  ! ixAmax^D=ixCmax^D+1;
530  ! fE(nghostcells^%3ixA^S,1)=-ws(nghostcells^%3ixA^S,3)*w(nghostcells^%3ixA^S,mom(2))
531  ! fE(nghostcells^%3ixA^S,2)= ws(nghostcells^%3ixA^S,3)*w(nghostcells^%3ixA^S,mom(1))
532  ! ixAmin^D=ixCmin^D+kr(2,^D);
533  ! ixAmax^D=ixCmax^D+kr(2,^D);
534  ! fE(nghostcells^%3ixC^S,1)=0.5d0*(fE(nghostcells^%3ixC^S,1)+fE(nghostcells^%3ixA^S,1))*&
535  ! qdt*s%dsC(nghostcells^%3ixC^S,1)
536  ! ixAmin^D=ixCmin^D+kr(1,^D);
537  ! ixAmax^D=ixCmax^D+kr(1,^D);
538  ! fE(nghostcells^%3ixC^S,2)=0.5d0*(fE(nghostcells^%3ixC^S,2)+fE(nghostcells^%3ixA^S,2))*&
539  ! qdt*s%dsC(nghostcells^%3ixC^S,2)
540  !end if
541 
542  !end associate
543 
544  end subroutine set_electric_field
545 
546  !> allow user to specify 'variables' left and right state at physical boundaries to control flux through the boundary surface
547  subroutine set_wlr(ixI^L,ixO^L,qt,wLC,wRC,wLp,wRp,s,idir)
549  integer, intent(in) :: ixI^L, ixO^L, idir
550  double precision, intent(in) :: qt
551  double precision, intent(inout) :: wLC(ixI^S,1:nw), wRC(ixI^S,1:nw)
552  double precision, intent(inout) :: wLp(ixI^S,1:nw), wRp(ixI^S,1:nw)
553  type(state) :: s
554 
555  !if(s%is_physical_boundary(3).and.idir==2) then
556  ! wLp(ixOmin2^%2ixO^S,mom(1))=1.d0
557  ! wRp(ixOmin2^%2ixO^S,mom(1))=wRp(ixOmin2^%2ixO^S,mom(1))
558  ! wLC(ixOmin2^%2ixO^S,mom(1))=wLp(ixOmin2^%2ixO^S,mom(1))*wLp(ixOmin2^%2ixO^S,rho_)
559  ! wRC(ixOmin2^%2ixO^S,mom(1))=wRp(ixOmin2^%2ixO^S,mom(1))*wRp(ixOmin2^%2ixO^S,rho_)
560  !end if
561  end subroutine set_wlr
562 
563  subroutine set_surface(ixI^L,x,delx,exp_factor,del_exp_factor,exp_factor_primitive)
565  integer, intent(in) :: ixI^L
566  double precision, intent(in) :: delx(ixI^S,1:ndim), x(ixI^S,1:ndim)
567  double precision, intent(out) :: exp_factor(ixI^S), del_exp_factor(ixI^S)
568  double precision, intent(out) :: exp_factor_primitive(ixI^S)
569 
570  end subroutine set_surface
571 
572  subroutine set_field_w(igrid,ip,xf,wP,wL,numP,nwP,nwL,dL,forward,ftype,tcondi)
574  !use mod_point_searching
575 
576  integer, intent(in) :: igrid,ip,numP,nwP,nwL
577  double precision, intent(in) :: xf(numP,ndim)
578  double precision, intent(inout) :: wP(numP,nwP),wL(1+nwL)
579  double precision, intent(in) :: dL
580  logical, intent(in) :: forward
581  character(len=std_len), intent(in) :: ftype,tcondi
582 
583  !double precision :: xpp(1:ndim),wpp(1:nw)
584 
585  !! nwP=2,nwL=0. get rho/T at line
586  !if (tcondi=='user') then
587  ! xpp(1:ndim)=xf(ip,1:ndim)
588  ! call get_point_w_ingrid(igrid,xpp,wpp,'primitive')
589  ! wP(ip,1)=wpp(rho_)
590  ! wP(ip,2)=wpp(p_)/wpp(rho_)
591  !endif
592 
593  end subroutine set_field_w
594 
595  subroutine set_field(xfn,igrid,field,ftype)
597 
598  integer,intent(in) :: igrid
599  double precision, intent(in) :: xfn(ndim)
600  double precision, intent(inout) :: field(ndim)
601  character(len=std_len), intent(in) :: ftype
602 
603  !if (ftype='xdir') then
604  ! field(:)=zero
605  ! field(1)=1.d0
606  !endif
607 
608  end subroutine set_field
609 
610  end interface
611 
612 end module mod_usr_methods
use different threshold in special regions for AMR to reduce/increase resolution there where nothing/...
Add names for the auxiliary variables.
this subroutine can be used in convert, to add auxiliary variables to the converted output file,...
Check arbitrary particle conditions or modifications.
flag=-1 : Treat all cells active, omit deactivation (onentry, default) flag=0 : Treat as normal domai...
Limit "dt" further if necessary, e.g. due to the special source terms. The getdt_courant (CFL conditi...
Calculation anormal pressure for hd & energy=.False.
initialize vector potential on cell edges for magnetic field
internal boundary, user defined This subroutine can be used to artificially overwrite ALL conservativ...
Associate fields to particle.
User-defined particle movement.
Calculate the 3d drag force of gas onto dust.
Calculate the time step associated with the usr drag force.
Calculate gravitational acceleration in each dimension.
Calculation anormal viscosity depending on space.
for processing after the advance (PIC-MHD, e.g.)
for processing after the advance (PIC-MHD, e.g.)
This subroutine is called at the beginning of each time step by each processor. No communication is s...
this subroutine is ONLY to be used for computing auxiliary variables which happen to be non-local (li...
Enforce additional refinement or coarsening One can use the coordinate info in x and/or time qt=t_n a...
Here one can add a steady (time-independent) equi vars.
Calculate w(iw)=w(iw)+qdt*SOURCE[wCT,qtC,x] within ixO for all indices iw=iwmin......
Set user defined, anisotropic opacity for use in diffusion coeff, heating and cooling,...
special boundary types, users must assign conservative variables in boundaries
Allow user to use their own data-postprocessing procedures.
Set user defined FLD diffusion coefficient.
Set user defined FLD flux limiter, lambda.
Special boundary type for radiation hydrodynamics module, only used to set the boundary conditions fo...
Set user defined opacity for use in diffusion coeff, heating and cooling, and radiation force....
Set user defined opacity for use in diffusion coeff, heating and cooling, and radiation force.
Set the "eta" array for resistive MHD based on w or the "current" variable which has components betwe...
If defined, this routine is called before writing output, and it can set/modify the variables in the ...
adjust w when restart from dat file with different w variables
Update payload of particles.
this is the place to compute a local auxiliary variable to be used as refinement criterion for the Lo...
This module contains definitions of global parameters and variables and some generic functions/subrou...
Module with all the methods that users can customize in AMRVAC.
procedure(rfactor), pointer usr_rfactor
procedure(source), pointer usr_source
procedure(special_resistivity), pointer usr_special_resistivity
procedure(special_opacity), pointer usr_special_opacity
procedure(process_grid), pointer usr_process_grid
procedure(phys_visco), pointer usr_setvisco
procedure(particle_position), pointer usr_particle_position
procedure(check_particle), pointer usr_check_particle
procedure(p_no_args), pointer usr_improve_initial_condition
procedure(a_refine_threshold), pointer usr_refine_threshold
procedure(set_surface), pointer usr_set_surface
procedure(phys_dust_get_dt), pointer usr_dust_get_dt
procedure(phys_gravity), pointer usr_gravity
procedure(aux_output), pointer usr_aux_output
procedure(p_no_args), pointer usr_print_log
procedure(phys_dust_get_3d_dragforce), pointer usr_get_3d_dragforce
procedure(special_diffcoef), pointer usr_special_diffcoef
procedure(process_adv_grid), pointer usr_process_adv_grid
procedure(particle_analytic), pointer usr_particle_analytic
procedure(special_convert), pointer usr_special_convert
procedure(create_particles), pointer usr_create_particles
procedure(init_one_grid), pointer usr_init_one_grid
Initialize earch grid block data.
procedure(update_payload), pointer usr_update_payload
procedure(p_no_args), pointer usr_write_analysis
procedure(sub_modify_io), pointer usr_modify_output
procedure(p_no_args), pointer usr_before_main_loop
procedure(special_fluxlimiter), pointer usr_special_fluxlimiter
procedure(set_field_w), pointer usr_set_field_w
procedure(flag_grid), pointer usr_flag_grid
procedure(process_global), pointer usr_process_global
procedure(special_bc), pointer usr_special_bc
procedure(process_adv_global), pointer usr_process_adv_global
procedure(internal_bc), pointer usr_internal_bc
procedure(set_equi_vars), pointer usr_set_equi_vars
procedure(special_mg_bc), pointer usr_special_mg_bc
procedure(set_j0), pointer usr_set_j0
procedure(special_opacity_qdot), pointer usr_special_opacity_qdot
procedure(particle_fields), pointer usr_particle_fields
procedure(refine_grid), pointer usr_refine_grid
procedure(hd_pthermal), pointer usr_set_pthermal
procedure(init_vector_potential), pointer usr_init_vector_potential
procedure(p_no_args), pointer usr_set_parameters
Initialize the user's settings (after initializing amrvac)
procedure(var_for_errest), pointer usr_var_for_errest
procedure(set_b0), pointer usr_set_b0
procedure(set_electric_field), pointer usr_set_electric_field
procedure(special_aniso_opacity), pointer usr_special_aniso_opacity
procedure(set_wlr), pointer usr_set_wlr
procedure(transform_w), pointer usr_transform_w
procedure(after_refine), pointer usr_after_refine
procedure(get_dt), pointer usr_get_dt
procedure(set_field), pointer usr_set_field
procedure(add_aux_names), pointer usr_add_aux_names