-
Notifications
You must be signed in to change notification settings - Fork 65
/
obj_psv.c
217 lines (160 loc) · 7.5 KB
/
obj_psv.c
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
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
/* --------------------------------------------------------------------------
* Calculate objective function for PSV problem
*
*
* D. Koehn
* Kiel, 03.06.2016
*
* --------------------------------------------------------------------------*/
#include "fd.h"
float obj_psv(struct wavePSV *wavePSV, struct wavePSV_PML *wavePSV_PML, struct matPSV *matPSV, struct fwiPSV *fwiPSV, struct mpiPSV *mpiPSV,
struct seisPSV *seisPSV, struct seisPSVfwi *seisPSVfwi, struct acq *acq, float *hc, int nsrc, int nsrc_loc, int nsrc_glob, int ntr,
int ntr_glob, int ns, int itest, int iter, float **Ws, float **Wr, int hin, int *DTINV_help, float eps_scale, MPI_Request * req_send, MPI_Request * req_rec){
/* global variables */
extern int RUN_MULTIPLE_SHOTS, TESTSHOT_START, TESTSHOT_END, TESTSHOT_INCR, N_STREAMER, SEISMO, QUELLART, QUELLTYP, ORDER_SPIKE;
extern int TIME_FILT, INV_STF, ORDER, L, MYID, LNORM, READREC, QUELLTYPB;
extern float FC_SPIKE_2,FC_SPIKE_1, FC, FC_START;
/* local variables */
float L2sum, L2_all_shots, energy_all_shots, energy_tmp, L2_tmp;
int ntr_loc, nt, ishot, nshots;
FILE *FP;
/* initialization of L2 calculation */
(*seisPSVfwi).L2=0.0;
(*seisPSVfwi).energy=0.0;
L2_all_shots=0.0;
energy_all_shots=0.0;
/* no differentiation of elastic and viscoelastic modelling because the viscoelastic parameters did not change during the forward modelling */
matcopy_elastic_PSV((*matPSV).prho,(*matPSV).ppi,(*matPSV).pu);
MPI_Barrier(MPI_COMM_WORLD);
av_mue((*matPSV).pu,(*matPSV).puipjp,(*matPSV).prho);
av_rho((*matPSV).prho,(*matPSV).prip,(*matPSV).prjp);
/* Preparing memory variables for update_s (viscoelastic) */
if (L) prepare_update_s_visc_PSV((*matPSV).etajm,(*matPSV).etaip,(*matPSV).peta,(*matPSV).fipjp,(*matPSV).pu,(*matPSV).puipjp,(*matPSV).ppi,(*matPSV).prho,(*matPSV).ptaus,(*matPSV).ptaup,
(*matPSV).ptausipjp,(*matPSV).f,(*matPSV).g,(*matPSV).bip,(*matPSV).bjm,(*matPSV).cip,(*matPSV).cjm,(*matPSV).dip,(*matPSV).d,(*matPSV).e);
if (RUN_MULTIPLE_SHOTS) nshots=nsrc; else nshots=1;
for (ishot=TESTSHOT_START;ishot<=TESTSHOT_END;ishot=ishot+TESTSHOT_INCR){
if(MYID==0){
printf("\n=================================================================================================\n");
printf("\n ***** Starting simulation (test-forward model) no. %d for shot %d of %d (rel. step length %.8f) \n",itest,ishot,nshots,eps_scale);
printf("\n=================================================================================================\n\n");
}
if((N_STREAMER>0)||(READREC==2)){
if (SEISMO){
(*acq).recpos=receiver(FP, &ntr, ishot);
(*acq).recswitch = ivector(1,ntr);
(*acq).recpos_loc = splitrec((*acq).recpos,&ntr_loc, ntr, (*acq).recswitch);
ntr_glob=ntr;
ntr=ntr_loc;
}
/* Memory for seismic data */
alloc_seisPSV(ntr,ns,seisPSV);
/* Memory for full data seismograms */
alloc_seisPSVfull(seisPSV,ntr_glob);
/* Memory for FWI seismic data */
alloc_seisPSVfwi(ntr,ntr_glob,ns,seisPSVfwi);
}
for (nt=1;nt<=8;nt++) (*acq).srcpos1[nt][1]=(*acq).srcpos[nt][ishot];
/* set QUELLTYP for each shot */
QUELLTYP = (*acq).srcpos[8][ishot];
if (RUN_MULTIPLE_SHOTS){
/* find this single source positions on subdomains */
if (nsrc_loc>0) free_matrix((*acq).srcpos_loc,1,8,1,1);
(*acq).srcpos_loc=splitsrc((*acq).srcpos1,&nsrc_loc, 1);
}
else{
/* Distribute multiple source positions on subdomains */
(*acq).srcpos_loc = splitsrc((*acq).srcpos,&nsrc_loc, nsrc);
}
MPI_Barrier(MPI_COMM_WORLD);
/*==================================================================================
Starting simulation (forward model)
==================================================================================*/
/* calculate wavelet for each source point */
(*acq).signals=NULL;
(*acq).signals=wavelet((*acq).srcpos_loc,nsrc_loc,ishot);
if (nsrc_loc){if(QUELLART==6){
/* time domain filtering of the source signal */
apply_tdfilt((*acq).signals,nsrc_loc,ns,ORDER_SPIKE,FC_SPIKE_2,FC_SPIKE_1);
}
}
/* time domain filtering*/
if ((TIME_FILT)&&(INV_STF==0)){
/* time domain filtering of the source signal */
apply_tdfilt((*acq).signals,nsrc_loc,ns,ORDER,FC,FC_START);
}
/* solve forward problem */
psv(wavePSV,wavePSV_PML,matPSV,fwiPSV,mpiPSV,seisPSV,seisPSVfwi,acq,hc,ishot,nshots,nsrc_loc,ns,ntr,Ws,Wr,hin,DTINV_help,2,req_send,req_rec);
/* ===============================================
Calculate objective function and data residuals
=============================================== */
if (ntr > 0){
calc_res_PSV(seisPSV,seisPSVfwi,(*acq).recswitch,(*acq).recpos,(*acq).recpos_loc,ntr_glob,ntr,nsrc_glob,(*acq).srcpos,ishot,ns,iter,1);
}
if((N_STREAMER>0)||(READREC==2)){
if (SEISMO) free_imatrix((*acq).recpos,1,3,1,ntr_glob);
if ((ntr>0) && (SEISMO)){
free_imatrix((*acq).recpos_loc,1,3,1,ntr);
(*acq).recpos_loc = NULL;
switch (SEISMO){
case 1 : /* particle velocities only */
free_matrix((*seisPSV).sectionvx,1,ntr,1,ns);
free_matrix((*seisPSV).sectionvy,1,ntr,1,ns);
(*seisPSV).sectionvx=NULL;
(*seisPSV).sectionvy=NULL;
break;
case 2 : /* pressure only */
free_matrix((*seisPSV).sectionp,1,ntr,1,ns);
break;
case 3 : /* curl and div only */
free_matrix((*seisPSV).sectioncurl,1,ntr,1,ns);
free_matrix((*seisPSV).sectiondiv,1,ntr,1,ns);
break;
case 4 : /* everything */
free_matrix((*seisPSV).sectionvx,1,ntr,1,ns);
free_matrix((*seisPSV).sectionvy,1,ntr,1,ns);
free_matrix((*seisPSV).sectionp,1,ntr,1,ns);
free_matrix((*seisPSV).sectioncurl,1,ntr,1,ns);
free_matrix((*seisPSV).sectiondiv,1,ntr,1,ns);
break;
}
}
free_matrix((*seisPSVfwi).sectionread,1,ntr_glob,1,ns);
free_ivector((*acq).recswitch,1,ntr);
if((QUELLTYPB==1)||(QUELLTYPB==3)||(QUELLTYPB==5)||(QUELLTYPB==7)){
free_matrix((*seisPSVfwi).sectionvxdata,1,ntr,1,ns);
free_matrix((*seisPSVfwi).sectionvxdiff,1,ntr,1,ns);
free_matrix((*seisPSVfwi).sectionvxdiffold,1,ntr,1,ns);
}
if((QUELLTYPB==1)||(QUELLTYPB==2)||(QUELLTYPB==6)||(QUELLTYPB==7)){
free_matrix((*seisPSVfwi).sectionvydata,1,ntr,1,ns);
free_matrix((*seisPSVfwi).sectionvydiff,1,ntr,1,ns);
free_matrix((*seisPSVfwi).sectionvydiffold,1,ntr,1,ns);
}
if(QUELLTYPB>=4){
free_matrix((*seisPSVfwi).sectionpdata,1,ntr,1,ns);
free_matrix((*seisPSVfwi).sectionpdiff,1,ntr,1,ns);
free_matrix((*seisPSVfwi).sectionpdiffold,1,ntr,1,ns);
}
ntr=0;
ntr_glob=0;
}
nsrc_loc=0;
} /* end of loop over shots */
/* calculate L2 norm of all CPUs*/
L2sum = 0.0;
L2_tmp = (*seisPSVfwi).L2;
MPI_Allreduce(&L2_tmp,&L2sum,1,MPI_FLOAT,MPI_SUM,MPI_COMM_WORLD);
/* calculate L2 norm of all CPUs*/
energy_all_shots = 0.0;
energy_tmp = (*seisPSVfwi).energy;
MPI_Allreduce(&energy_tmp,&energy_all_shots,1,MPI_FLOAT,MPI_SUM,MPI_COMM_WORLD);
/* if(MYID==0){
printf("L2sum: %e\n", L2sum);
printf("energy_sum: %e\n\n", energy_all_shots);
}*/
if(LNORM==2){
L2sum = L2sum/energy_all_shots;
}
else{L2sum=L2sum;}
return L2sum;
}