The FARGO_THORIN code developer's guide
 All Data Structures Files Functions Variables Typedefs Macros
Macros | Functions | Variables
Pebbles.c File Reference

Contains functions reponsible for the pebble disk initialisation, evolution due to source terms and pebble accretion. More...

#include "fargo.h"
Include dependency graph for Pebbles.c:

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 PolarGridPebbleDensTemp
 
static PolarGridPebbleAccelrad
 
static PolarGridPebbleAcceltheta
 
static PolarGridPebbleSize
 
static PolarGridEtaFaceCentered
 
static PolarGridEtaCellCentered
 
static PolarGridAccretedMass
 
static real pebbulkdens
 

Detailed Description

Contains functions reponsible for the pebble disk initialisation, evolution due to source terms and pebble accretion.

Also controls several outputs.

Author
Ondřej Chrenko chren.nosp@m.ko@s.nosp@m.irrah.nosp@m..tro.nosp@m.ja.mf.nosp@m.f.cu.nosp@m.ni.cz

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.

LICENSE

Copyright (c) 2017 Ondřej Chrenko. See the LICENSE file of the distribution.

Definition in file Pebbles.c.

Macro Definition Documentation

#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().

Function Documentation

void AccretePebblesOntoPlanets ( PlanetarySystem sys,
PolarGrid Rho,
PolarGrid Energy,
PolarGrid Vtheta,
real  dt 
)
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().

Here is the caller graph for this function:

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().

Here is the call graph for this function:

Here is the caller graph for this function:

void CriticalCharTime ( PolarGrid Vrad,
PolarGrid Vtheta 
)

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().

Here is the call graph for this function:

Here is the caller graph for this function:

boolean DetectCrashPebbles ( )

Safety check for negative pebble densities.

Definition at line 753 of file Pebbles.c.

References DetectCrash(), and PebbleDens.

Referenced by AlgoGas().

Here is the call graph for this function:

Here is the caller graph for this function:

void EquilPebbleDisk ( PolarGrid Rho,
PolarGrid Vrad,
PolarGrid Vtheta 
)

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().

Here is the call graph for this function:

Here is the caller graph for this function:

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().

Here is the caller graph for this function:

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().

Here is the call graph for this function:

Here is the caller graph for this function:

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().

Here is the call graph for this function:

Here is the caller graph for this function:

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.

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().

Here is the call graph for this function:

Here is the caller graph for this function:

real IntegrateColumnMass ( real  a,
real  b,
real  Hpeb 
)

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().

Here is the call graph for this function:

Here is the caller graph for this function:

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().

Here is the caller graph for this function:

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().

Here is the call graph for this function:

Here is the caller graph for this function:

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().

Here is the caller graph for this function:

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().

Here is the call graph for this function:

Here is the caller graph for this function:

void SourceTermsPebbles ( PolarGrid Vrad,
PolarGrid Vtheta,
real  dt 
)

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().

Here is the call graph for this function:

Here is the caller graph for this function:

void SubStep1Pebbles ( PolarGrid Vrad,
PolarGrid Vtheta,
real  dt 
)

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().

Here is the caller graph for this function:

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().

Here is the call graph for this function:

Here is the caller graph for this function:

real Trapzd ( int  n,
real  a,
real  b,
real  H2 
)

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().

Here is the caller graph for this function:

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().

Here is the call graph for this function:

Here is the caller graph for this function:

Variable Documentation

PolarGrid* AccretedMass
static

Definition at line 29 of file Pebbles.c.

PolarGrid * EtaCellCentered
static

Definition at line 28 of file Pebbles.c.

PolarGrid* EtaFaceCentered
static

Definition at line 28 of file Pebbles.c.

boolean FastTransport

Definition at line 29 of file Interpret.c.

Referenced by CriticalCharTime().

PolarGrid * PebbleAccelrad
static

Definition at line 26 of file Pebbles.c.

PolarGrid * PebbleAcceltheta
static

Definition at line 26 of file Pebbles.c.

PolarGrid* PebbleDensTemp
static

Definition at line 26 of file Pebbles.c.

PolarGrid* PebbleSize
static

Definition at line 27 of file Pebbles.c.

real pebbulkdens
static

Definition at line 31 of file Pebbles.c.

Referenced by InitPebbleArrays(), InitPebblesViaFlux(), and PebbleStokesNumbers().

boolean Restart

Definition at line 14 of file main.c.

Referenced by main().