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