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

Last change on this file since 427 was 427, checked in by dubos, 8 years ago

New kinetic energy (diagnostic)

File size: 22.3 KB
Line 
1MODULE geometry
2  USE field_mod
3  IMPLICIT NONE
4
5  TYPE t_geometry
6    TYPE(t_field),POINTER :: centroid(:)
7    TYPE(t_field),POINTER :: xyz_i(:)
8    TYPE(t_field),POINTER :: xyz_e(:)
9    TYPE(t_field),POINTER :: xyz_v(:)
10    TYPE(t_field),POINTER :: lon_i(:)
11    TYPE(t_field),POINTER :: lon_e(:)
12    TYPE(t_field),POINTER :: lat_i(:)
13    TYPE(t_field),POINTER :: lat_e(:)
14    TYPE(t_field),POINTER :: ep_e(:)
15    TYPE(t_field),POINTER :: et_e(:)
16    TYPE(t_field),POINTER :: elon_i(:)
17    TYPE(t_field),POINTER :: elat_i(:)
18    TYPE(t_field),POINTER :: elon_e(:)
19    TYPE(t_field),POINTER :: elat_e(:)
20    TYPE(t_field),POINTER :: Ai(:)
21    TYPE(t_field),POINTER :: Av(:)
22    TYPE(t_field),POINTER :: de(:)
23    TYPE(t_field),POINTER :: le(:)
24    TYPE(t_field),POINTER :: le_de(:) ! le/de, 0. if de=0.
25    TYPE(t_field),POINTER :: Riv(:)
26    TYPE(t_field),POINTER :: S1(:)
27    TYPE(t_field),POINTER :: S2(:)
28    TYPE(t_field),POINTER :: Riv2(:)
29    TYPE(t_field),POINTER :: ne(:)
30    TYPE(t_field),POINTER :: Wee(:)
31    TYPE(t_field),POINTER :: bi(:)
32    TYPE(t_field),POINTER :: fv(:)
33  END TYPE t_geometry   
34
35  TYPE(t_geometry),SAVE,TARGET :: geom
36
37 
38  REAL(rstd),POINTER :: Ai(:)          ! area of a cell
39!$OMP THREADPRIVATE(Ai)
40  REAL(rstd),POINTER :: centroid(:,:)  ! coordinate of the centroid of the cell
41!$OMP THREADPRIVATE(centroid)
42  REAL(rstd),POINTER :: xyz_i(:,:)     ! coordinate of the center of the cell (voronoi)
43!$OMP THREADPRIVATE(xyz_i)
44  REAL(rstd),POINTER :: xyz_e(:,:)     ! coordinate of a wind point on the cell on a edge
45!$OMP THREADPRIVATE(xyz_e)
46  REAL(rstd),POINTER :: xyz_v(:,:)     ! coordinate of a vertex (center of the dual mesh)
47!$OMP THREADPRIVATE(xyz_v)
48  REAL(rstd),POINTER :: lon_i(:)       ! longitude of the center of the cell (voronoi)
49!$OMP THREADPRIVATE(lon_i)
50  REAL(rstd),POINTER :: lon_e(:)       ! longitude of a wind point on the cell on a edge
51!$OMP THREADPRIVATE(lon_e)
52  REAL(rstd),POINTER :: lat_i(:)       ! latitude of the center of the cell (voronoi)
53!$OMP THREADPRIVATE(lat_i)
54  REAL(rstd),POINTER :: lat_e(:)       ! latitude of a wind point on the cell on a edge
55!$OMP THREADPRIVATE(lat_e)
56  REAL(rstd),POINTER :: ep_e(:,:)      ! perpendicular unit vector of a edge (outsider)
57!$OMP THREADPRIVATE(ep_e)
58  REAL(rstd),POINTER :: et_e(:,:)      ! tangeantial unit vector of a edge
59!$OMP THREADPRIVATE(et_e)
60  REAL(rstd),POINTER :: elon_i(:,:)    ! unit longitude vector on the center
61!$OMP THREADPRIVATE(elon_i)
62  REAL(rstd),POINTER :: elat_i(:,:)    ! unit latitude vector on the center
63!$OMP THREADPRIVATE(elat_i)
64  REAL(rstd),POINTER :: elon_e(:,:)    ! unit longitude vector on a wind point
65!$OMP THREADPRIVATE(elon_e)
66  REAL(rstd),POINTER :: elat_e(:,:)    ! unit latitude vector on a wind point
67!$OMP THREADPRIVATE(elat_e)
68  REAL(rstd),POINTER :: Av(:)          ! area of dual mesk cell
69!$OMP THREADPRIVATE(Av)
70  REAL(rstd),POINTER :: de(:)          ! distance from a neighbour == lenght of an edge of the dual mesh
71!$OMP THREADPRIVATE(de)
72  REAL(rstd),POINTER :: le(:)          ! lenght of a edge
73!$OMP THREADPRIVATE(le)
74  REAL(rstd),POINTER :: le_de(:)       ! le/de
75!$OMP THREADPRIVATE(le_de)
76  REAL(rstd),POINTER :: S1(:,:)        ! area of sub-triangle
77!$OMP THREADPRIVATE(S1)
78  REAL(rstd),POINTER :: S2(:,:)        ! area of sub-tirangle
79!$OMP THREADPRIVATE(S2)
80  REAL(rstd),POINTER :: Riv(:,:)       ! weight
81!$OMP THREADPRIVATE(Riv)
82  REAL(rstd),POINTER :: Riv2(:,:)      ! weight
83!$OMP THREADPRIVATE(Riv2)
84  INTEGER,POINTER    :: ne(:,:)        ! convention for the way on the normal wind on an edge
85!$OMP THREADPRIVATE(ne)
86  REAL(rstd),POINTER :: Wee(:,:,:)     ! weight
87!$OMP THREADPRIVATE(Wee)
88  REAL(rstd),POINTER :: bi(:)          ! orographie
89!$OMP THREADPRIVATE(bi)
90  REAL(rstd),POINTER :: fv(:)          ! coriolis (evaluted on a vertex)
91!$OMP THREADPRIVATE(fv)
92     
93
94  INTEGER, PARAMETER :: ne_right=1
95  INTEGER, PARAMETER :: ne_rup=-1
96  INTEGER, PARAMETER :: ne_lup=1
97  INTEGER, PARAMETER :: ne_left=-1
98  INTEGER, PARAMETER :: ne_ldown=1
99  INTEGER, PARAMETER :: ne_rdown=-1
100
101CONTAINS
102
103  SUBROUTINE allocate_geometry
104  USE field_mod
105  IMPLICIT NONE
106 
107    CALL allocate_field(geom%Ai,field_t,type_real,name='Ai')
108
109    CALL allocate_field(geom%xyz_i,field_t,type_real,3)
110    CALL allocate_field(geom%lon_i,field_t,type_real)
111    CALL allocate_field(geom%lat_i,field_t,type_real)
112    CALL allocate_field(geom%elon_i,field_t,type_real,3)
113    CALL allocate_field(geom%elat_i,field_t,type_real,3)
114    CALL allocate_field(geom%centroid,field_t,type_real,3)
115
116    CALL allocate_field(geom%xyz_e,field_u,type_real,3)
117    CALL allocate_field(geom%lon_e,field_u,type_real)
118    CALL allocate_field(geom%lat_e,field_u,type_real)
119    CALL allocate_field(geom%elon_e,field_u,type_real,3)
120    CALL allocate_field(geom%elat_e,field_u,type_real,3)
121    CALL allocate_field(geom%ep_e,field_u,type_real,3)
122    CALL allocate_field(geom%et_e,field_u,type_real,3)
123
124    CALL allocate_field(geom%xyz_v,field_z,type_real,3)
125    CALL allocate_field(geom%de,field_u,type_real)
126    CALL allocate_field(geom%le,field_u,type_real)
127    CALL allocate_field(geom%le_de,field_u,type_real)
128    CALL allocate_field(geom%bi,field_t,type_real)
129    CALL allocate_field(geom%Av,field_z,type_real)
130    CALL allocate_field(geom%S1,field_t,type_real,6)
131    CALL allocate_field(geom%S2,field_t,type_real,6)
132    CALL allocate_field(geom%Riv,field_t,type_real,6)
133    CALL allocate_field(geom%Riv2,field_t,type_real,6)
134    CALL allocate_field(geom%ne,field_t,type_integer,6)
135    CALL allocate_field(geom%Wee,field_u,type_real,5,2)
136    CALL allocate_field(geom%bi,field_t,type_real)
137    CALL allocate_field(geom%fv,field_z,type_real)
138
139  END SUBROUTINE allocate_geometry
140 
141 
142  SUBROUTINE swap_geometry(ind)
143  USE field_mod
144  IMPLICIT NONE
145    INTEGER,INTENT(IN) :: ind
146!!$OMP MASTER
147    Ai=geom%Ai(ind)
148    xyz_i=geom%xyz_i(ind)
149    centroid=geom%centroid(ind)
150    xyz_e=geom%xyz_e(ind)
151    ep_e=geom%ep_e(ind)
152    et_e=geom%et_e(ind)
153    lon_i=geom%lon_i(ind)
154    lat_i=geom%lat_i(ind)
155    lon_e=geom%lon_e(ind)
156    lat_e=geom%lat_e(ind)
157    elon_i=geom%elon_i(ind)
158    elat_i=geom%elat_i(ind)
159    elon_e=geom%elon_e(ind)
160    elat_e=geom%elat_e(ind)
161    xyz_v=geom%xyz_v(ind)
162    de=geom%de(ind)
163    le=geom%le(ind)
164    le_de=geom%le_de(ind)
165    Av=geom%Av(ind)
166    S1=geom%S1(ind)
167    S2=geom%S2(ind)
168    Riv=geom%Riv(ind)
169    Riv2=geom%Riv2(ind)
170    ne=geom%ne(ind)
171    Wee=geom%Wee(ind)
172    bi=geom%bi(ind)
173    fv=geom%fv(ind)
174!!$OMP END MASTER
175!!$OMP BARRIER   
176  END SUBROUTINE swap_geometry
177
178  SUBROUTINE update_circumcenters 
179    USE domain_mod
180    USE dimensions
181    USE spherical_geom_mod
182    USE vector
183    USE transfert_mod
184    USE omp_para
185
186    IMPLICIT NONE
187    REAL(rstd) :: x1(3),x2(3)
188    REAL(rstd) :: vect(3,6)
189    REAL(rstd) :: centr(3)
190    INTEGER :: ind,i,j,n,k
191    TYPE(t_message),SAVE :: message0, message1
192    LOGICAL, SAVE :: first=.TRUE.
193!$OMP THREADPRIVATE(first)   
194   
195    IF (first) THEN
196      CALL init_message(geom%xyz_i, req_i0 ,message0)
197      CALL init_message(geom%xyz_i, req_i1 ,message1)
198      first=.FALSE.
199    ENDIF
200   
201    CALL transfert_message(geom%xyz_i,message0)
202    CALL transfert_message(geom%xyz_i,message1)
203   
204    DO ind=1,ndomain
205      IF (.NOT. assigned_domain(ind) .OR. .NOT. is_omp_level_master) CYCLE
206      CALL swap_dimensions(ind)
207      CALL swap_geometry(ind)
208      DO j=jj_begin,jj_end
209        DO i=ii_begin,ii_end
210          n=(j-1)*iim+i
211          DO k=0,5
212            x1(:) = xyz_i(n+t_pos(k+1),:)
213            x2(:) = xyz_i(n+t_pos(MOD(k+1,6)+1),:)
214            if (norm(x1-x2)<1e-16) x2(:) = xyz_i(n+t_pos(MOD(k+2,6)+1),:)
215            CALL circumcenter(xyz_i(n,:), x1, x2, xyz_v(n+z_pos(k+1),:))
216          ENDDO
217        ENDDO
218      ENDDO
219    ENDDO
220
221  END SUBROUTINE update_circumcenters
222
223  SUBROUTINE remap_schmidt_loc
224    USE spherical_geom_mod
225    USE getin_mod
226    USE omp_para
227    USE domain_mod
228    USE dimensions
229    IMPLICIT NONE
230    INTEGER :: ind,i,j,n
231    REAL(rstd) :: schmidt_factor, schmidt_lon, schmidt_lat
232
233    ! Schmidt transform parameters
234    schmidt_factor = 1.
235    CALL getin('schmidt_factor', schmidt_factor)
236    schmidt_factor =  schmidt_factor**2.
237   
238    schmidt_lon = 0.
239    CALL getin('schmidt_lon', schmidt_lon)
240    schmidt_lon = schmidt_lon * pi/180.
241
242    schmidt_lat = 45.
243    CALL getin('schmidt_lat', schmidt_lat)
244    schmidt_lat = schmidt_lat * pi/180.
245
246    DO ind=1,ndomain
247      IF (.NOT. assigned_domain(ind) .OR. .NOT. is_omp_level_master) CYCLE
248      CALL swap_dimensions(ind)
249      CALL swap_geometry(ind)
250      DO j=jj_begin,jj_end
251        DO i=ii_begin,ii_end
252          n=(j-1)*iim+i
253          CALL schmidt_transform(xyz_i(n,:), schmidt_factor, schmidt_lon, schmidt_lat)
254        ENDDO
255      ENDDO
256    ENDDO
257  END SUBROUTINE remap_schmidt_loc
258
259  SUBROUTINE optimize_geometry
260  USE metric
261  USE spherical_geom_mod
262  USE domain_mod
263  USE dimensions
264  USE transfert_mod
265  USE vector
266  USE getin_mod
267  USE omp_para
268  IMPLICIT NONE
269    INTEGER :: nb_it=0
270    TYPE(t_domain),POINTER :: d
271    INTEGER :: ind,it,i,j,n,k
272    REAL(rstd) :: x1(3),x2(3)
273    REAL(rstd) :: vect(3,6)
274    REAL(rstd) :: centr(3)
275    REAL(rstd) :: sum
276    LOGICAL    :: check
277   
278   
279    CALL getin('optim_it',nb_it)
280   
281    DO ind=1,ndomain
282      IF (.NOT. assigned_domain(ind)  .OR. .NOT. is_omp_level_master) CYCLE
283      d=>domain(ind)
284      CALL swap_dimensions(ind)
285      CALL swap_geometry(ind)
286      xyz_i(:,1) = 0 ; xyz_i(:,2) = 0 ;  xyz_i(:,3) = 1
287     
288      DO j=jj_begin,jj_end
289        DO i=ii_begin,ii_end
290          n=(j-1)*iim+i
291          xyz_i(n,:)=d%xyz(:,i,j) 
292        ENDDO
293      ENDDO
294    ENDDO
295   
296    CALL update_circumcenters
297
298    DO ind=1,ndomain
299      IF (.NOT. assigned_domain(ind)  .OR. .NOT. is_omp_level_master ) CYCLE
300      d=>domain(ind)
301      CALL swap_dimensions(ind)
302      CALL swap_geometry(ind)
303      DO j=jj_begin,jj_end
304        DO i=ii_begin,ii_end
305          n=(j-1)*iim+i
306          DO k=0,5
307            x1(:) = xyz_v(n+z_pos(k+1),:)
308            x2(:) = d%vertex(:,k,i,j) 
309            IF (norm(x1-x2)>1e-10) THEN
310              PRINT*,"vertex diff ",ind,i,j,k
311              PRINT*,x1
312              PRINT*,x2
313            ENDIF
314          ENDDO
315        ENDDO
316      ENDDO
317    ENDDO
318   
319   
320    DO it=1,nb_it
321      IF (MOD(it,100)==0) THEN
322        check=is_master
323      ELSE
324        check=.FALSE.
325      ENDIF
326     
327      sum=0
328      DO ind=1,ndomain
329      IF (.NOT. assigned_domain(ind)  .OR. .NOT. is_omp_level_master ) CYCLE
330        CALL swap_dimensions(ind)
331        CALL swap_geometry(ind)
332        DO j=jj_begin,jj_end
333          DO i=ii_begin,ii_end
334            n=(j-1)*iim+i
335            vect(:,1)=xyz_v(n+z_rup,:)
336            vect(:,2)=xyz_v(n+z_up,:)
337            vect(:,3)=xyz_v(n+z_lup,:)
338            vect(:,4)=xyz_v(n+z_ldown,:)
339            vect(:,5)=xyz_v(n+z_down,:)
340            vect(:,6)=xyz_v(n+z_rdown,:)
341            CALL compute_centroid(vect,6,centr)
342            IF (check) THEN
343              sum=MAX(sum,norm(xyz_i(n,:)-centr(:)))
344            ENDIF
345            xyz_i(n,:)=centr(:)
346          ENDDO
347        ENDDO
348       
349      ENDDO
350
351      IF (check) PRINT *,"it = ",it,"  diff centroid circumcenter ",sum
352
353     CALL update_circumcenters
354
355    ENDDO
356   
357  END SUBROUTINE optimize_geometry
358
359  SUBROUTINE update_domain 
360    ! copy position of generators and vertices back into domain(:)%xyz/vertex
361    ! so that XIOS/create_header_gen gets the right values
362    USE omp_para
363    USE dimensions
364    USE domain_mod 
365    USE transfert_mpi_mod
366
367    INTEGER :: ind,i,j,k,n
368    TYPE(t_domain),POINTER :: d
369    TYPE(t_field),POINTER :: xyz_glo(:), xyz_loc(:), vertex_glo(:), vertex_loc(:)
370    REAL(rstd), POINTER :: xyz(:,:), vertex(:,:)
371   
372    CALL allocate_field(xyz_loc, field_t, type_real, 3)
373    CALL allocate_field(vertex_loc, field_z, type_real, 3)
374
375    DO ind=1,ndomain
376       IF (.NOT. assigned_domain(ind)  .OR. .NOT. is_omp_level_master ) CYCLE
377       CALL swap_dimensions(ind)
378       CALL swap_geometry(ind)
379       xyz = xyz_loc(ind)
380       xyz(:,:) = xyz_i(:,:)
381       vertex = vertex_loc(ind)
382       vertex(:,:) = xyz_v(:,:)
383    END DO
384   
385    CALL allocate_field_glo(xyz_glo, field_t, type_real, 3)
386    CALL allocate_field_glo(vertex_glo, field_z, type_real, 3)
387    CALL gather_field(xyz_loc, xyz_glo)
388    CALL bcast_field(xyz_glo)
389    CALL gather_field(vertex_loc, vertex_glo)
390    CALL bcast_field(vertex_glo)
391    CALL deallocate_field(vertex_loc)
392    CALL deallocate_field(xyz_loc)
393
394    DO ind=1,ndomain_glo
395       d=>domain_glo(ind)
396       xyz = xyz_glo(ind)
397       vertex = vertex_glo(ind)
398       DO j=d%jj_begin,d%jj_end
399          DO i=d%ii_begin,d%ii_end         
400             n=(j-1)*d%iim+i
401             d%xyz(:,i,j)=xyz(n,:)
402             DO k=0,5
403                d%vertex(:,k,i,j) = vertex(n+d%z_pos(k+1),:)
404             END DO
405          END DO
406       END DO
407    END DO
408
409    CALL deallocate_field_glo(vertex_glo)
410    CALL deallocate_field_glo(xyz_glo)
411
412  END SUBROUTINE update_domain
413 
414  SUBROUTINE set_geometry
415  USE metric
416  USE vector
417  USE spherical_geom_mod
418  USE domain_mod
419  USE dimensions
420  USE transfert_mod
421  USE getin_mod
422  USE omp_para
423  IMPLICIT NONE
424
425    REAL(rstd) :: surf(6)
426    REAL(rstd) :: surf_v(6)
427    REAL(rstd) :: vect(3,6)
428    REAL(rstd) :: centr(3)
429    REAL(rstd) :: vet(3),vep(3), vertex(3)
430    INTEGER :: ind,i,j,k,n
431    TYPE(t_domain),POINTER :: d
432    REAL(rstd) ::  S12
433    REAL(rstd) :: w(6)
434    REAL(rstd) :: lon,lat
435    INTEGER :: ii_glo,jj_glo
436    REAL(rstd) :: S
437         
438     
439    CALL optimize_geometry
440    CALL remap_schmidt_loc
441    CALL update_circumcenters
442    ! copy position of generators and vertices back into domain(:)%xyz/vertex
443    ! so that XIOS gets the right values
444    CALL update_domain
445
446    DO ind=1,ndomain
447      IF (.NOT. assigned_domain(ind)  .OR. .NOT. is_omp_level_master ) CYCLE
448      d=>domain(ind)
449      CALL swap_dimensions(ind)
450      CALL swap_geometry(ind)
451      lon_i(:)=0 ; lat_i(:)=0
452      lon_e(:)=0 ; lat_e(:)=0
453      DO j=jj_begin-1,jj_end+1
454        DO i=ii_begin-1,ii_end+1
455          n=(j-1)*iim+i
456
457          DO k=0,5
458            ne(n,k+1)=d%ne(k,i,j)
459          ENDDO
460
461          vect(:,1)=xyz_v(n+z_rup,:)
462          vect(:,2)=xyz_v(n+z_up,:)
463          vect(:,3)=xyz_v(n+z_lup,:)
464          vect(:,4)=xyz_v(n+z_ldown,:)
465          vect(:,5)=xyz_v(n+z_down,:)
466          vect(:,6)=xyz_v(n+z_rdown,:)
467          CALL compute_centroid(vect,6,centr)
468          centroid(n,:)=centr(:)
469
470               
471          CALL xyz2lonlat(xyz_v(n+z_up,:),lon,lat)
472          fv(n+z_up)=2*sin(lat)*omega
473          CALL xyz2lonlat(xyz_v(n+z_down,:),lon,lat)
474          fv(n+z_down)=2*sin(lat)*omega
475         
476          bi(n)=0. 
477       
478          CALL dist_cart(xyz_i(n,:),xyz_i(n+t_right,:),de(n+u_right))
479          CALL dist_cart(xyz_i(n,:),xyz_i(n+t_lup,:),de(n+u_lup))
480          CALL dist_cart(xyz_i(n,:),xyz_i(n+t_ldown,:),de(n+u_ldown))
481         
482          CALL div_arc_bis(xyz_i(n,:),xyz_i(n+t_right,:),0.5,xyz_e(n+u_right,:))
483          CALL div_arc_bis(xyz_i(n,:),xyz_i(n+t_lup,:),0.5,xyz_e(n+u_lup,:))
484          CALL div_arc_bis(xyz_i(n,:),xyz_i(n+t_ldown,:),0.5,xyz_e(n+u_ldown,:))
485
486          CALL dist_cart(xyz_v(n+z_rdown,:), xyz_v(n+z_rup,:),le(n+u_right))
487          CALL dist_cart(xyz_v(n+z_up,:), xyz_v(n+z_lup,:),le(n+u_lup))
488          CALL dist_cart(xyz_v(n+z_ldown,:), xyz_v(n+z_down,:),le(n+u_ldown))
489
490          le_de(n+u_right)=le(n+u_right)/de(n+u_right) ! NaN possible but should be harmless
491          le_de(n+u_lup)  =le(n+u_lup)  /de(n+u_lup)
492          le_de(n+u_ldown)=le(n+u_ldown)/de(n+u_ldown)
493
494          Ai(n)=0
495          DO k=0,5
496            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))
497            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))
498            Ai(n)=Ai(n)+surf(k+1)
499            IF (i==ii_end .AND. j==jj_begin) THEN
500              IF (Ai(n)<1e20) THEN
501              ELSE
502                PRINT *,"PB !!",Ai(n),k,surf(k+1)
503                PRINT*,xyz_i(n,:),xyz_v(n+z_pos(MOD((k-1+6),6)+1),:),xyz_v(n+z_pos(k+1),:)
504              ENDIF
505            ENDIF
506          ENDDO
507
508          ! Sign convention : Ringler et al., JCP 2010, eq. 21 p. 3071
509          ! Normal component is along outgoing normal vector if ne=1
510
511          CALL cross_product2(xyz_v(n+z_rdown,:),xyz_v(n+z_rup,:),vep)
512          IF (norm(vep)>1e-30) THEN
513            vep(:)=vep(:)/norm(vep)                         ! Inward normal vector
514            CALL cross_product2(vep,xyz_e(n+u_right,:),vet) ! Counter-clockwise tangent vector
515            vet(:)=vet(:)/norm(vet)
516            ep_e(n+u_right,:)=-vep(:)*ne(n,right)
517            et_e(n+u_right,:)=vet(:)*ne(n,right)
518          ENDIF
519
520          CALL cross_product2(xyz_v(n+z_up,:),xyz_v(n+z_lup,:),vep)
521          IF (norm(vep)>1e-30) THEN
522            vep(:)=vep(:)/norm(vep)
523            CALL cross_product2(vep,xyz_e(n+u_lup,:),vet)
524            vet(:)=vet(:)/norm(vet)
525            ep_e(n+u_lup,:)=-vep(:)*ne(n,lup)
526            et_e(n+u_lup,:)=vet(:)*ne(n,lup)
527          ENDIF
528
529          CALL cross_product2(xyz_v(n+z_ldown,:),xyz_v(n+z_down,:),vep)
530          IF (norm(vep)>1e-30) THEN
531            vep(:)=vep(:)/norm(vep)
532            CALL cross_product2(vep,xyz_e(n+u_ldown,:),vet)
533            vet(:)=vet(:)/norm(vet)
534            ep_e(n+u_ldown,:)=-vep(:)*ne(n,ldown)
535            et_e(n+u_ldown,:)=vet(:)*ne(n,ldown)
536          ENDIF
537
538          CALL xyz2lonlat(xyz_i(n,:),lon,lat)
539          lon_i(n)=lon
540          lat_i(n)=lat
541          elon_i(n,1) = -sin(lon)
542          elon_i(n,2) = cos(lon)
543          elon_i(n,3) = 0
544          elat_i(n,1) = -cos(lon)*sin(lat)
545          elat_i(n,2) = -sin(lon)*sin(lat)
546          elat_i(n,3) = cos(lat)
547
548         
549          CALL xyz2lonlat(xyz_e(n+u_right,:),lon,lat)
550          lon_e(n+u_right)=lon
551          lat_e(n+u_right)=lat
552          elon_e(n+u_right,1) = -sin(lon)
553          elon_e(n+u_right,2) = cos(lon)
554          elon_e(n+u_right,3) = 0
555          elat_e(n+u_right,1) = -cos(lon)*sin(lat)
556          elat_e(n+u_right,2) = -sin(lon)*sin(lat)
557          elat_e(n+u_right,3) = cos(lat)
558         
559          CALL xyz2lonlat(xyz_e(n+u_lup,:),lon,lat)
560          lon_e(n+u_lup)=lon
561          lat_e(n+u_lup)=lat
562          elon_e(n+u_lup,1) = -sin(lon)
563          elon_e(n+u_lup,2) = cos(lon)
564          elon_e(n+u_lup,3) = 0
565          elat_e(n+u_lup,1) = -cos(lon)*sin(lat)
566          elat_e(n+u_lup,2) = -sin(lon)*sin(lat)
567          elat_e(n+u_lup,3) = cos(lat)
568 
569          CALL xyz2lonlat(xyz_e(n+u_ldown,:),lon,lat)
570          lon_e(n+u_ldown)=lon
571          lat_e(n+u_ldown)=lat
572          elon_e(n+u_ldown,1) = -sin(lon)
573          elon_e(n+u_ldown,2) = cos(lon)
574          elon_e(n+u_ldown,3) = 0
575          elat_e(n+u_ldown,1) = -cos(lon)*sin(lat)
576          elat_e(n+u_ldown,2) = -sin(lon)*sin(lat)
577          elat_e(n+u_ldown,3) = cos(lat)
578 
579         
580          DO k=0,5
581            CALL surf_triangle(xyz_i(n,:), xyz_v(n+z_pos(k+1),:), xyz_i(n+t_pos(k+1),:),S1(n,k+1) )
582            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(n,k+1) )
583            S12 = .5*(S1(n,k+1)+S2(n,k+1))
584            Riv(n,k+1)=S12/Ai(n)
585            Riv2(n,k+1)=S12/surf_v(k+1)
586          ENDDO
587
588          DO k=1,6
589            IF (ABS(surf_v(k))<1e-30) THEN
590              Riv(n,k)=0.
591            ENDIF
592          ENDDO
593         
594          Av(n+z_up)=surf_v(vup)+1e-100
595          Av(n+z_down)=surf_v(vdown)+1e-100
596       
597        ENDDO
598      ENDDO
599     
600      DO j=jj_begin,jj_end
601        DO i=ii_begin,ii_end
602          n=(j-1)*iim+i
603   
604          CALL compute_wee(n,right,w)
605          Wee(n+u_right,:,1)=w(1:5)
606
607          CALL compute_wee(n+t_right,left,w)
608          Wee(n+u_right,:,2)=w(1:5)
609
610
611          CALL compute_wee(n,lup,w)
612          Wee(n+u_lup,:,1)=w(1:5)
613
614          CALL compute_wee(n+t_lup,rdown,w)
615          Wee(n+u_lup,:,2)=w(1:5)
616
617
618          CALL compute_wee(n,ldown,w)
619          Wee(n+u_ldown,:,1)=w(1:5)
620
621          CALL compute_wee(n+t_ldown,rup,w)
622          Wee(n+u_ldown,:,2)=w(1:5)
623
624        ENDDO
625      ENDDO
626     
627      DO j=jj_begin,jj_end
628        DO i=ii_begin,ii_end
629          n=(j-1)*iim+i
630          ii_glo=d%ii_begin_glo-d%ii_begin+i
631          jj_glo=d%jj_begin_glo-d%jj_begin+j
632         
633          IF (ii_glo==1 .AND. jj_glo==1) THEN
634            le(n+u_ldown)=0
635            xyz_v(n+z_ldown,:)=xyz_v(n+z_down,:)
636                       
637          ENDIF
638
639          IF (ii_glo==iim_glo .AND. jj_glo==1) THEN
640            le(n+u_right)=0
641            xyz_v(n+z_rdown,:)=xyz_v(n+z_rup,:)
642          ENDIF
643
644          IF (ii_glo==iim_glo .AND. jj_glo==jjm_glo) THEN
645            le(n+u_rup)=0
646            xyz_v(n+z_rup,:)=xyz_v(n+z_up,:)
647          ENDIF
648
649          IF (ii_glo==1 .AND. jj_glo==jjm_glo) THEN
650            le(n+u_lup)=0
651            xyz_v(n+z_up,:)=xyz_v(n+z_lup,:)
652          ENDIF
653         
654        ENDDO
655      ENDDO
656
657      DO j=jj_begin-1,jj_end+1
658        DO i=ii_begin-1,ii_end+1
659          n=(j-1)*iim+i
660          xyz_i(n,:)=xyz_i(n,:) * radius
661          xyz_v(n+z_up,:)=xyz_v(n+z_up,:) * radius
662          xyz_v(n+z_down,:)=xyz_v(n+z_down,:) *radius
663          de(n+u_right)=de(n+u_right) * radius
664          de(n+u_lup)=de(n+u_lup)*radius
665          de(n+u_ldown)=de(n+u_ldown)*radius
666          xyz_e(n+u_right,:)=xyz_e(n+u_right,:)*radius
667          xyz_e(n+u_lup,:)=xyz_e(n+u_lup,:)*radius
668          xyz_e(n+u_ldown,:)=xyz_e(n+u_ldown,:)*radius
669          le(n+u_right)=le(n+u_right)*radius
670          le(n+u_lup)=le(n+u_lup)*radius
671          le(n+u_ldown)=le(n+u_ldown)*radius
672          Ai(n)=Ai(n)*radius**2
673          Av(n+z_up)=Av(n+z_up)*radius**2
674          Av(n+z_down)=Av(n+z_down)*radius**2
675        ENDDO
676      ENDDO
677                 
678    ENDDO
679
680    CALL transfert_request(geom%Ai,req_i1)
681    CALL transfert_request(geom%centroid,req_i1)
682
683!    CALL surf_triangle(d%xyz(:,ii_begin,jj_begin),d%xyz(:,ii_begin,jj_end),d%xyz(:,ii_end,jj_begin),S)
684 
685  END SUBROUTINE set_geometry
686 
687  SUBROUTINE compute_wee(n,pos,w)
688  IMPLICIT NONE
689    INTEGER,INTENT(IN) :: n
690    INTEGER,INTENT(IN) :: pos
691    REAL(rstd),INTENT(OUT) ::w(6)
692
693    REAL(rstd) :: ne_(0:5)
694    REAL(rstd) :: Riv_(6)
695    INTEGER    :: k
696   
697   
698    DO k=0,5
699      ne_(k)=ne(n,MOD(pos-1+k+6,6)+1) 
700      Riv_(k+1)=Riv(n,MOD(pos-1+k+6,6)+1)
701    ENDDO
702         
703    w(1)=-ne_(0)*ne_(1)*(Riv_(1)-0.5)
704    w(2)=-ne_(2)*(ne_(0)*Riv_(2)-w(1)*ne_(1))
705    w(3)=-ne_(3)*(ne_(0)*Riv_(3)-w(2)*ne_(2))
706    w(4)=-ne_(4)*(ne_(0)*Riv_(4)-w(3)*ne_(3))
707    w(5)=-ne_(5)*(ne_(0)*Riv_(5)-w(4)*ne_(4))
708    w(6)=ne_(0)*ne_(5)*(Riv_(6)-0.5)
709   
710!    IF ( ABS(w(5)-w(6))>1e-20) PRINT *, "pb pour wee : w(5)!=w(6)",sum(Riv_(:))
711
712   END SUBROUTINE compute_wee
713           
714
715 
716  SUBROUTINE compute_geometry
717  IMPLICIT NONE
718    CALL allocate_geometry
719    CALL set_geometry
720   
721  END SUBROUTINE compute_geometry
722 
723END MODULE geometry
Note: See TracBrowser for help on using the repository browser.