modelupdate.c 4.44 KB
Newer Older
Simone Butzer's avatar
Simone Butzer committed
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
/*------------------------------------------------------------------------
 * Copyright (C) 2015 For the list of authors, see file AUTHORS.
 *
 * This file is part of IFOS3D.
 * 
 * IFOS3D 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.
 * 
 * IFOS3D 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 IFOS3D. See file COPYING and/or 
 * <http://www.gnu.org/licenses/gpl-2.0.html>.
--------------------------------------------------------------------------*/

/*--------------------------------------------------------------------------
 * update model for next iteration
 * S. Butzer 2013
 * ------------------------------------------------------------------------*/

#include "fd.h"
void modelupdate(int nx, int ny, int nz, float ***gradvp, float ***gradvs, float ***gradrho, float  ***  rho, float ***  pi, float ***  u, float **bfgsmod1, float step, float *beta, int it_group){

	extern FILE *FP;
	extern int NX, NY, NZ, LBFGS;

	int j,i,k,l,l1;
	float vpnew,vsnew,rhonew;
	float max[3],buf[3],max1[3],dummy[3];
	float vp,vs;
	int bfgsnum=5,w;
	float scale1=1.0, scale2=1.0;
	float vp0=6200.0, vs0=3600.0, rho0=2800.0;
	
	buf[0]=0.0;buf[1]=0.0;buf[2]=0.0;
	max[0]=0.0;max[1]=0.0;max[2]=0.0;
	dummy[0]=0.0; dummy[1]=0.0; dummy[2]=0.0;
	max1[0]=0.0;max1[1]=0.0;max1[2]=0.0;
	
	/*if(iteration<20)scale1=0.7;
	if(iteration>40)scale2=0.7;
	if(iteration>60)scale2=0.4;*/
	
	if(LBFGS) {w=it_group%bfgsnum;
		  if(w==0) w=bfgsnum;}

	 fprintf(FP,"Modelupdatefunction \n");
	
	/*find and exchange model maxima*/
	for (j=1;j<=ny;j++){
		for (i=1;i<=nx;i++){
			for (k=1;k<=nz;k++){
			vp=0.0;
			vs=0.0;
			vp=sqrt(pi[j][i][k]/rho[j][i][k]);
			vs=sqrt(u[j][i][k]/rho[j][i][k]);
			
			if (vp>max[0]) max[0]=vp;
			if (vs>max[1]) max[1]=vs;
			if (rho[j][i][k]>max[2]) max[2]=rho[j][i][k];
			}
		}
	}
	
	MPI_Barrier(MPI_COMM_WORLD);
	MPI_Allreduce(&max,&buf,3,MPI_FLOAT,MPI_MAX,MPI_COMM_WORLD);
	MPI_Barrier(MPI_COMM_WORLD);
	max[0]=buf[0];
	max[1]=buf[1];
	max[2]=buf[2];
	fprintf(FP,"rhomax=%4.2f, vpmax=%4.2f, vsmax=%4.2f \n", max[2], max[0], max[1]);
	
	/*if(iteration==1){vs0=vs0*vs0;
	vp0=vp0*vp0;}
	 else {vs0=1.0;
	 vp0=1.0;}*/

        if(!LBFGS){
		  /*Normalise gradients to maximum*/
		  for (j=1;j<=ny;j++){
			  for (i=1;i<=nx;i++){
				  for (k=1;k<=nz;k++){
				    
				  if (fabs(gradvp[j][i][k])>max1[0]) max1[0]=fabs(gradvp[j][i][k]);
				  if (fabs(gradvs[j][i][k])>max1[1]) max1[1]=fabs(gradvs[j][i][k]);
				  if (fabs(gradrho[j][i][k])>max1[2]) max1[2]=fabs(gradrho[j][i][k]);
				    
				  }
			  }
		  }
		  MPI_Barrier(MPI_COMM_WORLD);
		  MPI_Allreduce(&max1,&dummy,3,MPI_FLOAT,MPI_MAX,MPI_COMM_WORLD);
		  MPI_Barrier(MPI_COMM_WORLD);
		  
		  max1[0]=dummy[0];
		  max1[1]=dummy[1];
		  max1[2]=dummy[2];
		  if(max1[2]==0.0) max1[2]=1.0;
		  fprintf(FP,"max1[0]=%e,max1[1]=%e  \n",max1[0],max1[1]);
		  
		  if (max1[0]>0.0 && max1[1]>0.0){
			  for (j=1;j<=ny;j++){
				  for (i=1;i<=nx;i++){
					  for (k=1;k<=nz;k++){
						  gradvp[j][i][k]= gradvp[j][i][k]/max1[0];
						  gradvs[j][i][k]= gradvs[j][i][k]/max1[1];
						  gradrho[j][i][k]=gradrho[j][i][k]/max1[2];
					  }
				  }
			  }
		  }
		  else fprintf(FP,"max1*max2*max3 =0 \n");
	}
	/*if(it_group==1){scale1=4.0; scale2=4.0;}*/
	l=0;
	l1=NX*NY*NZ;
	for (j=1;j<=ny;j++){
		for (i=1;i<=nx;i++){
			for (k=1;k<=nz;k++){

			
			/*update model*/
			vpnew=0.0;
			/*vpnew=sqrt(pi[j][i][k]/rho[j][i][k])+max[0]*step*gradconvp[j][i][k];*/
			vpnew=sqrt(pi[j][i][k]/rho[j][i][k])+step*scale1*gradvp[j][i][k]*vp0;
			
			vsnew=0.0;
			/*vsnew=sqrt(u[j][i][k]/rho[j][i][k])+max[1]*step*gradconvs[j][i][k];*/
			vsnew=sqrt(u[j][i][k]/rho[j][i][k])+step*scale2*gradvs[j][i][k]*vs0;
		      
			rhonew=0.0;
			rhonew=rho[j][i][k]+0.0*step*gradrho[j][i][k]*rho0;
			
			    if(LBFGS){
			/*save normalised model differences for LBFGS*/
			 l++; l1++;
			bfgsmod1[w][l]=(vpnew-sqrt(pi[j][i][k]/rho[j][i][k]))/vp0;
			bfgsmod1[w][l1]=(vsnew-sqrt(u[j][i][k]/rho[j][i][k]))/vs0;
			/*bfgsmod3[w][l]=(rhonew-rho[j][i][k])/rho0;*/
			}
			  
			u[j][i][k]=rhonew*vsnew*vsnew;
			pi[j][i][k]=rhonew*vpnew*vpnew;
			rho[j][i][k]=rhonew;
			
			}
		}
	}
		
}