Changeset 73 for trunk/SOURCES/resol_adv_diff_2D-sept2009.f90
- Timestamp:
- 06/22/16 15:43:40 (8 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/SOURCES/resol_adv_diff_2D-sept2009.f90
r65 r73 101 101 subroutine resol_adv_diff_2D_vect(Dfx,Dfy,advx,advy,H_presc,i_Hpresc,bil,vieuxH,newH) 102 102 103 !$ USE OMP_LIB 103 104 implicit none 104 105 … … 148 149 if (itracebug.eq.1) call tracebug(' Entree dans routine resolution_diffusion') 149 150 150 151 !$OMP PARALLEL 152 !$OMP WORKSHARE 151 153 ! calcul de bdx et hdx 152 154 hdx(:,:)=dx1*(vieuxH(:,:)-eoshift(vieuxH(:,:),shift=-1,boundary=0.,dim=1)) … … 164 166 Frelax(:,:)=0. 165 167 DeltaH(:,:)=0. 166 168 !$OMP END WORKSHARE 167 169 168 170 ! schema spatial 169 171 170 172 if (upwind.eq.1.) then !schema amont 171 173 !$OMP WORKSHARE 172 174 where (Advx(:,:).ge.0.) 173 175 c_west(:,:)=1. … … 185 187 c_north(:,:)=1. 186 188 end where 187 189 !$OMP END WORKSHARE 188 190 else if (upwind.lt.1.) then ! schema centre 191 !$OMP WORKSHARE 189 192 c_west(:,:)=0.5 190 193 c_east(:,:)=0.5 191 194 c_south(:,:)=0.5 192 195 c_north(:,:)=0.5 196 !$OMP END WORKSHARE 193 197 end if 194 198 199 195 200 ! attribution des elements des diagonales 201 !$OMP DO PRIVATE(frdx,frdy,fraxw,fraxe,frays,frayn) 196 202 do j=2,ny-1 197 203 do i=2,nx-1 … … 254 260 end do 255 261 end do 256 257 262 !$OMP END DO 263 264 !$OMP WORKSHARE 258 265 where (i_hpresc(:,:) .eq.1) ! thickness prescribed 259 266 frelax(:,:) = H_presc(:,:) … … 264 271 erelax(:,:) = 0. 265 272 end where 273 !$OMP END WORKSHARE 274 !$OMP END PARALLEL 266 275 267 276 stopp = .false. 268 277 ntour=0 269 278 270 271 272 279 relax_loop: do while(.not.stopp) 273 280 ntour=ntour+1 274 275 281 !$OMP PARALLEL 282 !$OMP DO PRIVATE(reste) 276 283 do j=2,ny-1 277 284 do i=2,nx-1 … … 281 288 + crelax(i,j)*newH(i,j))- frelax(i,j) 282 289 283 if (ntour.eq.1) debug_3D(i,j,49)=(((arelax(i,j)*newH(i-1,j) +drelax(i,j)*newH(i,j-1)) &284 + (brelax(i,j)*newH(i+1,j) + erelax(i,j)*newH(i,j+1))) &285 + crelax(i,j)*newH(i,j))290 !if (ntour.eq.1) debug_3D(i,j,49)=(((arelax(i,j)*newH(i-1,j) +drelax(i,j)*newH(i,j-1)) & 291 ! + (brelax(i,j)*newH(i+1,j) + erelax(i,j)*newH(i,j+1))) & 292 ! + crelax(i,j)*newH(i,j)) 286 293 287 294 … … 290 297 end do 291 298 end do 292 293 debug_3D(:,:,50)=arelax(:,:) 294 debug_3D(:,:,51)=brelax(:,:) 295 debug_3D(:,:,52)=crelax(:,:) 296 debug_3D(:,:,53)=drelax(:,:) 297 debug_3D(:,:,54)=erelax(:,:) 298 debug_3D(:,:,55)=frelax(:,:) 299 300 301 299 !$OMP END DO 300 301 !debug_3D(:,:,50)=arelax(:,:) 302 !debug_3D(:,:,51)=brelax(:,:) 303 !debug_3D(:,:,52)=crelax(:,:) 304 !debug_3D(:,:,53)=drelax(:,:) 305 !debug_3D(:,:,54)=erelax(:,:) 306 !debug_3D(:,:,55)=frelax(:,:) 307 308 309 !$OMP WORKSHARE 302 310 newH(:,:)=newH(:,:)-deltaH(:,:) 303 304 311 !$OMP END WORKSHARE 312 !$OMP END PARALLEL 305 313 306 314 … … 310 318 delh=0 311 319 312 320 !$OMP PARALLEL 321 !$OMP DO REDUCTION(+:delh) 313 322 do j=2,ny-1 314 323 do i=2,nx-1 … … 316 325 end do 317 326 end do 327 !$OMP END DO 328 !$OMP END PARALLEL 318 329 319 330 if (delh.gt.0.) then … … 330 341 331 342 ! thickness at the upwind node 332 do j = 1, ny 333 do i = 2, nx 334 debug_3D(i,j,92) = c_west(i,j) * newH(i-1,j) + c_east(i,j) * newH(i,j) 335 end do 336 end do 337 do j = 2, ny 338 do i = 1, nx 339 debug_3D(i,j,93) = c_south(i,j) * newH(i,j-1) + c_north(i,j) * newH(i,j) 340 end do 341 end do 342 343 if (itracebug.eq.1) call tracebug(' fin routine resolution_diffusion') 343 !do j = 1, ny 344 ! do i = 2, nx 345 ! debug_3D(i,j,92) = c_west(i,j) * newH(i-1,j) + c_east(i,j) * newH(i,j) 346 ! end do 347 !end do 348 !do j = 2, ny 349 ! do i = 1, nx 350 ! debug_3D(i,j,93) = c_south(i,j) * newH(i,j-1) + c_north(i,j) * newH(i,j) 351 ! end do 352 !end do 353 354 if (itracebug.eq.1) call tracebug(' fin routine resolution_diffusion') 355 356 344 357 return 345 358
Note: See TracChangeset
for help on using the changeset viewer.