This repository was archived by the owner on Jan 17, 2024. It is now read-only.
-
Notifications
You must be signed in to change notification settings - Fork 0
Expand file tree
/
Copy pathwavempi.c
More file actions
136 lines (106 loc) · 3.66 KB
/
wavempi.c
File metadata and controls
136 lines (106 loc) · 3.66 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
#include <mpi.h>
#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <unistd.h>
#include "globals.h"
#include "parsempi.h"
#include "utils.h"
#include "tic.h"
#include "velocity_model.h"
#include "ricker.h"
#include "wavefield.h"
#include "wavefieldmpi.h"
#include "simulation.h"
#include "simulationmpi.h"
#include "domain_comm.h"
int main(int argc, char *argv[])
{
args *ap = arg_parse(argc, argv);
if(isstable_from_param(ap->frequency, ap->dx, ap->dz, ap->vel, ap->vel, ap->dt) != SIMUL_OK)
{
fprintf(stderr, "Unstable model\n");
exit(EXIT_FAILURE);
}
/* Create simulation */
simulation_params *simulation = simulation__create();
simulation->time = ap->time;
simulation->sample = ap->sample;
simulation->dt = ap->dt;
simulation->steps = simulation->time / simulation->dt + 1;
simulation->ntrec = simulation->sample/simulation->dt;
/* Pre calculate Laplacian parameters */
laplacian_params *lp = wavefield__laplacian_params_(ap->dx, ap->dz, ap->order, simulation->dt);
/***************************************************
*
* MPI Stuff
*
***************************************************/
MPI_Init(&argc,&argv);
int world_size;
MPI_Comm_size(MPI_COMM_WORLD, &world_size);
int world_rank;
MPI_Comm_rank(MPI_COMM_WORLD, &world_rank);
size_t number_slices_x = 2;
size_t number_slices_z = world_size/number_slices_x;
if (world_rank == 0 && world_size % number_slices_x != 0) {
fprintf(stderr, "Right now only multiples of 2 workers are supported\n");
MPI_Abort(MPI_COMM_WORLD, 1);
}
domain_info *domain = malloc(sizeof *domain);
decompose_domain(world_size, world_rank,
number_slices_x, number_slices_z,
lp->border_size, lp->border_size,
ap->nx, ap->nz,
ap->sx, ap->sz,
domain);
ricker_wavelet *wavelet;
ricker_source *ricker;
// We need to reposition the source to the local domain sx and sz
if(domain->has_source) {
wavelet = ricker__read_from_file( ap->sourcefile, ap->frequency, 0);
ricker = ricker__model(wavelet, ap->sx, ap->sz, 0);
ricker->x = domain->lsx;
ricker->z = domain->lsz;
}
// Create a model that is part of full model based on local domain grid
velocity_model *model = velocity_model__create_submodel_from_file(
ap->modelfile,
ap->nx, ap->nz, ap->dx, ap->dz,
domain->nxa, domain->nza, domain->nxb, domain->nzb);
wavefield *P = wavefield__create(domain->lnx, domain->lnz);
/* Iterate timesteps */
for(size_t it = 0; it < simulation->steps; it++)
{
if(domain->has_source)
simulation__inject_source(P, model, ricker, simulation, it);
wavefield__laplacian(P, model, lp);
wavefieldmpi__perfect_match_layer_by_domain(P, model, lp, domain);
wavefieldmpi__share_border_by_domain(P, lp, domain);
wavefield__swap(P);
simulationmpi__write(it, P, simulation, domain, ticprt);
}
/* Memory cleanup */
fflush(stdout);
wavefield__destroy(P);
velocity_model__destroy(model);
if(domain->has_source) {
ricker__destroy_source(ricker);
}
free(domain);
MPI_Finalize();
wavefield__destroy_laplacian_params(lp);
simulation__destroy(simulation);
free(ap);
// TODO: THIS
simulationmpi__merge_wavefields(domain);
/* Finishing */
if(ticprt) {
tic();
}
if(verbose)
fprintf(stderr, "xmovie n1=%zu n2=%zu d1=%lf d2=%lf clip=0.5 loop=2\n",
ap->nx, ap->nz, ap->dx, ap->dz);
fflush(stderr);
return 0;
}