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

Last change on this file since 47 was 47, checked in by dubos, 12 years ago

Fixed sign error in geometry.f90

File size: 17.9 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  USE ioipsl
118  IMPLICIT NONE
119    INTEGER :: nb_it=0
120    TYPE(t_domain),POINTER :: d
121    INTEGER :: ind,it,i,j,n,k
122    REAL(rstd) :: x1(3),x2(3)
123    REAL(rstd) :: vect(3,6)
124    REAL(rstd) :: centr(3)
125    REAL(rstd) :: sum
126    LOGICAL    :: check
127   
128   
129    CALL getin('optim_it',nb_it)
130   
131    DO ind=1,ndomain
132      d=>domain(ind)
133      CALL swap_dimensions(ind)
134      CALL swap_geometry(ind)
135      DO j=jj_begin,jj_end
136        DO i=ii_begin,ii_end
137          n=(j-1)*iim+i
138          xyz_i(n,:)=d%xyz(:,i,j) 
139        ENDDO
140      ENDDO
141    ENDDO
142   
143    CALL transfert_request(geom%xyz_i,req_i1)
144   
145    DO ind=1,ndomain
146      CALL swap_dimensions(ind)
147      CALL swap_geometry(ind)
148      DO j=jj_begin,jj_end
149        DO i=ii_begin,ii_end
150          n=(j-1)*iim+i
151          DO k=0,5
152            x1(:) = xyz_i(n+t_pos(k+1),:)
153            x2(:) = xyz_i(n+t_pos(MOD(k+1,6)+1),:)
154            if (norm(x1-x2)<1e-16) x2(:) = xyz_i(n+t_pos(MOD(k+2,6)+1),:)
155            CALL circumcenter(xyz_i(n,:), x1, x2, xyz_v(n+z_pos(k+1),:))
156          ENDDO
157        ENDDO
158      ENDDO
159    ENDDO
160   
161    DO ind=1,ndomain
162      d=>domain(ind)
163      CALL swap_dimensions(ind)
164      CALL swap_geometry(ind)
165      DO j=jj_begin,jj_end
166        DO i=ii_begin,ii_end
167          n=(j-1)*iim+i
168          DO k=0,5
169            x1(:) = xyz_v(n+z_pos(k+1),:)
170            x2(:) = d%vertex(:,k,i,j) 
171            IF (norm(x1-x2)>1e-10) THEN
172              PRINT*,"vertex diff ",ind,i,j,k
173              PRINT*,x1
174              PRINT*,x2
175            ENDIF
176          ENDDO
177        ENDDO
178      ENDDO
179    ENDDO
180   
181   
182    DO it=1,nb_it
183      IF (MOD(it,100)==0) THEN
184        check=.TRUE.
185      ELSE
186        check=.FALSE.
187      ENDIF
188     
189      sum=0
190      DO ind=1,ndomain
191        CALL swap_dimensions(ind)
192        CALL swap_geometry(ind)
193        DO j=jj_begin,jj_end
194          DO i=ii_begin,ii_end
195            n=(j-1)*iim+i
196            vect(:,1)=xyz_v(n+z_rup,:)
197            vect(:,2)=xyz_v(n+z_up,:)
198            vect(:,3)=xyz_v(n+z_lup,:)
199            vect(:,4)=xyz_v(n+z_ldown,:)
200            vect(:,5)=xyz_v(n+z_down,:)
201            vect(:,6)=xyz_v(n+z_rdown,:)
202            CALL compute_centroid(vect,6,centr)
203            IF (check) THEN
204              sum=MAX(sum,norm(xyz_i(n,:)-centr(:)))
205            ENDIF
206            xyz_i(n,:)=centr(:)
207          ENDDO
208        ENDDO
209!        i=ii_begin ; j=jj_begin ; n=(j-1)*iim+i ; xyz_i(n,:)=domain(ind)%xyz(:,i,j)
210!        i=ii_begin ; j=jj_end   ; n=(j-1)*iim+i ; xyz_i(n,:)=domain(ind)%xyz(:,i,j)
211!        i=ii_end   ; j=jj_begin ; n=(j-1)*iim+i ;   xyz_i(n,:)=domain(ind)%xyz(:,i,j)
212!        PRINT *,"Pb ?? : "
213!        PRINT *, xyz_i(n,:), domain(ind)%xyz(:,i,j), norm(xyz_i(n,:)- domain(ind)%xyz(:,i,j))
214!        i=ii_end   ; j=jj_end   ; n=(j-1)*iim+i ; xyz_i(n,:)=domain(ind)%xyz(:,i,j)
215       
216      ENDDO
217
218      IF (check) THEN
219!        sum=sum/(ndomain*ii_nb*jj_nb)
220        PRINT *,"it = ",it,"  diff centroid circumcenter ",sum
221      ENDIF
222     
223     CALL transfert_request(geom%xyz_i,req_i1)
224
225    DO ind=1,ndomain
226      CALL swap_dimensions(ind)
227      CALL swap_geometry(ind)
228      DO j=jj_begin,jj_end
229        DO i=ii_begin,ii_end
230          n=(j-1)*iim+i
231          DO k=0,5
232            x1(:) = xyz_i(n+t_pos(k+1),:)
233            x2(:) = xyz_i(n+t_pos(MOD(k+1,6)+1),:)
234            if (norm(x1-x2)<1e-16) x2(:) = xyz_i(n+t_pos(MOD(k+2,6)+1),:)
235            CALL circumcenter(xyz_i(n,:), x1, x2, xyz_v(n+z_pos(k+1),:))
236          ENDDO
237        ENDDO
238      ENDDO
239    ENDDO
240
241     
242    ENDDO
243   
244  END SUBROUTINE optimize_geometry
245 
246   
247  SUBROUTINE set_geometry
248  USE metric
249  USE vector
250  USE spherical_geom_mod
251  USE domain_mod
252  USE dimensions
253  USE transfert_mod
254  IMPLICIT NONE
255
256    REAL(rstd) :: surf(6)
257    REAL(rstd) :: surf_v(6)
258    REAL(rstd) :: vect(3,6)
259    REAL(rstd) :: centr(3)
260    REAL(rstd) :: vet(3),vep(3)
261    INTEGER :: ind,i,j,k,n
262    TYPE(t_domain),POINTER :: d
263    REAL(rstd) ::  S
264    REAL(rstd) :: w(6)
265    REAL(rstd) :: lon,lat
266    INTEGER :: ii_glo,jj_glo
267    REAL(rstd) :: S1,S2
268         
269     
270    CALL optimize_geometry
271   
272    DO ind=1,ndomain
273      d=>domain(ind)
274      CALL swap_dimensions(ind)
275      CALL swap_geometry(ind)
276      DO j=jj_begin-1,jj_end+1
277        DO i=ii_begin-1,ii_end+1
278          n=(j-1)*iim+i
279
280          DO k=0,5
281            ne(n,k+1)=d%ne(k,i,j)
282          ENDDO
283                 
284!          xyz_i(n,:)=d%xyz(:,i,j)
285!          xyz_v(n+z_up,:)=d%vertex(:,vup-1,i,j)
286!          xyz_v(n+z_down,:)=d%vertex(:,vdown-1,i,j)
287
288          vect(:,1)=xyz_v(n+z_rup,:)
289          vect(:,2)=xyz_v(n+z_up,:)
290          vect(:,3)=xyz_v(n+z_lup,:)
291          vect(:,4)=xyz_v(n+z_ldown,:)
292          vect(:,5)=xyz_v(n+z_down,:)
293          vect(:,6)=xyz_v(n+z_rdown,:)
294          CALL compute_centroid(vect,6,centr)
295          centroid(n,:)=centr(:)
296
297               
298          CALL xyz2lonlat(xyz_v(n+z_up,:),lon,lat)
299          fv(n+z_up)=2*sin(lat)*omega
300          CALL xyz2lonlat(xyz_v(n+z_down,:),lon,lat)
301          fv(n+z_down)=2*sin(lat)*omega
302         
303          bi(n)=0. 
304       
305!          CALL dist_cart(d%xyz(:,i,j),d%neighbour(:,right-1,i,j),de(n+u_right))
306!          CALL dist_cart(d%xyz(:,i,j),d%neighbour(:,lup-1,i,j),de(n+u_lup))
307!          CALL dist_cart(d%xyz(:,i,j),d%neighbour(:,ldown-1,i,j),de(n+u_ldown))
308
309          CALL dist_cart(xyz_i(n,:),xyz_i(n+t_right,:),de(n+u_right))
310          CALL dist_cart(xyz_i(n,:),xyz_i(n+t_lup,:),de(n+u_lup))
311          CALL dist_cart(xyz_i(n,:),xyz_i(n+t_ldown,:),de(n+u_ldown))
312         
313!          CALL div_arc_bis(d%xyz(:,i,j),d%neighbour(:,right-1,i,j),0.5,xyz_e(n+u_right,:))
314!          CALL div_arc_bis(d%xyz(:,i,j),d%neighbour(:,lup-1,i,j),0.5,xyz_e(n+u_lup,:))
315!          CALL div_arc_bis(d%xyz(:,i,j),d%neighbour(:,ldown-1,i,j),0.5,xyz_e(n+u_ldown,:))
316
317          CALL div_arc_bis(xyz_i(n,:),xyz_i(n+t_right,:),0.5,xyz_e(n+u_right,:))
318          CALL div_arc_bis(xyz_i(n,:),xyz_i(n+t_lup,:),0.5,xyz_e(n+u_lup,:))
319          CALL div_arc_bis(xyz_i(n,:),xyz_i(n+t_ldown,:),0.5,xyz_e(n+u_ldown,:))
320         
321!          CALL dist_cart(d%vertex(:,vrdown-1,i,j),d%vertex(:,vrup-1,i,j),le(n+u_right))
322!          CALL dist_cart(d%vertex(:,vup-1,i,j),d%vertex(:,vlup-1,i,j),le(n+u_lup))
323!          CALL dist_cart(d%vertex(:,vldown-1,i,j),d%vertex(:,vdown-1,i,j),le(n+u_ldown))
324
325          CALL dist_cart(xyz_v(n+z_rdown,:), xyz_v(n+z_rup,:),le(n+u_right))
326          CALL dist_cart(xyz_v(n+z_up,:), xyz_v(n+z_lup,:),le(n+u_lup))
327          CALL dist_cart(xyz_v(n+z_ldown,:), xyz_v(n+z_down,:),le(n+u_ldown))
328         
329          Ai(n)=0
330          DO k=0,5
331!            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))
332!            CALL surf_triangle(d%xyz(:,i,j),d%vertex(:,MOD((k-1+6),6),i,j),d%vertex(:,k,i,j),surf(k+1))
333            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))
334            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))
335            Ai(n)=Ai(n)+surf(k+1)
336            IF (i==ii_end .AND. j==jj_begin) THEN
337              IF (Ai(n)<1e20) THEN
338              ELSE
339                PRINT *,"PB !!",Ai(n),k,surf(k+1)
340                PRINT*,xyz_i(n,:),xyz_v(n+z_pos(MOD((k-1+6),6)+1),:),xyz_v(n+z_pos(k+1),:)
341              ENDIF
342            ENDIF
343          ENDDO
344
345          ! Sign convention : Ringler et al., JCP 2010, eq. 21 p. 3071
346          ! Normal component is along outgoing normal vector if ne=1
347
348          CALL cross_product2(xyz_v(n+z_rdown,:),xyz_v(n+z_rup,:),vep)
349          IF (norm(vep)>1e-30) THEN
350            vep(:)=vep(:)/norm(vep)                         ! Inward normal vector
351            CALL cross_product2(vep,xyz_e(n+u_right,:),vet) ! Counter-clockwise tangent vector
352            vet(:)=vet(:)/norm(vet)
353            ep_e(n+u_right,:)=-vep(:)*ne(n,right)
354            et_e(n+u_right,:)=vet(:)*ne(n,right)
355          ENDIF
356
357          CALL cross_product2(xyz_v(n+z_up,:),xyz_v(n+z_lup,:),vep)
358          IF (norm(vep)>1e-30) THEN
359            vep(:)=vep(:)/norm(vep)
360            CALL cross_product2(vep,xyz_e(n+u_lup,:),vet)
361            vet(:)=vet(:)/norm(vet)
362            ep_e(n+u_lup,:)=-vep(:)*ne(n,lup)
363            et_e(n+u_lup,:)=vet(:)*ne(n,lup)
364          ENDIF
365
366          CALL cross_product2(xyz_v(n+z_ldown,:),xyz_v(n+z_down,:),vep)
367          IF (norm(vep)>1e-30) THEN
368            vep(:)=vep(:)/norm(vep)
369            CALL cross_product2(vep,xyz_e(n+u_ldown,:),vet)
370            vet(:)=vet(:)/norm(vet)
371            ep_e(n+u_ldown,:)=-vep(:)*ne(n,ldown)
372            et_e(n+u_ldown,:)=vet(:)*ne(n,ldown)
373          ENDIF
374
375          CALL xyz2lonlat(xyz_i(n,:),lon,lat)
376          elon_i(n,1) = -sin(lon)
377          elon_i(n,2) = cos(lon)
378          elon_i(n,3) = 0
379          elat_i(n,1) = -cos(lon)*sin(lat)
380          elat_i(n,2) = -sin(lon)*sin(lat)
381          elat_i(n,3) = cos(lat)
382
383         
384          CALL xyz2lonlat(xyz_e(n+u_right,:),lon,lat)
385          elon_e(n+u_right,1) = -sin(lon)
386          elon_e(n+u_right,2) = cos(lon)
387          elon_e(n+u_right,3) = 0
388          elat_e(n+u_right,1) = -cos(lon)*sin(lat)
389          elat_e(n+u_right,2) = -sin(lon)*sin(lat)
390          elat_e(n+u_right,3) = cos(lat)
391         
392          CALL xyz2lonlat(xyz_e(n+u_lup,:),lon,lat)
393          elon_e(n+u_lup,1) = -sin(lon)
394          elon_e(n+u_lup,2) = cos(lon)
395          elon_e(n+u_lup,3) = 0
396          elat_e(n+u_lup,1) = -cos(lon)*sin(lat)
397          elat_e(n+u_lup,2) = -sin(lon)*sin(lat)
398          elat_e(n+u_lup,3) = cos(lat)
399 
400          CALL xyz2lonlat(xyz_e(n+u_ldown,:),lon,lat)
401          elon_e(n+u_ldown,1) = -sin(lon)
402          elon_e(n+u_ldown,2) = cos(lon)
403          elon_e(n+u_ldown,3) = 0
404          elat_e(n+u_ldown,1) = -cos(lon)*sin(lat)
405          elat_e(n+u_ldown,2) = -sin(lon)*sin(lat)
406          elat_e(n+u_ldown,3) = cos(lat)
407 
408         
409          DO k=0,5
410!            CALL surf_triangle(d%xyz(:,i,j),d%vertex(:,k,i,j),d%neighbour(:,k,i,j),S1)
411!            CALL surf_triangle(d%xyz(:,i,j),d%vertex(:,k,i,j),d%neighbour(:,MOD(k+1+6,6),i,j),S2)
412            CALL surf_triangle(xyz_i(n,:), xyz_v(n+z_pos(k+1),:), xyz_i(n+t_pos(k+1),:),S1)
413            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)
414            Riv(n,k+1)=0.5*(S1+S2)/Ai(n)
415            Riv2(n,k+1)=0.5*(S1+S2)/surf_v(k+1)
416          ENDDO
417
418          DO k=1,6
419            IF (ABS(surf_v(k))<1e-30) THEN
420              Riv(n,k)=0.
421            ENDIF
422          ENDDO
423         
424          Av(n+z_up)=surf_v(vup)+1e-100
425          Av(n+z_down)=surf_v(vdown)+1e-100
426       
427        ENDDO
428      ENDDO
429     
430      DO j=jj_begin,jj_end
431        DO i=ii_begin,ii_end
432          n=(j-1)*iim+i
433   
434          CALL compute_wee(n,right,w)
435          Wee(n+u_right,:,1)=w(1:5)
436
437          CALL compute_wee(n+t_right,left,w)
438          Wee(n+u_right,:,2)=w(1:5)
439
440
441          CALL compute_wee(n,lup,w)
442          Wee(n+u_lup,:,1)=w(1:5)
443
444          CALL compute_wee(n+t_lup,rdown,w)
445          Wee(n+u_lup,:,2)=w(1:5)
446
447
448          CALL compute_wee(n,ldown,w)
449          Wee(n+u_ldown,:,1)=w(1:5)
450
451          CALL compute_wee(n+t_ldown,rup,w)
452          Wee(n+u_ldown,:,2)=w(1:5)
453
454        ENDDO
455      ENDDO
456     
457      DO j=jj_begin,jj_end
458        DO i=ii_begin,ii_end
459          n=(j-1)*iim+i
460          ii_glo=d%ii_begin_glo-d%ii_begin+i
461          jj_glo=d%jj_begin_glo-d%jj_begin+j
462         
463          IF (ii_glo==1 .AND. jj_glo==1) THEN
464            le(n+u_ldown)=0
465            xyz_v(n+z_ldown,:)=xyz_v(n+z_down,:)
466                       
467          ENDIF
468
469          IF (ii_glo==iim_glo .AND. jj_glo==1) THEN
470            le(n+u_right)=0
471            xyz_v(n+z_rdown,:)=xyz_v(n+z_rup,:)
472          ENDIF
473
474          IF (ii_glo==iim_glo .AND. jj_glo==jjm_glo) THEN
475            le(n+u_rup)=0
476            xyz_v(n+z_rup,:)=xyz_v(n+z_up,:)
477          ENDIF
478
479          IF (ii_glo==1 .AND. jj_glo==jjm_glo) THEN
480            le(n+u_lup)=0
481            xyz_v(n+z_up,:)=xyz_v(n+z_lup,:)
482          ENDIF
483         
484        ENDDO
485      ENDDO
486
487      DO j=jj_begin-1,jj_end+1
488        DO i=ii_begin-1,ii_end+1
489          n=(j-1)*iim+i
490          xyz_i(n,:)=xyz_i(n,:) * radius
491          xyz_v(n+z_up,:)=xyz_v(n+z_up,:) * radius
492          xyz_v(n+z_down,:)=xyz_v(n+z_down,:) *radius
493          de(n+u_right)=de(n+u_right) * radius
494          de(n+u_lup)=de(n+u_lup)*radius
495          de(n+u_ldown)=de(n+u_ldown)*radius
496          xyz_e(n+u_right,:)=xyz_e(n+u_right,:)*radius
497          xyz_e(n+u_lup,:)=xyz_e(n+u_lup,:)*radius
498          xyz_e(n+u_ldown,:)=xyz_e(n+u_ldown,:)*radius
499          le(n+u_right)=le(n+u_right)*radius
500          le(n+u_lup)=le(n+u_lup)*radius
501          le(n+u_ldown)=le(n+u_ldown)*radius
502          Ai(n)=Ai(n)*radius**2
503          Av(n+z_up)=Av(n+z_up)*radius**2
504          Av(n+z_down)=Av(n+z_down)*radius**2
505        ENDDO
506      ENDDO
507                 
508    ENDDO
509    CALL transfert_request(geom%Ai,req_i1)
510    CALL transfert_request(geom%centroid,req_i1)
511    CALL surf_triangle(d%xyz(:,ii_begin,jj_begin),d%xyz(:,ii_begin,jj_end),d%xyz(:,ii_end,jj_begin),S)
512!    PRINT *,"Surf triangle : ",S*20/(4*Pi)
513 
514  END SUBROUTINE set_geometry
515 
516  SUBROUTINE compute_wee(n,pos,w)
517  IMPLICIT NONE
518    INTEGER,INTENT(IN) :: n
519    INTEGER,INTENT(IN) :: pos
520    REAL(rstd),INTENT(OUT) ::w(6)
521
522    REAL(rstd) :: ne_(0:5)
523    REAL(rstd) :: Riv_(6)
524    INTEGER    :: k
525   
526   
527    DO k=0,5
528      ne_(k)=ne(n,MOD(pos-1+k+6,6)+1) 
529      Riv_(k+1)=Riv(n,MOD(pos-1+k+6,6)+1)
530    ENDDO
531         
532    w(1)=-ne_(0)*ne_(1)*(Riv_(1)-0.5)
533    w(2)=-ne_(2)*(ne_(0)*Riv_(2)-w(1)*ne_(1))
534    w(3)=-ne_(3)*(ne_(0)*Riv_(3)-w(2)*ne_(2))
535    w(4)=-ne_(4)*(ne_(0)*Riv_(4)-w(3)*ne_(3))
536    w(5)=-ne_(5)*(ne_(0)*Riv_(5)-w(4)*ne_(4))
537    w(6)=ne_(0)*ne_(5)*(Riv_(6)-0.5)
538   
539!    IF ( ABS(w(5)-w(6))>1e-20) PRINT *, "pb pour wee : w(5)!=w(6)",sum(Riv_(:))
540
541   END SUBROUTINE compute_wee
542           
543
544 
545  SUBROUTINE compute_geometry
546  IMPLICIT NONE
547    CALL allocate_geometry
548    CALL set_geometry
549   
550  END SUBROUTINE compute_geometry
551 
552END MODULE geometry
Note: See TracBrowser for help on using the repository browser.