[MITgcm-devel] Re: downslope pkg
Michael Schodlok
Michael.P.Schodlok at jpl.nasa.gov
Wed Sep 10 13:40:10 EDT 2008
Hi Martin, Jean-Michel
thanks guys,
there is a figure of RHOAnoma and UVEL on
http://ecco2.jpl.nasa.gov/data4/MPS_output/pub/MITgcm/
UVEL refers to Jean-Michel earlier email, no non zero velocities
in-land
and RHOAnoma to Martins
regardless which packages is applied the solid ground is zero as well
(as in the ice shelf). I'll check out the new code.
happy lawn moving.
Michael
> Micheal,
>
> please have a look at density (diagnostics variable "RHOAnoma"). This
> variable should now be zero in the ice shelf (above the ocean, this is
> required for shelfice) and non-zero in the solid ground (below the
> ocean, necessary for down_slope). Non-zero "theta" and "salt" in "dry
> cells" (hFacC=0) within the ice shelf have no effect.
>
> Martin
>
> On 9 Sep 2008, at 15:34, Jean-Michel Campin wrote:
>
>> Hi Michael,
>>
>> I will check again that it's OK to use both pkgs, iceshelf and
>> down-slope.
>> regarding "the three dimensional diagnostic",
>> hFacC is the real mask, and you should not worry about any value of
>> theta
>> in-land, it's just there for the algorithm to work
>> (this will be fixed one day, but busy right now).
>> So, do you see, e.g, non zero velocities in-land ? or other than
>> theta,salt,rho fields that appear to change in-land ?
>> Thanks,
>> Jean-Michel
>>
>> PS: can you fill the subscription page to the development list:
>> http://mitgcm.org/mailman/listinfo/mitgcm-devel
>> so that we will just approve it.
>>
>> On Mon, Sep 08, 2008 at 07:27:46PM -0700, Michael Schodlok wrote:
>>>
>>> Hi Jean-Michel,
>>>
>>> first I thought the combination of packages iceshelf and down-slope was
>>> running without
>>> a problem as the STDOUT files didn't reveal anything abnormal.
>>> But, now it seems that the mask is lost somewhere when averaging to
>>> three dimensional
>>> fields. While hFacC is showing land and water grid points, the three
>>> dimensional diagnostic
>>> fields seem to have only the mask of one layer (the coastline)
>>> throughout the water column,
>>> i.e., in case of iceshelves included, the cavities are outlined very
>>> nicely.
>>> Could you please check it out:
>>> I've put three directories with different runs on:
>>> http://ecco2.jpl.nasa.gov/data4/MPS_output/pub/MITgcm/
>>>
>>> run_ohne: without downslope and without iceshelves
>>> run_down: with downslope and without iceshelves
>>> run_ice: with downslope and with iceshelves.
>>> the directories include the data files, hFacC and THETA (first month
>>> average).
>>> the model domain is 220x220x50;
>>>
>>> thanks
>>> Michael
>>>
>>>
>>>> Hi Jean-Michel,
>>>> the masking is absolutely necessary for the shelfice package in the
>>>> current implementation (who gave this package this stoopid name, some
>>>> German idiot!), because otherwise the model integrates some ridiculous
>>>> density (for Theta=Salt=0, or even worse from some other arbitrary
>>>> Theta and Salt) from the top for phiHyd. When useShelfice=.true. you
>>>> are supposed to provide a reference pressure at the bottom of the ice
>>>> shelf and integrating phiHyd from the very top (k=1) is not necessary.
>>>> JMD95P and MDJFW should work with shelfice, too, for this reason. I
>>>> don't know about the downSlope package. That one was not around, when
>>>> I implemented shelfice.
>>>> Without the masking, one would have to use if statements in the
>>>> ij-loops for the phiHyd intergration to determine, where the real
>>>> intergration of phiHyd is supposed to start; not good. Alternatively,
>>>> one could only mask alphaRho in the ice shelf, like this:
>>>> IF ( useShelfIce ) THEN
>>>> DO j=jMin,jMax
>>>> DO i=iMin,iMax
>>>> IF ( kTopC(I,J,bi,bj) .GE. k ) alphaRho(i,j) =
>>>> alphaRho(i,j)*maskC(i,j,k,bi,bj)
>>>> ENDDO
>>>> ENDDO
>>>> ENDIF
>>>> Do you have a better idea?
>>>>
>>>> Martin
>>>>
>>>> On 5 Sep 2008, at 18:31, Jean-Michel Campin wrote:
>>>>
>>>>> Hi Michael,
>>>>>
>>>>> I found a problem with shelfice & down_slope pkgs:
>>>>> looks like useShelfIce is masking the density (in CALC_PHI_HYD)
>>>>> so that we will have a problem when using an EOS directly function
>>>>> of pressure (e.g.: JMD95P). I looked to
>>>>> MITgcm_contrib/high_res_cube/data.hr
>>>>> and it's not the case (EOS=JMD95Z).
>>>>> Can you confirm that you are not using an EOS directly function
>>>>> of pressure ?
>>>>>
>>>>> Otherwise, how is the test going ?
>>>>>
>>>>> I cc this to the devel list, so that Martin will know why I will
>>>>> bother him to remove this masking with useShelfIce ...
>>>>>
>>>>> Cheers,
>>>>> Jean-Michel
>>>>> _______________________________________________
>>>>> MITgcm-devel mailing list
>>>>> MITgcm-devel at mitgcm.org
>>>>> http://mitgcm.org/mailman/listinfo/mitgcm-devel
>>>>
>>>
>>>
>>> ----
>>> --
>>>
>>> *************************************************
>>>
>>> Michael Schodlok <schodlok at jpl.nasa.gov>
>>> Jet Propulsion Lab/
>>> California Institute of Technology
>>> MS 300-323
>>> 4800 Oak Grove Dr
>>> Pasadena CA 91109-8099 USA
>>> Tel: +1-818-354-4015 Fax: +1-818-393-6720
>>>
>>> *************************************************
>>>
>>>
>>>
>> _______________________________________________
>> MITgcm-devel mailing list
>> MITgcm-devel at mitgcm.org
>> http://mitgcm.org/mailman/listinfo/mitgcm-devel
>
--
--
--
*************************************************
Michael Schodlok <schodlok at jpl.nasa.gov>
Jet Propulsion Lab/
California Institute of Technology
MS 300-323
4800 Oak Grove Dr
Pasadena
CA 91109-8099
USA
Tel: +1-818-354-4015
Fax: +1-818-393-6720
*************************************************
More information about the MITgcm-devel
mailing list