PML_pro.c 14.8 KB
Newer Older
Tilman Steinweg's avatar
Tilman Steinweg committed
1
/*-----------------------------------------------------------------------------------------
Florian Wittkamp's avatar
Florian Wittkamp committed
2
 * Copyright (C) 2016  For the list of authors, see file AUTHORS.
Tilman Steinweg's avatar
Tilman Steinweg committed
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
218
219
220
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
271
272
273
274
275
276
277
278
279
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
332
333
334
335
336
337
338
339
340
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
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
 *
 * This file is part of DENISE.
 * 
 * DENISE is free software: you can redistribute it and/or modify
 * it under the terms of the GNU General Public License as published by
 * the Free Software Foundation, version 2.0 of the License only.
 * 
 * DENISE is distributed in the hope that it will be useful,
 * 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
 * along with DENISE. See file COPYING and/or <http://www.gnu.org/licenses/gpl-2.0.html>.
-----------------------------------------------------------------------------------------*/

/*
* Define damping profiles for CPML boundary condition
* This C-PML implementation is adapted from the 2nd order isotropic CPML code by Dimitri Komatitsch and based in part on formulas given in Roden and Gedney (2000). 
* Additionally the code is based on the following references:
* 
* @ARTICLE{KoMa07,
* author = {Dimitri Komatitsch and Roland Martin},
* title = {An unsplit convolutional {P}erfectly {M}atched {L}ayer improved
*          at grazing incidence for the seismic wave equation},
* journal = {Geophysics},
* year = {2007},
* volume = {72},
* number = {5},
* pages = {SM155-SM167},
* doi = {10.1190/1.2757586}}
*
* @ARTICLE{MaKoEz08,
* author = {Roland Martin and Dimitri Komatitsch and Abdela\^aziz Ezziani},
* title = {An unsplit convolutional perfectly matched layer improved at grazing
* incidence for seismic wave equation in poroelastic media},
* journal = {Geophysics},
* year = {2008},
* volume = {73},
* pages = {T51-T61},
* number = {4},
* doi = {10.1190/1.2939484}}
*
* @ARTICLE{MaKo09,
* author = {Roland Martin and Dimitri Komatitsch},
* title = {An unsplit convolutional perfectly matched layer technique improved
* at grazing incidence for the viscoelastic wave equation},
* journal = {Geophysical Journal International},
* year = {2009},
* volume = {179},
* pages = {333-344},
* number = {1},
* doi = {10.1111/j.1365-246X.2009.04278.x}}
*
* @ARTICLE{MaKoGe08,
* author = {Roland Martin and Dimitri Komatitsch and Stephen D. Gedney},
* title = {A variational formulation of a stabilized unsplit convolutional perfectly
* matched layer for the isotropic or anisotropic seismic wave equation},
* journal = {Computer Modeling in Engineering and Sciences},
* year = {2008},
* volume = {37},
* pages = {274-304},
* number = {3}}
*
* The original CPML technique for Maxwell's equations is described in:
*
* @ARTICLE{RoGe00,
* author = {J. A. Roden and S. D. Gedney},
* title = {Convolution {PML} ({CPML}): {A}n Efficient {FDTD} Implementation
*          of the {CFS}-{PML} for Arbitrary Media},
* journal = {Microwave and Optical Technology Letters},
* year = {2000},
* volume = {27},
* number = {5},
* pages = {334-339},
* doi = {10.1002/1098-2760(20001205)27:5<334::AID-MOP14>3.0.CO;2-A}}
*/

#include "fd.h"

void PML_pro(float * d_x, float * K_x, float * alpha_prime_x, float * a_x, float * b_x, 
            float * d_x_half, float * K_x_half, float * alpha_prime_x_half, float * a_x_half, float * b_x_half,
            float * d_y, float * K_y, float * alpha_prime_y, float * a_y, float * b_y, 
            float * d_y_half, float * K_y_half, float * alpha_prime_y_half, float * a_y_half, float * b_y_half)
{

	/* extern variables */

	extern float DH, DAMPING, DT, FPML;
	extern int FREE_SURF, NX, NY, BOUNDARY, NXG, NYG;
	extern int NPROCX, NPROCY, MYID, POS[3], FW;
	extern FILE *FP;
	
	/* local variables */
	int i, j, ifw, ii, jj, xb, yb, xe, ye, h;
	int USE_PML_XMIN, USE_PML_XMAX, USE_PML_YMIN, USE_PML_YMAX;
	float amp, a, *coeff, fc, cp;
      
      extern float npower, k_max_PML;
      /*const float npower = 2.0;*/ /*  power to compute d0 profile */
      /*const float k_max_PML = 1.0;*/ /* from Gedney page 8.11 */
      const float alpha_max_PML = 2.0 * PI * (FPML/2.0); /* from festa and Vilotte */

      float thickness_PML_x, thickness_PML_y, xoriginleft, xoriginright, yoriginbottom, yorigintop;
      float Rcoef , d0_x, d0_y, xval, yval, abscissa_in_PML, abscissa_normalized;

      char profile[STRING_SIZE];

      /* define profile of absorption in PML region */

      /* thickness of the PML layer in meters */
      thickness_PML_x = FW*DH;
      thickness_PML_y = FW*DH;

      /* reflection coefficient (INRIA report section 6.1) */
      Rcoef = 0.001;

      /* compute d0 from INRIA report section 6.1 */
      d0_x = - (npower + 1) * DAMPING * log(Rcoef) / (2.0 * thickness_PML_x);
      d0_y = - (npower + 1) * DAMPING * log(Rcoef) / (2.0 * thickness_PML_y);
	
      /* damping in the X direction */
      /* -------------------------- */

      /* origin of the PML layer (position of right edge minus thickness, in meters) */
      xoriginleft = thickness_PML_x;
      xoriginright = (NXG-1) * DH - thickness_PML_x;

      /* left boundary */

      for (i=1;i<=FW;i++){
         
          K_x[i] = 1.0;
          K_x_half[i] = 1.0;
          xval = DH * (i-1);

            /* define damping profile at the grid points */
            abscissa_in_PML = xoriginleft - xval;
      
            if(abscissa_in_PML >= 0.0){
               abscissa_normalized = abscissa_in_PML / thickness_PML_x;
               d_x[i] = d0_x * pow(abscissa_normalized,npower);

            /* this taken from Gedney page 8.2 */
               K_x[i] = 1.0 + (k_max_PML - 1.0) * pow(abscissa_normalized,npower);
               alpha_prime_x[i] = alpha_max_PML * (1.0 - abscissa_normalized);
            }

            /* define damping profile at half the grid points */
            abscissa_in_PML = xoriginleft - (xval + DH/2.0);

            if(abscissa_in_PML >= 0.0){
               abscissa_normalized = abscissa_in_PML / thickness_PML_x;
               d_x_half[i] = d0_x * pow(abscissa_normalized,npower);

               /* this taken from Gedney page 8.2 */
               K_x_half[i] = 1.0 + (k_max_PML - 1.0) * pow(abscissa_normalized,npower);
               alpha_prime_x_half[i] = alpha_max_PML * (1.0 - abscissa_normalized);
            }

       /* just in case, for -5 at the end */
       if(alpha_prime_x[i] < 0.0){ alpha_prime_x[i] = 0.0;}
       if(alpha_prime_x_half[i] < 0.0) {alpha_prime_x_half[i] = 0.0;}

       b_x[i] = exp(- (d_x[i] / K_x[i] + alpha_prime_x[i]) * DT);
       b_x_half[i] = exp(- (d_x_half[i] / K_x_half[i] + alpha_prime_x_half[i]) * DT);

       /* avoid division by zero outside the PML */
       if(abs(d_x[i]) > 1.0e-6){ a_x[i] = d_x[i] * (b_x[i] - 1.0) / (K_x[i] * (d_x[i] + K_x[i] * alpha_prime_x[i]));}
       if(abs(d_x_half[i]) > 1.0e-6){ a_x_half[i] = d_x_half[i] * (b_x_half[i] - 1.0) / (K_x_half[i] * (d_x_half[i] + K_x_half[i] * alpha_prime_x_half[i]));}

      } /* end of left boundary */ 

      /* right boundary */
           
      for (i=NXG-FW+1;i<=NXG;i++){
             
             h=i-NXG+2*FW;
             
             K_x[h] = 1.0;
	     K_x_half[h] = 1.0;
	     xval = DH * (i-1);                                 
              
            /* define damping profile at the grid points */
            abscissa_in_PML = xval - xoriginright;
      
            if(abscissa_in_PML >= 0.0){
               abscissa_normalized = abscissa_in_PML / thickness_PML_x;
               d_x[h] = d0_x * pow(abscissa_normalized,npower);

            /* this taken from Gedney page 8.2 */
               K_x[h] = 1.0 + (k_max_PML - 1.0) * pow(abscissa_normalized,npower);
               alpha_prime_x[h] = alpha_max_PML * (1.0 - abscissa_normalized);
            }

            /* define damping profile at half the grid points */
            abscissa_in_PML = xval + DH/2.0 - xoriginright;

            if(abscissa_in_PML >= 0.0){
               abscissa_normalized = abscissa_in_PML / thickness_PML_x;
               d_x_half[h] = d0_x * pow(abscissa_normalized,npower);

               /* this taken from Gedney page 8.2 */
               K_x_half[h] = 1.0 + (k_max_PML - 1.0) * pow(abscissa_normalized,npower);
               alpha_prime_x_half[h] = alpha_max_PML * (1.0 - abscissa_normalized);
            }

            /* just in case, for -5 at the end */
               if(alpha_prime_x[h] < 0.0){ alpha_prime_x[h] = 0.0;}
               if(alpha_prime_x_half[h] < 0.0) {alpha_prime_x_half[h] = 0.0;}

               b_x[h] = exp(- (d_x[h] / K_x[h] + alpha_prime_x[h]) * DT);
               b_x_half[h] = exp(- (d_x_half[h] / K_x_half[h] + alpha_prime_x_half[h]) * DT);

            /* avoid division by zero outside the PML */
               if(abs(d_x[h]) > 1.0e-6){ a_x[h] = d_x[h] * (b_x[h] - 1.0) / (K_x[h] * (d_x[h] + K_x[h] * alpha_prime_x[h]));}
               if(abs(d_x_half[h]) > 1.0e-6){ a_x_half[h] = d_x_half[h] * (b_x_half[h] - 1.0) / (K_x_half[h] * (d_x_half[h] + K_x_half[h] * alpha_prime_x_half[h]));}

       } /* end of right boundary */   



	
      /* damping in the Y direction */
      /* -------------------------- */

      /* origin of the PML layer (position of right edge minus thickness, in meters) */
      yoriginbottom = thickness_PML_y;
      yorigintop = (NYG-1) * DH - thickness_PML_y;

      for (i=1;i<=FW;i++){
          
          K_y[i] = 1.0;
          K_y_half[i] = 1.0;
          yval = DH * (i-1);

          /* left boundary */

            /* define damping profile at the grid points */
            abscissa_in_PML = yoriginbottom - yval;
      
            if(abscissa_in_PML >= 0.0){
               abscissa_normalized = abscissa_in_PML / thickness_PML_y;
               d_y[i] = d0_y * pow(abscissa_normalized,npower);

            /* this taken from Gedney page 8.2 */
               K_y[i] = 1.0 + (k_max_PML - 1.0) * pow(abscissa_normalized,npower);
               alpha_prime_y[i] = alpha_max_PML * (1.0 - abscissa_normalized);
            }

            /* define damping profile at half the grid points */
            abscissa_in_PML = yoriginbottom - (yval + DH/2.0);

            if(abscissa_in_PML >= 0.0){
               abscissa_normalized = abscissa_in_PML / thickness_PML_y;
               d_y_half[i] = d0_y * pow(abscissa_normalized,npower);

               /* this taken from Gedney page 8.2 */
               K_y_half[i] = 1.0 + (k_max_PML - 1.0) * pow(abscissa_normalized,npower);
               alpha_prime_y_half[i] = alpha_max_PML * (1.0 - abscissa_normalized);
            }
          b_y[i] = exp(- (d_y[i] / K_y[i] + alpha_prime_y[i]) * DT);
          b_y_half[i] = exp(- (d_y_half[i] / K_y_half[i] + alpha_prime_y_half[i]) * DT);

          /* avoid division by zero outside the PML */
          if(abs(d_y[i]) > 1.0e-6){ a_y[i] = d_y[i] * (b_y[i] - 1.0) / (K_y[i] * (d_y[i] + K_y[i] * alpha_prime_y[i]));}
          if(abs(d_y_half[i]) > 1.0e-6){ a_y_half[i] = d_y_half[i] * (b_y_half[i] - 1.0) / (K_y_half[i] * (d_y_half[i] + K_y_half[i] * alpha_prime_y_half[i]));}

      } /* end of left boundary */ 

      /* top boundary */
      for (i=NYG-FW+1;i<=NYG;i++){  
           
            h=i-NYG+2*FW;
            
            K_y[h] = 1.0;
            K_y_half[h] = 1.0;
            yval = DH * (i-1);

            /* define damping profile at the grid points */
            abscissa_in_PML = yval - yorigintop;
      
            if(abscissa_in_PML >= 0.0){
               abscissa_normalized = abscissa_in_PML / thickness_PML_y;
               d_y[h] = d0_y * pow(abscissa_normalized,npower);

            /* this taken from Gedney page 8.2 */
               K_y[h] = 1.0 + (k_max_PML - 1.0) * pow(abscissa_normalized,npower);
               alpha_prime_y[h] = alpha_max_PML * (1.0 - abscissa_normalized);
            }

            /* define damping profile at half the grid points */
            abscissa_in_PML = yval + DH/2.0 - yorigintop;

            if(abscissa_in_PML >= 0.0){
               abscissa_normalized = abscissa_in_PML / thickness_PML_y;
               d_y_half[h] = d0_x * pow(abscissa_normalized,npower);

               /* this taken from Gedney page 8.2 */
               K_y_half[h] = 1.0 + (k_max_PML - 1.0) * pow(abscissa_normalized,npower);
               alpha_prime_y_half[h] = alpha_max_PML * (1.0 - abscissa_normalized);
            }

          b_y[h] = exp(- (d_y[h] / K_y[h] + alpha_prime_y[h]) * DT);
          b_y_half[h] = exp(- (d_y_half[h] / K_y_half[h] + alpha_prime_y_half[h]) * DT);

          /* avoid division by zero outside the PML */
          if(abs(d_y[h]) > 1.0e-6){ a_y[h] = d_y[h] * (b_y[h] - 1.0) / (K_y[h] * (d_y[h] + K_y[h] * alpha_prime_y[h]));}
          if(abs(d_y_half[h]) > 1.0e-6){ a_y_half[h] = d_y_half[h] * (b_y_half[h] - 1.0) / (K_y_half[h] * (d_y_half[h] + K_y_half[h] * alpha_prime_y_half[h]));}
      
       } /* end of top boundary */   



      
      
      /*if(MYID==0){
      sprintf(profile,"Kx.dat");
      FP=fopen(profile,"w");
      
    
      for (i=1;i<=2*FW+1;i=i++){
       fprintf(FP,"%e \n",K_y[i]);
      }
      fclose(FP);
      
      sprintf(profile,"Kx_half.dat");                                  
      FP=fopen(profile,"w");  
                    

      for (i=1;i<=2*FW+1;i++){
	fprintf(FP,"%e \n",K_y_half[i]);                          
      }
      fclose(FP); 

      sprintf(profile,"ax.dat");                                  
      FP=fopen(profile,"w");  
                    
     
      for (i=1;i<=2*FW+1;i++){
	fprintf(FP,"%e \n",a_y[i]);                          
      }
      fclose(FP); 
      
      sprintf(profile,"bx.dat");                                  
      FP=fopen(profile,"w");  
                    
   
      for (i=1;i<=2*FW+1;i++){
	fprintf(FP,"%e \n",b_y[i]);                          
      }
      fclose(FP);
              
      }*/
      

	/* compute coefficients for left and right grid boundaries (x-direction) */
	/*if ((!BOUNDARY) && (POS[1]==0))
	{
		yb=1; ye=NY; 
		for (i=1;i<=ifw;i++)
		{
			if ((POS[2]==0) && (!(FREE_SURF))) yb=i;
			if (POS[2]==NPROCY-1) ye=NY-i+1;
			for (j=yb;j<=ye;j++)
				absorb_coeff[j][i]=coeff[i];
		}
	}
			
	if ((!BOUNDARY) && (POS[1]==NPROCX-1))
	{
		yb=1; ye=NY;
		for (i=1;i<=ifw;i++){
			ii=NX-i+1;
			if ((POS[2]==0) && (!(FREE_SURF))) yb=i;
			if (POS[2]==NPROCY-1) ye=NY-i+1;
			for (j=yb;j<=ye;j++)
				absorb_coeff[j][ii]=coeff[i];
		}
	}*/
	

	/* compute coefficients for top and bottom grid boundaries (y-direction) */

	/*if ((POS[2]==0) && (!(FREE_SURF)))
	{
		xb=1; xe=NX;
		for (j=1;j<=ifw;j++)
		{
			if ((!BOUNDARY) && (POS[1]==0)) xb=j;
			if ((!BOUNDARY) && (POS[1]==NPROCX-1)) xe=NX-j+1;
			for (i=xb;i<=xe;i++)
				absorb_coeff[j][i]=coeff[j];
		}
	}

	if (POS[2]==NPROCY-1)
	{
		xb=1; xe=NX;
		for (j=1;j<=ifw;j++)
		{
			jj=NY-j+1;
			if ((!BOUNDARY) && (POS[1]==0)) xb=j;
			if ((!BOUNDARY) && (POS[1]==NPROCX-1)) xe=NX-j+1;
			for (i=xb;i<=xe;i++)
				absorb_coeff[jj][i]=coeff[j];
		}
	}*/


/*	sprintf(modfile,"absorb_coeff.bin");

	writemod(modfile,absorb_coeff,3); 

	MPI_Barrier(MPI_COMM_WORLD);

	if (MYID==0) mergemod(modfile,3); 
*/

}