update_s_elastic_PML_SH.c 17.4 KB
Newer Older
1
/*-----------------------------------------------------------------------------------------
Florian Wittkamp's avatar
Florian Wittkamp committed
2
 * Copyright (C) 2016  For the list of authors, see file AUTHORS.
3
 *
Florian Wittkamp's avatar
Florian Wittkamp committed
4
 * This file is part of IFOS.
5
 *
Florian Wittkamp's avatar
Florian Wittkamp committed
6
 * IFOS is free software: you can redistribute it and/or modify
7
8
9
 * it under the terms of the GNU General Public License as published by
 * the Free Software Foundation, version 2.0 of the License only.
 *
Florian Wittkamp's avatar
Florian Wittkamp committed
10
 * IFOS is distributed in the hope that it will be useful,
11
12
13
14
15
 * but WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 * GNU General Public License for more details.
 *
 * You should have received a copy of the GNU General Public License
Florian Wittkamp's avatar
Florian Wittkamp committed
16
 * along with IFOS. See file COPYING and/or <http://www.gnu.org/licenses/gpl-2.0.html>.
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
 -----------------------------------------------------------------------------------------*/

/*------------------------------------------------------------------------
 *   updating stress components at gridpoints [nx1...nx2][ny1...ny2]
 *   by a staggered grid finite difference scheme of arbitrary (FDORDER) order accuracy in space
 *   and second order accuracy in time
 *
 *   SH-Version
 *
 *  ----------------------------------------------------------------------*/

#include "fd.h"

void update_s_elastic_PML_SH(int nx1, int nx2, int ny1, int ny2, float **  vz, float **   sxz, float **   syz, float ** uxz, float ** uyz, float *hc,  int infoout,float * K_x, float * a_x, float * b_x, float * K_x_half, float * a_x_half, float * b_x_half,
                             float * K_y, float * a_y, float * b_y, float * K_y_half, float * a_y_half, float * b_y_half,float ** psi_vzx, float ** psi_vzy,float ** uipjp,float ** u,float ** rho){
    
    int i,j, fdoh, h, h1;
    float fipjp,f = 0.0;
    float vzx, vzy;
    
    float  dhi;
    extern float DT, DH;
39
    extern int MYID, FDORDER, FW, L,PARAMETERIZATION;
40
41
42
    extern int FREE_SURF, BOUNDARY;
    extern int NPROCX, NPROCY, POS[3];
    extern FILE *FP;
Florian Wittkamp's avatar
Florian Wittkamp committed
43
    double time1 = 0.0, time2;
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
    
    
    dhi = DT/DH;
    fdoh = FDORDER/2;
    
    
    if (infoout && (MYID==0)){
        time1=MPI_Wtime();
        fprintf(FP,"\n **Message from update_s_elastic_SH (printed by PE %d):\n",MYID);
        fprintf(FP," Updating stress components ...");
    }
    
    
    
    switch (FDORDER){
            
        case 2:
            for (j=ny1;j<=ny2;j++){
                for (i=nx1;i<=nx2;i++){
                    vzx = (  hc[1]*(vz[j][i+1]-vz[j][i]))*dhi;
                    
                    vzy = (  hc[1]*(vz[j][i]-vz[j-1][i]))*dhi;
                    
                    
                    /* left boundary */
                    if((!BOUNDARY) && (POS[1]==0) && (i<=FW)){
                        psi_vzx[j][i] = b_x_half[i] * psi_vzx[j][i] + a_x_half[i] * vzx;
                        vzx = vzx / K_x_half[i] + psi_vzx[j][i];
                    }
                    
                    /* right boundary */
                    if((!BOUNDARY) && (POS[1]==NPROCX-1) && (i>=nx2-FW+1)){
                        h1 = (i-nx2+2*FW);
                        h = i;
                        psi_vzx[j][h1] = b_x_half[h1] * psi_vzx[j][h1] + a_x_half[h1] * vzx;
                        vzx = vzx / K_x_half[h1] + psi_vzx[j][h1];
                    }
                    
                    /* top boundary */
                    if((POS[2]==0) && (!(FREE_SURF)) && (j<=FW)){
                        psi_vzy[j][i] = b_y[j] * psi_vzy[j][i] + a_y[j] * vzy;
                        vzy = vzy / K_y[j] + psi_vzy[j][i];
                    }
                    
                    /* bottom boundary */
                    if((POS[2]==NPROCY-1) && (j>=ny2-FW+1)){
                        h1 = (j-ny2+2*FW);
                        h = j;
                        psi_vzy[h1][i] = b_y[h1] * psi_vzy[h1][i] + a_y[h1] * vzy;
                        vzy = vzy / K_y[h1] + psi_vzy[h1][i];
                    }
                    
                    fipjp=uipjp[j][i];
                    
                    /* lambda - mu relationship*/
99
100
                    if (PARAMETERIZATION==3) f = u[j][i];
                    if (PARAMETERIZATION==1) f = rho[j][i] * u[j][i] * u[j][i];
101
102
103
104
105
                    
                    /* updating components of the stress tensor, partially */
                    sxz[j][i]+=(fipjp*vzx);
                    syz[j][i]+=(f*vzy);
                    
106
107
                    uxz[j][i]=(fipjp*vzx)/DT;
                    uyz[j][i]=(f*vzy)/DT;
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
                    
                }
            }
            
            break;
            
        case 4:
            for (j=ny1;j<=ny2;j++){
                for (i=nx1;i<=nx2;i++){
                    vzx = (  hc[1]*(vz[j][i+1]-vz[j][i])
                           + hc[2]*(vz[j][i+2]-vz[j][i-1]))*dhi;
                    
                    vzy = (  hc[1]*(vz[j][i]-vz[j-1][i])
                           + hc[2]*(vz[j+1][i]-vz[j-2][i]))*dhi;
                    
                    
                    /* left boundary */
                    if((!BOUNDARY) && (POS[1]==0) && (i<=FW)){
                        psi_vzx[j][i] = b_x_half[i] * psi_vzx[j][i] + a_x_half[i] * vzx;
                        vzx = vzx / K_x_half[i] + psi_vzx[j][i];
                    }
                    
                    /* right boundary */
                    if((!BOUNDARY) && (POS[1]==NPROCX-1) && (i>=nx2-FW+1)){
                        h1 = (i-nx2+2*FW);
                        h = i;
                        psi_vzx[j][h1] = b_x_half[h1] * psi_vzx[j][h1] + a_x_half[h1] * vzx;
                        vzx = vzx / K_x_half[h1] + psi_vzx[j][h1];
                    }
                    
                    /* top boundary */
                    if((POS[2]==0) && (!(FREE_SURF)) && (j<=FW)){
                        psi_vzy[j][i] = b_y[j] * psi_vzy[j][i] + a_y[j] * vzy;
                        vzy = vzy / K_y[j] + psi_vzy[j][i];
                    }
                    
                    /* bottom boundary */
                    if((POS[2]==NPROCY-1) && (j>=ny2-FW+1)){
                        h1 = (j-ny2+2*FW);
                        h = j;
                        psi_vzy[h1][i] = b_y[h1] * psi_vzy[h1][i] + a_y[h1] * vzy;
                        vzy = vzy / K_y[h1] + psi_vzy[h1][i];
                    }
                    
                    fipjp=uipjp[j][i];
                    
                    /* lambda - mu relationship*/
155
156
                    if (PARAMETERIZATION==3) f = u[j][i];
                    if (PARAMETERIZATION==1) f = rho[j][i] * u[j][i] * u[j][i];
157
158
159
160
161
                    
                    /* updating components of the stress tensor, partially */
                    sxz[j][i]+=(fipjp*vzx);
                    syz[j][i]+=(f*vzy);
                    
162
163
                    uxz[j][i]=(fipjp*vzx)/DT;
                    uyz[j][i]=(f*vzy)/DT;
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
                    
                }
            }
            break;
            
        case 6:
            for (j=ny1;j<=ny2;j++){
                for (i=nx1;i<=nx2;i++){
                    vzx = (  hc[1]*(vz[j][i+1]-vz[j][i])
                           + hc[2]*(vz[j][i+2]-vz[j][i-1])
                           + hc[3]*(vz[j][i+3]-vz[j][i-2]))*dhi;
                    
                    vzy = (  hc[1]*(vz[j][i]-vz[j-1][i])
                           + hc[2]*(vz[j+1][i]-vz[j-2][i])
                           + hc[3]*(vz[j+2][i]-vz[j-3][i]))*dhi;
                    
                    
                    /* left boundary */
                    if((!BOUNDARY) && (POS[1]==0) && (i<=FW)){
                        psi_vzx[j][i] = b_x_half[i] * psi_vzx[j][i] + a_x_half[i] * vzx;
                        vzx = vzx / K_x_half[i] + psi_vzx[j][i];
                    }
                    
                    /* right boundary */
                    if((!BOUNDARY) && (POS[1]==NPROCX-1) && (i>=nx2-FW+1)){
                        h1 = (i-nx2+2*FW);
                        h = i;
                        psi_vzx[j][h1] = b_x_half[h1] * psi_vzx[j][h1] + a_x_half[h1] * vzx;
                        vzx = vzx / K_x_half[h1] + psi_vzx[j][h1];
                    }
                    
                    /* top boundary */
                    if((POS[2]==0) && (!(FREE_SURF)) && (j<=FW)){
                        psi_vzy[j][i] = b_y[j] * psi_vzy[j][i] + a_y[j] * vzy;
                        vzy = vzy / K_y[j] + psi_vzy[j][i];
                    }
                    
                    /* bottom boundary */
                    if((POS[2]==NPROCY-1) && (j>=ny2-FW+1)){
                        h1 = (j-ny2+2*FW);
                        h = j;
                        psi_vzy[h1][i] = b_y[h1] * psi_vzy[h1][i] + a_y[h1] * vzy;
                        vzy = vzy / K_y[h1] + psi_vzy[h1][i];
                    }
                    
                    fipjp=uipjp[j][i];
                    
                    /* lambda - mu relationship*/
212
213
                    if (PARAMETERIZATION==3) f = u[j][i];
                    if (PARAMETERIZATION==1) f = rho[j][i] * u[j][i] * u[j][i];
214
215
216
217
218
                    
                    /* updating components of the stress tensor, partially */
                    sxz[j][i]+=(fipjp*vzx);
                    syz[j][i]+=(f*vzy);
                    
219
220
                    uxz[j][i]=(fipjp*vzx)/DT;
                    uyz[j][i]=(f*vzy)/DT;
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
                    
                }
            }
            break;
            
        case 8:
            for (j=ny1;j<=ny2;j++){
                for (i=nx1;i<=nx2;i++){
                    vzx = (  hc[1]*(vz[j][i+1]-vz[j][i])
                           + hc[2]*(vz[j][i+2]-vz[j][i-1])
                           + hc[3]*(vz[j][i+3]-vz[j][i-2])
                           + hc[4]*(vz[j][i+4]-vz[j][i-3]))*dhi;
                    
                    vzy = (  hc[1]*(vz[j][i]-vz[j-1][i])
                           + hc[2]*(vz[j+1][i]-vz[j-2][i])
                           + hc[3]*(vz[j+2][i]-vz[j-3][i])
                           + hc[4]*(vz[j+3][i]-vz[j-4][i]))*dhi;
                    
                    
                    /* left boundary */
                    if((!BOUNDARY) && (POS[1]==0) && (i<=FW)){
                        psi_vzx[j][i] = b_x_half[i] * psi_vzx[j][i] + a_x_half[i] * vzx;
                        vzx = vzx / K_x_half[i] + psi_vzx[j][i];
                    }
                    
                    /* right boundary */
                    if((!BOUNDARY) && (POS[1]==NPROCX-1) && (i>=nx2-FW+1)){
                        h1 = (i-nx2+2*FW);
                        h = i;
                        psi_vzx[j][h1] = b_x_half[h1] * psi_vzx[j][h1] + a_x_half[h1] * vzx;
                        vzx = vzx / K_x_half[h1] + psi_vzx[j][h1];
                    }
                    
                    /* top boundary */
                    if((POS[2]==0) && (!(FREE_SURF)) && (j<=FW)){
                        psi_vzy[j][i] = b_y[j] * psi_vzy[j][i] + a_y[j] * vzy;
                        vzy = vzy / K_y[j] + psi_vzy[j][i];
                    }
                    
                    /* bottom boundary */
                    if((POS[2]==NPROCY-1) && (j>=ny2-FW+1)){
                        h1 = (j-ny2+2*FW);
                        h = j;
                        psi_vzy[h1][i] = b_y[h1] * psi_vzy[h1][i] + a_y[h1] * vzy;
                        vzy = vzy / K_y[h1] + psi_vzy[h1][i];
                    }
                    
                    fipjp=uipjp[j][i];
                    
                    /* lambda - mu relationship*/
271
272
                    if (PARAMETERIZATION==3) f = u[j][i];
                    if (PARAMETERIZATION==1) f = rho[j][i] * u[j][i] * u[j][i];
273
274
275
276
277
                    
                    /* updating components of the stress tensor, partially */
                    sxz[j][i]+=(fipjp*vzx);
                    syz[j][i]+=(f*vzy);
                    
278
279
                    uxz[j][i]=(fipjp*vzx)/DT;
                    uyz[j][i]=(f*vzy)/DT;
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
                    
                }
            }
            break;
            
        case 10:
            for (j=ny1;j<=ny2;j++){
                for (i=nx1;i<=nx2;i++){
                    vzx = (  hc[1]*(vz[j][i+1]-vz[j][i])
                           + hc[2]*(vz[j][i+2]-vz[j][i-1])
                           + hc[3]*(vz[j][i+3]-vz[j][i-2])
                           + hc[4]*(vz[j][i+4]-vz[j][i-3])
                           + hc[5]*(vz[j][i+5]-vz[j][i-4]))*dhi;
                    
                    vzy = (  hc[1]*(vz[j][i]-vz[j-1][i])
                           + hc[2]*(vz[j+1][i]-vz[j-2][i])
                           + hc[3]*(vz[j+2][i]-vz[j-3][i])
                           + hc[4]*(vz[j+3][i]-vz[j-4][i])
                           + hc[5]*(vz[j+4][i]-vz[j-5][i]))*dhi;
                    
                    
                    /* left boundary */
                    if((!BOUNDARY) && (POS[1]==0) && (i<=FW)){
                        psi_vzx[j][i] = b_x_half[i] * psi_vzx[j][i] + a_x_half[i] * vzx;
                        vzx = vzx / K_x_half[i] + psi_vzx[j][i];
                    }
                    
                    /* right boundary */
                    if((!BOUNDARY) && (POS[1]==NPROCX-1) && (i>=nx2-FW+1)){
                        h1 = (i-nx2+2*FW);
                        h = i;
                        psi_vzx[j][h1] = b_x_half[h1] * psi_vzx[j][h1] + a_x_half[h1] * vzx;
                        vzx = vzx / K_x_half[h1] + psi_vzx[j][h1];
                    }
                    
                    /* top boundary */
                    if((POS[2]==0) && (!(FREE_SURF)) && (j<=FW)){
                        psi_vzy[j][i] = b_y[j] * psi_vzy[j][i] + a_y[j] * vzy;
                        vzy = vzy / K_y[j] + psi_vzy[j][i];
                    }
                    
                    /* bottom boundary */
                    if((POS[2]==NPROCY-1) && (j>=ny2-FW+1)){
                        h1 = (j-ny2+2*FW);
                        h = j;
                        psi_vzy[h1][i] = b_y[h1] * psi_vzy[h1][i] + a_y[h1] * vzy;
                        vzy = vzy / K_y[h1] + psi_vzy[h1][i];
                    }
                    
                    fipjp=uipjp[j][i];
                    
                    /* lambda - mu relationship*/
332
333
                    if (PARAMETERIZATION==3) f = u[j][i];
                    if (PARAMETERIZATION==1) f = rho[j][i] * u[j][i] * u[j][i];
334
335
336
337
338
                    
                    /* updating components of the stress tensor, partially */
                    sxz[j][i]+=(fipjp*vzx);
                    syz[j][i]+=(f*vzy);
                    
339
340
                    uxz[j][i]=(fipjp*vzx)/DT;
                    uyz[j][i]=(f*vzy)/DT;
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
                    
                }
            }
            break;
            
        case 12:
            for (j=ny1;j<=ny2;j++){
                for (i=nx1;i<=nx2;i++){
                    vzx = (  hc[1]*(vz[j][i+1]-vz[j][i])
                           + hc[2]*(vz[j][i+2]-vz[j][i-1])
                           + hc[3]*(vz[j][i+3]-vz[j][i-2])
                           + hc[4]*(vz[j][i+4]-vz[j][i-3])
                           + hc[5]*(vz[j][i+5]-vz[j][i-4])
                           + hc[6]*(vz[j][i+6]-vz[j][i-5]))*dhi;
                    
                    vzy = (  hc[1]*(vz[j][i]-vz[j-1][i])
                           + hc[2]*(vz[j+1][i]-vz[j-2][i])
                           + hc[3]*(vz[j+2][i]-vz[j-3][i])
                           + hc[4]*(vz[j+3][i]-vz[j-4][i])
                           + hc[5]*(vz[j+4][i]-vz[j-5][i])
                           + hc[6]*(vz[j+5][i]-vz[j-6][i]))*dhi;
                    
                    
                    /* left boundary */
                    if((!BOUNDARY) && (POS[1]==0) && (i<=FW)){
                        psi_vzx[j][i] = b_x_half[i] * psi_vzx[j][i] + a_x_half[i] * vzx;
                        vzx = vzx / K_x_half[i] + psi_vzx[j][i];
                    }
                    
                    /* right boundary */
                    if((!BOUNDARY) && (POS[1]==NPROCX-1) && (i>=nx2-FW+1)){
                        h1 = (i-nx2+2*FW);
                        h = i;
                        psi_vzx[j][h1] = b_x_half[h1] * psi_vzx[j][h1] + a_x_half[h1] * vzx;
                        vzx = vzx / K_x_half[h1] + psi_vzx[j][h1];
                    }
                    
                    /* top boundary */
                    if((POS[2]==0) && (!(FREE_SURF)) && (j<=FW)){
                        psi_vzy[j][i] = b_y[j] * psi_vzy[j][i] + a_y[j] * vzy;
                        vzy = vzy / K_y[j] + psi_vzy[j][i];
                    }
                    
                    /* bottom boundary */
                    if((POS[2]==NPROCY-1) && (j>=ny2-FW+1)){
                        h1 = (j-ny2+2*FW);
                        h = j;
                        psi_vzy[h1][i] = b_y[h1] * psi_vzy[h1][i] + a_y[h1] * vzy;
                        vzy = vzy / K_y[h1] + psi_vzy[h1][i];
                    }
                    
                    fipjp=uipjp[j][i];
                    
                    /* lambda - mu relationship*/
395
396
                    if (PARAMETERIZATION==3) f = u[j][i];
                    if (PARAMETERIZATION==1) f = rho[j][i] * u[j][i] * u[j][i];
397
398
399
400
401
                    
                    /* updating components of the stress tensor, partially */
                    sxz[j][i]+=(fipjp*vzx);
                    syz[j][i]+=(f*vzy);
                    
402
403
                    uxz[j][i]=(fipjp*vzx)/DT;
                    uyz[j][i]=(f*vzy)/DT;
404
405
406
407
408
409
410
411
412
413
414
415
                    
                }
            }
            break;
    } /* end of switch(FDORDER) */
    
    
    if (infoout && (MYID==0)){
        time2=MPI_Wtime();
        fprintf(FP," finished (real time: %4.2f s).\n",time2-time1);
    }
}