Parent Directory | Revision Log
Define macros of the preprocessor CPP_IIM, CPP_JJM, CPP_LLM so we can control the resolution from the compilation command, and automate compilation for several resolutions. In module yoethf_m, transform variables into named constants. So we do not need procedure yoethf any longer. Bug fix in program test_inter_barxy, missing calls to fyhyp and fxhyp, and definition of rlatu. Remove variable iecri of module conf_gcm_m. The files dyn_hist*.nc are written every time step. We are simplifying the output system, pending replacement by a whole new system. Modify possible value of vert_sampling from "param" to "strato_custom", following LMDZ. Default values of corresponding namelist variables are now the values used for LMDZ CMIP6.
1 | module filtreg_hemisph_m |
2 | |
3 | implicit none |
4 | |
5 | contains |
6 | |
7 | subroutine filtreg_hemisph(champ, sdd, matri) |
8 | |
9 | USE dimensions, ONLY: iim |
10 | |
11 | REAL, intent(inout):: champ(:, :, :) ! (iim + 1, :, :) |
12 | REAL, intent(in):: sdd(:) ! (iim) xprim[uv]^{\pm 1/2} |
13 | |
14 | real, intent(in):: matri(:, :, :) ! (iim, iim, :) |
15 | ! filtering matrix, last dimension is latitude |
16 | |
17 | ! Local: |
18 | integer l, j |
19 | |
20 | !----------------------------------------------------------------- |
21 | |
22 | forall (j = 1:size(champ, 2), l = 1:size(champ, 3)) & |
23 | champ(:iim, j, l) = champ(:iim, j, l) & |
24 | + matmul(matri(:, :, j), champ(:iim, j, l) * sdd) / sdd |
25 | champ(iim + 1, :, :) = champ(1, :, :) |
26 | |
27 | end subroutine filtreg_hemisph |
28 | |
29 | end module filtreg_hemisph_m |
ViewVC Help | |
Powered by ViewVC 1.1.21 |