Possible WET_DRY bug

Bug reports, work arounds and fixes

Moderators: arango, robertson

Post Reply
Message
Author
ianthomas
Posts: 2
Joined: Fri May 09, 2008 3:21 pm
Location: Knowtra Ltd
Contact:

Possible WET_DRY bug

#1 Unread post by ianthomas »

I think I have found a bug in the WET_DRY code, but it is also possible that I am doing something wrong instead (!) so I didn't want to fill out an erroneous Trac ticket.

I have a coastal model set up with land/sea masking, forced by tidal harmonics on both elevation and currents at the open boundaries. The peak tidal amplitude is just over 1 m at the coast. If I set my minimum water depth to 2 m, everything runs fine without WET_DRY. If I set the minimum water depth to 1 m then I enable WET_DRY otherwise I'll get negative water depths. I am using the default DCRIT of 0.1 m.

I understand from https://www.myroms.org/wiki/index.php/WET_DRY that if the water depth is less than DCRIT, i.e.

Code: Select all

zeta(i,j) + h(i,j) < Dcrit
then WET_DRY increases zeta so that this is no longer the case, i.e.

Code: Select all

zeta(i,j) = Dcrit - h(i,j)
Hence in my case, if zeta goes more negative than -0.9 m, then the WET_DRY code should set zeta to -0.9 m. This doesn't seem to be happening. I often get away with values of zeta between -1 and -0.9, but eventually it gets more negative than -1 such that I have negative water depth which results in the model blowing up of course.

Looking through the code for occurrences of Dcrit, it seems that Nonlinear/step2d_LF_AM3.h is probably the file of interest. Lines 2423-2433 check if zeta is less than Dcrit-h and set the wet/dry mask accordingly:

Code: Select all

      DO j=JstrV-1,Jend
        DO i=IstrU-1,Iend
          wetdry(i,j)=1.0_r8
          IF (zwrk(i,j).le.(Dcrit(ng)-h(i,j))) THEN
            wetdry(i,j)=0.0_r8
          END IF
#  ifdef MASKING
          wetdry(i,j)=wetdry(i,j)*rmask(i,j)
#  endif
        END DO
      END DO
I find that if I insert the following before line 2428 inside the IF loop:

Code: Select all

      zeta(i,j,knew) = (Dcrit(ng) - h(i,j))*rmask(i,j)
then wet/dry does what I expect it to and my problem goes away. Now I'm not suggesting that this is a correct solution as I am nowhere near understanding the wet/dry code, but it seems to work for me.

Could someone who understands the wet/dry code take a look at this and see if it is a real problem or if I am doing something wrong? I can supply my model input files if that helps.

jcwarner
Posts: 1200
Joined: Wed Dec 31, 2003 6:16 pm
Location: USGS, USA

Re: Possible WET_DRY bug

#2 Unread post by jcwarner »

I have made several modifications to the wetdry feature, and need to finalize these and send them to the main trunk. We should not be using 'zwrk', and there are several other issues that are too lengthy here to describe. I will try to get this done in the next few weeks. For now, if your fix is working then go with that. But it really does not seem to be the correct solution.

Sorry for the delay on this.

ianthomas
Posts: 2
Joined: Fri May 09, 2008 3:21 pm
Location: Knowtra Ltd
Contact:

Re: Possible WET_DRY bug

#3 Unread post by ianthomas »

Thanks for the quick reply. As you suggest, I'll stick with my temporary fix until the correct solution is available.

Post Reply