/[lmdze]/trunk/libf/dyn3d/Dissipation/divgrad2.f90
ViewVC logotype

Annotation of /trunk/libf/dyn3d/Dissipation/divgrad2.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 65 - (hide annotations)
Thu Sep 20 09:57:03 2012 UTC (11 years, 8 months ago) by guez
File size: 1252 byte(s)
Removed unused procedure "divgrad".

In procedure "dissip", save memory by using intermediary arrays "gdx"
and "gdy" several times instead of additional array "grx" and "gry".

In procedure "inidissip", write "dtdiss * teta*" instead of "teta*".

In "comvert", change name of s_sampling from "LMD5" to "tropo" and
from "strato2" to "strato".

1 guez 65 module divgrad2_m
2 guez 3
3 guez 65 IMPLICIT NONE
4 guez 3
5 guez 65 contains
6 guez 3
7 guez 65 SUBROUTINE divgrad2(klevel, h, deltapres, lh, divgra, cdivh)
8 guez 3
9 guez 65 ! From LMDZ4/libf/dyn3d/divgrad2.F, version 1.1.1.1 2004/05/19 12:53:06
10     ! P. Le Van
11    
12     ! Calcul de div(grad) de (pext * h)
13    
14     USE comgeom, ONLY: cuvscvgam2, cvuscugam2, unsair_gam2, unsapolnga2, &
15     unsapolsga2
16     USE laplacien_m, ONLY: laplacien
17     USE paramet_m, ONLY: ip1jmp1
18    
19     INTEGER, intent(in):: klevel
20     REAL, intent(in):: h(ip1jmp1, klevel), deltapres(ip1jmp1, klevel)
21     integer, intent(in):: lh
22     REAL, intent(out):: divgra(ip1jmp1, klevel)
23     real, intent(in):: cdivh
24    
25     ! Variables locales
26     REAL sqrtps(ip1jmp1, klevel)
27     INTEGER iter
28    
29     !-----------------------------------------------------------------
30    
31     divgra = h
32     CALL laplacien(klevel, divgra)
33     sqrtps = SQRT(deltapres)
34     divgra = divgra * sqrtps
35    
36     ! ItĂ©ration de l'opĂ©rateur laplacien_gam
37     DO iter = 1, lh - 2
38     CALL laplacien_gam(klevel, cuvscvgam2, cvuscugam2, unsair_gam2, &
39     unsapolnga2, unsapolsga2, divgra, divgra)
40     ENDDO
41    
42     divgra = divgra * sqrtps
43     CALL laplacien(klevel, divgra)
44     divgra = (-1.)**lh * cdivh * divgra / deltapres
45    
46     END SUBROUTINE divgrad2
47    
48     end module divgrad2_m

  ViewVC Help
Powered by ViewVC 1.1.21