Original WRF subgrid support version from John Michalakes without fire
[wrffire.git] / wrfv2_fire / external / RSL / RSL / handle_spec3.c
blobe642b2c3a1bce1e6d4722bf073f2ed816c295c47
1 /***********************************************************************
3 COPYRIGHT
5 The following is a notice of limited availability of the code and
6 Government license and disclaimer which must be included in the
7 prologue of the code and in all source listings of the code.
9 Copyright notice
10 (c) 1977 University of Chicago
12 Permission is hereby granted to use, reproduce, prepare
13 derivative works, and to redistribute to others at no charge. If
14 you distribute a copy or copies of the Software, or you modify a
15 copy or copies of the Software or any portion of it, thus forming
16 a work based on the Software and make and/or distribute copies of
17 such work, you must meet the following conditions:
19 a) If you make a copy of the Software (modified or verbatim)
20 it must include the copyright notice and Government
21 license and disclaimer.
23 b) You must cause the modified Software to carry prominent
24 notices stating that you changed specified portions of
25 the Software.
27 This software was authored by:
29 Argonne National Laboratory
30 J. Michalakes: (630) 252-6646; email: michalak@mcs.anl.gov
31 Mathematics and Computer Science Division
32 Argonne National Laboratory, Argonne, IL 60439
34 ARGONNE NATIONAL LABORATORY (ANL), WITH FACILITIES IN THE STATES
35 OF ILLINOIS AND IDAHO, IS OWNED BY THE UNITED STATES GOVERNMENT,
36 AND OPERATED BY THE UNIVERSITY OF CHICAGO UNDER PROVISION OF A
37 CONTRACT WITH THE DEPARTMENT OF ENERGY.
39 GOVERNMENT LICENSE AND DISCLAIMER
41 This computer code material was prepared, in part, as an account
42 of work sponsored by an agency of the United States Government.
43 The Government is granted for itself and others acting on its
44 behalf a paid-up, nonexclusive, irrevocable worldwide license in
45 this data to reproduce, prepare derivative works, distribute
46 copies to the public, perform publicly and display publicly, and
47 to permit others to do so. NEITHER THE UNITED STATES GOVERNMENT
48 NOR ANY AGENCY THEREOF, NOR THE UNIVERSITY OF CHICAGO, NOR ANY OF
49 THEIR EMPLOYEES, MAKES ANY WARRANTY, EXPRESS OR IMPLIED, OR
50 ASSUMES ANY LEGAL LIABILITY OR RESPONSIBILITY FOR THE ACCURACY,
51 COMPLETENESS, OR USEFULNESS OF ANY INFORMATION, APPARATUS,
52 PRODUCT, OR PROCESS DISCLOSED, OR REPRESENTS THAT ITS USE WOULD
53 NOT INFRINGE PRIVATELY OWNED RIGHTS.
55 ***************************************************************************/
57 #include <stdio.h>
58 #include <stdlib.h>
59 #include "rsl.h"
60 #include "which_boundary.h"
62 /* this routine is specific to MM5. It's purpose is to read in boundary
63 data from the boundary file and these records have unusual shapes. There
64 are also multiple fields per record. For not it easier to do this this
65 was, with a special routine, than to incorporate the necessary generality
66 into rsl. This particular implementation of this is especially
67 dreadful -- there is ample room for improvement as time permits. */
69 /* IMPORTANT -- note that the three-d variables have k (level) as their
70 second, not third, dimension in MM5. EG: ueb(MIX,MKX,NSPGD) 2/1/94 */
72 /* rev: 2/3/94, added code to pass and receive extra boundary points to
73 neighboring processors along a 4 pt stencil */
75 /* rev: 7/14/94, changed to a 12 pt stencil to make sure we have enough
76 points for the dot and the cross grids. Also modified rsl_mm_io.c */
78 /* rev: 9/8/94, fixed memory leak wherein pbuf was only freed when
79 rsl_myproc was equal to mdest. Change to free unconditionally
80 for each of the boundaries. */
82 /* rev: 1/10/95, added test to make sure domain is always c.d. This
83 function will not work properly on a nest and should never be
84 called on one. MM5 does not read in boundary data on a nest. */
86 int
87 handle_special3( req, which_boundary, rbuf, buf )
88 rsl_read_req_t * req ;
89 int which_boundary ;
90 char *rbuf ;
91 char **buf ;
93 int dim, i, k, ig, jg, d ;
94 int maj, min, majlen, minlen ;
95 int nelem ;
96 int P ;
97 int bwdth ;
98 int mlen, mtag, mdest ;
99 int psendto[ RSL_MAXPROC ] ; /* size of messages to each processor */
100 rsl_point_t *domain ;
101 int i_am_monitor ;
102 int typelen ;
104 RSL_C_IAMMONITOR( &i_am_monitor ) ;
106 bwdth = req->speciala ;
108 if ( which_boundary == WHICH_BDY_NORTH || which_boundary == WHICH_BDY_SOUTH )
109 nelem = bwdth * req->glen[1] * ( req->ndim==3?req->glen[2]:1 ) ;
110 else if ( which_boundary == WHICH_BDY_WEST || which_boundary == WHICH_BDY_EAST )
111 nelem = bwdth * req->glen[0] * ( req->ndim==3?req->glen[2]:1 ) ;
112 else
113 RSL_TEST_ERR(1,"handle_special3 bad which_boundary") ;
115 typelen = elemsize( req->type ) ;
117 /* figure out sizes of buffers needed */
118 /*caller will free */
119 *buf = RSL_MALLOC( char, nelem*typelen+100) ; /* 100 is safety */
121 d = req->domain ;
122 RSL_TEST_ERR( d != 0,
123 "attempt to read boundary data on domain other than c.d.") ;
124 majlen = domain_info[d].len_n ;
125 minlen = domain_info[d].len_m ;
126 domain = domain_info[d].domain ;
128 for ( i = 0 ; i < rsl_nproc_all ; i++ ) /* 95/02/22 */
130 psendto[i] = 0 ;
132 for ( jg = 0 ; jg < majlen ; jg++ )
134 for ( ig = 0 ; ig < minlen ; ig++ )
136 if ( jg < bwdth ) /* west 1 */
138 psendto[domain[INDEX_2(jg,ig,minlen)].P] = 1 ;
140 if ( jg >= majlen - bwdth && jg < majlen ) /* east 2 */
142 psendto[domain[INDEX_2(jg,ig,minlen)].P] = 1 ;
144 if ( ig >= minlen - bwdth && ig < minlen ) /* north 8 */
146 psendto[domain[INDEX_2(jg,ig,minlen)].P] = 1 ;
148 if ( ig < bwdth ) /* south 4 */
150 psendto[domain[INDEX_2(jg,ig,minlen)].P] = 1 ;
155 if ( i_am_monitor )
157 for ( P = 0 ; P < rsl_nproc_all ; P++ ) /* 95/02/22 */
159 if ( psendto[P] != 0 )
161 mdest = rsl_c_comp2phys_proc( P ) ;
162 if ( ! mdest == rsl_myproc )
164 mtag = MTYPE_FROMTO( MSG_SPECIAL1_RESPONSE, rsl_myproc, mdest ) ;
165 mlen = nelem*typelen ;
166 RSL_SEND( rbuf, mlen, mtag, mdest ) ;
170 bcopy( rbuf, *buf, nelem*typelen ) ;
172 else
174 if ( psendto[rsl_c_phys2comp_proc(rsl_myproc)] )
176 mdest = RSL_C_MONITOR_PROC () ;
177 mtag = MTYPE_FROMTO( MSG_SPECIAL1_RESPONSE, mdest, rsl_myproc ) ;
178 mlen = nelem*typelen ;
179 RSL_RECV( *buf, mlen, mtag ) ;
183 if ( psendto[rsl_c_phys2comp_proc(rsl_myproc)] )
184 return( nelem*typelen ) ; /* whether to expect message */
185 else
186 return(0) ;