SUBROUTINE obs_cost (ng, model) ! !git $Id$ !svn $Id: obs_cost.F 1151 2023-02-09 03:08:53Z arango $ !================================================== Hernan G. Arango === ! Copyright (c) 2002-2023 The ROMS/TOMS Group Andrew M. Moore ! ! Licensed under a MIT/X style license ! ! See License_ROMS.md ! !======================================================================= ! ! ! This routine computes the data penalty function directly in during ! ! runs of the nonlinear model: ! ! ! ! Jdata = transpose(H X - Xo) * O^(-1) * (H X - Xo) ! ! ! ! H : observation operator (linearized if incremental) ! ! Xo : observations vector ! ! H X : representer model at observation points ! ! O : observations error covariance ! ! ! !======================================================================= ! USE mod_param USE mod_parallel USE mod_fourdvar USE mod_scalars ! implicit none ! ! Imported variable declarations. ! integer, intent(in) :: ng, model ! ! Local variable declarations. ! integer :: NSUB, iobs, ivar real(r8) :: cff, cff1 real(r8), dimension(0:NobsVar(ng)) :: my_ObsCost ! !----------------------------------------------------------------------- ! Compute observation misfit cost function (ObsCost). !----------------------------------------------------------------------- ! ! Compute nonlinear model data penalty function. ! IF (model.eq.iNLM) THEN DO ivar=0,NobsVar(ng) my_ObsCost(ivar)=0.0_r8 END DO DO iobs=NstrObs(ng),NendObs(ng) ivar=ObsType2State(ObsType(iobs)) IF ((ivar.gt.0).and.(ObsScale(iobs).gt.0.0_r8).and. & & (ObsErr(iobs).ne.0.0_r8)) THEN cff=ObsScale(iobs)*(NLmodVal(iobs)-ObsVal(iobs))**2/ & & ObsErr(iobs) my_ObsCost(0)=my_ObsCost(0)+cff my_ObsCost(ivar)=my_ObsCost(ivar)+cff END IF END DO END IF ! !----------------------------------------------------------------------- ! Load global values. Notice that there is not need for a global ! reduction here since all the threads have the same copy of all ! the vectors used. !----------------------------------------------------------------------- ! NSUB=1 ! distributed-memory !$OMP CRITICAL (COST_FUN) tile_count=tile_count+1 IF (tile_count.eq.NSUB) THEN tile_count=0 IF (model.eq.iNLM) THEN DO ivar=0,NobsVar(ng) FOURDVAR(ng)%NLPenalty(ivar)=my_ObsCost(ivar)+ & & FOURDVAR(ng)%NLPenalty(ivar) END DO END IF END IF !$OMP END CRITICAL (COST_FUN) RETURN END SUBROUTINE obs_cost