1 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
2 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
3 | !----------------------------------------------------------------------- |
---|
4 | ! CVS m_Accumulator.F90,v 1.37 2008-05-15 13:20:10 larson Exp |
---|
5 | ! CVS MCT_2_8_0 |
---|
6 | !BOP ------------------------------------------------------------------- |
---|
7 | ! |
---|
8 | ! !MODULE: m_Accumulator - Time Averaging/Accumlation Buffer |
---|
9 | ! |
---|
10 | ! !DESCRIPTION: |
---|
11 | ! |
---|
12 | ! An {\em accumulator} is a data class used for computing running sums |
---|
13 | ! and/or time averages of {\tt AttrVect} class data. |
---|
14 | ! The period of time over which data are accumulated/averaged is the |
---|
15 | ! {\em accumulation cycle}, which is defined by the total number |
---|
16 | ! of accumulation steps (the component {\tt Accumulator\%num\_steps}). When |
---|
17 | ! the accumulation routine {\tt accumulate\_} is invoked, the number |
---|
18 | ! of accumulation cycle steps (the component |
---|
19 | ! {\tt Accumulator\%steps\_done})is incremented, and compared with |
---|
20 | ! the number of steps in the accumulation cycle to determine if the |
---|
21 | ! accumulation cycle has been completed. The accumulation buffers |
---|
22 | ! of the {\tt Accumulator} are stored in an {\tt AttrVect} (namely |
---|
23 | ! the component {\tt Accumulator\%data}), which allows the user to |
---|
24 | ! define the number of variables and their names at run-time. |
---|
25 | ! Finally, one can define for each field |
---|
26 | ! being accumulated the specific accumulation {\em action}. Currently, |
---|
27 | ! there are two options: Time Averaging and Time Summation. The |
---|
28 | ! user chooses the specific action by setting an integer action |
---|
29 | ! flag for each attribute being accumulated. The supported options |
---|
30 | ! are defined by the public data member constants {\tt MCT\_SUM} and |
---|
31 | ! {\tt MCT\_AVG}. |
---|
32 | ! \\ |
---|
33 | ! This module also supports a simple usage of accumulator where all |
---|
34 | ! the actions are SUM ({\tt inits\_} and {\tt initavs\_}) and the user |
---|
35 | ! must call {\tt average\_} to calculate the average from the current |
---|
36 | ! value of {\tt Accumulator\%steps\_done}. {\tt Accumulator\%num\_steps} |
---|
37 | ! is ignored in this case. |
---|
38 | ! |
---|
39 | ! !INTERFACE: |
---|
40 | |
---|
41 | module m_Accumulator |
---|
42 | ! |
---|
43 | ! !USES: |
---|
44 | ! |
---|
45 | use m_List, only : List |
---|
46 | use m_AttrVect, only : AttrVect |
---|
47 | use m_realkinds,only : SP,DP,FP |
---|
48 | |
---|
49 | implicit none |
---|
50 | |
---|
51 | private ! except |
---|
52 | |
---|
53 | ! !PUBLIC TYPES: |
---|
54 | |
---|
55 | public :: Accumulator ! The class data structure |
---|
56 | |
---|
57 | Type Accumulator |
---|
58 | #ifdef SEQUENCE |
---|
59 | sequence |
---|
60 | #endif |
---|
61 | integer :: num_steps ! total number of accumulation steps |
---|
62 | integer :: steps_done ! number of accumulation steps performed |
---|
63 | integer, pointer, dimension(:) :: iAction ! index of integer actions |
---|
64 | integer, pointer, dimension(:) :: rAction ! index of real actions |
---|
65 | type(AttrVect) :: data ! accumulated sum field storage |
---|
66 | End Type Accumulator |
---|
67 | |
---|
68 | ! !PUBLIC MEMBER FUNCTIONS: |
---|
69 | ! |
---|
70 | public :: init ! creation method |
---|
71 | public :: initp ! partial creation method (MCT USE ONLY) |
---|
72 | public :: clean ! destruction method |
---|
73 | public :: initialized ! check if initialized |
---|
74 | public :: lsize ! local length of the data arrays |
---|
75 | public :: NumSteps ! number of steps in a cycle |
---|
76 | public :: StepsDone ! number of steps completed in the |
---|
77 | ! current cycle |
---|
78 | public :: nIAttr ! number of integer fields |
---|
79 | public :: nRAttr ! number of real fields |
---|
80 | public :: indexIA ! index the integer fields |
---|
81 | public :: indexRA ! index the real fields |
---|
82 | public :: getIList ! Return tag from INTEGER |
---|
83 | ! attribute list |
---|
84 | public :: getRList ! Return tag from REAL attribute |
---|
85 | ! list |
---|
86 | public :: exportIAttr ! Return INTEGER attribute as a vector |
---|
87 | public :: exportRAttr ! Return REAL attribute as a vector |
---|
88 | public :: importIAttr ! Insert INTEGER vector as attribute |
---|
89 | public :: importRAttr ! Insert REAL vector as attribute |
---|
90 | public :: zero ! Clear an accumulator |
---|
91 | public :: SharedAttrIndexList ! Returns the number of shared |
---|
92 | ! attributes, and lists of the |
---|
93 | ! respective locations of these |
---|
94 | ! shared attributes |
---|
95 | public :: accumulate ! Add AttrVect data into an Accumulator |
---|
96 | public :: average ! Calculate an average in an Accumulator |
---|
97 | |
---|
98 | ! Definition of interfaces for the methods for the Accumulator: |
---|
99 | |
---|
100 | interface init ; module procedure & |
---|
101 | init_, & |
---|
102 | inits_, & |
---|
103 | initv_, & |
---|
104 | initavs_ |
---|
105 | end interface |
---|
106 | interface initp ; module procedure initp_ ; end interface |
---|
107 | interface clean ; module procedure clean_ ; end interface |
---|
108 | interface initialized; module procedure initialized_ ; end interface |
---|
109 | interface lsize ; module procedure lsize_ ; end interface |
---|
110 | interface NumSteps ; module procedure NumSteps_ ; end interface |
---|
111 | interface StepsDone ; module procedure StepsDone_ ; end interface |
---|
112 | interface nIAttr ; module procedure nIAttr_ ; end interface |
---|
113 | interface nRAttr ; module procedure nRAttr_ ; end interface |
---|
114 | interface indexIA; module procedure indexIA_; end interface |
---|
115 | interface indexRA; module procedure indexRA_; end interface |
---|
116 | interface getIList; module procedure getIList_; end interface |
---|
117 | interface getRList; module procedure getRList_; end interface |
---|
118 | interface exportIAttr ; module procedure exportIAttr_ ; end interface |
---|
119 | interface exportRAttr ; module procedure & |
---|
120 | exportRAttrSP_, & |
---|
121 | exportRAttrDP_ |
---|
122 | end interface |
---|
123 | interface importIAttr ; module procedure importIAttr_ ; end interface |
---|
124 | interface importRAttr ; module procedure & |
---|
125 | importRAttrSP_, & |
---|
126 | importRAttrDP_ |
---|
127 | end interface |
---|
128 | interface zero ; module procedure zero_ ; end interface |
---|
129 | interface SharedAttrIndexList ; module procedure & |
---|
130 | aCaCSharedAttrIndexList_, & |
---|
131 | aVaCSharedAttrIndexList_ |
---|
132 | end interface |
---|
133 | interface accumulate ; module procedure accumulate_ ; end interface |
---|
134 | interface average ; module procedure average_ ; end interface |
---|
135 | |
---|
136 | ! !PUBLIC DATA MEMBERS: |
---|
137 | ! |
---|
138 | public :: MCT_SUM |
---|
139 | public :: MCT_AVG |
---|
140 | |
---|
141 | integer, parameter :: MCT_SUM = 1 |
---|
142 | integer, parameter :: MCT_AVG = 2 |
---|
143 | |
---|
144 | ! !REVISION HISTORY: |
---|
145 | ! 7Sep00 - Jay Larson <larson@mcs.anl.gov> - initial prototype |
---|
146 | ! 7Feb01 - Jay Larson <larson@mcs.anl.gov> - Public interfaces |
---|
147 | ! to getIList() and getRList(). |
---|
148 | ! 9Aug01 - E.T. Ong <eong@mcs.anl.gov> - added initialized and |
---|
149 | ! initp_ routines. Added 'action' in Accumulator type. |
---|
150 | ! 6May02 - Jay Larson <larson@mcs.anl.gov> - added import/export |
---|
151 | ! routines. |
---|
152 | ! 26Aug02 - E.T. Ong <eong@mcs.anl.gov> - thourough code revision; |
---|
153 | ! no added routines |
---|
154 | ! 10Jan08 - R. Jacob <jacob@mcs.anl.gov> - add simple accumulator |
---|
155 | ! use support and check documentation. |
---|
156 | !EOP ___________________________________________________________________ |
---|
157 | |
---|
158 | character(len=*),parameter :: myname='MCT::m_Accumulator' |
---|
159 | |
---|
160 | contains |
---|
161 | |
---|
162 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
163 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
164 | !BOP ------------------------------------------------------------------- |
---|
165 | ! |
---|
166 | ! !IROUTINE: init_ - Initialize an Accumulator and its Registers |
---|
167 | ! |
---|
168 | ! !DESCRIPTION: |
---|
169 | ! This routine allocates space for the output {\tt Accumulator} argument |
---|
170 | ! {\tt aC}, and at a minimum sets the number of time steps in an |
---|
171 | ! accumulation cycle (defined by the input {\tt INTEGER} argument |
---|
172 | ! {\tt num\_steps}), and the {\em length} of the {\tt Accumulator} |
---|
173 | ! register buffer (defined by the input {\tt INTEGER} argument {\tt |
---|
174 | ! lsize}). If one wishes to accumulate integer fields, the list of |
---|
175 | ! these fields is defined by the input {\tt CHARACTER} argument |
---|
176 | ! {\tt iList}, which is specified as a colon-delimited set of |
---|
177 | ! substrings (further information regarding this is available in the |
---|
178 | ! routine {\tt init\_()} of the module {\tt m\_AttrVect}). If no |
---|
179 | ! value of {\tt iList} is supplied, no integer attribute accumulation |
---|
180 | ! buffers will be allocated. The accumulation action on each of the |
---|
181 | ! integer attributes can be defined by supplying the input {\tt INTEGER} |
---|
182 | ! array argument {\tt iAction(:)} (whose length must correspond to the |
---|
183 | ! number of items in {\tt iList}). The values of the elements of |
---|
184 | ! {\tt iAction(:)} must be one of the values among the public data |
---|
185 | ! members defined in the declaration section of this module. If the |
---|
186 | ! integer attributes are to be accumulated (i.e. one supplies {\tt iList}), |
---|
187 | ! but {\tt iAction(:)} is not specified, the default action for all |
---|
188 | ! integer accumulation operations will be summation. The input arguments |
---|
189 | ! {\tt rList} and {\tt rAction(:)} define the names of the real variables |
---|
190 | ! to be accumulated and the accumulation action for each. The arguments |
---|
191 | ! {\tt rList} and {\tt rAction(:)} are related to each other the same |
---|
192 | ! way as {\tt iList} and {\tt iAction(:)}. Finally, the user can |
---|
193 | ! manually set the number of completed steps in an accumulation cycle |
---|
194 | ! (e.g. for restart purposes) by supplying a value for the optional |
---|
195 | ! input {\tt INTEGER} argument {\tt steps\_done}. |
---|
196 | ! |
---|
197 | ! !INTERFACE: |
---|
198 | |
---|
199 | subroutine init_(aC, iList, iAction, rList, rAction, lsize, & |
---|
200 | num_steps,steps_done) |
---|
201 | ! |
---|
202 | ! !USES: |
---|
203 | ! |
---|
204 | use m_AttrVect, only : AttrVect_init => init |
---|
205 | use m_AttrVect, only : AttrVect_zero => zero |
---|
206 | |
---|
207 | use m_List, only: List |
---|
208 | use m_List, only: List_nullify => nullify |
---|
209 | use m_List, only: List_init => init |
---|
210 | use m_List, only: List_nitem => nitem |
---|
211 | use m_List, only: List_clean => clean |
---|
212 | |
---|
213 | use m_stdio |
---|
214 | use m_die |
---|
215 | |
---|
216 | implicit none |
---|
217 | |
---|
218 | ! !INPUT PARAMETERS: |
---|
219 | ! |
---|
220 | character(len=*), optional, intent(in) :: iList |
---|
221 | integer, dimension(:), optional, intent(in) :: iAction |
---|
222 | character(len=*), optional, intent(in) :: rList |
---|
223 | integer, dimension(:), optional, intent(in) :: rAction |
---|
224 | integer, intent(in) :: lsize |
---|
225 | integer, intent(in) :: num_steps |
---|
226 | integer, optional, intent(in) :: steps_done |
---|
227 | |
---|
228 | ! !OUTPUT PARAMETERS: |
---|
229 | ! |
---|
230 | type(Accumulator), intent(out) :: aC |
---|
231 | |
---|
232 | ! !REVISION HISTORY: |
---|
233 | ! 11Sep00 - Jay Larson <larson@mcs.anl.gov> - initial prototype |
---|
234 | ! 27JUL01 - E.T. Ong <eong@mcs.anl.gov> - added iAction, rAction, |
---|
235 | ! niAction, and nrAction to accumulator type. Also defined |
---|
236 | ! MCT_SUM and MCT_AVG for accumulator module. |
---|
237 | !EOP ___________________________________________________________________ |
---|
238 | ! |
---|
239 | character(len=*),parameter :: myname_=myname//'::init_' |
---|
240 | integer :: my_steps_done, nIAttr, nRAttr, ierr |
---|
241 | integer, dimension(:), pointer :: my_iAction, my_rAction |
---|
242 | logical :: status |
---|
243 | type(List) :: temp_iList, temp_rList |
---|
244 | |
---|
245 | nullify(my_iAction) |
---|
246 | nullify(my_rAction) |
---|
247 | |
---|
248 | call List_nullify(temp_iList) |
---|
249 | call List_nullify(temp_rList) |
---|
250 | |
---|
251 | ! Argument consistency checks: |
---|
252 | |
---|
253 | ! 1) Terminate with error message if optional argument iAction (rAction) |
---|
254 | ! is supplied but optional argument iList (rList) is not. |
---|
255 | |
---|
256 | if(present(iAction) .and. (.not. present(iList))) then |
---|
257 | write(stderr,'(2a)') myname_,'::FATAL--Argument iAction supplied but action iList absent!' |
---|
258 | call die(myname_) |
---|
259 | endif |
---|
260 | |
---|
261 | if(present(rAction) .and. (.not. present(rList))) then |
---|
262 | write(stderr,'(2a)') myname_,'::FATAL--Argument rAction supplied but action rList absent!' |
---|
263 | call die(myname_) |
---|
264 | endif |
---|
265 | |
---|
266 | ! 2) For iList and rList, generate temporary List data structures to facilitate |
---|
267 | ! attribute counting. |
---|
268 | |
---|
269 | if(present(iList)) then ! create temp_iList |
---|
270 | call List_init(temp_iList, iList) |
---|
271 | nIAttr = List_nitem(temp_iList) |
---|
272 | endif |
---|
273 | |
---|
274 | if(present(rList)) then ! create temp_iList |
---|
275 | call List_init(temp_rList, rList) |
---|
276 | nRAttr = List_nitem(temp_rList) |
---|
277 | endif |
---|
278 | |
---|
279 | ! 3) Terminate with error message if optional arguments iAction (rAction) |
---|
280 | ! and iList (rList) are supplied but the size of iAction (rAction) does not |
---|
281 | ! match the number of items in iList (rList). |
---|
282 | |
---|
283 | if(present(iAction) .and. present(iList)) then |
---|
284 | if(size(iAction) /= nIAttr) then |
---|
285 | write(stderr,'(2a,2(a,i8))') myname_, & |
---|
286 | '::FATAL--Size mismatch between iAction and iList! ', & |
---|
287 | 'size(iAction)=',size(iAction),', ','No. items in iList=',nIAttr |
---|
288 | call die(myname_) |
---|
289 | endif |
---|
290 | endif |
---|
291 | |
---|
292 | if(present(rAction) .and. present(rList)) then |
---|
293 | if(size(rAction) /= nRAttr) then |
---|
294 | write(stderr,'(2a,2(a,i8))') myname_, & |
---|
295 | '::FATAL--Size mismatch between rAction and rList! ', & |
---|
296 | 'size(rAction)=',size(rAction),', ','No items in rList=',nRAttr |
---|
297 | call die(myname_) |
---|
298 | endif |
---|
299 | endif |
---|
300 | |
---|
301 | ! Initialize the Accumulator components. |
---|
302 | |
---|
303 | ! steps_done: |
---|
304 | |
---|
305 | if(present(steps_done)) then |
---|
306 | my_steps_done = steps_done |
---|
307 | else |
---|
308 | my_steps_done = 0 |
---|
309 | endif |
---|
310 | |
---|
311 | ! my_iAction (if iList is present) |
---|
312 | |
---|
313 | if(present(iList)) then ! set up my_iAction |
---|
314 | |
---|
315 | allocate(my_iAction(nIAttr), stat=ierr) |
---|
316 | if(ierr /= 0) then |
---|
317 | write(stderr,'(2a,i8)') myname_, & |
---|
318 | '::FATAL: allocate(my_iAction) failed with ierr=',ierr |
---|
319 | call die(myname_) |
---|
320 | endif |
---|
321 | |
---|
322 | if(present(iAction)) then ! use its values |
---|
323 | my_iAction = iAction |
---|
324 | else ! go with default summation by assigning value MCT_SUM |
---|
325 | my_iAction = MCT_SUM |
---|
326 | endif |
---|
327 | |
---|
328 | endif |
---|
329 | |
---|
330 | ! my_rAction (if rList is present) |
---|
331 | |
---|
332 | if(present(rList)) then ! set up my_rAction |
---|
333 | |
---|
334 | allocate(my_rAction(nRAttr), stat=ierr) |
---|
335 | if(ierr /= 0) then |
---|
336 | write(stderr,'(2a,i8)') myname_, & |
---|
337 | '::FATAL: allocate(my_rAction) failed with ierr=',ierr |
---|
338 | call die(myname_) |
---|
339 | endif |
---|
340 | |
---|
341 | if(present(rAction)) then ! use its values |
---|
342 | my_rAction = rAction |
---|
343 | else ! go with default summation by assigning value MCT_SUM |
---|
344 | my_rAction = MCT_SUM |
---|
345 | endif |
---|
346 | |
---|
347 | endif |
---|
348 | |
---|
349 | ! Build the Accumulator aC minus its data component: |
---|
350 | |
---|
351 | if(present(iList) .and. present(rList)) then ! Both REAL and INTEGER registers |
---|
352 | |
---|
353 | call initp_(aC,my_iAction,my_rAction,num_steps,my_steps_done) |
---|
354 | |
---|
355 | deallocate(my_iAction, my_rAction, stat=ierr) |
---|
356 | if(ierr /= 0) then |
---|
357 | write(stderr,'(2a,i8)') myname_, & |
---|
358 | '::FATAL: deallocate(my_iAction, my_rAction) failed with ierr=',ierr |
---|
359 | call die(myname_) |
---|
360 | endif |
---|
361 | |
---|
362 | else ! Either only REAL or only INTEGER registers in aC |
---|
363 | |
---|
364 | if(present(iList)) then ! Only INTEGER REGISTERS |
---|
365 | |
---|
366 | call initp_(aC=aC, iAction=my_iAction, num_steps=num_steps, & |
---|
367 | steps_done=my_steps_done) |
---|
368 | |
---|
369 | deallocate(my_iAction, stat=ierr) |
---|
370 | if(ierr /= 0) then |
---|
371 | write(stderr,'(2a,i8)') myname_, & |
---|
372 | '::FATAL: deallocate(my_iAction) failed with ierr=',ierr |
---|
373 | call die(myname_) |
---|
374 | endif |
---|
375 | |
---|
376 | endif |
---|
377 | |
---|
378 | if(present(rList)) then ! Only REAL REGISTERS |
---|
379 | |
---|
380 | call initp_(aC=aC, rAction=my_rAction, num_steps=num_steps, & |
---|
381 | steps_done=my_steps_done) |
---|
382 | |
---|
383 | deallocate(my_rAction, stat=ierr) |
---|
384 | if(ierr /= 0) then |
---|
385 | write(stderr,'(2a,i8)') myname_, & |
---|
386 | '::FATAL: deallocate(my_rAction) failed with ierr=',ierr |
---|
387 | call die(myname_) |
---|
388 | endif |
---|
389 | |
---|
390 | endif |
---|
391 | |
---|
392 | endif |
---|
393 | |
---|
394 | ! Initialize the AttrVect data component for aC: |
---|
395 | |
---|
396 | if(present(iList) .and. present(rList)) then |
---|
397 | call AttrVect_init(aC%data,iList,rList,lsize) |
---|
398 | else |
---|
399 | if(present(iList)) then |
---|
400 | call AttrVect_init(aV=aC%data,iList=iList,lsize=lsize) |
---|
401 | endif |
---|
402 | if(present(rList)) then |
---|
403 | call AttrVect_init(aV=aC%data,rList=rList,lsize=lsize) |
---|
404 | endif |
---|
405 | endif |
---|
406 | |
---|
407 | call AttrVect_zero(aC%data) |
---|
408 | |
---|
409 | ! Clean up |
---|
410 | |
---|
411 | if(present(iList)) call List_clean(temp_iList) |
---|
412 | if(present(rList)) call List_clean(temp_rList) |
---|
413 | |
---|
414 | ! Check that aC has been properly initialized |
---|
415 | |
---|
416 | status = initialized_(aC=aC,die_flag=.true.,source_name=myname_) |
---|
417 | |
---|
418 | end subroutine init_ |
---|
419 | |
---|
420 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
421 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
422 | !BOP ------------------------------------------------------------------- |
---|
423 | ! |
---|
424 | ! !IROUTINE: inits_ - Initialize a simple Accumulator and its Registers |
---|
425 | ! |
---|
426 | ! !DESCRIPTION: |
---|
427 | ! This routine allocates space for the output simple {\tt Accumulator} argument |
---|
428 | ! {\tt aC}, and sets the {\em length} of the {\tt Accumulator} |
---|
429 | ! register buffer (defined by the input {\tt INTEGER} argument {\tt |
---|
430 | ! lsize}). If one wishes to accumulate integer fields, the list of |
---|
431 | ! these fields is defined by the input {\tt CHARACTER} argument |
---|
432 | ! {\tt iList}, which is specified as a colon-delimited set of |
---|
433 | ! substrings (further information regarding this is available in the |
---|
434 | ! routine {\tt init\_()} of the module {\tt m\_AttrVect}). If no |
---|
435 | ! value of {\tt iList} is supplied, no integer attribute accumulation |
---|
436 | ! buffers will be allocated. The input argument {\tt rList} define |
---|
437 | ! the names of the real variables to be accumulated. Finally, the user can |
---|
438 | ! manually set the number of completed steps in an accumulation cycle |
---|
439 | ! (e.g. for restart purposes) by supplying a value for the optional |
---|
440 | ! input {\tt INTEGER} argument {\tt steps\_done}. |
---|
441 | ! Its default value is zero. |
---|
442 | ! |
---|
443 | ! In a simple accumulator, the action is always SUM. |
---|
444 | ! |
---|
445 | ! |
---|
446 | ! !INTERFACE: |
---|
447 | |
---|
448 | subroutine inits_(aC, iList, rList, lsize,steps_done) |
---|
449 | ! |
---|
450 | ! !USES: |
---|
451 | ! |
---|
452 | use m_List, only : List_init => init |
---|
453 | use m_List, only : List_clean => clean |
---|
454 | use m_List, only : List_nitem => nitem |
---|
455 | use m_AttrVect, only : AttrVect_init => init |
---|
456 | use m_AttrVect, only : AttrVect_zero => zero |
---|
457 | use m_die |
---|
458 | |
---|
459 | implicit none |
---|
460 | |
---|
461 | ! !INPUT PARAMETERS: |
---|
462 | ! |
---|
463 | character(len=*), optional, intent(in) :: iList |
---|
464 | character(len=*), optional, intent(in) :: rList |
---|
465 | integer, intent(in) :: lsize |
---|
466 | integer, optional, intent(in) :: steps_done |
---|
467 | |
---|
468 | ! !OUTPUT PARAMETERS: |
---|
469 | ! |
---|
470 | type(Accumulator), intent(out) :: aC |
---|
471 | |
---|
472 | ! !REVISION HISTORY: |
---|
473 | ! 10Jan08 - R. Jacob <jacob@mcs.anlgov> - initial version based on init_ |
---|
474 | ! |
---|
475 | !EOP ___________________________________________________________________ |
---|
476 | ! |
---|
477 | character(len=*),parameter :: myname_=myname//'::inits_' |
---|
478 | type(List) :: tmplist |
---|
479 | integer :: my_steps_done,ier,i,actsize |
---|
480 | logical :: status |
---|
481 | |
---|
482 | ! Initialize the Accumulator components. |
---|
483 | |
---|
484 | if(present(steps_done)) then |
---|
485 | my_steps_done = steps_done |
---|
486 | else |
---|
487 | my_steps_done = 0 |
---|
488 | endif |
---|
489 | |
---|
490 | aC%num_steps = -1 ! special value for simple aC |
---|
491 | aC%steps_done = my_steps_done |
---|
492 | |
---|
493 | nullify(aC%iAction,aC%rAction) |
---|
494 | |
---|
495 | if(present(iList)) then |
---|
496 | call List_init(tmplist,iList) |
---|
497 | actsize=List_nitem(tmplist) |
---|
498 | allocate(aC%iAction(actsize),stat=ier) |
---|
499 | if(ier /= 0) call die(myname_,"iAction allocate",ier) |
---|
500 | do i=1,lsize |
---|
501 | aC%iAction=MCT_SUM |
---|
502 | enddo |
---|
503 | call List_clean(tmplist) |
---|
504 | endif |
---|
505 | |
---|
506 | if(present(rList)) then |
---|
507 | call List_init(tmplist,rList) |
---|
508 | actsize=List_nitem(tmpList) |
---|
509 | allocate(aC%rAction(actsize),stat=ier) |
---|
510 | if(ier /= 0) call die(myname_,"rAction allocate",ier) |
---|
511 | do i=1,lsize |
---|
512 | aC%rAction=MCT_SUM |
---|
513 | enddo |
---|
514 | call List_clean(tmplist) |
---|
515 | endif |
---|
516 | |
---|
517 | ! Initialize the AttrVect component aC: |
---|
518 | |
---|
519 | if(present(iList) .and. present(rList)) then |
---|
520 | call AttrVect_init(aC%data,iList,rList,lsize) |
---|
521 | else |
---|
522 | if(present(iList)) then |
---|
523 | call AttrVect_init(aV=aC%data,iList=iList,lsize=lsize) |
---|
524 | endif |
---|
525 | if(present(rList)) then |
---|
526 | call AttrVect_init(aV=aC%data,rList=rList,lsize=lsize) |
---|
527 | endif |
---|
528 | endif |
---|
529 | |
---|
530 | call AttrVect_zero(aC%data) |
---|
531 | |
---|
532 | ! Check that aC has been properly initialized |
---|
533 | |
---|
534 | status = initialized_(aC=aC,die_flag=.true.,source_name=myname_) |
---|
535 | |
---|
536 | end subroutine inits_ |
---|
537 | |
---|
538 | |
---|
539 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
540 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
541 | !BOP ------------------------------------------------------------------- |
---|
542 | ! |
---|
543 | ! !IROUTINE: initp_ - Initialize an Accumulator but not its Registers |
---|
544 | ! |
---|
545 | ! !DESCRIPTION: |
---|
546 | ! This routine is an internal service routine for use by the other |
---|
547 | ! initialization routines in this module. It sets up some---but not |
---|
548 | ! all---of the components of the output {\tt Accumulator} argument |
---|
549 | ! {\tt aC}. This routine can set up the following components of |
---|
550 | ! {\tt aC}: |
---|
551 | ! \begin{enumerate} |
---|
552 | ! \item {\tt aC\%iAction}, the array of accumlation actions for the |
---|
553 | ! integer attributes of {\tt aC} (if the input {\tt INTEGER} array |
---|
554 | ! argument {\tt iAction(:)} is supplied); |
---|
555 | ! \item {\tt aC\%rAction}, the array of accumlation actions for the |
---|
556 | ! real attributes of {\tt aC} (if the input {\tt INTEGER} array |
---|
557 | ! argument {\tt rAction(:)} is supplied); |
---|
558 | ! \item {\tt aC\%num\_steps}, the number of steps in an accumulation |
---|
559 | ! cycle (if the input {\tt INTEGER} argument {\tt num\_steps} is |
---|
560 | ! supplied); and |
---|
561 | ! \item {\tt aC\%steps\_done}, the number of steps completed so far |
---|
562 | ! in an accumulation cycle (if the input {\tt INTEGER} argument |
---|
563 | ! {\tt steps\_done} is supplied). |
---|
564 | ! \end{enumerate} |
---|
565 | ! |
---|
566 | ! !INTERFACE: |
---|
567 | |
---|
568 | subroutine initp_(aC, iAction, rAction, num_steps, steps_done) |
---|
569 | |
---|
570 | ! |
---|
571 | ! !USES: |
---|
572 | ! |
---|
573 | use m_die |
---|
574 | |
---|
575 | implicit none |
---|
576 | |
---|
577 | ! !INPUT PARAMETERS: |
---|
578 | ! |
---|
579 | integer, dimension(:), optional, intent(in) :: iAction |
---|
580 | integer, dimension(:), optional, intent(in) :: rAction |
---|
581 | integer, intent(in) :: num_steps |
---|
582 | integer, optional, intent(in) :: steps_done |
---|
583 | |
---|
584 | ! !OUTPUT PARAMETERS: |
---|
585 | ! |
---|
586 | type(Accumulator), intent(out) :: aC |
---|
587 | |
---|
588 | ! !REVISION HISTORY: |
---|
589 | ! 11Sep00 - Jay Larson <larson@mcs.anl.gov> - initial prototype |
---|
590 | ! 27JUL01 - E.T. Ong <eong@mcs.anl.gov> - added iAction, rAction, |
---|
591 | ! niAction, and nrAction to accumulator type. Also defined |
---|
592 | ! MCT_SUM and MCT_AVG for accumulator module. |
---|
593 | !EOP ___________________________________________________________________ |
---|
594 | ! |
---|
595 | character(len=*),parameter :: myname_=myname//'::initp_' |
---|
596 | integer :: i,ier |
---|
597 | integer :: steps_completed |
---|
598 | |
---|
599 | ! if the argument steps_done is not present, assume |
---|
600 | ! the accumulator is starting at step zero, that is, |
---|
601 | ! set steps_completed to zero |
---|
602 | |
---|
603 | steps_completed = 0 |
---|
604 | if(present(steps_done)) steps_completed = steps_done |
---|
605 | |
---|
606 | ! Set the stepping info: |
---|
607 | |
---|
608 | aC%num_steps = num_steps |
---|
609 | aC%steps_done = steps_completed |
---|
610 | |
---|
611 | |
---|
612 | ! Assign iAction and niAction components |
---|
613 | |
---|
614 | nullify(aC%iAction,aC%rAction) |
---|
615 | |
---|
616 | if(present(iAction)) then |
---|
617 | |
---|
618 | if(size(iAction)>0) then |
---|
619 | |
---|
620 | allocate(aC%iAction(size(iAction)),stat=ier) |
---|
621 | if(ier /= 0) call die(myname_,"iAction allocate",ier) |
---|
622 | |
---|
623 | do i=1,size(iAction) |
---|
624 | aC%iAction(i) = iAction(i) |
---|
625 | enddo |
---|
626 | |
---|
627 | endif |
---|
628 | |
---|
629 | endif |
---|
630 | |
---|
631 | if(present(rAction)) then |
---|
632 | |
---|
633 | if(size(rAction)>0) then |
---|
634 | |
---|
635 | allocate(aC%rAction(size(rAction)),stat=ier) |
---|
636 | if(ier /= 0) call die(myname_,"iAction allocate",ier) |
---|
637 | |
---|
638 | do i=1,size(rAction) |
---|
639 | aC%rAction(i) = rAction(i) |
---|
640 | enddo |
---|
641 | |
---|
642 | endif |
---|
643 | |
---|
644 | endif |
---|
645 | |
---|
646 | end subroutine initp_ |
---|
647 | |
---|
648 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
649 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
650 | !BOP ------------------------------------------------------------------- |
---|
651 | ! |
---|
652 | ! !IROUTINE: initv_ - Initialize One Accumulator using Another |
---|
653 | ! |
---|
654 | ! !DESCRIPTION: |
---|
655 | ! This routine takes the integer and real attribute information (including |
---|
656 | ! accumulation action settings for each attribute) from a previously |
---|
657 | ! initialized {\tt Accumulator} (the input argument {\tt bC}), and uses |
---|
658 | ! it to create another {\tt Accumulator} (the output argument {\tt aC}). |
---|
659 | ! In the absence of the {\tt INTEGER} input arguments {\tt lsize}, |
---|
660 | ! {\tt num\_steps}, and {\tt steps\_done}, {\tt aC} will inherit from |
---|
661 | ! {\tt bC} its length, the number of steps in its accumulation cycle, and |
---|
662 | ! the number of steps completed in its present accumulation cycle, |
---|
663 | ! respectively. |
---|
664 | ! |
---|
665 | ! !INTERFACE: |
---|
666 | |
---|
667 | subroutine initv_(aC, bC, lsize, num_steps, steps_done) |
---|
668 | ! |
---|
669 | ! !USES: |
---|
670 | ! |
---|
671 | use m_List, only : List |
---|
672 | use m_List, only : ListExportToChar => exportToChar |
---|
673 | use m_List, only : List_copy => copy |
---|
674 | use m_List, only : List_allocated => allocated |
---|
675 | use m_List, only : List_clean => clean |
---|
676 | use m_die |
---|
677 | |
---|
678 | implicit none |
---|
679 | |
---|
680 | ! !INPUT PARAMETERS: |
---|
681 | ! |
---|
682 | type(Accumulator), intent(in) :: bC |
---|
683 | integer, optional, intent(in) :: lsize |
---|
684 | integer, optional, intent(in) :: num_steps |
---|
685 | integer, optional, intent(in) :: steps_done |
---|
686 | |
---|
687 | ! !OUTPUT PARAMETERS: |
---|
688 | ! |
---|
689 | type(Accumulator), intent(out) :: aC |
---|
690 | |
---|
691 | ! !REVISION HISTORY: |
---|
692 | ! 11Sep00 - Jay Larson <larson@mcs.anl.gov> - initial prototype |
---|
693 | ! 17May01 - R. Jacob <jacob@mcs.anl.gov> - change string_get to |
---|
694 | ! list_get |
---|
695 | ! 27JUL01 - E.T. Ong <eong@mcs.anl.gov> - added iaction,raction |
---|
696 | ! compatibility |
---|
697 | ! 2Aug02 - J. Larson <larson@mcs.anl.gov> made argument num_steps |
---|
698 | ! optional |
---|
699 | !EOP ___________________________________________________________________ |
---|
700 | |
---|
701 | character(len=*),parameter :: myname_=myname//'::initv_' |
---|
702 | |
---|
703 | type(List) :: temp_iList, temp_rList |
---|
704 | integer :: myNumSteps, myStepsDone |
---|
705 | integer :: aC_lsize |
---|
706 | integer :: niActions, nrActions |
---|
707 | integer, dimension(:), allocatable :: iActionArray, rActionArray |
---|
708 | integer :: i,ier |
---|
709 | logical :: status |
---|
710 | |
---|
711 | ! Check that bC has been initialized |
---|
712 | |
---|
713 | status = initialized(aC=bC,die_flag=.true.,source_name=myname_) |
---|
714 | |
---|
715 | ! If the argument steps_done is present, set myStepsDone |
---|
716 | ! to this value; otherwise, set it to zero |
---|
717 | |
---|
718 | if(present(num_steps)) then ! set it manually |
---|
719 | myNumSteps = num_steps |
---|
720 | else ! inherit it from bC |
---|
721 | myNumSteps = bC%num_steps |
---|
722 | endif |
---|
723 | |
---|
724 | ! If the argument steps_done is present, set myStepsDone |
---|
725 | ! to this value; otherwise, set it to zero |
---|
726 | |
---|
727 | if(present(steps_done)) then ! set it manually |
---|
728 | myStepsDone= steps_done |
---|
729 | else ! inherit it from bC |
---|
730 | myStepsDone = bC%steps_done |
---|
731 | endif |
---|
732 | |
---|
733 | ! If the argument lsize is present, |
---|
734 | ! set aC_lsize to this value; otherwise, set it to the lsize of bC |
---|
735 | |
---|
736 | if(present(lsize)) then ! set it manually |
---|
737 | aC_lsize = lsize |
---|
738 | else ! inherit it from bC |
---|
739 | aC_lsize = lsize_(bC) |
---|
740 | endif |
---|
741 | |
---|
742 | ! Convert the two Lists to two Strings |
---|
743 | |
---|
744 | niActions = 0 |
---|
745 | nrActions = 0 |
---|
746 | |
---|
747 | if(List_allocated(bC%data%iList)) then |
---|
748 | call List_copy(temp_iList,bC%data%iList) |
---|
749 | niActions = nIAttr_(bC) |
---|
750 | endif |
---|
751 | |
---|
752 | if(List_allocated(bC%data%rList)) then |
---|
753 | call List_copy(temp_rList,bC%data%rList) |
---|
754 | nrActions = nRAttr_(bC) |
---|
755 | endif |
---|
756 | |
---|
757 | ! Convert the pointers to arrays |
---|
758 | |
---|
759 | allocate(iActionArray(niActions),rActionArray(nrActions),stat=ier) |
---|
760 | if(ier /= 0) call die(myname_,"iActionArray/rActionArray allocate",ier) |
---|
761 | |
---|
762 | if( niActions>0 ) then |
---|
763 | do i=1,niActions |
---|
764 | iActionArray(i)=bC%iAction(i) |
---|
765 | enddo |
---|
766 | endif |
---|
767 | |
---|
768 | if( nrActions>0 ) then |
---|
769 | do i=1,nrActions |
---|
770 | rActionArray(i)=bC%rAction(i) |
---|
771 | enddo |
---|
772 | endif |
---|
773 | |
---|
774 | ! Call init with present arguments |
---|
775 | |
---|
776 | if( (niActions>0) .and. (nrActions>0) ) then |
---|
777 | |
---|
778 | call init_(aC, iList=ListExportToChar(temp_iList), & |
---|
779 | iAction=iActionArray, & |
---|
780 | rList=ListExportToChar(temp_rList), & |
---|
781 | rAction=rActionArray, & |
---|
782 | lsize=aC_lsize, & |
---|
783 | num_steps=myNumSteps, & |
---|
784 | steps_done=myStepsDone) |
---|
785 | |
---|
786 | else |
---|
787 | |
---|
788 | if( niActions>0 ) then |
---|
789 | |
---|
790 | call init_(aC, iList=ListExportToChar(temp_iList), & |
---|
791 | iAction=iActionArray, & |
---|
792 | lsize=aC_lsize, & |
---|
793 | num_steps=myNumSteps, & |
---|
794 | steps_done=myStepsDone) |
---|
795 | |
---|
796 | endif |
---|
797 | |
---|
798 | if( nrActions>0 ) then |
---|
799 | |
---|
800 | call init_(aC, rList=ListExportToChar(temp_rList), & |
---|
801 | rAction=rActionArray, & |
---|
802 | lsize=aC_lsize, & |
---|
803 | num_steps=myNumSteps, & |
---|
804 | steps_done=myStepsDone) |
---|
805 | endif |
---|
806 | |
---|
807 | endif |
---|
808 | |
---|
809 | if(List_allocated(bC%data%iList)) call List_clean(temp_iList) |
---|
810 | if(List_allocated(bC%data%rList)) call List_clean(temp_rList) |
---|
811 | |
---|
812 | deallocate(iActionArray,rActionArray,stat=ier) |
---|
813 | if(ier /= 0) call die(myname_,"iActionArray/rActionArray deallocate",ier) |
---|
814 | |
---|
815 | ! Check that aC as been properly initialized |
---|
816 | |
---|
817 | status = initialized(aC=aC,die_flag=.true.,source_name=myname_) |
---|
818 | |
---|
819 | end subroutine initv_ |
---|
820 | |
---|
821 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
822 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
823 | !BOP ------------------------------------------------------------------- |
---|
824 | ! |
---|
825 | ! !IROUTINE: initavs_ - Initialize a simple Accumulator from an AttributeVector |
---|
826 | ! |
---|
827 | ! !DESCRIPTION: |
---|
828 | ! This routine takes the integer and real attribute information (including |
---|
829 | ! from a previously initialized {\tt AttributeVector} (the input argument {\tt aV}), and uses |
---|
830 | ! it to create a simple (sum only) {\tt Accumulator} (the output argument {\tt aC}). |
---|
831 | ! In the absence of the {\tt INTEGER} input argument {\tt lsize}, |
---|
832 | ! {\tt aC} will inherit from {\tt Av} its length. In the absence of the |
---|
833 | ! optional INTEGER argument, {\tt steps\_done} will be set to zero. |
---|
834 | ! |
---|
835 | ! !INTERFACE: |
---|
836 | |
---|
837 | subroutine initavs_(aC, aV, acsize, steps_done) |
---|
838 | ! |
---|
839 | ! !USES: |
---|
840 | ! |
---|
841 | use m_AttrVect, only: AttrVect_lsize => lsize |
---|
842 | use m_AttrVect, only: AttrVect_nIAttr => nIAttr |
---|
843 | use m_AttrVect, only: AttrVect_nRAttr => nRAttr |
---|
844 | use m_AttrVect, only: AttrVect_exIL2c => exportIListToChar |
---|
845 | use m_AttrVect, only: AttrVect_exRL2c => exportRListToChar |
---|
846 | use m_die |
---|
847 | |
---|
848 | implicit none |
---|
849 | |
---|
850 | ! !INPUT PARAMETERS: |
---|
851 | ! |
---|
852 | type(AttrVect), intent(in) :: aV |
---|
853 | integer, optional, intent(in) :: acsize |
---|
854 | integer, optional, intent(in) :: steps_done |
---|
855 | |
---|
856 | ! !OUTPUT PARAMETERS: |
---|
857 | ! |
---|
858 | type(Accumulator), intent(out) :: aC |
---|
859 | |
---|
860 | ! !REVISION HISTORY: |
---|
861 | ! 10Jan08 - R. Jacob <jacob@mcs.anl.gov> - initial version based on initv_ |
---|
862 | !EOP ___________________________________________________________________ |
---|
863 | |
---|
864 | character(len=*),parameter :: myname_=myname//'::initavs_' |
---|
865 | |
---|
866 | integer :: myNumSteps, myStepsDone |
---|
867 | integer :: aC_lsize |
---|
868 | integer :: i,ier |
---|
869 | integer :: nIatt,nRatt |
---|
870 | logical :: status |
---|
871 | |
---|
872 | |
---|
873 | ! If the argument steps_done is present, set myStepsDone |
---|
874 | ! to this value; otherwise, set it to zero |
---|
875 | |
---|
876 | if(present(steps_done)) then ! set it manually |
---|
877 | myStepsDone= steps_done |
---|
878 | else ! set it to zero |
---|
879 | myStepsDone = 0 |
---|
880 | endif |
---|
881 | |
---|
882 | ! If the argument acsize is present, |
---|
883 | ! set aC_lsize to this value; otherwise, set it to the lsize of bC |
---|
884 | |
---|
885 | if(present(acsize)) then ! set it manually |
---|
886 | aC_lsize = acsize |
---|
887 | else ! inherit it from bC |
---|
888 | aC_lsize = AttrVect_lsize(aV) |
---|
889 | endif |
---|
890 | nIatt=AttrVect_nIAttr(aV) |
---|
891 | nRatt=AttrVect_nRAttr(aV) |
---|
892 | |
---|
893 | if((nIAtt>0) .and. (nRatt>0)) then |
---|
894 | call inits_(aC,AttrVect_exIL2c(aV),AttrVect_exRL2c(aV), & |
---|
895 | aC_lsize,myStepsDone) |
---|
896 | else |
---|
897 | if(nIatt>0) then |
---|
898 | call inits_(aC,iList=AttrVect_exIL2c(aV),lsize=aC_lsize, & |
---|
899 | steps_done=myStepsDone) |
---|
900 | endif |
---|
901 | if(nRatt>0) then |
---|
902 | call inits_(aC,rList=AttrVect_exRL2c(aV),lsize=aC_lsize, & |
---|
903 | steps_done=myStepsDone) |
---|
904 | endif |
---|
905 | endif |
---|
906 | |
---|
907 | |
---|
908 | ! Check that aC as been properly initialized |
---|
909 | |
---|
910 | status = initialized(aC=aC,die_flag=.true.,source_name=myname_) |
---|
911 | |
---|
912 | end subroutine initavs_ |
---|
913 | |
---|
914 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
915 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
916 | !BOP ------------------------------------------------------------------- |
---|
917 | ! |
---|
918 | ! !IROUTINE: clean_ - Destroy an Accumulator |
---|
919 | ! |
---|
920 | ! !DESCRIPTION: |
---|
921 | ! This routine deallocates all allocated memory structures associated |
---|
922 | ! with the input/output {\tt Accumulator} argument {\tt aC}. The |
---|
923 | ! success (failure) of this operation is signified by the zero (non-zero) |
---|
924 | ! value of the optional {\tt INTEGER} output argument {\tt stat}. If |
---|
925 | ! {\tt clean\_()} is invoked with {\tt stat} present, it is the user's |
---|
926 | ! obligation to check this return code and act accordingly. If {\tt stat} |
---|
927 | ! is not supplied and any of the deallocation operations fail, this |
---|
928 | ! routine will terminate execution with an error statement. |
---|
929 | ! |
---|
930 | ! !INTERFACE: |
---|
931 | |
---|
932 | subroutine clean_(aC, stat) |
---|
933 | ! |
---|
934 | ! !USES: |
---|
935 | ! |
---|
936 | use m_mall |
---|
937 | use m_stdio |
---|
938 | use m_die |
---|
939 | use m_AttrVect, only : AttrVect_clean => clean |
---|
940 | |
---|
941 | implicit none |
---|
942 | |
---|
943 | ! !INPUT/OUTPUT PARAMETERS: |
---|
944 | ! |
---|
945 | type(Accumulator), intent(inout) :: aC |
---|
946 | |
---|
947 | ! !OUTPUT PARAMETERS: |
---|
948 | ! |
---|
949 | integer, optional, intent(out) :: stat |
---|
950 | |
---|
951 | ! !REVISION HISTORY: |
---|
952 | ! 11Sep00 - Jay Larson <larson@mcs.anl.gov> - initial prototype |
---|
953 | ! 27JUL01 - E.T. Ong <eong@mcs.anl.gov> - deallocate pointers iAction |
---|
954 | ! and rAction. |
---|
955 | ! 1Mar02 - E.T. Ong <eong@mcs.anl.gov> removed the die to prevent |
---|
956 | ! crashes and added stat argument. |
---|
957 | !EOP ___________________________________________________________________ |
---|
958 | |
---|
959 | character(len=*),parameter :: myname_=myname//'::clean_' |
---|
960 | integer :: ier |
---|
961 | |
---|
962 | if(present(stat)) then |
---|
963 | stat=0 |
---|
964 | call AttrVect_clean(aC%data,stat) |
---|
965 | else |
---|
966 | call AttrVect_clean(aC%data) |
---|
967 | endif |
---|
968 | |
---|
969 | if( associated(aC%iAction) ) then |
---|
970 | |
---|
971 | deallocate(aC%iAction,stat=ier) |
---|
972 | |
---|
973 | if(ier /= 0) then |
---|
974 | if(present(stat)) then |
---|
975 | stat=ier |
---|
976 | else |
---|
977 | call warn(myname_,'deallocate(aC%iAction)',ier) |
---|
978 | endif |
---|
979 | endif |
---|
980 | |
---|
981 | endif |
---|
982 | |
---|
983 | if( associated(aC%rAction) ) then |
---|
984 | |
---|
985 | deallocate(aC%rAction,stat=ier) |
---|
986 | |
---|
987 | if(ier /= 0) then |
---|
988 | if(present(stat)) then |
---|
989 | stat=ier |
---|
990 | else |
---|
991 | call warn(myname_,'deallocate(aC%rAction)',ier) |
---|
992 | endif |
---|
993 | endif |
---|
994 | |
---|
995 | endif |
---|
996 | |
---|
997 | end subroutine clean_ |
---|
998 | |
---|
999 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1000 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1001 | !BOP ------------------------------------------------------------------- |
---|
1002 | ! |
---|
1003 | ! !IROUTINE: initialized_ - Check if an Accumulator is Initialized |
---|
1004 | ! |
---|
1005 | ! !DESCRIPTION: |
---|
1006 | ! This logical function returns a value of {\tt .TRUE.} if the input |
---|
1007 | ! {\tt Accumulator} argument {\tt aC} is initialized correctly. The |
---|
1008 | ! term "correctly initialized" means there is internal consistency |
---|
1009 | ! between the number of integer and real attributes in {\tt aC}, and |
---|
1010 | ! their respective data structures for accumulation registers, and |
---|
1011 | ! accumulation action flags. The optional {\tt LOGICAL} input argument |
---|
1012 | ! {\tt die\_flag} if present, can result in messages written to |
---|
1013 | ! {\tt stderr}: |
---|
1014 | ! \begin {itemize} |
---|
1015 | ! \item if {\tt die\_flag} is true and {\tt aC} is correctly initialized, |
---|
1016 | ! and |
---|
1017 | ! \item if {\tt die\_flag} is false and {\tt aC} is incorrectly |
---|
1018 | ! initialized. |
---|
1019 | ! \end{itemize} |
---|
1020 | ! Otherwise, inconsistencies in how {\tt aC} is set up will result in |
---|
1021 | ! termination with an error message. |
---|
1022 | ! The optional {\tt CHARACTER} input argument {\tt source\_name} allows |
---|
1023 | ! the user to, in the event of error, generate traceback information |
---|
1024 | ! (e.g., the name of the routine that invoked this one). |
---|
1025 | ! |
---|
1026 | ! !INTERFACE: |
---|
1027 | |
---|
1028 | logical function initialized_(aC, die_flag, source_name) |
---|
1029 | ! |
---|
1030 | ! !USES: |
---|
1031 | ! |
---|
1032 | |
---|
1033 | use m_stdio |
---|
1034 | use m_die |
---|
1035 | use m_List, only : List |
---|
1036 | use m_List, only : List_allocated => allocated |
---|
1037 | |
---|
1038 | use m_AttrVect, only : AttrVect |
---|
1039 | use m_AttrVect, only : Attr_nIAttr => nIAttr |
---|
1040 | use m_AttrVect, only : Attr_nRAttr => nRAttr |
---|
1041 | |
---|
1042 | implicit none |
---|
1043 | |
---|
1044 | ! !INPUT PARAMETERS: |
---|
1045 | ! |
---|
1046 | type(Accumulator), intent(in) :: aC |
---|
1047 | logical, optional, intent(in) :: die_flag |
---|
1048 | character(len=*), optional, intent(in) :: source_name |
---|
1049 | |
---|
1050 | ! !REVISION HISTORY: |
---|
1051 | ! 7AUG01 - E.T. Ong <eong@mcs.anl.gov> - initital prototype |
---|
1052 | ! |
---|
1053 | !EOP ___________________________________________________________________ |
---|
1054 | |
---|
1055 | character(len=*),parameter :: myname_=myname//'::initialized_' |
---|
1056 | integer :: i |
---|
1057 | logical :: kill |
---|
1058 | logical :: aC_associated |
---|
1059 | |
---|
1060 | if(present(die_flag)) then |
---|
1061 | kill = .true. |
---|
1062 | else |
---|
1063 | kill = .false. |
---|
1064 | endif |
---|
1065 | |
---|
1066 | ! Initial value |
---|
1067 | initialized_ = .true. |
---|
1068 | aC_associated = .true. |
---|
1069 | |
---|
1070 | ! Check the association status of pointers in aC |
---|
1071 | |
---|
1072 | if( associated(aC%iAction) .or. associated(aC%rAction) ) then |
---|
1073 | aC_associated = .true. |
---|
1074 | else |
---|
1075 | initialized_ = .false. |
---|
1076 | aC_associated = .false. |
---|
1077 | if(kill) then |
---|
1078 | if(present(source_name)) write(stderr,*) source_name, myname_, & |
---|
1079 | ":: ERROR, Neither aC%iAction nor aC%rAction are associated" |
---|
1080 | call die(myname_,"Neither aC%iAction nor aC%rAction are associated") |
---|
1081 | endif |
---|
1082 | endif |
---|
1083 | |
---|
1084 | if( List_allocated(aC%data%iList) .or. List_allocated(aC%data%rList) ) then |
---|
1085 | aC_associated = .true. |
---|
1086 | else |
---|
1087 | initialized_ = .false. |
---|
1088 | aC_associated = .false. |
---|
1089 | if(kill) then |
---|
1090 | if(present(source_name)) write(stderr,*) source_name, myname_, & |
---|
1091 | ":: ERROR, Neither aC%data%iList nor aC%data%rList are allocated" |
---|
1092 | call die(myname_,"Neither aC%data%iList nor aC%data%rList are allocated") |
---|
1093 | endif |
---|
1094 | endif |
---|
1095 | |
---|
1096 | ! Make sure iAction and rAction sizes are greater than zero |
---|
1097 | |
---|
1098 | if(associated(aC%iAction)) then |
---|
1099 | if(size(aC%iAction)<=0) then |
---|
1100 | initialized_ = .false. |
---|
1101 | aC_associated = .false. |
---|
1102 | if(kill) then |
---|
1103 | if(present(source_name)) write(stderr,*) source_name, myname_, & |
---|
1104 | ":: ERROR, size(aC%iAction<=0), size = ", size(aC%iAction) |
---|
1105 | call die(myname_,"size(aC%iAction<=0), size = ", size(aC%iAction)) |
---|
1106 | endif |
---|
1107 | endif |
---|
1108 | endif |
---|
1109 | |
---|
1110 | if(associated(aC%rAction)) then |
---|
1111 | if(size(aC%rAction)<=0) then |
---|
1112 | initialized_ = .false. |
---|
1113 | aC_associated = .false. |
---|
1114 | if(kill) then |
---|
1115 | if(present(source_name)) write(stderr,*) source_name, myname_, & |
---|
1116 | ":: ERROR, size(aC%rAction<=0), size = ", size(aC%rAction) |
---|
1117 | call die(myname_,"size(aC%rAction<=0), size = ", size(aC%rAction)) |
---|
1118 | endif |
---|
1119 | endif |
---|
1120 | endif |
---|
1121 | |
---|
1122 | ! More sanity checking... |
---|
1123 | |
---|
1124 | if( aC_associated ) then |
---|
1125 | |
---|
1126 | if( (Attr_nIAttr(aC%data) == 0) .and. (Attr_nRAttr(aC%data) == 0) ) then |
---|
1127 | initialized_ = .false. |
---|
1128 | if(kill) then |
---|
1129 | if(present(source_name)) write(stderr,*) source_name, myname_, & |
---|
1130 | ":: ERROR, No attributes found in aC%data" |
---|
1131 | call die(myname_,"No attributes found in aC%data") |
---|
1132 | endif |
---|
1133 | endif |
---|
1134 | |
---|
1135 | if(Attr_nIAttr(aC%data) > 0) then |
---|
1136 | |
---|
1137 | if( size(aC%iAction) /= Attr_nIAttr(aC%data) ) then |
---|
1138 | initialized_ = .false. |
---|
1139 | if(kill) then |
---|
1140 | if(present(source_name)) write(stderr,*) source_name, myname_, & |
---|
1141 | ":: ERROR, size(aC%iAction) /= nIAttr(aC%data)" |
---|
1142 | call die(myname_,"size(aC%iAction) /= nIAttr(aC%data)") |
---|
1143 | endif |
---|
1144 | endif |
---|
1145 | |
---|
1146 | do i=1,Attr_nIAttr(aC%data) |
---|
1147 | if( (aC%iAction(i) /= MCT_SUM) .and. & |
---|
1148 | (aC%iAction(i) /= MCT_AVG) ) then |
---|
1149 | initialized_ = .false. |
---|
1150 | if(kill) then |
---|
1151 | if(present(source_name)) write(stderr,*) source_name, & |
---|
1152 | myname_, ":: ERROR, Invalid value found in aC%iAction" |
---|
1153 | call die(myname_,"Invalid value found in aC%iAction", & |
---|
1154 | aC%iAction(i)) |
---|
1155 | endif |
---|
1156 | endif |
---|
1157 | enddo |
---|
1158 | |
---|
1159 | endif ! if(Attr_nIAttr(aC%data) > 0) |
---|
1160 | |
---|
1161 | if(Attr_nRAttr(aC%data) > 0) then |
---|
1162 | |
---|
1163 | if( size(aC%rAction) /= Attr_nRAttr(aC%data) ) then |
---|
1164 | initialized_ = .false. |
---|
1165 | if(kill) then |
---|
1166 | if(present(source_name)) write(stderr,*) source_name, & |
---|
1167 | myname_, ":: ERROR, size(aC%rAction) /= nRAttr(aC%data)" |
---|
1168 | call die(myname_,"size(aC%rAction) /= nRAttr(aC%data)") |
---|
1169 | endif |
---|
1170 | endif |
---|
1171 | |
---|
1172 | do i=1,Attr_nRAttr(aC%data) |
---|
1173 | if( (aC%rAction(i) /= MCT_SUM) .and. & |
---|
1174 | (aC%rAction(i) /= MCT_AVG) ) then |
---|
1175 | initialized_ = .false. |
---|
1176 | if(kill) then |
---|
1177 | if(present(source_name)) write(stderr,*) source_name, & |
---|
1178 | myname_, ":: ERROR, Invalid value found in aC%rAction", & |
---|
1179 | aC%rAction(i) |
---|
1180 | call die(myname_,"Invalid value found in aC%rAction", & |
---|
1181 | aC%iAction(i)) |
---|
1182 | endif |
---|
1183 | endif |
---|
1184 | enddo |
---|
1185 | |
---|
1186 | endif ! if(Attr_nRAttr(aC%data) > 0) |
---|
1187 | |
---|
1188 | endif ! if (aC_associated) |
---|
1189 | |
---|
1190 | end function initialized_ |
---|
1191 | |
---|
1192 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1193 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1194 | !BOP ------------------------------------------------------------------- |
---|
1195 | ! |
---|
1196 | ! !IROUTINE: lsize_ - Length of an Accumulator |
---|
1197 | ! |
---|
1198 | ! !DESCRIPTION: |
---|
1199 | ! This {\tt INTEGER} query function returns the number of data points |
---|
1200 | ! for which the input {\tt Accumulator} argument {\tt aC} is performing |
---|
1201 | ! accumulation. This value corresponds to the length of the {\tt AttrVect} |
---|
1202 | ! component {\tt aC\%data} that stores the accumulation registers. |
---|
1203 | ! |
---|
1204 | ! !INTERFACE: |
---|
1205 | |
---|
1206 | integer function lsize_(aC) |
---|
1207 | ! |
---|
1208 | ! !USES: |
---|
1209 | ! |
---|
1210 | use m_AttrVect, only : AttrVect_lsize => lsize |
---|
1211 | |
---|
1212 | implicit none |
---|
1213 | |
---|
1214 | ! !INPUT PARAMETERS: |
---|
1215 | ! |
---|
1216 | type(Accumulator), intent(in) :: aC |
---|
1217 | |
---|
1218 | ! !REVISION HISTORY: |
---|
1219 | ! 12Sep00 - Jay Larson <larson@mcs.anl.gov> - initial prototype |
---|
1220 | !EOP ___________________________________________________________________ |
---|
1221 | |
---|
1222 | character(len=*),parameter :: myname_=myname//'::lsize_' |
---|
1223 | |
---|
1224 | |
---|
1225 | ! The function AttrVect_lsize is called to return |
---|
1226 | ! its local size data |
---|
1227 | |
---|
1228 | lsize_=AttrVect_lsize(aC%data) |
---|
1229 | |
---|
1230 | end function lsize_ |
---|
1231 | |
---|
1232 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1233 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1234 | !BOP ------------------------------------------------------------------- |
---|
1235 | ! |
---|
1236 | ! !IROUTINE: NumSteps_ - Number of Accumulation Cycle Time Steps |
---|
1237 | ! |
---|
1238 | ! !DESCRIPTION: |
---|
1239 | ! This {\tt INTEGER} query function returns the number of time steps in an |
---|
1240 | ! accumulation cycle for the input {\tt Accumulator} argument {\tt aC}. |
---|
1241 | ! |
---|
1242 | ! !INTERFACE: |
---|
1243 | |
---|
1244 | integer function NumSteps_(aC) |
---|
1245 | ! |
---|
1246 | ! !USES: |
---|
1247 | ! |
---|
1248 | use m_die, only : die |
---|
1249 | use m_stdio, only : stderr |
---|
1250 | |
---|
1251 | implicit none |
---|
1252 | |
---|
1253 | ! !INPUT PARAMETERS: |
---|
1254 | ! |
---|
1255 | type(Accumulator), intent(in) :: aC |
---|
1256 | |
---|
1257 | ! !REVISION HISTORY: |
---|
1258 | ! 7Aug02 - Jay Larson <larson@mcs.anl.gov> - initial prototype |
---|
1259 | !EOP ___________________________________________________________________ |
---|
1260 | |
---|
1261 | character(len=*),parameter :: myname_=myname//'::NumSteps_' |
---|
1262 | |
---|
1263 | integer :: myNumSteps |
---|
1264 | |
---|
1265 | |
---|
1266 | ! Retrieve the number of cycle steps from aC: |
---|
1267 | |
---|
1268 | myNumSteps = aC%num_steps |
---|
1269 | |
---|
1270 | if(myNumSteps <= 0) then |
---|
1271 | write(stderr,'(2a,i8)') myname_, & |
---|
1272 | ':: FATAL--illegal number of steps in an accumulation cycle = ',& |
---|
1273 | myNumSteps |
---|
1274 | call die(myname_) |
---|
1275 | endif |
---|
1276 | |
---|
1277 | NumSteps_ = myNumSteps |
---|
1278 | |
---|
1279 | end function NumSteps_ |
---|
1280 | |
---|
1281 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1282 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1283 | !BOP ------------------------------------------------------------------- |
---|
1284 | ! |
---|
1285 | ! !IROUTINE: StepsDone_ - Number of Completed Steps in the Current Cycle |
---|
1286 | ! |
---|
1287 | ! !DESCRIPTION: |
---|
1288 | ! This {\tt INTEGER} query function returns the of time steps that have |
---|
1289 | ! been completed in the current accumulation cycle for the input |
---|
1290 | ! {\tt Accumulator} argument {\tt aC}. |
---|
1291 | ! |
---|
1292 | ! !INTERFACE: |
---|
1293 | |
---|
1294 | integer function StepsDone_(aC) |
---|
1295 | ! |
---|
1296 | ! !USES: |
---|
1297 | ! |
---|
1298 | use m_die, only : die |
---|
1299 | use m_stdio, only : stderr |
---|
1300 | |
---|
1301 | implicit none |
---|
1302 | |
---|
1303 | ! !INPUT PARAMETERS: |
---|
1304 | ! |
---|
1305 | type(Accumulator), intent(in) :: aC |
---|
1306 | |
---|
1307 | ! !REVISION HISTORY: |
---|
1308 | ! 7Aug02 - Jay Larson <larson@mcs.anl.gov> - initial prototype |
---|
1309 | !EOP ___________________________________________________________________ |
---|
1310 | |
---|
1311 | character(len=*),parameter :: myname_=myname//'::StepsDone_' |
---|
1312 | |
---|
1313 | integer :: myStepsDone |
---|
1314 | |
---|
1315 | ! Retrieve the number of completed steps from aC: |
---|
1316 | |
---|
1317 | myStepsDone = aC%steps_done |
---|
1318 | |
---|
1319 | if(myStepsDone < 0) then |
---|
1320 | write(stderr,'(2a,i8)') myname_, & |
---|
1321 | ':: FATAL--illegal number of completed steps = ',& |
---|
1322 | myStepsDone |
---|
1323 | call die(myname_) |
---|
1324 | endif |
---|
1325 | |
---|
1326 | StepsDone_ = myStepsDone |
---|
1327 | |
---|
1328 | end function StepsDone_ |
---|
1329 | |
---|
1330 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1331 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1332 | !BOP ------------------------------------------------------------------- |
---|
1333 | ! |
---|
1334 | ! !IROUTINE: nIAttr_ - Return the Number of INTEGER Attributes |
---|
1335 | ! |
---|
1336 | ! !DESCRIPTION: |
---|
1337 | ! This {\tt INTEGER} query function returns the number of integer |
---|
1338 | ! attributes that are stored in the input {\tt Accumulator} argument |
---|
1339 | ! {\tt aC}. This value is equal to the number of integer attributes |
---|
1340 | ! in the {\tt AttrVect} component {\tt aC\%data} that stores the |
---|
1341 | ! accumulation registers. |
---|
1342 | ! |
---|
1343 | ! !INTERFACE: |
---|
1344 | |
---|
1345 | integer function nIAttr_(aC) |
---|
1346 | ! |
---|
1347 | ! !USES: |
---|
1348 | ! |
---|
1349 | use m_AttrVect, only : AttrVect_nIAttr => nIAttr |
---|
1350 | |
---|
1351 | implicit none |
---|
1352 | |
---|
1353 | ! !INPUT PARAMETERS: |
---|
1354 | ! |
---|
1355 | type(Accumulator),intent(in) :: aC |
---|
1356 | |
---|
1357 | ! !REVISION HISTORY: |
---|
1358 | ! 12Sep00 - Jay Larson <larson@mcs.anl.gov> - initial prototype |
---|
1359 | !EOP ___________________________________________________________________ |
---|
1360 | |
---|
1361 | character(len=*),parameter :: myname_=myname//'::nIAttr_' |
---|
1362 | |
---|
1363 | ! The function AttrVect_nIAttr is called to return the |
---|
1364 | ! number of integer fields |
---|
1365 | |
---|
1366 | nIAttr_=AttrVect_nIAttr(aC%data) |
---|
1367 | |
---|
1368 | end function nIAttr_ |
---|
1369 | |
---|
1370 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1371 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1372 | !BOP ------------------------------------------------------------------- |
---|
1373 | ! |
---|
1374 | ! !IROUTINE: nRAttr_ - number of REAL fields stored in the Accumulator. |
---|
1375 | ! |
---|
1376 | ! !DESCRIPTION: |
---|
1377 | ! This {\tt INTEGER} query function returns the number of real |
---|
1378 | ! attributes that are stored in the input {\tt Accumulator} argument |
---|
1379 | ! {\tt aC}. This value is equal to the number of real attributes |
---|
1380 | ! in the {\tt AttrVect} component {\tt aC\%data} that stores the |
---|
1381 | ! accumulation registers. |
---|
1382 | ! |
---|
1383 | ! !INTERFACE: |
---|
1384 | |
---|
1385 | integer function nRAttr_(aC) |
---|
1386 | ! |
---|
1387 | ! !USES: |
---|
1388 | ! |
---|
1389 | use m_AttrVect, only : AttrVect_nRAttr => nRAttr |
---|
1390 | |
---|
1391 | implicit none |
---|
1392 | |
---|
1393 | ! !INPUT PARAMETERS: |
---|
1394 | ! |
---|
1395 | type(Accumulator),intent(in) :: aC |
---|
1396 | |
---|
1397 | ! !REVISION HISTORY: |
---|
1398 | ! 12Sep00 - Jay Larson <larson@mcs.anl.gov> - initial prototype |
---|
1399 | !EOP ___________________________________________________________________ |
---|
1400 | |
---|
1401 | character(len=*),parameter :: myname_=myname//'::nRAttr_' |
---|
1402 | |
---|
1403 | ! The function AttrVect_nRAttr is called to return the |
---|
1404 | ! number of real fields |
---|
1405 | |
---|
1406 | nRAttr_=AttrVect_nRAttr(aC%data) |
---|
1407 | |
---|
1408 | end function nRAttr_ |
---|
1409 | |
---|
1410 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1411 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1412 | !BOP ------------------------------------------------------------------- |
---|
1413 | ! |
---|
1414 | ! !IROUTINE: getIList_ - Retrieve a Numbered INTEGER Attribute Name |
---|
1415 | ! |
---|
1416 | ! !DESCRIPTION: |
---|
1417 | ! This routine returns as a {\tt String} (see the mpeu module |
---|
1418 | ! {\tt m\_String} for information) the name of the {\tt ith} item in |
---|
1419 | ! the integer registers of the {\tt Accumulator} argument {\tt aC}. |
---|
1420 | ! |
---|
1421 | ! !INTERFACE: |
---|
1422 | |
---|
1423 | subroutine getIList_(item, ith, aC) |
---|
1424 | ! |
---|
1425 | ! !USES: |
---|
1426 | ! |
---|
1427 | use m_AttrVect, only : AttrVect_getIList => getIList |
---|
1428 | use m_String, only : String |
---|
1429 | |
---|
1430 | implicit none |
---|
1431 | |
---|
1432 | ! !INPUT PARAMETERS: |
---|
1433 | ! |
---|
1434 | integer, intent(in) :: ith |
---|
1435 | type(Accumulator), intent(in) :: aC |
---|
1436 | |
---|
1437 | ! !OUTPUT PARAMETERS: |
---|
1438 | ! |
---|
1439 | type(String), intent(out) :: item |
---|
1440 | |
---|
1441 | ! !REVISION HISTORY: |
---|
1442 | ! 12Sep00 - Jay Larson <larson@mcs.anl.gov> - initial prototype |
---|
1443 | !EOP ___________________________________________________________________ |
---|
1444 | |
---|
1445 | character(len=*),parameter :: myname_=myname//'::getIList_' |
---|
1446 | |
---|
1447 | call AttrVect_getIList(item,ith,aC%data) |
---|
1448 | |
---|
1449 | end subroutine getIList_ |
---|
1450 | |
---|
1451 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1452 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1453 | !BOP ------------------------------------------------------------------- |
---|
1454 | ! |
---|
1455 | ! !IROUTINE: getRList_ - Retrieve a Numbered REAL Attribute Name |
---|
1456 | ! |
---|
1457 | ! !DESCRIPTION: |
---|
1458 | ! This routine returns as a {\tt String} (see the mpeu module |
---|
1459 | ! {\tt m\_String} for information) the name of the {\tt ith} item in |
---|
1460 | ! the real registers of the {\tt Accumulator} argument {\tt aC}. |
---|
1461 | ! |
---|
1462 | ! !INTERFACE: |
---|
1463 | |
---|
1464 | subroutine getRList_(item, ith, aC) |
---|
1465 | ! |
---|
1466 | ! !USES: |
---|
1467 | ! |
---|
1468 | use m_AttrVect, only : AttrVect_getRList => getRList |
---|
1469 | use m_String, only : String |
---|
1470 | |
---|
1471 | implicit none |
---|
1472 | |
---|
1473 | ! !INPUT PARAMETERS: |
---|
1474 | ! |
---|
1475 | integer, intent(in) :: ith |
---|
1476 | type(Accumulator),intent(in) :: aC |
---|
1477 | |
---|
1478 | ! !OUTPUT PARAMETERS: |
---|
1479 | ! |
---|
1480 | type(String), intent(out) :: item |
---|
1481 | |
---|
1482 | ! !REVISION HISTORY: |
---|
1483 | ! 12Sep00 - Jay Larson <larson@mcs.anl.gov> - initial prototype |
---|
1484 | !EOP ___________________________________________________________________ |
---|
1485 | |
---|
1486 | character(len=*),parameter :: myname_=myname//'::getRList_' |
---|
1487 | |
---|
1488 | call AttrVect_getRList(item,ith,aC%data) |
---|
1489 | |
---|
1490 | end subroutine getRList_ |
---|
1491 | |
---|
1492 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1493 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1494 | !BOP ------------------------------------------------------------------- |
---|
1495 | ! |
---|
1496 | ! !IROUTINE: indexIA_ - Index an INTEGER Attribute |
---|
1497 | ! |
---|
1498 | ! !DESCRIPTION: |
---|
1499 | ! This {\tt INTEGER} query function returns the index in the integer |
---|
1500 | ! accumulation register buffer of the {\tt Accumulator} argument {\tt aC} |
---|
1501 | ! the attribute named by the {\tt CHARACTER} argument {\tt item}. That |
---|
1502 | ! is, all the accumulator running tallies for the attribute named |
---|
1503 | ! {\tt item} reside in |
---|
1504 | !\begin{verbatim} |
---|
1505 | ! aC%data%iAttr(indexIA_(aC,item),:). |
---|
1506 | !\end{verbatim} |
---|
1507 | ! The user may request traceback information (e.g., the name of the |
---|
1508 | ! routine from which this one is called) by providing values for either |
---|
1509 | ! of the optional {\tt CHARACTER} arguments {\tt perrWith} or {\tt dieWith} |
---|
1510 | ! In the event {\tt indexIA\_()} can not find {\tt item} in {\tt aC}, |
---|
1511 | ! the routine behaves as follows: |
---|
1512 | ! \begin{enumerate} |
---|
1513 | ! \item if neither {\tt perrWith} nor {\tt dieWith} are present, |
---|
1514 | ! {\tt indexIA\_()} returns a value of zero; |
---|
1515 | ! \item if {\tt perrWith} is present, but {\tt dieWith} is not, an error |
---|
1516 | ! message is written to {\tt stderr} incorporating user-supplied traceback |
---|
1517 | ! information stored in the argument {\tt perrWith}; |
---|
1518 | ! \item if {\tt dieWith} is present, execution terminates with an error |
---|
1519 | ! message written to {\tt stderr} that incorporates user-supplied traceback |
---|
1520 | ! information stored in the argument {\tt dieWith}. |
---|
1521 | ! \end{enumerate} |
---|
1522 | ! !INTERFACE: |
---|
1523 | |
---|
1524 | integer function indexIA_(aC, item, perrWith, dieWith) |
---|
1525 | ! |
---|
1526 | ! !USES: |
---|
1527 | ! |
---|
1528 | use m_AttrVect, only : AttrVect_indexIA => indexIA |
---|
1529 | use m_die, only : die |
---|
1530 | use m_stdio,only : stderr |
---|
1531 | |
---|
1532 | implicit none |
---|
1533 | |
---|
1534 | ! !INPUT PARAMETERS: |
---|
1535 | ! |
---|
1536 | type(Accumulator), intent(in) :: aC |
---|
1537 | character(len=*), intent(in) :: item |
---|
1538 | character(len=*), optional, intent(in) :: perrWith |
---|
1539 | character(len=*), optional, intent(in) :: dieWith |
---|
1540 | |
---|
1541 | ! !REVISION HISTORY: |
---|
1542 | ! 14Sep00 - Jay Larson <larson@mcs.anl.gov> - initial prototype |
---|
1543 | !EOP ___________________________________________________________________ |
---|
1544 | |
---|
1545 | character(len=*),parameter :: myname_=myname//'::indexIA_' |
---|
1546 | |
---|
1547 | indexIA_=AttrVect_indexIA(aC%data,item) |
---|
1548 | |
---|
1549 | if(indexIA_==0) then |
---|
1550 | if(.not.present(dieWith)) then |
---|
1551 | if(present(perrWith)) write(stderr,'(4a)') perrWith, & |
---|
1552 | '" indexIA_() error, not found "',trim(item),'"' |
---|
1553 | else |
---|
1554 | write(stderr,'(4a)') dieWith, & |
---|
1555 | '" indexIA_() error, not found "',trim(item),'"' |
---|
1556 | call die(dieWith) |
---|
1557 | endif |
---|
1558 | endif |
---|
1559 | |
---|
1560 | end function indexIA_ |
---|
1561 | |
---|
1562 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1563 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1564 | !BOP ------------------------------------------------------------------- |
---|
1565 | ! |
---|
1566 | ! !IROUTINE: indexRA_ - index the Accumulator real attribute list. |
---|
1567 | ! |
---|
1568 | ! !DESCRIPTION: |
---|
1569 | ! This {\tt INTEGER} query function returns the index in the real |
---|
1570 | ! accumulation register buffer of the {\tt Accumulator} argument {\tt aC} |
---|
1571 | ! the attribute named by the {\tt CHARACTER} argument {\tt item}. That |
---|
1572 | ! is, all the accumulator running tallies for the attribute named |
---|
1573 | ! {\tt item} reside in |
---|
1574 | !\begin{verbatim} |
---|
1575 | ! aC%data%rAttr(indexRA_(aC,item),:). |
---|
1576 | !\end{verbatim} |
---|
1577 | ! The user may request traceback information (e.g., the name of the |
---|
1578 | ! routine from which this one is called) by providing values for either |
---|
1579 | ! of the optional {\tt CHARACTER} arguments {\tt perrWith} or {\tt dieWith} |
---|
1580 | ! In the event {\tt indexRA\_()} can not find {\tt item} in {\tt aC}, |
---|
1581 | ! the routine behaves as follows: |
---|
1582 | ! \begin{enumerate} |
---|
1583 | ! \item if neither {\tt perrWith} nor {\tt dieWith} are present, |
---|
1584 | ! {\tt indexRA\_()} returns a value of zero; |
---|
1585 | ! \item if {\tt perrWith} is present, but {\tt dieWith} is not, an error |
---|
1586 | ! message is written to {\tt stderr} incorporating user-supplied traceback |
---|
1587 | ! information stored in the argument {\tt perrWith}; |
---|
1588 | ! \item if {\tt dieWith} is present, execution terminates with an error |
---|
1589 | ! message written to {\tt stderr} that incorporates user-supplied traceback |
---|
1590 | ! information stored in the argument {\tt dieWith}. |
---|
1591 | ! \end{enumerate} |
---|
1592 | ! |
---|
1593 | ! !INTERFACE: |
---|
1594 | |
---|
1595 | integer function indexRA_(aC, item, perrWith, dieWith) |
---|
1596 | ! |
---|
1597 | ! !USES: |
---|
1598 | ! |
---|
1599 | use m_AttrVect, only : AttrVect_indexRA => indexRA |
---|
1600 | use m_die, only : die |
---|
1601 | use m_stdio,only : stderr |
---|
1602 | |
---|
1603 | implicit none |
---|
1604 | |
---|
1605 | ! !INPUT PARAMETERS: |
---|
1606 | ! |
---|
1607 | type(Accumulator), intent(in) :: aC |
---|
1608 | character(len=*), intent(in) :: item |
---|
1609 | character(len=*), optional, intent(in) :: perrWith |
---|
1610 | character(len=*), optional, intent(in) :: dieWith |
---|
1611 | |
---|
1612 | ! !REVISION HISTORY: |
---|
1613 | ! 14Sep00 - Jay Larson <larson@mcs.anl.gov> - initial prototype |
---|
1614 | !EOP ___________________________________________________________________ |
---|
1615 | |
---|
1616 | character(len=*),parameter :: myname_=myname//'::indexRA_' |
---|
1617 | |
---|
1618 | indexRA_=AttrVect_indexRA(aC%data,item) |
---|
1619 | |
---|
1620 | if(indexRA_==0) then |
---|
1621 | if(.not.present(dieWith)) then |
---|
1622 | if(present(perrWith)) write(stderr,'(4a)') perrWith, & |
---|
1623 | '" indexRA_() error, not found "',trim(item),'"' |
---|
1624 | else |
---|
1625 | write(stderr,'(4a)') dieWith, & |
---|
1626 | '" indexRA_() error, not found "',trim(item),'"' |
---|
1627 | call die(dieWith) |
---|
1628 | endif |
---|
1629 | endif |
---|
1630 | |
---|
1631 | end function indexRA_ |
---|
1632 | |
---|
1633 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1634 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1635 | !BOP ------------------------------------------------------------------- |
---|
1636 | ! |
---|
1637 | ! !IROUTINE: exportIAttr_ - Export INTEGER Attribute to a Vector |
---|
1638 | ! |
---|
1639 | ! !DESCRIPTION: |
---|
1640 | ! This routine extracts from the input {\tt Accumulator} argument |
---|
1641 | ! {\tt aC} the integer attribute corresponding to the tag defined in |
---|
1642 | ! the input {\tt CHARACTER} argument {\tt AttrTag}, and returns it in |
---|
1643 | ! the {\tt INTEGER} output array {\tt outVect}, and its length in the |
---|
1644 | ! output {\tt INTEGER} argument {\tt lsize}. |
---|
1645 | ! |
---|
1646 | ! {\bf N.B.:} This routine will fail if the {\tt AttrTag} is not in |
---|
1647 | ! the {\tt Accumulator} {\tt List} component {\tt aC\%data\%iList}. |
---|
1648 | ! |
---|
1649 | ! {\bf N.B.:} The flexibility of this routine regarding the pointer |
---|
1650 | ! association status of the output argument {\tt outVect} means the |
---|
1651 | ! user must invoke this routine with care. If the user wishes this |
---|
1652 | ! routine to fill a pre-allocated array, then obviously this array |
---|
1653 | ! must be allocated prior to calling this routine. If the user wishes |
---|
1654 | ! that the routine {\em create} the output argument array {\tt outVect}, |
---|
1655 | ! then the user must ensure this pointer is not allocated (i.e. the user |
---|
1656 | ! must nullify this pointer) at the time this routine is invoked. |
---|
1657 | ! |
---|
1658 | ! {\bf N.B.:} If the user has relied on this routine to allocate memory |
---|
1659 | ! associated with the pointer {\tt outVect}, then the user is responsible |
---|
1660 | ! for deallocating this array once it is no longer needed. Failure to |
---|
1661 | ! do so will result in a memory leak. |
---|
1662 | ! |
---|
1663 | ! !INTERFACE: |
---|
1664 | |
---|
1665 | subroutine exportIAttr_(aC, AttrTag, outVect, lsize) |
---|
1666 | ! |
---|
1667 | ! !USES: |
---|
1668 | ! |
---|
1669 | use m_die |
---|
1670 | use m_stdio |
---|
1671 | |
---|
1672 | use m_AttrVect, only : AttrVect_exportIAttr => exportIAttr |
---|
1673 | |
---|
1674 | implicit none |
---|
1675 | |
---|
1676 | ! !INPUT PARAMETERS: |
---|
1677 | |
---|
1678 | type(Accumulator), intent(in) :: aC |
---|
1679 | character(len=*), intent(in) :: AttrTag |
---|
1680 | |
---|
1681 | ! !OUTPUT PARAMETERS: |
---|
1682 | |
---|
1683 | integer, dimension(:), pointer :: outVect |
---|
1684 | integer, optional, intent(out) :: lsize |
---|
1685 | |
---|
1686 | ! !REVISION HISTORY: |
---|
1687 | |
---|
1688 | ! 6May02 - J.W. Larson <larson@mcs.anl.gov> - initial prototype. |
---|
1689 | ! |
---|
1690 | !EOP ___________________________________________________________________ |
---|
1691 | |
---|
1692 | character(len=*),parameter :: myname_=myname//'::exportIAttr_' |
---|
1693 | |
---|
1694 | ! Export the data (inheritance from AttrVect) |
---|
1695 | if(present(lsize)) then |
---|
1696 | call AttrVect_exportIAttr(aC%data, AttrTag, outVect, lsize) |
---|
1697 | else |
---|
1698 | call AttrVect_exportIAttr(aC%data, AttrTag, outVect) |
---|
1699 | endif |
---|
1700 | |
---|
1701 | end subroutine exportIAttr_ |
---|
1702 | |
---|
1703 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1704 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1705 | !BOP ------------------------------------------------------------------- |
---|
1706 | ! |
---|
1707 | ! !IROUTINE: exportRAttrSP_ - Export REAL Attribute to a Vector |
---|
1708 | ! |
---|
1709 | ! !DESCRIPTION: |
---|
1710 | ! This routine extracts from the input {\tt Accumulator} argument |
---|
1711 | ! {\tt aC} the real attribute corresponding to the tag defined in |
---|
1712 | ! the input {\tt CHARACTER} argument {\tt AttrTag}, and returns it in |
---|
1713 | ! the {\tt REAL} output array {\tt outVect}, and its length in the |
---|
1714 | ! output {\tt INTEGER} argument {\tt lsize}. |
---|
1715 | ! |
---|
1716 | ! {\bf N.B.:} This routine will fail if the {\tt AttrTag} is not in |
---|
1717 | ! the {\tt Accumulator} {\tt List} component {\tt aC\%data\%iList}. |
---|
1718 | ! |
---|
1719 | ! {\bf N.B.:} The flexibility of this routine regarding the pointer |
---|
1720 | ! association status of the output argument {\tt outVect} means the |
---|
1721 | ! user must invoke this routine with care. If the user wishes this |
---|
1722 | ! routine to fill a pre-allocated array, then obviously this array |
---|
1723 | ! must be allocated prior to calling this routine. If the user wishes |
---|
1724 | ! that the routine {\em create} the output argument array {\tt outVect}, |
---|
1725 | ! then the user must ensure this pointer is not allocated (i.e. the user |
---|
1726 | ! must nullify this pointer) at the time this routine is invoked. |
---|
1727 | ! |
---|
1728 | ! {\bf N.B.:} If the user has relied on this routine to allocate memory |
---|
1729 | ! associated with the pointer {\tt outVect}, then the user is responsible |
---|
1730 | ! for deallocating this array once it is no longer needed. Failure to |
---|
1731 | ! do so will result in a memory leak. |
---|
1732 | ! |
---|
1733 | ! !INTERFACE: |
---|
1734 | |
---|
1735 | subroutine exportRAttrSP_(aC, AttrTag, outVect, lsize) |
---|
1736 | ! |
---|
1737 | ! !USES: |
---|
1738 | ! |
---|
1739 | use m_die |
---|
1740 | use m_stdio |
---|
1741 | |
---|
1742 | use m_AttrVect, only : AttrVect_exportRAttr => exportRAttr |
---|
1743 | |
---|
1744 | implicit none |
---|
1745 | |
---|
1746 | ! !INPUT PARAMETERS: |
---|
1747 | |
---|
1748 | type(Accumulator), intent(in) :: aC |
---|
1749 | character(len=*), intent(in) :: AttrTag |
---|
1750 | |
---|
1751 | ! !OUTPUT PARAMETERS: |
---|
1752 | |
---|
1753 | real(SP), dimension(:), pointer :: outVect |
---|
1754 | integer, optional, intent(out) :: lsize |
---|
1755 | |
---|
1756 | ! !REVISION HISTORY: |
---|
1757 | ! 6May02 - J.W. Larson <larson@mcs.anl.gov> - initial prototype. |
---|
1758 | ! |
---|
1759 | !EOP ___________________________________________________________________ |
---|
1760 | |
---|
1761 | character(len=*),parameter :: myname_=myname//'::exportRAttrSP_' |
---|
1762 | |
---|
1763 | ! Export the data (inheritance from AttrVect) |
---|
1764 | |
---|
1765 | if(present(lsize)) then |
---|
1766 | call AttrVect_exportRAttr(aC%data, AttrTag, outVect, lsize) |
---|
1767 | else |
---|
1768 | call AttrVect_exportRAttr(aC%data, AttrTag, outVect) |
---|
1769 | endif |
---|
1770 | |
---|
1771 | end subroutine exportRAttrSP_ |
---|
1772 | |
---|
1773 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1774 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1775 | ! ---------------------------------------------------------------------- |
---|
1776 | ! |
---|
1777 | ! !IROUTINE: exportRAttrDP_ - Export REAL Attribute to a Vector |
---|
1778 | ! |
---|
1779 | ! !DESCRIPTION: |
---|
1780 | ! Double precision version of exportRAttrSP_ |
---|
1781 | ! |
---|
1782 | ! !INTERFACE: |
---|
1783 | |
---|
1784 | subroutine exportRAttrDP_(aC, AttrTag, outVect, lsize) |
---|
1785 | ! |
---|
1786 | ! !USES: |
---|
1787 | ! |
---|
1788 | use m_die |
---|
1789 | use m_stdio |
---|
1790 | |
---|
1791 | use m_AttrVect, only : AttrVect_exportRAttr => exportRAttr |
---|
1792 | |
---|
1793 | implicit none |
---|
1794 | |
---|
1795 | ! !INPUT PARAMETERS: |
---|
1796 | |
---|
1797 | type(Accumulator), intent(in) :: aC |
---|
1798 | character(len=*), intent(in) :: AttrTag |
---|
1799 | |
---|
1800 | ! !OUTPUT PARAMETERS: |
---|
1801 | |
---|
1802 | real(DP), dimension(:), pointer :: outVect |
---|
1803 | integer, optional, intent(out) :: lsize |
---|
1804 | |
---|
1805 | ! !REVISION HISTORY: |
---|
1806 | ! 6May02 - J.W. Larson <larson@mcs.anl.gov> - initial prototype. |
---|
1807 | ! |
---|
1808 | ! ______________________________________________________________________ |
---|
1809 | |
---|
1810 | character(len=*),parameter :: myname_=myname//'::exportRAttrDP_' |
---|
1811 | |
---|
1812 | ! Export the data (inheritance from AttrVect) |
---|
1813 | |
---|
1814 | if(present(lsize)) then |
---|
1815 | call AttrVect_exportRAttr(aC%data, AttrTag, outVect, lsize) |
---|
1816 | else |
---|
1817 | call AttrVect_exportRAttr(aC%data, AttrTag, outVect) |
---|
1818 | endif |
---|
1819 | |
---|
1820 | end subroutine exportRAttrDP_ |
---|
1821 | |
---|
1822 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1823 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1824 | !BOP ------------------------------------------------------------------- |
---|
1825 | ! |
---|
1826 | ! !IROUTINE: importIAttr_ - Import INTEGER Attribute from a Vector |
---|
1827 | ! |
---|
1828 | ! !DESCRIPTION: |
---|
1829 | ! This routine imports data provided in the input {\tt INTEGER} vector |
---|
1830 | ! {\tt inVect} into the {\tt Accumulator} argument {\tt aC}, storing |
---|
1831 | ! it as the integer attribute corresponding to the tag defined in |
---|
1832 | ! the input {\tt CHARACTER} argument {\tt AttrTag}. The input |
---|
1833 | ! {\tt INTEGER} argument {\tt lsize} is used to ensure there is |
---|
1834 | ! sufficient space in the {\tt Accumulator} to store the data. |
---|
1835 | ! |
---|
1836 | ! {\bf N.B.:} This routine will fail if the {\tt AttrTag} is not in |
---|
1837 | ! the {\tt Accumulator} {\tt List} component {\tt aC\%data\%rList}. |
---|
1838 | ! |
---|
1839 | ! !INTERFACE: |
---|
1840 | |
---|
1841 | subroutine importIAttr_(aC, AttrTag, inVect, lsize) |
---|
1842 | ! |
---|
1843 | ! !USES: |
---|
1844 | ! |
---|
1845 | use m_die |
---|
1846 | use m_stdio , only : stderr |
---|
1847 | |
---|
1848 | use m_AttrVect, only : AttrVect_importIAttr => importIAttr |
---|
1849 | |
---|
1850 | implicit none |
---|
1851 | |
---|
1852 | ! !INPUT PARAMETERS: |
---|
1853 | |
---|
1854 | character(len=*), intent(in) :: AttrTag |
---|
1855 | integer, dimension(:), pointer :: inVect |
---|
1856 | integer, intent(in) :: lsize |
---|
1857 | |
---|
1858 | ! !INPUT/OUTPUT PARAMETERS: |
---|
1859 | |
---|
1860 | type(Accumulator), intent(inout) :: aC |
---|
1861 | |
---|
1862 | ! !REVISION HISTORY: |
---|
1863 | ! 6May02 - J.W. Larson <larson@mcs.anl.gov> - initial prototype. |
---|
1864 | !EOP ___________________________________________________________________ |
---|
1865 | |
---|
1866 | character(len=*),parameter :: myname_=myname//'::importIAttr_' |
---|
1867 | |
---|
1868 | ! Argument Check: |
---|
1869 | |
---|
1870 | if(lsize > lsize_(aC)) then |
---|
1871 | write(stderr,*) myname_,':: ERROR, lsize > lsize_(aC).', & |
---|
1872 | 'lsize = ',lsize,'lsize_(aC) = ',lsize_(ac) |
---|
1873 | call die(myname_) |
---|
1874 | endif |
---|
1875 | |
---|
1876 | ! Import the data (inheritance from AttrVect) |
---|
1877 | |
---|
1878 | call AttrVect_importIAttr(aC%data, AttrTag, inVect, lsize) |
---|
1879 | |
---|
1880 | end subroutine importIAttr_ |
---|
1881 | |
---|
1882 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1883 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1884 | !BOP ------------------------------------------------------------------- |
---|
1885 | ! |
---|
1886 | ! !IROUTINE: importRAttrSP_ - Import REAL Attribute from a Vector |
---|
1887 | ! |
---|
1888 | ! !DESCRIPTION: |
---|
1889 | ! This routine imports data provided in the input {\tt REAL} vector |
---|
1890 | ! {\tt inVect} into the {\tt Accumulator} argument {\tt aC}, storing |
---|
1891 | ! it as the real attribute corresponding to the tag defined in |
---|
1892 | ! the input {\tt CHARACTER} argument {\tt AttrTag}. The input |
---|
1893 | ! {\tt INTEGER} argument {\tt lsize} is used to ensure there is |
---|
1894 | ! sufficient space in the {\tt Accumulator} to store the data. |
---|
1895 | ! |
---|
1896 | ! {\bf N.B.:} This routine will fail if the {\tt AttrTag} is not in |
---|
1897 | ! the {\tt Accumulator} {\tt List} component {\tt aC\%data\%rList}. |
---|
1898 | ! |
---|
1899 | ! !INTERFACE: |
---|
1900 | |
---|
1901 | subroutine importRAttrSP_(aC, AttrTag, inVect, lsize) |
---|
1902 | ! |
---|
1903 | ! !USES: |
---|
1904 | ! |
---|
1905 | use m_die |
---|
1906 | use m_stdio , only : stderr |
---|
1907 | |
---|
1908 | use m_AttrVect, only : AttrVect_importRAttr => importRAttr |
---|
1909 | |
---|
1910 | implicit none |
---|
1911 | |
---|
1912 | ! !INPUT PARAMETERS: |
---|
1913 | |
---|
1914 | character(len=*), intent(in) :: AttrTag |
---|
1915 | real(SP), dimension(:), pointer :: inVect |
---|
1916 | integer, intent(in) :: lsize |
---|
1917 | |
---|
1918 | ! !INPUT/OUTPUT PARAMETERS: |
---|
1919 | |
---|
1920 | type(Accumulator), intent(inout) :: aC |
---|
1921 | |
---|
1922 | ! !REVISION HISTORY: |
---|
1923 | ! 6May02 - J.W. Larson <larson@mcs.anl.gov> - initial prototype. |
---|
1924 | !EOP ___________________________________________________________________ |
---|
1925 | |
---|
1926 | character(len=*),parameter :: myname_=myname//'::importRAttrSP_' |
---|
1927 | |
---|
1928 | ! Argument Check: |
---|
1929 | |
---|
1930 | if(lsize > lsize_(aC)) then |
---|
1931 | write(stderr,*) myname_,':: ERROR, lsize > lsize_(aC).', & |
---|
1932 | 'lsize = ',lsize,'lsize_(aC) = ',lsize_(ac) |
---|
1933 | call die(myname_) |
---|
1934 | endif |
---|
1935 | |
---|
1936 | ! Import the data (inheritance from AttrVect) |
---|
1937 | |
---|
1938 | call AttrVect_importRAttr(aC%data, AttrTag, inVect, lsize) |
---|
1939 | |
---|
1940 | end subroutine importRAttrSP_ |
---|
1941 | |
---|
1942 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1943 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1944 | ! ---------------------------------------------------------------------- |
---|
1945 | ! |
---|
1946 | ! !IROUTINE: importRAttrDP_ - Import REAL Attribute from a Vector |
---|
1947 | ! |
---|
1948 | ! !DESCRIPTION: |
---|
1949 | ! Double precision version of importRAttrSP_ |
---|
1950 | ! |
---|
1951 | ! !INTERFACE: |
---|
1952 | |
---|
1953 | subroutine importRAttrDP_(aC, AttrTag, inVect, lsize) |
---|
1954 | ! |
---|
1955 | ! !USES: |
---|
1956 | ! |
---|
1957 | use m_die |
---|
1958 | use m_stdio , only : stderr |
---|
1959 | |
---|
1960 | use m_AttrVect, only : AttrVect_importRAttr => importRAttr |
---|
1961 | |
---|
1962 | implicit none |
---|
1963 | |
---|
1964 | ! !INPUT PARAMETERS: |
---|
1965 | |
---|
1966 | character(len=*), intent(in) :: AttrTag |
---|
1967 | real(DP), dimension(:), pointer :: inVect |
---|
1968 | integer, intent(in) :: lsize |
---|
1969 | |
---|
1970 | ! !INPUT/OUTPUT PARAMETERS: |
---|
1971 | |
---|
1972 | type(Accumulator), intent(inout) :: aC |
---|
1973 | |
---|
1974 | ! !REVISION HISTORY: |
---|
1975 | ! 6May02 - J.W. Larson <larson@mcs.anl.gov> - initial prototype. |
---|
1976 | ! ______________________________________________________________________ |
---|
1977 | |
---|
1978 | character(len=*),parameter :: myname_=myname//'::importRAttrDP_' |
---|
1979 | |
---|
1980 | ! Argument Check: |
---|
1981 | |
---|
1982 | if(lsize > lsize_(aC)) then |
---|
1983 | write(stderr,*) myname_,':: ERROR, lsize > lsize_(aC).', & |
---|
1984 | 'lsize = ',lsize,'lsize_(aC) = ',lsize_(ac) |
---|
1985 | call die(myname_) |
---|
1986 | endif |
---|
1987 | |
---|
1988 | ! Import the data (inheritance from AttrVect) |
---|
1989 | |
---|
1990 | call AttrVect_importRAttr(aC%data, AttrTag, inVect, lsize) |
---|
1991 | |
---|
1992 | end subroutine importRAttrDP_ |
---|
1993 | |
---|
1994 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
1995 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
1996 | !BOP ------------------------------------------------------------------- |
---|
1997 | ! |
---|
1998 | ! !IROUTINE: zero_ - Zero an Accumulator |
---|
1999 | ! |
---|
2000 | ! !DESCRIPTION: |
---|
2001 | ! This subroutine clears the the {\tt Accumulator} argument {\tt aC}. |
---|
2002 | ! This is accomplished by setting the number of completed steps in the |
---|
2003 | ! accumulation cycle to zero, and zeroing out all of the accumlation |
---|
2004 | ! registers. |
---|
2005 | ! |
---|
2006 | ! !INTERFACE: |
---|
2007 | |
---|
2008 | subroutine zero_(aC) |
---|
2009 | ! |
---|
2010 | ! !USES: |
---|
2011 | ! |
---|
2012 | use m_AttrVect, only : AttrVect_zero => zero |
---|
2013 | |
---|
2014 | implicit none |
---|
2015 | |
---|
2016 | ! !INPUT/OUTPUT PARAMETERS: |
---|
2017 | ! |
---|
2018 | type(Accumulator), intent(inout) :: aC |
---|
2019 | |
---|
2020 | ! !REVISION HISTORY: |
---|
2021 | ! 7Aug02 - Jay Larson <larson@mcs.anl.gov> - initial prototype |
---|
2022 | !EOP ___________________________________________________________________ |
---|
2023 | |
---|
2024 | character(len=*),parameter :: myname_=myname//'::zero_' |
---|
2025 | |
---|
2026 | ! Set number of completed cycle steps to zero: |
---|
2027 | |
---|
2028 | aC%steps_done = 0 |
---|
2029 | |
---|
2030 | ! Zero out the accumulation registers: |
---|
2031 | |
---|
2032 | call AttrVect_zero(aC%data) |
---|
2033 | |
---|
2034 | end subroutine zero_ |
---|
2035 | |
---|
2036 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
2037 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
2038 | !BOP ------------------------------------------------------------------- |
---|
2039 | ! |
---|
2040 | ! !IROUTINE: aCaCSharedAttrIndexList_ - Cross-index Two Accumulators |
---|
2041 | ! |
---|
2042 | ! !DESCRIPTION: {\tt aCaCSharedAttrIndexList\_()} takes a pair of |
---|
2043 | ! user-supplied {\tt Accumulator} variables {\tt aC1} and {\tt aC2}, |
---|
2044 | ! and for choice of either {\tt REAL} or {\tt INTEGER} attributes (as |
---|
2045 | ! specified literally in the input {\tt CHARACTER} argument {\tt attrib}) |
---|
2046 | ! returns the number of shared attributes {\tt NumShared}, and arrays of |
---|
2047 | ! indices {\tt Indices1} and {\tt Indices2} to their storage locations |
---|
2048 | ! in {\tt aC1} and {\tt aC2}, respectively. |
---|
2049 | ! |
---|
2050 | ! {\bf N.B.:} This routine returns two allocated arrays---{\tt Indices1(:)} |
---|
2051 | ! and {\tt Indices2(:)}---which must be deallocated once the user no longer |
---|
2052 | ! needs them. Failure to do this will create a memory leak. |
---|
2053 | ! |
---|
2054 | ! !INTERFACE: |
---|
2055 | |
---|
2056 | subroutine aCaCSharedAttrIndexList_(aC1, aC2, attrib, NumShared, & |
---|
2057 | Indices1, Indices2) |
---|
2058 | |
---|
2059 | ! |
---|
2060 | ! !USES: |
---|
2061 | ! |
---|
2062 | use m_stdio |
---|
2063 | use m_die, only : MP_perr_die, die, warn |
---|
2064 | |
---|
2065 | use m_List, only : GetSharedListIndices |
---|
2066 | |
---|
2067 | implicit none |
---|
2068 | |
---|
2069 | ! !INPUT PARAMETERS: |
---|
2070 | ! |
---|
2071 | type(Accumulator), intent(in) :: aC1 |
---|
2072 | type(Accumulator), intent(in) :: aC2 |
---|
2073 | character*7, intent(in) :: attrib |
---|
2074 | |
---|
2075 | ! !OUTPUT PARAMETERS: |
---|
2076 | ! |
---|
2077 | integer, intent(out) :: NumShared |
---|
2078 | integer,dimension(:), pointer :: Indices1 |
---|
2079 | integer,dimension(:), pointer :: Indices2 |
---|
2080 | |
---|
2081 | ! !REVISION HISTORY: |
---|
2082 | ! 7Feb01 - J.W. Larson <larson@mcs.anl.gov> - initial version |
---|
2083 | !EOP ___________________________________________________________________ |
---|
2084 | |
---|
2085 | character(len=*),parameter :: myname_=myname//'::aCaCSharedAttrIndexList_' |
---|
2086 | |
---|
2087 | integer :: ierr |
---|
2088 | |
---|
2089 | ! Based on the value of the argument attrib, pass the |
---|
2090 | ! appropriate pair of Lists for comparison... |
---|
2091 | |
---|
2092 | select case(trim(attrib)) |
---|
2093 | case('REAL','real') |
---|
2094 | call GetSharedListIndices(aC1%data%rList, aC2%data%rList, NumShared, & |
---|
2095 | Indices1, Indices2) |
---|
2096 | case('INTEGER','integer') |
---|
2097 | call GetSharedListIndices(aC1%data%iList, aC2%data%iList, NumShared, & |
---|
2098 | Indices1, Indices2) |
---|
2099 | case default |
---|
2100 | write(stderr,'(4a)') myname_,":: value of argument attrib=",attrib, & |
---|
2101 | " not recognized. Allowed values: REAL, real, INTEGER, integer" |
---|
2102 | ierr = 1 |
---|
2103 | call die(myname_, 'invalid value for attrib', ierr) |
---|
2104 | end select |
---|
2105 | |
---|
2106 | end subroutine aCaCSharedAttrIndexList_ |
---|
2107 | |
---|
2108 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
2109 | ! Math and Computer Science Division, Argonne National Laboratory ! |
---|
2110 | !BOP ------------------------------------------------------------------- |
---|
2111 | ! |
---|
2112 | ! !IROUTINE: aVaCSharedAttrIndexList_ - Cross-index with an AttrVect |
---|
2113 | ! |
---|
2114 | ! !DESCRIPTION: {\tt aVaCSharedAttrIndexList\_()} a user-supplied |
---|
2115 | ! {\tt AttrVect} variable {\tt aV} and an {\tt Accumulator} variable |
---|
2116 | ! {\tt aC}, and for choice of either {\tt REAL} or {\tt INTEGER} |
---|
2117 | ! attributes (as ! specified literally in the input {\tt CHARACTER} |
---|
2118 | ! argument {\tt attrib}) returns the number of shared attributes |
---|
2119 | ! {\tt NumShared}, and arrays of indices {\tt Indices1} and {\tt Indices2} |
---|
2120 | ! to their storage locations in {\tt aV} and {\tt aC}, respectively. |
---|
2121 | ! |
---|
2122 | ! {\bf N.B.:} This routine returns two allocated arrays---{\tt Indices1(:)} |
---|
2123 | ! and {\tt Indices2(:)}---which must be deallocated once the user no longer |
---|
2124 | ! needs them. Failure to do this will create a memory leak. |
---|
2125 | ! |
---|
2126 | ! !INTERFACE: |
---|
2127 | |
---|
2128 | subroutine aVaCSharedAttrIndexList_(aV, aC, attrib, NumShared, & |
---|
2129 | Indices1, Indices2) |
---|
2130 | |
---|
2131 | ! |
---|
2132 | ! !USES: |
---|
2133 | ! |
---|
2134 | use m_stdio |
---|
2135 | use m_die, only : MP_perr_die, die, warn |
---|
2136 | |
---|
2137 | use m_AttrVect, only : AttrVect |
---|
2138 | |
---|
2139 | use m_List, only : GetSharedListIndices |
---|
2140 | |
---|
2141 | |
---|
2142 | implicit none |
---|
2143 | |
---|
2144 | ! !INPUT PARAMETERS: |
---|
2145 | ! |
---|
2146 | type(AttrVect), intent(in) :: aV |
---|
2147 | type(Accumulator), intent(in) :: aC |
---|
2148 | character(len=*), intent(in) :: attrib |
---|
2149 | |
---|
2150 | ! !OUTPUT PARAMETERS: |
---|
2151 | ! |
---|
2152 | integer, intent(out) :: NumShared |
---|
2153 | integer,dimension(:), pointer :: Indices1 |
---|
2154 | integer,dimension(:), pointer :: Indices2 |
---|
2155 | |
---|
2156 | ! !REVISION HISTORY: |
---|
2157 | ! 7Feb01 - J.W. Larson <larson@mcs.anl.gov> - initial version |
---|
2158 | !EOP ___________________________________________________________________ |
---|
2159 | |
---|
2160 | character(len=*),parameter :: myname_=myname//'::aVaCSharedAttrIndexList_' |
---|
2161 | |
---|
2162 | integer :: ierr |
---|
2163 | |
---|
2164 | ! Based on the value of the argument attrib, pass the |
---|
2165 | ! appropriate pair of Lists for comparison... |
---|
2166 | |
---|
2167 | select case(trim(attrib)) |
---|
2168 | case('REAL','real') |
---|
2169 | call GetSharedListIndices(aV%rList, aC%data%rList, NumShared, & |
---|
2170 | Indices1, Indices2) |
---|
2171 | case('INTEGER','integer') |
---|
2172 | call GetSharedListIndices(aV%iList, aC%data%iList, NumShared, & |
---|
2173 | Indices1, Indices2) |
---|
2174 | case default |
---|
2175 | write(stderr,'(4a)') myname_,":: value of argument attrib=",attrib, & |
---|
2176 | " not recognized. Allowed values: REAL, real, INTEGER, integer" |
---|
2177 | ierr = 1 |
---|
2178 | call die(myname_, 'invalid value for attrib', ierr) |
---|
2179 | end select |
---|
2180 | |
---|
2181 | end subroutine aVaCSharedAttrIndexList_ |
---|
2182 | |
---|
2183 | !BOP ------------------------------------------------------------------- |
---|
2184 | ! |
---|
2185 | ! !IROUTINE: accumulate_--Acumulate from an AttrVect to an Accumulator. |
---|
2186 | ! |
---|
2187 | ! !DESCRIPTION: |
---|
2188 | ! This routine performs time {\em accumlation} of data present in an |
---|
2189 | ! MCT field data {\tt AttrVect} variable {\tt aV} and combines it with |
---|
2190 | ! the running tallies stored in the MCT {\tt Accumulator} variable {\tt aC}. |
---|
2191 | ! This routine automatically identifies which |
---|
2192 | ! fields are held in common by {\tt aV} and {\tt aC} and uses the |
---|
2193 | ! accumulation action information stored in {\tt aC} to decide how |
---|
2194 | ! each field in {\tt aV} is to be combined into its corresponding |
---|
2195 | ! running tally in {\tt aC}. The accumulation operations currently |
---|
2196 | ! supported are: |
---|
2197 | ! \begin {itemize} |
---|
2198 | ! \item {\tt MCT\_SUM}: Add the current values in the {\tt Av} to the current values in {\tt Ac}. |
---|
2199 | ! \item {\tt MCT\_AVG}: Same as {\tt MCT\_SUM} except when {\tt steps\_done} is equal |
---|
2200 | ! to {\tt num\_steps} then perform one more sum and replaced with average. |
---|
2201 | ! \end {itemize} |
---|
2202 | ! |
---|
2203 | ! This routine also automatically increments the counter in {\tt aC} |
---|
2204 | ! signifying the number of steps completed in the accumulation cycle. |
---|
2205 | ! |
---|
2206 | ! NOTE: The user must reset (zero) the {\tt Accumulator} after the average |
---|
2207 | ! has been formed or the next call to {\tt accumulate} will add to the average. |
---|
2208 | ! |
---|
2209 | ! !INTERFACE: |
---|
2210 | |
---|
2211 | subroutine accumulate_(aV, aC) |
---|
2212 | |
---|
2213 | ! |
---|
2214 | ! !USES: |
---|
2215 | ! |
---|
2216 | use m_stdio, only : stdout,stderr |
---|
2217 | use m_die, only : die |
---|
2218 | |
---|
2219 | use m_AttrVect, only : AttrVect |
---|
2220 | use m_AttrVect, only : AttrVect_lsize => lsize |
---|
2221 | use m_AttrVect, only : AttrVect_nIAttr => nIAttr |
---|
2222 | use m_AttrVect, only : AttrVect_nRAttr => nRAttr |
---|
2223 | use m_AttrVect, only : AttrVect_indexRA => indexRA |
---|
2224 | use m_AttrVect, only : AttrVect_indexIA => indexIA |
---|
2225 | |
---|
2226 | implicit none |
---|
2227 | |
---|
2228 | ! !INPUT PARAMETERS: |
---|
2229 | ! |
---|
2230 | type(AttrVect), intent(in) :: aV ! Input AttrVect |
---|
2231 | |
---|
2232 | ! !INPUT/OUTPUT PARAMETERS: |
---|
2233 | ! |
---|
2234 | type(Accumulator), intent(inout) :: aC ! Output Accumulator |
---|
2235 | |
---|
2236 | ! !REVISION HISTORY: |
---|
2237 | ! 18Sep00 - J.W. Larson <larson@mcs.anl.gov> -- initial version. |
---|
2238 | ! 7Feb01 - J.W. Larson <larson@mcs.anl.gov> -- General version. |
---|
2239 | ! 10Jun01 - E.T. Ong -- fixed divide-by-zero problem in integer |
---|
2240 | ! attribute accumulation. |
---|
2241 | ! 27Jul01 - E.T. Ong <eong@mcs.anl.gov> -- removed action argument. |
---|
2242 | ! Make compatible with new Accumulator type. |
---|
2243 | !EOP ___________________________________________________________________ |
---|
2244 | |
---|
2245 | character(len=*),parameter :: myname_=myname//'::accumulate_' |
---|
2246 | |
---|
2247 | ! Overlapping attribute index number |
---|
2248 | integer :: num_indices |
---|
2249 | |
---|
2250 | ! Overlapping attribute index storage arrays: |
---|
2251 | integer, dimension(:), pointer :: aCindices, aVindices |
---|
2252 | integer :: aCindex, aVindex |
---|
2253 | |
---|
2254 | ! Error flag and loop indices |
---|
2255 | integer :: ierr, l, n |
---|
2256 | |
---|
2257 | ! Averaging time-weighting factor: |
---|
2258 | real(FP) :: step_weight |
---|
2259 | integer :: num_steps |
---|
2260 | |
---|
2261 | ! Character variable used as a data type flag: |
---|
2262 | character*7 :: data_flag |
---|
2263 | |
---|
2264 | ! Sanity check of arguments: |
---|
2265 | |
---|
2266 | if(lsize_(aC) /= AttrVect_lsize(aV)) then |
---|
2267 | write(stderr,'(2a,i8,a,i8)') myname_, & |
---|
2268 | ':: Mismatched Accumulator/AttrVect lengths. AttrVect_lsize(aV) = ',& |
---|
2269 | AttrVect_lsize(aV), 'lsize_(aC) = ',lsize_(aC) |
---|
2270 | call die(myname_) |
---|
2271 | endif |
---|
2272 | |
---|
2273 | if(aC%num_steps == 0) then |
---|
2274 | write(stderr,'(2a)') myname,':: FATAL--Zero steps in accumulation cycle.' |
---|
2275 | call die(myname_) |
---|
2276 | endif |
---|
2277 | |
---|
2278 | ! Set num_steps from aC: |
---|
2279 | |
---|
2280 | num_steps = aC%num_steps |
---|
2281 | |
---|
2282 | ! Accumulation of REAL attribute data: |
---|
2283 | |
---|
2284 | if( associated(aC%rAction) ) then ! if summing or avergaging reals... |
---|
2285 | |
---|
2286 | ! Accumulate only if fields are present |
---|
2287 | |
---|
2288 | data_flag = 'REAL' |
---|
2289 | call aVaCSharedAttrIndexList_(aV, aC, data_flag, num_indices, & |
---|
2290 | aVindices, aCindices) |
---|
2291 | |
---|
2292 | if(num_indices > 0) then |
---|
2293 | do n=1,num_indices |
---|
2294 | aVindex = aVindices(n) |
---|
2295 | aCindex = aCindices(n) |
---|
2296 | |
---|
2297 | ! Accumulate if the action is MCT_SUM or MCT_AVG |
---|
2298 | if( (aC%rAction(aCindex) == MCT_SUM).or. & |
---|
2299 | (aC%rAction(aCindex) == MCT_AVG) ) then |
---|
2300 | do l=1,AttrVect_lsize(aV) |
---|
2301 | aC%data%rAttr(aCindex,l) = aC%data%rAttr(aCindex,l) + & |
---|
2302 | aV%rAttr(aVindex,l) |
---|
2303 | end do |
---|
2304 | endif |
---|
2305 | end do |
---|
2306 | |
---|
2307 | deallocate(aVindices, aCindices, stat=ierr) |
---|
2308 | if(ierr /= 0) then |
---|
2309 | write(stderr,'(2a,i8)') myname_, & |
---|
2310 | ':: Error in first deallocate(aVindices...), ierr = ',ierr |
---|
2311 | call die(myname_) |
---|
2312 | endif |
---|
2313 | |
---|
2314 | endif ! if(num_indices > 0) |
---|
2315 | |
---|
2316 | endif ! if( associated(aC%rAction) ) |
---|
2317 | |
---|
2318 | |
---|
2319 | ! Accumulation of INTEGER attribute data: |
---|
2320 | |
---|
2321 | if( associated(aC%iAction) ) then ! if summing or avergaging ints... |
---|
2322 | |
---|
2323 | ! Accumulate only if fields are present |
---|
2324 | |
---|
2325 | |
---|
2326 | data_flag = 'INTEGER' |
---|
2327 | call aVaCSharedAttrIndexList_(aV, aC, data_flag, num_indices, & |
---|
2328 | aVindices, aCindices) |
---|
2329 | |
---|
2330 | if(num_indices > 0) then |
---|
2331 | |
---|
2332 | do n=1,num_indices |
---|
2333 | aVindex = aVindices(n) |
---|
2334 | aCindex = aCindices(n) |
---|
2335 | |
---|
2336 | ! Accumulate if the action is MCT_SUM or MCT_AVG |
---|
2337 | if( (aC%iAction(aCindex) == MCT_SUM) .or. & |
---|
2338 | (aC%iAction(aCindex) == MCT_AVG) ) then |
---|
2339 | do l=1,AttrVect_lsize(aV) |
---|
2340 | aC%data%iAttr(aCindex,l) = aC%data%iAttr(aCindex,l) + & |
---|
2341 | aV%iAttr(aVindex,l) |
---|
2342 | end do |
---|
2343 | endif |
---|
2344 | end do |
---|
2345 | |
---|
2346 | deallocate(aVindices, aCindices, stat=ierr) |
---|
2347 | if(ierr /= 0) then |
---|
2348 | write(stderr,'(2a,i8)') myname_, & |
---|
2349 | ':: Error in second deallocate(aVindices...), ierr = ',ierr |
---|
2350 | call die(myname_) |
---|
2351 | endif |
---|
2352 | |
---|
2353 | endif ! if(num_indices > 0) |
---|
2354 | |
---|
2355 | endif ! if( associated(aC%iAction) ) |
---|
2356 | |
---|
2357 | ! Increment aC%steps_done: |
---|
2358 | |
---|
2359 | aC%steps_done = aC%steps_done + 1 |
---|
2360 | |
---|
2361 | ! If we are at the end of an averaging period, compute the |
---|
2362 | ! average (if desired). |
---|
2363 | |
---|
2364 | if(aC%steps_done == num_steps) then |
---|
2365 | |
---|
2366 | step_weight = 1.0_FP / REAL(num_steps,FP) |
---|
2367 | do n=1,nRAttr_(aC) |
---|
2368 | if( aC%rAction(n) == MCT_AVG ) then |
---|
2369 | do l=1,lsize_(aC) |
---|
2370 | aC%data%rAttr(n,l) = step_weight * aC%data%rAttr(n,l) |
---|
2371 | enddo |
---|
2372 | endif |
---|
2373 | enddo |
---|
2374 | |
---|
2375 | do n=1,nIAttr_(aC) |
---|
2376 | if( aC%iAction(n) == MCT_AVG ) then |
---|
2377 | do l=1,lsize_(aC) |
---|
2378 | aC%data%iAttr(n,l) = aC%data%iAttr(n,l) / num_steps |
---|
2379 | enddo |
---|
2380 | endif |
---|
2381 | enddo |
---|
2382 | |
---|
2383 | endif |
---|
2384 | |
---|
2385 | end subroutine accumulate_ |
---|
2386 | |
---|
2387 | !BOP ------------------------------------------------------------------- |
---|
2388 | ! |
---|
2389 | ! !IROUTINE: average_ -- Force an average to be taken on an Accumulator |
---|
2390 | ! |
---|
2391 | ! !DESCRIPTION: |
---|
2392 | ! This routine will compute the average of the current values in an |
---|
2393 | ! {\tt Accumulator} using the current value of {\tt steps\_done} |
---|
2394 | ! in the {\tt Accumulator} |
---|
2395 | ! |
---|
2396 | ! !INTERFACE: |
---|
2397 | |
---|
2398 | subroutine average_(aC) |
---|
2399 | |
---|
2400 | ! |
---|
2401 | ! !USES: |
---|
2402 | ! |
---|
2403 | use m_stdio, only : stdout,stderr |
---|
2404 | use m_die, only : die |
---|
2405 | |
---|
2406 | use m_AttrVect, only : AttrVect |
---|
2407 | use m_AttrVect, only : AttrVect_lsize => lsize |
---|
2408 | use m_AttrVect, only : AttrVect_nIAttr => nIAttr |
---|
2409 | use m_AttrVect, only : AttrVect_nRAttr => nRAttr |
---|
2410 | use m_AttrVect, only : AttrVect_indexRA => indexRA |
---|
2411 | use m_AttrVect, only : AttrVect_indexIA => indexIA |
---|
2412 | |
---|
2413 | implicit none |
---|
2414 | |
---|
2415 | ! !INPUT/OUTPUT PARAMETERS: |
---|
2416 | ! |
---|
2417 | type(Accumulator), intent(inout) :: aC ! Output Accumulator |
---|
2418 | |
---|
2419 | ! !REVISION HISTORY: |
---|
2420 | ! 11Jan08 - R.Jacob <jacob@mcs.anl.gov> -- initial version based on accumulate_ |
---|
2421 | !EOP ___________________________________________________________________ |
---|
2422 | |
---|
2423 | character(len=*),parameter :: myname_=myname//'::average_' |
---|
2424 | |
---|
2425 | ! Overlapping attribute index number |
---|
2426 | integer :: num_indices |
---|
2427 | |
---|
2428 | ! Overlapping attribute index storage arrays: |
---|
2429 | integer, dimension(:), pointer :: aCindices, aVindices |
---|
2430 | integer :: aCindex, aVindex |
---|
2431 | |
---|
2432 | ! Error flag and loop indices |
---|
2433 | integer :: ierr, l, n |
---|
2434 | |
---|
2435 | ! Averaging time-weighting factor: |
---|
2436 | real(FP) :: step_weight |
---|
2437 | integer :: steps_done |
---|
2438 | |
---|
2439 | |
---|
2440 | if(aC%num_steps == 0) then |
---|
2441 | write(stderr,'(2a)') myname_,':: FATAL--Zero steps in accumulation cycle.' |
---|
2442 | call die(myname_) |
---|
2443 | endif |
---|
2444 | |
---|
2445 | if(aC%steps_done == 0) then |
---|
2446 | write(stderr,'(2a)') myname_,':: FATAL--Zero steps completed in accumulation cycle.' |
---|
2447 | call die(myname_) |
---|
2448 | endif |
---|
2449 | |
---|
2450 | ! Set num_steps from aC: |
---|
2451 | |
---|
2452 | steps_done = aC%steps_done |
---|
2453 | |
---|
2454 | |
---|
2455 | step_weight = 1.0_FP / REAL(steps_done,FP) |
---|
2456 | do n=1,nRAttr_(aC) |
---|
2457 | do l=1,lsize_(aC) |
---|
2458 | aC%data%rAttr(n,l) = step_weight * aC%data%rAttr(n,l) |
---|
2459 | enddo |
---|
2460 | enddo |
---|
2461 | |
---|
2462 | do n=1,nIAttr_(aC) |
---|
2463 | do l=1,lsize_(aC) |
---|
2464 | aC%data%iAttr(n,l) = aC%data%iAttr(n,l) / steps_done |
---|
2465 | enddo |
---|
2466 | enddo |
---|
2467 | |
---|
2468 | |
---|
2469 | end subroutine average_ |
---|
2470 | |
---|
2471 | end module m_Accumulator |
---|