!-------------------------------------- LICENCE BEGIN ------------------------------------ !Environment Canada - Atmospheric Science and Technology License/Disclaimer, ! version 3; Last Modified: May 7, 2008. !This is free but copyrighted software; you can use/redistribute/modify it under the terms !of the Environment Canada - Atmospheric Science and Technology License/Disclaimer !version 3 or (at your option) any later version that should be found at: !http://collaboration.cmc.ec.gc.ca/science/rpn.comm/license.html ! !This software is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; !without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. !See the above mentioned License/Disclaimer for more details. !You should have received a copy of the License/Disclaimer along with this software; !if not, you can write to: EC-RPN COMM Group, 2121 TransCanada, suite 500, Dorval (Quebec), !CANADA, H9P 1J3; or send e-mail to service.rpn@ec.gc.ca !-------------------------------------- LICENCE END -------------------------------------- ***s/r set_poic - preparation of projection matrix in the east-west, * C grid model * #include "model_macros_f.h"*
subroutine set_poic ( F_eval_8, F_evec_8, F_xp0_8, F_xp2_8, 2,2 % F_npts, NSTOR ) * #include "impnone.cdk"
* integer F_npts, NSTOR real*8 F_eval_8(NSTOR), F_evec_8(NSTOR,NSTOR) real*8 F_xp0_8(NSTOR,3), F_xp2_8(NSTOR,3) * *author * jean cote - sept 1995 - from setpois * *revision * v2_00 - Desgagne/Lee - initial MPI version (from setpoic v1_03) * v3_00 - Desgagne & Lee - Lam configuration * *object * See above id * *arguments * Name I/O Description *---------------------------------------------------------------- * F_eval_8 O - eigenvalue vector * F_evec_8 O - eigenvector matrix * F_xp0_8 I - alongx projector diagonals * F_xp2_8 I - alongx second derivatives diagonals * F_npts I - number of points to operate on * *implicits #include "dcst.cdk"
#include "glb_ld.cdk"
#include "glb_pil.cdk"
#include "lun.cdk"
* integer i, j, fni0 real*8 a_8(F_npts,F_npts),b_8(F_npts,F_npts) real*8 d_8(3*F_npts-1), r_8(F_npts) real*8 zero,one,two,pdfaz parameter( zero = 0.0 , one = 1.0, two = 2.0 ) * -------------------------------------------------------------------- * * put input arguments in real*8 arrays * do j=1,NSTOR do i=1,NSTOR F_evec_8(i,j)=zero enddo F_eval_8(j)=zero enddo do i=1,F_npts do j=1,F_npts a_8(i,j) = zero b_8(i,j) = zero end do end do do i = 1, F_npts-1 a_8(i,i+1) = F_xp2_8(i+Lam_pil_w,3) a_8(i,i ) = F_xp2_8(i+Lam_pil_w,2) a_8(i+1,i) = a_8(i,i+1) b_8(i,i+1) = F_xp0_8(i+Lam_pil_w,3) b_8(i,i ) = F_xp0_8(i+Lam_pil_w,2) b_8(i+1,i) = b_8(i,i+1) end do a_8(F_npts,F_npts) = F_xp2_8(NSTOR-Lam_pil_e,2) b_8(F_npts,F_npts) = F_xp0_8(NSTOR-Lam_pil_e,2) if (.not. G_lam) then a_8(F_npts,1 ) = F_xp2_8(F_npts,3) a_8(1,F_npts ) = a_8(F_npts,1) b_8(F_npts,1 ) = F_xp0_8(F_npts,3) b_8(1,F_npts ) = b_8(F_npts,1) endif * if (Lun_debug_L) then print *,'F_npts=',F_npts,' G_ni=',G_ni,'NSTOR=',NSTOR print *,'set_poic: G_xg_8(G_ni-5)=',G_xg_8(G_ni-5),'G_xg_8(6)=',G_xg_8(6) print *,'set_poic: G_yg_8(G_nj-5)=',G_yg_8(G_nj-5),'G_yg_8(6)=',G_yg_8(6) endif if (G_lam) then pdfaz = ONE / sqrt( (G_xg_8(G_ni-glb_pil_e)-G_xg_8(glb_pil_w) ) ) if (Lun_debug_L) then print *,'pdfaz=ONE / sqrt( (G_xg_8(G_ni-glb_pil_e)-G_xg_8(glb_pil_w) ) )' endif call prepoic
( r_8, a_8, b_8, d_8, pdfaz, -ONE, % 'N', fni0, F_npts, F_npts ) else pdfaz = ONE / sqrt( TWO*Dcst_pi_8 ) call prepoic
( r_8, a_8, b_8, d_8, pdfaz, -ONE, % 'P', fni0, F_npts, NSTOR ) endif if (Lun_debug_L) then print *,'set_poic: pdfaz=',pdfaz endif * * put real*8 results in output arguments arrays * NOTE: for non-LAM case, the loops are F_npts=NSTOR * do j = 1,F_npts do i = 1,F_npts F_evec_8(i+Lam_pil_w,j+Lam_pil_s) = a_8(i,j) enddo enddo * do i = 1, F_npts F_eval_8(i+Lam_pil_w) = r_8(i) end do * return end