Back to home page

MITgcm

 
 

    


File indexing completed on 2025-08-05 05:09:15 UTC

view on githubraw file Latest commit 13ce79fe on 2025-08-04 21:05:34 UTC
24462d2fa8 Patr*0001 #include "PROFILES_OPTIONS.h"
57c22ecc45 Jean*0002 #include "AD_CONFIG.h"
6a770e0a24 Patr*0003 
c9bf163375 Ivan*0004 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
13ce79fe94 Ivan*0005 C     BOP
                0006 C     !ROUTINE: PROFILES_INIT_FIXED
6a770e0a24 Patr*0007 
13ce79fe94 Ivan*0008 C     !INTERFACE:
                0009       SUBROUTINE PROFILES_INIT_FIXED( myThid )
                0010 
                0011 C     !DESCRIPTION:
                0012 C     Initialization for netcdf profiles data
6a770e0a24 Patr*0013 
13ce79fe94 Ivan*0014 C     !USES:
                0015       IMPLICIT NONE
                0016 C     == Global variables ===
6a770e0a24 Patr*0017 #include "SIZE.h"
                0018 #include "EEPARAMS.h"
                0019 #include "PARAMS.h"
                0020 #include "GRID.h"
                0021 #include "DYNVARS.h"
d28c90138c Patr*0022 #ifdef ALLOW_CAL
6a770e0a24 Patr*0023 #include "cal.h"
d28c90138c Patr*0024 #endif
24462d2fa8 Patr*0025 #ifdef ALLOW_PROFILES
6328b73337 Gael*0026 # include "PROFILES_SIZE.h"
6e4c90fea3 Patr*0027 # include "profiles.h"
                0028 # include "netcdf.inc"
                0029 #endif
39ce977435 Gael*0030 
13ce79fe94 Ivan*0031 C     !INPUT/OUTPUT PARAMETERS:
                0032 C     myThid: my thread ID number
c9bf163375 Ivan*0033       INTEGER myThid
13ce79fe94 Ivan*0034 C     EOP
39ce977435 Gael*0035 
13ce79fe94 Ivan*0036 #ifdef ALLOW_PROFILES
                0037 C     !FUNCTIONS:
                0038       INTEGER  ILNBLNK
c9bf163375 Ivan*0039       EXTERNAL ILNBLNK
13ce79fe94 Ivan*0040       INTEGER  MDS_RECLEN
c9bf163375 Ivan*0041       EXTERNAL MDS_RECLEN
6a770e0a24 Patr*0042 
13ce79fe94 Ivan*0043 C     !LOCAL VARIABLES:
                0044 C     i,j,bi,bj    :: loop indices
                0045 C     kLev         :: vertical level index
                0046 C     prof_num, num_file, num_var
                0047 C                  :: loop indices for profiles, files, and variables
                0048 C     iG,jG        :: global indices
                0049 C     ic, chunkProf:: chunk-related indices
                0050 C     chunk        :: chunk size
                0051 C     ProfNo_tile  :: counter
c9bf163375 Ivan*0052       CHARACTER*(MAX_LEN_MBUF) msgBuf
13ce79fe94 Ivan*0053       INTEGER iUnit
                0054       INTEGER kLev, prof_num
                0055       INTEGER chunk, chunkProf, ic, recLen
                0056       INTEGER i,j,bi,bj,iG,jG,num_file,num_var,ProfNo_tile
c9bf163375 Ivan*0057       INTEGER stopProfiles
13ce79fe94 Ivan*0058       INTEGER fid, dimId, varId1, varId1a, varId1b
                0059       INTEGER varId2,varId3
6a770e0a24 Patr*0060       _RL tmpyymmdd(1000),tmphhmmss(1000),diffsecs
5042c05de8 Gael*0061       _RL yymmddMin,yymmddMax
                0062       _RL hhmmssMin,hhmmssMax
                0063 
c9bf163375 Ivan*0064       INTEGER tmpdate(4),tmpdiff(4),profIsInRunTime
13ce79fe94 Ivan*0065       _RL tmp_lon, tmp_lon2(1000), tmp_lat2(1000)
                0066       _RL lon_cur, lat_cur
39ce977435 Gael*0067       _RL lon_1, lon_2, lat_1, lat_2
2767dff983 Jean*0068       _RL lon_tmp1, lon_tmp2
39ce977435 Gael*0069       _RL lat_fac, lon_fac
c9bf163375 Ivan*0070       INTEGER prof_i, prof_j
13ce79fe94 Ivan*0071       INTEGER vec_start(2), vec_count(2), profno_div1000
c9bf163375 Ivan*0072       CHARACTER*(MAX_LEN_FNAM) profilesfile, fnamedatanc
                0073       CHARACTER*(MAX_LEN_FNAM) fnameequinc
13ce79fe94 Ivan*0074 # ifdef ALLOW_ADJOINT_RUN
                0075       CHARACTER*(MAX_LEN_FNAM) adfnameequinc
                0076 # endif
                0077 # ifdef ALLOW_TANGENTLINEAR_RUN
                0078       CHARACTER*(MAX_LEN_FNAM) tlfnameequinc
                0079 # endif
                0080       INTEGER IL, JL, KL, err, nerr
                0081       LOGICAL exst
                0082 
                0083       INTEGER varId_intp1, varId_intp2, varId_intp11, varId_intp22
                0084       INTEGER varId_intp3, varId_intp4, varId_intp5, iq, iINTERP
ba63501b4c Gael*0085       _RL tmp_i(1000,NUM_INTERP_POINTS)
                0086       _RL tmp_j(1000,NUM_INTERP_POINTS)
                0087       _RL tmp_weights(1000,NUM_INTERP_POINTS),tmp_sum_weights
                0088       _RL tmp_xC11(1000),tmp_yC11(1000)
                0089       _RL tmp_xCNINJ(1000),tmp_yCNINJ(1000)
c9bf163375 Ivan*0090       INTEGER stopGenericGrid
ba63501b4c Gael*0091       Real*8 xy_buffer_r8(0:sNx+1,0:sNy+1)
c9bf163375 Ivan*0092       INTEGER vec_start2(2), vec_count2(2)
                0093       INTEGER hh, ProfNo_hh
13ce79fe94 Ivan*0094 # ifdef ALLOW_PROFILES_SAMPLESPLIT_COST
                0095       INTEGER m, kC, kCMax
                0096       INTEGER varId4
6b2230d510 Ou W*0097       _RL tmp_avgbin(1000)
13ce79fe94 Ivan*0098 # endif /* ALLOW_PROFILES_SAMPLESPLIT_COST */
                0099 
                0100       iUnit = standardMessageUnit
                0101       WRITE(msgBuf,'(A)') ' '
                0102       CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                0103       WRITE(msgBuf,'(A)')
                0104      &     '// ======================================================='
                0105       CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                0106       WRITE(msgBuf,'(A)')
                0107      &     '// insitu profiles model sampling >>> START <<<'
                0108       CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                0109       WRITE(msgBuf,'(A)')
                0110      &     '// ======================================================='
                0111       CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                0112       WRITE(msgBuf,'(A)') ' '
                0113       CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                0114 
                0115       stopProfiles = 0
                0116       stopGenericGrid = 0
                0117 
                0118       IF ( (.NOT.profilesDoGenGrid) .AND.
f0e4bffe35 Gael*0119      &     (.NOT.usingSphericalPolarGrid .OR. rotateGrid) ) THEN
13ce79fe94 Ivan*0120        WRITE(msgBuf,'(2A)') 'PROFILES_INIT_FIXED: ',
                0121      &      'profilesDoGenGrid=.true. is required'
                0122        CALL PRINT_ERROR( msgBuf , myThid )
                0123        WRITE(msgBuf,'(2A)') 'PROFILES_INIT_FIXED: ',
                0124      &      'unless usingSphericalGrid=.TRUE. and rotateGrid=.FALSE.'
                0125        CALL PRINT_ERROR( msgBuf , myThid )
                0126        CALL ALL_PROC_DIE( myThid )
                0127        STOP 'ABNORMAL END: S/R PROFILES_INIT_FIXED'
f0e4bffe35 Gael*0128       ENDIF
                0129 
13ce79fe94 Ivan*0130       WRITE(msgBuf,'(A)') ' '
                0131       CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                0132       WRITE(msgBuf,'(A)') 'general packages parameters :'
                0133 
                0134       JL = ILNBLNK( profilesDir )
c9bf163375 Ivan*0135       IF (JL.NE.0) THEN
13ce79fe94 Ivan*0136        WRITE(msgBuf,'(2A)') '  profilesDir ', profilesDir(1:JL)
c9bf163375 Ivan*0137       ELSE
13ce79fe94 Ivan*0138        WRITE(msgBuf,'(2A)') '  profilesDir ','./'
c9bf163375 Ivan*0139       ENDIF
                0140 
13ce79fe94 Ivan*0141       CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                0142       WRITE(msgBuf,'(A,L5)') '  profilesDoGenGrid  ', profilesDoGenGrid
                0143       CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                0144       WRITE(msgBuf,'(A,L5)') '  profilesDoNcOutput ', profilesDoNcOutput
                0145       CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                0146       WRITE(msgBuf,'(A)') ' '
                0147       CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
b00d6c1700 Gael*0148 
13ce79fe94 Ivan*0149       _BEGIN_MASTER( myThid )
71a5587721 Gael*0150 
13ce79fe94 Ivan*0151       DO bj = 1, nSy
                0152        DO bi = 1, nSx
                0153         profiles_curfile_buff(bi,bj) = 0
5042c05de8 Gael*0154         yymmddMin=modelstartdate(1)
                0155         yymmddMax=modelenddate(1)
                0156         hhmmssMin=modelstartdate(2)
                0157         hhmmssMax=modelenddate(2)
6e4c90fea3 Patr*0158 
13ce79fe94 Ivan*0159         DO kLev = 1, NLEVELMAX
                0160          DO ic = 1, 1000
                0161           DO num_var = 1, NVARMAX
                0162            profiles_data_buff  (kLev,ic,num_var,bi,bj) = 0. _d 0
                0163            profiles_weight_buff(kLev,ic,num_var,bi,bj) = 0. _d 0
c9bf163375 Ivan*0164           ENDDO
                0165          ENDDO
                0166         ENDDO
d5aa75d39a Jean*0167 
13ce79fe94 Ivan*0168         DO num_file = 1, NFILESPROFMAX
                0169          ProfNo_hh=0
                0170          profilesfile=' '
                0171          IL = ILNBLNK( profilesfiles(num_file) )
                0172          IF (IL.NE.0) THEN
                0173           WRITE(profilesfile,'(A)') profilesfiles(num_file)(1:IL)
                0174           WRITE(msgBuf,'(A)') ' '
                0175           CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                0176           WRITE(msgBuf,'(A,I3,2A)') 'profiles file #', num_file,
                0177      &         ' is ', profilesfile(1:IL)
                0178           CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                0179          ENDIF
6a770e0a24 Patr*0180 
13ce79fe94 Ivan*0181          IL = ILNBLNK( profilesfile )
                0182          IF (IL.NE.0) THEN
6a770e0a24 Patr*0183 C===========================================================
13ce79fe94 Ivan*0184 C     Open data files and read information
6a770e0a24 Patr*0185 C===========================================================
13ce79fe94 Ivan*0186           WRITE(fnamedatanc,'(2A)') profilesfile(1:IL),'.nc'
                0187           err = NF_OPEN( fnamedatanc, NF_NOWRITE,
                0188      &         fiddata(num_file,bi,bj) )
                0189           CALL PROFILES_NF_ERROR(
                0190      &         'INIT_FIXED: NF_OPEN fiddata',err,bi,bj,myThid )
                0191 
                0192 C     1) Read the number of profiles and available dimensions:
                0193           fid = fiddata(num_file,bi,bj)
                0194 
                0195           err = NF_INQ_DIMID( fid, 'iPROF', dimId )
                0196           CALL PROFILES_NF_ERROR(
                0197      &         'INIT_FIXED: NF_INQ_DIMID iPROF',err,bi,bj,myThid )
                0198           err = NF_INQ_DIMLEN( fid,
                0199      &         dimId, ProfNo(num_file,bi,bj) )
                0200           CALL PROFILES_NF_ERROR(
                0201      &         'INIT_FIXED: NF_INQ_DIMLEN ProfNo',
                0202      &         err,bi,bj,myThid )
                0203 
                0204           err = NF_INQ_DIMID( fid, 'iDEPTH', dimId )
                0205           CALL PROFILES_NF_ERROR(
                0206      &         'INIT_FIXED: NF_INQ_DIMID iDEPTH',
                0207      &         err,bi,bj,myThid )
                0208           IF (err.NE.NF_NOERR) THEN
                0209            err = NF_INQ_DIMID( fid, 'Z', dimId )
                0210            CALL PROFILES_NF_ERROR(
                0211      &          'INIT_FIXED: NF_INQ_DIMID Z',err,bi,bj,myThid )
6a770e0a24 Patr*0212 
13ce79fe94 Ivan*0213           ENDIF
6a770e0a24 Patr*0214 
13ce79fe94 Ivan*0215           err = NF_INQ_DIMLEN( fid, dimId,
                0216      &         ProfDepthNo(num_file,bi,bj) )
                0217           CALL PROFILES_NF_ERROR(
                0218      &         'INIT_FIXED: NF_INQ_DIMLEN ProfDepthNo',
                0219      &         err,bi,bj,myThid )
f0e4bffe35 Gael*0220 
13ce79fe94 Ivan*0221           err = NF_INQ_DIMID( fid, 'iINTERP', dimId )
                0222           IF (err.EQ.NF_NOERR) THEN
                0223            err = NF_INQ_DIMLEN( fid, dimId, iINTERP )
                0224            CALL PROFILES_NF_ERROR(
                0225      &          'INIT_FIXED: NF_INQ_DIMLEN iINTERP',
                0226      &          err,bi,bj,myThid )
b632e3ba1b Gael*0227 
13ce79fe94 Ivan*0228           ELSE
6a770e0a24 Patr*0229 
13ce79fe94 Ivan*0230            IF (debugLevel .GE. debLevA) THEN
                0231             CALL PROFILES_NF_ERROR(
                0232      &           'INIT_FIXED: NF_INQ_DIMID iINTERP',
                0233      &           err,bi,bj,myThid )
                0234             WRITE(msgBuf,'(3A,I3)')
                0235      &           'S/R PROFILES_INIT_FIXED: ',
                0236      &           'no iINTERP dim in data file using iINTERP ',
                0237      &           '= NUM_INTERP_POINTS =', NUM_INTERP_POINTS
                0238             CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
ba63501b4c Gael*0239 
13ce79fe94 Ivan*0240            ENDIF
ba63501b4c Gael*0241 
13ce79fe94 Ivan*0242            iINTERP = NUM_INTERP_POINTS
6a770e0a24 Patr*0243 
13ce79fe94 Ivan*0244           ENDIF
6a770e0a24 Patr*0245 
13ce79fe94 Ivan*0246           WRITE(msgBuf,'(2(A,I4))')
                0247      &         '  current tile is bi,bj                      =',
                0248      &         bi,',',bj
                0249           CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                0250           WRITE(msgBuf,'(A,I9)')
                0251      &         '  # of depth levels in file                  =',
                0252      &         ProfDepthNo(num_file,bi,bj)
                0253           CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                0254           WRITE(msgBuf,'(A,I9)')
                0255      &         '  # of profiles in file                      =',
                0256      &         ProfNo(num_file,bi,bj)
                0257           CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                0258 
                0259 C     2) Read dates and positions:
                0260           err = NF_INQ_VARID( fid,'prof_depth', varId1a )
                0261           CALL PROFILES_NF_ERROR(
                0262      &         'INIT_FIXED: NF_INQ_VARID prof_depth',
                0263      &         err,bi,bj,myThid )
                0264           IF (err.NE.NF_NOERR) THEN
                0265 C     LEGACY: Try old variable name: depth
                0266            err = NF_INQ_VARID( fid,'depth', varId1a )
                0267            CALL PROFILES_NF_ERROR(
                0268      &          'INIT_FIXED: NF_INQ_VARID depth',
                0269      &          err,bi,bj,myThid )
6a770e0a24 Patr*0270 
13ce79fe94 Ivan*0271           ENDIF
ba63501b4c Gael*0272 
13ce79fe94 Ivan*0273           IF (err.NE.NF_NOERR) THEN
                0274 C     If neither is found, then stop
                0275            IL = ILNBLNK( profilesfile )
                0276            WRITE(msgBuf,'(4A)')
                0277      &          'PROFILES_INIT_FIXED: file ', profilesfile(1:IL),
                0278      &          '.nc is not in the pkg/profiles format',
                0279      &          ' (no prof_depth etc.)'
                0280            CALL PRINT_ERROR( msgBuf, myThid )
6a770e0a24 Patr*0281 
13ce79fe94 Ivan*0282            stopProfiles = 1
                0283           ENDIF
39ce977435 Gael*0284 
13ce79fe94 Ivan*0285           DO kLev = 1, ProfDepthNo(num_file,bi,bj)
                0286            err = NF_GET_VAR1_DOUBLE( fid, varId1a, kLev,
                0287      &          prof_depth(num_file,kLev,bi,bj) )
                0288            CALL PROFILES_NF_ERROR(
                0289      &          'INIT_FIXED: NF_VAR1_DOUBLE prof_depth',
                0290      &          err,bi,bj,myThid )
                0291           ENDDO
39ce977435 Gael*0292 
13ce79fe94 Ivan*0293 C     Get time, lon, lat varIds
                0294           err = NF_INQ_VARID( fid,'prof_YYYYMMDD', varId1a )
                0295           nerr = err
                0296           CALL PROFILES_NF_ERROR(
                0297      &         'INIT_FIXED: NF_INQ_VARID prof_YYYYMMDD',
                0298      &         err,bi,bj,myThid )
                0299           err = NF_INQ_VARID( fid,'prof_HHMMSS', varId1b )
                0300           nerr = nerr + err
                0301           CALL PROFILES_NF_ERROR(
                0302      &         'INIT_FIXED: NF_INQ_VARID prof_HHMMSS',
                0303      &         err,bi,bj,myThid )
                0304           err = NF_INQ_VARID( fid,'prof_lon', varId2 )
                0305           nerr = nerr + err
                0306           CALL PROFILES_NF_ERROR(
                0307      &         'INIT_FIXED: NF_INQ_VARID prof_lon',
                0308      &         err,bi,bj,myThid )
                0309           err = NF_INQ_VARID( fid,'prof_lat', varId3 )
                0310           nerr = nerr + err
                0311           CALL PROFILES_NF_ERROR(
                0312      &         'INIT_FIXED: NF_INQ_VARID prof_lat',
                0313      &         err,bi,bj,myThid )
                0314 # ifdef ALLOW_PROFILES_SAMPLESPLIT_COST
                0315           err = NF_INQ_VARID( fid,'prof_bin_id_a', varId4 )
                0316           nerr = nerr + err
                0317           CALL PROFILES_NF_ERROR(
                0318      &         'INIT_FIXED: NF_INQ_VARID prof_bin_id_a',
                0319      &         err,bi,bj,myThid )
                0320 # endif /* ALLOW_PROFILES_SAMPLESPLIT_COST */
                0321 
                0322           IF (nerr.NE.NF_NOERR) THEN
                0323            IL = ILNBLNK( profilesfile )
                0324            WRITE(msgBuf,'(4A)')
                0325      &          'PROFILES_INIT_FIXED: file ', profilesfile(1:IL),
                0326      &          '.nc not in pkg/profiles format',
                0327      &          ' (no prof_YYYYMMDD etc.)'
                0328            CALL PRINT_ERROR( msgBuf, myThid )
                0329 
                0330            stopProfiles = 1
                0331           ENDIF
39ce977435 Gael*0332 
13ce79fe94 Ivan*0333           IF (profilesDoGenGrid) THEN
                0334 C     3) Read interpolation information (grid points, coeffs, etc.)
                0335            err = NF_INQ_VARID( fid,'prof_interp_XC11', varId_intp1 )
                0336            nerr = err
                0337            CALL PROFILES_NF_ERROR(
                0338      &          'INIT_FIXED: NF_INQ_VARID prof_interp_XC11',
                0339      &          err,bi,bj,myThid )
                0340            err = NF_INQ_VARID( fid,'prof_interp_YC11', varId_intp2 )
                0341            nerr = nerr + err
                0342            CALL PROFILES_NF_ERROR(
                0343      &          'INIT_FIXED: NF_INQ_VARID prof_interp_YC11',
                0344      &          err,bi,bj,myThid )
                0345            err = NF_INQ_VARID( fid,'prof_interp_XCNINJ', varId_intp11 )
                0346            nerr = nerr + err
                0347            CALL PROFILES_NF_ERROR(
                0348      &          'INIT_FIXED: NF_INQ_VARID prof_interp_XCNINJ',
                0349      &          err,bi,bj,myThid )
                0350            err = NF_INQ_VARID( fid,'prof_interp_YCNINJ', varId_intp22 )
                0351            nerr = nerr + err
                0352            CALL PROFILES_NF_ERROR(
                0353      &          'INIT_FIXED: NF_INQ_VARID prof_interp_YCNINJ',
                0354      &          err,bi,bj,myThid )
                0355            err = NF_INQ_VARID( fid,'prof_interp_weights', varId_intp3 )
                0356            nerr = nerr + err
                0357            CALL PROFILES_NF_ERROR(
                0358      &          'INIT_FIXED: NF_INQ_VARID prof_interp_weights',
                0359      &          err,bi,bj,myThid )
                0360            err = NF_INQ_VARID( fid,'prof_interp_i', varId_intp4 )
                0361            nerr = nerr + err
                0362            CALL PROFILES_NF_ERROR(
                0363      &          'INIT_FIXED: NF_INQ_VARID prof_interp_i',
                0364      &          err,bi,bj,myThid )
                0365            err = NF_INQ_VARID( fid,'prof_interp_j', varId_intp5 )
                0366            nerr = nerr + err
                0367            CALL PROFILES_NF_ERROR(
                0368      &          'INIT_FIXED: NF_INQ_VARID prof_interp_j',
                0369      &          err,bi,bj,myThid )
                0370 
                0371            IF (nerr.NE.NF_NOERR) THEN
                0372             IL = ILNBLNK( profilesfile )
                0373             WRITE(msgBuf,'(4A)')
                0374      &           'PROFILES_INIT_FIXED: file ', profilesfile(1:IL),
                0375      &           '.nc missing interpolation information',
                0376      &           ' (profilesDoGenGrid)'
                0377             CALL PRINT_ERROR( msgBuf, myThid )
                0378 
                0379             stopGenericGrid = 2
                0380            ENDIF
                0381           ENDIF !IF (profilesDoGenGrid)
                0382 
                0383 C     4) Default values
                0384           DO prof_num = 1, NOBSGLOB
                0385            prof_time     (num_file,prof_num,bi,bj) = -999. _d 0
                0386            prof_lon      (num_file,prof_num,bi,bj) = -999. _d 0
                0387            prof_lat      (num_file,prof_num,bi,bj) = -999. _d 0
                0388            prof_ind_glob (num_file,prof_num,bi,bj) = 0
                0389 # ifdef ALLOW_PROFILES_SAMPLESPLIT_COST
                0390            prof_ind_avgbin(num_file,prof_num,bi,bj) = -999
                0391 # endif /* ALLOW_PROFILES_SAMPLESPLIT_COST */
                0392 
                0393            DO iq = 1, NUM_INTERP_POINTS
                0394             prof_interp_i(num_file,prof_num,iq,bi,bj) = 1
                0395             prof_interp_j(num_file,prof_num,iq,bi,bj) = 1
                0396             prof_interp_weights(num_file,prof_num,iq,bi,bj) = 0. _d 0
                0397            ENDDO
39ce977435 Gael*0398 
13ce79fe94 Ivan*0399            prof_interp_xC11  (num_file,prof_num,bi,bj) = -999. _d 0
                0400            prof_interp_yC11  (num_file,prof_num,bi,bj) = -999. _d 0
                0401            prof_interp_xCNINJ(num_file,prof_num,bi,bj) = -999. _d 0
                0402            prof_interp_yCNINJ(num_file,prof_num,bi,bj) = -999. _d 0
                0403           ENDDO !DO prof_num
                0404 
                0405 C     5) Main loop: look for profiles in this bi,bj tile
                0406           ProfNo_tile = 0
                0407           profNo_div1000 = MAX(0,INT(ProfNo(num_file,bi,bj)/1000))
                0408 
                0409           DO chunkProf = 1, profno_div1000+1
                0410 C     5.1) Read a chunk
                0411            chunk = 1000*(chunkProf-1)
                0412 
                0413            IF (MIN(ProfNo(num_file,bi,bj), 1000*chunkProf).GE.
                0414      &          1+chunk) THEN
                0415             vec_start(1) = 1
                0416             vec_start(2) = 1+chunk
                0417             vec_count(1) = 1
                0418             vec_count(2) = MIN(1000, ProfNo(num_file,bi,bj)-chunk)
                0419 
                0420             IF ( (vec_count(2).LE.0) .OR.
                0421      &           (vec_count(2).GT.1000) .OR.
                0422      &           (vec_start(2).LE.0) .OR.
                0423      &           (vec_count(2)+vec_start(2)-1.GT.
                0424      &           ProfNo(num_file,bi,bj)) ) THEN
                0425              IL  = ILNBLNK( profilesfile )
                0426              WRITE(msgBuf,'(3A)')
                0427      &            'PROFILES_INIT_FIXED: file ',profilesfile(1:IL),
                0428      &            '.nc was not read properly (case 1).'
                0429              CALL PRINT_ERROR( msgBuf, myThid )
                0430 
                0431              stopProfiles = 1
                0432             ENDIF
                0433 
                0434             err = NF_GET_VARA_DOUBLE( fid, varId1a,
                0435      &           vec_start(2),vec_count(2), tmpyymmdd )
                0436             nerr = err
                0437             CALL PROFILES_NF_ERROR(
                0438      &           'INIT_FIXED: NF_GET_VARA_DOUBLE tmpyymmdd',
                0439      &           err,bi,bj,myThid )
                0440             err = NF_GET_VARA_DOUBLE( fid, varId1b,
                0441      &           vec_start(2), vec_count(2), tmphhmmss )
                0442             nerr = nerr + err
                0443             CALL PROFILES_NF_ERROR(
                0444      &           'INIT_FIXED: NF_GET_VARA_DOUBLE tmphhmmss',
                0445      &           err,bi,bj,myThid )
                0446             err = NF_GET_VARA_DOUBLE( fid, varId2,
                0447      &           vec_start(2),vec_count(2), tmp_lon2 )
                0448             nerr = nerr + err
                0449             CALL PROFILES_NF_ERROR(
                0450      &           'INIT_FIXED: NF_GET_VARA_DOUBLE tmp_lon2',
                0451      &           err,bi,bj,myThid )
                0452             err = NF_GET_VARA_DOUBLE( fid, varId3,
                0453      &           vec_start(2),vec_count(2), tmp_lat2 )
                0454             nerr = nerr + err
                0455             CALL PROFILES_NF_ERROR(
                0456      &           'INIT_FIXED: NF_GET_VARA_DOUBLE tmp_lat2',
                0457      &           err,bi,bj,myThid )
                0458 # ifdef ALLOW_PROFILES_SAMPLESPLIT_COST
                0459             err = NF_GET_VARA_DOUBLE( fid, varId4,
                0460      &           vec_start(2),vec_count(2), tmp_avgbin )
                0461             nerr = nerr + err
                0462             CALL PROFILES_NF_ERROR(
                0463      &           'INIT_FIXED: NF_GET_VARA_DOUBLE tmp_avgbin',
                0464      &           err,bi,bj,myThid )
                0465 # endif /* ALLOW_PROFILES_SAMPLESPLIT_COST */
                0466 
                0467             IF (nerr.NE.NF_NOERR) THEN
                0468              WRITE(msgBuf,'(3A)')
                0469      &            'PROFILES_INIT_FIXED: file ', profilesfile(1:IL),
                0470      &            '.nc get_vara_double (case 2).'
                0471              CALL PRINT_ERROR( msgBuf, myThid )
                0472 
                0473              stopProfiles = 1
                0474             ENDIF
                0475 
                0476 C     If profilesDoGenGrid, then read interpolation coeffs and indices
                0477             IF (profilesDoGenGrid) THEN
                0478              err = NF_GET_VARA_DOUBLE( fid, varId_intp1,
                0479      &            vec_start(2),vec_count(2), tmp_xC11 )
                0480              CALL PROFILES_NF_ERROR(
                0481      &            'INIT_FIXED: NF_GET_VARA_DOUBLE tmp_xC11',
                0482      &            err,bi,bj,myThid )
                0483              err = NF_GET_VARA_DOUBLE( fid, varId_intp2,
                0484      &            vec_start(2),vec_count(2), tmp_yC11 )
                0485              CALL PROFILES_NF_ERROR(
                0486      &            'INIT_FIXED: NF_GET_VARA_DOUBLE tmp_yC11',
                0487      &            err,bi,bj,myThid )
                0488              err = NF_GET_VARA_DOUBLE( fid, varId_intp11,
                0489      &            vec_start(2),vec_count(2), tmp_xCNINJ )
                0490              CALL PROFILES_NF_ERROR(
                0491      &            'INIT_FIXED: NF_GET_VARA_DOUBLE tmp_xCNINJ',
                0492      &            err,bi,bj,myThid )
                0493              err = NF_GET_VARA_DOUBLE( fid, varId_intp22,
                0494      &            vec_start(2),vec_count(2), tmp_yCNINJ )
                0495              CALL PROFILES_NF_ERROR(
                0496      &            'INIT_FIXED: NF_GET_VARA_DOUBLE tmp_yCNINJ',
                0497      &            err,bi,bj,myThid )
                0498 
                0499              DO iq = 1, iINTERP
                0500               vec_start2(1) = iq
                0501               vec_start2(2) = 1+chunk
                0502               vec_count2(1) = 1
                0503               vec_count2(2) = MIN(1000,ProfNo(num_file,bi,bj)-chunk)
                0504 
                0505               err = NF_GET_VARA_DOUBLE( fid, varId_intp3,
                0506      &             vec_start2,vec_count2, tmp_weights(1,iq) )
                0507               CALL PROFILES_NF_ERROR(
                0508      &             'INIT_FIXED: NF_GET_VARA_DOUBLE tmp_weights',
                0509      &             err,bi,bj,myThid )
                0510               err = NF_GET_VARA_DOUBLE( fid, varId_intp4,
                0511      &             vec_start2,vec_count2, tmp_i(1,iq) )
                0512               CALL PROFILES_NF_ERROR(
                0513      &             'INIT_FIXED: NF_GET_VARA_DOUBLE tmp_i',
                0514      &             err,bi,bj,myThid )
                0515               err = NF_GET_VARA_DOUBLE(fid, varId_intp5,
                0516      &             vec_start2,vec_count2, tmp_j(1,iq) )
                0517               CALL PROFILES_NF_ERROR(
                0518      &             'INIT_FIXED: NF_GET_VARA_DOUBLE tmp_j',
                0519      &             err,bi,bj,myThid )
                0520              ENDDO
                0521             ENDIF !IF (profilesDoGenGrid)
                0522 
                0523 C     5.2) Loop through this chunk
                0524             DO ic = 1, MIN(1000,ProfNo(num_file,bi,bj)-chunk)
                0525              IF (stopProfiles.EQ.0) THEN
                0526 
                0527               profIsInRunTime = 1
                0528 
                0529               IF (( (tmpyymmdd(ic).GT.yymmddMin).OR.
                0530      &             ((tmpyymmdd(ic).EQ.yymmddMin) .AND.
                0531      &              (tmphhmmss(ic).GT.hhmmssMin)) ) .AND.
                0532      &            ( (tmpyymmdd(ic).LT.yymmddMax).OR.
                0533      &             ((tmpyymmdd(ic).EQ.yymmddMax) .AND.
                0534      &              (tmphhmmss(ic).LT.hhmmssMax)) )) THEN
                0535                hh = INT(tmphhmmss(ic))/10000
                0536                IF (hh.LT.hoursPerDay) THEN
                0537                 profIsInRunTime = 1
                0538 
                0539                 CALL CAL_FULLDATE( INT(tmpyymmdd(ic)),
                0540      &               INT(tmphhmmss(ic)),tmpdate,myThid )
                0541                 CALL CAL_TIMEPASSED( modelstartdate,tmpdate,
                0542      &               tmpdiff,myThid )
                0543                 CALL CAL_TOSECONDS(
                0544      &               tmpdiff,diffsecs,myThid )
                0545 
                0546                 diffsecs = diffsecs+nIter0*deltaTClock
                0547 
                0548                ELSE
                0549 C     If tmp hhmmss is out of range then disregard profile
                0550                 profIsInRunTime = 0
                0551                 diffsecs  = -deltaTClock
                0552                 ProfNo_hh = ProfNo_hh+1
                0553                ENDIF
                0554 
                0555               ELSE
                0556 
                0557                profIsInRunTime = 0
                0558                diffsecs = -deltaTClock
                0559 
                0560               ENDIF !IF (( (tmpyymmdd(ic)
                0561 
                0562 C     5.2a) Determine profiles in current tile domain (lat-lon grid case)
                0563               IF ( (.NOT.profilesDoGenGrid) .AND.
                0564      &             (profIsInRunTime.EQ.1) ) THEN
                0565                IF (xC(sNx+1,1,bi,bj).LT.xC(1,1,bi,bj)) THEN
                0566                 tmp_lon = xC(sNx+1,1,bi,bj)+360. _d 0
                0567                ELSE
                0568                 tmp_lon = xC(sNx+1,1,bi,bj)
                0569                ENDIF
                0570 
                0571                IF ( (xC(1,1,bi,bj).LE.tmp_lon2(ic)) .AND.
                0572      &              (tmp_lon.GT.tmp_lon2(ic)) .AND.
                0573      &              (yC(1,1,bi,bj).LE.tmp_lat2(ic)) .AND.
                0574      &              (yC(1,sNy+1,bi,bj).GT.tmp_lat2(ic)) ) THEN
                0575                 lon_cur = tmp_lon2(ic)
                0576                 lat_cur = tmp_lat2(ic)
                0577                ELSEIF ( (xC(sNx+1,1,bi,bj).LT.xC(1,1,bi,bj))
                0578      &               .AND. (xC(1,1,bi,bj).LE. tmp_lon2(ic)+360. _d 0)
                0579      &               .AND. (tmp_lon.GT.tmp_lon2(ic)+360. _d 0)
                0580      &               .AND. (yC(1,1,bi,bj).LE.tmp_lat2(ic))
                0581      &               .AND. (yC(1,sNy+1,bi,bj).GT.tmp_lat2(ic)) ) THEN
                0582                 lon_cur = tmp_lon2(ic)+360. _d 0
                0583                 lat_cur = tmp_lat2(ic)
                0584                ELSE
                0585                 profIsInRunTime = 0
                0586                ENDIF
                0587 
                0588 C     Determine value of i,j to the south-ouest of data point
                0589                prof_i = -10
                0590                prof_j = -10
                0591                lon_1  = -10
                0592                lon_2  = -10
                0593                lat_1  = -10
                0594                lat_2  = -10
                0595 
                0596                IF (profIsInRunTime.EQ.1) THEN
                0597                 DO j = 1, sNy+1
                0598                  DO i = 1, sNx+1
                0599 C     Value of j, south of the data point:
                0600                   IF ( (yC(i,j,bi,bj).LE.lat_cur) .AND.
                0601      &                 (yC(i,j+1,bi,bj).GT.lat_cur) ) THEN
                0602                    prof_j = j
                0603                    lat_1  = yC(i,j,bi,bj)
                0604                    lat_2  = yC(i,j+1,bi,bj)
                0605                   ENDIF
                0606 
                0607 C     Value of i, west of the data point:
                0608                   IF (xC(i+1,j,bi,bj).LT.xC(1,j,bi,bj)) THEN
                0609                    lon_tmp2 = xC(i+1,j,bi,bj)+360
                0610                   ELSE
                0611                    lon_tmp2 = xC(i+1,j,bi,bj)
                0612                   ENDIF
                0613                   IF (xC(i,j,bi,bj).LT.xC(1,j,bi,bj)) THEN
                0614                    lon_tmp1 = xC(i,j,bi,bj)+360
                0615                   ELSE
                0616                    lon_tmp1 = xC(i,j,bi,bj)
                0617                   ENDIF
                0618 
                0619                   IF ( lon_tmp1.LE.lon_cur .AND.
                0620      &                 lon_tmp2.GT.lon_cur ) THEN
                0621                    prof_i = i
                0622                    lon_1  = lon_tmp1
                0623                    lon_2  = lon_tmp2
                0624                   ENDIF
                0625 
                0626                  ENDDO !DO i
                0627                 ENDDO !DO j
                0628                ENDIF !IF (profIsRunTime.EQ.1)
                0629 
                0630                IF ((prof_i.EQ.-10).OR.(prof_j.EQ.-10)) THEN
                0631                 profIsInRunTime = 0
                0632                ENDIF
                0633 
                0634                IF (profIsInRunTime.EQ.1) THEN
                0635 C     If yes then store prof_time and longitude and latitude:
                0636                 ProfNo_tile = ProfNo_tile+1
                0637                 prof_time(num_file,ProfNo_tile,bi,bj)     =  diffsecs
                0638                 prof_lon(num_file,ProfNo_tile,bi,bj)      = lon_cur
                0639                 prof_lat(num_file,ProfNo_tile,bi,bj)      = lat_cur
                0640                 prof_ind_glob(num_file,ProfNo_tile,bi,bj) = ic+chunk
                0641 # ifdef ALLOW_PROFILES_SAMPLESPLIT_COST
                0642                 prof_ind_avgbin(num_file,ProfNo_tile,bi,bj)
                0643      &               = tmp_avgbin(ic)
                0644 # endif /* ALLOW_PROFILES_SAMPLESPLIT_COST */
                0645 C     Then store interpolation coeffs and indices
                0646                 lon_fac = (lon_cur-lon_1)/(lon_2-lon_1)
                0647                 lat_fac = (lat_cur-lat_1)/(lat_2-lat_1)
                0648 
                0649                 prof_interp_weights(num_file,ProfNo_tile,1,bi,bj)
                0650      &               = (1-lon_fac)*(1-lat_fac)
                0651                 prof_interp_i(num_file,ProfNo_tile,1,bi,bj) = prof_i
                0652                 prof_interp_j(num_file,ProfNo_tile,1,bi,bj) = prof_j
                0653 
                0654                 prof_interp_weights(num_file,ProfNo_tile,2,bi,bj)
                0655      &               = lon_fac*(1-lat_fac)
                0656                 prof_interp_i(num_file,ProfNo_tile,2,bi,bj) = prof_i+1
                0657                 prof_interp_j(num_file,ProfNo_tile,2,bi,bj) = prof_j
                0658 
                0659                 prof_interp_weights(num_file,ProfNo_tile,3,bi,bj)
                0660      &               = (1-lon_fac)*lat_fac
                0661                 prof_interp_i(num_file,ProfNo_tile,3,bi,bj) = prof_i
                0662                 prof_interp_j(num_file,ProfNo_tile,3,bi,bj) = prof_j+1
                0663 
                0664                 prof_interp_weights(num_file,ProfNo_tile,4,bi,bj)
                0665      &               = lon_fac*lat_fac
                0666                 prof_interp_i(num_file,ProfNo_tile,4,bi,bj) = prof_i+1
                0667                 prof_interp_j(num_file,ProfNo_tile,4,bi,bj) = prof_j+1
                0668                ENDIF !IF (profIsRunTime.EQ.1)
                0669 
                0670 C     5.2a) Determine profiles in current tile domain (gen grid case)
                0671               ELSEIF (profIsInRunTime.EQ.1) THEN
                0672 
                0673                IF (stopGenericGrid.EQ.0) THEN
                0674                 IF (ABS(tmp_xC11(ic)-xC(1,1,bi,bj)).LT.0.0001 _d 0
                0675      &         .AND.ABS(tmp_yC11(ic)-yC(1,1,bi,bj)).LT.0.0001 _d 0
                0676      &         .AND.ABS(tmp_xCNINJ(ic)-xC(sNx,sNy,bi,bj)).LT.0.0001 _d 0
                0677      &         .AND.ABS(tmp_yCNINJ(ic)-yC(sNx,sNy,bi,bj)).LT.0.0001 _d 0
                0678      &         .AND.profIsInRunTime.EQ.1 ) THEN
                0679 C     If yes then store prof_time and interpolation coeffs and indices:
                0680                  ProfNo_tile = ProfNo_tile+1
                0681                  prof_time(num_file,ProfNo_tile,bi,bj) = diffsecs
                0682 # ifdef ALLOW_PROFILES_SAMPLESPLIT_COST
                0683                  prof_ind_avgbin(num_file,rofNo_tile,bi,bj)
                0684      &                = tmp_avgbin(ic)
                0685 # endif /* ALLOW_PROFILES_SAMPLESPLIT_COST */
                0686                  prof_interp_xC11(num_file,ProfNo_tile,bi,bj)
                0687      &                = tmp_xC11(ic)
                0688                  prof_interp_yC11(num_file,ProfNo_tile,bi,bj)
                0689      &                = tmp_yC11(ic)
                0690                  prof_interp_xCNINJ(num_file,ProfNo_tile,bi,bj)
                0691      &                = tmp_xCNINJ(ic)
                0692                  prof_interp_yCNINJ(num_file,ProfNo_tile,bi,bj)
                0693      &                = tmp_yCNINJ(ic)
                0694                  tmp_sum_weights = 0. _d 0
                0695 
                0696                  DO iq = 1, iINTERP
                0697 
                0698                   prof_interp_weights(num_file,ProfNo_tile,iq,bi,bj)
                0699      &                 = tmp_weights(ic,iq)
                0700                   prof_interp_i(num_file,ProfNo_tile,iq,bi,bj)
                0701      &                 = INT(tmp_i(ic,iq))
                0702                   prof_interp_j(num_file,ProfNo_tile,iq,bi,bj)
                0703      &                 = INT(tmp_j(ic,iq))
                0704                   tmp_sum_weights = tmp_sum_weights + tmp_weights(ic,iq)
                0705 C     More test of the inputs: is the offline-computed
                0706 C     Interpolation information consistent (self and with grid)
                0707                   IF ( tmp_i(ic,iq).LT.0     .OR. tmp_j(ic,iq).LT.0 .OR.
                0708      &                 tmp_i(ic,iq).GT.sNx+1 .OR. tmp_j(ic,iq).GT.sNy+1
                0709      &                 ) THEN
                0710                    WRITE(msgBuf,'(4A)')
                0711      &                  'PROFILES_INIT_FIXED: file ',
                0712      &                  profilesfile(1:IL),
                0713      &                  '.nc has inconsistent interp ',
                0714      &                  '(profilesDoGenGrid; out of tile)'
                0715                    CALL PRINT_ERROR( msgBuf, myThid )
                0716 
                0717                    stopGenericGrid = 1
                0718                   ENDIF
                0719 
                0720 # ifdef ALLOW_PROFILES_EXCLUDE_CORNERS
                0721                   IF (tmp_weights(ic,iq).NE.0. _d 0) THEN
                0722 C     We should not reach the corners just outside of the current tile.
                0723 c             IF ( (tmp_i(ic,iq).EQ.0.AND.tmp_j(ic,iq).EQ.0)
                0724 c    &        .OR. (tmp_i(ic,iq).EQ.sNx+1.AND.tmp_j(ic,iq).EQ.sNy+1)
                0725 c    &        .OR. (tmp_i(ic,iq).EQ.0.AND.tmp_j(ic,iq).EQ.sNy+1)
                0726 c    &        .OR. (tmp_i(ic,iq).EQ.sNx+1.AND.tmp_j(ic,iq).EQ.0) ) THEN
                0727 C     This should be the same, fewer evaluations, and easier to read?
                0728                    IF((tmp_i(ic,iq).EQ.0.OR.tmp_i(ic,iq).EQ.sNx+1) .AND.
                0729      &                (tmp_j(ic,iq).EQ.0.OR.tmp_j(ic,iq).EQ.sNy+1)
                0730      &                ) THEN
                0731                     WRITE(msgBuf,'(5A)')
                0732      &                   'PROFILES_INIT_FIXED: file ',
                0733      &                   profilesfile(1:IL),
                0734      &                   '.nc has inconsistent interp',
                0735      &                   '(profilesDoGenGrid; overlapping',
                0736      &                   'corner)'
                0737                     CALL PRINT_ERROR( msgBuf, myThid )
                0738 
                0739                     stopGenericGrid = 1
                0740                    ENDIF
                0741                   ENDIF
                0742 # endif /* ALLOW_PROFILES_EXCLUDE_CORNERS */
                0743 
                0744                   IF ( (tmp_weights(ic,iq).LT.0. _d 0) .OR.
                0745      &                 (tmp_weights(ic,iq).GT.1. _d 0) ) THEN
                0746                    WRITE(msgBuf,'(5A)')
                0747      &                  'PROFILES_INIT_FIXED: file ',
                0748      &                  profilesfile(1:IL),
                0749      &                  '.nc has inconsistent interp',
                0750      &                  'weights (profilesDoGenGrid; ',
                0751      &                  'sum oustide 0-1)'
                0752                    CALL PRINT_ERROR( msgBuf, myThid)
                0753 
                0754                    stopGenericGrid = 1
                0755                   ENDIF
                0756 
                0757                  ENDDO !DO iq
                0758 
                0759                  IF (ABS(tmp_sum_weights -1. _d 0).GT. 0.0001 _d 0) THEN
                0760                   WRITE(msgBuf,'(4A)')
                0761      &                 'PROFILES_INIT_FIXED: file ', profilesfile(1:IL),
                0762      &                 '.nc has inconsistent interp weights',
                0763      &                 ' (profilesDoGenGrid; dont add to 1)'
                0764                   CALL PRINT_ERROR( msgBuf, myThid)
                0765 
                0766                   stopGenericGrid = 1
                0767                  ENDIF
                0768 
                0769                  prof_ind_glob(num_file,ProfNo_tile,bi,bj) = ic + chunk
                0770 
                0771                 ENDIF !IF ( ( ABS(tmp_xC11(ic)...
                0772                ENDIF !IF (stopGenericGrid.EQ.0)
                0773               ENDIF!IF (.NOT.profilesDoGenGrid ... ELSEIF ...
                0774 
                0775 C     Check that maximum size was not reached:
                0776               IF (ProfNo_tile.GE.NOBSGLOB) THEN
                0777                WRITE(msgBuf,'(3A)')
                0778      &              'PROFILES_INIT_FIXED: file ', profilesfile(1:IL),
                0779      &              '.nc not read properly (increase NOBSGLOB).'
                0780                CALL PRINT_ERROR( msgBuf, myThid)
                0781 
                0782                stopProfiles = 1
39ce977435 Gael*0783 
13ce79fe94 Ivan*0784               ENDIF
39ce977435 Gael*0785 
13ce79fe94 Ivan*0786              ENDIF !IF (stopProfiles.EQ.0)
                0787             ENDDO !DO ic
39ce977435 Gael*0788 
13ce79fe94 Ivan*0789            ENDIF !IF (MIN(ProfNo(num_file,bi,bj)
                0790           ENDDO !DO chunkProf
39ce977435 Gael*0791 
13ce79fe94 Ivan*0792           ProfNo(num_file,bi,bj) = ProfNo_tile
39ce977435 Gael*0793 
13ce79fe94 Ivan*0794           WRITE(msgBuf,'(A,I9)')
                0795      &         '  # of profiles with erroneous HHMMSS values =',
                0796      &         ProfNo_hh
                0797           CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
ba63501b4c Gael*0798 
13ce79fe94 Ivan*0799           WRITE(msgBuf,'(A,I9)')
                0800      &         '  # of profiles within tile and time period  =',
                0801      &         ProfNo(num_file,bi,bj)
                0802           CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
ba63501b4c Gael*0803 
13ce79fe94 Ivan*0804 C     6) Available variables in the data set
                0805           DO num_var = 1, NVARMAX
                0806            prof_num_var_cur(num_file,num_var,bi,bj) = 0
                0807           ENDDO
                0808           prof_num_var_tot(num_file,bi,bj) = 0
                0809 
                0810           DO num_var = 1, NVARMAX
                0811            JL = ILNBLNK( prof_names(num_file,num_var) )
                0812            err = NF_INQ_VARID( fid,
                0813      &          prof_names(num_file,num_var)(1:JL), varId1 )
                0814 
                0815            IF (err.EQ.NF_NOERR) THEN
                0816             vec_quantities(num_file,num_var,bi,bj) = .TRUE.
                0817 
                0818             prof_num_var_tot(num_file,bi,bj) =
                0819      &           prof_num_var_tot(num_file,bi,bj) + 1
                0820             prof_num_var_cur(num_file,num_var,bi,bj) =
                0821      &           prof_num_var_tot(num_file,bi,bj)
                0822            ELSE
                0823             CALL PROFILES_NF_ERROR(
                0824      &           'INIT_FIXED: NF_INQ_VARID prof_names = '//
                0825      &           prof_names(num_file,num_var)(1:JL),err,bi,bj,myThid )
                0826             IF (debugLevel .GE. debLevA) THEN
                0827              WRITE(msgBuf,'(3A)') 'S/R PROFILES_INIT_FIXED: no ',
                0828      &            prof_names(num_file,num_var)(1:JL),
                0829      &            ', setting corresponding vec_quantities = F'
                0830              CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                0831             ENDIF
                0832             vec_quantities(num_file,num_var,bi,bj) = .FALSE.
                0833            ENDIF
                0834 
                0835            IF (vec_quantities(num_file,num_var,bi,bj)) THEN
                0836             KL = ILNBLNK( prof_names(num_file,num_var) )
                0837             JL = ILNBLNK( prof_namesmod(num_file,num_var) )
                0838 
                0839             IF (prof_namesmod(num_file,num_var).EQ.'pTracer') THEN
                0840              WRITE(msgBuf,'(A,I3,5A,I3)') '  variable #',
                0841      &            num_var,' is ' ,
                0842      &            prof_names(num_file,num_var)(1:KL),' and ',
                0843      &            prof_namesmod(num_file,num_var)(1:JL),' #',
                0844      &            prof_itracer(num_file,num_var)
                0845             ELSE
                0846              WRITE(msgBuf,'(A,I3,4A)') '  variable #',
                0847      &            num_var,' is            ' ,
                0848      &            prof_names(num_file,num_var)(1:KL),' and ',
                0849      &            prof_namesmod(num_file,num_var)(1:JL)
                0850             ENDIF
                0851             CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                0852            ENDIF !IF (vec_quantities
                0853           ENDDO !DO num_var
6a770e0a24 Patr*0854 
                0855 C===========================================================
13ce79fe94 Ivan*0856 C     Create files for model counterparts to observations
6a770e0a24 Patr*0857 C===========================================================
                0858 
13ce79fe94 Ivan*0859           IF (ProfNo(num_file,bi,bj).GT.0) THEN
                0860            iG=bi+(myXGlobalLo-1)/sNx ! Kludge until unstructered tiles
                0861            jG=bj+(myYGlobalLo-1)/sNy ! Kludge until unstructered tiles
                0862            JL = ILNBLNK( profilesDir )
                0863 
                0864            IF (profilesDoNcOutput) THEN
                0865             WRITE(fnameequinc,'(3A,2(I3.3,A))')
                0866      &           profilesDir(1:JL),profilesfile(1:IL),
                0867      &           '.',iG,'.',jG,'.equi.nc'
                0868 # ifdef ALLOW_ADJOINT_RUN
                0869             WRITE(adfnameequinc,'(4A,2(I3.3,A))')
                0870      &           profilesDir(1:JL),'ad',profilesfile(1:IL),
                0871      &           '.',iG,'.',jG,'.equi.nc'
                0872 # endif /* ALLOW_ADJOINT_RUN */
                0873 # ifdef ALLOW_TANGENTLINEAR_RUN
                0874             WRITE(tlfnameequinc,'(4A,2(I3.3,A))')
                0875      &           profilesDir(1:JL),'tl',profilesfile(1:IL),
                0876      &           '.',iG,'.',jG,'.equi.nc'
                0877 # endif /* ALLOW_TANGENTLINEAR_RUN */
                0878 
                0879             JL = ILNBLNK( fnameequinc )
                0880             INQUIRE( FILE = fnameequinc(1:JL), EXIST = exst )
                0881             IF (.NOT.exst) THEN
                0882              CALL PROFILES_INIT_NCFILE( num_file,
                0883      &            fiddata(num_file,bi,bj), fnameequinc(1:JL),
                0884      &            fidforward(num_file,bi,bj),
                0885      &            ProfNo(num_file,bi,bj),
                0886      &            ProfDepthNo(num_file,bi,bj),
                0887      &            bi, bj, myThid )
                0888             ELSE
                0889 C     Obtain existing NetCDF file id
                0890              err = NF_OPEN( fnameequinc(1:JL), NF_WRITE,
                0891      &            fidforward(num_file,bi,bj) )
                0892              CALL PROFILES_NF_ERROR( 'INIT_FIXED: NF_OPEN fidforward',
                0893      &            err,bi,bj,myThid )
                0894             ENDIF
                0895 
                0896 # ifdef ALLOW_ADJOINT_RUN
                0897             JL = ILNBLNK( adfnameequinc )
                0898             INQUIRE( FILE = adfnameequinc(1:JL), EXIST = exst )
                0899             IF (.NOT.exst) THEN
                0900              CALL PROFILES_INIT_NCFILE( num_file,
                0901      &            fiddata(num_file,bi,bj), adfnameequinc(1:JL),
                0902      &            fidadjoint(num_file,bi,bj),
                0903      &            ProfNo(num_file,bi,bj),
                0904      &            ProfDepthNo(num_file,bi,bj),
                0905      &            bi, bj, myThid )
                0906             ELSE
                0907 C     Obtain existing NetCDF file id
                0908              err = NF_OPEN( adfnameequinc(1:JL), NF_WRITE,
                0909      &            fidadjoint(num_file,bi,bj) )
                0910              CALL PROFILES_NF_ERROR( 'INIT_FIXED: NF_OPEN fidadjoint',
                0911      &            err,bi,bj,myThid )
                0912             ENDIF
                0913 
                0914 # endif /* ALLOW_ADJOINT_RUN */
                0915 # ifdef ALLOW_TANGENTLINEAR_RUN
                0916             JL = ILNBLNK( tlfnameequinc )
                0917             INQUIRE( FILE = tlfnameequinc(1:JL), EXIST = exst )
                0918             IF (.NOT.exst) THEN
                0919              CALL PROFILES_INIT_NCFILE( num_file,
                0920      &            fiddata(num_file,bi,bj), tlfnameequinc(1:JL),
                0921      &            fidtangent(num_file,bi,bj),
                0922      &            ProfNo(num_file,bi,bj),
                0923      &            ProfDepthNo(num_file,bi,bj),
                0924      &            bi, bj, myThid )
                0925             ELSE
                0926 C     Obtain existing NetCDF file id
                0927              err = NF_OPEN( tlfnameequinc(1:JL), NF_WRITE,
                0928      &            fidtangent(num_file,bi,bj) )
                0929              CALL PROFILES_NF_ERROR( 'INIT_FIXED: NF_OPEN fidtangent',
                0930      &            err,bi,bj,myThid )
                0931             ENDIF
                0932 
                0933 # endif /* ALLOW_TANGENTLINEAR_RUN */
                0934            ELSE !IF (profilesDoNcOutput
                0935             WRITE(fnameequinc,'(3A,2(I3.3,A))')
                0936      &           profilesDir(1:JL),profilesfile(1:IL),
                0937      &           '.',iG,'.',jG,'.equi.data'
                0938 
                0939 # ifdef ALLOW_ADJOINT_RUN
                0940             WRITE(adfnameequinc,'(4A,2(I3.3,A))')
                0941      &           profilesDir(1:JL),'ad',profilesfile(1:IL),
                0942      &           '.',iG,'.',jG,'.equi.data'
                0943 
                0944 # endif /* ALLOW_ADJOINT_RUN */
                0945 # ifdef ALLOW_TANGENTLINEAR_RUN
                0946             WRITE(tlfnameequinc,'(4A,2(I3.3,A))')
                0947      &           profilesDir(1:JL),'tl',profilesfile(1:IL),
                0948      &           '.',iG,'.',jG,'.equi.data'
                0949 
                0950 # endif /* ALLOW_TANGENTLINEAR_RUN */
                0951 
                0952             JL = ILNBLNK( fnameequinc )
                0953             INQUIRE( FILE = fnameequinc(1:JL), EXIST = exst )
                0954 # ifdef PROFILES_USE_MDSFINDUNITS
                0955             CALL MDSFINDUNIT( fidforward(num_file,bi,bj), myThid )
                0956 # else
                0957             CALL PROFILES_FINDUNIT( fidforward(num_file,bi,bj), myThid )
                0958 # endif
                0959             IF (.NOT.exst) THEN
                0960              CALL PROFILES_INIT_NCFILE( num_file,
                0961      &            fiddata(num_file,bi,bj), fnameequinc(1:JL),
                0962      &            fidforward(num_file,bi,bj),
                0963      &            ProfNo(num_file,bi,bj),
                0964      &            ProfDepthNo(num_file,bi,bj),
                0965      &            bi, bj, myThid )
                0966             ELSE
                0967              OPEN( fidforward(num_file,bi,bj),
                0968      &            FILE = fnameequinc(1:JL), FORM = 'unformatted',
                0969      &            STATUS = 'unknown', ACCESS = 'direct',
                0970      &            RECL = (ProfDepthNo(num_file,bi,bj)+1)*WORDLENGTH*2 )
                0971             ENDIF
                0972 
                0973 # ifdef ALLOW_ADJOINT_RUN
                0974             JL = ILNBLNK( adfnameequinc )
                0975             INQUIRE( FILE = adfnameequinc(1:JL), EXIST = exst )
                0976 #  ifdef PROFILES_USE_MDSFINDUNITS
                0977             CALL MDSFINDUNIT( fidadjoint(num_file,bi,bj), myThid )
                0978 #  else
                0979             CALL PROFILES_FINDUNIT( fidadjoint(num_file,bi,bj), myThid )
                0980 #  endif
                0981             IF (.NOT.exst) THEN
                0982              CALL PROFILES_INIT_NCFILE( num_file,
                0983      &            fiddata(num_file,bi,bj), adfnameequinc(1:JL),
                0984      &            fidadjoint(num_file,bi,bj),
                0985      &            ProfNo(num_file,bi,bj),
                0986      &            ProfDepthNo(num_file,bi,bj),
                0987      &            bi, bj, myThid )
                0988             ELSE
                0989              OPEN( fidadjoint(num_file,bi,bj),
                0990      &            FILE = adfnameequinc(1:JL), FORM = 'unformatted',
                0991      &            STATUS = 'unknown', ACCESS = 'direct',
                0992      &            RECL = (ProfDepthNo(num_file,bi,bj)+1)*WORDLENGTH*2 )
                0993             ENDIF
                0994 
                0995 # endif /* ALLOW_ADJOINT_RUN */
                0996 # ifdef ALLOW_TANGENTLINEAR_RUN
                0997             JL = ILNBLNK( tlfnameequinc )
                0998             INQUIRE( FILE = tlfnameequinc(1:JL), EXIST = exst )
                0999 #  ifdef PROFILES_USE_MDSFINDUNITS
                1000             CALL MDSFINDUNIT( fidtangent(num_file,bi,bj), myThid )
                1001 #  else
                1002             CALL PROFILES_FINDUNIT( fidtangent(num_file,bi,bj), myThid )
                1003 #  endif
                1004             IF (.NOT.exst) THEN
                1005              CALL PROFILES_INIT_NCFILE( num_file,
                1006      &            fiddata(num_file,bi,bj), tlfnameequinc(1:JL),
                1007      &            fidtangent(num_file,bi,bj),
                1008      &            ProfNo(num_file,bi,bj),
                1009      &            ProfDepthNo(num_file,bi,bj),
                1010      &            bi, bj, myThid )
                1011             ELSE
                1012              OPEN( fidtangent(num_file,bi,bj),
                1013      &            FILE = tlfnameequinc(1:JL), FORM = 'unformatted',
                1014      &            STATUS = 'unknown', ACCESS = 'direct',
                1015      &            RECL = (ProfDepthNo(num_file,bi,bj)+1)*WORDLENGTH*2 )
                1016             ENDIF
                1017 
                1018 # endif /* ALLOW_TANGENTLINEAR_RUN */
                1019 
                1020            ENDIF !IF (profilesDoNcOutput
                1021 
                1022           ENDIF !IF (ProfNo(num_file,bi,bj).GT.0)
                1023 
                1024 C     ===========================================================
                1025          ELSE !IF (IL.NE.0)
                1026 
                1027           ProfNo(num_file,bi,bj) = 0
                1028 
                1029           DO num_var = 1, NVARMAX
                1030            prof_num_var_cur(num_file,num_var,bi,bj) = 0
                1031            vec_quantities(num_file,num_var,bi,bj) = .FALSE.
                1032           ENDDO
6a770e0a24 Patr*1033 
13ce79fe94 Ivan*1034           prof_num_var_tot(num_file,bi,bj) = 0
                1035 
                1036           DO prof_num = 1, NOBSGLOB
                1037            prof_time      (num_file,prof_num,bi,bj) = -999. _d 0
                1038            prof_lon       (num_file,prof_num,bi,bj) = -999. _d 0
                1039            prof_lat       (num_file,prof_num,bi,bj) = -999. _d 0
                1040            prof_ind_glob  (num_file,prof_num,bi,bj) = 0
                1041 # ifdef ALLOW_PROFILES_SAMPLESPLIT_COST
                1042            prof_ind_avgbin(num_file,prof_num,bi,bj) = -999
                1043 # endif /* ALLOW_PROFILES_SAMPLESPLIT_COST */
                1044 
                1045            DO iq = 1, NUM_INTERP_POINTS
                1046             prof_interp_i(num_file,prof_num,iq,bi,bj) = 1
                1047             prof_interp_j(num_file,prof_num,iq,bi,bj) = 1
                1048             prof_interp_weights(num_file,prof_num,iq,bi,bj) = 0. _d 0
                1049            ENDDO
6a770e0a24 Patr*1050 
13ce79fe94 Ivan*1051            prof_interp_xC11  (num_file,prof_num,bi,bj) = -999. _d 0
                1052            prof_interp_yC11  (num_file,prof_num,bi,bj) = -999. _d 0
                1053            prof_interp_xCNINJ(num_file,prof_num,bi,bj) = -999. _d 0
                1054            prof_interp_yCNINJ(num_file,prof_num,bi,bj) = -999. _d 0
                1055           ENDDO !DO prof_num
6a770e0a24 Patr*1056 
13ce79fe94 Ivan*1057          ENDIF !IF (IL.NE.0)
6a770e0a24 Patr*1058 
13ce79fe94 Ivan*1059         ENDDO !DO num_file
6a770e0a24 Patr*1060 
13ce79fe94 Ivan*1061 # ifdef ALLOW_PROFILES_SAMPLESPLIT_COST
                1062 C     Find unique depth levels from all profile datasets
                1063 C     initialize prof_depth_comb
                1064         IF (bi.EQ.1.AND.bj.EQ.1) THEN
6b2230d510 Ou W*1065          NLEVELCOMB = 0
                1066          NLEVELCOMBRL = NLEVELCOMB
13ce79fe94 Ivan*1067         ENDIF
6b2230d510 Ou W*1068 
13ce79fe94 Ivan*1069         DO kCMax = 1, NLEVELCOMBMAX
                1070          prof_depth_comb(kCMax,bi,bj) = -999. _d 0
                1071         ENDDO
6b2230d510 Ou W*1072 
13ce79fe94 Ivan*1073         m = 1
                1074         DO num_file = 1, NFILESPROFMAX
                1075          DO kLev = 1, ProfDepthNo(num_file,bi,bj)
                1076           IF (m.EQ.1) THEN
                1077            prof_depth_comb(m,bi,bj) = prof_depth(num_file,kLev,bi,bj)
6b2230d510 Ou W*1078            m = m + 1
13ce79fe94 Ivan*1079 C     Sort
                1080           ELSE !IF (m.EQ.1)
                1081            DO kCMax = 1, NLEVELCOMBMAX-1
                1082             IF (prof_depth_comb(kCMax,bi,bj) .NE. -999. _d 0) THEN
                1083              IF ( prof_depth(num_file,kLev,bi,bj).LT.
                1084      &            prof_depth_comb(kCMax,bi,bj).AND. kCMax.EQ.1 ) THEN
                1085 
                1086               prof_depth_comb(NLEVELCOMBMAX,bi,bj) =
                1087      &             prof_depth_comb(kCMax,bi,bj)
                1088 
                1089               prof_depth_comb(kCMax,bi,bj) =
                1090      &             prof_depth(num_file,kLev,bi,bj)
                1091 
                1092               DO kC = NLEVELCOMBMAX-1, kCMax+2, -1
                1093                prof_depth_comb(kC,bi,bj) = prof_depth_comb(kC-1,bi,bj)
                1094               ENDDO !DO kC
                1095 
                1096               prof_depth_comb(kCMax+1,bi,bj) =
                1097      &             prof_depth_comb(NLEVELCOMBMAX,bi,bj)
                1098 
                1099              ELSEIF ( prof_depth(num_file,kLev,bi,bj).GT.
                1100      &             prof_depth_comb(kCMax,bi,bj) .AND.
                1101      &             prof_depth(num_file,kLev,bi,bj).LT.
                1102      &             prof_depth_comb(kCMax+1,bi,bj) ) THEN
                1103 
                1104               prof_depth_comb(NLEVELCOMBMAX,bi,bj) =
                1105      &             prof_depth_comb(kCMax+1,bi,bj)
                1106 
                1107               prof_depth_comb(kCMax+1,bi,bj) =
                1108      &            prof_depth(num_file, kLev,bi,bj)
                1109 
                1110               DO kC = NLEVELCOMBMAX-1, kCMax+3, -1
                1111                prof_depth_comb(kC,bi,bj) = prof_depth_comb(kC-1,bi,bj)
                1112               ENDDO !DO kC
                1113 
                1114               prof_depth_comb(kCMax+2,bi,bj)=
                1115      &             prof_depth_comb(NLEVELCOMBMAX,bi,bj)
                1116 
                1117              ELSEIF ( prof_depth(num_file,kLev,bi,bj).GT.
                1118      &             prof_depth_comb(kCMax,bi,bj) .AND.
                1119      &             prof_depth_comb(kCMax+1,bi,bj).EQ. -999. _d 0 ) THEN
                1120 
                1121               prof_depth_comb(kCMax+1,bi,bj) =
                1122      &             prof_depth(num_file, kLev,bi,bj)
6b2230d510 Ou W*1123 
13ce79fe94 Ivan*1124              ENDIF !IF (prof_depth(num_file,kLev,bi,bj).LT.
                1125             ENDIF !IF (prof_depth_comb(kCMax,bi,bj).NE.-999.
                1126            ENDDO !DO kCMax
                1127 
                1128           ENDIF !IF (m.EQ.1)
                1129 
                1130           IF (m.GE.NLEVELCOMBMAX-2) THEN
                1131            WRITE(msgBuf,'(A)') 'increase NLEVELCOMBMAX'
                1132            CALL PRINT_ERROR( msgBuf, myThid )
c9bf163375 Ivan*1133           ENDIF
13ce79fe94 Ivan*1134 
                1135          ENDDO ! DO kLev
                1136         ENDDO ! DO num_file
                1137 
                1138         prof_depth_comb(NLEVELCOMBMAX,bi,bj) = -999. _d 0
                1139 
                1140 C     Diagnostics output
                1141         DO kCMax = 1, NLEVELCOMBMAX
                1142          IF ( prof_depth_comb(kCMax,bi,bj) .GE. 0. _d 0
                1143      &        .AND. NLEVELCOMB.LT.kCMax ) THEN
                1144           NLEVELCOMB = kCMax
                1145 
                1146           IF (kCMax.GE.NLEVELCOMBMAX-2) THEN
                1147            WRITE(msgBuf,'(A,2I6)')
                1148      &          'increase NLEVELCOMBMAX: kCMax, NLEVELCOMBMA  ',
                1149      &          kCMax, NLEVELCOMBMAX
                1150            CALL PRINT_ERROR( msgBuf, myThid )
                1151           ENDIF
                1152          ENDIF !IF ( prof_depth_comb(m,bi,bj) .GE. 0. _d 0
                1153         ENDDO !DO kCMax
                1154 
                1155         WRITE(msgBuf,'(A,I6,D20.5)') 'NLEVELCOMB = ', NLEVELCOMB
                1156         CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                1157 # endif /* ALLOW_PROFILES_SAMPLESPLIT_COST */
6b2230d510 Ou W*1158 
6a770e0a24 Patr*1159 C===========================================================
13ce79fe94 Ivan*1160 C     Error analysis:
ba63501b4c Gael*1161 C===========================================================
                1162 
13ce79fe94 Ivan*1163 C     1) Provide interpolation information
                1164         IF (stopGenericGrid.EQ.2) THEN
ba63501b4c Gael*1165          iG=bi+(myXGlobalLo-1)/sNx ! Kludge until unstructered tiles
                1166          jG=bj+(myYGlobalLo-1)/sNy ! Kludge until unstructered tiles
13ce79fe94 Ivan*1167 C     XC grid
                1168          CALL MDSFINDUNIT( fid, myThid )
                1169          WRITE(fnameequinc,'(A,2(I3.3,A),2(I4.4,A))')
                1170      &        'profilesXCincl1PointOverlap.',iG,'.',jG,'.',
                1171      &        sNx,'.',sNy,'.data'
                1172 
                1173          recLen = MDS_RECLEN( 64,(sNx+2)*(sNy+2),myThid )
                1174 
                1175          WRITE(iUnit,'(2A,/,A)')
                1176      &        'PROFILES_INIT_FIXED: creating grid from profiles; ',
                1177      &        'file:', fnameequinc
                1178          JL  = ILNBLNK( fnameequinc )
                1179          OPEN( fid, FILE = fnameequinc(1:JL), FORM = 'unformatted',
                1180      &        STATUS = 'unknown', ACCESS = 'direct', RECL = recLen )
                1181 
                1182          DO j = 0, sNy+1
                1183           DO i = 0, sNx+1
                1184            xy_buffer_r8(i,j) = xC(i,j,bi,bj)
                1185           ENDDO
ba63501b4c Gael*1186          ENDDO
13ce79fe94 Ivan*1187 
                1188 # ifdef _BYTESWAPIO
                1189          CALL MDS_BYTESWAPR8( (sNx+2)*(sNy+2),xy_buffer_r8 )
                1190 # endif
                1191          WRITE(fid,rec=1) xy_buffer_r8
                1192          CLOSE(fid)
                1193 C     YC grid
                1194          CALL MDSFINDUNIT( fid, myThid )
                1195          WRITE(fnameequinc,'(A,2(I3.3,A),2(I4.4,A))')
                1196      &        'profilesYCincl1PointOverlap.',iG,'.',jG,'.',
                1197      &        sNx,'.',sNy,'.data'
                1198 
                1199          recLen = MDS_RECLEN( 64,(sNx+2)*(sNy+2),myThid )
                1200 
                1201          WRITE(iUnit,'(2A,/,A)')
                1202      &        'PROFILES_INIT_FIXED: creating grid from profiles; ',
                1203      &        'file:', fnameequinc
                1204          JL  = ILNBLNK( fnameequinc )
                1205          OPEN( fid, FILE = fnameequinc(1:JL), FORM = 'unformatted',
                1206      &        STATUS = 'unknown', ACCESS = 'direct', RECL = recLen )
                1207 
                1208          DO j = 0, sNy+1
                1209           DO i = 0, sNx+1
                1210            xy_buffer_r8(i,j) = yC(i,j,bi,bj)
                1211           ENDDO
ba63501b4c Gael*1212          ENDDO
                1213 
13ce79fe94 Ivan*1214 # ifdef _BYTESWAPIO
                1215          CALL MDS_BYTESWAPR8( (sNx+2)*(sNy+2),xy_buffer_r8 )
                1216 # endif
                1217          WRITE(fid,rec=1) xy_buffer_r8
                1218          CLOSE(fid)
ba63501b4c Gael*1219 
13ce79fe94 Ivan*1220          WRITE(msgBuf,'(2A)')
                1221      &        'PROFILES_INIT_FIXED : when using profilesDoGenGrid ',
                1222      &        'you have to provide interpolation coeffs etc. '
                1223          CALL PRINT_ERROR( msgBuf, myThid )
71a5587721 Gael*1224 
13ce79fe94 Ivan*1225          WRITE(msgBuf,'(2A)')
                1226      &        'and some of your nc files dont have them. ',
                1227      &        'You could use profiles_prep_mygrid.m and/or'
                1228          CALL PRINT_ERROR( msgBuf, myThid )
                1229 
                1230          WRITE(msgBuf,'(A)')
                1231      &        'use the grid info in profiles*incl1PointOverlap*data'
                1232          CALL PRINT_ERROR( msgBuf, myThid )
                1233 
                1234          stopProfiles = 1
                1235         ENDIF !IF (stopGenericGrid.EQ.2)
                1236 
                1237        ENDDO !DO bi
                1238       ENDDO !DO bj
6b2230d510 Ou W*1239 
c9bf163375 Ivan*1240       _END_MASTER( myThid )
ba63501b4c Gael*1241       _BARRIER
                1242 
13ce79fe94 Ivan*1243 # ifdef ALLOW_PROFILES_SAMPLESPLIT_COST
                1244       NLEVELCOMBRL = NLEVELCOMB
                1245       _GLOBAL_MAX_RL( NLEVELCOMBRL, myThid )
                1246       NLEVELCOMB = NLEVELCOMBRL
                1247 # endif
                1248 
                1249 C     2) Stop after other kind of errors
                1250       CALL GLOBAL_SUM_INT( stopProfiles, myThid )
c9bf163375 Ivan*1251       IF ( stopProfiles.GE.1) THEN
13ce79fe94 Ivan*1252        CALL ALL_PROC_DIE( myThid )
                1253        STOP 'ABNORMAL END: S/R PROFILES_INIT_FIXED'
c9bf163375 Ivan*1254       ENDIF
39ce977435 Gael*1255 
13ce79fe94 Ivan*1256       CALL GLOBAL_SUM_INT( stopGenericGrid, myThid )
c9bf163375 Ivan*1257       IF ( stopGenericGrid.GE.1) THEN
13ce79fe94 Ivan*1258        CALL ALL_PROC_DIE( myThid )
                1259        STOP 'ABNORMAL END: S/R PROFILES_INIT_FIXED'
c9bf163375 Ivan*1260       ENDIF
ba63501b4c Gael*1261 
13ce79fe94 Ivan*1262       WRITE(msgBuf,'(A)') ' '
                1263       CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                1264       WRITE(msgBuf,'(A)')
                1265      &     '// ======================================================='
                1266       CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                1267       WRITE(msgBuf,'(A)')
                1268      &     '// insitu profiles model sampling >>> END <<<'
                1269       CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                1270       WRITE(msgBuf,'(A)')
                1271      &     '// ======================================================='
                1272       CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
                1273       WRITE(msgBuf,'(A)') ' '
                1274       CALL PRINT_MESSAGE( msgBuf, iUnit, SQUEEZE_RIGHT, myThid )
c9bf163375 Ivan*1275 #endif /* ALLOW_PROFILES */
6a770e0a24 Patr*1276 
d5aa75d39a Jean*1277       RETURN
6a770e0a24 Patr*1278       END