The FARGO_THORIN code developer's guide
|
Contains functions reponsible for the pebble disk initialisation, evolution due to source terms and pebble accretion. More...
#include "fargo.h"
Go to the source code of this file.
Macros | |
#define | TRAPEZMAX 35 |
#define | TRAPEZEPS 1.0e-7 |
Functions | |
void | InitPebbleArrays () |
Initialise polar arrays associated with the pebble disk. More... | |
void | EquilPebbleDisk (PolarGrid *Rho, PolarGrid *Vrad, PolarGrid *Vtheta) |
Sets up a pebble disk in a coagulation-drift equilibrium. More... | |
void | RestartPebbleDisk (PolarGrid *Rho, int index) |
Reads arrays necessary to restart the pebble disk. More... | |
void | InitPebblesViaFlux (PolarGrid *Rho, PolarGrid *Vrad) |
Imposing the radial mass flux of pebbles in a coagulation-drift equilibrium, initialises their surface density, flow velocity and local Stokes numbers. More... | |
void | BckpFieldsForBC () |
Backs up the initial state of the pebble disk to impose damping boundary conditions later. More... | |
void | EtaPressureSupport (PolarGrid *Vtheta) |
Calculates tha gas rotation parameter eta. More... | |
void | PebbleStokesNumbers (PolarGrid *Rho) |
Calculates the local Stokes number using the dominant pebble size, local gas density and parametric pebble material density. More... | |
real | Trapzd (int n, real a, real b, real H2) |
Trapezoidal rule to integrate the column mass of pebbles located within the accretion radius. More... | |
real | IntegrateColumnMass (real a, real b, real Hpeb) |
Primitive algorithm to find the mass in a vertical column of pebbles overlapping the accretion radius. More... | |
void | AccretePebblesOntoPlanets (PlanetarySystem *sys, PolarGrid *Rho, PolarGrid *Energy, PolarGrid *Vtheta, real dt) |
Finds the amount of pebbles to be transfered from the pebble disk onto the planets. More... | |
void | CorrectPebblesVtheta (real domega) |
Corrects the azimuthal flow velocity of pebbles to keep up with the frame rotation. More... | |
void | SourceTermsPebbles (PolarGrid *Vrad, PolarGrid *Vtheta, real dt) |
Calculates the source terms acting on pebbles. More... | |
void | SubStep1Pebbles (PolarGrid *Vrad, PolarGrid *Vtheta, real dt) |
Applies a semi-implicit method to evolve the pebbles dynamically. More... | |
void | ParticleDiffusion (PolarGrid *Rho) |
Applies the particle diffusion term acting on pebbles. More... | |
void | EvolvePebbleDisk (real dt) |
Calls the transport routines and applies the boundary conditions for pebbles. More... | |
void | SynchronizePebbleDisc () |
Synchronises pebble fluid hydrodynamic quantities among the overlapping grid zones. More... | |
void | CriticalCharTime (PolarGrid *Vrad, PolarGrid *Vtheta) |
Restricts the time step using the CFL condition for the pebble fluid. More... | |
void | WritePebbles (int index) |
Outputs the pebble fluid arrays. More... | |
boolean | DetectCrashPebbles () |
Safety check for negative pebble densities. More... | |
void | ParametricAccretion (PlanetarySystem *sys, real dt) |
Writes filtering factors if needed. More... | |
Variables | |
boolean | Restart |
boolean | FastTransport |
static PolarGrid * | PebbleDensTemp |
static PolarGrid * | PebbleAccelrad |
static PolarGrid * | PebbleAcceltheta |
static PolarGrid * | PebbleSize |
static PolarGrid * | EtaFaceCentered |
static PolarGrid * | EtaCellCentered |
static PolarGrid * | AccretedMass |
static real | pebbulkdens |
Contains functions reponsible for the pebble disk initialisation, evolution due to source terms and pebble accretion.
Also controls several outputs.
The pebble disk equilibrium model is inspired by Lambrechts & Johansen (2014). The evolution follows a standard two-fluid approximation with linear drag coupling and particle diffusion.
Copyright (c) 2017 Ondřej Chrenko. See the LICENSE file of the distribution.
Definition in file Pebbles.c.
#define TRAPEZEPS 1.0e-7 |
Definition at line 22 of file Pebbles.c.
Referenced by IntegrateColumnMass().
#define TRAPEZMAX 35 |
Definition at line 21 of file Pebbles.c.
Referenced by IntegrateColumnMass().
void AccretePebblesOntoPlanets | ( | PlanetarySystem * | sys, |
PolarGrid * | Rho, | ||
PolarGrid * | Energy, | ||
PolarGrid * | Vtheta, | ||
real | dt | ||
) |
Finds the amount of pebbles to be transfered from the pebble disk onto the planets.
Definition at line 312 of file Pebbles.c.
References AccretHeating, CellAbscissa, CellOrdinate, DT, polargrid::Field, GetGlobalIFrac(), HEATINGDELAY, heatsrc, heatsrc_index, heatsrc_max, IntegrateColumnMass(), InvRmed, Max_or_active, MPI_Allreduce(), MPI_COMM_WORLD, MPI_DOUBLE, MPI_INT, mpi_make1Dprofile(), MPI_SUM, polargrid::Nrad, polargrid::Nsec, OmegaFrame, OmegaInv, PEBBLEALPHA, PebbleDens, PebbleVrad, PebbleVtheta, PhysicalTime, PI, PLANETARYDENSITY, RHO2CGS, Rinf, Rmed, Rsup, SoundSpeed, SQRT2PI_INV, SQRT_ADIABIND_INV, StokesNumber, Surf, vt1D, VXplanet, VYplanet, Xplanet, Yplanet, and Zero_or_active.
Referenced by AlgoGas().
void BckpFieldsForBC | ( | ) |
Backs up the initial state of the pebble disk to impose damping boundary conditions later.
Definition at line 168 of file Pebbles.c.
References polargrid::Field, polargrid::Nrad, polargrid::Nsec, OmegaFrame, PebbleDens, PebbleVrad, PebbleVtheta, PebDensInit, PebVradInit, PebVthetaInit, and Rmed.
Referenced by RestartPebbleDisk().
void CorrectPebblesVtheta | ( | real | domega | ) |
Corrects the azimuthal flow velocity of pebbles to keep up with the frame rotation.
Definition at line 516 of file Pebbles.c.
References CorrectVtheta(), and PebbleVtheta.
Referenced by AlgoGas().
Restricts the time step using the CFL condition for the pebble fluid.
Definition at line 702 of file Pebbles.c.
References FastTransport, polargrid::Field, invdtpeb_sq, max2(), Max_or_active, polargrid::Nsec, One_or_active, PebbleVrad, PebbleVtheta, PI, Rinf, Rmed, and Rsup.
Referenced by AlgoGas().
boolean DetectCrashPebbles | ( | ) |
Safety check for negative pebble densities.
Definition at line 753 of file Pebbles.c.
References DetectCrash(), and PebbleDens.
Referenced by AlgoGas().
Sets up a pebble disk in a coagulation-drift equilibrium.
Definition at line 58 of file Pebbles.c.
References EtaPressureSupport(), and InitPebblesViaFlux().
Referenced by Initialization().
void EtaPressureSupport | ( | PolarGrid * | Vtheta | ) |
Calculates tha gas rotation parameter eta.
Definition at line 195 of file Pebbles.c.
References polargrid::Field, polargrid::Nrad, OmegaFrame, and Rmed.
Referenced by EquilPebbleDisk().
void EvolvePebbleDisk | ( | real | dt | ) |
Calls the transport routines and applies the boundary conditions for pebbles.
Definition at line 675 of file Pebbles.c.
References DampPebbles(), PebbleDens, PebbleVrad, PebbleVtheta, and TransportPebbles().
Referenced by AlgoGas().
void InitPebbleArrays | ( | ) |
Initialise polar arrays associated with the pebble disk.
Definition at line 37 of file Pebbles.c.
References CreatePolarGrid(), DragForceRad, DragForceTheta, GasAccelrad, GasAcceltheta, NRAD, NSEC, PEBBLEBULKDENS, PebbleDens, PebbleVrad, PebbleVtheta, pebbulkdens, RHO2CGS, and StokesNumber.
Referenced by InitEuler().
Imposing the radial mass flux of pebbles in a coagulation-drift equilibrium, initialises their surface density, flow velocity and local Stokes numbers.
Definition at line 80 of file Pebbles.c.
References ADIABIND, CPU_Master, CPU_Number, polargrid::Field, FLUX2CU, Merge, OmegaFrame, OUTPUTDIR, PEBBLECOAGULATION, PebbleDens, PEBBLEFLUX, PebbleVrad, PebbleVtheta, pebbulkdens, PebDensInit, PebVradInit, PebVthetaInit, PI, Rmed, StokesNumber, and WriteDiskPolar().
Referenced by EquilPebbleDisk().
Primitive algorithm to find the mass in a vertical column of pebbles overlapping the accretion radius.
Definition at line 292 of file Pebbles.c.
References masterprint(), prs_exit(), TRAPEZEPS, TRAPEZMAX, and Trapzd().
Referenced by AccretePebblesOntoPlanets().
void ParametricAccretion | ( | PlanetarySystem * | sys, |
real | dt | ||
) |
Writes filtering factors if needed.
Mass accretion onto planets is provided by a parametric prescription using a given mass doubling time.
Definition at line 779 of file Pebbles.c.
References AccretHeating, DT, HEATINGDELAY, heatsrc, heatsrc_index, heatsrc_max, polargrid::Nrad, polargrid::Nsec, PARAMETRICACCRETION, PhysicalTime, PI, PLANETARYDENSITY, RHO2CGS, Rinf, Rsup, SoundSpeed, Surf, Xplanet, and Yplanet.
Referenced by AlgoGas().
void ParticleDiffusion | ( | PolarGrid * | Rho | ) |
Applies the particle diffusion term acting on pebbles.
See Eq. (35) in Chrenko et al. (2017)
Definition at line 629 of file Pebbles.c.
References polargrid::Field, FViscosity(), InvDiffRmed, PebbleDens, PebbleVrad, PebbleVtheta, PI, Rinf, Rmed, and SCHMIDTNUMBER.
Referenced by AlgoGas().
void PebbleStokesNumbers | ( | PolarGrid * | Rho | ) |
Calculates the local Stokes number using the dominant pebble size, local gas density and parametric pebble material density.
Definition at line 240 of file Pebbles.c.
References polargrid::Field, polargrid::Nrad, pebbulkdens, SQRT_ADIABIND_INV, and StokesNumber.
Referenced by AlgoGas().
void RestartPebbleDisk | ( | PolarGrid * | Rho, |
int | index | ||
) |
Reads arrays necessary to restart the pebble disk.
Definition at line 66 of file Pebbles.c.
References BckpFieldsForBC(), PebbleDens, PebbleVrad, PebbleVtheta, and ReadfromFile().
Referenced by Initialization().
Calculates the source terms acting on pebbles.
Definition at line 523 of file Pebbles.c.
References BackReaction, DampPebbles(), DragForceRad, DragForceTheta, polargrid::Field, InvRinf, OmegaFrame, OmegaInv, PebbleDens, PebbleGravAccelRad, PebbleGravAccelTheta, PebbleVrad, PebbleVtheta, Rinf, and StokesNumber.
Referenced by AlgoGas().
Applies a semi-implicit method to evolve the pebbles dynamically.
See Appendix C in Chrenko et al. (2017)
Definition at line 590 of file Pebbles.c.
References polargrid::Field, GasAccelrad, GasAcceltheta, OmegaInv, PebbleVrad, PebbleVtheta, and StokesNumber.
Referenced by AlgoGas().
void SynchronizePebbleDisc | ( | ) |
Synchronises pebble fluid hydrodynamic quantities among the overlapping grid zones.
Definition at line 685 of file Pebbles.c.
References CPU_Number, CPUOVERLAP, polargrid::Field, polargrid::Nrad, PebbleDens, PebbleVrad, PebbleVtheta, and SynchronizeOverlapFields().
Referenced by AlgoGas().
Trapezoidal rule to integrate the column mass of pebbles located within the accretion radius.
Adopted from Numerical Recipes.
Definition at line 265 of file Pebbles.c.
References a.
Referenced by IntegrateColumnMass().
void WritePebbles | ( | int | index | ) |
Outputs the pebble fluid arrays.
Definition at line 743 of file Pebbles.c.
References PebbleDens, PebbleVrad, PebbleVtheta, Write_Eta, and WriteDiskPolar().
Referenced by SendOutput().
boolean FastTransport |
Definition at line 29 of file Interpret.c.
Referenced by CriticalCharTime().
|
static |
Definition at line 31 of file Pebbles.c.
Referenced by InitPebbleArrays(), InitPebblesViaFlux(), and PebbleStokesNumbers().