void write_grid_file(void) { /* The following subroutine saves the depth, reduced gravity of */ /* each interface, the potential density of each layer, and the */ /* Coriolis parameter. A variety of metric terms are written to a */ /* separate file. */ char filename[40]; /* The file name of the save file. */ char filepath[120]; /* The path (dir/file) to the file. */ int i, j, cdfid, timeid; size_t err = 1; struct varcdfinfo varinfo[10]; /* vardesc is a structure defined in HIM_io.h. The elements of */ /* this structure, in order, are: (1) the variable name for the NetCDF*/ /* file; (2) the variable's long name; (3) a character indicating the */ /* horizontal grid, which may be '1' (column), 'h', 'q', 'u', or 'v', */ /* for the corresponding C-grid variable; (4) a character indicating */ /* the vertical grid, which may be 'L' (layer), 'i' (interface), */ /* '2' (mixed-layers), or '1' (no vertical location); (5) a character */ /* indicating the time levels of the field, which may be 's' (snap- */ /* shot), 'a' (average between snapshots), 'm' (monthly average), or */ /* '1' (no time variation); (6) the variable's units; and (7) a */ /* character indicating the size in memory to write, which may be */ /* 'd' (8-byte) or 'f' (4-byte). */ vardesc vars[4] = { {"D","Basin Depth",'h','1','1',"meter", 'd'}, {"g","Reduced gravity",'1','L','1',"meter second-2", 'd'}, {"R","Target Potential Density",'1','L','1',"kilogram meter-3", 'd'}, {"f","Coriolis Parameter",'q','1','1',"second-1", 'd'} }; sprintf(filename,"D.%d.%d.%d",NXTOT,NYTOT,NZ); strcpy(filepath, directory); strcat(filepath, filename); create_file(filepath, vars, 4, &cdfid, &timeid, varinfo); err *= write_field(cdfid, vars[0], varinfo[0], 0, D[0]); err *= write_field(cdfid, vars[1], varinfo[1], 0, g); err *= write_field(cdfid, vars[2], varinfo[2], 0, Rlay); err *= write_field(cdfid, vars[3], varinfo[3], 0, f[0]); close_file(&cdfid); if (err == 0) printf("Problems saving general parameters.\n"); { /* ### REVISIT THIS PART OF THE ROUTINE AFTER METRICS... ### */ double out[NYMEM][NXMEM]; /* An array for output. */ extern double hmask[NYMEM][NXMEM]; vardesc vars2[10]={ {"geolatb","latitude at q points",'q','1','1',"degree",'d'}, {"geolonb","longitude at q points",'q','1','1',"degree",'d'}, {"wet", "land or ocean?", 'h','1','1',"none",'d'}, {"geolat", "latitude at h points", 'h','1','1',"degree",'d'}, {"geolon","longitude at h points",'h','1','1',"degree",'d'}, {"dxh","Zonal grid spacing at h points",'h','1','1',"m",'d'}, {"dxq","Zonal grid spacing at q points",'q','1','1',"m",'d'}, {"dyh","Meridional grid spacing at h points",'h','1','1',"m",'d'}, {"dyq","Meridional grid spacing at q points",'q','1','1',"m",'d'}, {"Ah","Area of h cells",'h','1','1',"m2",'d'}, }; sprintf(filename,"grid.%d.%d",NXTOT,NYTOT); strcpy(filepath, directory); strcat(filepath, filename); create_file(filepath, vars2, 10, &cdfid, &timeid, varinfo); err *= write_field(cdfid, vars2[0], varinfo[0], 0, geolatq[0]); err *= write_field(cdfid, vars2[1], varinfo[1], 0, geolonq[0]); err *= write_field(cdfid, vars2[2], varinfo[2], 0, hmask[0]); err *= write_field(cdfid, vars2[3], varinfo[3], 0, geolath[0]); err *= write_field(cdfid, vars2[4], varinfo[4], 0, geolonh[0]); for (j=0;j<NYMEM;j++) for (i=0;i<NXMEM;i++) out[j][i]=DXh(j,i); err *= write_field(cdfid, vars2[5], varinfo[5], 0, out[0]); for (j=0;j<NYMEM;j++) for (i=0;i<NXMEM;i++) out[j][i]=DXq(j,i); err *= write_field(cdfid, vars2[6], varinfo[6], 0, out[0]); for (j=0;j<NYMEM;j++) for (i=0;i<NXMEM;i++) out[j][i]=DYh(j,i); err *= write_field(cdfid, vars2[7], varinfo[7], 0, out[0]); for (j=0;j<NYMEM;j++) for (i=0;i<NXMEM;i++) out[j][i]=DYq(j,i); err *= write_field(cdfid, vars2[8], varinfo[8], 0, out[0]); for (j=0;j<NYMEM;j++) for (i=0;i<NXMEM;i++) out[j][i]=DXDYh(j,i); err *= write_field(cdfid, vars2[9], varinfo[9], 0, out[0]); close_file(&cdfid); if (err == 0) printf("Problems saving latitude and longitude and wet.\n"); } }
void tracer(int itts) { /* This subroutine time steps the tracer concentration. */ /* A positive definite scheme is used. */ double minslope; /* The maximum concentration slope per */ /* grid point consistent with mono- */ /* tonicity, in conc. (nondim.). */ double ***hvol; /* The cell volume of an h-element */ double slope[NXMEM+NYMEM][NTR]; /* The concentration slope per grid */ /* point in units of concentration (nondim.). */ double fluxtr[NXMEM+NYMEM][NTR];/* The flux of tracer across a */ /* boundary, in m3 * conc. (nondim.). */ double ***uhr; /* The remaining zonal and meridional */ double ***vhr; /* thickness fluxes, in m3.*/ double uhh[NXMEM]; /* uhh and vhh are the reduced fluxes */ double vhh[NYMEM]; /* during the current iteration, in m3.d */ double hup, hlos; /* hup is the upwind volume, hlos is the */ /* part of that volume that might be lost */ /* due to advection out the other side of */ /* the grid box, both in m3. */ double ts2; double landvolfill; /* An arbitrary? nonzero cell volume, m3. */ double ***ear; double ***ebr; double ***wdh; double bet[NXMEM]; /* bet and gam are variables used by the */ double gam[NZ][NXMEM]; /* tridiagonal solver. */ double hnew0[NXMEM]; /* The original topmost layer thickness, */ #if defined AGE2 || defined AGE3 // extern double hnew[NZ][NXMEM][NYMEM]; extern double ***hnew; #else double ***hnew; #endif double hlst1, Ihnew; double hlst[NYMEM]; // double MLMIN = EPSILON; /* min depth for ML */ double MLMIN = 4.25; double BLMIN = 0.20; #ifdef ENTRAIN double nts = dt/DT; /* number of timesteps (#day*86400/3600seconds) */ #endif int i, j, k, m, ii, pstage; int itt; double fract1; double fract2; # ifdef WRTTS double wrts; # endif hvol = alloc3d(NZ,NXMEM,NYMEM); if(hvol == NULL) { fprintf(stderr,"not enough memory for hvol!\n"); } uhr = alloc3d(NZ,NXMEM,NYMEM); if(uhr == NULL) { fprintf(stderr,"not enough memory for uhr!\n"); } vhr = alloc3d(NZ,NXMEM,NYMEM); if(vhr == NULL) { fprintf(stderr,"not enough memory for vhr!\n"); } ear = alloc3d(NZ,NXMEM,NYMEM); if(ear == NULL) { fprintf(stderr,"not enough memory for ear!\n"); } ebr = alloc3d(NZ,NXMEM,NYMEM); if(ebr == NULL) { fprintf(stderr,"not enough memory for ebr!\n"); } wdh = alloc3d(NZ,NXMEM,NYMEM); if(wdh == NULL) { fprintf(stderr,"not enough memory for wdh!\n"); } #if !defined AGE2 && !defined AGE3 hnew = alloc3d(NZ,NXMEM,NYMEM); if(hnew == NULL) { fprintf(stderr,"not enough memory for hnew!\n"); } #endif landvolfill = EPSILON*1000000.0; /* This is arbitrary. */ /* zonal re-entrance */ #pragma omp parallel { #pragma omp for private(j,k) for (j=0;j<=NYMEM-1;j++) { for (k=0;k<NZ;k++) { uhtm[k][nx+1][j] = uhtm[k][2][j]; uhtm[k][nx+2][j] = uhtm[k][3][j]; uhtm[k][0][j] = uhtm[k][nx-1][j]; uhtm[k][1][j] = uhtm[k][nx][j]; vhtm[k][nx+1][j] = vhtm[k][2][j]; vhtm[k][nx+2][j] = vhtm[k][3][j]; vhtm[k][0][j] = vhtm[k][nx-1][j]; vhtm[k][1][j] = vhtm[k][nx][j]; } for (k=0;k<NZ+1;k++) { wd[k][nx+1][j] = wd[k][2][j]; wd[k][nx+2][j] = wd[k][3][j]; wd[k][0][j] = wd[k][nx-1][j]; wd[k][1][j] = wd[k][nx][j]; } } /* meridional re-entrance */ #pragma omp for private(i,k,ii) for (i=2;i<=nx;i++) { ii = 363 - i; for (k=0;k<NZ;k++) { uhtm[k][ii][ny+1] = (-1)*uhtm[k][i][ny]; uhtm[k][ii][ny+2] = (-1)*uhtm[k][i][ny-1]; vhtm[k][ii][ny+1] = (-1)*vhtm[k][i][ny]; vhtm[k][ii][ny+2] = (-1)*vhtm[k][i][ny-1]; } for (k=0;k<NZ+1;k++) { wd[k][ii][ny+1] = wd[k][i][ny]; wd[k][ii][ny+2] = wd[k][i][ny-1]; } } #pragma omp for private(i,j,k) for (k=0;k<NZ;k++) { /* Put the thickness fluxes into uhr and vhr. */ for (j=0;j<=ny+2;j++) { for (i=0;i<=nx+2;i++) { uhr[k][i][j] = uhtm[k][i][j]*dt; vhr[k][i][j] = vhtm[k][i][j]*dt; if (h[k][i][j] < EPSILON) { h[k][i][j] = 1.0*EPSILON; } /* This line calculates the cell volume */ hvol[k][i][j] = DXDYh(i,j)*h[k][i][j]; hnew[k][i][j] = h[k][i][j]; } } } /* calculate the diapycnal velocities at the interfaces */ /* if we read in the ea, eb and eaml variables */ /* Otherwise we read in wd directly */ #ifdef ENTRAIN #pragma omp for private(i,j) for (i=X0;i<=nx+1;i++) for (j=Y0;j<=ny;j++) wd[0][i][j] = nts*eaml[i][j]; #pragma omp for private(i,j,k) for (k=1;k<NZ;k++) { for (i=X0;i<=nx+1;i++) for (j=Y0;j<=ny;j++) wd[k][i][j] = nts*(ea[k][i][j] - eb[k-1][i][j]); } #endif } // omp #define STANDARD_ADVECTION //#undef STANDARD_ADVECTION #ifdef STANDARD_ADVECTION /* pstage=1; print_tr(pstage); */ /* beginning of itt loop */ for (itt = 0; itt < NUM_ADV_ITER; itt++) { /* big loop over k */ //ompfail #pragma omp parallel { #pragma omp for private(i,j,k,m,minslope,slope,uhh,vhh,fluxtr,hup,hlos,ts2,hlst,hlst1,Ihnew) for (k=0;k<NZ;k++) { /* To insure positive definiteness of the thickness at each */ /* iteration, the mass fluxes out of each layer are checked each */ /* time. This may not be very efficient, but it should be reliable. */ /* ============================================================ */ /* first advect zonally */ /* ============================================================ */ #ifndef ADV1D for (j=Y1;j<=ny;j++) { /* Calculate the i-direction profiles (slopes) of each tracer that */ /* is being advected. */ //#pragma omp for private(i,m,minslope) for (i=X0;i<=nx+1;i++) { for (m=0;m<NTR;m++) { minslope = 4.0*((fabs(tr[m][k][i+1][j]-tr[m][k][i][j]) < fabs(tr[m][k][i][j]-tr[m][k][i-1][j])) ? (tr[m][k][i+1][j]-tr[m][k][i][j]) : (tr[m][k][i][j]-tr[m][k][i-1][j])); slope[i][m] = umask[i][j]*umask[i-1][j] * (((tr[m][k][i+1][j]-tr[m][k][i][j]) * (tr[m][k][i][j]-tr[m][k][i-1][j]) < 0.0) ? 0.0 : ((fabs(tr[m][k][i+1][j]-tr[m][k][i-1][j])<fabs(minslope)) ? 0.5*(tr[m][k][i+1][j]-tr[m][k][i-1][j]) : 0.5*minslope)); } } //#pragma omp barrier /* Calculate the i-direction fluxes of each tracer, using as much */ /* the minimum of the remaining mass flux (uhr) and the half the mass */ /* in the cell plus whatever part of its half of the mass flux that */ /* the flux through the other side does not require. */ //#pragma omp for private(i,m,hup,hlos,ts2) for (i=X0;i<=nx;i++) { if (uhr[k][i][j] == 0.0) { uhh[i] = 0.0; for (m=0;m<NTR;m++) fluxtr[i][m] = 0.0; } else if (uhr[k][i][j] < 0.0) { if (k==0 || k==1 ) { hup = (hvol[k][i+1][j]-DXDYh(i+1,j)*MLMIN); } else { hup = (hvol[k][i+1][j]-DXDYh(i+1,j)*EPSILON); } hlos = D_MAX(0.0,uhr[k][i+1][j]); if (((hup + uhr[k][i][j] - hlos) < 0.0) && ((0.5*hup + uhr[k][i][j]) < 0.0)) { uhh[i] = D_MIN(-0.5*hup,-hup+hlos); } else uhh[i] = uhr[k][i][j]; ts2 = 0.5*(1.0 + uhh[i]/hvol[k][i+1][j]); for (m=0;m<NTR;m++) { fluxtr[i][m] = uhh[i]*(tr[m][k][i+1][j] - slope[i+1][m]*ts2); } } else { if (k==0 || k==1 ) { hup = (hvol[k][i][j]-DXDYh(i,j)*MLMIN); } else { hup = (hvol[k][i][j]-DXDYh(i,j)*EPSILON); } hlos = D_MAX(0.0,-uhr[k][i-1][j]); if (((hup - uhr[k][i][j] - hlos) < 0.0) && ((0.5*hup - uhr[k][i][j]) < 0.0)) { uhh[i] = D_MAX(0.5*hup,hup-hlos); } else uhh[i] = uhr[k][i][j]; ts2 = 0.5*(1.0 - uhh[i]/hvol[k][i][j]); for (m=0;m<NTR;m++) { fluxtr[i][m] = uhh[i]*(tr[m][k][i][j] + slope[i][m]*ts2); } } } //#pragma omp barrier /* Calculate new tracer concentration in each cell after accounting */ /* for the i-direction fluxes. */ uhr[k][X0][j] -= uhh[X0]; // #pragma omp barrier //#pragma omp for private(i,m,hlst1,Ihnew) for (i=X1;i<=nx;i++) { if ((uhh[i] != 0.0) || (uhh[i-1] != 0.0)) { uhr[k][i][j] -= uhh[i]; hlst1 = hvol[k][i][j]; hvol[k][i][j] -= (uhh[i] - uhh[i-1]); Ihnew = 1.0 / hvol[k][i][j]; for (m=0;m<NTR;m++) { tr[m][k][i][j] *= hlst1; tr[m][k][i][j] = (tr[m][k][i][j] - (fluxtr[i][m]-fluxtr[i-1][m])) * Ihnew; } } } // #pragma omp barrier } /* j loop */ #endif /* ============================================================ */ /* now advect meridionally */ /* ============================================================ */ #ifndef ADV1D for (i=X1;i<=nx;i++) { /* Calculate the j-direction profiles (slopes) of each tracer that */ /* is being advected. */ //#pragma omp for private(j,m,minslope) for (j=Y0;j<=ny+1;j++) { for (m=0;m<NTR;m++) { minslope = 4.0*((fabs(tr[m][k][i][j+1]-tr[m][k][i][j]) < fabs(tr[m][k][i][j]-tr[m][k][i][j-1])) ? (tr[m][k][i][j+1]-tr[m][k][i][j]) : (tr[m][k][i][j]-tr[m][k][i][j-1])); slope[j][m] = vmask[i][j] * vmask[i][j-1] * (((tr[m][k][i][j+1]-tr[m][k][i][j]) * (tr[m][k][i][j]-tr[m][k][i][j-1]) < 0.0) ? 0.0 : ((fabs(tr[m][k][i][j+1]-tr[m][k][i][j-1])<fabs(minslope)) ? 0.5*(tr[m][k][i][j+1]-tr[m][k][i][j-1]) : 0.5*minslope)); } } // #pragma omp barrier /* Calculate the j-direction fluxes of each tracer, using as much */ /* the minimum of the remaining mass flux (vhr) and the half the mass */ /* in the cell plus whatever part of its half of the mass flux that */ /* the flux through the other side does not require. */ //#pragma omp for private(j,m,hup,hlos,ts2) for (j=Y0;j<=ny;j++) { if (vhr[k][i][j] == 0.0) { vhh[j] = 0.0; for (m=0;m<NTR;m++) fluxtr[j][m] = 0.0; } else if (vhr[k][i][j] < 0.0) { if (k==0 || k==1 ) { hup = (hvol[k][i][j+1]-DXDYh(i,j+1)*MLMIN); } else { hup = (hvol[k][i][j+1]-DXDYh(i,j+1)*EPSILON); } hlos = D_MAX(0.0,vhr[k][i][j+1]); if (((hup + vhr[k][i][j] - hlos) < 0.0) && ((0.5*hup + vhr[k][i][j]) < 0.0)) { vhh[j] = D_MIN(-0.5*hup,-hup+hlos); } else vhh[j] = vhr[k][i][j]; ts2 = 0.5*(1.0 + vhh[j]/(hvol[k][i][j+1])); for (m=0;m<NTR;m++) { fluxtr[j][m] = vhh[j]*(tr[m][k][i][j+1] - slope[j+1][m]*ts2); } } else { if (k==0 || k==1 ) { hup = (hvol[k][i][j]-DXDYh(i,j)*MLMIN); } else { hup = (hvol[k][i][j]-DXDYh(i,j)*EPSILON); } hlos = D_MAX(0.0,-vhr[k][i][j-1]); if (((hup - vhr[k][i][j] - hlos) < 0.0) && ((0.5*hup - vhr[k][i][j]) < 0.0)) { vhh[j] = D_MAX(0.5*hup,hup-hlos); } else vhh[j] = vhr[k][i][j]; ts2 = 0.5*(1.0 - vhh[j] / (hvol[k][i][j])); for (m=0;m<NTR;m++) { fluxtr[j][m] = vhh[j]*(tr[m][k][i][j] + slope[j][m]*ts2); } } } // #pragma omp barrier /* Calculate new tracer concentration in each cell after accounting */ /* for the j-direction fluxes. */ vhr[k][i][Y0] -= vhh[Y0]; // #pragma omp barrier //#pragma omp for private(j,m,Ihnew) for (j=Y1;j<=ny;j++) { if ((vhh[j] != 0.0) || (vhh[j-1] != 0.0)) { hlst[j] = hvol[k][i][j]; hvol[k][i][j] -= (vhh[j] - vhh[j-1]); Ihnew = 1.0 / hvol[k][i][j]; vhr[k][i][j] -= vhh[j]; for (m=0;m<NTR;m++) { tr[m][k][i][j] *= hlst[j]; tr[m][k][i][j] = (tr[m][k][i][j] - fluxtr[j][m] + fluxtr[j-1][m]) * Ihnew; } } } // #pragma omp barrier } /* i loop */ #endif } /* end of big loop over k */ /* calculate new thickness field - to be used for vertical */ /* tracer advection from updated volumes (vol + fluxes) */ #pragma omp for private(i,j,k) for (k=0; k<=NZ-1; k++) { for (i=X1; i<=nx; i++) { for (j=Y1; j<=ny; j++) { hnew[k][i][j] = hvol[k][i][j]/DXDYh(i,j); if (hnew[k][i][j] < EPSILON) { hnew[k][i][j] = EPSILON; } } } } /* ============================================================ */ /* now advect vertically */ /* ============================================================ */ #pragma omp for private(i,j,hup,hlos) for (j=Y1; j<=ny; j++) { for (i=X1; i<=nx; i++) { /* work from top to bottom - by interfaces - interface k is the */ /* interface between layer k and layer k-1. net flux at this */ /* interface is wd[[k][i][j]= ea[k][i][j] and eb[k-1][i][j] */ /* k=0 */ if (wd[0][i][j] == 0.0) { wdh[0][i][j] = 0.0; } else if (wd[0][i][j] < 0.0) { hup = hnew[0][i][j] - MLMIN; hlos = D_MAX(0.0, wd[1][i][j]); if (((hup + wd[0][i][j] - hlos) < 0.0) && ((0.5*hup + wd[0][i][j]) < 0.0)) { wdh[0][i][j] = D_MIN(-0.5*hup,-hup+hlos); } else wdh[0][i][j] = wd[0][i][j]; } else { wdh[0][i][j] = wd[0][i][j]; } } } #pragma omp for private(i,j,hup,hlos) for (j=Y1; j<=ny; j++) { for (i=X1; i<=nx; i++) { /* k=1 */ if (wd[1][i][j] == 0.0) { wdh[1][i][j] = 0.0; } else if (wd[1][i][j] > 0.0) { hup = hnew[0][i][j] - MLMIN; hlos = D_MAX(0.0, -wd[0][i][j]); if (((hup - wd[1][i][j] - hlos) < 0.0) && ((0.5*hup - wd[1][i][j]) < 0.0)) { wdh[1][i][j] = D_MAX(0.5*hup,hup-hlos); } else wdh[1][i][j] = wd[1][i][j]; } else { hup = hnew[1][i][j] - MLMIN; hlos = D_MAX(0.0,wd[2][i][j]); if (((hup + wd[1][i][j] - hlos) < 0.0) && ((0.5*hup + wd[1][i][j]) < 0.0)) { wdh[1][i][j] = D_MIN(-0.5*hup,-hup+hlos); } else wdh[1][i][j] = wd[1][i][j]; } } } #pragma omp for private(i,j,hup,hlos) for (j=Y1; j<=ny; j++) { for (i=X1; i<=nx; i++) { /* k=2 */ if (wd[2][i][j] == 0.0) { wdh[2][i][j] = 0.0; } else if (wd[2][i][j] > 0.0) { hup = hnew[1][i][j] - MLMIN; hlos = D_MAX(0.0, -wd[1][i][j]); if (((hup - wd[2][i][j] - hlos) < 0.0) && ((0.5*hup - wd[2][i][j]) < 0.0)) { wdh[2][i][j] = D_MAX(0.5*hup,hup-hlos); } else wdh[2][i][j] = wd[2][i][j]; } else { hup = hnew[2][i][j] - EPSILON; hlos = D_MAX(0.0,wd[3][i][j]); if (((hup + wd[2][i][j] - hlos) < 0.0) && ((0.5*hup + wd[2][i][j]) < 0.0)) { wdh[2][i][j] = D_MIN(-0.5*hup,-hup+hlos); } else wdh[2][i][j] = wd[2][i][j]; } } } /* k=3 --> NZ-1 */ #pragma omp for private(i,j,k,hup,hlos) for (k=3; k<=NZ-1; k++) { for (j=Y1; j<=ny; j++) { for (i=X1; i<=nx; i++) { if (wd[k][i][j] == 0.0) { wdh[k][i][j] = 0.0; } else if (wd[k][i][j] > 0.0) { hup = hnew[k-1][i][j] - EPSILON; hlos = D_MAX(0.0, -wd[k-1][i][j]); if (((hup - wd[k][i][j] - hlos) < 0.0) && ((0.5*hup - wd[k][i][j]) < 0.0)) { wdh[k][i][j] = D_MAX(0.5*hup,hup-hlos); } else wdh[k][i][j] = wd[k][i][j]; } else { hup = hnew[k][i][j] - EPSILON; if (k != NZ-1) { hlos = D_MAX(0.0,wd[k+1][i][j]); } else { hlos = 0.0; } if (((hup + wd[k][i][j] - hlos) < 0.0) && ((0.5*hup + wd[k][i][j]) < 0.0)) { wdh[k][i][j] = D_MIN(-0.5*hup,-hup+hlos); } else wdh[k][i][j] = wd[k][i][j]; } } /* k */ } /* j */ } /* i */ #pragma omp for private(i,j) for (i=X1; i<=nx; i++) for (j=Y1; j<=ny; j++) { ear[0][i][j] = wdh[0][i][j]; /* added by Curtis - bottom ebr wasn't set anywhere else */ ebr[NZ-1][i][j] = 0; } #pragma omp for private(i,j,k) for (k=1;k<=NZ-1;k++) { for (i=X1; i<=nx; i++) { for (j=Y1; j<=ny; j++) { ear[k][i][j] = 0.5 * (fabs(wdh[k][i][j]) + wdh[k][i][j]); ebr[k-1][i][j] = 0.5 * (fabs(wdh[k][i][j]) - wdh[k][i][j]); } } } #pragma omp for private(i,j,k,m,hnew0,bet,gam) for (j=Y1; j<=ny; j++) { for (i=X1; i<=nx; i++) { hnew0[i] = hnew[0][i][j]; bet[i]=1.0/(hnew[0][i][j] + ebr[0][i][j] + wdh[0][i][j]); for (m=0;m<NTR;m++) tr[m][0][i][j] = bet[i]*(hnew0[i]*tr[m][0][i][j]); } for (k=1;k<=NZ-1;k++) { for (i=X1;i<=nx;i++) { gam[k][i] = ebr[k-1][i][j] * bet[i]; bet[i]=1.0/(hnew[k][i][j] + ebr[k][i][j] + (1.0-gam[k][i])*ear[k][i][j]); for (m=0;m<NTR;m++) tr[m][k][i][j] = bet[i] * (hnew[k][i][j]*tr[m][k][i][j] + ear[k][i][j]*(tr[m][k-1][i][j]) ); } } for (m=0;m<NTR;m++) for (k=NZ-2;k>=0;k--) { for (i=X1;i<=nx;i++) { tr[m][k][i][j] += gam[k+1][i]*tr[m][k+1][i][j]; } } } /*j*/ /* update hvol with diapycnal fluxes */ #pragma omp for private(i,j,k) for (k=0;k<NZ-1;k++) { for (i=X1; i<=nx; i++) for (j=Y1; j<=ny; j++) hnew[k][i][j] += (wdh[k][i][j] - wdh[k+1][i][j]); } #pragma omp for private(i,j) for (i=X1; i<=nx; i++) for (j=Y1; j<=ny; j++) hnew[NZ-1][i][j] += wdh[NZ-1][i][j]; #pragma omp for private(i,j,k) for (k=0;k<=NZ-1;k++) for (i=X1; i<=nx; i++) for (j=Y1; j<=ny; j++) { if (hnew[k][i][j] < EPSILON) hnew[k][i][j] = EPSILON; hvol[k][i][j] = DXDYh(i,j)*hnew[k][i][j]; if ( wd[k][i][j] > 0.0 && ( wdh[k][i][j] > wd[k][i][j] )) printf("case 1 wdh[k]\n"); else if ( wd[k][i][j] < 0.0 && ( wdh[k][i][j] < wd[k][i][j] )) printf("case 2 wdh[k]\n"); wd[k][i][j] -= wdh[k][i][j]; } #else /* STANDARD_ADVECTION */ /* big loop over k */ // printf("phos(%d,%d,%d)=%g,uhtm=%g\n",0,190,26,tr[mPHOSPHATE][0][190][26], // uhtm[0][190][26]); // exit(1); //yanxu: note these are null cycles, so I comment them out //yanxu for (k=0;k<NZ;k++) { // first advect zonally //yanxu for (j=Y1;j<=ny;j++) { //yanxu for (i=X0;i<=nx+1;i++) { //yanxu for (m=0;m<NTR;m++) { // fluxtr[i][m] = uhtm[i]*(tr[m][k][i][j]); //yanxu } //yanxu } //yanxu } // now advect meridionally // null cycles again, yanxu //yanxu for (i=X1;i<=nx;i++) { //yanxu } //yanxu } /* end of big loop over k */ /* calculate new thickness field - to be used for vertical */ /* tracer advection from updated volumes (vol + fluxes) */ #pragma omp for private(i,j,k) for (k=0; k<=NZ-1; k++) for (i=X1; i<=nx; i++) for (j=Y1; j<=ny; j++) { hnew[k][i][j] = hvol[k][i][j]/DXDYh(i,j); if (hnew[k][i][j] < EPSILON) hnew[k][i][j] = EPSILON; } // now advect vertically //yanxu: null cycles //yanxu for (j=Y1; j<=ny; j++) { //yanxu for (i=X1; i<=nx; i++) { //yanxu //yanxu } //yanxu } #endif /* STANDARD_ADVECTION */ /* zonal re-entrance */ #pragma omp for private(j,k,m) for (j=0;j<NYMEM;j++) { for (k=0;k<NZ;k++) { uhr[k][nx+1][j] = uhr[k][2][j]; uhr[k][nx+2][j] = uhr[k][3][j]; uhr[k][0][j] = uhr[k][nx-1][j]; uhr[k][1][j] = uhr[k][nx][j]; vhr[k][nx+1][j] = vhr[k][2][j]; vhr[k][nx+2][j] = vhr[k][3][j]; vhr[k][0][j] = vhr[k][nx-1][j]; vhr[k][1][j] = vhr[k][nx][j]; hvol[k][nx+1][j] = hvol[k][2][j]; hvol[k][nx+2][j] = hvol[k][3][j]; hvol[k][0][j] = hvol[k][nx-1][j]; hvol[k][1][j] = hvol[k][nx][j]; for (m=0;m<NTR;m++) { tr[m][k][nx+1][j] = tr[m][k][2][j]; tr[m][k][nx+2][j] = tr[m][k][3][j]; tr[m][k][0][j] = tr[m][k][nx-1][j]; tr[m][k][1][j] = tr[m][k][nx][j]; } } } /* meridional re-entrance */ /* meridional re-entrance */ #pragma omp for private(i,k,ii,m) for (i=2;i<=nx;i++) { ii = 363 - i; for (k=0;k<NZ;k++) { uhr[k][ii][ny+1] = (-1)*uhr[k][i][ny]; uhr[k][ii][ny+2] = (-1)*uhr[k][i][ny-1]; vhr[k][ii][ny+1] = (-1)*vhr[k][i][ny]; vhr[k][ii][ny+2] = (-1)*vhr[k][i][ny-1]; hvol[k][ii][ny+1] = hvol[k][i][ny]; hvol[k][ii][ny+2] = hvol[k][i][ny-1]; for (m=0;m<NTR;m++) { tr[m][k][ii][ny+1] = tr[m][k][i][ny]; tr[m][k][ii][ny+2] = tr[m][k][i][ny-1]; } } } } // omp #ifdef STANDARD_ADVECTION # ifdef WRTTS printf("itt = %i\n",itt); wrts = (double)(itts)+(double)(itt)/11.+0.0001; write_ts(wrts); # endif // WRTTS //**************************************************************************************** } /* end of temp itt iteration loop */ //**************************************************************************************** #endif fract1 = (double)(NTSTEP-itts) / (double)NTSTEP; fract2 = 1.0 - fract1; //#pragma omp parallel for private(i,j,k) schedule(dynamic) for (k=0;k<=NZ-1;k++) { for (i=X1; i<=nx; i++) { for (j=Y1; j<=ny; j++) { # ifdef USE_CALC_H h[k][i][j] = hnew[k][i][j]; if (h[k][i][j] < 0.0) printf("tracadv l 796 - h[%d][%d][%d] = %g\n", k,i,j, h[k][i][j]); # else h[k][i][j] = fract1*hstart[k][i][j] + fract2*hend[k][i][j]; //BX h[k][i][j] = hend[k][i][j]; # endif #ifdef HTEST htest[k][i][j] = hnew[k][i][j]; printf("htest(%d,%d,%d)=%g,hend=%g\n", k,i,j, // htest[k][i][j] = h[k][i][j]; #endif } } } //HF // zonal re-entrance //#pragma omp parallel for private(j,k) schedule(dynamic) for (k=0;k<NZ;k++) { for (j=0;j<=NYMEM-1;j++) { h[k][nx+1][j] = h[k][2][j]; h[k][nx+2][j] = h[k][3][j]; h[k][0][j] = h[k][nx-1][j]; h[k][1][j] = h[k][nx][j]; } } // meridional re-entrance //#pragma omp parallel for private(i,k,ii) schedule(dynamic) for (i=2;i<=nx;i++) { ii = 363 - i; for (k=0;k<NZ;k++) { h[k][ii][ny+1] = h[k][i][ny]; h[k][ii][ny+2] = h[k][i][ny-1]; } } //HF-e # ifdef WRTTS printf("End of tracadv\n"); wrts = (double)(itts)+(double)(itt)/11.+0.0005; write_ts(wrts); # endif // WRTTS #ifdef DIFFUSE_TRACER if ((KD>0.0) || (KDML>0.0)) { diffuse_tracer(); } #endif pstage=4; print_tr(pstage); if ((KHTR>0.0)) { tracer_hordiff(); } pstage=5; print_tr(pstage); free3d(hvol, NZ); free3d(uhr, NZ); free3d(vhr, NZ); free3d(ear, NZ); free3d(ebr, NZ); free3d(wdh, NZ); #if !defined AGE2 && !defined AGE3 free3d(hnew, NZ); #endif }
void set_metrics(void) { int i,j; extern double areagr[NXMEM][NYMEM]; double Iareagr[NXMEM][NYMEM]; /* Calculate the values of the metric terms that might be used */ /* and save them in arrays. */ #ifdef CARTESIAN /* On a cartesian grid, the various DX... and DY... macros all */ /* point to the same scalars. */ i = 0; j = 0; DXh(i,j) = RE * LENLON * M_PI / (180.0 * NXTOT); DYh(i,j) = RE * LENLAT * M_PI / (180.0 * NYTOT); IDXh(i,j) = 1.0 / DXh(i,j); IDYh(i,j) = 1.0 / DYh(i,j); DXDYh(i,j) = DXh(i,j) * DYh(i,j); IDXDYh(i,j) = IDXh(i,j) * IDYh(i,j); for (j=Y0-1;j<=ny+2;j++) latq[j] = LOWLAT + LENLAT*(double)(j-Y0+Y0abs)/(double)NYTOT; for (j=Y0;j<=ny;j++) lath[j] = LOWLAT + LENLAT*((double)(j-Y0+Y0abs)-0.5)/(double)NYTOT; for (i=X0-1;i<=nx+2;i++) lonq[i] = WESTLON + LENLON*(double)(i-X0+X0abs)/(double)NXTOT; for (i=X0;i<=nx;i++) lonh[i] = WESTLON + LENLON*((double)(i-X0+X0abs)-0.5)/(double)NXTOT; # if defined(PARALLEL_Y) && !defined(PARALLEL_IO) && defined(NETCDF_OUTPUT) for (j=Y0;j<=NYTOT+Y0;j++) latq_tot[j] = LOWLAT + LENLAT*(double)(j-Y0)/(double)NYTOT; for (j=Y0;j<=NYTOT+Y0;j++) lath_tot[j] = LOWLAT + LENLAT*((double)(j-Y0)-0.5)/(double)NYTOT; # endif # if defined(PARALLEL_X) && !defined(PARALLEL_IO) && defined(NETCDF_OUTPUT) for (i=X0;i<=NXTOT+X0;i++) lonq_tot[i] = WESTLON + LENLON*(double)(i-X0)/(double)NXTOT; for (i=X0;i<=NXTOT+X0;i++) lonh_tot[i] = WESTLON + LENLON*((double)(i-X0)-0.5)/(double)NXTOT; # endif #else /* All of the metric terms should be defined over the domain from */ /* X0-1 to nx+2. Outside of the physical domain, both the metrics */ /* and their inverses may be set to zero. */ /* Any points that are outside of the computational domain should */ /* have their values set to zero _BEFORE_ setting the other metric */ /* terms, because these macros may or may not expand to 2-dimensional */ /* arrays. */ /* The metric terms within the computational domain are set here. */ # ifdef ISOTROPIC { double C0, I_C0, yq, yh, jd; double fnRef, yRef; /* fnRef is the value of the integral of */ /* 1/(dx*cos(lat)) from the equator to a */ /* reference latitude, while yRef is the */ /* j-index of that reference latitude. */ int itt1, itt2; C0 = M_PI*((double) LENLON / (double) (180*NXTOT)); I_C0 = 1.0 / C0; /* With an isotropic grid, the north-south extent of the domain, */ /* the east-west extent, and the number of grid points in each */ /* direction are _not_ independent. Here the north-south extent */ /* will be determined to fit the east-west extent and the number of */ /* grid points. The grid is perfectly isotropic. */ # ifdef NORTHREFERENCE /* The following 2 lines fixes the refererence latitude at the */ /* northern edge of the domain, LOWLAT+LENLAT at j=NYTOT+Y0. */ yRef = Y0abs - Y0 - NYTOT; fnRef = I_C0 * INTSECY(((LOWLAT+LENLAT)*M_PI/180.0)); # else /* The following line sets the refererence latitude LOWLAT at j=Y0. */ yRef = Y0abs - Y0; fnRef = I_C0 * INTSECY((LOWLAT*M_PI/180.0)); # endif /* Everything else should pretty much work. */ yq = LOWLAT*M_PI/180.0; /* If the model is in parallel in the Y-direction, do the same set of */ /* calculations which would occur on a single processor. */ for (j=Y0-1-Y0abs;j<Y0-1;j++) { jd = fnRef + (double) (j + yRef) - 0.5; yh = find_root(jd,yq,&itt1); jd = fnRef + (double) (j + yRef); yq = find_root(jd,yh,&itt2); # if defined(PARALLEL_Y) && !defined(PARALLEL_IO) && defined(NETCDF_OUTPUT) latq_tot[j+Y0abs] = yq*180.0/M_PI; lath_tot[j+Y0abs] = yh*180.0/M_PI; # endif } for (j=Y0-1;j<=ny+2;j++) { jd = fnRef + (double) (j + yRef) - 0.5; yh = find_root(jd,yq,&itt1); jd = fnRef + (double) (j + yRef); yq = find_root(jd,yh,&itt2); latq[j] = yq*180.0/M_PI; lath[j] = yh*180.0/M_PI; # if defined(PARALLEL_Y) && !defined(PARALLEL_IO) && defined(NETCDF_OUTPUT) latq_tot[j+Y0abs] = yq*180.0/M_PI; lath_tot[j+Y0abs] = yh*180.0/M_PI; # endif for (i=X0-1;i<=nx+2;i++) { DXq(i,j) = cos(yq) * (RE * C0); DYq(i,j) = DXq(i,j); DXv(i,j) = DXq(i,j); DYv(i,j) = DYq(i,j); DXh(i,j) = cos(yh) * (RE * C0); DYh(i,j) = DXh(i,j); DXu(i,j) = DXh(i,j); DYu(i,j) = DYh(i,j); } } # if defined(PARALLEL_Y) && !defined(PARALLEL_IO) && defined(NETCDF_OUTPUT) for (j=ny+3;j<=NYTOT+Y0-Y0abs;j++) { jd = fnRef + (double) (j + yRef) - 0.5; yh = find_root(jd,yq,&itt1); jd = fnRef + (double) (j + yRef); yq = find_root(jd,yh,&itt2); latq_tot[j+Y0abs] = yq*180.0/M_PI; lath_tot[j+Y0abs] = yh*180.0/M_PI; } # endif } for (i=X0-1;i<=nx+2;i++) lonq[i] = WESTLON + LENLON*(double)(i-X0+X0abs)/(double)NXTOT; for (i=X0;i<=nx;i++) lonh[i] = WESTLON + LENLON*((double)(i-X0+X0abs)-0.5)/(double)NXTOT; # if defined(PARALLEL_X) && !defined(PARALLEL_IO) && defined(NETCDF_OUTPUT) for (i=X0;i<=NXTOT+X0;i++) lonq_tot[i] = WESTLON + LENLON*(double)(i-X0)/(double)NXTOT; for (i=X0;i<=NXTOT+X0;i++) lonh_tot[i] = WESTLON + LENLON*((double)(i-X0)-0.5)/(double)NXTOT; # endif # else /* ISOTROPIC */ /* This code implements latitude/longitude coordinates on a sphere. */ //BX-a // for (j=0;j<=NYTOT-1;j++) { //HF for (j=Y0-1;j<=ny+2;j++) { for (j=2;j<NYTOT+2;j++) { latq[j] = qlat[j]; lath[j] = hlat[j]; } //BX-e /* HF for (i=0;i<=NXTOT-1;i++) { for (j=0;j<=NYTOT-1;j++) { DXq(i+2,j+2) = dxq[j][i]; DYq(i+2,j+2) = dyq[j][i]; DXv(i+2,j+2) = dxv[j][i]; DYv(i+2,j+2) = dyv[j][i]; DXh(i+2,j+2) = dxh[j][i]; DYh(i+2,j+2) = dyh[j][i]; DXu(i+2,j+2) = dxu[j][i]; DYu(i+2,j+2) = dyu[j][i]; } } */ //BX-a for (i=X0-1;i<=nx+2;i++) lonq[i] = WESTLON + LENLON*(double)(i-X0+X0abs)/(double)NXTOT; for (i=X0;i<=nx;i++) lonh[i] = WESTLON + LENLON*((double)(i-X0+X0abs)-0.5)/(double)NXTOT; //BX-e #endif /* The remaining code should not be changed. */ for (i=X0-1;i<=nx+2;i++) { for (j=Y0-1;j<=ny+2;j++) { IDXh(i,j) = (DXh(i,j) > 0.0) ? (1.0 / DXh(i,j)) : 0.0; IDXu(i,j) = (DXu(i,j) > 0.0) ? (1.0 / DXu(i,j)) : 0.0; IDXv(i,j) = (DXv(i,j) > 0.0) ? (1.0 / DXv(i,j)) : 0.0; IDXq(i,j) = (DXq(i,j) > 0.0) ? (1.0 / DXq(i,j)) : 0.0; IDYh(i,j) = (DYh(i,j) > 0.0) ? (1.0 / DYh(i,j)) : 0.0; IDYu(i,j) = (DYu(i,j) > 0.0) ? (1.0 / DYu(i,j)) : 0.0; IDYv(i,j) = (DYv(i,j) > 0.0) ? (1.0 / DYv(i,j)) : 0.0; IDYq(i,j) = (DYq(i,j) > 0.0) ? (1.0 / DYq(i,j)) : 0.0; //HF alternate def: DXDYh(i,j) = DXh(i,j) * DYh(i,j); //HF alternate def: IDXDYh(i,j) = IDXh(i,j) * IDYh(i,j); DXDYq(i,j) = DXq(i,j) * DYq(i,j); IDXDYq(i,j) = IDXq(i,j) * IDYq(i,j); Iareagr[i][j] = (areagr[i][j] > 0.0) ? (1.0 / areagr[i][j]) : 0.0; DXDYh(i,j) = areagr[i][j]; IDXDYh(i,j) = Iareagr[i][j]; } } #endif /* CARTESIAN */ }