| 1 |
C $Header: /u/gcmpack/MITgcm/pkg/streamice/streamice_petsc_numerate.F,v 1.3 2014/03/26 17:49:51 dgoldberg Exp $ |
| 2 |
C $Name: $ |
| 3 |
|
| 4 |
#include "STREAMICE_OPTIONS.h" |
| 5 |
|
| 6 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
| 7 |
|
| 8 |
CBOP |
| 9 |
SUBROUTINE STREAMICE_PETSC_NUMERATE (myThid) |
| 10 |
|
| 11 |
C /============================================================\ |
| 12 |
C | SUBROUTINE | |
| 13 |
C | o | |
| 14 |
C |============================================================| |
| 15 |
C | | |
| 16 |
C \============================================================/ |
| 17 |
IMPLICIT NONE |
| 18 |
|
| 19 |
C === Global variables === |
| 20 |
#include "SIZE.h" |
| 21 |
#include "GRID.h" |
| 22 |
#include "EEPARAMS.h" |
| 23 |
#include "PARAMS.h" |
| 24 |
#include "STREAMICE.h" |
| 25 |
#ifdef ALLOW_PETSC |
| 26 |
# ifdef ALLOW_USE_MPI |
| 27 |
# include "EESUPPORT.h" |
| 28 |
# endif |
| 29 |
#endif |
| 30 |
|
| 31 |
INTEGER myThid |
| 32 |
|
| 33 |
#ifdef ALLOW_STREAMICE |
| 34 |
|
| 35 |
INTEGER i, j, bi, bj, ki, kj |
| 36 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
| 37 |
#ifdef ALLOW_USE_MPI |
| 38 |
integer mpiRC, mpiMyWid |
| 39 |
#endif |
| 40 |
#ifdef ALLOW_PETSC |
| 41 |
_RS DoFCount |
| 42 |
integer n_dofs_proc_loc (0:nPx*nPy-1) |
| 43 |
integer n_dofs_cum_sum (0:nPx*nPy-1) |
| 44 |
|
| 45 |
|
| 46 |
DO bj = myByLo(myThid), myByHi(myThid) |
| 47 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
| 48 |
DO j=1,sNy |
| 49 |
DO i=1,sNx |
| 50 |
streamice_petsc_dofs_u (i,j,bi,bj) = -2.0 |
| 51 |
streamice_petsc_dofs_v (i,j,bi,bj) = -2.0 |
| 52 |
ENDDO |
| 53 |
ENDDO |
| 54 |
ENDDO |
| 55 |
ENDDO |
| 56 |
|
| 57 |
DoFCount = -1.0 |
| 58 |
DO bj = myByLo(myThid), myByHi(myThid) |
| 59 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
| 60 |
DO j=1,sNy |
| 61 |
DO i=1,sNx |
| 62 |
|
| 63 |
|
| 64 |
C DOFS ARE NUMBERED AS FOLLOWS ON PROCESSOR DOMAIN: |
| 65 |
C grid is stepped through in order bj, bi, j, i |
| 66 |
C 1) if umask(i,j,bi,bj)==1, the counter is updated by 1; |
| 67 |
C streamice_petsc_dofs_u is assigned the counter; |
| 68 |
C o/w streamice_petsc_dofs_u is assigned -1 |
| 69 |
C 2) if vmask(i,j,bi,bj)==1, the counter is updated by 1; |
| 70 |
C streamice_petsc_dofs_v is assigned the counter; |
| 71 |
C o/w streamice_petsc_dofs_v is assigned -1 |
| 72 |
C NOTE THESE NUMBERING ARRAYS ARE USED TO CONSTRUCT PETSC VECTORS AND MATRIX |
| 73 |
|
| 74 |
if (STREAMICE_umask (i,j,bi,bj).eq.1.0) THEN |
| 75 |
DoFCount = DoFCount + 1.0 |
| 76 |
streamice_petsc_dofs_u (i,j,bi,bj) = DoFCount |
| 77 |
else |
| 78 |
streamice_petsc_dofs_u (i,j,bi,bj) = -1.0 |
| 79 |
endif |
| 80 |
|
| 81 |
if (STREAMICE_vmask (i,j,bi,bj).eq.1.0) THEN |
| 82 |
DoFCount = DoFCount + 1.0 |
| 83 |
streamice_petsc_dofs_v (i,j,bi,bj) = DoFCount |
| 84 |
else |
| 85 |
streamice_petsc_dofs_v (i,j,bi,bj) = -1.0 |
| 86 |
endif |
| 87 |
|
| 88 |
ENDDO |
| 89 |
ENDDO |
| 90 |
ENDDO |
| 91 |
ENDDO |
| 92 |
|
| 93 |
print *, "DOF_COUNT", dofcount |
| 94 |
|
| 95 |
#ifdef ALLOW_USE_MPI |
| 96 |
|
| 97 |
DO i=0,nPx*nPy-1 |
| 98 |
n_dofs_proc_loc (i) = 0 |
| 99 |
ENDDO |
| 100 |
|
| 101 |
CALL MPI_COMM_RANK( MPI_COMM_WORLD, mpiMyWId, mpiRC ) |
| 102 |
|
| 103 |
n_dofs_proc_loc (mpiMyWId) = INT(DoFCount)+1 |
| 104 |
|
| 105 |
CALL MPI_Allreduce(n_dofs_proc_loc,n_dofs_process,nPx*nPy, |
| 106 |
& MPI_INTEGER, MPI_SUM,MPI_COMM_MODEL,mpiRC) |
| 107 |
|
| 108 |
n_dofs_cum_sum(0) = 0 |
| 109 |
|
| 110 |
DO i=1,nPx*nPy-1 |
| 111 |
n_dofs_cum_sum(i) = n_dofs_cum_sum(i-1)+ |
| 112 |
& n_dofs_process(i-1) |
| 113 |
ENDDO |
| 114 |
|
| 115 |
#else /* ALLOW_USE_MPI */ |
| 116 |
|
| 117 |
n_dofs_process (0) = INT(DoFCount)+1 |
| 118 |
n_dofs_cum_sum (0) = INT(DoFCount)+1 |
| 119 |
|
| 120 |
#endif /* ALLOW_USE_MPI */ |
| 121 |
|
| 122 |
DO bj = myByLo(myThid), myByHi(myThid) |
| 123 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
| 124 |
DO j=1,sNy |
| 125 |
DO i=1,sNx |
| 126 |
IF (streamice_petsc_dofs_u(i,j,bi,bj).ge.0 ) THEN |
| 127 |
streamice_petsc_dofs_u(i,j,bi,bj) = |
| 128 |
& streamice_petsc_dofs_u(i,j,bi,bj) + |
| 129 |
& n_dofs_cum_sum(mpimywid) |
| 130 |
ENDIF |
| 131 |
IF (streamice_petsc_dofs_v(i,j,bi,bj).ge.0 ) THEN |
| 132 |
streamice_petsc_dofs_v(i,j,bi,bj) = |
| 133 |
& streamice_petsc_dofs_v(i,j,bi,bj) + |
| 134 |
& n_dofs_cum_sum(mpimywid) |
| 135 |
ENDIF |
| 136 |
ENDDO |
| 137 |
ENDDO |
| 138 |
ENDDO |
| 139 |
ENDDO |
| 140 |
|
| 141 |
_EXCH_XY_RS(streamice_petsc_dofs_u,myThid) |
| 142 |
_EXCH_XY_RS(streamice_petsc_dofs_v,myThid) |
| 143 |
|
| 144 |
|
| 145 |
#endif /* ALLOW_PETSC */ |
| 146 |
|
| 147 |
|
| 148 |
#endif |
| 149 |
RETURN |
| 150 |
END |