[MITgcm-support] Excluding Tiles that are All Land
Matthew Mazloff
mmazloff at ucsd.edu
Thu Sep 12 00:51:15 EDT 2019
Hi Senja
It can be done.
First turn on the exch2 package by compiling with exch2 in packages.conf
In SIZE.h change so nPx is the number of cores you are using and nPy = 1
put a data.exch2 in your run directory.
Here is an example
________________________
# EXCH2 Package: Wrapper-2 User Choice
#--------------------
# preDefTopol :: pre-defined Topology selector:
# :: = 0 : topology defined from processing "data.exch2";
# :: = 1 : simple, single facet topology;
# :: = 2 : customized topology (w2_set_myown_facets)
# :: = 3 : 6-facet Cube (3 face-dims: nRed, nGreen, nBlue).
# dimsFacets :: facet pair of dimensions (n1x,n1y, n2x,n2y ...)
# facetEdgeLink :: Face-Edge connectivity map:
# facetEdgeLink(i,j)=XX.1 : face(j)-edge(i) (i=1,2,3,4 <==> N,S,E,W)
# is connected to Northern edge of face "XX" ; similarly,
# = XX.2 : to Southern.E, XX.3 = Eastern.E, XX.4 = Western.E of face "XX"
# blankList :: List of "blank" tiles
# W2_mapIO :: global map IO selector (-1 = old type ; 0 = 1 long line in X
# :: 1 = compact, mostly in Y dir)
# W2_printMsg :: option for information messages printing
# :: <0 : write to log file ; =0 : minimum print ;
# :: =1 : no duplicated print ; =2 : all processes do print
#--------------------
&W2_EXCH2_PARM01
W2_printMsg= 0,
# preDefTopol= 0,
W2_mapIO = 1,
dimsFacets = NX, NY,
&
________________________
where you need to replace NX and NY by your true domain size.
Nx = sNx*nSx*nPx, Ny = sNy*nSy*nPy,
Now run a timestep. And make sure you run where you get all your STDOUT. (i.e. #undef SINGLE_DISK_IO) You can turn this back on later if you want.
Once you have all your STDOUT. Run these commands
% GET THE BLANK TILES
grep Empty STDOUT.* > trash1.txt
%CROP BEGINNG OF THIS TEXT
sed -r 's .{47} ' trash1.txt > trash2.txt
%AND FIX THE END
sed 's/ (bi,bj= 1 1 )/, /g' trash2.txt > trash3.txt
%NOW GRAB THE TILE NUMBERS
cat trash3.txt | cut -c 1-8 > trash4.txt
%MAKE ALL ONE LINE
cat trash4.txt | tr '\n' ' ' > blanks.txt
%count how many blanks
wc -l trash4.txt
The output of this is how many blank tiles you have.
Now copy the blank list in blanks.txt into data.exch2.
For example if you have 19 blank tiles you would write in data.exch2:
blankList(1:19)=1,2,3,4,5,6,7,8,9,10,19,24,28,29,30,31,32,33,34
Then you would recompile with nPx being 19 smaller. (e.g. if it was 119, its now 100. You now only need 100 cores)
When you now call this executable ask for this new number of cores. It will run the exact same as before, but needs less cores
Hope this helps. I may have forgotten something so don’t hesitate to let me know if something isn’t working
Matt
> On Sep 11, 2019, at 7:51 AM, Senja Walberg <senja.w at gmail.com> wrote:
>
> I've read in some of the documentation that it might be possible to exclude processors from needing to be allocated to tiles which do not have any fluid-containing cells in them.
>
> Does anyone know if this can be easily done or have an example which does this?
>
> Regards,
> Senja
> _______________________________________________
> MITgcm-support mailing list
> MITgcm-support at mitgcm.org
> http://mailman.mitgcm.org/mailman/listinfo/mitgcm-support
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://mailman.mitgcm.org/pipermail/mitgcm-support/attachments/20190911/b423893f/attachment-0001.html>
More information about the MITgcm-support
mailing list