traceio.h File Reference

#include "utils.h"
#include "usempi.h"
#include "cubic.h"
#include "iwave_fopen.h"
#include "su.h"
#include "header.h"
#include "segy.h"

Go to the source code of this file.

Classes

struct  s_offsegy
struct  s_tracegeom

Defines

#define IWAVE_USE_FMGR
#define MAX_TRACES   500000
#define MAX_RECS   1000
#define SRC_TOL   0.001
#define W_EOF   -10

Typedefs

typedef s_offsegy offsegy
typedef s_tracegeom tracegeom

Functions

void mygethdval (segy *tr, const char *ch, Value *val)
void myputhdval (segy *tr, const char *ch, Value *val)
int construct_tracegeom (tracegeom *tg, const char *fin, float dt, float tol, FILE *stream)
 Trace geometry constructor.
int init_tracegeom (tracegeom *tg, int irec, RPNT og, IPNT n, RPNT d, RPNT o, IPNT axord, int order, int ndim, int initbuf, FILE *stream)
 Initialization for trace geometry object.
void destroy_tracegeom (tracegeom *tg)
 destructor: frees all dynamically allocated memory, closes files.
void setnull_tracegeom (tracegeom *tg)
 default initialization (as distinct from construction)- sets tracelength to zero, all other internals to zero or other sensible default values.
void fprint_tracegeom (tracegeom const *tg, FILE *fp)
 diagnostic output to stream (fp).
void print_tracegeom (tracegeom const *tg)
 diagnostic output to stdout.
void sampletraces (tracegeom *tg, int order, int load, int it, IPNT allocstrides, IPNT allocorigs, IPNT strides, IPNT origs, ireal *field, ireal mult)
 Sampling function.
void tapermutetraces (tracegeom *tg, int it, float dx, int width, int wtime, float mslope, float mzotime, float mwidth)
int writetraces (tracegeom const *tg, RPNT d, RPNT og, FILE *stream)
 Trace output.
void calc_group (int *first, int *last, int nrec)
void indexrange (tracegeom const *tg, IPNT gmin, IPNT gmax)
 returns index range of receiver array.


Define Documentation

#define IWAVE_USE_FMGR

Definition at line 5 of file traceio.h.

#define MAX_TRACES   500000

Definition at line 26 of file traceio.h.

#define MAX_RECS   1000

Definition at line 28 of file traceio.h.

#define SRC_TOL   0.001

Definition at line 30 of file traceio.h.

#define W_EOF   -10

Definition at line 32 of file traceio.h.


Typedef Documentation

typedef struct s_offsegy offsegy

typedef struct s_tracegeom tracegeom

tracegeom

Trace Geometry. Consists of struct defining sampling etc. and storing essential trace header and i/o status information, together with constructor, initializer, destructor, sampler, writer, and printer functions. The struct stores information for only one shot record at a time. Constructor connects the object to header (trace input) and data (trace output) files, which may be the same. Initializer populates the struct with data from the next shot record, read from the header file; this data includes all essential header words and may include trace samples. Sampler interacts with a grid, extracting samples at trace locations and storing them in the trace buffer data member of the struct. Writer dumps the current trace buffer to the data file. These functions rely on an internal traceserver object, which hides details of parallelism and data distribution.


Function Documentation

void mygethdval ( segy tr,
const char *  ch,
Value val 
)

void myputhdval ( segy tr,
const char *  ch,
Value val 
)

int construct_tracegeom ( tracegeom tg,
const char *  fin,
float  dt,
float  tol,
FILE *  stream 
)

Trace geometry constructor.

Calling this function

Parameters:
[out] tg - (tracegeom *) trace geometry object to be constructed
[in] fin - (char *) name of header (input) file.
[in] dt - (float) time step for internal time loop, used to reserve buffer space
[in] tol - (float) tolerance for distinguishing source positions with l-infty norm; should be computed in calling unit as fraction of maximum space step of grid to be sampled, possibly using cpp macro SRC_TOL defined in the header file.
[in] stream - (FILE *) output stream for diagnostic information.
Returns:
0 on successful completion, else error code as in base/include/utils.h.

int init_tracegeom ( tracegeom tg,
int  irec,
RPNT  og,
IPNT  n,
RPNT  d,
RPNT  o,
IPNT  axord,
int  order,
int  ndim,
int  initbuf,
FILE *  stream 
)

Initialization for trace geometry object.

Takes spatial geometry and start and end time from input file. Time sample rate specified externally (arg dt), determines buffer allocation. Time unit is ms. It is ASSUMED that any common spatial data (eg. scalel, scalco) is uniform, and may be read from any trace.

The data to be sampled repeatedly to form traces are assumed to functions on a uniform rectangular grid. The construction of samples differentiates between a LOCAL grid, meaning the grid on which the data to be sampled are being updated by the process, and a GLOBAL grid, of which the LOCAL grid is a subgrid. Index tuples are defined by the GLOBAL grid. The sampling function (below) is supplied with the index tuple of the LOCAL grid origin, which enables conversion from global indices to local offset.

This function does several tasks, using the tracegeom struct to record its results:

The routine may be called several times to initialize any number of time series on the same internal time grid. This internal (stepping) time grid is fixed by the requirement tha (1) the step be given by the input argument dt, and (2) t=0 is a sample point.

arguments:

Parameters:
[out] tg (tracegeom *) trace geometry object to be initialized
[in] og (RPNT) coordinates of GLOBAL grid origin (grid point with all indices = 0) relative to physical coordinate system. The physical coordinate system refers to the physical domain, and is the coordinate system with respect to which the source and receiver coordinates are defined. These numbers (lengths, not indices!) give the coordinates of the global computational grid origin in the physical system.
[in] n (IPNT) number of samples on each axis for the LOCAL subgrid in which traces are to be sampled.
[in] o (RPNT) coordinates of LOCAL subgrid origin.
[in] d (RPNT) spatial steps (common to GLOBAL grid and LOCAL subgrid)
[in] dt (float) time step in ms. input argument, NOT read from file - computed by the simulator, actual time step between time levels of simulation. Recorded in tracegeom object as tg->dt.
[in] ndim (int) dimension of GLOBAL and LOCAL grids.
[in] usernt,usert0 (int, float) optional specification of number output samples and time of first output sample. If set (i.e. usernt > 0), then
  1. output sample rate (tg->dtout) set to simulator-supplied rate (dt);
  2. number of output samples (tg->ntout) set to usernt;
  3. number of simulation time steps (tg->nt) set to usernt;
  4. time of first output sample (tg->t0out) is set to usert0, possibly adjusted so that t=0 is a sample;
  5. time of first simulation sample (tg->t0) also set to usert0.
[this option is mostly useful in convergence tests and the like]. If NOT set (i.e. usernt <= 0) then
  1. output sample rate (tg->dtout) read from file, and used together with cubic interpolation (see Cubic Spline Interpolation and Adjoint Interpolation) to write traces;
  2. number of output samples (tg->ntout) read from file (usernt ignored);
  3. time of first output sample (tg->t0out) read from file (usert0 ignored);
  4. number of simulation time steps (tg->nt) computed from dt, tg->dtout, tg->t0out, and tg->ntout;
  5. time of first simulation sample (tg->t0) computed from tg->t0out and dt so that simulation grid includes t=0 as a sample point.
Thus in the second case (usernt NOT set) the simulator output is interpolated onto the time sampling of the file used to initialized the tracegeom.
[in] initbuf (int) flag for data initiation. If set, data traces are read and transferred to the internal simulation (tg->dt,tg->nt,tg->t0) grid by cubic spline interpolation (initbuf > 0) or adjoint cubic interpolation (initbut < 0) (see Cubic Spline Interpolation and Adjoint Interpolation). To avoid temp storage, this option implemented by another pass through file.
ADDED 22.10.08: AXIS ORDERING

Parameters:
[in] axord (int) The axord argument codes axis order as follows:
  1. for 1D, who cares - z index is zero in any case, so check that axord[0]=0;
  2. for 2D, axord[0]=index of z, axord[1]=index of x - check that values are legit;
  3. for 3D, axord[0]=index of z, axord[1]=index of x, axord[2]=index of y - check that values are legit.
[in] order (int) interpolation order for non-grid source and receiver positions - current legal values are 0 and 1
[in] stream (FILE *) verbose output stream
Returns:
0 on successful completion, else error code as in base/include/utils.h.

void destroy_tracegeom ( tracegeom tg  ) 

destructor: frees all dynamically allocated memory, closes files.

Should be called as part of "default constructor"

void setnull_tracegeom ( tracegeom tg  ) 

default initialization (as distinct from construction)- sets tracelength to zero, all other internals to zero or other sensible default values.

Leaves files open, does not reset whole-data-set dependent quantities. Should be called to initialized structure for new record.

void fprint_tracegeom ( tracegeom const *  tg,
FILE *  fp 
)

diagnostic output to stream (fp).

void print_tracegeom ( tracegeom const *  tg  ) 

diagnostic output to stdout.

void sampletraces ( tracegeom tg,
int  order,
int  load,
int  it,
IPNT  allocstrides,
IPNT  allocorigs,
IPNT  strides,
IPNT  origs,
ireal field,
ireal  mult 
)

Sampling function.

To be called at every time step. Index arrays defined relative to global grid. Requires both local (to each process) strides and origins of allocated grid, to compute offset into array to be sampled, and strides and origins of computational grid (array of points to be updated by each processor. Only those summands in interpolation formulas computed corresponding to these computational grid points - since computational grids do not overlap, each interpolation term is computed by at most one processor. If no processor computes a term, then point contributing belongs to no computational grid, i.e. lies on Dirichlet for the corresponding field. TODO: this idea certainly works for the pressure, but needs to be checked for sampling velocity in the near-boundary layer.

Version 1.5 addendum: the final argument is a scale factor, applied to every sample. For load (adjoint to save), an additional factor of cell volume ratio is applied - that is, time step / volume of spatial cell.

Parameters:
[out] tg (tracegeom) trace geometry - output for save, input for load. defines sampling of grid, samples saved to buffer tg.buf in save mode, loaded from tg.buf in load mode.
[in] order (int) order of interpolation rule.- either trilinear (order=1) or nearest neighbor (order=0 are supported;
[in] load (int) if set, transfer data from traces to grid by adj interp; else, transfer data from grid to traces by interp
[in] it (int) time step index for sampling - function is no-op if it < 0 or it > tg.nt-1;
[in] allocstrides (IPNT) number of samples on each local ALLOCATED grid axis;
[in] allocorigs (IPNT) index tuple of local ALLOCATED grid origin relative to global grid;
[in] strides (IPNT) number of samples on each local COMPUTATIONAL grid axis (updated points);
[in] origs (IPNT) index tuple of local COMPUTATIONAL grid origin relative to global grid;
[in] field (float *) (input for save, output for load) buffer containing grid data to be sampled.
[in] mult (ireal) scale factor for samples - for load (adjoint of save), scaled by ratio of cell volumes

void tapermutetraces ( tracegeom tg,
int  it,
float  dx,
int  width,
int  wtime,
float  mslope,
float  mzotime,
float  mwidth 
)

int writetraces ( tracegeom const *  tg,
RPNT  d,
RPNT  og,
FILE *  stream 
)

Trace output.

Sets up internal segy for writing, by transferring ntout, dtout, scalel, scalco, etc. to header fields. writes out traces. Uses axis steps and origin coordinates combined with index and relative coordinate data for source and receiver locations to compute trace headers. On return, fname will be SU file in working directory.

Parameters:
[in] tg (tracegeom) header and trace data to be sampled;
[in] d (RPNT) axis steps in grid;
[in] og (RPNT) coordinate tuple for origin of GLOBAL grid;
[in] stream (FILE *) verbose output stream
Returns:
0 on successful write, else error code as in base/include/utils.h.

void calc_group ( int *  first,
int *  last,
int  nrec 
)

void indexrange ( tracegeom const *  tg,
IPNT  gmin,
IPNT  gmax 
)

returns index range of receiver array.

Coordinate order is nominal, that is, not reinterpreted via axord but same as internal grid order. NOTE: gmin, gmax overwritten, not updated!!!

Parameters:
[in] tg - header and trace data to be sampled;
[out] gmin - min index array of receivers
[out] gmax - max index array of receivers


Generated on 5 Jan 2017 for IWAVETRACE by  doxygen 1.4.7