source: codes/icosagcm/trunk/src/geometry.f90 @ 15

Last change on this file since 15 was 15, checked in by ymipsl, 12 years ago

Update on 3D dynamic

YM

File size: 17.6 KB
Line 
1MODULE geometry
2  USE field_mod
3 
4  TYPE t_geometry
5    TYPE(t_field),POINTER :: xyz_i(:)
6    TYPE(t_field),POINTER :: centroid(:)
7    TYPE(t_field),POINTER :: xyz_e(:)
8    TYPE(t_field),POINTER :: xyz_v(:)
9    TYPE(t_field),POINTER :: ep_e(:)
10    TYPE(t_field),POINTER :: et_e(:)
11    TYPE(t_field),POINTER :: elon_i(:)
12    TYPE(t_field),POINTER :: elat_i(:)
13    TYPE(t_field),POINTER :: elon_e(:)
14    TYPE(t_field),POINTER :: elat_e(:)
15    TYPE(t_field),POINTER :: Ai(:)
16    TYPE(t_field),POINTER :: Av(:)
17    TYPE(t_field),POINTER :: de(:)
18    TYPE(t_field),POINTER :: le(:)
19    TYPE(t_field),POINTER :: Riv(:)
20    TYPE(t_field),POINTER :: Riv2(:)
21    TYPE(t_field),POINTER :: ne(:)
22    TYPE(t_field),POINTER :: Wee(:)
23    TYPE(t_field),POINTER :: bi(:)
24    TYPE(t_field),POINTER :: fv(:)
25  END TYPE t_geometry   
26
27  TYPE(t_geometry),TARGET :: geom
28 
29  REAL(rstd),POINTER :: Ai(:)          ! area of a cell
30  REAL(rstd),POINTER :: xyz_i(:,:)     ! coordinate of the center of the cell (voronoi)
31  REAL(rstd),POINTER :: centroid(:,:)  ! coordinate of the centroid of the cell
32  REAL(rstd),POINTER :: xyz_e(:,:)     ! coordinate of a wind point on the cell on a edge
33  REAL(rstd),POINTER :: ep_e(:,:)      ! perpendicular unit vector of a edge (outsider)
34  REAL(rstd),POINTER :: et_e(:,:)      ! tangeantial unit vector of a edge
35  REAL(rstd),POINTER :: elon_i(:,:)    ! unit longitude vector on the center
36  REAL(rstd),POINTER :: elat_i(:,:)    ! unit latitude vector on the center
37  REAL(rstd),POINTER :: elon_e(:,:)    ! unit longitude vector on a wind point
38  REAL(rstd),POINTER :: elat_e(:,:)    ! unit latitude vector on a wind point
39  REAL(rstd),POINTER :: xyz_v(:,:)     ! coordinate of a vertex (center of the dual mesh)
40  REAL(rstd),POINTER :: Av(:)          ! area of dual mesk cell
41  REAL(rstd),POINTER :: de(:)          ! distance from a neighbour == lenght of an edge of the dual mesh
42  REAL(rstd),POINTER :: le(:)          ! lenght of a edge
43  REAL(rstd),POINTER :: Riv(:,:)       ! weight
44  REAL(rstd),POINTER :: Riv2(:,:)      ! weight
45  INTEGER,POINTER    :: ne(:,:)        ! convention for the way on the normal wind on an edge
46  REAL(rstd),POINTER :: Wee(:,:,:)     ! weight
47  REAL(rstd),POINTER :: bi(:)          ! orographie
48  REAL(rstd),POINTER :: fv(:)          ! coriolis (evaluted on a vertex)
49     
50
51
52CONTAINS
53
54  SUBROUTINE allocate_geometry
55  USE field_mod
56  IMPLICIT NONE
57 
58    CALL allocate_field(geom%Ai,field_t,type_real)
59    CALL allocate_field(geom%xyz_i,field_t,type_real,3)
60    CALL allocate_field(geom%centroid,field_t,type_real,3)
61    CALL allocate_field(geom%xyz_e,field_u,type_real,3)
62    CALL allocate_field(geom%ep_e,field_u,type_real,3)
63    CALL allocate_field(geom%et_e,field_u,type_real,3)
64    CALL allocate_field(geom%elon_i,field_t,type_real,3)
65    CALL allocate_field(geom%elat_i,field_t,type_real,3)
66    CALL allocate_field(geom%elon_e,field_u,type_real,3)
67    CALL allocate_field(geom%elat_e,field_u,type_real,3)
68    CALL allocate_field(geom%xyz_v,field_z,type_real,3)
69    CALL allocate_field(geom%de,field_u,type_real)
70    CALL allocate_field(geom%le,field_u,type_real)
71    CALL allocate_field(geom%bi,field_t,type_real)
72    CALL allocate_field(geom%Av,field_z,type_real)
73    CALL allocate_field(geom%Riv,field_t,type_real,6)
74    CALL allocate_field(geom%Riv2,field_t,type_real,6)
75    CALL allocate_field(geom%ne,field_t,type_integer,6)
76    CALL allocate_field(geom%Wee,field_u,type_real,5,2)
77    CALL allocate_field(geom%bi,field_t,type_real)
78    CALL allocate_field(geom%fv,field_z,type_real)
79
80  END SUBROUTINE allocate_geometry
81 
82 
83  SUBROUTINE swap_geometry(ind)
84  USE field_mod
85  IMPLICIT NONE
86    INTEGER,INTENT(IN) :: ind
87    Ai=geom%Ai(ind)
88    xyz_i=geom%xyz_i(ind)
89    centroid=geom%centroid(ind)
90    xyz_e=geom%xyz_e(ind)
91    ep_e=geom%ep_e(ind)
92    et_e=geom%et_e(ind)
93    elon_i=geom%elon_i(ind)
94    elat_i=geom%elat_i(ind)
95    elon_e=geom%elon_e(ind)
96    elat_e=geom%elat_e(ind)
97    xyz_v=geom%xyz_v(ind)
98    de=geom%de(ind)
99    le=geom%le(ind)
100    Av=geom%Av(ind)
101    Riv=geom%Riv(ind)
102    Riv2=geom%Riv2(ind)
103    ne=geom%ne(ind)
104    Wee=geom%Wee(ind)
105    bi=geom%bi(ind)
106    fv=geom%fv(ind)
107   
108  END SUBROUTINE swap_geometry
109 
110  SUBROUTINE optimize_geometry
111  USE metric
112  USE spherical_geom_mod
113  USE domain_mod
114  USE dimensions
115  USE transfert_mod
116  USE vector
117  IMPLICIT NONE
118    INTEGER,PARAMETER :: nb_it=3000
119    TYPE(t_domain),POINTER :: d
120    INTEGER :: ind,it,i,j,n,k
121    REAL(rstd) :: x1(3),x2(3)
122    REAL(rstd) :: vect(3,6)
123    REAL(rstd) :: centr(3)
124    REAL(rstd) :: sum
125    LOGICAL    :: check
126   
127    DO ind=1,ndomain
128      d=>domain(ind)
129      CALL swap_dimensions(ind)
130      CALL swap_geometry(ind)
131      DO j=jj_begin,jj_end
132        DO i=ii_begin,ii_end
133          n=(j-1)*iim+i
134          xyz_i(n,:)=d%xyz(:,i,j) 
135        ENDDO
136      ENDDO
137    ENDDO
138   
139    CALL transfert_request(geom%xyz_i,req_i1)
140   
141    DO ind=1,ndomain
142      CALL swap_dimensions(ind)
143      CALL swap_geometry(ind)
144      DO j=jj_begin,jj_end
145        DO i=ii_begin,ii_end
146          n=(j-1)*iim+i
147          DO k=0,5
148            x1(:) = xyz_i(n+t_pos(k+1),:)
149            x2(:) = xyz_i(n+t_pos(MOD(k+1,6)+1),:)
150            if (norm(x1-x2)<1e-16) x2(:) = xyz_i(n+t_pos(MOD(k+2,6)+1),:)
151            CALL circumcenter(xyz_i(n,:), x1, x2, xyz_v(n+z_pos(k+1),:))
152          ENDDO
153        ENDDO
154      ENDDO
155    ENDDO
156   
157    DO ind=1,ndomain
158      d=>domain(ind)
159      CALL swap_dimensions(ind)
160      CALL swap_geometry(ind)
161      DO j=jj_begin,jj_end
162        DO i=ii_begin,ii_end
163          n=(j-1)*iim+i
164          DO k=0,5
165            x1(:) = xyz_v(n+z_pos(k+1),:)
166            x2(:) = d%vertex(:,k,i,j) 
167            IF (norm(x1-x2)>1e-10) THEN
168              PRINT*,"vertex diff ",ind,i,j,k
169              PRINT*,x1
170              PRINT*,x2
171            ENDIF
172          ENDDO
173        ENDDO
174      ENDDO
175    ENDDO
176   
177   
178    DO it=1,nb_it
179      IF (MOD(it,100)==0) THEN
180        check=.TRUE.
181      ELSE
182        check=.FALSE.
183      ENDIF
184     
185      sum=0
186      DO ind=1,ndomain
187        CALL swap_dimensions(ind)
188        CALL swap_geometry(ind)
189        DO j=jj_begin,jj_end
190          DO i=ii_begin,ii_end
191            n=(j-1)*iim+i
192            vect(:,1)=xyz_v(n+z_rup,:)
193            vect(:,2)=xyz_v(n+z_up,:)
194            vect(:,3)=xyz_v(n+z_lup,:)
195            vect(:,4)=xyz_v(n+z_ldown,:)
196            vect(:,5)=xyz_v(n+z_down,:)
197            vect(:,6)=xyz_v(n+z_rdown,:)
198            CALL compute_centroid(vect,6,centr)
199            IF (check) THEN
200              sum=MAX(sum,norm(xyz_i(n,:)-centr(:)))
201            ENDIF
202            xyz_i(n,:)=centr(:)
203          ENDDO
204        ENDDO
205!        i=ii_begin ; j=jj_begin ; n=(j-1)*iim+i ; xyz_i(n,:)=domain(ind)%xyz(:,i,j)
206!        i=ii_begin ; j=jj_end   ; n=(j-1)*iim+i ; xyz_i(n,:)=domain(ind)%xyz(:,i,j)
207!        i=ii_end   ; j=jj_begin ; n=(j-1)*iim+i ;   xyz_i(n,:)=domain(ind)%xyz(:,i,j)
208!        PRINT *,"Pb ?? : "
209!        PRINT *, xyz_i(n,:), domain(ind)%xyz(:,i,j), norm(xyz_i(n,:)- domain(ind)%xyz(:,i,j))
210!        i=ii_end   ; j=jj_end   ; n=(j-1)*iim+i ; xyz_i(n,:)=domain(ind)%xyz(:,i,j)
211       
212      ENDDO
213
214      IF (check) THEN
215!        sum=sum/(ndomain*ii_nb*jj_nb)
216        PRINT *,"it = ",it,"  diff centroid circumcenter ",sum
217      ENDIF
218     
219     CALL transfert_request(geom%xyz_i,req_i1)
220
221    DO ind=1,ndomain
222      CALL swap_dimensions(ind)
223      CALL swap_geometry(ind)
224      DO j=jj_begin,jj_end
225        DO i=ii_begin,ii_end
226          n=(j-1)*iim+i
227          DO k=0,5
228            x1(:) = xyz_i(n+t_pos(k+1),:)
229            x2(:) = xyz_i(n+t_pos(MOD(k+1,6)+1),:)
230            if (norm(x1-x2)<1e-16) x2(:) = xyz_i(n+t_pos(MOD(k+2,6)+1),:)
231            CALL circumcenter(xyz_i(n,:), x1, x2, xyz_v(n+z_pos(k+1),:))
232          ENDDO
233        ENDDO
234      ENDDO
235    ENDDO
236
237     
238    ENDDO
239   
240  END SUBROUTINE optimize_geometry
241 
242   
243  SUBROUTINE set_geometry
244  USE metric
245  USE vector
246  USE spherical_geom_mod
247  USE domain_mod
248  USE dimensions
249  USE transfert_mod
250  IMPLICIT NONE
251
252    REAL(rstd) :: surf(6)
253    REAL(rstd) :: surf_v(6)
254    REAL(rstd) :: vect(3,6)
255    REAL(rstd) :: centr(6)
256    REAL(rstd) :: vet(3),vep(3)
257    INTEGER :: ind,i,j,k,n
258    TYPE(t_domain),POINTER :: d
259    REAL(rstd) ::  S
260    REAL(rstd) :: w(6)
261    REAL(rstd) :: lon,lat
262    INTEGER :: ii_glo,jj_glo
263    REAL(rstd) :: S1,S2
264         
265     
266    CALL optimize_geometry
267   
268    DO ind=1,ndomain
269      d=>domain(ind)
270      CALL swap_dimensions(ind)
271      CALL swap_geometry(ind)
272      DO j=jj_begin-1,jj_end+1
273        DO i=ii_begin-1,ii_end+1
274          n=(j-1)*iim+i
275
276          DO k=0,5
277            ne(n,k+1)=d%ne(k,i,j)
278          ENDDO
279                 
280!          xyz_i(n,:)=d%xyz(:,i,j)
281!          xyz_v(n+z_up,:)=d%vertex(:,vup-1,i,j)
282!          xyz_v(n+z_down,:)=d%vertex(:,vdown-1,i,j)
283
284          vect(:,1)=xyz_v(n+z_rup,:)
285          vect(:,2)=xyz_v(n+z_up,:)
286          vect(:,3)=xyz_v(n+z_lup,:)
287          vect(:,4)=xyz_v(n+z_ldown,:)
288          vect(:,5)=xyz_v(n+z_down,:)
289          vect(:,6)=xyz_v(n+z_rdown,:)
290          CALL compute_centroid(vect,6,centr)
291          centroid(n,:)=centr(:)
292
293               
294          CALL xyz2lonlat(xyz_v(n+z_up,:),lon,lat)
295          fv(n+z_up)=2*sin(lat)*omega
296          CALL xyz2lonlat(xyz_v(n+z_down,:),lon,lat)
297          fv(n+z_down)=2*sin(lat)*omega
298         
299          bi(n)=0. 
300       
301!          CALL dist_cart(d%xyz(:,i,j),d%neighbour(:,right-1,i,j),de(n+u_right))
302!          CALL dist_cart(d%xyz(:,i,j),d%neighbour(:,lup-1,i,j),de(n+u_lup))
303!          CALL dist_cart(d%xyz(:,i,j),d%neighbour(:,ldown-1,i,j),de(n+u_ldown))
304
305          CALL dist_cart(xyz_i(n,:),xyz_i(n+t_right,:),de(n+u_right))
306          CALL dist_cart(xyz_i(n,:),xyz_i(n+t_lup,:),de(n+u_lup))
307          CALL dist_cart(xyz_i(n,:),xyz_i(n+t_ldown,:),de(n+u_ldown))
308         
309!          CALL div_arc_bis(d%xyz(:,i,j),d%neighbour(:,right-1,i,j),0.5,xyz_e(n+u_right,:))
310!          CALL div_arc_bis(d%xyz(:,i,j),d%neighbour(:,lup-1,i,j),0.5,xyz_e(n+u_lup,:))
311!          CALL div_arc_bis(d%xyz(:,i,j),d%neighbour(:,ldown-1,i,j),0.5,xyz_e(n+u_ldown,:))
312
313          CALL div_arc_bis(xyz_i(n,:),xyz_i(n+t_right,:),0.5,xyz_e(n+u_right,:))
314          CALL div_arc_bis(xyz_i(n,:),xyz_i(n+t_lup,:),0.5,xyz_e(n+u_lup,:))
315          CALL div_arc_bis(xyz_i(n,:),xyz_i(n+t_ldown,:),0.5,xyz_e(n+u_ldown,:))
316         
317!          CALL dist_cart(d%vertex(:,vrdown-1,i,j),d%vertex(:,vrup-1,i,j),le(n+u_right))
318!          CALL dist_cart(d%vertex(:,vup-1,i,j),d%vertex(:,vlup-1,i,j),le(n+u_lup))
319!          CALL dist_cart(d%vertex(:,vldown-1,i,j),d%vertex(:,vdown-1,i,j),le(n+u_ldown))
320
321          CALL dist_cart(xyz_v(n+z_rdown,:), xyz_v(n+z_rup,:),le(n+u_right))
322          CALL dist_cart(xyz_v(n+z_up,:), xyz_v(n+z_lup,:),le(n+u_lup))
323          CALL dist_cart(xyz_v(n+z_ldown,:), xyz_v(n+z_down,:),le(n+u_ldown))
324         
325          Ai(n)=0
326          DO k=0,5
327!            CALL surf_triangle(d%xyz(:,i,j),d%neighbour(:,k,i,j),d%neighbour(:,MOD((k+1+6),6),i,j),surf_v(k+1))
328!            CALL surf_triangle(d%xyz(:,i,j),d%vertex(:,MOD((k-1+6),6),i,j),d%vertex(:,k,i,j),surf(k+1))
329            CALL surf_triangle(xyz_i(n,:),xyz_i(n+t_pos(k+1),:),xyz_i(n+t_pos(MOD((k+1+6),6)+1),:),surf_v(k+1))
330            CALL surf_triangle(xyz_i(n,:),xyz_v(n+z_pos(MOD((k-1+6),6)+1),:),xyz_v(n+z_pos(k+1),:),surf(k+1))
331            Ai(n)=Ai(n)+surf(k+1)
332            IF (i==ii_end .AND. j==jj_begin) THEN
333              IF (Ai(n)<1e20) THEN
334              ELSE
335                PRINT *,"PB !!",Ai(n),k,surf(k+1)
336                PRINT*,xyz_i(n,:),xyz_v(n+z_pos(MOD((k-1+6),6)+1),:),xyz_v(n+z_pos(k+1),:)
337              ENDIF
338            ENDIF
339          ENDDO
340
341
342          CALL cross_product2(xyz_v(n+z_rdown,:),xyz_v(n+z_rup,:),vep)
343          IF (norm(vep)>1e-30) THEN
344            vep(:)=vep(:)/norm(vep)
345            CALL cross_product2(vep,xyz_e(n+u_right,:),vet)
346            vet(:)=vet(:)/norm(vet)
347            ep_e(n+u_right,:)=vep(:)*ne(n,right)
348            et_e(n+u_right,:)=vet(:)*ne(n,right)
349          ENDIF
350
351          CALL cross_product2(xyz_v(n+z_up,:),xyz_v(n+z_lup,:),vep)
352          IF (norm(vep)>1e-30) THEN
353            vep(:)=vep(:)/norm(vep)
354            CALL cross_product2(vep,xyz_e(n+u_lup,:),vet)
355            vet(:)=vet(:)/norm(vet)
356            ep_e(n+u_lup,:)=vep(:)*ne(n,lup)
357            et_e(n+u_lup,:)=vet(:)*ne(n,lup)
358          ENDIF
359
360          CALL cross_product2(xyz_v(n+z_ldown,:),xyz_v(n+z_down,:),vep)
361          IF (norm(vep)>1e-30) THEN
362            vep(:)=vep(:)/norm(vep)
363            CALL cross_product2(vep,xyz_e(n+u_ldown,:),vet)
364            vet(:)=vet(:)/norm(vet)
365            ep_e(n+u_ldown,:)=vep(:)*ne(n,ldown)
366            et_e(n+u_ldown,:)=vet(:)*ne(n,ldown)
367          ENDIF
368
369          CALL xyz2lonlat(xyz_i(n,:),lon,lat)
370          elon_i(n,1) = -sin(lon)
371          elon_i(n,2) = cos(lon)
372          elon_i(n,3) = 0
373          elat_i(n,1) = -cos(lon)*sin(lat)
374          elat_i(n,2) = -sin(lon)*sin(lat)
375          elat_i(n,3) = cos(lat)
376
377         
378          CALL xyz2lonlat(xyz_e(n+u_right,:),lon,lat)
379          elon_e(n+u_right,1) = -sin(lon)
380          elon_e(n+u_right,2) = cos(lon)
381          elon_e(n+u_right,3) = 0
382          elat_e(n+u_right,1) = -cos(lon)*sin(lat)
383          elat_e(n+u_right,2) = -sin(lon)*sin(lat)
384          elat_e(n+u_right,3) = cos(lat)
385         
386          CALL xyz2lonlat(xyz_e(n+u_lup,:),lon,lat)
387          elon_e(n+u_lup,1) = -sin(lon)
388          elon_e(n+u_lup,2) = cos(lon)
389          elon_e(n+u_lup,3) = 0
390          elat_e(n+u_lup,1) = -cos(lon)*sin(lat)
391          elat_e(n+u_lup,2) = -sin(lon)*sin(lat)
392          elat_e(n+u_lup,3) = cos(lat)
393 
394          CALL xyz2lonlat(xyz_e(n+u_ldown,:),lon,lat)
395          elon_e(n+u_ldown,1) = -sin(lon)
396          elon_e(n+u_ldown,2) = cos(lon)
397          elon_e(n+u_ldown,3) = 0
398          elat_e(n+u_ldown,1) = -cos(lon)*sin(lat)
399          elat_e(n+u_ldown,2) = -sin(lon)*sin(lat)
400          elat_e(n+u_ldown,3) = cos(lat)
401 
402         
403          DO k=0,5
404!            CALL surf_triangle(d%xyz(:,i,j),d%vertex(:,k,i,j),d%neighbour(:,k,i,j),S1)
405!            CALL surf_triangle(d%xyz(:,i,j),d%vertex(:,k,i,j),d%neighbour(:,MOD(k+1+6,6),i,j),S2)
406            CALL surf_triangle(xyz_i(n,:), xyz_v(n+z_pos(k+1),:), xyz_i(n+t_pos(k+1),:),S1)
407            CALL surf_triangle(xyz_i(n,:), xyz_v(n+z_pos(k+1),:), xyz_i(n+t_pos(MOD(k+1+6,6)+1),:),S2)
408            Riv(n,k+1)=0.5*(S1+S2)/Ai(n)
409            Riv2(n,k+1)=0.5*(S1+S2)/surf_v(k+1)
410          ENDDO
411
412          DO k=1,6
413            IF (ABS(surf_v(k))<1e-30) THEN
414              Riv(n,k)=0.
415            ENDIF
416          ENDDO
417         
418          Av(n+z_up)=surf_v(vup)+1e-100
419          Av(n+z_down)=surf_v(vdown)+1e-100
420       
421        ENDDO
422      ENDDO
423     
424      DO j=jj_begin,jj_end
425        DO i=ii_begin,ii_end
426          n=(j-1)*iim+i
427   
428          CALL compute_wee(n,right,w)
429          Wee(n+u_right,:,1)=w(1:5)
430
431          CALL compute_wee(n+t_right,left,w)
432          Wee(n+u_right,:,2)=w(1:5)
433
434
435          CALL compute_wee(n,lup,w)
436          Wee(n+u_lup,:,1)=w(1:5)
437
438          CALL compute_wee(n+t_lup,rdown,w)
439          Wee(n+u_lup,:,2)=w(1:5)
440
441
442          CALL compute_wee(n,ldown,w)
443          Wee(n+u_ldown,:,1)=w(1:5)
444
445          CALL compute_wee(n+t_ldown,rup,w)
446          Wee(n+u_ldown,:,2)=w(1:5)
447
448        ENDDO
449      ENDDO
450     
451      DO j=jj_begin,jj_end
452        DO i=ii_begin,ii_end
453          n=(j-1)*iim+i
454          ii_glo=d%ii_begin_glo-d%ii_begin+i
455          jj_glo=d%jj_begin_glo-d%jj_begin+j
456         
457          IF (ii_glo==1 .AND. jj_glo==1) THEN
458            le(n+u_ldown)=0
459            xyz_v(n+z_ldown,:)=xyz_v(n+z_down,:)
460                       
461          ENDIF
462
463          IF (ii_glo==iim_glo .AND. jj_glo==1) THEN
464            le(n+u_right)=0
465            xyz_v(n+z_rdown,:)=xyz_v(n+z_rup,:)
466          ENDIF
467
468          IF (ii_glo==iim_glo .AND. jj_glo==jjm_glo) THEN
469            le(n+u_rup)=0
470            xyz_v(n+z_rup,:)=xyz_v(n+z_up,:)
471          ENDIF
472
473          IF (ii_glo==1 .AND. jj_glo==jjm_glo) THEN
474            le(n+u_lup)=0
475            xyz_v(n+z_up,:)=xyz_v(n+z_lup,:)
476          ENDIF
477         
478        ENDDO
479      ENDDO
480
481      DO j=jj_begin-1,jj_end+1
482        DO i=ii_begin-1,ii_end+1
483          n=(j-1)*iim+i
484          xyz_i(n,:)=xyz_i(n,:) * radius
485          xyz_v(n+z_up,:)=xyz_v(n+z_up,:) * radius
486          xyz_v(n+z_down,:)=xyz_v(n+z_down,:) *radius
487          de(n+u_right)=de(n+u_right) * radius
488          de(n+u_lup)=de(n+u_lup)*radius
489          de(n+u_ldown)=de(n+u_ldown)*radius
490          xyz_e(n+u_right,:)=xyz_e(n+u_right,:)*radius
491          xyz_e(n+u_lup,:)=xyz_e(n+u_lup,:)*radius
492          xyz_e(n+u_ldown,:)=xyz_e(n+u_ldown,:)*radius
493          le(n+u_right)=le(n+u_right)*radius
494          le(n+u_lup)=le(n+u_lup)*radius
495          le(n+u_ldown)=le(n+u_ldown)*radius
496          Ai(n)=Ai(n)*radius**2
497          Av(n+z_up)=Av(n+z_up)*radius**2
498          Av(n+z_down)=Av(n+z_down)*radius**2
499        ENDDO
500      ENDDO
501                 
502    ENDDO
503    CALL transfert_request(geom%Ai,req_i1)
504    CALL transfert_request(geom%centroid,req_i1)
505    CALL surf_triangle(d%xyz(:,ii_begin,jj_begin),d%xyz(:,ii_begin,jj_end),d%xyz(:,ii_end,jj_begin),S)
506!    PRINT *,"Surf triangle : ",S*20/(4*Pi)
507 
508  END SUBROUTINE set_geometry
509 
510  SUBROUTINE compute_wee(n,pos,w)
511  IMPLICIT NONE
512    INTEGER,INTENT(IN) :: n
513    INTEGER,INTENT(IN) :: pos
514    REAL(rstd),INTENT(OUT) ::w(6)
515
516    REAL(rstd) :: ne_(0:5)
517    REAL(rstd) :: Riv_(6)
518    INTEGER    :: k
519   
520   
521    DO k=0,5
522      ne_(k)=ne(n,MOD(pos-1+k+6,6)+1) 
523      Riv_(k+1)=Riv(n,MOD(pos-1+k+6,6)+1)
524    ENDDO
525         
526    w(1)=-ne_(0)*ne_(1)*(Riv_(1)-0.5)
527    w(2)=-ne_(2)*(ne_(0)*Riv_(2)-w(1)*ne_(1))
528    w(3)=-ne_(3)*(ne_(0)*Riv_(3)-w(2)*ne_(2))
529    w(4)=-ne_(4)*(ne_(0)*Riv_(4)-w(3)*ne_(3))
530    w(5)=-ne_(5)*(ne_(0)*Riv_(5)-w(4)*ne_(4))
531    w(6)=ne_(0)*ne_(5)*(Riv_(6)-0.5)
532   
533!    IF ( ABS(w(5)-w(6))>1e-20) PRINT *, "pb pour wee : w(5)!=w(6)",sum(Riv_(:))
534
535   END SUBROUTINE compute_wee
536           
537
538 
539  SUBROUTINE compute_geometry
540  IMPLICIT NONE
541    CALL allocate_geometry
542    CALL set_geometry
543   
544  END SUBROUTINE compute_geometry
545 
546END MODULE geometry
Note: See TracBrowser for help on using the repository browser.