MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
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 ! allow user to specify adiabatic index and gamma dependent on space
112 procedure(set_adiab), pointer :: usr_set_adiab => null()
113 procedure(set_adiab), pointer :: usr_set_gamma => null()
114
115 abstract interface
116
117 subroutine p_no_args()
118 end subroutine p_no_args
119
120 !> Initialize one grid
121 subroutine init_one_grid(ixI^L,ixO^L,w,x)
123 integer, intent(in) :: ixI^L, ixO^L
124 double precision, intent(in) :: x(ixI^S,1:ndim)
125 double precision, intent(inout) :: w(ixI^S,1:nw)
126 end subroutine init_one_grid
127
128 !> special boundary types, users must assign conservative
129 !> variables in boundaries
130 subroutine special_bc(qt,ixI^L,ixO^L,iB,w,x)
132 !> Shape of input arrays
133 integer, intent(in) :: ixI^L
134 !> Region where boundary values have to be set
135 integer, intent(in) :: ixO^L
136 !> Integer indicating direction of boundary
137 integer, intent(in) :: iB
138 double precision, intent(in) :: qt, x(ixI^S,1:ndim)
139 double precision, intent(inout) :: w(ixI^S,1:nw)
140 end subroutine special_bc
141
142 !> Special boundary type for radiation hydrodynamics module, only used to
143 !> set the boundary conditions for the radiation energy.
144 subroutine special_mg_bc(iB)
146 integer, intent(in) :: iB
147 end subroutine special_mg_bc
148
149 !> internal boundary, user defined
150 !> This subroutine can be used to artificially overwrite ALL conservative
151 !> variables in a user-selected region of the mesh, and thereby act as
152 !> an internal boundary region. It is called just before external (ghost cell)
153 !> boundary regions will be set by the BC selection. Here, you could e.g.
154 !> want to introduce an extra variable (nwextra, to be distinguished from nwaux)
155 !> which can be used to identify the internal boundary region location.
156 !> Its effect should always be local as it acts on the mesh.
157 subroutine internal_bc(level,qt,ixI^L,ixO^L,w,x)
159 integer, intent(in) :: ixI^L,ixO^L,level
160 double precision, intent(in) :: qt
161 double precision, intent(inout) :: w(ixI^S,1:nw)
162 double precision, intent(in) :: x(ixI^S,1:ndim)
163 end subroutine internal_bc
164
165 !> this subroutine is ONLY to be used for computing auxiliary variables
166 !> which happen to be non-local (like div v), and are in no way used for
167 !> flux computations. As auxiliaries, they are also not advanced
168 subroutine process_grid(igrid,level,ixI^L,ixO^L,qt,w,x)
170 integer, intent(in) :: igrid,level,ixI^L,ixO^L
171 double precision, intent(in) :: qt,x(ixI^S,1:ndim)
172 double precision, intent(inout) :: w(ixI^S,1:nw)
173 end subroutine process_grid
174
175 !> If defined, this routine is called before writing output, and it can
176 !> set/modify the variables in the w array.
177 subroutine sub_modify_io(ixI^L,ixO^L,qt,w,x)
179 integer, intent(in) :: ixI^L,ixO^L
180 double precision, intent(in) :: qt,x(ixI^S,1:ndim)
181 double precision, intent(inout) :: w(ixI^S,1:nw)
182 end subroutine sub_modify_io
183
184 !> This subroutine is called at the beginning of each time step
185 !> by each processor. No communication is specified, so the user
186 !> has to implement MPI routines if information has to be shared
187 subroutine process_global(iit,qt)
189 integer, intent(in) :: iit
190 double precision, intent(in) :: qt
191 end subroutine process_global
192
193 !> for processing after the advance (PIC-MHD, e.g.)
194 subroutine process_adv_grid(igrid,level,ixI^L,ixO^L,qt,w,x)
196 integer, intent(in) :: igrid,level,ixI^L,ixO^L
197 double precision, intent(in) :: qt,x(ixI^S,1:ndim)
198 double precision, intent(inout) :: w(ixI^S,1:nw)
199 end subroutine process_adv_grid
200
201 !> for processing after the advance (PIC-MHD, e.g.)
202 subroutine process_adv_global(iit,qt)
204 integer, intent(in) :: iit
205 double precision, intent(in) :: qt
206 end subroutine process_adv_global
207
208 !> this subroutine can be used in convert, to add auxiliary variables to the
209 !> converted output file, for further analysis using tecplot, paraview, ....
210 !> these auxiliary values need to be stored in the nw+1:nw+nwauxio slots
211 !
212 !> the array normconv can be filled in the (nw+1:nw+nwauxio) range with
213 !> corresponding normalization values (default value 1)
214 subroutine aux_output(ixI^L,ixO^L,w,x,normconv)
216 integer, intent(in) :: ixI^L,ixO^L
217 double precision, intent(in) :: x(ixI^S,1:ndim)
218 double precision :: w(ixI^S,nw+nwauxio)
219 double precision :: normconv(0:nw+nwauxio)
220 end subroutine aux_output
221
222 !> Add names for the auxiliary variables
223 subroutine add_aux_names(varnames)
225 character(len=*) :: varnames
226 end subroutine add_aux_names
227
228 !> Calculate w(iw)=w(iw)+qdt*SOURCE[wCT,qtC,x] within ixO for all indices
229 !> iw=iwmin...iwmax. wCT is at time qCT
230 subroutine source(qdt,ixI^L,ixO^L,iw^LIM,qtC,wCT,qt,w,x)
232 integer, intent(in) :: ixI^L, ixO^L, iw^LIM
233 double precision, intent(in) :: qdt, qtC, qt
234 double precision, intent(in) :: wCT(ixI^S,1:nw), x(ixI^S,1:ndim)
235 double precision, intent(inout) :: w(ixI^S,1:nw)
236 end subroutine source
237
238 !> Limit "dt" further if necessary, e.g. due to the special source terms.
239 !> The getdt_courant (CFL condition) and the getdt subroutine in the AMRVACPHYS
240 !> module have already been called.
241 subroutine get_dt(w,ixI^L,ixO^L,dtnew,dx^D,x)
243 integer, intent(in) :: ixI^L, ixO^L
244 double precision, intent(in) :: dx^D, x(ixI^S,1:ndim)
245 double precision, intent(in) :: w(ixI^S,1:nw)
246 double precision, intent(inout) :: dtnew
247 end subroutine get_dt
248
249 !> Calculate gravitational acceleration in each dimension
250 subroutine phys_gravity(ixI^L,ixO^L,wCT,x,gravity_field)
252 integer, intent(in) :: ixI^L, ixO^L
253 double precision, intent(in) :: x(ixI^S,1:ndim)
254 double precision, intent(in) :: wCT(ixI^S,1:nw)
255 double precision, intent(out) :: gravity_field(ixI^S,ndim)
256 end subroutine phys_gravity
257
258 !> Calculate the 3d drag force of gas onto dust
259 subroutine phys_dust_get_3d_dragforce(ixI^L, ixO^L, w, x, fdrag, ptherm, vgas,dust_n_species)
261 integer, intent(in) :: ixI^L, ixO^L, dust_n_species
262 double precision, intent(in) :: x(ixI^S, 1:ndim)
263 double precision, intent(in) :: w(ixI^S, 1:nw)
264 double precision, intent(out) :: &
265 fdrag(ixI^S, 1:ndir, 1:dust_n_species)
266 double precision, intent(in) :: ptherm(ixI^S), vgas(ixI^S, ndir)
267 end subroutine phys_dust_get_3d_dragforce
268
269 !> Calculate the time step associated with the usr drag force
270 subroutine phys_dust_get_dt(w, ixI^L, ixO^L, dtdust, dx^D, x, dust_n_species)
272 integer, intent(in) :: ixI^L, ixO^L, dust_n_species
273 double precision, intent(in) :: dx^D, x(ixI^S,1:ndim)
274 double precision, intent(in) :: w(ixI^S,1:nw)
275 double precision, intent(inout) :: dtdust(1:dust_n_species)
276 end subroutine phys_dust_get_dt
277
278 !>Calculation anormal viscosity depending on space
279 subroutine phys_visco(ixI^L,ixO^L,x,w,mu)
281 integer, intent(in) :: ixI^L, ixO^L
282 double precision, intent(in) :: x(ixI^S,1:ndim)
283 double precision, intent(in) :: w(ixI^S,1:nw)
284 double precision, intent(out) :: mu(ixI^S)
285 end subroutine phys_visco
286
287 !>Calculation anormal pressure for hd & energy=.False.
288 subroutine hd_pthermal(w,x,ixI^L,ixO^L,pth)
290 integer, intent(in) :: ixI^L, ixO^L
291 double precision, intent(in) :: x(ixI^S,1:ndim)
292 double precision, intent(in) :: w(ixI^S,1:nw)
293 double precision, intent(out) :: pth(ixI^S)
294 end subroutine hd_pthermal
295
296 !>Calculation R factor for ideal gas law with partial ionization
297 subroutine rfactor(w,x,ixI^L,ixO^L,pth)
299 integer, intent(in) :: ixI^L, ixO^L
300 double precision, intent(in) :: x(ixI^S,1:ndim)
301 double precision, intent(in) :: w(ixI^S,1:nw)
302 double precision, intent(out) :: pth(ixI^S)
303 end subroutine rfactor
304
305 !>set adiabatic index
306 subroutine set_adiab(w,x,ixI^L,ixO^L,adiab)
308 integer, intent(in) :: ixI^L, ixO^L
309 double precision, intent(in) :: x(ixI^S,1:ndim)
310 double precision, intent(in) :: w(ixI^S,1:nw)
311 double precision, intent(out) :: adiab(ixO^S)
312 end subroutine set_adiab
313
314 !> Set the "eta" array for resistive MHD based on w or the
315 !> "current" variable which has components between idirmin and 3.
316 subroutine special_resistivity(w,ixI^L,ixO^L,idirmin,x,current,eta)
318 integer, intent(in) :: ixI^L, ixO^L, idirmin
319 double precision, intent(in) :: w(ixI^S,nw), x(ixI^S,1:ndim)
320 double precision :: current(ixI^S,7-2*ndir:3), eta(ixI^S)
321 end subroutine special_resistivity
322
323
324 !> Set user defined opacity for use in diffusion coeff, heating and cooling, and radiation force
325 subroutine special_opacity(ixI^L,ixO^L,w,x,kappa)
327 integer, intent(in) :: ixI^L, ixO^L
328 double precision, intent(in) :: w(ixI^S,1:nw), x(ixI^S,1:ndim)
329 double precision, intent(out):: kappa(ixO^S)
330 end subroutine special_opacity
331
332 !> Set user defined, anisotropic opacity for use in diffusion coeff, heating and cooling, and radiation force
333 subroutine special_aniso_opacity(ixI^L,ixO^L,w,x,kappa,idir)
335 integer, intent(in) :: ixI^L, ixO^L, idir
336 double precision, intent(in) :: w(ixI^S,1:nw), x(ixI^S,1:ndim)
337 double precision, intent(out):: kappa(ixO^S)
338 end subroutine special_aniso_opacity
339
340 !> Set user defined opacity for use in diffusion coeff, heating and cooling, and radiation force. Overwrites special_opacity
341 subroutine special_opacity_qdot(ixI^L,ixO^L,w,x,kappa)
343 integer, intent(in) :: ixI^L, ixO^L
344 double precision, intent(in) :: w(ixI^S,1:nw), x(ixI^S,1:ndim)
345 double precision, intent(out):: kappa(ixO^S)
346 end subroutine special_opacity_qdot
347
348 !> Set user defined FLD flux limiter, lambda
349 subroutine special_fluxlimiter(ixI^L,ixO^L,w,x,fld_lambda,fld_R)
351 integer, intent(in) :: ixI^L, ixO^L
352 double precision, intent(in) :: w(ixI^S,1:nw), x(ixI^S,1:ndim)
353 double precision, intent(out):: fld_lambda(ixI^S),fld_R(ixI^S)
354 end subroutine special_fluxlimiter
355
356 !> Set user defined FLD diffusion coefficient
357 subroutine special_diffcoef(w, wCT, x, ixI^L, ixO^L)
359 integer, intent(in) :: ixI^L, ixO^L
360 double precision, intent(inout) :: w(ixI^S, 1:nw)
361 double precision, intent(in) :: wCT(ixI^S, 1:nw)
362 double precision, intent(in) :: x(ixI^S, 1:ndim)
363 end subroutine special_diffcoef
364
365 !> Enforce additional refinement or coarsening
366 !> One can use the coordinate info in x and/or time qt=t_n and w(t_n) values w.
367 !> you must set consistent values for integers refine/coarsen:
368 !> refine = -1 enforce to not refine
369 !> refine = 0 doesn't enforce anything
370 !> refine = 1 enforce refinement
371 !> coarsen = -1 enforce to not coarsen
372 !> coarsen = 0 doesn't enforce anything
373 !> coarsen = 1 enforce coarsen
374 !> e.g. refine for negative first coordinate x < 0 as
375 !> if (any(x(ix^S,1) < zero)) refine=1
376 subroutine refine_grid(igrid,level,ixI^L,ixO^L,qt,w,x,refine,coarsen)
378 integer, intent(in) :: igrid, level, ixI^L, ixO^L
379 double precision, intent(in) :: qt, w(ixI^S,1:nw), x(ixI^S,1:ndim)
380 integer, intent(inout) :: refine, coarsen
381 end subroutine refine_grid
382
383 !> this is the place to compute a local auxiliary variable to be used
384 !> as refinement criterion for the Lohner error estimator only
385 !> -->it is then requiring and iflag>nw
386 !> note that ixO=ixI=ixG, hence the term local (gradients need special attention!)
387 subroutine var_for_errest(ixI^L,ixO^L,iflag,w,x,var)
389 integer, intent(in) :: ixI^L,ixO^L,iflag
390 double precision, intent(in) :: w(ixI^S,1:nw), x(ixI^S,1:ndim)
391 double precision, intent(out) :: var(ixI^S)
392 end subroutine var_for_errest
393
394 !> Here one can add a steady (time-independent) potential background field
395 subroutine set_b0(ixI^L,ixO^L,x,wB0)
397 integer, intent(in) :: ixI^L,ixO^L
398 double precision, intent(in) :: x(ixI^S,1:ndim)
399 double precision, intent(inout) :: wB0(ixI^S,1:ndir)
400 end subroutine set_b0
401
402 !> Here one can add a time-independent background current density
403 subroutine set_j0(ixI^L,ixO^L,x,wJ0)
405 integer, intent(in) :: ixI^L,ixO^L
406 double precision, intent(in) :: x(ixI^S,1:ndim)
407 double precision, intent(inout) :: wJ0(ixI^S,7-2*ndir:ndir)
408 end subroutine set_j0
409
410 !> Here one can add a steady (time-independent) equi vars
411 subroutine set_equi_vars(ixI^L,ixO^L,x,w0)
413 integer, intent(in) :: ixI^L,ixO^L
414 double precision, intent(in) :: x(ixI^S,1:ndim)
415 double precision, intent(inout) :: w0(ixI^S,1:number_equi_vars)
416 end subroutine set_equi_vars
417
418 !> adjust w when restart from dat file with different w variables
419 subroutine transform_w(ixI^L,ixO^L,nw_in,w_in,x,w_out)
421 integer, intent(in) :: ixI^L, ixO^L, nw_in
422 double precision, intent(in) :: w_in(ixI^S,1:nw_in)
423 double precision, intent(in) :: x(ixI^S, 1:ndim)
424 double precision, intent(out) :: w_out(ixI^S,1:nw)
425 end subroutine transform_w
426
427 !> use different threshold in special regions for AMR to
428 !> reduce/increase resolution there where nothing/something interesting happens.
429 subroutine a_refine_threshold(wlocal,xlocal,threshold,qt,level)
431 double precision, intent(in) :: wlocal(1:nw),xlocal(1:ndim),qt
432 double precision, intent(inout) :: threshold
433 integer, intent(in) :: level
434 end subroutine a_refine_threshold
435
436 !> Allow user to use their own data-postprocessing procedures
437 subroutine special_convert(qunitconvert)
439 integer, intent(in) :: qunitconvert
440 character(len=20) :: userconvert_type
441 end subroutine special_convert
442
443 !> flag=-1 : Treat all cells active, omit deactivation (onentry, default)
444 !> flag=0 : Treat as normal domain
445 !> flag=1 : Treat as passive, but reduce by safety belt
446 !> flag=2 : Always treat as passive
447 subroutine flag_grid(qt,ixI^L,ixO^L,w,x,flag)
449 integer, intent(in) :: ixI^L, ixO^L
450 integer, intent(inout) :: flag
451 double precision, intent(in) :: qt
452 double precision, intent(inout) :: w(ixI^S,1:nw)
453 double precision, intent(in) :: x(ixI^S,1:ndim)
454 end subroutine flag_grid
455
456 !> Update payload of particles
457 subroutine update_payload(igrid,x,u,q,m,mypayload,mynpayload,particle_time)
459 integer, intent(in) :: igrid,mynpayload
460 double precision, intent(in) :: x(1:ndir),u(1:ndir),q,m,particle_time
461 double precision, intent(out) :: mypayload(mynpayload)
462 end subroutine update_payload
463
464 !> Create particles
465 subroutine create_particles(n_particles, x, v, q, m, follow)
466 integer, intent(in) :: n_particles
467 double precision, intent(out) :: x(3, n_particles)
468 double precision, intent(out) :: v(3, n_particles)
469 double precision, intent(out) :: q(n_particles)
470 double precision, intent(out) :: m(n_particles)
471 logical, intent(out) :: follow(n_particles)
472 end subroutine create_particles
473
474 !> Check arbitrary particle conditions or modifications
475 subroutine check_particle(igrid,x,v,q,m,follow,check)
477 integer, intent(in) :: igrid
478 double precision, intent(inout) :: x(1:ndir)
479 double precision, intent(inout) :: v(1:ndir),q,m
480 logical, intent(inout) :: follow
481 logical, intent(out) :: check
482 end subroutine check_particle
483
484 !> Associate fields to particle
485 subroutine particle_fields(w, x, E, B)
487 double precision, intent(in) :: w(ixG^T,1:nw)
488 double precision, intent(in) :: x(ixG^T,1:ndim)
489 double precision, intent(out) :: E(ixG^T, ndir)
490 double precision, intent(out) :: B(ixG^T, ndir)
491 end subroutine particle_fields
492
493 subroutine particle_analytic(ix, x, tloc, vec)
495 integer, intent(in) :: ix(ndir) !< Indices in gridvars
496 double precision, intent(in) :: x(ndir)
497 double precision, intent(in) :: tloc
498 double precision, intent(out) :: vec(ndir)
499 end subroutine particle_analytic
500
501 !> User-defined particle movement
502 subroutine particle_position(x, n, tloc, tlocnew)
504 integer, intent(in) :: n
505 double precision, intent(inout) :: x(3)
506 double precision, intent(in) :: tloc, tlocnew
507 end subroutine particle_position
508
509 subroutine after_refine(n_coarsen, n_refine)
510 integer, intent(in) :: n_coarsen
511 integer, intent(in) :: n_refine
512 end subroutine after_refine
513
514 !> initialize vector potential on cell edges for magnetic field
515 subroutine init_vector_potential(ixI^L, ixC^L, xC, A, idir)
517
518 integer, intent(in) :: ixI^L, ixC^L, idir
519 double precision, intent(in) :: xC(ixI^S,1:ndim)
520 double precision, intent(out) :: A(ixI^S)
521
522 end subroutine init_vector_potential
523
524 ! allow user to change inductive electric field, especially for boundary driven applications
525 subroutine set_electric_field(ixI^L,ixO^L,qt,qdt,fE,s)
527 integer, intent(in) :: ixI^L, ixO^L
528 double precision, intent(in) :: qt, qdt
529 type(state) :: s
530 double precision, intent(inout) :: fE(ixI^S,sdim:3)
531
532 !integer :: ixC^L,ixA^L
533 ! For example, to set inductive electric field at bottom boundary in a 3D box for induction equation
534 ! v and b are from observational data for data-driven application
535
536 !associate(w=>s%w,ws=>s%ws)
537
538 !if(s%is_physical_boundary(5)) then
539 ! ixCmin^D=ixOmin^D-1;
540 ! ixCmax^D=ixOmax^D;
541 ! ixAmin^D=ixCmin^D;
542 ! ixAmax^D=ixCmax^D+1;
543 ! fE(nghostcells^%3ixA^S,1)=-ws(nghostcells^%3ixA^S,3)*w(nghostcells^%3ixA^S,mom(2))
544 ! fE(nghostcells^%3ixA^S,2)= ws(nghostcells^%3ixA^S,3)*w(nghostcells^%3ixA^S,mom(1))
545 ! ixAmin^D=ixCmin^D+kr(2,^D);
546 ! ixAmax^D=ixCmax^D+kr(2,^D);
547 ! fE(nghostcells^%3ixC^S,1)=0.5d0*(fE(nghostcells^%3ixC^S,1)+fE(nghostcells^%3ixA^S,1))*&
548 ! qdt*s%dsC(nghostcells^%3ixC^S,1)
549 ! ixAmin^D=ixCmin^D+kr(1,^D);
550 ! ixAmax^D=ixCmax^D+kr(1,^D);
551 ! fE(nghostcells^%3ixC^S,2)=0.5d0*(fE(nghostcells^%3ixC^S,2)+fE(nghostcells^%3ixA^S,2))*&
552 ! qdt*s%dsC(nghostcells^%3ixC^S,2)
553 !end if
554
555 !end associate
556
557 end subroutine set_electric_field
558
559 !> allow user to specify 'variables' left and right state at physical boundaries to control flux through the boundary surface
560 subroutine set_wlr(ixI^L,ixO^L,qt,wLC,wRC,wLp,wRp,s,idir)
562 integer, intent(in) :: ixI^L, ixO^L, idir
563 double precision, intent(in) :: qt
564 double precision, intent(inout) :: wLC(ixI^S,1:nw), wRC(ixI^S,1:nw)
565 double precision, intent(inout) :: wLp(ixI^S,1:nw), wRp(ixI^S,1:nw)
566 type(state) :: s
567
568 !if(s%is_physical_boundary(3).and.idir==2) then
569 ! wLp(ixOmin2^%2ixO^S,mom(1))=1.d0
570 ! wRp(ixOmin2^%2ixO^S,mom(1))=wRp(ixOmin2^%2ixO^S,mom(1))
571 ! wLC(ixOmin2^%2ixO^S,mom(1))=wLp(ixOmin2^%2ixO^S,mom(1))*wLp(ixOmin2^%2ixO^S,rho_)
572 ! wRC(ixOmin2^%2ixO^S,mom(1))=wRp(ixOmin2^%2ixO^S,mom(1))*wRp(ixOmin2^%2ixO^S,rho_)
573 !end if
574 end subroutine set_wlr
575
576 subroutine set_surface(ixI^L,x,delx,exp_factor,del_exp_factor,exp_factor_primitive)
578 integer, intent(in) :: ixI^L
579 double precision, intent(in) :: delx(ixI^S,1:ndim), x(ixI^S,1:ndim)
580 double precision, intent(out) :: exp_factor(ixI^S), del_exp_factor(ixI^S)
581 double precision, intent(out) :: exp_factor_primitive(ixI^S)
582
583 end subroutine set_surface
584
585 subroutine set_field_w(igrid,ip,xf,wP,wL,numP,nwP,nwL,dL,forward,ftype,tcondi)
587 !use mod_point_searching
588
589 integer, intent(in) :: igrid,ip,numP,nwP,nwL
590 double precision, intent(in) :: xf(numP,ndim)
591 double precision, intent(inout) :: wP(numP,nwP),wL(1+nwL)
592 double precision, intent(in) :: dL
593 logical, intent(in) :: forward
594 character(len=std_len), intent(in) :: ftype,tcondi
595
596 !double precision :: xpp(1:ndim),wpp(1:nw)
597
598 !! nwP=2,nwL=0. get rho/T at line
599 !if (tcondi=='user') then
600 ! xpp(1:ndim)=xf(ip,1:ndim)
601 ! call get_point_w_ingrid(igrid,xpp,wpp,'primitive')
602 ! wP(ip,1)=wpp(rho_)
603 ! wP(ip,2)=wpp(p_)/wpp(rho_)
604 !endif
605
606 end subroutine set_field_w
607
608 subroutine set_field(xfn,igrid,field,ftype)
610
611 integer,intent(in) :: igrid
612 double precision, intent(in) :: xfn(ndim)
613 double precision, intent(inout) :: field(ndim)
614 character(len=std_len), intent(in) :: ftype
615
616 !if (ftype='xdir') then
617 ! field(:)=zero
618 ! field(1)=1.d0
619 !endif
620
621 end subroutine set_field
622
623 end interface
624
625end 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...
Calculation R factor for ideal gas law with partial ionization.
Here one can add a steady (time-independent) potential background field.
Here one can add a steady (time-independent) equi vars.
Here one can add a time-independent background current density.
allow user to specify 'variables' left and right state at physical boundaries to control flux through...
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(set_adiab), pointer usr_set_adiab
procedure(special_opacity), pointer usr_special_opacity
procedure(set_adiab), pointer usr_set_gamma
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