dfnTrans
Code for Particle Tracking simulations in 3D DFN
|
Go to the source code of this file.
Data Structures | |
struct | material |
struct | vertex |
struct | element |
struct | contam |
Macros | |
#define | pi 3.14159265359 |
Functions | |
void | ReadInit () |
void | ReadDataFiles () |
void | AdjacentCells (int ln, int i, int j, int k) |
void | Convertto2d () |
void | Convertto3d () |
void | DarcyVelocity () |
struct matr | MatrixProducts (double normxarea[][2], int number) |
struct lb | DefineBoundaryAngle (int i, unsigned int edge_1, unsigned int edge_2, int f1, int coorf) |
void | VelocityInteriorNode (double normxarea[][2], int i, int number, unsigned int indj[max_neighb], int vi) |
void | VelocityExteriorNode (double normxarea[][2], int i, int number, unsigned int indj[max_neighb], struct lb lbound, int vi) |
void | CheckNewCell () |
void | ParticleTrack () |
struct intcoef | CalculateWeights (int nn1, int nn2, int nn3) |
void | SearchNeighborCells (int nn1, int nn2, int nn3) |
int | InsideCell (unsigned int numc) |
void | NeighborCells (int k) |
void | PredictorStep () |
void | CorrectorStep () |
void | DefineTimeStep () |
int | CheckDistance () |
void | AcrossIntersection (int prevcell, int int1, int int2, int mixing_rule) |
void | ChangeFracture (int cell_win) |
struct posit3d | CalculatePosition3D () |
int | InitCell () |
int | InitPos () |
void | Moving2Center (int nnp, int cellnumber) |
int | Moving2NextCell (int stuck, int k) |
void | BoundaryCells () |
int | InitParticles_np (int k_current, int firstn, int lastn, int parts_fracture, int first_ind, int last_ind) |
int | InitParticles_eq (int k_current, int firstn, int lastn, double parts_dist, int first_ind, int last_ind) |
int | CornerVelocity (int i, int m1, int m2, int m3, int s1, int s2, int s3) |
void | ReadBoundaryNodes () |
FILE * | OpenFile (char filen[120], char fileopt[2]) |
void | ReadFEHMfile (int nedges) |
void | ReadPFLOTRANfile (int nedges) |
void | WritingInit () |
void | Velocity3D () |
double | CalculateCurrentDT () |
int | Yindex (int nodenum, int np) |
int | Xindex (int nodenum, int np) |
int | CompleteMixingRandomSampling (double products[4], double speedsq[4], int indj, int int1, int indk) |
int | StreamlineRandomSampling (double products[4], double speedsq[4], int indj, int int1, int indk, int neighborcellind[4], int neighborfracind[4], int prevfrac, int prevcell) |
void | OutputVelocities () |
int | XindexC (int nodenum, int ii) |
int | YindexC (int nodenum, int ii) |
double | DefineAngle (double u1, double u2, double v1, double v2) |
void | HalfPolygonVelocity (int i, int k, int fractn, int indc, unsigned int fractj[max_neighb]) |
struct posit3d | CalculateVelocity3D () |
void | BoundaryLine (int n1, int n2, int n3) |
void | CheckGrid () |
int | BVelocityDirection (int b1, int b2) |
double | InOutFlowCell (int indcell, int int1, double nposx, double nposy) |
void | Moving2NextCellBound (int prevcell) |
struct inpfile | Control_File (char fileobject[], int ctr) |
struct inpfile | Control_Data (char fileobject[], int ctr) |
void | ParticleOutput (int currentt, int frac_p) |
struct inpfile | Control_Param (char fileobject[], int ctr) |
void | FlowInWeight (int numbpart) |
int | InitParticles_ones (int k_current, double inter_p[][4], int fracture, int parts_fracture, int ii, double thirdcoor, int zonenumb_in, int first_ind, int last_ind) |
void | Coordinations2D () |
void | ReadAperture () |
void | InitInMatrix () |
double | TimeFromMatrix (double pdist) |
void | FinalPosition () |
struct lagrangian | CalculateLagrangian (double xcurrent, double ycurrent, double zcurrent, double xprev, double yprev, double zprev) |
void | OutputMarPlumDisp (int currentnum, char path[125]) |
int | String_Compare (char string1[], char string2[]) |
struct inpfile | Control_File_Optional (char fileobject[], int ctr) |
double | TimeDomainRW (double time_advect) |
int | InitParticles_flux (int k_current, int firstn, int lastn, double weight_p) |
int | InitInWell (int nodepart) |
Variables | |
char | maindir [125] |
double | porosity |
double | density |
unsigned long int | timesteps |
double | thickness |
double | saturation |
double | timeunit |
char | controlfile [120] |
double | totalFluxIn |
unsigned int | pflotran |
unsigned int | fehm |
unsigned int | nnodes |
unsigned int | ncells |
unsigned int | nfract |
unsigned int | max_neighb |
unsigned int | npart |
unsigned int | np |
unsigned int | nzone_in |
unsigned int * | nodezonein |
unsigned int * | nodezoneout |
unsigned int | flag_w |
struct material * | fracture |
struct vertex * | node |
struct contam * | particle |
struct element * | cell |
void AcrossIntersection | ( | int | prevcell, |
int | int1, | ||
int | int2, | ||
int | mixing_rule | ||
) |
Particle is moving through intersection line
Definition at line 1784 of file TrackingPart.c.
References CalculateWeights(), cell, vertex::cells, ChangeFracture(), CompleteMixingRandomSampling(), vertex::coord_xy, vertex::fracture, element::fracture, contam::intcell, mixing_rule, node, element::node_ind, np, particle, contam::position, contam::prev_pos, StreamlineRandomSampling(), element::veloc_ind, vertex::velocity, intcoef::weights, XindexC(), and YindexC().
Referenced by CheckDistance().
void AdjacentCells | ( | int | ln, |
int | in, | ||
int | jn, | ||
int | kn | ||
) |
The function defines adjacent triangular cells for each node
Definition at line 485 of file ReadGridInit.c.
References cell, vertex::cells, vertex::fracts, vertex::fracture, element::fracture, fracture, node, and vertex::numneighb.
Referenced by ReadDataFiles().
void BoundaryCells | ( | ) |
The function checks velocities at boundary nodes and fixes pathological cases.
In the rare event, when two velocity vectors are directed towards each other along the same cell edge, a particle will be stuck there. To avoid this situation, one of velocities is redirected along the edge with node associated with lower pressure value.
Definition at line 447 of file VelocityReconstruction.c.
References cell, vertex::cells, vertex::coord_xy, CornerVelocity(), DefineAngle(), fracture, vertex::indnodes, ncells, node, element::node_ind, vertex::numneighb, pi, vertex::pressure, element::veloc_ind, vertex::velocity, XindexC(), and YindexC().
Referenced by DarcyVelocity().
void BoundaryLine | ( | int | n1, |
int | n2, | ||
int | n3 | ||
) |
int BVelocityDirection | ( | int | b1, |
int | b2 | ||
) |
double CalculateCurrentDT | ( | ) |
Functions returns particle instanteneous time step
Definition at line 2239 of file TrackingPart.c.
References contam::cell, cell, node, element::node_ind, np, particle, vertex::timestep, element::veloc_ind, and contam::weight.
Referenced by CalculateWeights(), CheckDistance(), CorrectorStep(), and PredictorStep().
struct lagrangian CalculateLagrangian | ( | double | xcurrent, |
double | ycurrent, | ||
double | zcurrent, | ||
double | xprev, | ||
double | yprev, | ||
double | zprev | ||
) |
Function calculates Lagrangian variables: tau and beta.
Definition at line 2634 of file TrackingPart.c.
References contam::cell, cell, vertex::cells, vertex::coord_xy, vertex::indnodes, node, element::node_ind, np, vertex::numneighb, particle, contam::position, contam::prev_pos, contam::time, Xindex(), and Yindex().
Referenced by ParticleTrack().
struct posit3d CalculatePosition3D | ( | ) |
Function calculates 3D coordinates of current particle's position at 2D fracture plane
Definition at line 138 of file RotateFracture.c.
References CalculatePosition3D(), contam::cell, cell, posit3d::cord3, element::fracture, contam::fracture, fracture, np, particle, contam::position, and material::rot2mat.
Referenced by CalculatePosition3D(), ParticleOutput(), and ParticleTrack().
struct posit3d CalculateVelocity3D | ( | ) |
The function converts particle's 2D velocity vector to 3D velocity vector
Definition at line 190 of file RotateFracture.c.
References posit3d::cord3, vertex::fracture, fracture, maindir, nnodes, node, OpenFile(), material::rot3mat, and vertex::velocity.
Referenced by ParticleOutput(), and ParticleTrack().
struct intcoef CalculateWeights | ( | int | nn1, |
int | nn2, | ||
int | nn3 | ||
) |
Function calculates interpolation weights for velocity instanteneous particles velocity and time step.
Definition at line 1111 of file TrackingPart.c.
References CalculateCurrentDT(), CalculateWeights(), contam::cell, cell, FLAG_OUT, contam::fracture, node, element::node_ind, np, vertex::numneighb, particle, contam::position, contam::prev_pos, SearchNeighborCells(), contam::velocity, contam::weight, and intcoef::weights.
Referenced by AcrossIntersection(), CalculateWeights(), InOutFlowCell(), and InsideCell().
void ChangeFracture | ( | int | cell_win | ) |
This function recalculates particles coordinations at intersection lines. Particles XY coordinations at one fracture are recalculated to 3D positions and then new XY coordinations of an intersecting fracture are defined.
Definition at line 138 of file RotateFracture.c.
Referenced by AcrossIntersection(), and InOutFlowCell().
int CheckDistance | ( | ) |
Function checks the distance from particle to intersection line when particles is located at the intersection triangular cell
Definition at line 1474 of file TrackingPart.c.
References AcrossIntersection(), CalculateCurrentDT(), contam::cell, cell, vertex::cells, CheckNewCell(), vertex::coord_xy, vertex::fracture, fracture, vertex::indnodes, InOutFlowCell(), mixing_rule, no_out, node, element::node_ind, np, vertex::numneighb, particle, ParticleOutput(), contam::position, PredictorStep(), contam::prev_pos, t, t_adv, t_adv0, contam::t_adv_diff, contam::t_diff, tdrw, tdrw_o, contam::time, timediff, TimeDomainRW(), Xindex(), and Yindex().
Referenced by ParticleTrack().
void CheckGrid | ( | ) |
The function checks the grid: looking for nodes that are defined as internal or external, but belong to two fractures and should be defined as interface nodes
Definition at line 965 of file ReadGridInit.c.
References fracture, nnodes, node, and vertex::numneighb.
Referenced by main().
void CheckNewCell | ( | ) |
Function performs a check, did particle move to a new triangular cell during last time step or stayed at the same cell
Definition at line 1111 of file TrackingPart.c.
Referenced by CheckDistance(), and ParticleTrack().
int CompleteMixingRandomSampling | ( | double | products[4], |
double | speedsq[4], | ||
int | indj, | ||
int | int1, | ||
int | indk | ||
) |
Particle motion at a fracture intersection is determined by the streamline routing rule Arg 1: Cross product to define outgoing and incoming cells. Vector contains value for each cell at the intersection Arg 2: Vector of the speed squared for each cell at the intersection Arg 3: index of current cell in int1 node list Arg 4: int1 is the index of the node of the incoming cell that lies on the intersection Arg 5: indk gives the position in the list of the 4 intersection cells that is previous cell Return: The Exit cell at the intersection
Definition at line 2089 of file TrackingPart.c.
References vertex::cells, and node.
Referenced by AcrossIntersection().
struct inpfile Control_Data | ( | char | fileobject[], |
int | ctr | ||
) |
The function reads control file with input parameters to dfnTrans; returns input parameter value. If the parameter is not defined in the control file, the program is terminated.
Definition at line 1089 of file ReadGridInit.c.
Referenced by InitPos(), main(), and ReadBoundaryNodes().
struct inpfile Control_File | ( | char | fileobject[], |
int | ctr | ||
) |
The function reads control file with input parameters to dfnTrans; returns the file name and/or input parameter value. If the parameter is not defined in the control file, the program is terminated.
Definition at line 965 of file ReadGridInit.c.
Referenced by InitInMatrix(), InitPos(), main(), OutputMarPlumDisp(), ParticleTrack(), ReadAperture(), ReadBoundaryNodes(), ReadDataFiles(), ReadFEHMfile(), ReadInit(), and ReadPFLOTRANfile().
struct inpfile Control_File_Optional | ( | char | fileobject[], |
int | ctr | ||
) |
The function reads control file with input parameters to dfnTrans; returns the file name and/or input parameter value. This function is called for the optional parameters only. If the parameter is not defined in the control file, the default value is used.
Definition at line 965 of file ReadGridInit.c.
Referenced by InitPos(), ParticleTrack(), and ReadInit().
struct inpfile Control_Param | ( | char | fileobject[], |
int | ctr | ||
) |
The function reads control file with input parameters to dfnTrans; returns input parameter value. If the parameter is not defined in the control file, the program is terminated.
Definition at line 1089 of file ReadGridInit.c.
Referenced by InitPos(), main(), ParticleTrack(), and TimeFromMatrix().
void Convertto2d | ( | ) |
The function uses rotation matrix and fracture's normal vector to rotate fracture from its position in 3D to XY plane
if node belongs to intersection, belongs to two fractures
Definition at line 17 of file RotateFracture.c.
References vertex::coord, vertex::coord_xy, vertex::fracture, fracture, nfract, nnodes, node, material::nvect_xy, pi, material::rot2mat, material::rot3mat, and material::theta.
Referenced by main().
void Convertto3d | ( | ) |
The function uses rotation matrix to rotate fracture from its position in XY plane to 3D domain
Definition at line 100 of file RotateFracture.c.
References fracture, nfract, material::nvect_z, pi, material::rot3mat, and material::theta.
Referenced by main().
void Coordinations2D | ( | ) |
The function outputs 2D coordinations of nodes: every fracture into separate file. Those files are used as input to gstat for length correlation of aperture.
Definition at line 280 of file RotateFracture.c.
References vertex::coord_xy, fracture, nfract, nnodes, node, and OpenFile().
int CornerVelocity | ( | int | i, |
int | m1, | ||
int | m2, | ||
int | m3, | ||
int | s1, | ||
int | s2, | ||
int | s3 | ||
) |
The Function identifies pathologycal cases at corner of fracture boundary. The angle between velocity on corner node and velocities/edges on surrounding boundary nodes are calculated. Then, if velocity is pointing outside of fracture, turns the velocity direction along the fracture edge.
Definition at line 658 of file VelocityReconstruction.c.
References vertex::coord_xy, DefineAngle(), node, pi, vertex::velocity, XindexC(), and YindexC().
Referenced by BoundaryCells().
void CorrectorStep | ( | ) |
Corrector step in Predictor-Corrector technique. Function calculates new particle position using calculated velocity in Predictor step.
Definition at line 1408 of file TrackingPart.c.
References CalculateCurrentDT(), contam::cell, cell, element::node_ind, np, particle, contam::position, contam::prev_pos, element::veloc_ind, and contam::velocity.
Referenced by ParticleTrack().
void DarcyVelocity | ( | ) |
Function performs Darcy's velocity reconstruction (main function of Darcy velocities reconstruction procedure) using linear least square algorithm. Flow solver provides flow fluxes on each edge of control volume cells. Velocity on each control volume cell center (each node) is reconstructed using flow fluxes.
Velocity of interior and interior-interface nodes is reconstructed according to Eq.5 (Painter,2011); exterior node' velocity according to Eq.7 (Painter, 2011)
Function makes a loop over all nodes and calls functions for velocity reconstruction depending of type of the node (external, internal, internal-interface, external interface)
Definition at line 29 of file VelocityReconstruction.c.
References BoundaryCells(), vertex::coord_xy, DefineBoundaryAngle(), fehm, vertex::fracture, fracture, HalfPolygonVelocity(), vertex::indnodes, max_neighb, nnodes, node, vertex::numneighb, pflotran, vertex::velocity, VelocityExteriorNode(), and VelocityInteriorNode().
Referenced by main().
double DefineAngle | ( | double | u1, |
double | u2, | ||
double | v1, | ||
double | v2 | ||
) |
Function defines an angle between two edges of voronoy cell
Definition at line 347 of file VelocityReconstruction.c.
References lb::angle.
Referenced by BoundaryCells(), CornerVelocity(), and ParticleOutput().
struct lb DefineBoundaryAngle | ( | int | i, |
unsigned int | edge_1, | ||
unsigned int | edge_2, | ||
int | f1, | ||
int | coorf | ||
) |
The function defines angle between two edges of a boundary cell and return the result. The angle is used in flow velocity reconstruction on fracture edge.
Definition at line 214 of file VelocityReconstruction.c.
References fracture, maindir, nfract, nnodes, node, and vertex::velocity.
Referenced by DarcyVelocity(), and HalfPolygonVelocity().
void DefineTimeStep | ( | ) |
Function performs a loop through all the nodes in the mesh and defines a time step at each node. The time step of particles will be interpolated from time steps defined at each node.
Definition at line 1447 of file TrackingPart.c.
References nnodes, node, vertex::timestep, and vertex::velocity.
Referenced by main().
void FinalPosition | ( | ) |
Function calculates particles final position at out-flow boundary
Definition at line 2634 of file TrackingPart.c.
Referenced by ParticleTrack().
void FlowInWeight | ( | int | numberpart | ) |
Function calculates particle's in-flow flux weight. Used in option #1, #2, #3.
Definition at line 747 of file InitialPartPositions.c.
References contam::cell, cell, contam::fl_weight, InitCell(), node, element::node_ind, np, vertex::numneighb, particle, and contam::weight.
Referenced by ParticleTrack().
void HalfPolygonVelocity | ( | int | i, |
int | k, | ||
int | fractn, | ||
int | indc, | ||
unsigned int | fractj[max_neighb] | ||
) |
Velocity reconstruction on an intersection node. In this case, the control volume cell is devided onto two polygons by intersection line. The velocity reconstruction is performed on each part of control volume cell at intersection.
Definition at line 806 of file VelocityReconstruction.c.
References vertex::area, cell, vertex::cells, vertex::coord_xy, DefineBoundaryAngle(), fehm, fracture, vertex::indnodes, max_neighb, node, pflotran, vertex::type, vertex::typeN, element::veloc_ind, VelocityExteriorNode(), and VelocityInteriorNode().
Referenced by DarcyVelocity().
int InitCell | ( | ) |
Function performs a search to find cell Id where the particle was initially placed.
Definition at line 542 of file InitialPartPositions.c.
References vertex::cells, InsideCell(), node, nodezonein, vertex::numneighb, and nzone_in.
Referenced by FlowInWeight(), and ParticleTrack().
void InitInMatrix | ( | ) |
Function performs data reading from files for Option #5, init_matrix. In this option particles are placed at random positions in matrix, the time to move from initial position to nearest fracture is calculated in function TimeFromMatrix.
Input files ParticleInitCoordR.dat and ClosestNodeR.inp provide particles initial positions at matrix and te closes node ID at DFN mesh, respectively.
Definition at line 832 of file InitialPartPositions.c.
References contam::cell, Control_File(), vertex::coord, vertex::coord_xy, inpfile::filename, contam::fl_weight, vertex::fracture, contam::fracture, contam::intcell, node, npart, OpenFile(), particle, contam::position, contam::pressure, contam::prev_pos, contam::time, TimeFromMatrix(), contam::velocity, and contam::weight.
Referenced by InitPos().
int InitInWell | ( | int | nodepart | ) |
Definition at line 1028 of file InitialPartPositions.c.
References vertex::aperture, contam::cell, vertex::coord_xy, contam::fl_weight, vertex::fracture, contam::fracture, node, nodezonein, nzone_in, particle, contam::position, contam::pressure, contam::time, and contam::velocity.
Referenced by InitPos().
int InitParticles_eq | ( | int | k_current, |
int | firstn, | ||
int | lastn, | ||
double | parts_dist, | ||
int | first_ind, | ||
int | last_ind | ||
) |
Function defines particles initial positions at a single fracture edge, using calculated before distance between particles. Option #2, init_eqd.
Definition at line 623 of file InitialPartPositions.c.
References vertex::coord_xy, contam::fl_weight, vertex::fracture, contam::fracture, contam::intcell, node, npart, particle, contam::position, contam::pressure, contam::time, and contam::velocity.
Referenced by InitPos().
int InitParticles_flux | ( | int | k_current, |
int | first_ind, | ||
int | last_ind, | ||
double | weight_p | ||
) |
Function defines particle's initial positions at single fracture edge in Option #6. Particles are placed according to input flux weights. In this case, all the particles have the same flux weight, but number of particles per fracture edge depends on in-flow flux of this fracture.
Definition at line 939 of file InitialPartPositions.c.
References contam::cell, vertex::coord_xy, density, contam::fl_weight, vertex::flux, vertex::fracture, contam::fracture, contam::intcell, node, nodezonein, np, npart, vertex::numneighb, particle, contam::position, contam::pressure, contam::time, and contam::velocity.
Referenced by InitPos().
int InitParticles_np | ( | int | k_current, |
int | firstnd, | ||
int | lastnd, | ||
int | parts_fracture, | ||
int | first_ind, | ||
int | last_ind | ||
) |
Function defines particle's initial positions on a single fracture edge. Option #1, init_nf
Definition at line 578 of file InitialPartPositions.c.
References vertex::coord_xy, contam::fl_weight, vertex::fracture, contam::fracture, contam::intcell, node, nodezonein, npart, particle, contam::position, contam::pressure, contam::time, and contam::velocity.
Referenced by InitPos().
int InitParticles_ones | ( | int | k_current, |
double | inter_p[][4], | ||
int | fracture_n, | ||
int | parts_fracture, | ||
int | ii, | ||
double | thirdcoor, | ||
int | zonenumb_in, | ||
int | first_ind, | ||
int | last_ind | ||
) |
Function defines particles initial positions in option #3, where user defines region at in-flow boundary face for particles.
Definition at line 668 of file InitialPartPositions.c.
References contam::fl_weight, contam::fracture, fracture, contam::intcell, npart, particle, contam::position, contam::pressure, material::rot2mat, contam::time, and contam::velocity.
Referenced by InitPos().
int InitPos | ( | ) |
Function defines the required option of particles initial positions defined at input control file; calculates number of particles, allocates memory.
First option init_nf: equal number of particles on every boundary edge regardless of edge length
Second option init_eqd: calculate total length of boundary edges; define the distance between particles and place particles equidistant from each other on all edges
Third option init_oneregion: user specifies a region and all particles start from the fracture edges that located inside the region
Fourth option init_random: particles will be set randomly over all fractures surfaces; in this case the particles are weighted acoording to the aperture
5th option init_matrix: particles are set randomly in rock matrix, their time from initial positions is calculated, then particles move through fractures to out-flow boundary. This option requires preprocessing; the script RandomPositGener.py will provide necessary input files
Sixth option init_flux: place particles according to input fluxes. Each particle has the same flux weight. The number of particles are placed propportionally to inflow flux on each cell and fracture edge.
Definition at line 21 of file InitialPartPositions.c.
References vertex::aperture, contam::cell, cell, Control_Data(), Control_File(), Control_File_Optional(), Control_Param(), vertex::coord, inpfile::filename, contam::fl_weight, inpfile::flag, flag_w, vertex::fracture, element::fracture, contam::fracture, fracture, InitInMatrix(), InitInWell(), InitParticles_eq(), InitParticles_flux(), InitParticles_np(), InitParticles_ones(), InsideCell(), Moving2Center(), ncells, node, element::node_ind, nodezonein, npart, nzone_in, inpfile::param, particle, contam::pressure, contam::time, totalFluxIn, and contam::velocity.
Referenced by ParticleTrack().
double InOutFlowCell | ( | int | indcell, |
int | int1, | ||
double | nposx, | ||
double | nposy | ||
) |
Function defines if velocities on cell vertices pointing in or out of intersection line
Definition at line 1730 of file TrackingPart.c.
References CalculateWeights(), contam::cell, cell, ChangeFracture(), vertex::coord_xy, contam::fracture, node, element::node_ind, np, particle, contam::position, element::veloc_ind, vertex::velocity, intcoef::weights, XindexC(), and YindexC().
Referenced by CheckDistance().
int InsideCell | ( | unsigned int | numc | ) |
Function checks if particle is in the cell (numc is cell ID)
Definition at line 1324 of file TrackingPart.c.
References CalculateWeights(), contam::cell, cell, contam::intcell, node, element::node_ind, np, particle, contam::weight, and intcoef::weights.
Referenced by InitCell(), InitPos(), Moving2Center(), and NeighborCells().
struct matr MatrixProducts | ( | double | normxarea[][2], |
int | number | ||
) |
Function performs matrix dot product (GTG)-1 and returns the result.
Definition at line 214 of file VelocityReconstruction.c.
Referenced by VelocityExteriorNode(), and VelocityInteriorNode().
void Moving2Center | ( | int | nnp, |
int | cellnumber | ||
) |
Function moves particle to the center of the same cell
Definition at line 2175 of file TrackingPart.c.
References cell, vertex::coord_xy, InsideCell(), node, element::node_ind, particle, contam::position, Xindex(), and Yindex().
Referenced by InitPos(), Moving2NextCell(), and Moving2NextCellBound().
int Moving2NextCell | ( | int | stuck, |
int | k | ||
) |
Functions performs the movement of particle from one cell to the center of neighbouring cell.
Definition at line 2198 of file TrackingPart.c.
References contam::cell, cell, vertex::cells, fracture, Moving2Center(), node, np, and particle.
void Moving2NextCellBound | ( | int | prevcell | ) |
In the pathological rare case, when particle is out of fracture, the function is called and it's moving particle to internal cell
Definition at line 2286 of file TrackingPart.c.
References contam::cell, cell, vertex::cells, fracture, Moving2Center(), node, element::node_ind, np, and particle.
Referenced by ParticleTrack().
void NeighborCells | ( | int | k | ) |
Function checks neighboring cells to find a particle
Definition at line 1426 of file TrackingPart.c.
References vertex::cells, fracture, InsideCell(), node, np, and particle.
Referenced by SearchNeighborCells().
FILE* OpenFile | ( | char | filen[120], |
char | fileopt[2] | ||
) |
The function opens file for reading or writing. If error - the program is terminated.
Definition at line 664 of file ReadGridInit.c.
Referenced by CalculateVelocity3D(), Coordinations2D(), InitInMatrix(), OutputMarPlumDisp(), ParticleOutput(), ParticleTrack(), ReadAperture(), ReadBoundaryNodes(), ReadDataFiles(), ReadFEHMfile(), ReadInit(), ReadPFLOTRANfile(), and WritingInit().
void OutputMarPlumDisp | ( | int | currentnum, |
char | path[125] | ||
) |
This function is used when particles trajectories should be output in format that MARFA and/or PLUMECALC codes will be able to read and process.
Definition at line 19 of file output.c.
References cell, Control_File(), inpfile::filename, maindir, marfa, OpenFile(), plumec, timeunit, and traj_o.
Referenced by ParticleTrack().
void OutputVelocities | ( | ) |
This function is called to output velocities on nodes in the order of fractture IDs. The output can be used to visualise flow field of 2D reconstructed Darcy velocities in 3D simulation domain
Definition at line 214 of file VelocityReconstruction.c.
void ParticleOutput | ( | int | currentt, |
int | fract_p | ||
) |
The function of particles trajectories outputs. Function is called at every intersection and outputs to file at each segment of particles trajectory: from intersection to intersection. The curvature of the trajectory is defined and dictate number of time steps for outputs (unless user requested every time step output).
Definition at line 2367 of file TrackingPart.c.
References vertex::aperture, avs_o, tempout::betap, CalculatePosition3D(), CalculateVelocity3D(), contam::cell, cell, tempout::cellp, posit3d::cord3, curv_o, DefineAngle(), contam::fracture, fracture, tempout::fracturep, tempout::length_t, maindir, no_out, node, element::node_ind, nodeID, np, OpenFile(), particle, pi, contam::position, tempout::position2d, tempout::position3d, contam::pressure, tempout::pressure, t, tempdata, tfile, contam::time, timecounter, tempout::timep, tempout::times, traj_o, and tempout::velocity3d.
Referenced by CheckDistance(), and ParticleTrack().
void ParticleTrack | ( | ) |
The main driving function of particles tracking procedure.
Definition at line 66 of file TrackingPart.c.
References all_out, vertex::aperture, avs_o, tempout::betap, lagrangian::betta, CalculateLagrangian(), CalculatePosition3D(), CalculateVelocity3D(), contam::cell, cell, tempout::cellp, CheckDistance(), CheckNewCell(), Control_File(), Control_File_Optional(), Control_Param(), vertex::coord, posit3d::cord3, CorrectorStep(), curv_o, disp_o, inpfile::filename, FinalPosition(), contam::fl_weight, inpfile::flag, FLAG_OUT, flag_w, FlowInWeight(), frac_o, vertex::fracture, contam::fracture, fracture, tempout::fracturep, InitCell(), InitPos(), lagrangian::initx, lagrangian::inity, lagrangian::initz, lagvariable, tempout::length_t, maindir, marfa, mixing_rule, Moving2NextCellBound(), nfract, no_out, node, element::node_ind, nodeID, nodezonein, nodezoneout, np, OpenFile(), OutputMarPlumDisp(), inpfile::param, particle, ParticleOutput(), plumec, contam::position, tempout::position2d, tempout::position3d, PredictorStep(), contam::pressure, tempout::pressure, String_Compare(), t, t_adv, t_adv0, contam::t_adv_diff, contam::t_diff, lagrangian::tau, tdrw, tdrw_diffcoeff, tdrw_lambda, tdrw_limited, tdrw_o, tdrw_porosity, tempdata, tfile, contam::time, timecounter, timediff, TimeDomainRW(), tempout::timep, tempout::times, timesteps, timeunit, traj_o, and tempout::velocity3d.
Referenced by main().
void PredictorStep | ( | ) |
Predictor step in Predictor-Corrector technique. Function calculates new velocities and new particle position.
Definition at line 1383 of file TrackingPart.c.
References CalculateCurrentDT(), contam::cell, cell, node, element::node_ind, np, particle, contam::position, vertex::pressure, contam::pressure, contam::prev_pos, element::veloc_ind, contam::velocity, vertex::velocity, and contam::weight.
Referenced by CheckDistance(), and ParticleTrack().
void ReadAperture | ( | ) |
The function reads apertures of fractures in DFN. Aperture can be defined for each fracture or for every node/cell, representing internal heterogeneity.
Definition at line 1203 of file ReadGridInit.c.
References vertex::aperture, Control_File(), inpfile::filename, vertex::fracture, fracture, nfract, nnodes, node, and OpenFile().
Referenced by ReadDataFiles().
void ReadBoundaryNodes | ( | ) |
The functions reads node's IDs assigned to in-flow and out-flow boundaries. the total input flow flux is calculated.
Definition at line 526 of file ReadGridInit.c.
References Control_Data(), Control_File(), density, inpfile::filename, inpfile::flag, vertex::flux, maindir, node, nodezonein, nodezoneout, vertex::numneighb, nzone_in, OpenFile(), totalFluxIn, and vertex::typeN.
Referenced by main().
void ReadDataFiles | ( | ) |
The function reads DFN mesh from inp and stor files
if no aperture file specified, the aperture of all fractures will be equal to thickness value
Definition at line 155 of file ReadGridInit.c.
References AdjacentCells(), vertex::aperture, vertex::area, cell, vertex::cells, Control_File(), fehm, inpfile::filename, material::firstcell, material::firstnode, vertex::flux, vertex::fracts, vertex::fracture, element::fracture, fracture, vertex::indnodes, material::lastnode, max_neighb, ncells, nfract, nnodes, node, material::numbcells, vertex::numneighb, OpenFile(), pflotran, ReadAperture(), ReadFEHMfile(), ReadPFLOTRANfile(), thickness, vertex::type, vertex::typeN, and element::veloc_ind.
Referenced by main().
void ReadFEHMfile | ( | int | nedges | ) |
The function opens and reads FEHM outputs; read in flow fluxes and cell volumes.
Definition at line 796 of file ReadGridInit.c.
References Control_File(), inpfile::filename, max_neighb, nnodes, node, OpenFile(), and vertex::pressure.
Referenced by ReadDataFiles().
void ReadInit | ( | ) |
The function reads total number of nodes, triangular cells, fractures in DFN mesh. The memory is allocated for data structures: NODE, CELL, FRACTURE
Definition at line 16 of file ReadGridInit.c.
References vertex::area, cell, vertex::cells, Control_File(), Control_File_Optional(), inpfile::filename, inpfile::flag, vertex::flux, vertex::fracts, fracture, vertex::indnodes, max_neighb, ncells, nfract, nnodes, node, material::nvect_xy, material::nvect_z, OpenFile(), material::theta, and vertex::type.
Referenced by main().
void ReadPFLOTRANfile | ( | int | nedges | ) |
The function opens and reads PFLOTRAN files, read in flow fluxes, areas, pressure.
Definition at line 677 of file ReadGridInit.c.
References vertex::aperture, vertex::area, Control_File(), density, inpfile::filename, inpfile::flag, vertex::flux, nnodes, node, vertex::numneighb, OpenFile(), vertex::pressure, and vertex::pvolume.
Referenced by ReadDataFiles().
void SearchNeighborCells | ( | int | nn1, |
int | nn2, | ||
int | nn3 | ||
) |
Function performs a search of neighbouring cells of current particles position.
Definition at line 1285 of file TrackingPart.c.
References cell, NeighborCells(), np, and particle.
Referenced by CalculateWeights().
int StreamlineRandomSampling | ( | double | products[4], |
double | speedsq[4], | ||
int | indj, | ||
int | int1, | ||
int | indk, | ||
int | neighborcellind[4], | ||
int | neighborfracind[4], | ||
int | prevfrac, | ||
int | prevcell | ||
) |
Particle motion at a fracture intersection is determined by the streamline routing rule Arg 1: Cross product to define outgoing and incoming cells. Vector contains value for each cell at the intersection Arg 2: Vector of the speed squared for each cell at the intersection Arg 3: index of current cell in int1 node list Arg 4: int1 is the index of the node of the incoming cell that lies on the intersection Arg 5: indk gives the position in the list of the 4 intersection cells that is previous cell Arg 6: Vector of the indicies for the 4 neighboring cells Arg 7: Vector of the fracture index for the 4 neighboring cells Arg 8: Fracture a particle is coming from Arg 9: Cell index a particle is coming from Return: The Exit cell at the intersection
Definition at line 1901 of file TrackingPart.c.
References cell, vertex::cells, fracture, and node.
Referenced by AcrossIntersection().
int String_Compare | ( | char | string1[], |
char | string2[] | ||
) |
The finction performs a string comparison, used in the control file reading.
Definition at line 1089 of file ReadGridInit.c.
Referenced by main(), and ParticleTrack().
double TimeDomainRW | ( | double | time_advect | ) |
Time Domain Random Walk (TDRW) procedure to account for matrix diffusion. Returns a diffusion time of particle per fracture. This function is called at each intersection.
Details on the method can be found in
Hyman, Jeffrey D., Harihar Rajaram, Shriram Srinivasan, Nataliia Makedonska, Satish Karra, Hari Viswanathan, and Gowri Srinivasan. "Matrix Diffusion in Fractured Media: New Insights Into Power Law Scaling of Breakthrough Curves." Geophysical Research Letters 46, no. 23 (2019): 13785-13795.
Definition at line 2769 of file TrackingPart.c.
References vertex::aperture, contam::cell, cell, material::firstnode, contam::fracture, fracture, node, element::node_ind, np, particle, pi, tdrw_diffcoeff, tdrw_lambda, tdrw_limited, tdrw_porosity, and timediff.
Referenced by CheckDistance(), and ParticleTrack().
double TimeFromMatrix | ( | double | pdist | ) |
Option #5. Estimation of travel time of particles moving from matrix to the closest fracture.
Definition at line 913 of file InitialPartPositions.c.
References Control_Param(), inpfile::param, pi, and timeunit.
Referenced by InitInMatrix().
void Velocity3D | ( | ) |
Recalculates 2D velocities at XY fracture plane to 3D velocties in the simulation domain. This procedure is not used for particle tracking, but can be used for velocity field visualization.
Definition at line 190 of file RotateFracture.c.
Referenced by main().
void VelocityExteriorNode | ( | double | norm_xarea[][2], |
int | i, | ||
int | number, | ||
unsigned int | indj[max_neighb], | ||
struct lb | lbound, | ||
int | vi | ||
) |
The function reconstructs Darcy velocity on exterior (boundary) node
Definition at line 391 of file VelocityReconstruction.c.
References density, vertex::flux, lb::length_b, matr::matinv, matr::matrinvG, MatrixProducts(), node, lb::norm, porosity, timeunit, and vertex::velocity.
Referenced by DarcyVelocity(), and HalfPolygonVelocity().
void VelocityInteriorNode | ( | double | normx_area[][2], |
int | i, | ||
int | number, | ||
unsigned int | indj[max_neighb], | ||
int | vi | ||
) |
The function reconstructs Darcy velocity on interior cell center, at interior node
Definition at line 356 of file VelocityReconstruction.c.
References density, vertex::flux, matr::matinv, matr::matrinvG, MatrixProducts(), node, porosity, timeunit, and vertex::velocity.
Referenced by DarcyVelocity(), and HalfPolygonVelocity().
void WritingInit | ( | ) |
Functions write the data structure of nodes, cells and fractures into files in ASCII format with detail explanations. This output is optional and can be helpful in code debugging.
Definition at line 898 of file ReadGridInit.c.
References cell, inpfile::filename, fracture, maindir, nfract, nnodes, node, material::numbcells, vertex::numneighb, and OpenFile().
Referenced by main().
int Xindex | ( | int | nodenum, |
int | nnp | ||
) |
Functions returns the index of X coordination of intersection node
Definition at line 2254 of file TrackingPart.c.
References fracture, node, and particle.
Referenced by CalculateLagrangian(), CheckDistance(), and Moving2Center().
int XindexC | ( | int | nodenum, |
int | ii | ||
) |
Functions returns the index of X coordination of intersection node (cells)
Definition at line 776 of file VelocityReconstruction.c.
References cell, fracture, and node.
Referenced by AcrossIntersection(), BoundaryCells(), CornerVelocity(), and InOutFlowCell().
int Yindex | ( | int | nodenum, |
int | nnp | ||
) |
Functions returns the index of Y coordination of intersection node
Definition at line 2270 of file TrackingPart.c.
References fracture, node, and particle.
Referenced by CalculateLagrangian(), CheckDistance(), and Moving2Center().
int YindexC | ( | int | nodenum, |
int | ii | ||
) |
Functions returns the index of Y coordination of intersection node (cells)
Definition at line 791 of file VelocityReconstruction.c.
References cell, fracture, and node.
Referenced by AcrossIntersection(), BoundaryCells(), CornerVelocity(), and InOutFlowCell().
|
extern |
DYNAMIC ARRAY OF TRIANGULAR CELLS in DFN mesh
Definition at line 42 of file main.c.
Referenced by AcrossIntersection(), AdjacentCells(), BoundaryCells(), CalculateCurrentDT(), CalculateLagrangian(), CalculatePosition3D(), CalculateWeights(), CheckDistance(), CorrectorStep(), FlowInWeight(), HalfPolygonVelocity(), InitPos(), InOutFlowCell(), InsideCell(), main(), Moving2Center(), Moving2NextCell(), Moving2NextCellBound(), OutputMarPlumDisp(), ParticleOutput(), ParticleTrack(), PredictorStep(), ReadDataFiles(), ReadInit(), SearchNeighborCells(), StreamlineRandomSampling(), TimeDomainRW(), WritingInit(), XindexC(), and YindexC().
|
extern |
|
extern |
Flow density, defined by user
Definition at line 47 of file main.c.
Referenced by InitParticles_flux(), main(), ReadBoundaryNodes(), ReadPFLOTRANfile(), VelocityExteriorNode(), and VelocityInteriorNode().
|
extern |
fehm=1 when FEHM flow solver is used
Definition at line 28 of file main.c.
Referenced by DarcyVelocity(), HalfPolygonVelocity(), and ReadDataFiles().
|
extern |
Definition at line 38 of file main.c.
Referenced by InitPos(), and ParticleTrack().
|
extern |
DYNAMIC ARRAY OF FRACTURES
Definition at line 39 of file main.c.
Referenced by AdjacentCells(), BoundaryCells(), CalculatePosition3D(), CalculateVelocity3D(), CheckDistance(), CheckGrid(), Convertto2d(), Convertto3d(), Coordinations2D(), DarcyVelocity(), DefineBoundaryAngle(), HalfPolygonVelocity(), InitParticles_ones(), InitPos(), main(), Moving2NextCell(), Moving2NextCellBound(), NeighborCells(), ParticleOutput(), ParticleTrack(), ReadAperture(), ReadDataFiles(), ReadInit(), StreamlineRandomSampling(), TimeDomainRW(), WritingInit(), Xindex(), XindexC(), Yindex(), and YindexC().
|
extern |
The directory/path for particle tracking outputs, defined by user
Definition at line 43 of file main.c.
Referenced by CalculateVelocity3D(), DefineBoundaryAngle(), main(), OutputMarPlumDisp(), ParticleOutput(), ParticleTrack(), ReadBoundaryNodes(), and WritingInit().
|
extern |
maximum number of nodes neighbours/connections
Definition at line 32 of file main.c.
Referenced by DarcyVelocity(), HalfPolygonVelocity(), ReadDataFiles(), ReadFEHMfile(), and ReadInit().
|
extern |
total number of triangular elements in the DFN mesh
Definition at line 30 of file main.c.
Referenced by BoundaryCells(), InitPos(), ReadDataFiles(), and ReadInit().
|
extern |
total number of fractures in the DFN mesh
Definition at line 31 of file main.c.
Referenced by Convertto2d(), Convertto3d(), Coordinations2D(), DefineBoundaryAngle(), ParticleTrack(), ReadAperture(), ReadDataFiles(), ReadInit(), and WritingInit().
|
extern |
total number of nodes in the DFN mesh
Definition at line 29 of file main.c.
Referenced by CalculateVelocity3D(), CheckGrid(), Convertto2d(), Coordinations2D(), DarcyVelocity(), DefineBoundaryAngle(), DefineTimeStep(), ReadAperture(), ReadDataFiles(), ReadFEHMfile(), ReadInit(), ReadPFLOTRANfile(), and WritingInit().
|
extern |
DYNAMIC ARRAY OF NODES in DFN mesh
Definition at line 40 of file main.c.
Referenced by AcrossIntersection(), AdjacentCells(), BoundaryCells(), CalculateCurrentDT(), CalculateLagrangian(), CalculateVelocity3D(), CalculateWeights(), CheckDistance(), CheckGrid(), CompleteMixingRandomSampling(), Convertto2d(), Coordinations2D(), CornerVelocity(), DarcyVelocity(), DefineBoundaryAngle(), DefineTimeStep(), FlowInWeight(), HalfPolygonVelocity(), InitCell(), InitInMatrix(), InitInWell(), InitParticles_eq(), InitParticles_flux(), InitParticles_np(), InitPos(), InOutFlowCell(), InsideCell(), main(), Moving2Center(), Moving2NextCell(), Moving2NextCellBound(), NeighborCells(), ParticleOutput(), ParticleTrack(), PredictorStep(), ReadAperture(), ReadBoundaryNodes(), ReadDataFiles(), ReadFEHMfile(), ReadInit(), ReadPFLOTRANfile(), StreamlineRandomSampling(), TimeDomainRW(), VelocityExteriorNode(), VelocityInteriorNode(), WritingInit(), Xindex(), XindexC(), Yindex(), and YindexC().
|
extern |
pointer to the dynamic array with a list of in-flow boundary nodes
Definition at line 36 of file main.c.
Referenced by InitCell(), InitInWell(), InitParticles_flux(), InitParticles_np(), InitPos(), main(), ParticleTrack(), and ReadBoundaryNodes().
|
extern |
pointer to the dynamic array with a list of out-flow boundary nodes
Definition at line 37 of file main.c.
Referenced by main(), ParticleTrack(), and ReadBoundaryNodes().
|
extern |
index of current particle, idex in particle's loop
Definition at line 34 of file main.c.
Referenced by AcrossIntersection(), CalculateCurrentDT(), CalculateLagrangian(), CalculatePosition3D(), CalculateWeights(), CheckDistance(), CorrectorStep(), FlowInWeight(), InitParticles_flux(), InOutFlowCell(), InsideCell(), Moving2NextCell(), Moving2NextCellBound(), NeighborCells(), ParticleOutput(), ParticleTrack(), PredictorStep(), SearchNeighborCells(), and TimeDomainRW().
|
extern |
initial number of particles set up in the simulation
Definition at line 33 of file main.c.
Referenced by InitInMatrix(), InitParticles_eq(), InitParticles_flux(), InitParticles_np(), InitParticles_ones(), and InitPos().
|
extern |
number of nodes in in-flow boundary face/zone
Definition at line 35 of file main.c.
Referenced by InitCell(), InitInWell(), InitPos(), and ReadBoundaryNodes().
|
extern |
DYNAMIC ARRAY OF PARTICLES
Definition at line 41 of file main.c.
Referenced by AcrossIntersection(), CalculateCurrentDT(), CalculateLagrangian(), CalculatePosition3D(), CalculateWeights(), CheckDistance(), CorrectorStep(), FlowInWeight(), InitInMatrix(), InitInWell(), InitParticles_eq(), InitParticles_flux(), InitParticles_np(), InitParticles_ones(), InitPos(), InOutFlowCell(), InsideCell(), main(), Moving2Center(), Moving2NextCell(), Moving2NextCellBound(), NeighborCells(), ParticleOutput(), ParticleTrack(), PredictorStep(), SearchNeighborCells(), TimeDomainRW(), Xindex(), and Yindex().
|
extern |
pflotran=1 when PFLOTRAN flow solver is used
Definition at line 27 of file main.c.
Referenced by DarcyVelocity(), HalfPolygonVelocity(), and ReadDataFiles().
|
extern |
Fracture porosity, used in velocity reconstructions, defined by user
Definition at line 46 of file main.c.
Referenced by main(), VelocityExteriorNode(), and VelocityInteriorNode().
|
extern |
|
extern |
One value for all fractures aperture (defined by user), used in case when fracture aperture is not provided by user
Definition at line 49 of file main.c.
Referenced by main(), and ReadDataFiles().
|
extern |
Max number of time steps used for each particles movements, defined by user
Definition at line 48 of file main.c.
Referenced by main(), and ParticleTrack().
|
extern |
Time unit multiplier, converts calculated time/velocities according to required time units
Definition at line 51 of file main.c.
Referenced by main(), OutputMarPlumDisp(), ParticleTrack(), TimeFromMatrix(), VelocityExteriorNode(), and VelocityInteriorNode().
|
extern |
Total Flow flux on in-flow boundary, calculated in the code
Definition at line 52 of file main.c.
Referenced by InitPos(), and ReadBoundaryNodes().