TABLE OF CONTENTS


ABINIT/prttagm_images [ Functions ]

[ Top ] [ Functions ]

NAME

 prttagm_images

FUNCTION

 Extension to prttagm to include the printing of
 images information, in those cases the same variable
 is printed several times for each dataset

 Cases where images information are relevant includes
 xcart, xred, acell, fcart.

 INPUT
 (see prttagm.F90)

OUTPUT

  (only writing)

PARENTS

      outvar_a_h,outvar_i_n,outvar_o_z

CHILDREN

      appdig,prttagm,write_var_netcdf

SOURCE

 29 #if defined HAVE_CONFIG_H
 30 #include "config.h"
 31 #endif
 32 
 33 #include "abi_common.h"
 34 
 35 
 36 subroutine prttagm_images(dprarr_images,iout,jdtset_,length,&
 37 & marr,narrm,ncid,ndtset_alloc,token,typevarphys,&
 38 & mxnimage,nimage,ndtset,prtimg,strimg,firstchar,forceprint)
 39 
 40  use defs_basis
 41  use m_profiling_abi
 42 
 43  use m_nctk,   only : write_var_netcdf
 44 
 45 !This section has been created automatically by the script Abilint (TD).
 46 !Do not modify the following lines by hand.
 47 #undef ABI_FUNC
 48 #define ABI_FUNC 'prttagm_images'
 49  use interfaces_32_util
 50  use interfaces_57_iovars, except_this_one => prttagm_images
 51 !End of the abilint section
 52 
 53  implicit none
 54 
 55 !Arguments ------------------------------------
 56 !scalars
 57  integer,intent(in) :: iout,length,marr,ndtset_alloc,ncid
 58  integer,intent(in) :: mxnimage,ndtset
 59  integer,intent(in),optional :: forceprint
 60  character(len=*),intent(in) :: token
 61  character(len=3),intent(in) :: typevarphys
 62  character(len=1),intent(in),optional :: firstchar
 63 !arrays
 64  integer,intent(in) :: prtimg(mxnimage,0:ndtset_alloc)
 65  integer,intent(in) :: jdtset_(0:ndtset_alloc)
 66  integer,intent(in) :: nimage(0:ndtset_alloc)
 67  character(len=8),intent(in) :: strimg(mxnimage)
 68  integer,intent(in) :: narrm(0:ndtset_alloc)
 69  real(dp),intent(in) :: dprarr_images(marr,mxnimage,0:ndtset_alloc)
 70 
 71 !Local variables-------------------------------
 72  integer :: iarr,idtset,iimage,jdtset,multi_narr,narr
 73  integer :: intarr_images(marr,mxnimage,0:ndtset_alloc)
 74  integer,allocatable :: intarr(:,:)
 75  real(dp), allocatable :: dprarr(:,:)
 76  logical :: print_out,print_netcdf,test_multiimages
 77  character(len=1) :: first_column
 78  character(len=4) :: appen
 79  character(len=16) :: keywd
 80  character(len=50) :: full_format
 81  character(len=*), parameter :: format_1  ='",a16,t22,'
 82  character(len=*), parameter :: format_1a ='",a16,a,t22,'
 83  character(len=*), parameter :: format_2  ='",t22,'
 84  character(len=*), parameter :: long_dpr  ='3es18.10)'
 85 !character(len=*), parameter :: format01160 ="(1x,a16,1x,(t22,3es18.10)) "
 86 !character(len=*), parameter :: format01160a="(1x,a16,a,1x,(t22,3es18.10)) "
 87 
 88 ! *************************************************************************
 89 
 90  test_multiimages=.false.
 91  do idtset=1,ndtset_alloc
 92    if(nimage(idtset)>1)then
 93      do iarr=1,narrm(idtset)
 94        if(sum(abs( dprarr_images(iarr,2:nimage(idtset),idtset)- &
 95 &       dprarr_images(iarr,1              ,idtset)))>tol12)then
 96          test_multiimages=.true.
 97        end if
 98      end do
 99    end if
100  end do
101 
102  if(nimage(0)==0)test_multiimages=.true.
103 
104 !DEBUG
105 !if(trim(token)=='vel')then
106 !write(ab_out,*)' test_multiimages=',test_multiimages
107 !endif
108 !ENDDEBUG
109 
110  if(.not.test_multiimages)then
111 
112    narr=narrm(1)
113    ABI_ALLOCATE(intarr,(marr,0:ndtset_alloc))
114    ABI_ALLOCATE(dprarr,(marr,0:ndtset_alloc))
115    do idtset=0,ndtset_alloc
116      dprarr(1:narrm(idtset),idtset)=dprarr_images(1:narrm(idtset),1,idtset)
117 
118 !    DEBUG
119 !    if(trim(token)=='vel')then
120 !    write(ab_out,*)' idtset,narrm(idtset),dprarr(1:narrm(idtset),idtset)=',&
121 !    &    idtset,narrm(idtset),dprarr(1:narrm(idtset),idtset)
122 !    endif
123 !    ENDDEBUG
124 
125    end do
126    multi_narr=0
127    if(ndtset_alloc>1)then
128      do idtset=1,ndtset_alloc
129        if(narrm(1)/=narrm(idtset))multi_narr=1
130      end do
131    end if
132 !  if(narrm(0)==0)multi_narr=1
133 !  DEBUG
134 !  if(trim(token)=='fcart')then
135 !  write(std_out,*)' will call prttagm with fcart '
136 !  write(std_out,*)' narrm(0:ndtset_alloc)=',narrm(0:ndtset_alloc)
137 !  write(std_out,*)' multi_narr=',multi_narr
138 !  do idtset=0,ndtset_alloc
139 !  write(std_out,*)' dprarr_images(1:narrm(idtset),1,idtset)=',dprarr_images(1:narrm(idtset),1,idtset)
140 !  enddo
141 !  endif
142 !  ENDDEBUG
143    if (present(firstchar).and.present(forceprint)) then
144      call prttagm(dprarr,intarr,iout,jdtset_,length,marr,narr,&
145 &     narrm,ncid,ndtset_alloc,token,typevarphys,multi_narr,&
146 &     firstchar=firstchar,forceprint=forceprint)
147    else if (present(firstchar)) then
148      call prttagm(dprarr,intarr,iout,jdtset_,length,marr,narr,&
149 &     narrm,ncid,ndtset_alloc,token,typevarphys,multi_narr,&
150 &     firstchar=firstchar)
151    else if (present(forceprint)) then
152      call prttagm(dprarr,intarr,iout,jdtset_,length,marr,narr,&
153 &     narrm,ncid,ndtset_alloc,token,typevarphys,multi_narr,&
154 &     forceprint=forceprint)
155    else
156      call prttagm(dprarr,intarr,iout,jdtset_,length,marr,narr,&
157 &     narrm,ncid,ndtset_alloc,token,typevarphys,multi_narr)
158    end if
159    ABI_DEALLOCATE(intarr)
160    ABI_DEALLOCATE(dprarr)
161 
162  else
163 
164    first_column=' ';if (present(firstchar)) first_column=firstchar
165 
166    do idtset=1,ndtset_alloc
167 
168      if (narrm(idtset)>0)then
169        do iimage=1,nimage(idtset)
170 
171          print_out=.true.
172          if (prtimg(iimage,idtset)==0) print_out=.false.
173          if (nimage(0)>=nimage(idtset)) then
174            if (sum(abs(dprarr_images(1:narrm(idtset),iimage,idtset) &
175 &           -dprarr_images(1:narrm(idtset),iimage,0)))<tol12) print_out=.false.
176          end if
177          print_netcdf=print_out
178 
179          if (present(forceprint)) then
180            if (forceprint==1.or.forceprint==3) print_out=.true.
181            if (forceprint==1.or.forceprint==2) print_netcdf=.true.
182          end if
183 
184          if (print_out.or.print_netcdf.or.(ncid<0))then
185            keywd=token//trim(strimg(iimage))
186 
187            if(ndtset>0)then
188              jdtset=jdtset_(idtset)
189              call appdig(jdtset,'',appen)
190              if (print_out) then
191                full_format='("'//first_column//trim(format_1a)//'("'// &
192 &               first_column//trim(format_2)//trim(long_dpr)//")"
193                write(iout,full_format) &
194 &               trim(keywd),appen,dprarr_images(1:narrm(idtset),iimage,idtset)
195              end if
196 #ifdef HAVE_NETCDF
197              if (print_netcdf) then
198                call write_var_netcdf(intarr_images(1:narrm(idtset),iimage,idtset),&
199 &               dprarr_images(1:narrm(idtset),iimage,idtset),&
200 &               marr,narrm(idtset),ncid,'DPR',trim(keywd)//appen)
201              end if
202 #endif
203            else
204 
205              if (print_out) then
206                full_format='("'//first_column//trim(format_1)//'("'// &
207 &               first_column//trim(format_2)//trim(long_dpr)//")"
208                write(iout,full_format) &
209 &               trim(keywd),dprarr_images(1:narrm(idtset),iimage,idtset)
210              end if
211 #ifdef HAVE_NETCDF
212              if (print_netcdf) then
213                call write_var_netcdf(intarr_images(1:narrm(idtset),iimage,idtset),&
214 &               dprarr_images(1:narrm(idtset),iimage,idtset),&
215 &               marr,narrm(idtset),abs(ncid),'DPR',trim(keywd))
216              end if
217 #endif
218 
219            end if
220          end if
221        end do
222      end if
223    end do
224 
225  end if
226 
227 end subroutine prttagm_images